xref: /petsc/include/petscmat.h (revision 8e7ba810fe86bd58b4e28290b119b414003d9d51)
12eac72dbSBarry Smith /*
22eac72dbSBarry Smith      Include file for the matrix component of PETSc
32eac72dbSBarry Smith */
40a835dfdSSatish Balay #ifndef __PETSCMAT_H
50a835dfdSSatish Balay #define __PETSCMAT_H
62c8e378dSBarry Smith #include <petscvec.h>
72eac72dbSBarry Smith 
8d9274352SBarry Smith /*S
98f6c3df8SBarry Smith      Mat - Abstract PETSc matrix object used to manage all linear operators in PETSc, even those without
108f6c3df8SBarry Smith            an explicit sparse representation (such as matrix-free operators)
112eac72dbSBarry Smith 
12d91e6319SBarry Smith    Level: beginner
13d91e6319SBarry Smith 
14d9274352SBarry Smith   Concepts: matrix; linear operator
15d9274352SBarry Smith 
168f6c3df8SBarry Smith .seealso:  MatCreate(), MatType, MatSetType(), MatDestroy()
17d9274352SBarry Smith S*/
18d9274352SBarry Smith typedef struct _p_Mat*           Mat;
19d9274352SBarry Smith 
2076bdecfbSBarry Smith /*J
218f6c3df8SBarry Smith     MatType - String with the name of a PETSc matrix type
22d9274352SBarry Smith 
23d9274352SBarry Smith    Level: beginner
24d9274352SBarry Smith 
258f6c3df8SBarry Smith .seealso: MatSetType(), Mat, MatSolverPackage, MatRegister()
2676bdecfbSBarry Smith J*/
2719fd82e9SBarry Smith typedef const char* MatType;
28273d9f13SBarry Smith #define MATSAME            "same"
295a11e1b2SBarry Smith #define MATMAIJ            "maij"
30273d9f13SBarry Smith #define MATSEQMAIJ         "seqmaij"
31273d9f13SBarry Smith #define MATMPIMAIJ         "mpimaij"
32273d9f13SBarry Smith #define MATIS              "is"
335a11e1b2SBarry Smith #define MATAIJ             "aij"
34273d9f13SBarry Smith #define MATSEQAIJ          "seqaij"
357d6a0e61SBarry Smith #define MATSEQAIJPTHREAD   "seqaijpthread"
36bf2c1783SBarry Smith #define MATAIJPTHREAD      "aijpthread"
37273d9f13SBarry Smith #define MATMPIAIJ          "mpiaij"
385a11e1b2SBarry Smith #define MATAIJCRL          "aijcrl"
395a11e1b2SBarry Smith #define MATSEQAIJCRL       "seqaijcrl"
405a11e1b2SBarry Smith #define MATMPIAIJCRL       "mpiaijcrl"
418154be41SBarry Smith #define MATAIJCUSP         "aijcusp"
428154be41SBarry Smith #define MATSEQAIJCUSP      "seqaijcusp"
438154be41SBarry Smith #define MATMPIAIJCUSP      "mpiaijcusp"
449ae82921SPaul Mullowney #define MATAIJCUSPARSE     "aijcusparse"
459ae82921SPaul Mullowney #define MATSEQAIJCUSPARSE  "seqaijcusparse"
469ae82921SPaul Mullowney #define MATMPIAIJCUSPARSE  "mpiaijcusparse"
47d67ff14aSKarl Rupp #define MATAIJVIENNACL     "aijviennacl"
48d67ff14aSKarl Rupp #define MATSEQAIJVIENNACL  "seqaijviennacl"
49d67ff14aSKarl Rupp #define MATMPIAIJVIENNACL  "mpiaijviennacl"
505a11e1b2SBarry Smith #define MATAIJPERM         "aijperm"
515a11e1b2SBarry Smith #define MATSEQAIJPERM      "seqaijperm"
525a11e1b2SBarry Smith #define MATMPIAIJPERM      "mpiaijperm"
53273d9f13SBarry Smith #define MATSHELL           "shell"
545a11e1b2SBarry Smith #define MATDENSE           "dense"
55209238afSKris Buschelman #define MATSEQDENSE        "seqdense"
56273d9f13SBarry Smith #define MATMPIDENSE        "mpidense"
57db31f6deSJed Brown #define MATELEMENTAL       "elemental"
585a11e1b2SBarry Smith #define MATBAIJ            "baij"
59273d9f13SBarry Smith #define MATSEQBAIJ         "seqbaij"
60273d9f13SBarry Smith #define MATMPIBAIJ         "mpibaij"
61273d9f13SBarry Smith #define MATMPIADJ          "mpiadj"
625a11e1b2SBarry Smith #define MATSBAIJ           "sbaij"
63273d9f13SBarry Smith #define MATSEQSBAIJ        "seqsbaij"
64273d9f13SBarry Smith #define MATMPISBAIJ        "mpisbaij"
65c0cdd4a1SDahai Guo #define MATSEQBSTRM        "seqbstrm"
66c0cdd4a1SDahai Guo #define MATMPIBSTRM        "mpibstrm"
67c0cdd4a1SDahai Guo #define MATBSTRM           "bstrm"
68aa5a9175SDahai Guo #define MATSEQSBSTRM       "seqsbstrm"
69aa5a9175SDahai Guo #define MATMPISBSTRM       "mpisbstrm"
70aa5a9175SDahai Guo #define MATSBSTRM          "sbstrm"
71cebc7f6cSBarry Smith #define MATDAAD            "daad"
72cebc7f6cSBarry Smith #define MATMFFD            "mffd"
73c8a8475eSBarry Smith #define MATNORMAL          "normal"
74c5e4d11fSDmitry Karpeev #define MATNORMALHERMITIAN "normalh"
75ab92ecdeSBarry Smith #define MATLRC             "lrc"
762a6744ebSBarry Smith #define MATSCATTER         "scatter"
77421e10b8SBarry Smith #define MATBLOCKMAT        "blockmat"
78793850ffSBarry Smith #define MATCOMPOSITE       "composite"
791f8c7532SHong Zhang #define MATFFT             "fft"
801f8c7532SHong Zhang #define MATFFTW            "fftw"
81e133240eSMatthew G Knepley #define MATSEQCUFFT        "seqcufft"
82557cca28SSatish Balay #define MATTRANSPOSEMAT    "transpose"
8372ca8751SBarry Smith #define MATSCHURCOMPLEMENT "schurcomplement"
841d6018f0SLisandro Dalcin #define MATPYTHON          "python"
85f91d8e95SBarry Smith #define MATHYPRESTRUCT     "hyprestruct"
86a9e6138eSGlenn Hammond #define MATHYPRESSTRUCT    "hypresstruct"
87ee1cef2cSJed Brown #define MATSUBMATRIX       "submatrix"
882c0dbf93SJed Brown #define MATLOCALREF        "localref"
89d8588912SDave May #define MATNEST            "nest"
90773941d6SBarry Smith 
9176bdecfbSBarry Smith /*J
92c7393fdbSBarry Smith     MatSolverPackage - String with the name of a PETSc matrix solver type.
939989ab13SBarry Smith 
94c2b89b5dSBarry Smith     For example: "petsc" indicates what PETSc provides, "superlu_dist" the parallel SuperLU_DIST package etc
959989ab13SBarry Smith 
969989ab13SBarry Smith    Level: beginner
979989ab13SBarry Smith 
985c9eb25fSBarry Smith .seealso: MatGetFactor(), Mat, MatSetType(), MatType
9976bdecfbSBarry Smith J*/
100c7393fdbSBarry Smith #define MatSolverPackage char*
1012692d6eeSBarry Smith #define MATSOLVERSUPERLU      "superlu"
1022692d6eeSBarry Smith #define MATSOLVERSUPERLU_DIST "superlu_dist"
1032692d6eeSBarry Smith #define MATSOLVERUMFPACK      "umfpack"
10420db9a53SJed Brown #define MATSOLVERCHOLMOD      "cholmod"
1052692d6eeSBarry Smith #define MATSOLVERESSL         "essl"
1062692d6eeSBarry Smith #define MATSOLVERLUSOL        "lusol"
1072692d6eeSBarry Smith #define MATSOLVERMUMPS        "mumps"
108d615f992SSatish Balay #define MATSOLVERMKL_PARDISO  "mkl_pardiso"
109d305a81bSVasiliy Kozyrev #define MATSOLVERMKL_CPARDISO "mkl_cpardiso"
1102692d6eeSBarry Smith #define MATSOLVERPASTIX       "pastix"
1112692d6eeSBarry Smith #define MATSOLVERMATLAB       "matlab"
1122692d6eeSBarry Smith #define MATSOLVERPETSC        "petsc"
1132692d6eeSBarry Smith #define MATSOLVERBAS          "bas"
1149ae82921SPaul Mullowney #define MATSOLVERCUSPARSE     "cusparse"
115aa5a9175SDahai Guo #define MATSOLVERBSTRM        "bstrm"
116aa5a9175SDahai Guo #define MATSOLVERSBSTRM       "sbstrm"
11715767789SHong Zhang #define MATSOLVERELEMENTAL    "elemental"
11817f1a0eaSHong Zhang #define MATSOLVERCLIQUE       "clique"
11969e15a41SShri Abhyankar #define MATSOLVERKLU          "klu"
120c0cdd4a1SDahai Guo 
121b24902e0SBarry Smith /*E
122b24902e0SBarry Smith     MatFactorType - indicates what type of factorization is requested
123b24902e0SBarry Smith 
124b24902e0SBarry Smith     Level: beginner
125b24902e0SBarry Smith 
126af0996ceSBarry Smith    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
127b24902e0SBarry Smith 
128c7393fdbSBarry Smith .seealso: MatSolverPackage, MatGetFactor()
129b24902e0SBarry Smith E*/
130599ef60dSHong Zhang typedef enum {MAT_FACTOR_NONE, MAT_FACTOR_LU, MAT_FACTOR_CHOLESKY, MAT_FACTOR_ILU, MAT_FACTOR_ICC,MAT_FACTOR_ILUDT} MatFactorType;
131014dd563SJed Brown PETSC_EXTERN const char *const MatFactorTypes[];
132e92e720dSBarry Smith 
133014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactor(Mat,const MatSolverPackage,MatFactorType,Mat*);
134014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactorAvailable(Mat,const MatSolverPackage,MatFactorType,PetscBool *);
135014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFactorGetSolverPackage(Mat,const MatSolverPackage*);
136014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactorType(Mat,MatFactorType*);
13718713533SBarry Smith PETSC_EXTERN PetscErrorCode MatSolverPackageRegister(const MatSolverPackage,const MatType,MatFactorType,PetscErrorCode(*)(Mat,MatFactorType,Mat*));
13842c9c57cSBarry Smith PETSC_EXTERN PetscErrorCode MatSolverPackageGet(const MatSolverPackage,const MatType,MatFactorType,PetscBool*,PetscBool*,PetscErrorCode (**)(Mat,MatFactorType,Mat*));
1399989ab13SBarry Smith 
140c06d978dSMatthew Knepley /* Logging support */
1410700a824SBarry Smith #define    MAT_FILE_CLASSID 1211216    /* used to indicate matrices in binary files */
142014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_CLASSID;
143335efc43SPeter Brune PETSC_EXTERN PetscClassId MAT_COLORING_CLASSID;
144014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_FDCOLORING_CLASSID;
145014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_TRANSPOSECOLORING_CLASSID;
146014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_PARTITIONING_CLASSID;
147014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_COARSEN_CLASSID;
148014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_NULLSPACE_CLASSID;
149014dd563SJed Brown PETSC_EXTERN PetscClassId MATMFFD_CLASSID;
150c06d978dSMatthew Knepley 
151ceb03754SKris Buschelman /*E
152ceb03754SKris Buschelman     MatReuse - Indicates if matrices obtained from a previous call to MatGetSubMatrices()
153d6eff37eSBarry Smith      or MatGetSubMatrix() are to be reused to store the new matrix values. For MatConvert() is used to indicate
154d6eff37eSBarry Smith      that the input matrix is to be replaced with the converted matrix.
155ceb03754SKris Buschelman 
156ceb03754SKris Buschelman     Level: beginner
157ceb03754SKris Buschelman 
158af0996ceSBarry Smith    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
159ceb03754SKris Buschelman 
1600c451bc4SBarry Smith .seealso: MatGetSubMatrices(), MatGetSubMatrix(), MatDestroyMatrices(), MatConvert()
161ceb03754SKris Buschelman E*/
162dfe085dbSJed Brown typedef enum {MAT_INITIAL_MATRIX,MAT_REUSE_MATRIX,MAT_IGNORE_MATRIX} MatReuse;
163ceb03754SKris Buschelman 
1645494a064SHong Zhang /*E
1655494a064SHong Zhang     MatGetSubMatrixOption - Indicates if matrices obtained from a call to MatGetSubMatrices()
166829201f2SHong Zhang      include the matrix values. Currently it is only used by MatGetSeqNonzerostructure().
1675494a064SHong Zhang 
1685494a064SHong Zhang     Level: beginner
1695494a064SHong Zhang 
170829201f2SHong Zhang .seealso: MatGetSeqNonzerostructure()
1715494a064SHong Zhang E*/
1725494a064SHong Zhang typedef enum {MAT_DO_NOT_GET_VALUES,MAT_GET_VALUES} MatGetSubMatrixOption;
1735494a064SHong Zhang 
174607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatInitializePackage(void);
175c06d978dSMatthew Knepley 
176014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreate(MPI_Comm,Mat*);
177014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetSizes(Mat,PetscInt,PetscInt,PetscInt,PetscInt);
17819fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatSetType(Mat,MatType);
179014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetFromOptions(Mat);
180685405a1SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatViewFromOptions(Mat A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);}
181bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatRegister(const char[],PetscErrorCode(*)(Mat));
182014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRegisterBaseName(const char[],const char[],const char[]);
183014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetOptionsPrefix(Mat,const char[]);
184014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAppendOptionsPrefix(Mat,const char[]);
185014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOptionsPrefix(Mat,const char*[]);
186422a814eSBarry Smith PETSC_EXTERN PetscErrorCode MatSetErrorIfFPE(Mat,PetscBool);
187f69a0ea3SMatthew Knepley 
188140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatList;
189140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatColoringList;
190140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatPartitioningList;
191140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatCoarsenList;
19228988994SBarry Smith 
1933b224e63SBarry Smith /*E
194aa6c7ce3SBarry Smith     MatStructure - Indicates if two matrices have the same nonzero structure
1953b224e63SBarry Smith 
1963b224e63SBarry Smith     Level: beginner
1973b224e63SBarry Smith 
198af0996ceSBarry Smith    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
1993b224e63SBarry Smith 
200aa6c7ce3SBarry Smith .seealso: MatCopy(), MatAXPY()
2013b224e63SBarry Smith E*/
202aa6c7ce3SBarry Smith typedef enum {DIFFERENT_NONZERO_PATTERN,SUBSET_NONZERO_PATTERN,SAME_NONZERO_PATTERN} MatStructure;
2033b224e63SBarry Smith 
204014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqDense(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*);
205014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateDense(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*);
206014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
207014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
208014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *);
209014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithSplitArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],PetscInt[],PetscInt[],PetscScalar[],Mat*);
2108d7a6e47SBarry Smith 
211014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
212014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
213014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat*);
214d21a29f3SJed Brown 
215014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAdj(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscInt[],Mat*);
216014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
217d21a29f3SJed Brown 
218014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
219014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPISBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *);
22038f409ebSLisandro Dalcin PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
221014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
2224cce697cSJed Brown PETSC_EXTERN PetscErrorCode MatXAIJSetPreallocation(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscInt[],const PetscInt[]);
223dfb205c3SBarry Smith 
224014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateShell(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,void *,Mat*);
225014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateNormal(Mat,Mat*);
226c5e4d11fSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatCreateNormalHermitian(Mat,Mat*);
227014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateLRC(Mat,Mat,Mat,Mat*);
228c5e4d11fSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatCreateIS(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,ISLocalToGlobalMapping,ISLocalToGlobalMapping,Mat*);
229014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
230014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
231c0cdd4a1SDahai Guo 
232014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
233014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
234014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqSBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
235014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPISBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
236c0cdd4a1SDahai Guo 
237014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateScatter(MPI_Comm,VecScatter,Mat*);
238014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScatterSetVecScatter(Mat,VecScatter);
239014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScatterGetVecScatter(Mat,VecScatter*);
240014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateBlockMat(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,Mat*);
241014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeAddMat(Mat,Mat);
242014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeMerge(Mat);
243014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateComposite(MPI_Comm,PetscInt,const Mat*,Mat*);
2446d7c1e57SBarry Smith typedef enum {MAT_COMPOSITE_ADDITIVE,MAT_COMPOSITE_MULTIPLICATIVE} MatCompositeType;
245014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeSetType(Mat,MatCompositeType);
2466d7c1e57SBarry Smith 
24719fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatCreateFFT(MPI_Comm,PetscInt,const PetscInt[],MatType,Mat*);
248014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqCUFFT(MPI_Comm,PetscInt,const PetscInt[],Mat*);
249dedccee8SHong Zhang 
250014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateTranspose(Mat,Mat*);
251d0de2241SAndrew Spott PETSC_EXTERN PetscErrorCode MatCreateHermitianTranspose(Mat,Mat*);
252014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSubMatrix(Mat,IS,IS,Mat*);
253014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSubMatrixUpdate(Mat,Mat,IS,IS);
254014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateLocalRef(Mat,IS,IS,Mat*);
2551472f72bSBarry Smith 
256014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPythonSetType(Mat,const char[]);
2571d6018f0SLisandro Dalcin 
258014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetUp(Mat);
259014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroy(Mat*);
260e56f5c9eSBarry Smith PETSC_EXTERN PetscErrorCode MatGetNonzeroState(Mat,PetscObjectState*);
26121c89e3eSBarry Smith 
262014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatConjugate(Mat);
263014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRealPart(Mat);
264014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatImaginaryPart(Mat);
265014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetDiagonalBlock(Mat,Mat*);
266014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetTrace(Mat,PetscScalar*);
267713ccfa9SJed Brown PETSC_EXTERN PetscErrorCode MatInvertBlockDiagonal(Mat,const PetscScalar **);
26899cafbc1SBarry Smith 
2698ed539a5SBarry Smith /* ------------------------------------------------------------*/
270014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
271014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
272014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesRow(Mat,PetscInt,const PetscScalar[]);
273014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesRowLocal(Mat,PetscInt,const PetscScalar[]);
274014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBatch(Mat,PetscInt,PetscInt,PetscInt[],const PetscScalar[]);
27573a71a0fSBarry Smith PETSC_EXTERN PetscErrorCode MatSetRandom(Mat,PetscRandom);
27684cb2905SBarry Smith 
2772ef4de8bSBarry Smith /*S
2782ef4de8bSBarry Smith      MatStencil - Data structure (C struct) for storing information about a single row or
279be479b30SJed Brown         column of a matrix as indexed on an associated grid.
280be479b30SJed Brown 
281be479b30SJed Brown    Fortran usage is different, see MatSetValuesStencil() for details.
2822ef4de8bSBarry Smith 
2832ef4de8bSBarry Smith    Level: beginner
2842ef4de8bSBarry Smith 
2852ef4de8bSBarry Smith   Concepts: matrix; linear operator
2862ef4de8bSBarry Smith 
287be479b30SJed Brown .seealso:  MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockedStencil()
2882ef4de8bSBarry Smith S*/
289435da068SBarry Smith typedef struct {
290c1ac3661SBarry Smith   PetscInt k,j,i,c;
291435da068SBarry Smith } MatStencil;
2922ef4de8bSBarry Smith 
293014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
294014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
295014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt);
296435da068SBarry Smith 
297014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetColoring(Mat,ISColoring);
298014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesAdifor(Mat,PetscInt,void*);
2993a7fca6bSBarry Smith 
300d91e6319SBarry Smith /*E
301d91e6319SBarry Smith     MatAssemblyType - Indicates if the matrix is now to be used, or if you plan
302d91e6319SBarry Smith      to continue to add values to it
303d91e6319SBarry Smith 
304d91e6319SBarry Smith     Level: beginner
305d91e6319SBarry Smith 
306d91e6319SBarry Smith .seealso: MatAssemblyBegin(), MatAssemblyEnd()
307d91e6319SBarry Smith E*/
3086d4a8577SBarry Smith typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType;
309014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssemblyBegin(Mat,MatAssemblyType);
310014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssemblyEnd(Mat,MatAssemblyType);
311014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssembled(Mat,PetscBool *);
3124f9c727eSBarry Smith 
3131d73ed98SBarry Smith 
31430de9b25SBarry Smith 
315d91e6319SBarry Smith /*E
316d91e6319SBarry Smith     MatOption - Options that may be set for a matrix and its behavior or storage
317d91e6319SBarry Smith 
318d91e6319SBarry Smith     Level: beginner
319d91e6319SBarry Smith 
320af0996ceSBarry Smith    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
321d91e6319SBarry Smith 
322382164b0SBarry Smith    Developer Notes: Entries that are negative need not be called collectively by all processes.
323382164b0SBarry Smith 
324d91e6319SBarry Smith .seealso: MatSetOption()
325d91e6319SBarry Smith E*/
326c5e4d11fSDmitry Karpeev typedef enum {MAT_OPTION_MIN = -3,
327c5e4d11fSDmitry Karpeev               MAT_UNUSED_NONZERO_LOCATION_ERR = -2,
32892d4d306SBarry Smith               MAT_ROW_ORIENTED = -1,
329382164b0SBarry Smith               MAT_SYMMETRIC = 1,
330382164b0SBarry Smith               MAT_STRUCTURALLY_SYMMETRIC = 2,
331382164b0SBarry Smith               MAT_NEW_DIAGONALS = 3,
332382164b0SBarry Smith               MAT_IGNORE_OFF_PROC_ENTRIES = 4,
333382164b0SBarry Smith               MAT_USE_HASH_TABLE = 5,
334382164b0SBarry Smith               MAT_KEEP_NONZERO_PATTERN = 6,
335382164b0SBarry Smith               MAT_IGNORE_ZERO_ENTRIES = 7,
336382164b0SBarry Smith               MAT_USE_INODES = 8,
337382164b0SBarry Smith               MAT_HERMITIAN = 9,
338382164b0SBarry Smith               MAT_SYMMETRY_ETERNAL = 10,
339c5e4d11fSDmitry Karpeev               MAT_NEW_NONZERO_LOCATION_ERR = 11,
340382164b0SBarry Smith               MAT_IGNORE_LOWER_TRIANGULAR = 12,
341382164b0SBarry Smith               MAT_ERROR_LOWER_TRIANGULAR = 13,
342382164b0SBarry Smith               MAT_GETROW_UPPERTRIANGULAR = 14,
343382164b0SBarry Smith               MAT_SPD = 15,
34492d4d306SBarry Smith               MAT_NO_OFF_PROC_ZERO_ROWS = 16,
34592d4d306SBarry Smith               MAT_NO_OFF_PROC_ENTRIES = 17,
34692d4d306SBarry Smith               MAT_NEW_NONZERO_LOCATIONS = 18,
347c5e4d11fSDmitry Karpeev               MAT_NEW_NONZERO_ALLOCATION_ERR = 19,
348c5e4d11fSDmitry Karpeev               MAT_SUBSET_OFF_PROC_ENTRIES = 20,
349c5e4d11fSDmitry Karpeev               MAT_OPTION_MAX = 21} MatOption;
350e057df02SPaul Mullowney 
351014dd563SJed Brown PETSC_EXTERN const char *MatOptions[];
352014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetOption(Mat,MatOption,PetscBool);
3537d68702bSBarry Smith PETSC_EXTERN PetscErrorCode MatGetOption(Mat,MatOption,PetscBool*);
35419fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatGetType(Mat,MatType*);
35584cb2905SBarry Smith 
356014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]);
357014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
358014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
359014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowUpperTriangular(Mat);
360014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRowUpperTriangular(Mat);
361014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
362014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
363014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnVector(Mat,Vec,PetscInt);
3648c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJGetArray(Mat,PetscScalar *[]);
3658c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJRestoreArray(Mat,PetscScalar *[]);
36621e72a00SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJGetMaxRowNonzeros(Mat,PetscInt*);
367bd04181cSBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJSetValuesLocalFast(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
3688c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatDenseGetArray(Mat,PetscScalar *[]);
3698c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatDenseRestoreArray(Mat,PetscScalar *[]);
370014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBlockSize(Mat,PetscInt *);
371014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSize(Mat,PetscInt);
372014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBlockSizes(Mat,PetscInt *,PetscInt *);
373014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSizes(Mat,PetscInt,PetscInt);
37433d57670SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSizesFromMats(Mat,Mat,Mat);
375014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNThreads(Mat,PetscInt);
376014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNThreads(Mat,PetscInt*);
3771620fd73SBarry Smith 
378014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMult(Mat,Vec,Vec);
379014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultDiagonalBlock(Mat,Vec,Vec);
380014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultAdd(Mat,Vec,Vec,Vec);
381014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTranspose(Mat,Vec,Vec);
382014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultHermitianTranspose(Mat,Vec,Vec);
383014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsTranspose(Mat,Mat,PetscReal,PetscBool *);
384014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitianTranspose(Mat,Mat,PetscReal,PetscBool *);
385014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeAdd(Mat,Vec,Vec,Vec);
386014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultHermitianTransposeAdd(Mat,Vec,Vec,Vec);
387014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultConstrained(Mat,Vec,Vec);
388014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeConstrained(Mat,Vec,Vec);
389014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatSolve(Mat,Mat,Mat);
390f9426fe0SMark Adams PETSC_EXTERN PetscErrorCode MatResidual(Mat,Vec,Vec,Vec);
391f5edf698SHong Zhang 
392d91e6319SBarry Smith /*E
393d91e6319SBarry Smith     MatDuplicateOption - Indicates if a duplicated sparse matrix should have
394d91e6319SBarry Smith   its numerical values copied over or just its nonzero structure.
395d91e6319SBarry Smith 
396d91e6319SBarry Smith     Level: beginner
397d91e6319SBarry Smith 
398af0996ceSBarry Smith    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
399d91e6319SBarry Smith 
40070dcbbb9SBarry Smith $   MAT_SHARE_NONZERO_PATTERN - the i and j arrays in the new matrix will be shared with the original matrix
40170dcbbb9SBarry Smith $                               this also triggers the MAT_DO_NOT_COPY_VALUES option. This is used when you
40270dcbbb9SBarry Smith $                               have several matrices with the same nonzero pattern.
40370dcbbb9SBarry Smith 
404d91e6319SBarry Smith .seealso: MatDuplicate()
405d91e6319SBarry Smith E*/
40670dcbbb9SBarry Smith typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES,MAT_SHARE_NONZERO_PATTERN} MatDuplicateOption;
4072e8a6d31SBarry Smith 
40819fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatConvert(Mat,MatType,MatReuse,Mat*);
409014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDuplicate(Mat,MatDuplicateOption,Mat*);
41094a9d846SBarry Smith 
411cb5b572fSBarry Smith 
412014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCopy(Mat,Mat,MatStructure);
413014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatView(Mat,PetscViewer);
414014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsSymmetric(Mat,PetscReal,PetscBool *);
415014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsStructurallySymmetric(Mat,PetscBool *);
416014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitian(Mat,PetscReal,PetscBool *);
417014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsSymmetricKnown(Mat,PetscBool *,PetscBool *);
418014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitianKnown(Mat,PetscBool *,PetscBool *);
419014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMissingDiagonal(Mat,PetscBool  *,PetscInt *);
420014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLoad(Mat, PetscViewer);
4217b80b807SBarry Smith 
4221a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool  *);
4231a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool  *);
4241a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool  *);
4251a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool  *);
426d4fbbf0eSBarry Smith 
427d91e6319SBarry Smith /*S
428d91e6319SBarry Smith      MatInfo - Context of matrix information, used with MatGetInfo()
429d91e6319SBarry Smith 
430d91e6319SBarry Smith    In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE
431d91e6319SBarry Smith 
432d91e6319SBarry Smith    Level: intermediate
433d91e6319SBarry Smith 
434d91e6319SBarry Smith   Concepts: matrix^nonzero information
435d91e6319SBarry Smith 
436d9274352SBarry Smith .seealso:  MatGetInfo(), MatInfoType
437d91e6319SBarry Smith S*/
4384e220ebcSLois Curfman McInnes typedef struct {
439b0a32e0cSBarry Smith   PetscLogDouble block_size;                         /* block size */
440b0a32e0cSBarry Smith   PetscLogDouble nz_allocated,nz_used,nz_unneeded;   /* number of nonzeros */
441b0a32e0cSBarry Smith   PetscLogDouble memory;                             /* memory allocated */
442b0a32e0cSBarry Smith   PetscLogDouble assemblies;                         /* number of matrix assemblies called */
443b0a32e0cSBarry Smith   PetscLogDouble mallocs;                            /* number of mallocs during MatSetValues() */
444b0a32e0cSBarry Smith   PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */
445b0a32e0cSBarry Smith   PetscLogDouble factor_mallocs;                     /* number of mallocs during factorization */
4464e220ebcSLois Curfman McInnes } MatInfo;
4474e220ebcSLois Curfman McInnes 
448d9274352SBarry Smith /*E
449d9274352SBarry Smith     MatInfoType - Indicates if you want information about the local part of the matrix,
450d9274352SBarry Smith      the entire parallel matrix or the maximum over all the local parts.
451d9274352SBarry Smith 
452d9274352SBarry Smith     Level: beginner
453d9274352SBarry Smith 
454af0996ceSBarry Smith    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
455d9274352SBarry Smith 
456d9274352SBarry Smith .seealso: MatGetInfo(), MatInfo
457d9274352SBarry Smith E*/
4587b80b807SBarry Smith typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType;
459014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetInfo(Mat,MatInfoType,MatInfo*);
460014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetDiagonal(Mat,Vec);
461014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMax(Mat,Vec,PetscInt[]);
462014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMin(Mat,Vec,PetscInt[]);
463014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMaxAbs(Mat,Vec,PetscInt[]);
464014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMinAbs(Mat,Vec,PetscInt[]);
465014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowSum(Mat,Vec);
466014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTranspose(Mat,MatReuse,Mat*);
467014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatHermitianTranspose(Mat,MatReuse,Mat*);
468014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPermute(Mat,IS,IS,Mat *);
469014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalScale(Mat,Vec,Vec);
470014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalSet(Mat,Vec,InsertMode);
471014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatEqual(Mat,Mat,PetscBool *);
472014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultEqual(Mat,Mat,PetscInt,PetscBool *);
473014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultAddEqual(Mat,Mat,PetscInt,PetscBool *);
474014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeEqual(Mat,Mat,PetscInt,PetscBool *);
475014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscBool *);
4767b80b807SBarry Smith 
477014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNorm(Mat,NormType,PetscReal *);
478014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnNorms(Mat,NormType,PetscReal *);
479014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroEntries(Mat);
480014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
481014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsIS(Mat,IS,PetscScalar,Vec,Vec);
482014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsStencil(Mat,PetscInt,const MatStencil [],PetscScalar,Vec,Vec);
483014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsStencil(Mat,PetscInt,const MatStencil[],PetscScalar,Vec,Vec);
484014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumns(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
485014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsIS(Mat,IS,PetscScalar,Vec,Vec);
4867b80b807SBarry Smith 
487014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSize(Mat,PetscInt*,PetscInt*);
488014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalSize(Mat,PetscInt*,PetscInt*);
489014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRange(Mat,PetscInt*,PetscInt*);
490014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRanges(Mat,const PetscInt**);
491014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRangeColumn(Mat,PetscInt*,PetscInt*);
492014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRangesColumn(Mat,const PetscInt**);
493566876eaSJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipIS(Mat,IS*,IS*);
4945ef9f2a5SBarry Smith 
495014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
49653dd7562SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatGetSubMatricesMPI(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
497014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroyMatrices(PetscInt,Mat *[]);
498014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSubMatrix(Mat,IS,IS,MatReuse,Mat *);
499014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalSubMatrix(Mat,IS,IS,Mat*);
500014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreLocalSubMatrix(Mat,IS,IS,Mat*);
501014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSeqNonzeroStructure(Mat,Mat*);
502014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroySeqNonzeroStructure(Mat*);
5037b80b807SBarry Smith 
504014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJ(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*);
505014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJSymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*);
506014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJNumeric(Mat,Mat);
507014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMat(Mat,MatReuse,Mat*);
508014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*);
509014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,Mat*);
510014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetGhosts(Mat, PetscInt *,const PetscInt *[]);
5118efafbd8SBarry Smith 
512014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt);
513aa1e27eaSFande Kong PETSC_EXTERN PetscErrorCode MatIncreaseOverlapSplit(Mat mat,PetscInt n,IS is[],PetscInt ov);
514d2b2a242SBarry Smith PETSC_EXTERN PetscErrorCode MatMPIAIJSetUseScalableIncreaseOverlap(Mat,PetscBool);
5157b80b807SBarry Smith 
516014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*);
517014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMultSymbolic(Mat,Mat,PetscReal,Mat*);
518014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMultNumeric(Mat,Mat,Mat);
51922440eb1SKris Buschelman 
5207bab7c10SHong Zhang PETSC_EXTERN PetscErrorCode MatMatMatMult(Mat,Mat,Mat,MatReuse,PetscReal,Mat*);
5217bab7c10SHong Zhang PETSC_EXTERN PetscErrorCode MatMatMatMultSymbolic(Mat,Mat,Mat,PetscReal,Mat*);
5227bab7c10SHong Zhang PETSC_EXTERN PetscErrorCode MatMatMatMultNumeric(Mat,Mat,Mat,Mat);
5237bab7c10SHong Zhang 
524014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*);
525014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAPSymbolic(Mat,Mat,PetscReal,Mat*);
526014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAPNumeric(Mat,Mat,Mat);
527014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARt(Mat,Mat,MatReuse,PetscReal,Mat*);
528014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARtSymbolic(Mat,Mat,PetscReal,Mat*);
529014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARtNumeric(Mat,Mat,Mat);
53022440eb1SKris Buschelman 
531014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeMatMult(Mat,Mat,MatReuse,PetscReal,Mat*);
532014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposetMatMultSymbolic(Mat,Mat,PetscReal,Mat*);
533014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposetMatMultNumeric(Mat,Mat,Mat);
534014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMult(Mat,Mat,MatReuse,PetscReal,Mat*);
535014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMultSymbolic(Mat,Mat,PetscReal,Mat*);
536014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMultNumeric(Mat,Mat,Mat);
537bc011b1eSHong Zhang 
538014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAXPY(Mat,PetscScalar,Mat,MatStructure);
539014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAYPX(Mat,PetscScalar,Mat,MatStructure);
5401c741599SHong Zhang 
541014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScale(Mat,PetscScalar);
542014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShift(Mat,PetscScalar);
5437b80b807SBarry Smith 
544014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping);
545014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*);
546a93ff8c4SPeter Brune PETSC_EXTERN PetscErrorCode MatGetLayouts(Mat,PetscLayout*,PetscLayout*);
547014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
548014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsLocalIS(Mat,IS,PetscScalar,Vec,Vec);
549014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
550014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocalIS(Mat,IS,PetscScalar,Vec,Vec);
551014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
552014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
553052efed2SBarry Smith 
554014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStashSetInitialSize(Mat,PetscInt,PetscInt);
555014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
55690f02eecSBarry Smith 
557014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInterpolate(Mat,Vec,Vec);
558014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInterpolateAdd(Mat,Vec,Vec,Vec);
559014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestrict(Mat,Vec,Vec);
5602a7a6963SBarry Smith PETSC_EXTERN PetscErrorCode MatCreateVecs(Mat,Vec*,Vec*);
56184cc2a5aSBarry Smith PETSC_DEPRECATED("Use MatCreateVecs()") PETSC_STATIC_INLINE PetscErrorCode MatGetVecs(Mat mat,Vec *x,Vec *y) {return MatCreateVecs(mat,x,y);}
56253cd1579SHong Zhang PETSC_EXTERN PetscErrorCode MatCreateRedundantMatrix(Mat,PetscInt,MPI_Comm,MatReuse,Mat*);
563014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetMultiProcBlock(Mat,MPI_Comm,MatReuse,Mat*);
564014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFindZeroDiagonals(Mat,IS*);
5653a062f41SBarry Smith PETSC_EXTERN PetscErrorCode MatFindOffBlockDiagonalEntries(Mat,IS*);
5669c8f2541SHong Zhang PETSC_EXTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
567bd481603SBarry Smith 
568bd481603SBarry Smith /*MC
5691620fd73SBarry Smith    MatSetValue - Set a single entry into a matrix.
5701620fd73SBarry Smith 
5711620fd73SBarry Smith    Not collective
5721620fd73SBarry Smith 
573a9834a7dSSatish Balay    Synopsis:
574a9834a7dSSatish Balay      #include <petscmat.h>
575a9834a7dSSatish Balay      PetscErrorCode MatSetValue(Mat m,PetscInt row,PetscInt col,PetscScalar value,InsertMode mode)
576a9834a7dSSatish Balay 
5771620fd73SBarry Smith    Input Parameters:
5781620fd73SBarry Smith +  m - the matrix
5791620fd73SBarry Smith .  row - the row location of the entry
5801620fd73SBarry Smith .  col - the column location of the entry
5811620fd73SBarry Smith .  value - the value to insert
5821620fd73SBarry Smith -  mode - either INSERT_VALUES or ADD_VALUES
5831620fd73SBarry Smith 
5841620fd73SBarry Smith    Notes:
5851620fd73SBarry Smith    For efficiency one should use MatSetValues() and set several or many
5861620fd73SBarry Smith    values simultaneously if possible.
5871620fd73SBarry Smith 
5881620fd73SBarry Smith    Level: beginner
5891620fd73SBarry Smith 
5901620fd73SBarry Smith .seealso: MatSetValues(), MatSetValueLocal()
5911620fd73SBarry Smith M*/
5921620fd73SBarry 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);}
5931620fd73SBarry Smith 
5941620fd73SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatGetValue(Mat v,PetscInt i,PetscInt j,PetscScalar *va) {return MatGetValues(v,1,&i,1,&j,va);}
5951620fd73SBarry Smith 
5961620fd73SBarry 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);}
5971620fd73SBarry Smith 
5981620fd73SBarry Smith /*MC
5990d2a7ff2SSatish Balay    MatPreallocateInitialize - Begins the block of code that will count the number of nonzeros per
600bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
601bd481603SBarry Smith 
602bd481603SBarry Smith    Synopsis:
603aaa7dc30SBarry Smith    #include <petscmat.h>
604c1ac3661SBarry Smith    PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)
605bd481603SBarry Smith 
606bd481603SBarry Smith    Collective on MPI_Comm
607bd481603SBarry Smith 
608bd481603SBarry Smith    Input Parameters:
609bd481603SBarry Smith +  comm - the communicator that will share the eventually allocated matrix
610859fcb39SBarry Smith .  nrows - the number of LOCAL rows in the matrix
611859fcb39SBarry Smith -  ncols - the number of LOCAL columns in the matrix
612bd481603SBarry Smith 
613bd481603SBarry Smith    Output Parameters:
614bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
615bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
616bd481603SBarry Smith 
617bd481603SBarry Smith    Level: intermediate
618bd481603SBarry Smith 
619bd481603SBarry Smith    Notes:
620a7f22e61SSatish Balay     See Users-Manual: ch_performance for more details.
621bd481603SBarry Smith 
6221d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
623bd481603SBarry Smith 
624bd481603SBarry Smith    Use MatPreallocateInitializeSymmetric() for symmetric matrices (MPISBAIJ matrices)
625bd481603SBarry Smith 
6261620fd73SBarry Smith    This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize().
6271620fd73SBarry Smith 
628bd481603SBarry Smith   Concepts: preallocation^Matrix
629bd481603SBarry Smith 
630d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(),
631d6e23781SBarry Smith           MatPreallocateInitializeSymmetric(), MatPreallocateSymmetricSetLocalBlock()
632bd481603SBarry Smith M*/
633c4f061fbSSatish Balay #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \
6347c922b88SBarry Smith { \
635a89bc211SBarry Smith   PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__start,__end; \
6361795a4d1SJed Brown   _4_ierr = PetscCalloc2(__nrows,&dnz,__nrows,&onz);CHKERRQ(_4_ierr); \
6371795a4d1SJed Brown   __start = 0; __end = __start;                                         \
638c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ctmp;\
639a89bc211SBarry Smith   _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows;
6407c922b88SBarry Smith 
641bd481603SBarry Smith /*MC
642bd481603SBarry Smith    MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
643bd481603SBarry Smith        inserted using a local number of the rows and columns
644bd481603SBarry Smith 
645bd481603SBarry Smith    Synopsis:
646aaa7dc30SBarry Smith    #include <petscmat.h>
647c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
648bd481603SBarry Smith 
649bd481603SBarry Smith    Not Collective
650bd481603SBarry Smith 
651bd481603SBarry Smith    Input Parameters:
652784ac674SJed Brown +  map - the row mapping from local numbering to global numbering
653bd481603SBarry Smith .  nrows - the number of rows indicated
6541d73ed98SBarry Smith .  rows - the indices of the rows
655784ac674SJed Brown .  cmap - the column mapping from local to global numbering
656bd481603SBarry Smith .  ncols - the number of columns in the matrix
657bd481603SBarry Smith .  cols - the columns indicated
658bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
659bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
660bd481603SBarry Smith 
661bd481603SBarry Smith    Level: intermediate
662bd481603SBarry Smith 
663bd481603SBarry Smith    Notes:
664a7f22e61SSatish Balay     See Users-Manual: ch_performance for more details.
665bd481603SBarry Smith 
6661d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
667bd481603SBarry Smith 
668bd481603SBarry Smith   Concepts: preallocation^Matrix
669bd481603SBarry Smith 
670d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(),
671d6e23781SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock()
672bd481603SBarry Smith M*/
673784ac674SJed Brown #define MatPreallocateSetLocal(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \
674c4f061fbSSatish Balay {\
675c1ac3661SBarry Smith   PetscInt __l;\
676784ac674SJed Brown   _4_ierr = ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\
677784ac674SJed Brown   _4_ierr = ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\
678c4f061fbSSatish Balay   for (__l=0;__l<nrows;__l++) {\
679ef66eb69SBarry Smith     _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
680c4f061fbSSatish Balay   }\
681c4f061fbSSatish Balay }
682c4f061fbSSatish Balay 
683bd481603SBarry Smith /*MC
684d6e23781SBarry Smith    MatPreallocateSetLocalBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be
685bd481603SBarry Smith        inserted using a local number of the rows and columns
686bd481603SBarry Smith 
687bd481603SBarry Smith    Synopsis:
688aaa7dc30SBarry Smith    #include <petscmat.h>
689d6e23781SBarry Smith    PetscErrorCode MatPreallocateSetLocalBlock(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
690d6e23781SBarry Smith 
691d6e23781SBarry Smith    Not Collective
692d6e23781SBarry Smith 
693d6e23781SBarry Smith    Input Parameters:
694d6e23781SBarry Smith +  map - the row mapping from local numbering to global numbering
695d6e23781SBarry Smith .  nrows - the number of rows indicated
696d6e23781SBarry Smith .  rows - the indices of the rows
697d6e23781SBarry Smith .  cmap - the column mapping from local to global numbering
698d6e23781SBarry Smith .  ncols - the number of columns in the matrix
699d6e23781SBarry Smith .  cols - the columns indicated
700d6e23781SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
701d6e23781SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
702d6e23781SBarry Smith 
703d6e23781SBarry Smith    Level: intermediate
704d6e23781SBarry Smith 
705d6e23781SBarry Smith    Notes:
706d6e23781SBarry Smith     See Users-Manual: ch_performance for more details.
707d6e23781SBarry Smith 
708d6e23781SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
709d6e23781SBarry Smith 
710d6e23781SBarry Smith   Concepts: preallocation^Matrix
711d6e23781SBarry Smith 
712d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(),
713d6e23781SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock()
714d6e23781SBarry Smith M*/
715d6e23781SBarry Smith #define MatPreallocateSetLocalBlock(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \
716d6e23781SBarry Smith {\
717d6e23781SBarry Smith   PetscInt __l;\
718d6e23781SBarry Smith   _4_ierr = ISLocalToGlobalMappingApplyBlock(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\
719d6e23781SBarry Smith   _4_ierr = ISLocalToGlobalMappingApplyBlock(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\
720d6e23781SBarry Smith   for (__l=0;__l<nrows;__l++) {\
721d6e23781SBarry Smith     _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
722d6e23781SBarry Smith   }\
723d6e23781SBarry Smith }
724d6e23781SBarry Smith 
725d6e23781SBarry Smith /*MC
726d6e23781SBarry Smith    MatPreallocateSymmetricSetLocalBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be
727d6e23781SBarry Smith        inserted using a local number of the rows and columns
728d6e23781SBarry Smith 
729d6e23781SBarry Smith    Synopsis:
730d6e23781SBarry Smith    #include <petscmat.h>
731d6e23781SBarry Smith    PetscErrorCode MatPreallocateSymmetricSetLocalBlock(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
732bd481603SBarry Smith 
733bd481603SBarry Smith    Not Collective
734bd481603SBarry Smith 
735bd481603SBarry Smith    Input Parameters:
736bd481603SBarry Smith +  map - the mapping between local numbering and global numbering
737bd481603SBarry Smith .  nrows - the number of rows indicated
7381d73ed98SBarry Smith .  rows - the indices of the rows
739bd481603SBarry Smith .  ncols - the number of columns in the matrix
740bd481603SBarry Smith .  cols - the columns indicated
741bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
742bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
743bd481603SBarry Smith 
744bd481603SBarry Smith    Level: intermediate
745bd481603SBarry Smith 
746bd481603SBarry Smith    Notes:
747a7f22e61SSatish Balay     See Users-Manual: ch_performance for more details.
748bd481603SBarry Smith 
749bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
750bd481603SBarry Smith 
751bd481603SBarry Smith   Concepts: preallocation^Matrix
752bd481603SBarry Smith 
753d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateInitialize(),
754d6e23781SBarry Smith           MatPreallocateInitialize(),  MatPreallocateSetLocal()
755bd481603SBarry Smith M*/
756d6e23781SBarry Smith #define MatPreallocateSymmetricSetLocalBlock(map,nrows,rows,ncols,cols,dnz,onz) 0;\
757d3d32019SBarry Smith {\
758c1ac3661SBarry Smith   PetscInt __l;\
759d6e23781SBarry Smith   _4_ierr = ISLocalToGlobalMappingApplyBlock(map,nrows,rows,rows);CHKERRQ(_4_ierr);\
760d6e23781SBarry Smith   _4_ierr = ISLocalToGlobalMappingApplyBlock(map,ncols,cols,cols);CHKERRQ(_4_ierr);\
761d3d32019SBarry Smith   for (__l=0;__l<nrows;__l++) {\
762d6e23781SBarry Smith     _4_ierr = MatPreallocateSymmetricSetBlock((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
763d3d32019SBarry Smith   }\
764d3d32019SBarry Smith }
765bd481603SBarry Smith /*MC
766bd481603SBarry Smith    MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
767bd481603SBarry Smith        inserted using a local number of the rows and columns
768bd481603SBarry Smith 
769bd481603SBarry Smith    Synopsis:
770aaa7dc30SBarry Smith    #include <petscmat.h>
771c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
772bd481603SBarry Smith 
773bd481603SBarry Smith    Not Collective
774bd481603SBarry Smith 
775bd481603SBarry Smith    Input Parameters:
77664075487SBarry Smith +  row - the row
777bd481603SBarry Smith .  ncols - the number of columns in the matrix
778a50f8bf6SHong Zhang -  cols - the columns indicated
779a50f8bf6SHong Zhang 
780a50f8bf6SHong Zhang    Output Parameters:
781a50f8bf6SHong Zhang +  dnz - the array that will be passed to the matrix preallocation routines
782bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
783bd481603SBarry Smith 
784bd481603SBarry Smith    Level: intermediate
785bd481603SBarry Smith 
786bd481603SBarry Smith    Notes:
787a7f22e61SSatish Balay     See Users-Manual: ch_performance for more details.
788bd481603SBarry Smith 
789bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
790bd481603SBarry Smith 
7911620fd73SBarry Smith    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
7921620fd73SBarry Smith 
793bd481603SBarry Smith   Concepts: preallocation^Matrix
794bd481603SBarry Smith 
795d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(),
796d6e23781SBarry Smith           MatPreallocateInitialize(), MatPreallocateSetLocal()
797bd481603SBarry Smith M*/
798c4f061fbSSatish Balay #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\
799c1ac3661SBarry Smith { PetscInt __i; \
800e32f2f54SBarry 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);\
801e32f2f54SBarry 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);\
8027c922b88SBarry Smith   for (__i=0; __i<nc; __i++) {\
80364075487SBarry Smith     if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \
8047cc688d7SBarry Smith     else dnz[row - __rstart]++;\
8057c922b88SBarry Smith   }\
8067c922b88SBarry Smith }
8077c922b88SBarry Smith 
808bd481603SBarry Smith /*MC
809d6e23781SBarry Smith    MatPreallocateSymmetricSetBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be
810bd481603SBarry Smith        inserted using a local number of the rows and columns
811bd481603SBarry Smith 
812bd481603SBarry Smith    Synopsis:
813aaa7dc30SBarry Smith    #include <petscmat.h>
814d6e23781SBarry Smith    PetscErrorCode MatPreallocateSymmetricSetBlock(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
815bd481603SBarry Smith 
816bd481603SBarry Smith    Not Collective
817bd481603SBarry Smith 
818bd481603SBarry Smith    Input Parameters:
819bd481603SBarry Smith +  nrows - the number of rows indicated
8201d73ed98SBarry Smith .  rows - the indices of the rows
821bd481603SBarry Smith .  ncols - the number of columns in the matrix
822bd481603SBarry Smith .  cols - the columns indicated
823bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
824bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
825bd481603SBarry Smith 
826bd481603SBarry Smith    Level: intermediate
827bd481603SBarry Smith 
828bd481603SBarry Smith    Notes:
829a7f22e61SSatish Balay     See Users-Manual: ch_performance for more details.
830bd481603SBarry Smith 
831bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
832bd481603SBarry Smith 
8331620fd73SBarry Smith    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
8341620fd73SBarry Smith 
835bd481603SBarry Smith   Concepts: preallocation^Matrix
836bd481603SBarry Smith 
837d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(),  MatPreallocateInitialize(),
838d6e23781SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocal()
839bd481603SBarry Smith M*/
840d6e23781SBarry Smith #define MatPreallocateSymmetricSetBlock(row,nc,cols,dnz,onz) 0;\
841c1ac3661SBarry Smith { PetscInt __i; \
842d3d32019SBarry Smith   for (__i=0; __i<nc; __i++) {\
843d3d32019SBarry Smith     if (cols[__i] >= __end) onz[row - __rstart]++; \
844d3d32019SBarry Smith     else if (cols[__i] >= row) dnz[row - __rstart]++;\
845d3d32019SBarry Smith   }\
846d3d32019SBarry Smith }
847d3d32019SBarry Smith 
848bd481603SBarry Smith /*MC
84916371a99SBarry Smith    MatPreallocateLocation -  An alternative to MatPreallocationSet() that puts the nonzero locations into the matrix if it exists
85016371a99SBarry Smith 
85116371a99SBarry Smith    Synopsis:
852aaa7dc30SBarry Smith    #include <petscmat.h>
85316371a99SBarry Smith    PetscErrorCode MatPreallocateLocations(Mat A,PetscInt row,PetscInt ncols,PetscInt *cols,PetscInt *dnz,PetscInt *onz)
85416371a99SBarry Smith 
85516371a99SBarry Smith    Not Collective
85616371a99SBarry Smith 
85716371a99SBarry Smith    Input Parameters:
85816371a99SBarry Smith .  A - matrix
85916371a99SBarry Smith .  row - row where values exist (must be local to this process)
86016371a99SBarry Smith .  ncols - number of columns
86116371a99SBarry Smith .  cols - columns with nonzeros
86216371a99SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
86316371a99SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
86416371a99SBarry Smith 
86516371a99SBarry Smith    Level: intermediate
86616371a99SBarry Smith 
86716371a99SBarry Smith    Notes:
868a7f22e61SSatish Balay     See Users-Manual: ch_performance for more details.
86916371a99SBarry Smith 
87016371a99SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
87116371a99SBarry Smith 
87216371a99SBarry Smith    This is a MACRO not a function because it uses a bunch of variables private to the MatPreallocation.... routines.
87316371a99SBarry Smith 
87416371a99SBarry Smith   Concepts: preallocation^Matrix
87516371a99SBarry Smith 
876d6e23781SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(),
877d6e23781SBarry Smith           MatPreallocateSymmetricSetLocalBlock()
87816371a99SBarry Smith M*/
8790298fd71SBarry 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);}
88016371a99SBarry Smith 
88116371a99SBarry Smith 
88216371a99SBarry Smith /*MC
8830d2a7ff2SSatish Balay    MatPreallocateFinalize - Ends the block of code that will count the number of nonzeros per
884bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
885bd481603SBarry Smith 
886bd481603SBarry Smith    Synopsis:
887aaa7dc30SBarry Smith    #include <petscmat.h>
888c1ac3661SBarry Smith    PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz)
889bd481603SBarry Smith 
890bd481603SBarry Smith    Collective on MPI_Comm
891bd481603SBarry Smith 
892bd481603SBarry Smith    Input Parameters:
89316371a99SBarry Smith +  dnz - the array that was be passed to the matrix preallocation routines
894bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
895bd481603SBarry Smith 
896bd481603SBarry Smith    Level: intermediate
897bd481603SBarry Smith 
898bd481603SBarry Smith    Notes:
899a7f22e61SSatish Balay     See Users-Manual: ch_performance for more details.
900bd481603SBarry Smith 
901bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
902bd481603SBarry Smith 
9031620fd73SBarry Smith    This is a MACRO not a function because it closes the { started in MatPreallocateInitialize().
9041620fd73SBarry Smith 
905bd481603SBarry Smith   Concepts: preallocation^Matrix
906bd481603SBarry Smith 
907d6e23781SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(),
908d6e23781SBarry Smith           MatPreallocateSymmetricSetLocalBlock()
909bd481603SBarry Smith M*/
910a89bc211SBarry Smith #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree2(dnz,onz);CHKERRQ(_4_ierr);}
9117c922b88SBarry Smith 
9127b80b807SBarry Smith /* Routines unique to particular data structures */
913014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellGetContext(Mat,void *);
914698d4c6aSKris Buschelman 
915014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInodeAdjustForInodes(Mat,IS*,IS*);
916014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *);
917698d4c6aSKris Buschelman 
918014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetColumnIndices(Mat,PetscInt[]);
919014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetColumnIndices(Mat,PetscInt[]);
920014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
921014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
922014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
923014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJFromTriple(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*,PetscInt,PetscBool);
9247b80b807SBarry Smith 
925a96a251dSBarry Smith #define MAT_SKIP_ALLOCATION -4
926a96a251dSBarry Smith 
927014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
928014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
929014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]);
930273d9f13SBarry Smith 
931014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
932014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
933014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
934014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []);
935014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
936014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]);
937014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
938014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]);
939014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIDenseSetPreallocation(Mat,PetscScalar[]);
940014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqDenseSetPreallocation(Mat,PetscScalar[]);
9419230625dSJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,const PetscInt*[]);
9429230625dSJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,const PetscInt*[]);
943014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAdjCreateNonemptySubcommMat(Mat,Mat*);
944273d9f13SBarry Smith 
9452e1947a5SStefano Zampini PETSC_EXTERN PetscErrorCode MatISSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
946014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqDenseSetLDA(Mat,PetscInt);
947014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDenseGetLocalMatrix(Mat,Mat*);
9481b807ce4Svictorle 
949014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStoreValues(Mat);
950014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRetrieveValues(Mat);
9512e8a6d31SBarry Smith 
952014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDAADSetCtx(Mat,void*);
9533a7fca6bSBarry Smith 
954014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFindNonzeroRows(Mat,IS*);
9557b80b807SBarry Smith /*
9567b80b807SBarry Smith   These routines are not usually accessed directly, rather solving is
95794b7f48cSBarry Smith   done through the KSP and PC interfaces.
9587b80b807SBarry Smith */
9597b80b807SBarry Smith 
96076bdecfbSBarry Smith /*J
9618f6c3df8SBarry Smith     MatOrderingType - String with the name of a PETSc matrix ordering
962d9274352SBarry Smith 
963d9274352SBarry Smith    Level: beginner
964d9274352SBarry Smith 
965d9274352SBarry Smith .seealso: MatGetOrdering()
96676bdecfbSBarry Smith J*/
96719fd82e9SBarry Smith typedef const char* MatOrderingType;
9682692d6eeSBarry Smith #define MATORDERINGNATURAL     "natural"
9692692d6eeSBarry Smith #define MATORDERINGND          "nd"
9702692d6eeSBarry Smith #define MATORDERING1WD         "1wd"
9712692d6eeSBarry Smith #define MATORDERINGRCM         "rcm"
9722692d6eeSBarry Smith #define MATORDERINGQMD         "qmd"
9732692d6eeSBarry Smith #define MATORDERINGROWLENGTH   "rowlength"
974510184a7SJed Brown #define MATORDERINGWBM         "wbm"
975fc1bef75SJed Brown #define MATORDERINGSPECTRAL    "spectral"
976312e372aSBarry Smith #define MATORDERINGAMD         "amd"            /* only works if UMFPACK is installed with PETSc */
977b12f92e5SBarry Smith 
97819fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatGetOrdering(Mat,MatOrderingType,IS*,IS*);
979140e18c1SBarry Smith PETSC_EXTERN PetscErrorCode MatGetOrderingList(PetscFunctionList*);
980bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatOrderingRegister(const char[],PetscErrorCode(*)(Mat,MatOrderingType,IS*,IS*));
981140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatOrderingList;
982d4fbbf0eSBarry Smith 
983014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS);
984fc1bef75SJed Brown PETSC_EXTERN PetscErrorCode MatCreateLaplacian(Mat,PetscReal,PetscBool,Mat*);
985a2ce50c7SBarry Smith 
986d91e6319SBarry Smith /*S
987d90ac83dSHong Zhang     MatFactorShiftType - Numeric Shift.
988d90ac83dSHong Zhang 
989d90ac83dSHong Zhang    Level: beginner
990d90ac83dSHong Zhang 
991d90ac83dSHong Zhang S*/
992d90ac83dSHong Zhang typedef enum {MAT_SHIFT_NONE,MAT_SHIFT_NONZERO,MAT_SHIFT_POSITIVE_DEFINITE,MAT_SHIFT_INBLOCKS} MatFactorShiftType;
9936a6fc655SJed Brown PETSC_EXTERN const char *const MatFactorShiftTypes[];
9945e9742b9SJed Brown PETSC_EXTERN const char *const MatFactorShiftTypesDetail[];
995d90ac83dSHong Zhang 
996d90ac83dSHong Zhang /*S
997422a814eSBarry Smith    MatFactorInfo - Data passed into the matrix factorization routines, and information about the resulting factorization
998b00f7748SHong Zhang 
99961cad860SBarry Smith    In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE, that is use
100061cad860SBarry Smith $     MatFactorInfo  info(MAT_FACTORINFO_SIZE)
1001b00f7748SHong Zhang 
100215e8a5b3SHong Zhang    Notes: These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC.
1003b00f7748SHong Zhang 
100461cad860SBarry Smith       You can use MatFactorInfoInitialize() to set default values.
100561cad860SBarry Smith 
1006b00f7748SHong Zhang    Level: developer
1007b00f7748SHong Zhang 
1008d7d82daaSBarry Smith .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(),
1009d7d82daaSBarry Smith           MatFactorInfoInitialize()
1010b00f7748SHong Zhang 
1011b00f7748SHong Zhang S*/
1012b00f7748SHong Zhang typedef struct {
101315e8a5b3SHong Zhang   PetscReal     diagonal_fill;  /* force diagonal to fill in if initially not filled */
101485317021SBarry Smith   PetscReal     usedt;
101515e8a5b3SHong Zhang   PetscReal     dt;             /* drop tolerance */
1016b00f7748SHong Zhang   PetscReal     dtcol;          /* tolerance for pivoting */
101715e8a5b3SHong Zhang   PetscReal     dtcount;        /* maximum nonzeros to be allowed per row */
101867e28bfeSBarry Smith   PetscReal     fill;           /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */
1019348344bbSBarry Smith   PetscReal     levels;         /* ICC/ILU(levels) */
1020bcd9e38bSBarry Smith   PetscReal     pivotinblocks;  /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0
1021bcd9e38bSBarry Smith                                    factorization may be faster if do not pivot */
102215e8a5b3SHong Zhang   PetscReal     zeropivot;      /* pivot is called zero if less than this */
1023f4db908eSBarry Smith   PetscReal     shifttype;      /* type of shift added to matrix factor to prevent zero pivots */
1024f4db908eSBarry Smith   PetscReal     shiftamount;     /* how large the shift is */
102515e8a5b3SHong Zhang } MatFactorInfo;
1026ffa6d0a5SLois Curfman McInnes 
1027014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFactorInfoInitialize(MatFactorInfo*);
1028014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCholeskyFactor(Mat,IS,const MatFactorInfo*);
1029014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCholeskyFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
1030014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCholeskyFactorNumeric(Mat,Mat,const MatFactorInfo*);
1031014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLUFactor(Mat,IS,IS,const MatFactorInfo*);
1032014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatILUFactor(Mat,IS,IS,const MatFactorInfo*);
1033014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
1034014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatILUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
1035014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatICCFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
1036014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatICCFactor(Mat,IS,const MatFactorInfo*);
1037014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLUFactorNumeric(Mat,Mat,const MatFactorInfo*);
1038014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*);
1039014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolve(Mat,Vec,Vec);
1040014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatForwardSolve(Mat,Vec,Vec);
1041014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatBackwardSolve(Mat,Vec,Vec);
1042014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveAdd(Mat,Vec,Vec,Vec);
1043014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveTranspose(Mat,Vec,Vec);
1044014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveTransposeAdd(Mat,Vec,Vec,Vec);
1045014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolves(Mat,Vecs,Vecs);
1046*8e7ba810SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorSetSchurIS(Mat,IS);
1047014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetUnfactored(Mat);
1048bb5a7306SBarry Smith 
1049d91e6319SBarry Smith /*E
1050d91e6319SBarry Smith     MatSORType - What type of (S)SOR to perform
1051bb1eb677SSatish Balay 
1052d91e6319SBarry Smith     Level: beginner
1053d91e6319SBarry Smith 
1054d9274352SBarry Smith    May be bitwise ORd together
1055d9274352SBarry Smith 
1056af0996ceSBarry Smith    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
1057d91e6319SBarry Smith 
10584e7234bfSBarry Smith    MatSORType may be bitwise ORd together, so do not change the numbers
10594e7234bfSBarry Smith 
106041f059aeSBarry Smith .seealso: MatSOR()
1061d91e6319SBarry Smith E*/
1062ee50ffe9SBarry Smith typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3,
1063ee50ffe9SBarry Smith               SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8,
1064ee50ffe9SBarry Smith               SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16,
106584cb2905SBarry Smith               SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType;
1066014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSOR(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
10678ed539a5SBarry Smith 
1068d4fbbf0eSBarry Smith /*
1069639f9d9dSBarry Smith     These routines are for efficiently computing Jacobians via finite differences.
1070639f9d9dSBarry Smith */
1071b12f92e5SBarry Smith 
10727086a01eSPeter Brune /*S
10737086a01eSPeter Brune      MatColoring - Object for managing the coloring of matrices.
10747086a01eSPeter Brune 
10757086a01eSPeter Brune    Level: beginner
10767086a01eSPeter Brune 
10777086a01eSPeter Brune   Concepts: matrix, coloring
10787086a01eSPeter Brune 
10797086a01eSPeter Brune .seealso:  MatFDColoringCreate() ISColoring MatFDColoring
10807086a01eSPeter Brune S*/
10817086a01eSPeter Brune typedef struct _p_MatColoring* MatColoring;
108276bdecfbSBarry Smith /*J
10838f6c3df8SBarry Smith     MatColoringType - String with the name of a PETSc matrix coloring
1084d9274352SBarry Smith 
1085d9274352SBarry Smith    Level: beginner
1086d9274352SBarry Smith 
10877086a01eSPeter Brune .seealso: MatColoringSetType(), MatColoring
108876bdecfbSBarry Smith J*/
10897086a01eSPeter Brune 
109019fd82e9SBarry Smith typedef const  char*           MatColoringType;
10915567d71aSPeter Brune #define MATCOLORINGJP      "jp"
1092b9acb617SPeter Brune #define MATCOLORINGPOWER   "power"
10932692d6eeSBarry Smith #define MATCOLORINGNATURAL "natural"
10942692d6eeSBarry Smith #define MATCOLORINGSL      "sl"
10952692d6eeSBarry Smith #define MATCOLORINGLF      "lf"
10962692d6eeSBarry Smith #define MATCOLORINGID      "id"
10978121bdceSPeter Brune #define MATCOLORINGGREEDY  "greedy"
1098b12f92e5SBarry Smith 
10998ac6da0aSPeter Brune /*E
11008ac6da0aSPeter Brune    MatColoringWeightType - Type of weight scheme
11018ac6da0aSPeter Brune 
11028ac6da0aSPeter Brune     Not Collective
11038ac6da0aSPeter Brune 
11048ac6da0aSPeter Brune +   MAT_COLORING_RANDOM  - Random weights
11058ac6da0aSPeter Brune .   MAT_COLORING_LEXICAL - Lexical weighting based upon global numbering.
11068ac6da0aSPeter Brune -   MAT_COLORING_LF      - Last-first weighting.
11078ac6da0aSPeter Brune 
11088ac6da0aSPeter Brune     Level: intermediate
11098ac6da0aSPeter Brune 
1110af0996ceSBarry Smith    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
11118ac6da0aSPeter Brune 
11128ac6da0aSPeter Brune .seealso: MatCUSPSetFormat(), MatCUSPFormatOperation
11138ac6da0aSPeter Brune E*/
1114301c140bSPeter Brune typedef enum {MAT_COLORING_WEIGHT_RANDOM,MAT_COLORING_WEIGHT_LEXICAL,MAT_COLORING_WEIGHT_LF,MAT_COLORING_WEIGHT_SL} MatColoringWeightType;
11158ac6da0aSPeter Brune 
1116335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringCreate(Mat,MatColoring*);
1117d7f2a307SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringGetDegrees(Mat,PetscInt,PetscInt*);
1118335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringDestroy(MatColoring*);
1119335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringView(MatColoring,PetscViewer);
1120335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetType(MatColoring,MatColoringType);
1121335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetFromOptions(MatColoring);
1122335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetDistance(MatColoring,PetscInt);
1123335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringGetDistance(MatColoring,PetscInt*);
1124335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetMaxColors(MatColoring,PetscInt);
1125335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringGetMaxColors(MatColoring,PetscInt*);
1126335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringApply(MatColoring,ISColoring*);
1127335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringRegister(const char[],PetscErrorCode(*)(MatColoring));
1128014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*);
11298ac6da0aSPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetWeightType(MatColoring,MatColoringWeightType);
1130c29971fcSPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetWeights(MatColoring,PetscReal*,PetscInt*);
11318ac6da0aSPeter Brune PETSC_EXTERN PetscErrorCode MatColoringCreateWeights(MatColoring,PetscReal **,PetscInt **lperm);
1132639f9d9dSBarry Smith 
1133d9274352SBarry Smith /*S
1134d9274352SBarry Smith      MatFDColoring - Object for computing a sparse Jacobian via finite differences
1135d9274352SBarry Smith         and coloring
1136639f9d9dSBarry Smith 
1137d9274352SBarry Smith    Level: beginner
1138d9274352SBarry Smith 
1139d9274352SBarry Smith   Concepts: coloring, sparse Jacobian, finite differences
1140d9274352SBarry Smith 
1141d9274352SBarry Smith .seealso:  MatFDColoringCreate()
1142d9274352SBarry Smith S*/
1143e2a1c21fSSatish Balay typedef struct _p_MatFDColoring* MatFDColoring;
1144639f9d9dSBarry Smith 
1145014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringCreate(Mat,ISColoring,MatFDColoring *);
1146014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringDestroy(MatFDColoring*);
1147014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringView(MatFDColoring,PetscViewer);
1148014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*);
1149014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**);
1150014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal);
1151014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetFromOptions(MatFDColoring);
1152d1e9a80fSBarry Smith PETSC_EXTERN PetscErrorCode MatFDColoringApply(Mat,MatFDColoring,Vec,void *);
1153014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetF(MatFDColoring,Vec);
1154014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,PetscInt*[]);
1155f86b9fbaSHong Zhang PETSC_EXTERN PetscErrorCode MatFDColoringSetUp(Mat,ISColoring,MatFDColoring);
1156f86b9fbaSHong Zhang PETSC_EXTERN PetscErrorCode MatFDColoringSetBlockSize(MatFDColoring,PetscInt,PetscInt);
1157f86b9fbaSHong Zhang 
1158b1683b59SHong Zhang 
1159b1683b59SHong Zhang /*S
1160b9af6bddSHong Zhang      MatTransposeColoring - Object for computing a sparse matrix product C=A*B^T via coloring
1161b1683b59SHong Zhang 
1162b1683b59SHong Zhang    Level: beginner
1163b1683b59SHong Zhang 
1164b1683b59SHong Zhang   Concepts: coloring, sparse matrix product
1165b1683b59SHong Zhang 
1166b9af6bddSHong Zhang .seealso:  MatTransposeColoringCreate()
1167b1683b59SHong Zhang S*/
1168b9af6bddSHong Zhang typedef struct _p_MatTransposeColoring* MatTransposeColoring;
1169b1683b59SHong Zhang 
1170014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeColoringCreate(Mat,ISColoring,MatTransposeColoring *);
1171014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransColoringApplySpToDen(MatTransposeColoring,Mat,Mat);
1172014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransColoringApplyDenToSp(MatTransposeColoring,Mat,Mat);
1173014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeColoringDestroy(MatTransposeColoring*);
1174b1683b59SHong Zhang 
1175639f9d9dSBarry Smith /*
11760752156aSBarry Smith     These routines are for partitioning matrices: currently used only
11773eda8832SBarry Smith   for adjacency matrix, MatCreateMPIAdj().
11780752156aSBarry Smith */
1179ca161407SBarry Smith 
1180d9274352SBarry Smith /*S
1181d9274352SBarry Smith      MatPartitioning - Object for managing the partitioning of a matrix or graph
1182d9274352SBarry Smith 
1183d9274352SBarry Smith    Level: beginner
1184d9274352SBarry Smith 
1185d9274352SBarry Smith   Concepts: partitioning
1186d9274352SBarry Smith 
1187743a1026Svictorle .seealso:  MatPartitioningCreate(), MatPartitioningType
1188d9274352SBarry Smith S*/
118991e9ee9fSBarry Smith typedef struct _p_MatPartitioning* MatPartitioning;
1190d9274352SBarry Smith 
119176bdecfbSBarry Smith /*J
11928f6c3df8SBarry Smith     MatPartitioningType - String with the name of a PETSc matrix partitioning
1193d9274352SBarry Smith 
1194d9274352SBarry Smith    Level: beginner
11950d04baf8SBarry Smith dm
1196b547a13bSvictorle .seealso: MatPartitioningCreate(), MatPartitioning
119776bdecfbSBarry Smith J*/
119819fd82e9SBarry Smith typedef const char* MatPartitioningType;
11992692d6eeSBarry Smith #define MATPARTITIONINGCURRENT  "current"
1200aa1e27eaSFande Kong #define MATPARTITIONINGAVERAGE   "average"
12012692d6eeSBarry Smith #define MATPARTITIONINGSQUARE   "square"
12022692d6eeSBarry Smith #define MATPARTITIONINGPARMETIS "parmetis"
12032692d6eeSBarry Smith #define MATPARTITIONINGCHACO    "chaco"
12042692d6eeSBarry Smith #define MATPARTITIONINGPARTY    "party"
120550d91057SBarry Smith #define MATPARTITIONINGPTSCOTCH "ptscotch"
120688d2ac2bSFande Kong #define MATPARTITIONINGHIERARCH  "hierarch"
120771306c60Sjroman 
1208ca161407SBarry Smith 
1209014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningCreate(MPI_Comm,MatPartitioning*);
121019fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningSetType(MatPartitioning,MatPartitioningType);
1211014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetNParts(MatPartitioning,PetscInt);
1212014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetAdjacency(MatPartitioning,Mat);
1213014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]);
1214014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []);
1215014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningApply(MatPartitioning,IS*);
1216014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningDestroy(MatPartitioning*);
12172aabb6bbSBarry Smith 
1218bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningRegister(const char[],PetscErrorCode (*)(MatPartitioning));
121930de9b25SBarry Smith 
1220f1144a30SSatish Balay 
12212bad1931SBarry Smith 
1222014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningView(MatPartitioning,PetscViewer);
1223014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetFromOptions(MatPartitioning);
122419fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningGetType(MatPartitioning,MatPartitioningType*);
1225ca161407SBarry Smith 
122627538973SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetRepartition(MatPartitioning part);
1227014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetCoarseSequential(MatPartitioning);
1228014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningParmetisGetEdgeCut(MatPartitioning, PetscInt *);
12290752156aSBarry Smith 
1230b6956c12SJose Roman typedef enum { MP_CHACO_MULTILEVEL=1,MP_CHACO_SPECTRAL=2,MP_CHACO_LINEAR=4,MP_CHACO_RANDOM=5,MP_CHACO_SCATTERED=6 } MPChacoGlobalType;
12316a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoGlobalTypes[];
1232b6956c12SJose Roman typedef enum { MP_CHACO_KERNIGHAN=1,MP_CHACO_NONE=2 } MPChacoLocalType;
12336a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoLocalTypes[];
1234b6956c12SJose Roman typedef enum { MP_CHACO_LANCZOS=0,MP_CHACO_RQI=1 } MPChacoEigenType;
12356a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoEigenTypes[];
1236b6956c12SJose Roman 
1237014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetGlobal(MatPartitioning,MPChacoGlobalType);
1238014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetGlobal(MatPartitioning,MPChacoGlobalType*);
1239014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetLocal(MatPartitioning,MPChacoLocalType);
1240014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetLocal(MatPartitioning,MPChacoLocalType*);
1241014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal);
1242014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType);
1243014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenSolver(MatPartitioning,MPChacoEigenType*);
1244014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenTol(MatPartitioning,PetscReal);
1245014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenTol(MatPartitioning,PetscReal*);
1246014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenNumber(MatPartitioning,PetscInt);
1247014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenNumber(MatPartitioning,PetscInt*);
124871306c60Sjroman 
124971306c60Sjroman #define MP_PARTY_OPT "opt"
125071306c60Sjroman #define MP_PARTY_LIN "lin"
125171306c60Sjroman #define MP_PARTY_SCA "sca"
125271306c60Sjroman #define MP_PARTY_RAN "ran"
125371306c60Sjroman #define MP_PARTY_GBF "gbf"
125471306c60Sjroman #define MP_PARTY_GCF "gcf"
125571306c60Sjroman #define MP_PARTY_BUB "bub"
125671306c60Sjroman #define MP_PARTY_DEF "def"
1257014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetGlobal(MatPartitioning,const char*);
125871306c60Sjroman #define MP_PARTY_HELPFUL_SETS "hs"
125971306c60Sjroman #define MP_PARTY_KERNIGHAN_LIN "kl"
126071306c60Sjroman #define MP_PARTY_NONE "no"
1261014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetLocal(MatPartitioning,const char*);
1262014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal);
1263014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetBipart(MatPartitioning,PetscBool);
1264014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscBool);
126571306c60Sjroman 
126650d91057SBarry Smith typedef enum { MP_PTSCOTCH_QUALITY,MP_PTSCOTCH_SPEED,MP_PTSCOTCH_BALANCE,MP_PTSCOTCH_SAFETY,MP_PTSCOTCH_SCALABILITY } MPPTScotchStrategyType;
12676a6fc655SJed Brown PETSC_EXTERN const char *const MPPTScotchStrategyTypes[];
1268e0f1cffaSJose Roman 
1269014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetImbalance(MatPartitioning,PetscReal);
1270014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetImbalance(MatPartitioning,PetscReal*);
1271014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetStrategy(MatPartitioning,MPPTScotchStrategyType);
1272014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetStrategy(MatPartitioning,MPPTScotchStrategyType*);
127371306c60Sjroman 
1274b43b03e9SMark F. Adams /*
12756294fa2bSFande Kong  * hierarchical partitioning
12766294fa2bSFande Kong  */
12774e713f55SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalGetFineparts(MatPartitioning,IS*);
12784e713f55SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalGetCoarseparts(MatPartitioning,IS*);
12794e713f55SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalSetNcoarseparts(MatPartitioning,PetscInt);
12804e713f55SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalSetNfineparts(MatPartitioning, PetscInt);
12816294fa2bSFande Kong 
12826294fa2bSFande Kong /*
1283b43b03e9SMark F. Adams     These routines are for coarsening matrices:
1284b43b03e9SMark F. Adams */
1285b43b03e9SMark F. Adams 
1286b43b03e9SMark F. Adams /*S
1287b43b03e9SMark F. Adams      MatCoarsen - Object for managing the coarsening of a graph (symmetric matrix)
1288b43b03e9SMark F. Adams 
1289b43b03e9SMark F. Adams    Level: beginner
1290b43b03e9SMark F. Adams 
1291b43b03e9SMark F. Adams   Concepts: coarsen
1292b43b03e9SMark F. Adams 
1293b43b03e9SMark F. Adams .seealso:  MatCoarsenCreate), MatCoarsenType
1294b43b03e9SMark F. Adams S*/
1295b43b03e9SMark F. Adams typedef struct _p_MatCoarsen* MatCoarsen;
1296b43b03e9SMark F. Adams 
1297b43b03e9SMark F. Adams /*J
12988f6c3df8SBarry Smith     MatCoarsenType - String with the name of a PETSc matrix coarsen
1299b43b03e9SMark F. Adams 
1300b43b03e9SMark F. Adams    Level: beginner
13018f6c3df8SBarry Smith 
1302b43b03e9SMark F. Adams .seealso: MatCoarsenCreate(), MatCoarsen
1303b43b03e9SMark F. Adams J*/
130419fd82e9SBarry Smith typedef const char* MatCoarsenType;
1305b43b03e9SMark F. Adams #define MATCOARSENMIS  "mis"
13069057884aSMark F. Adams #define MATCOARSENHEM  "hem"
1307b43b03e9SMark F. Adams 
13080cbbd2e1SMark F. Adams /* linked list for aggregates */
130941b27cdeSMark F. Adams typedef struct _PetscCDIntNd{
131041b27cdeSMark F. Adams   struct _PetscCDIntNd *next;
13110cbbd2e1SMark F. Adams   PetscInt             gid;
131241b27cdeSMark F. Adams }PetscCDIntNd;
13130cbbd2e1SMark F. Adams 
13140cbbd2e1SMark F. Adams /* only used by node pool */
131541b27cdeSMark F. Adams typedef struct _PetscCDArrNd{
131641b27cdeSMark F. Adams   struct _PetscCDArrNd *next;
131741b27cdeSMark F. Adams   struct _PetscCDIntNd *array;
131841b27cdeSMark F. Adams }PetscCDArrNd;
13190cbbd2e1SMark F. Adams 
13200cbbd2e1SMark F. Adams typedef struct _PetscCoarsenData{
132182c86c8fSBarry Smith   PetscCDArrNd pool_list;  /* node pool */
132241b27cdeSMark F. Adams   PetscCDIntNd *new_node;
13230cbbd2e1SMark F. Adams   PetscInt     new_left;
13240cbbd2e1SMark F. Adams   PetscInt     chk_sz;
132541b27cdeSMark F. Adams   PetscCDIntNd *extra_nodes;
132682c86c8fSBarry Smith   PetscCDIntNd **array;  /* Array of lists */
13270cbbd2e1SMark F. Adams   PetscInt     size;
132882c86c8fSBarry Smith   Mat          mat;  /* cache a Mat for communication data */
13290cbbd2e1SMark F. Adams }PetscCoarsenData;
13300cbbd2e1SMark F. Adams 
1331014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenCreate(MPI_Comm,MatCoarsen*);
133219fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenSetType(MatCoarsen,MatCoarsenType);
1333014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetAdjacency(MatCoarsen,Mat);
1334014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetGreedyOrdering(MatCoarsen,const IS);
1335014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetStrictAggs(MatCoarsen,PetscBool);
1336014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenGetData( MatCoarsen, PetscCoarsenData ** );
1337014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenApply(MatCoarsen);
1338014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenDestroy(MatCoarsen*);
1339b43b03e9SMark F. Adams 
1340bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenRegister(const char[],PetscErrorCode (*)(MatCoarsen));
1341b43b03e9SMark F. Adams 
1342b43b03e9SMark F. Adams 
1343b43b03e9SMark F. Adams 
1344014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenView(MatCoarsen,PetscViewer);
1345014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetFromOptions(MatCoarsen);
134619fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenGetType(MatCoarsen,MatCoarsenType*);
1347685405a1SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatCoarsenViewFromOptions(MatCoarsen A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);}
1348b43b03e9SMark F. Adams 
1349014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMeshToVertexGraph(Mat,PetscInt,Mat*);
1350014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMeshToCellGraph(Mat,PetscInt,Mat*);
1351591294e4SBarry Smith 
13520752156aSBarry Smith /*
1353af0996ceSBarry Smith     If you add entries here you must also add them to petsc/finclude/petscmat.h
1354d4fbbf0eSBarry Smith */
13551c1c02c0SLois Curfman McInnes typedef enum { MATOP_SET_VALUES=0,
13561c1c02c0SLois Curfman McInnes                MATOP_GET_ROW=1,
13571c1c02c0SLois Curfman McInnes                MATOP_RESTORE_ROW=2,
13581c1c02c0SLois Curfman McInnes                MATOP_MULT=3,
13591c1c02c0SLois Curfman McInnes                MATOP_MULT_ADD=4,
13607c922b88SBarry Smith                MATOP_MULT_TRANSPOSE=5,
13617c922b88SBarry Smith                MATOP_MULT_TRANSPOSE_ADD=6,
13621c1c02c0SLois Curfman McInnes                MATOP_SOLVE=7,
13631c1c02c0SLois Curfman McInnes                MATOP_SOLVE_ADD=8,
13647c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE=9,
13657c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE_ADD=10,
13661c1c02c0SLois Curfman McInnes                MATOP_LUFACTOR=11,
13671c1c02c0SLois Curfman McInnes                MATOP_CHOLESKYFACTOR=12,
1368a714c06dSBarry Smith                MATOP_SOR=13,
13691c1c02c0SLois Curfman McInnes                MATOP_TRANSPOSE=14,
13701c1c02c0SLois Curfman McInnes                MATOP_GETINFO=15,
13711c1c02c0SLois Curfman McInnes                MATOP_EQUAL=16,
13721c1c02c0SLois Curfman McInnes                MATOP_GET_DIAGONAL=17,
13731c1c02c0SLois Curfman McInnes                MATOP_DIAGONAL_SCALE=18,
13741c1c02c0SLois Curfman McInnes                MATOP_NORM=19,
13751c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_BEGIN=20,
13761c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_END=21,
1377d519adbfSMatthew Knepley                MATOP_SET_OPTION=22,
1378d519adbfSMatthew Knepley                MATOP_ZERO_ENTRIES=23,
1379d519adbfSMatthew Knepley                MATOP_ZERO_ROWS=24,
1380d519adbfSMatthew Knepley                MATOP_LUFACTOR_SYMBOLIC=25,
1381d519adbfSMatthew Knepley                MATOP_LUFACTOR_NUMERIC=26,
1382d519adbfSMatthew Knepley                MATOP_CHOLESKY_FACTOR_SYMBOLIC=27,
1383d519adbfSMatthew Knepley                MATOP_CHOLESKY_FACTOR_NUMERIC=28,
1384d519adbfSMatthew Knepley                MATOP_SETUP_PREALLOCATION=29,
1385d519adbfSMatthew Knepley                MATOP_ILUFACTOR_SYMBOLIC=30,
1386d519adbfSMatthew Knepley                MATOP_ICCFACTOR_SYMBOLIC=31,
13879d39f69dSJed Brown                /* MATOP_PLACEHOLDER_32=32, */
13889d39f69dSJed Brown                /* MATOP_PLACEHOLDER_33=33, */
1389d519adbfSMatthew Knepley                MATOP_DUPLICATE=34,
1390d519adbfSMatthew Knepley                MATOP_FORWARD_SOLVE=35,
1391d519adbfSMatthew Knepley                MATOP_BACKWARD_SOLVE=36,
1392d519adbfSMatthew Knepley                MATOP_ILUFACTOR=37,
1393d519adbfSMatthew Knepley                MATOP_ICCFACTOR=38,
1394d519adbfSMatthew Knepley                MATOP_AXPY=39,
1395d519adbfSMatthew Knepley                MATOP_GET_SUBMATRICES=40,
1396d519adbfSMatthew Knepley                MATOP_INCREASE_OVERLAP=41,
1397d519adbfSMatthew Knepley                MATOP_GET_VALUES=42,
1398d519adbfSMatthew Knepley                MATOP_COPY=43,
1399d519adbfSMatthew Knepley                MATOP_GET_ROW_MAX=44,
1400d519adbfSMatthew Knepley                MATOP_SCALE=45,
1401d519adbfSMatthew Knepley                MATOP_SHIFT=46,
140235153367SBarry Smith                MATOP_DIAGONAL_SET=47,
14039d39f69dSJed Brown                MATOP_ZERO_ROWS_COLUMNS=48,
14049d39f69dSJed Brown                MATOP_SET_RANDOM=49,
1405d519adbfSMatthew Knepley                MATOP_GET_ROW_IJ=50,
1406d519adbfSMatthew Knepley                MATOP_RESTORE_ROW_IJ=51,
1407d519adbfSMatthew Knepley                MATOP_GET_COLUMN_IJ=52,
1408d519adbfSMatthew Knepley                MATOP_RESTORE_COLUMN_IJ=53,
1409d519adbfSMatthew Knepley                MATOP_FDCOLORING_CREATE=54,
1410d519adbfSMatthew Knepley                MATOP_COLORING_PATCH=55,
1411d519adbfSMatthew Knepley                MATOP_SET_UNFACTORED=56,
1412d519adbfSMatthew Knepley                MATOP_PERMUTE=57,
1413d519adbfSMatthew Knepley                MATOP_SET_VALUES_BLOCKED=58,
1414d519adbfSMatthew Knepley                MATOP_GET_SUBMATRIX=59,
1415d519adbfSMatthew Knepley                MATOP_DESTROY=60,
1416d519adbfSMatthew Knepley                MATOP_VIEW=61,
1417d519adbfSMatthew Knepley                MATOP_CONVERT_FROM=62,
14187bab7c10SHong Zhang                MATOP_MATMAT_MULT=63,
14197bab7c10SHong Zhang                MATOP_MATMAT_MULT_SYMBOLIC=64,
14207bab7c10SHong Zhang                MATOP_MATMAT_MULT_NUMERIC=65,
1421d519adbfSMatthew Knepley                MATOP_SET_LOCAL_TO_GLOBAL_MAP=66,
1422d519adbfSMatthew Knepley                MATOP_SET_VALUES_LOCAL=67,
1423d519adbfSMatthew Knepley                MATOP_ZERO_ROWS_LOCAL=68,
1424d519adbfSMatthew Knepley                MATOP_GET_ROW_MAX_ABS=69,
1425d519adbfSMatthew Knepley                MATOP_GET_ROW_MIN_ABS=70,
1426d519adbfSMatthew Knepley                MATOP_CONVERT=71,
1427d519adbfSMatthew Knepley                MATOP_SET_COLORING=72,
14289d39f69dSJed Brown                /* MATOP_PLACEHOLDER_73=73, */
1429d519adbfSMatthew Knepley                MATOP_SET_VALUES_ADIFOR=74,
1430d519adbfSMatthew Knepley                MATOP_FD_COLORING_APPLY=75,
1431d519adbfSMatthew Knepley                MATOP_SET_FROM_OPTIONS=76,
1432bda6c4cbSBarry Smith                MATOP_MULT_CONSTRAINED=77,
1433bda6c4cbSBarry Smith                MATOP_MULT_TRANSPOSE_CONSTRAIN=78,
14349d39f69dSJed Brown                MATOP_FIND_ZERO_DIAGONALS=79,
1435d519adbfSMatthew Knepley                MATOP_MULT_MULTIPLE=80,
1436d519adbfSMatthew Knepley                MATOP_SOLVE_MULTIPLE=81,
1437d519adbfSMatthew Knepley                MATOP_GET_INERTIA=82,
1438d519adbfSMatthew Knepley                MATOP_LOAD=83,
1439d519adbfSMatthew Knepley                MATOP_IS_SYMMETRIC=84,
1440d519adbfSMatthew Knepley                MATOP_IS_HERMITIAN=85,
1441d519adbfSMatthew Knepley                MATOP_IS_STRUCTURALLY_SYMMETRIC=86,
1442820469bcSHong Zhang                MATOP_SET_VALUES_BLOCKEDLOCAL=87,
1443d519adbfSMatthew Knepley                MATOP_GET_VECS=88,
1444d519adbfSMatthew Knepley                MATOP_MAT_MULT=89,
1445d519adbfSMatthew Knepley                MATOP_MAT_MULT_SYMBOLIC=90,
1446d519adbfSMatthew Knepley                MATOP_MAT_MULT_NUMERIC=91,
1447d519adbfSMatthew Knepley                MATOP_PTAP=92,
1448d519adbfSMatthew Knepley                MATOP_PTAP_SYMBOLIC=93,
1449d519adbfSMatthew Knepley                MATOP_PTAP_NUMERIC=94,
1450bda6c4cbSBarry Smith                MATOP_MAT_TRANSPOSE_MULT=95,
1451bda6c4cbSBarry Smith                MATOP_MAT_TRANSPOSE_MULT_SYMBO=96,
1452bda6c4cbSBarry Smith                MATOP_MAT_TRANSPOSE_MULT_NUMER=97,
14539d39f69dSJed Brown                /* MATOP_PLACEHOLDER_98=98, */
14549d39f69dSJed Brown                /* MATOP_PLACEHOLDER_99=99, */
14559d39f69dSJed Brown                /* MATOP_PLACEHOLDER_100=100, */
14569d39f69dSJed Brown                /* MATOP_PLACEHOLDER_101=101, */
1457d519adbfSMatthew Knepley                MATOP_CONJUGATE=102,
14589d39f69dSJed Brown                /* MATOP_PLACEHOLDER_103=103, */
1459d519adbfSMatthew Knepley                MATOP_SET_VALUES_ROW=104,
1460d519adbfSMatthew Knepley                MATOP_REAL_PART=105,
1461bda6c4cbSBarry Smith                MATOP_IMAGINARY_PART=106,
1462bda6c4cbSBarry Smith                MATOP_GET_ROW_UPPER_TRIANGULAR=107,
1463bda6c4cbSBarry Smith                MATOP_RESTORE_ROW_UPPER_TRIANG=108,
1464bda6c4cbSBarry Smith                MATOP_MAT_SOLVE=109,
14650e8d04f7SBarry Smith                MATOP_GET_REDUNDANT_MATRIX=110,
1466d519adbfSMatthew Knepley                MATOP_GET_ROW_MIN=111,
14670e8d04f7SBarry Smith                MATOP_GET_COLUMN_VECTOR=112,
1468d519adbfSMatthew Knepley                MATOP_MISSING_DIAGONAL=113,
14690e8d04f7SBarry Smith                MATOP_GET_SEQ_NONZERO_STRUCTUR=114,
147089c6957cSBarry Smith                MATOP_CREATE=115,
147189c6957cSBarry Smith                MATOP_GET_GHOSTS=116,
14720e8d04f7SBarry Smith                MATOP_GET_LOCAL_SUB_MATRIX=117,
14730e8d04f7SBarry Smith                MATOP_RESTORE_LOCALSUB_MATRIX=118,
1474eeffb40dSHong Zhang                MATOP_MULT_DIAGONAL_BLOCK=119,
14750e8d04f7SBarry Smith                MATOP_HERMITIAN_TRANSPOSE=120,
14760e8d04f7SBarry Smith                MATOP_MULT_HERMITIAN_TRANSPOSE=121,
14770e8d04f7SBarry Smith                MATOP_MULT_HERMITIAN_TRANS_ADD=122,
14780e8d04f7SBarry Smith                MATOP_GET_MULTI_PROC_BLOCK=123,
14799d39f69dSJed Brown                MATOP_FIND_NONZERO_ROWS=124,
14800e8d04f7SBarry Smith                MATOP_GET_COLUMN_NORMS=125,
14819d39f69dSJed Brown                MATOP_INVERT_BLOCK_DIAGONAL=126,
14829d39f69dSJed Brown                /* MATOP_PLACEHOLDER_127=127, */
14830e8d04f7SBarry Smith                MATOP_GET_SUB_MATRICES_PARALLE=128,
14842b8ad9a3SHong Zhang                MATOP_SET_VALUES_BATCH=129,
14850e8d04f7SBarry Smith                MATOP_TRANSPOSE_MAT_MULT=130,
14860e8d04f7SBarry Smith                MATOP_TRANSPOSE_MAT_MULT_SYMBO=131,
1487e9b602ebSSatish Balay                MATOP_TRANSPOSE_MAT_MULT_NUMER=132,
14880e8d04f7SBarry Smith                MATOP_TRANSPOSE_COLORING_CREAT=133,
14890e8d04f7SBarry Smith                MATOP_TRANS_COLORING_APPLY_SPT=134,
14900e8d04f7SBarry Smith                MATOP_TRANS_COLORING_APPLY_DEN=135,
14910e8d04f7SBarry Smith                MATOP_RART=136,
14920e8d04f7SBarry Smith                MATOP_RART_SYMBOLIC=137,
14930e8d04f7SBarry Smith                MATOP_RART_NUMERIC=138,
1494e09a3074SHong Zhang                MATOP_SET_BLOCK_SIZES=139,
1495f9426fe0SMark Adams                MATOP_AYPX=140,
14961919a2e2SJed Brown                MATOP_RESIDUAL=141,
14979c8f2541SHong Zhang                MATOP_FDCOLORING_SETUP=142,
14989c8f2541SHong Zhang                MATOP_MPICONCATENATESEQ=144
1499fae171e0SBarry Smith              } MatOperation;
1500014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatHasOperation(Mat,MatOperation,PetscBool *);
1501014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellSetOperation(Mat,MatOperation,void(*)(void));
1502014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellGetOperation(Mat,MatOperation,void(**)(void));
1503014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellSetContext(Mat,void*);
1504112a2221SBarry Smith 
150590ace30eSBarry Smith /*
150690ace30eSBarry Smith    Codes for matrices stored on disk. By default they are
150790ace30eSBarry Smith    stored in a universal format. By changing the format with
15087973ad04SBarry Smith    PetscViewerSetFormat(viewer,PETSC_VIEWER_NATIVE); the matrices will
150990ace30eSBarry Smith    be stored in a way natural for the matrix, for example dense matrices
151090ace30eSBarry Smith    would be stored as dense. Matrices stored this way may only be
1511f4403165SShri Abhyankar    read into matrices of the same type.
151290ace30eSBarry Smith */
151390ace30eSBarry Smith #define MATRIX_BINARY_FORMAT_DENSE -1
151490ace30eSBarry Smith 
1515014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetHashTableFactor(Mat,PetscReal);
1516014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatISGetLocalMat(Mat,Mat*);
1517014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatISSetLocalMat(Mat,Mat);
1518b7ce53b6SStefano Zampini PETSC_EXTERN PetscErrorCode MatISGetMPIXAIJ(Mat,MatReuse,Mat*);
15191f4e1ec7SBarry Smith 
1520d9274352SBarry Smith /*S
1521d9274352SBarry Smith      MatNullSpace - Object that removes a null space from a vector, i.e.
1522d9274352SBarry Smith          orthogonalizes the vector to a subsapce
1523d9274352SBarry Smith 
1524f7a9e4ceSBarry Smith    Level: advanced
1525d9274352SBarry Smith 
1526d9274352SBarry Smith   Concepts: matrix; linear operator, null space
1527d9274352SBarry Smith 
15286e1639daSBarry Smith   Users manual sections:
15296e1639daSBarry Smith .   sec_singular
15306e1639daSBarry Smith 
1531d9274352SBarry Smith .seealso:  MatNullSpaceCreate()
1532d9274352SBarry Smith S*/
153374637425SBarry Smith typedef struct _p_MatNullSpace* MatNullSpace;
1534d9274352SBarry Smith 
1535014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceCreate(MPI_Comm,PetscBool ,PetscInt,const Vec[],MatNullSpace*);
1536014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(MatNullSpace,Vec,void*),void*);
1537014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceDestroy(MatNullSpace*);
1538d0195637SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceRemove(MatNullSpace,Vec);
1539014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNullSpace(Mat, MatNullSpace *);
15405fa7ec2dSBarry Smith PETSC_EXTERN PetscErrorCode MatGetTransposeNullSpace(Mat, MatNullSpace *);
15415fa7ec2dSBarry Smith PETSC_EXTERN PetscErrorCode MatSetTransposeNullSpace(Mat,MatNullSpace);
1542014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNullSpace(Mat,MatNullSpace);
1543014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNearNullSpace(Mat,MatNullSpace);
1544014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNearNullSpace(Mat,MatNullSpace*);
1545014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceTest(MatNullSpace,Mat,PetscBool  *);
1546014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceView(MatNullSpace,PetscViewer);
1547014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceGetVecs(MatNullSpace,PetscBool*,PetscInt*,const Vec**);
1548014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceCreateRigidBody(Vec,MatNullSpace*);
154974637425SBarry Smith 
1550014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatReorderingSeqSBAIJ(Mat,IS);
1551014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetHashTableFactor(Mat,PetscReal);
1552014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetColumnIndices(Mat,PetscInt *);
1553014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJInvertBlockDiagonal(Mat);
15543f1d51d7SBarry Smith 
1555014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMAIJ(Mat,PetscInt,Mat*);
1556014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMAIJRedimension(Mat,PetscInt,Mat*);
1557014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMAIJGetAIJ(Mat,Mat*);
1558c4f061fbSSatish Balay 
1559014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatComputeExplicitOperator(Mat,Mat*);
1560b0a32e0cSBarry Smith 
1561014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalScaleLocal(Mat,Vec);
156204f1ad80SBarry Smith 
1563014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*);
1564014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetBase(Mat,Vec,Vec);
1565014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*);
1566014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*));
1567014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec));
1568014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt);
1569014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDResetHHistory(Mat);
1570014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctionError(Mat,PetscReal);
1571014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetPeriod(Mat,PetscInt);
1572014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDGetH(Mat,PetscScalar *);
1573014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetOptionsPrefix(Mat,const char[]);
1574014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*);
1575014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*);
1576e884886fSBarry Smith 
15776370053bSBarry Smith /*S
15786370053bSBarry Smith     MatMFFD - A data structured used to manage the computation of the h differencing parameter for matrix-free
15796370053bSBarry Smith               Jacobian vector products
1580e884886fSBarry Smith 
15816370053bSBarry Smith     Notes: MATMFFD is a specific MatType which uses the MatMFFD data structure
15826370053bSBarry Smith 
15836370053bSBarry Smith            MatMFFD*() methods actually take the Mat as their first argument. Not a MatMFFD data structure
15846370053bSBarry Smith 
15856370053bSBarry Smith     Level: developer
15866370053bSBarry Smith 
15876370053bSBarry Smith .seealso: MATMFFD, MatCreateMFFD(), MatMFFDSetFuction(), MatMFFDSetType(), MatMFFDRegister()
15886370053bSBarry Smith S*/
1589e884886fSBarry Smith typedef struct _p_MatMFFD* MatMFFD;
1590e884886fSBarry Smith 
159176bdecfbSBarry Smith /*J
1592e884886fSBarry Smith     MatMFFDType - algorithm used to compute the h used in computing matrix-vector products via differencing of the function
1593e884886fSBarry Smith 
1594e884886fSBarry Smith    Level: beginner
1595e884886fSBarry Smith 
1596e884886fSBarry Smith .seealso: MatMFFDSetType(), MatMFFDRegister()
159776bdecfbSBarry Smith J*/
159819fd82e9SBarry Smith typedef const char* MatMFFDType;
1599e884886fSBarry Smith #define MATMFFD_DS  "ds"
1600e884886fSBarry Smith #define MATMFFD_WP  "wp"
1601e884886fSBarry Smith 
160219fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDSetType(Mat,MatMFFDType);
1603bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDRegister(const char[],PetscErrorCode (*)(MatMFFD));
1604e884886fSBarry Smith 
1605014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDDSSetUmin(Mat,PetscReal);
1606014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDWPSetComputeNormU(Mat,PetscBool );
1607e884886fSBarry Smith 
1608014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *);
1609014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *);
16107dbadf16SMatthew Knepley 
161197969023SHong Zhang /*
161297969023SHong Zhang    PETSc interface to MUMPS
161397969023SHong Zhang */
161497969023SHong Zhang #ifdef PETSC_HAVE_MUMPS
1615014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMumpsSetIcntl(Mat,PetscInt,PetscInt);
1616bc6112feSHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetIcntl(Mat,PetscInt,PetscInt*);
1617b9e9ea26SSatish Balay PETSC_EXTERN PetscErrorCode MatMumpsSetCntl(Mat,PetscInt,PetscReal);
1618bc6112feSHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetCntl(Mat,PetscInt,PetscReal*);
1619bc6112feSHong Zhang 
1620ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetInfo(Mat,PetscInt,PetscInt*);
1621ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetInfog(Mat,PetscInt,PetscInt*);
1622ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetRinfo(Mat,PetscInt,PetscReal*);
1623ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetRinfog(Mat,PetscInt,PetscReal*);
16246444a565SStefano Zampini 
162559ac8732SStefano Zampini PETSC_EXTERN PetscErrorCode MatMumpsInvertSchurComplement(Mat);
162659ac8732SStefano Zampini PETSC_EXTERN PetscErrorCode MatMumpsCreateSchurComplement(Mat,Mat*);
16276444a565SStefano Zampini PETSC_EXTERN PetscErrorCode MatMumpsGetSchurComplement(Mat,Mat*);
162859ac8732SStefano Zampini PETSC_EXTERN PetscErrorCode MatMumpsRestoreSchurComplement(Mat,Mat*);
1629e807eca7SStefano Zampini PETSC_EXTERN PetscErrorCode MatMumpsSolveSchurComplement(Mat,Vec,Vec);
16307404bcfbSStefano Zampini PETSC_EXTERN PetscErrorCode MatMumpsSolveSchurComplementTranspose(Mat,Vec,Vec);
163197969023SHong Zhang #endif
163223a5497aSJed Brown 
1633d954db57SHong Zhang /*
1634d615f992SSatish Balay    PETSc interface to Mkl_Pardiso
1635b500a6b7SJose David Bermeol */
1636b500a6b7SJose David Bermeol #ifdef PETSC_HAVE_MKL_PARDISO
1637d615f992SSatish Balay PETSC_EXTERN PetscErrorCode MatMkl_PardisoSetCntl(Mat,PetscInt,PetscInt);
1638b500a6b7SJose David Bermeol #endif
1639b500a6b7SJose David Bermeol 
1640b500a6b7SJose David Bermeol /*
1641d305a81bSVasiliy Kozyrev    PETSc interface to Mkl_CPardiso
1642d305a81bSVasiliy Kozyrev */
1643d305a81bSVasiliy Kozyrev #ifdef PETSC_HAVE_MKL_CPARDISO
1644d305a81bSVasiliy Kozyrev PETSC_EXTERN PetscErrorCode MatMkl_CPardisoSetCntl(Mat,PetscInt,PetscInt);
1645d305a81bSVasiliy Kozyrev #endif
1646d305a81bSVasiliy Kozyrev 
1647d305a81bSVasiliy Kozyrev /*
1648d954db57SHong Zhang    PETSc interface to SUPERLU
1649d954db57SHong Zhang */
1650d954db57SHong Zhang #ifdef PETSC_HAVE_SUPERLU
1651014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSuperluSetILUDropTol(Mat,PetscReal);
1652d954db57SHong Zhang #endif
1653d954db57SHong Zhang 
1654aa372e3fSPaul Mullowney #ifdef PETSC_HAVE_CUDA
16553f9c0db1SPaul Mullowney /*E
1656e057df02SPaul Mullowney     MatCUSPARSEStorageFormat - indicates the storage format for CUSPARSE (GPU)
16572692e278SPaul Mullowney     matrices.
1658e057df02SPaul Mullowney 
1659e057df02SPaul Mullowney     Not Collective
1660e057df02SPaul Mullowney 
16618468deeeSKarl Rupp +   MAT_CUSPARSE_CSR - Compressed Sparse Row
16622692e278SPaul Mullowney .   MAT_CUSPARSE_ELL - Ellpack (requires CUDA 4.2 or later).
16632692e278SPaul Mullowney -   MAT_CUSPARSE_HYB - Hybrid, a combination of Ellpack and Coordinate format (requires CUDA 4.2 or later).
1664e057df02SPaul Mullowney 
1665e057df02SPaul Mullowney     Level: intermediate
1666e057df02SPaul Mullowney 
1667af0996ceSBarry Smith    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
1668e057df02SPaul Mullowney 
1669e057df02SPaul Mullowney .seealso: MatCUSPARSESetFormat(), MatCUSPARSEFormatOperation
1670e057df02SPaul Mullowney E*/
1671e057df02SPaul Mullowney 
16723f9c0db1SPaul Mullowney typedef enum {MAT_CUSPARSE_CSR, MAT_CUSPARSE_ELL, MAT_CUSPARSE_HYB} MatCUSPARSEStorageFormat;
1673e057df02SPaul Mullowney 
1674e057df02SPaul Mullowney /* these will be strings associated with enumerated type defined above */
1675e057df02SPaul Mullowney PETSC_EXTERN const char *const MatCUSPARSEStorageFormats[];
1676e057df02SPaul Mullowney 
16773f9c0db1SPaul Mullowney /*E
1678e057df02SPaul Mullowney     MatCUSPARSEFormatOperation - indicates the operation of CUSPARSE (GPU)
16792692e278SPaul Mullowney     matrices whose operation should use a particular storage format.
1680e057df02SPaul Mullowney 
1681e057df02SPaul Mullowney     Not Collective
1682e057df02SPaul Mullowney 
16838468deeeSKarl Rupp +   MAT_CUSPARSE_MULT_DIAG - sets the storage format for the diagonal matrix in the parallel MatMult
16848468deeeSKarl Rupp .   MAT_CUSPARSE_MULT_OFFDIAG - sets the storage format for the offdiagonal matrix in the parallel MatMult
16858468deeeSKarl Rupp .   MAT_CUSPARSE_MULT - sets the storage format for the entire matrix in the serial (single GPU) MatMult
16868468deeeSKarl Rupp -   MAT_CUSPARSE_ALL - sets the storage format for all CUSPARSE (GPU) matrices
1687e057df02SPaul Mullowney 
1688e057df02SPaul Mullowney     Level: intermediate
1689e057df02SPaul Mullowney 
1690e057df02SPaul Mullowney .seealso: MatCUSPARSESetFormat(), MatCUSPARSEStorageFormat
1691e057df02SPaul Mullowney E*/
169236d62e41SPaul Mullowney typedef enum {MAT_CUSPARSE_MULT_DIAG, MAT_CUSPARSE_MULT_OFFDIAG, MAT_CUSPARSE_MULT, MAT_CUSPARSE_ALL} MatCUSPARSEFormatOperation;
1693e057df02SPaul Mullowney 
169410e9db80SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
169510e9db80SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCreateAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
1696e057df02SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCUSPARSESetFormat(Mat,MatCUSPARSEFormatOperation,MatCUSPARSEStorageFormat);
16979ae82921SPaul Mullowney #endif
16989ae82921SPaul Mullowney 
169990c53902SBarry Smith #if defined(PETSC_HAVE_CUSP)
1700014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCUSP(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
1701014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateAIJCUSP(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
1702e057df02SPaul Mullowney 
17033f9c0db1SPaul Mullowney /*E
1704e057df02SPaul Mullowney     MatCUSPStorageFormat - indicates the storage format for CUSP (GPU)
170536d62e41SPaul Mullowney     matrices.
1706e057df02SPaul Mullowney 
1707e057df02SPaul Mullowney     Not Collective
1708e057df02SPaul Mullowney 
17098468deeeSKarl Rupp +   MAT_CUSP_CSR - Compressed Sparse Row
17108468deeeSKarl Rupp .   MAT_CUSP_DIA - Diagonal
17118468deeeSKarl Rupp -   MAT_CUSP_ELL - Ellpack
1712e057df02SPaul Mullowney 
1713e057df02SPaul Mullowney     Level: intermediate
1714e057df02SPaul Mullowney 
1715af0996ceSBarry Smith    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
1716e057df02SPaul Mullowney 
1717e057df02SPaul Mullowney .seealso: MatCUSPSetFormat(), MatCUSPFormatOperation
1718e057df02SPaul Mullowney E*/
17193f9c0db1SPaul Mullowney typedef enum {MAT_CUSP_CSR, MAT_CUSP_DIA, MAT_CUSP_ELL} MatCUSPStorageFormat;
1720e057df02SPaul Mullowney 
1721e057df02SPaul Mullowney /* these will be strings associated with enumerated type defined above */
1722e057df02SPaul Mullowney PETSC_EXTERN const char *const MatCUSPStorageFormats[];
1723e057df02SPaul Mullowney 
17243f9c0db1SPaul Mullowney /*E
1725e057df02SPaul Mullowney     MatCUSPFormatOperation - indicates the operation of CUSP (GPU)
172636d62e41SPaul Mullowney     matrices whose operation should use a particular storage format.
1727e057df02SPaul Mullowney 
1728e057df02SPaul Mullowney     Not Collective
1729e057df02SPaul Mullowney 
17308468deeeSKarl Rupp +   MAT_CUSP_MULT_DIAG - sets the storage format for the diagonal matrix in the parallel MatMult
17318468deeeSKarl Rupp .   MAT_CUSP_MULT_OFFDIAG - sets the storage format for the offdiagonal matrix in the parallel MatMult
17328468deeeSKarl Rupp .   MAT_CUSP_MULT - sets the storage format for the entire matrix in the serial (single GPU) MatMult
17338468deeeSKarl Rupp -   MAT_CUSP_ALL - sets the storage format for all CUSP (GPU) matrices
1734e057df02SPaul Mullowney 
1735e057df02SPaul Mullowney     Level: intermediate
1736e057df02SPaul Mullowney 
1737af0996ceSBarry Smith    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
1738e057df02SPaul Mullowney 
1739e057df02SPaul Mullowney .seealso: MatCUSPSetFormat(), MatCUSPStorageFormat
1740e057df02SPaul Mullowney E*/
174136d62e41SPaul Mullowney typedef enum {MAT_CUSP_MULT_DIAG, MAT_CUSP_MULT_OFFDIAG, MAT_CUSP_MULT, MAT_CUSP_ALL} MatCUSPFormatOperation;
1742e057df02SPaul Mullowney 
1743e057df02SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCUSPSetFormat(Mat,MatCUSPFormatOperation,MatCUSPStorageFormat);
1744e057df02SPaul Mullowney #endif
174590c53902SBarry Smith 
1746d67ff14aSKarl Rupp #if defined(PETSC_HAVE_VIENNACL)
1747d67ff14aSKarl Rupp PETSC_EXTERN PetscErrorCode MatCreateSeqAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
1748d67ff14aSKarl Rupp PETSC_EXTERN PetscErrorCode MatCreateAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
1749d67ff14aSKarl Rupp #endif
1750d67ff14aSKarl Rupp 
175154efbe56SHong Zhang /*
175254efbe56SHong Zhang    PETSc interface to FFTW
175354efbe56SHong Zhang */
175454efbe56SHong Zhang #if defined(PETSC_HAVE_FFTW)
1755014dd563SJed Brown PETSC_EXTERN PetscErrorCode VecScatterPetscToFFTW(Mat,Vec,Vec);
1756014dd563SJed Brown PETSC_EXTERN PetscErrorCode VecScatterFFTWToPetsc(Mat,Vec,Vec);
17572a7a6963SBarry Smith PETSC_EXTERN PetscErrorCode MatCreateVecsFFTW(Mat,Vec*,Vec*,Vec*);
175854efbe56SHong Zhang #endif
175954efbe56SHong Zhang 
1760382fd914SHong Zhang /*
1761382fd914SHong Zhang    PETSc interface to ELEMENTAL
1762382fd914SHong Zhang */
1763db31f6deSJed Brown #if defined(PETSC_HAVE_ELEMENTAL)
17642ef15aa2SHong Zhang #if defined(PETSC_USE_COMPLEX)
17652ef15aa2SHong Zhang typedef El::Complex<PetscReal> PetscElemScalar;
17662ef15aa2SHong Zhang #else
17672ef15aa2SHong Zhang typedef PetscScalar PetscElemScalar;
1768382fd914SHong Zhang #endif
1769607a6623SBarry Smith PETSC_EXTERN PetscErrorCode PetscElementalInitializePackage(void);
1770db31f6deSJed Brown PETSC_EXTERN PetscErrorCode PetscElementalFinalizePackage(void);
1771db31f6deSJed Brown #endif
1772db31f6deSJed Brown 
1773014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateNest(MPI_Comm,PetscInt,const IS[],PetscInt,const IS[],const Mat[],Mat*);
1774014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSize(Mat,PetscInt*,PetscInt*);
1775014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetISs(Mat,IS[],IS[]);
1776014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetLocalISs(Mat,IS[],IS[]);
1777014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSubMats(Mat,PetscInt*,PetscInt*,Mat***);
1778014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSubMat(Mat,PetscInt,PetscInt,Mat*);
177919fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatNestSetVecType(Mat,VecType);
1780014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestSetSubMats(Mat,PetscInt,const IS[],PetscInt,const IS[],const Mat[]);
1781014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestSetSubMat(Mat,PetscInt,PetscInt,Mat);
1782d8588912SDave May 
17834325cce7SMatthew G Knepley PETSC_EXTERN PetscErrorCode MatChop(Mat,PetscReal);
1784e0a58c10SMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatComputeBandwidth(Mat,PetscReal,PetscInt*);
17854325cce7SMatthew G Knepley 
1786b541e6a4SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSubdomainsCreateCoalesce(Mat,PetscInt,PetscInt*,IS**);
178753dd7562SDmitry Karpeev 
178823a5497aSJed Brown #endif
1789