xref: /petsc/include/petscmat.h (revision df750dc879e6b61f257e3b41356584ef8bf91d84)
12eac72dbSBarry Smith /*
22eac72dbSBarry Smith      Include file for the matrix component of PETSc
32eac72dbSBarry Smith */
40a835dfdSSatish Balay #ifndef __PETSCMAT_H
50a835dfdSSatish Balay #define __PETSCMAT_H
62c8e378dSBarry Smith #include <petscvec.h>
72eac72dbSBarry Smith 
8d9274352SBarry Smith /*S
98f6c3df8SBarry Smith      Mat - Abstract PETSc matrix object used to manage all linear operators in PETSc, even those without
108f6c3df8SBarry Smith            an explicit sparse representation (such as matrix-free operators)
112eac72dbSBarry Smith 
12d91e6319SBarry Smith    Level: beginner
13d91e6319SBarry Smith 
14d9274352SBarry Smith   Concepts: matrix; linear operator
15d9274352SBarry Smith 
168f6c3df8SBarry Smith .seealso:  MatCreate(), MatType, MatSetType(), MatDestroy()
17d9274352SBarry Smith S*/
18d9274352SBarry Smith typedef struct _p_Mat*           Mat;
19d9274352SBarry Smith 
2076bdecfbSBarry Smith /*J
218f6c3df8SBarry Smith     MatType - String with the name of a PETSc matrix type
22d9274352SBarry Smith 
23d9274352SBarry Smith    Level: beginner
24d9274352SBarry Smith 
258f6c3df8SBarry Smith .seealso: MatSetType(), Mat, MatSolverPackage, MatRegister()
2676bdecfbSBarry Smith J*/
2719fd82e9SBarry Smith typedef const char* MatType;
28273d9f13SBarry Smith #define MATSAME            "same"
295a11e1b2SBarry Smith #define MATMAIJ            "maij"
30273d9f13SBarry Smith #define MATSEQMAIJ         "seqmaij"
31273d9f13SBarry Smith #define MATMPIMAIJ         "mpimaij"
32273d9f13SBarry Smith #define MATIS              "is"
335a11e1b2SBarry Smith #define MATAIJ             "aij"
34273d9f13SBarry Smith #define MATSEQAIJ          "seqaij"
35273d9f13SBarry Smith #define MATMPIAIJ          "mpiaij"
365a11e1b2SBarry Smith #define MATAIJCRL          "aijcrl"
375a11e1b2SBarry Smith #define MATSEQAIJCRL       "seqaijcrl"
385a11e1b2SBarry Smith #define MATMPIAIJCRL       "mpiaijcrl"
398154be41SBarry Smith #define MATAIJCUSP         "aijcusp"
408154be41SBarry Smith #define MATSEQAIJCUSP      "seqaijcusp"
418154be41SBarry Smith #define MATMPIAIJCUSP      "mpiaijcusp"
429ae82921SPaul Mullowney #define MATAIJCUSPARSE     "aijcusparse"
439ae82921SPaul Mullowney #define MATSEQAIJCUSPARSE  "seqaijcusparse"
449ae82921SPaul Mullowney #define MATMPIAIJCUSPARSE  "mpiaijcusparse"
45d67ff14aSKarl Rupp #define MATAIJVIENNACL     "aijviennacl"
46d67ff14aSKarl Rupp #define MATSEQAIJVIENNACL  "seqaijviennacl"
47d67ff14aSKarl Rupp #define MATMPIAIJVIENNACL  "mpiaijviennacl"
485a11e1b2SBarry Smith #define MATAIJPERM         "aijperm"
495a11e1b2SBarry Smith #define MATSEQAIJPERM      "seqaijperm"
505a11e1b2SBarry Smith #define MATMPIAIJPERM      "mpiaijperm"
51273d9f13SBarry Smith #define MATSHELL           "shell"
525a11e1b2SBarry Smith #define MATDENSE           "dense"
53209238afSKris Buschelman #define MATSEQDENSE        "seqdense"
54273d9f13SBarry Smith #define MATMPIDENSE        "mpidense"
55db31f6deSJed Brown #define MATELEMENTAL       "elemental"
565a11e1b2SBarry Smith #define MATBAIJ            "baij"
57273d9f13SBarry Smith #define MATSEQBAIJ         "seqbaij"
58273d9f13SBarry Smith #define MATMPIBAIJ         "mpibaij"
59273d9f13SBarry Smith #define MATMPIADJ          "mpiadj"
605a11e1b2SBarry Smith #define MATSBAIJ           "sbaij"
61273d9f13SBarry Smith #define MATSEQSBAIJ        "seqsbaij"
62273d9f13SBarry Smith #define MATMPISBAIJ        "mpisbaij"
63cebc7f6cSBarry Smith #define MATDAAD            "daad"
64cebc7f6cSBarry Smith #define MATMFFD            "mffd"
65c8a8475eSBarry Smith #define MATNORMAL          "normal"
66c5e4d11fSDmitry Karpeev #define MATNORMALHERMITIAN "normalh"
67ab92ecdeSBarry Smith #define MATLRC             "lrc"
682a6744ebSBarry Smith #define MATSCATTER         "scatter"
69421e10b8SBarry Smith #define MATBLOCKMAT        "blockmat"
70793850ffSBarry Smith #define MATCOMPOSITE       "composite"
711f8c7532SHong Zhang #define MATFFT             "fft"
721f8c7532SHong Zhang #define MATFFTW            "fftw"
73e133240eSMatthew G Knepley #define MATSEQCUFFT        "seqcufft"
74557cca28SSatish Balay #define MATTRANSPOSEMAT    "transpose"
7572ca8751SBarry Smith #define MATSCHURCOMPLEMENT "schurcomplement"
761d6018f0SLisandro Dalcin #define MATPYTHON          "python"
7763c07aadSStefano Zampini #define MATHYPRE           "hypre"
78f91d8e95SBarry Smith #define MATHYPRESTRUCT     "hyprestruct"
79a9e6138eSGlenn Hammond #define MATHYPRESSTRUCT    "hypresstruct"
80ee1cef2cSJed Brown #define MATSUBMATRIX       "submatrix"
812c0dbf93SJed Brown #define MATLOCALREF        "localref"
82d8588912SDave May #define MATNEST            "nest"
83c094ef40SMatthew G. Knepley #define MATPREALLOCATOR    "preallocator"
84a3b2e22bSHong Zhang #define MATDUMMY           "dummy"
85773941d6SBarry Smith 
8676bdecfbSBarry Smith /*J
87c7393fdbSBarry Smith     MatSolverPackage - String with the name of a PETSc matrix solver type.
889989ab13SBarry Smith 
89c2b89b5dSBarry Smith     For example: "petsc" indicates what PETSc provides, "superlu_dist" the parallel SuperLU_DIST package etc
909989ab13SBarry Smith 
919989ab13SBarry Smith    Level: beginner
929989ab13SBarry Smith 
935c9eb25fSBarry Smith .seealso: MatGetFactor(), Mat, MatSetType(), MatType
9476bdecfbSBarry Smith J*/
95c7393fdbSBarry Smith #define MatSolverPackage char*
962692d6eeSBarry Smith #define MATSOLVERSUPERLU          "superlu"
972692d6eeSBarry Smith #define MATSOLVERSUPERLU_DIST     "superlu_dist"
9808f5efcfSPieter Ghysels #define MATSOLVERSTRUMPACK        "strumpack"
992692d6eeSBarry Smith #define MATSOLVERUMFPACK          "umfpack"
10020db9a53SJed Brown #define MATSOLVERCHOLMOD          "cholmod"
101d89f5e7aSBarry Smith #define MATSOLVERKLU              "klu"
102418810c4SBarry Smith #define MATSOLVERSPARSEELEMENTAL  "sparseelemental"
103d89f5e7aSBarry Smith #define MATSOLVERELEMENTAL        "elemental"
1042692d6eeSBarry Smith #define MATSOLVERESSL             "essl"
1052692d6eeSBarry Smith #define MATSOLVERLUSOL            "lusol"
1062692d6eeSBarry Smith #define MATSOLVERMUMPS            "mumps"
107d615f992SSatish Balay #define MATSOLVERMKL_PARDISO      "mkl_pardiso"
108d305a81bSVasiliy Kozyrev #define MATSOLVERMKL_CPARDISO     "mkl_cpardiso"
1092692d6eeSBarry Smith #define MATSOLVERPASTIX           "pastix"
1102692d6eeSBarry Smith #define MATSOLVERMATLAB           "matlab"
1112692d6eeSBarry Smith #define MATSOLVERPETSC            "petsc"
1122692d6eeSBarry Smith #define MATSOLVERBAS              "bas"
1139ae82921SPaul Mullowney #define MATSOLVERCUSPARSE         "cusparse"
114c0cdd4a1SDahai Guo 
115b24902e0SBarry Smith /*E
116b24902e0SBarry Smith     MatFactorType - indicates what type of factorization is requested
117b24902e0SBarry Smith 
118b24902e0SBarry Smith     Level: beginner
119b24902e0SBarry Smith 
120af0996ceSBarry Smith    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
121b24902e0SBarry Smith 
122c7393fdbSBarry Smith .seealso: MatSolverPackage, MatGetFactor()
123b24902e0SBarry Smith E*/
124599ef60dSHong Zhang typedef enum {MAT_FACTOR_NONE, MAT_FACTOR_LU, MAT_FACTOR_CHOLESKY, MAT_FACTOR_ILU, MAT_FACTOR_ICC,MAT_FACTOR_ILUDT} MatFactorType;
125014dd563SJed Brown PETSC_EXTERN const char *const MatFactorTypes[];
126e92e720dSBarry Smith 
127014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactor(Mat,const MatSolverPackage,MatFactorType,Mat*);
128014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactorAvailable(Mat,const MatSolverPackage,MatFactorType,PetscBool *);
129014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFactorGetSolverPackage(Mat,const MatSolverPackage*);
130014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactorType(Mat,MatFactorType*);
13118713533SBarry Smith PETSC_EXTERN PetscErrorCode MatSolverPackageRegister(const MatSolverPackage,const MatType,MatFactorType,PetscErrorCode(*)(Mat,MatFactorType,Mat*));
13242c9c57cSBarry Smith PETSC_EXTERN PetscErrorCode MatSolverPackageGet(const MatSolverPackage,const MatType,MatFactorType,PetscBool*,PetscBool*,PetscErrorCode (**)(Mat,MatFactorType,Mat*));
1339989ab13SBarry Smith 
134c06d978dSMatthew Knepley /* Logging support */
1350700a824SBarry Smith #define    MAT_FILE_CLASSID 1211216    /* used to indicate matrices in binary files */
136014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_CLASSID;
137335efc43SPeter Brune PETSC_EXTERN PetscClassId MAT_COLORING_CLASSID;
138014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_FDCOLORING_CLASSID;
139014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_TRANSPOSECOLORING_CLASSID;
140014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_PARTITIONING_CLASSID;
141014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_COARSEN_CLASSID;
142014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_NULLSPACE_CLASSID;
143014dd563SJed Brown PETSC_EXTERN PetscClassId MATMFFD_CLASSID;
144c06d978dSMatthew Knepley 
145ceb03754SKris Buschelman /*E
146cf37664fSBarry Smith     MatReuse - Indicates if matrices obtained from a previous call to MatGetSubMatrices(), MatGetSubMatrix(), MatConvert() or several other functions
147cf37664fSBarry Smith      are to be reused to store the new matrix values.
148cf37664fSBarry Smith 
149cf37664fSBarry Smith $  MAT_INITIAL_MATRIX - create a new matrix
150cf37664fSBarry Smith $  MAT_REUSE_MATRIX - reuse the matrix created with a previous call that used MAT_INITIAL_MATRIX
151cf37664fSBarry Smith $  MAT_INPLACE_MATRIX - replace the first input matrix with the new matrix (not applicable to all functions)
152cf37664fSBarry Smith $  MAT_IGNORE_MATRIX - do not create a new matrix or reuse a give matrix, just ignore that matrix argument (not applicable to all functions)
153ceb03754SKris Buschelman 
154ceb03754SKris Buschelman     Level: beginner
155ceb03754SKris Buschelman 
156af0996ceSBarry Smith    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
157ceb03754SKris Buschelman 
1580c451bc4SBarry Smith .seealso: MatGetSubMatrices(), MatGetSubMatrix(), MatDestroyMatrices(), MatConvert()
159ceb03754SKris Buschelman E*/
160511c6705SHong Zhang typedef enum {MAT_INITIAL_MATRIX,MAT_REUSE_MATRIX,MAT_IGNORE_MATRIX,MAT_INPLACE_MATRIX} MatReuse;
161ceb03754SKris Buschelman 
1625494a064SHong Zhang /*E
1635494a064SHong Zhang     MatGetSubMatrixOption - Indicates if matrices obtained from a call to MatGetSubMatrices()
164829201f2SHong Zhang      include the matrix values. Currently it is only used by MatGetSeqNonzerostructure().
1655494a064SHong Zhang 
1665494a064SHong Zhang     Level: beginner
1675494a064SHong Zhang 
168829201f2SHong Zhang .seealso: MatGetSeqNonzerostructure()
1695494a064SHong Zhang E*/
1705494a064SHong Zhang typedef enum {MAT_DO_NOT_GET_VALUES,MAT_GET_VALUES} MatGetSubMatrixOption;
1715494a064SHong Zhang 
172607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatInitializePackage(void);
173c06d978dSMatthew Knepley 
174014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreate(MPI_Comm,Mat*);
175014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetSizes(Mat,PetscInt,PetscInt,PetscInt,PetscInt);
17619fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatSetType(Mat,MatType);
177014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetFromOptions(Mat);
178685405a1SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatViewFromOptions(Mat A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);}
179bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatRegister(const char[],PetscErrorCode(*)(Mat));
180014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRegisterBaseName(const char[],const char[],const char[]);
181014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetOptionsPrefix(Mat,const char[]);
182014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAppendOptionsPrefix(Mat,const char[]);
183014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOptionsPrefix(Mat,const char*[]);
18484d44b13SHong Zhang PETSC_EXTERN PetscErrorCode MatSetErrorIfFailure(Mat,PetscBool);
185f69a0ea3SMatthew Knepley 
186140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatList;
187140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatColoringList;
188140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatPartitioningList;
18928988994SBarry Smith 
1903b224e63SBarry Smith /*E
191aa6c7ce3SBarry Smith     MatStructure - Indicates if two matrices have the same nonzero structure
1923b224e63SBarry Smith 
1933b224e63SBarry Smith     Level: beginner
1943b224e63SBarry Smith 
195af0996ceSBarry Smith    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
1963b224e63SBarry Smith 
197aa6c7ce3SBarry Smith .seealso: MatCopy(), MatAXPY()
1983b224e63SBarry Smith E*/
199aa6c7ce3SBarry Smith typedef enum {DIFFERENT_NONZERO_PATTERN,SUBSET_NONZERO_PATTERN,SAME_NONZERO_PATTERN} MatStructure;
2003b224e63SBarry Smith 
201014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqDense(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*);
202014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateDense(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*);
203014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
204014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
205014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *);
206014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithSplitArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],PetscInt[],PetscInt[],PetscScalar[],Mat*);
2078d7a6e47SBarry Smith 
208014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
209014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
210014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat*);
211d21a29f3SJed Brown 
212014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAdj(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscInt[],Mat*);
213014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
214d21a29f3SJed Brown 
215014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
216014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPISBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *);
21738f409ebSLisandro Dalcin PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
218014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
2194cce697cSJed Brown PETSC_EXTERN PetscErrorCode MatXAIJSetPreallocation(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscInt[],const PetscInt[]);
220dfb205c3SBarry Smith 
221014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateShell(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,void *,Mat*);
222014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateNormal(Mat,Mat*);
223c5e4d11fSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatCreateNormalHermitian(Mat,Mat*);
224986c4d72SJose E. Roman PETSC_EXTERN PetscErrorCode MatCreateLRC(Mat,Mat,Vec,Mat,Mat*);
225267ca84cSJose E. Roman PETSC_EXTERN PetscErrorCode MatLRCGetMats(Mat,Mat*,Mat*,Vec*,Mat*);
226c5e4d11fSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatCreateIS(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,ISLocalToGlobalMapping,ISLocalToGlobalMapping,Mat*);
227014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
228014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
229c0cdd4a1SDahai Guo 
230014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateScatter(MPI_Comm,VecScatter,Mat*);
231014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScatterSetVecScatter(Mat,VecScatter);
232014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScatterGetVecScatter(Mat,VecScatter*);
233014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateBlockMat(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,Mat*);
234014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeAddMat(Mat,Mat);
235014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeMerge(Mat);
236014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateComposite(MPI_Comm,PetscInt,const Mat*,Mat*);
2376d7c1e57SBarry Smith typedef enum {MAT_COMPOSITE_ADDITIVE,MAT_COMPOSITE_MULTIPLICATIVE} MatCompositeType;
238014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeSetType(Mat,MatCompositeType);
2396d7c1e57SBarry Smith 
24019fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatCreateFFT(MPI_Comm,PetscInt,const PetscInt[],MatType,Mat*);
241014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqCUFFT(MPI_Comm,PetscInt,const PetscInt[],Mat*);
242dedccee8SHong Zhang 
243014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateTranspose(Mat,Mat*);
2448060fb66Sstefano_zampini PETSC_EXTERN PetscErrorCode MatTransposeGetMat(Mat,Mat*);
245d0de2241SAndrew Spott PETSC_EXTERN PetscErrorCode MatCreateHermitianTranspose(Mat,Mat*);
246014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSubMatrix(Mat,IS,IS,Mat*);
247014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSubMatrixUpdate(Mat,Mat,IS,IS);
248014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateLocalRef(Mat,IS,IS,Mat*);
2491472f72bSBarry Smith 
250978814f1SStefano Zampini #if defined(PETSC_HAVE_HYPRE)
251d975228cSstefano_zampini PETSC_EXTERN PetscErrorCode MatHYPRESetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
252978814f1SStefano Zampini #endif
253978814f1SStefano Zampini 
254014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPythonSetType(Mat,const char[]);
2551d6018f0SLisandro Dalcin 
256014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetUp(Mat);
257014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroy(Mat*);
258e56f5c9eSBarry Smith PETSC_EXTERN PetscErrorCode MatGetNonzeroState(Mat,PetscObjectState*);
25921c89e3eSBarry Smith 
260014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatConjugate(Mat);
261014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRealPart(Mat);
262014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatImaginaryPart(Mat);
263014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetDiagonalBlock(Mat,Mat*);
264014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetTrace(Mat,PetscScalar*);
265713ccfa9SJed Brown PETSC_EXTERN PetscErrorCode MatInvertBlockDiagonal(Mat,const PetscScalar **);
26699cafbc1SBarry Smith 
2678ed539a5SBarry Smith /* ------------------------------------------------------------*/
268014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
269014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
270014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesRow(Mat,PetscInt,const PetscScalar[]);
271014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesRowLocal(Mat,PetscInt,const PetscScalar[]);
272014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBatch(Mat,PetscInt,PetscInt,PetscInt[],const PetscScalar[]);
27373a71a0fSBarry Smith PETSC_EXTERN PetscErrorCode MatSetRandom(Mat,PetscRandom);
27484cb2905SBarry Smith 
2752ef4de8bSBarry Smith /*S
2762ef4de8bSBarry Smith      MatStencil - Data structure (C struct) for storing information about a single row or
27762ef0227SBarry Smith         column of a matrix as indexed on an associated grid. These are arguments to MatSetStencil() and MatSetBlockStencil()
27862ef0227SBarry Smith 
27962ef0227SBarry Smith    The i,j, and k represent the logical coordinates over the entire grid (for 2 and 1 dimensional problems the k and j entries are ignored).
28062ef0227SBarry Smith    The c represents the the degrees of freedom at each grid point (the dof argument to DMDASetDOF()). If dof is 1 then this entry is ignored.
28162ef0227SBarry Smith 
28262ef0227SBarry Smith    For stencil access to vectors see DMDAVecGetArray(), DMDAVecGetArrayF90().
283be479b30SJed Brown 
284be479b30SJed Brown    Fortran usage is different, see MatSetValuesStencil() for details.
2852ef4de8bSBarry Smith 
2862ef4de8bSBarry Smith    Level: beginner
2872ef4de8bSBarry Smith 
2882ef4de8bSBarry Smith   Concepts: matrix; linear operator
2892ef4de8bSBarry Smith 
29062ef0227SBarry Smith .seealso:  MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockedStencil(), DMDAVecGetArray(), DMDAVecGetArrayF90()
2912ef4de8bSBarry Smith S*/
292435da068SBarry Smith typedef struct {
293c1ac3661SBarry Smith   PetscInt k,j,i,c;
294435da068SBarry Smith } MatStencil;
2952ef4de8bSBarry Smith 
296014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
297014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
298014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt);
299435da068SBarry Smith 
300d91e6319SBarry Smith /*E
301d91e6319SBarry Smith     MatAssemblyType - Indicates if the matrix is now to be used, or if you plan
30262ef0227SBarry Smith      to continue to add or insert 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,
349c10200c1SHong Zhang               MAT_SUBMAT_SINGLEIS = 21,
350c10200c1SHong Zhang               MAT_OPTION_MAX = 22} MatOption;
351e057df02SPaul Mullowney 
352014dd563SJed Brown PETSC_EXTERN const char *MatOptions[];
353014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetOption(Mat,MatOption,PetscBool);
3547d68702bSBarry Smith PETSC_EXTERN PetscErrorCode MatGetOption(Mat,MatOption,PetscBool*);
35519fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatGetType(Mat,MatType*);
35684cb2905SBarry Smith 
357014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]);
358014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
359014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
360014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowUpperTriangular(Mat);
361014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRowUpperTriangular(Mat);
362014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnVector(Mat,Vec,PetscInt);
3638c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJGetArray(Mat,PetscScalar *[]);
3648c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJRestoreArray(Mat,PetscScalar *[]);
36521e72a00SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJGetMaxRowNonzeros(Mat,PetscInt*);
366bd04181cSBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJSetValuesLocalFast(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
3678397e458SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqSBAIJGetArray(Mat,PetscScalar *[]);
3688397e458SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqSBAIJRestoreArray(Mat,PetscScalar *[]);
3698c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatDenseGetArray(Mat,PetscScalar *[]);
3708c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatDenseRestoreArray(Mat,PetscScalar *[]);
371014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBlockSize(Mat,PetscInt *);
372014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSize(Mat,PetscInt);
373014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBlockSizes(Mat,PetscInt *,PetscInt *);
374014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSizes(Mat,PetscInt,PetscInt);
37533d57670SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSizesFromMats(Mat,Mat,Mat);
3761620fd73SBarry Smith 
377014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMult(Mat,Vec,Vec);
378014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultDiagonalBlock(Mat,Vec,Vec);
379014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultAdd(Mat,Vec,Vec,Vec);
380014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTranspose(Mat,Vec,Vec);
381014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultHermitianTranspose(Mat,Vec,Vec);
382014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsTranspose(Mat,Mat,PetscReal,PetscBool *);
383014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitianTranspose(Mat,Mat,PetscReal,PetscBool *);
384014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeAdd(Mat,Vec,Vec,Vec);
385014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultHermitianTransposeAdd(Mat,Vec,Vec,Vec);
386014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultConstrained(Mat,Vec,Vec);
387014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeConstrained(Mat,Vec,Vec);
388014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatSolve(Mat,Mat,Mat);
389f9426fe0SMark Adams PETSC_EXTERN PetscErrorCode MatResidual(Mat,Vec,Vec,Vec);
390f5edf698SHong Zhang 
391d91e6319SBarry Smith /*E
392d91e6319SBarry Smith     MatDuplicateOption - Indicates if a duplicated sparse matrix should have
393d91e6319SBarry Smith   its numerical values copied over or just its nonzero structure.
394d91e6319SBarry Smith 
395d91e6319SBarry Smith     Level: beginner
396d91e6319SBarry Smith 
397af0996ceSBarry Smith    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
398d91e6319SBarry Smith 
39970dcbbb9SBarry Smith $   MAT_SHARE_NONZERO_PATTERN - the i and j arrays in the new matrix will be shared with the original matrix
40070dcbbb9SBarry Smith $                               this also triggers the MAT_DO_NOT_COPY_VALUES option. This is used when you
40170dcbbb9SBarry Smith $                               have several matrices with the same nonzero pattern.
40270dcbbb9SBarry Smith 
403d91e6319SBarry Smith .seealso: MatDuplicate()
404d91e6319SBarry Smith E*/
40570dcbbb9SBarry Smith typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES,MAT_SHARE_NONZERO_PATTERN} MatDuplicateOption;
4062e8a6d31SBarry Smith 
40719fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatConvert(Mat,MatType,MatReuse,Mat*);
408014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDuplicate(Mat,MatDuplicateOption,Mat*);
40994a9d846SBarry Smith 
410cb5b572fSBarry Smith 
411014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCopy(Mat,Mat,MatStructure);
412014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatView(Mat,PetscViewer);
413014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsSymmetric(Mat,PetscReal,PetscBool *);
414014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsStructurallySymmetric(Mat,PetscBool *);
415014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitian(Mat,PetscReal,PetscBool *);
416014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsSymmetricKnown(Mat,PetscBool *,PetscBool *);
417014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitianKnown(Mat,PetscBool *,PetscBool *);
418014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMissingDiagonal(Mat,PetscBool  *,PetscInt *);
419014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLoad(Mat, PetscViewer);
4207b80b807SBarry Smith 
4211a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool  *);
4221a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool  *);
4231a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool  *);
4241a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool  *);
425d4fbbf0eSBarry Smith 
426d91e6319SBarry Smith /*S
427d91e6319SBarry Smith      MatInfo - Context of matrix information, used with MatGetInfo()
428d91e6319SBarry Smith 
429d91e6319SBarry Smith    In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE
430d91e6319SBarry Smith 
431d91e6319SBarry Smith    Level: intermediate
432d91e6319SBarry Smith 
433d91e6319SBarry Smith   Concepts: matrix^nonzero information
434d91e6319SBarry Smith 
435d9274352SBarry Smith .seealso:  MatGetInfo(), MatInfoType
436d91e6319SBarry Smith S*/
4374e220ebcSLois Curfman McInnes typedef struct {
438b0a32e0cSBarry Smith   PetscLogDouble block_size;                         /* block size */
439b0a32e0cSBarry Smith   PetscLogDouble nz_allocated,nz_used,nz_unneeded;   /* number of nonzeros */
440b0a32e0cSBarry Smith   PetscLogDouble memory;                             /* memory allocated */
441b0a32e0cSBarry Smith   PetscLogDouble assemblies;                         /* number of matrix assemblies called */
442b0a32e0cSBarry Smith   PetscLogDouble mallocs;                            /* number of mallocs during MatSetValues() */
443b0a32e0cSBarry Smith   PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */
444b0a32e0cSBarry Smith   PetscLogDouble factor_mallocs;                     /* number of mallocs during factorization */
4454e220ebcSLois Curfman McInnes } MatInfo;
4464e220ebcSLois Curfman McInnes 
447d9274352SBarry Smith /*E
448d9274352SBarry Smith     MatInfoType - Indicates if you want information about the local part of the matrix,
449d9274352SBarry Smith      the entire parallel matrix or the maximum over all the local parts.
450d9274352SBarry Smith 
451d9274352SBarry Smith     Level: beginner
452d9274352SBarry Smith 
453af0996ceSBarry Smith    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
454d9274352SBarry Smith 
455d9274352SBarry Smith .seealso: MatGetInfo(), MatInfo
456d9274352SBarry Smith E*/
4577b80b807SBarry Smith typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType;
458014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetInfo(Mat,MatInfoType,MatInfo*);
459014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetDiagonal(Mat,Vec);
460014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMax(Mat,Vec,PetscInt[]);
461014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMin(Mat,Vec,PetscInt[]);
462014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMaxAbs(Mat,Vec,PetscInt[]);
463014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMinAbs(Mat,Vec,PetscInt[]);
464014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowSum(Mat,Vec);
465014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTranspose(Mat,MatReuse,Mat*);
466014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatHermitianTranspose(Mat,MatReuse,Mat*);
467014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPermute(Mat,IS,IS,Mat*);
468014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalScale(Mat,Vec,Vec);
469014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalSet(Mat,Vec,InsertMode);
470a52676ecSHong Zhang 
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*);
476a52676ecSHong Zhang PETSC_EXTERN PetscErrorCode MatMatMultEqual(Mat,Mat,Mat,PetscInt,PetscBool*);
477a52676ecSHong Zhang PETSC_EXTERN PetscErrorCode MatTransposeMatMultEqual(Mat,Mat,Mat,PetscInt,PetscBool*);
4787b80b807SBarry Smith 
479014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNorm(Mat,NormType,PetscReal*);
480014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnNorms(Mat,NormType,PetscReal*);
481014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroEntries(Mat);
482014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
483014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsIS(Mat,IS,PetscScalar,Vec,Vec);
484014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsStencil(Mat,PetscInt,const MatStencil [],PetscScalar,Vec,Vec);
485014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsStencil(Mat,PetscInt,const MatStencil[],PetscScalar,Vec,Vec);
486014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumns(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
487014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsIS(Mat,IS,PetscScalar,Vec,Vec);
4887b80b807SBarry Smith 
489014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSize(Mat,PetscInt*,PetscInt*);
490014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalSize(Mat,PetscInt*,PetscInt*);
491014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRange(Mat,PetscInt*,PetscInt*);
492014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRanges(Mat,const PetscInt**);
493014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRangeColumn(Mat,PetscInt*,PetscInt*);
494014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRangesColumn(Mat,const PetscInt**);
495566876eaSJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipIS(Mat,IS*,IS*);
4965ef9f2a5SBarry Smith 
497014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
49853dd7562SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatGetSubMatricesMPI(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
499014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroyMatrices(PetscInt,Mat *[]);
500*df750dc8SHong Zhang PETSC_EXTERN PetscErrorCode MatDestroySubMatrices(PetscInt,Mat *[]);
501014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSubMatrix(Mat,IS,IS,MatReuse,Mat *);
502014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalSubMatrix(Mat,IS,IS,Mat*);
503014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreLocalSubMatrix(Mat,IS,IS,Mat*);
504014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSeqNonzeroStructure(Mat,Mat*);
505014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroySeqNonzeroStructure(Mat*);
5067b80b807SBarry Smith 
507014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJ(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*);
508014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJSymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*);
509014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJNumeric(Mat,Mat);
510014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMat(Mat,MatReuse,Mat*);
511014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*);
512014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,Mat*);
513014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetGhosts(Mat, PetscInt *,const PetscInt *[]);
5148efafbd8SBarry Smith 
515014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt);
516aa1e27eaSFande Kong PETSC_EXTERN PetscErrorCode MatIncreaseOverlapSplit(Mat mat,PetscInt n,IS is[],PetscInt ov);
517d2b2a242SBarry Smith PETSC_EXTERN PetscErrorCode MatMPIAIJSetUseScalableIncreaseOverlap(Mat,PetscBool);
5187b80b807SBarry Smith 
519014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*);
520014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMultSymbolic(Mat,Mat,PetscReal,Mat*);
521014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMultNumeric(Mat,Mat,Mat);
52222440eb1SKris Buschelman 
5237bab7c10SHong Zhang PETSC_EXTERN PetscErrorCode MatMatMatMult(Mat,Mat,Mat,MatReuse,PetscReal,Mat*);
5242df6c5c3SPatrick Farrell PETSC_EXTERN PetscErrorCode MatGalerkin(Mat,Mat,Mat,MatReuse,PetscReal,Mat*);
5257bab7c10SHong Zhang 
526014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*);
527014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAPSymbolic(Mat,Mat,PetscReal,Mat*);
528014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAPNumeric(Mat,Mat,Mat);
529014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARt(Mat,Mat,MatReuse,PetscReal,Mat*);
530014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARtSymbolic(Mat,Mat,PetscReal,Mat*);
531014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARtNumeric(Mat,Mat,Mat);
53222440eb1SKris Buschelman 
533014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeMatMult(Mat,Mat,MatReuse,PetscReal,Mat*);
534014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposetMatMultSymbolic(Mat,Mat,PetscReal,Mat*);
535014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposetMatMultNumeric(Mat,Mat,Mat);
536014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMult(Mat,Mat,MatReuse,PetscReal,Mat*);
537014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMultSymbolic(Mat,Mat,PetscReal,Mat*);
538014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMultNumeric(Mat,Mat,Mat);
539bc011b1eSHong Zhang 
540014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAXPY(Mat,PetscScalar,Mat,MatStructure);
541014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAYPX(Mat,PetscScalar,Mat,MatStructure);
5421c741599SHong Zhang 
543014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScale(Mat,PetscScalar);
544014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShift(Mat,PetscScalar);
5457b80b807SBarry Smith 
546014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping);
547014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*);
548a93ff8c4SPeter Brune PETSC_EXTERN PetscErrorCode MatGetLayouts(Mat,PetscLayout*,PetscLayout*);
549014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
550014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsLocalIS(Mat,IS,PetscScalar,Vec,Vec);
551014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
552014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocalIS(Mat,IS,PetscScalar,Vec,Vec);
553014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
554014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
555052efed2SBarry Smith 
556014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStashSetInitialSize(Mat,PetscInt,PetscInt);
557014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
55890f02eecSBarry Smith 
559014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInterpolate(Mat,Vec,Vec);
560014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInterpolateAdd(Mat,Vec,Vec,Vec);
561014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestrict(Mat,Vec,Vec);
5622a7a6963SBarry Smith PETSC_EXTERN PetscErrorCode MatCreateVecs(Mat,Vec*,Vec*);
56384cc2a5aSBarry Smith PETSC_DEPRECATED("Use MatCreateVecs()") PETSC_STATIC_INLINE PetscErrorCode MatGetVecs(Mat mat,Vec *x,Vec *y) {return MatCreateVecs(mat,x,y);}
56453cd1579SHong Zhang PETSC_EXTERN PetscErrorCode MatCreateRedundantMatrix(Mat,PetscInt,MPI_Comm,MatReuse,Mat*);
565014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetMultiProcBlock(Mat,MPI_Comm,MatReuse,Mat*);
566014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFindZeroDiagonals(Mat,IS*);
5673a062f41SBarry Smith PETSC_EXTERN PetscErrorCode MatFindOffBlockDiagonalEntries(Mat,IS*);
5689c8f2541SHong Zhang PETSC_EXTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
569bd481603SBarry Smith 
570bd481603SBarry Smith /*MC
5711620fd73SBarry Smith    MatSetValue - Set a single entry into a matrix.
5721620fd73SBarry Smith 
5731620fd73SBarry Smith    Not collective
5741620fd73SBarry Smith 
575a9834a7dSSatish Balay    Synopsis:
576a9834a7dSSatish Balay      #include <petscmat.h>
577a9834a7dSSatish Balay      PetscErrorCode MatSetValue(Mat m,PetscInt row,PetscInt col,PetscScalar value,InsertMode mode)
578a9834a7dSSatish Balay 
5791620fd73SBarry Smith    Input Parameters:
5801620fd73SBarry Smith +  m - the matrix
5811620fd73SBarry Smith .  row - the row location of the entry
5821620fd73SBarry Smith .  col - the column location of the entry
5831620fd73SBarry Smith .  value - the value to insert
5841620fd73SBarry Smith -  mode - either INSERT_VALUES or ADD_VALUES
5851620fd73SBarry Smith 
5861620fd73SBarry Smith    Notes:
5871620fd73SBarry Smith    For efficiency one should use MatSetValues() and set several or many
5881620fd73SBarry Smith    values simultaneously if possible.
5891620fd73SBarry Smith 
5901620fd73SBarry Smith    Level: beginner
5911620fd73SBarry Smith 
5921620fd73SBarry Smith .seealso: MatSetValues(), MatSetValueLocal()
5931620fd73SBarry Smith M*/
5941620fd73SBarry 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);}
5951620fd73SBarry Smith 
5961620fd73SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatGetValue(Mat v,PetscInt i,PetscInt j,PetscScalar *va) {return MatGetValues(v,1,&i,1,&j,va);}
5971620fd73SBarry Smith 
5981620fd73SBarry 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);}
5991620fd73SBarry Smith 
6001620fd73SBarry Smith /*MC
6010d2a7ff2SSatish Balay    MatPreallocateInitialize - Begins the block of code that will count the number of nonzeros per
602bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
603bd481603SBarry Smith 
604bd481603SBarry Smith    Synopsis:
605aaa7dc30SBarry Smith    #include <petscmat.h>
606c1ac3661SBarry Smith    PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)
607bd481603SBarry Smith 
608bd481603SBarry Smith    Collective on MPI_Comm
609bd481603SBarry Smith 
610bd481603SBarry Smith    Input Parameters:
611bd481603SBarry Smith +  comm - the communicator that will share the eventually allocated matrix
612859fcb39SBarry Smith .  nrows - the number of LOCAL rows in the matrix
613859fcb39SBarry Smith -  ncols - the number of LOCAL columns in the matrix
614bd481603SBarry Smith 
615bd481603SBarry Smith    Output Parameters:
616bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
617bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
618bd481603SBarry Smith 
619bd481603SBarry Smith    Level: intermediate
620bd481603SBarry Smith 
621bd481603SBarry Smith    Notes:
622a7f22e61SSatish Balay     See Users-Manual: ch_performance for more details.
623bd481603SBarry Smith 
6241d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
625bd481603SBarry Smith 
626bd481603SBarry Smith    Use MatPreallocateInitializeSymmetric() for symmetric matrices (MPISBAIJ matrices)
627bd481603SBarry Smith 
6281620fd73SBarry Smith    This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize().
6291620fd73SBarry Smith 
630bd481603SBarry Smith   Concepts: preallocation^Matrix
631bd481603SBarry Smith 
632d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(),
633d6e23781SBarry Smith           MatPreallocateInitializeSymmetric(), MatPreallocateSymmetricSetLocalBlock()
634bd481603SBarry Smith M*/
635c4f061fbSSatish Balay #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \
6367c922b88SBarry Smith { \
637a89bc211SBarry Smith   PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__start,__end; \
6381795a4d1SJed Brown   _4_ierr = PetscCalloc2(__nrows,&dnz,__nrows,&onz);CHKERRQ(_4_ierr); \
6391795a4d1SJed Brown   __start = 0; __end = __start;                                         \
640c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ctmp;\
641a89bc211SBarry Smith   _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows;
6427c922b88SBarry Smith 
643bd481603SBarry Smith /*MC
644bd481603SBarry Smith    MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
645bd481603SBarry Smith        inserted using a local number of the rows and columns
646bd481603SBarry Smith 
647bd481603SBarry Smith    Synopsis:
648aaa7dc30SBarry Smith    #include <petscmat.h>
649c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
650bd481603SBarry Smith 
651bd481603SBarry Smith    Not Collective
652bd481603SBarry Smith 
653bd481603SBarry Smith    Input Parameters:
654784ac674SJed Brown +  map - the row mapping from local numbering to global numbering
655bd481603SBarry Smith .  nrows - the number of rows indicated
6561d73ed98SBarry Smith .  rows - the indices of the rows
657784ac674SJed Brown .  cmap - the column mapping from local to global numbering
658bd481603SBarry Smith .  ncols - the number of columns in the matrix
659bd481603SBarry Smith .  cols - the columns indicated
660bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
661bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
662bd481603SBarry Smith 
663bd481603SBarry Smith    Level: intermediate
664bd481603SBarry Smith 
665bd481603SBarry Smith    Notes:
666a7f22e61SSatish Balay     See Users-Manual: ch_performance for more details.
667bd481603SBarry Smith 
6681d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
669bd481603SBarry Smith 
670bd481603SBarry Smith   Concepts: preallocation^Matrix
671bd481603SBarry Smith 
672d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(),
673c1154cd5SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocalRemoveDups()
674bd481603SBarry Smith M*/
675784ac674SJed Brown #define MatPreallocateSetLocal(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \
676c4f061fbSSatish Balay {\
677c1ac3661SBarry Smith   PetscInt __l;\
678784ac674SJed Brown   _4_ierr = ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\
679784ac674SJed Brown   _4_ierr = ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\
680c4f061fbSSatish Balay   for (__l=0;__l<nrows;__l++) {\
681ef66eb69SBarry Smith     _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
682c4f061fbSSatish Balay   }\
683c4f061fbSSatish Balay }
684c4f061fbSSatish Balay 
685bd481603SBarry Smith /*MC
686c1154cd5SBarry Smith    MatPreallocateSetLocalRemoveDups - Indicates the locations (rows and columns) in the matrix where nonzeros will be
687c1154cd5SBarry Smith        inserted using a local number of the rows and columns. This version removes any duplicate columns in cols
688c1154cd5SBarry Smith 
689c1154cd5SBarry Smith    Synopsis:
690c1154cd5SBarry Smith    #include <petscmat.h>
691c1154cd5SBarry Smith    PetscErrorCode MatPreallocateSetLocalRemoveDups(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
692c1154cd5SBarry Smith 
693c1154cd5SBarry Smith    Not Collective
694c1154cd5SBarry Smith 
695c1154cd5SBarry Smith    Input Parameters:
696c1154cd5SBarry Smith +  map - the row mapping from local numbering to global numbering
697c1154cd5SBarry Smith .  nrows - the number of rows indicated
698c1154cd5SBarry Smith .  rows - the indices of the rows (these values are mapped to the global values)
699c1154cd5SBarry Smith .  cmap - the column mapping from local to global numbering
700c1154cd5SBarry Smith .  ncols - the number of columns in the matrix   (this value will be changed if duplicate columns are found)
701c1154cd5SBarry Smith .  cols - the columns indicated (these values are mapped to the global values, they are then sorted and duplicates removed)
702c1154cd5SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
703c1154cd5SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
704c1154cd5SBarry Smith 
705c1154cd5SBarry Smith    Level: intermediate
706c1154cd5SBarry Smith 
707c1154cd5SBarry Smith    Notes:
708c1154cd5SBarry Smith     See Users-Manual: ch_performance for more details.
709c1154cd5SBarry Smith 
710c1154cd5SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
711c1154cd5SBarry Smith 
712c1154cd5SBarry Smith   Concepts: preallocation^Matrix
713c1154cd5SBarry Smith 
714c1154cd5SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(),
715c1154cd5SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocal()
716c1154cd5SBarry Smith M*/
717c1154cd5SBarry Smith #define MatPreallocateSetLocalRemoveDups(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \
718c1154cd5SBarry Smith {\
719c1154cd5SBarry Smith   PetscInt __l;\
720c1154cd5SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\
721c1154cd5SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\
722c1154cd5SBarry Smith   _4_ierr = PetscSortRemoveDupsInt(&ncols,cols);CHKERRQ(_4_ierr);\
723c1154cd5SBarry Smith   for (__l=0;__l<nrows;__l++) {\
724c1154cd5SBarry Smith     _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
725c1154cd5SBarry Smith   }\
726c1154cd5SBarry Smith }
727c1154cd5SBarry Smith 
728c1154cd5SBarry Smith /*MC
729d6e23781SBarry Smith    MatPreallocateSetLocalBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be
730bd481603SBarry Smith        inserted using a local number of the rows and columns
731bd481603SBarry Smith 
732bd481603SBarry Smith    Synopsis:
733aaa7dc30SBarry Smith    #include <petscmat.h>
734d6e23781SBarry Smith    PetscErrorCode MatPreallocateSetLocalBlock(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
735d6e23781SBarry Smith 
736d6e23781SBarry Smith    Not Collective
737d6e23781SBarry Smith 
738d6e23781SBarry Smith    Input Parameters:
739d6e23781SBarry Smith +  map - the row mapping from local numbering to global numbering
740d6e23781SBarry Smith .  nrows - the number of rows indicated
741d6e23781SBarry Smith .  rows - the indices of the rows
742d6e23781SBarry Smith .  cmap - the column mapping from local to global numbering
743d6e23781SBarry Smith .  ncols - the number of columns in the matrix
744d6e23781SBarry Smith .  cols - the columns indicated
745d6e23781SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
746d6e23781SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
747d6e23781SBarry Smith 
748d6e23781SBarry Smith    Level: intermediate
749d6e23781SBarry Smith 
750d6e23781SBarry Smith    Notes:
751d6e23781SBarry Smith     See Users-Manual: ch_performance for more details.
752d6e23781SBarry Smith 
753d6e23781SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
754d6e23781SBarry Smith 
755d6e23781SBarry Smith   Concepts: preallocation^Matrix
756d6e23781SBarry Smith 
757d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(),
758d6e23781SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock()
759d6e23781SBarry Smith M*/
760d6e23781SBarry Smith #define MatPreallocateSetLocalBlock(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \
761d6e23781SBarry Smith {\
762d6e23781SBarry Smith   PetscInt __l;\
763d6e23781SBarry Smith   _4_ierr = ISLocalToGlobalMappingApplyBlock(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\
764d6e23781SBarry Smith   _4_ierr = ISLocalToGlobalMappingApplyBlock(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\
765d6e23781SBarry Smith   for (__l=0;__l<nrows;__l++) {\
766d6e23781SBarry Smith     _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
767d6e23781SBarry Smith   }\
768d6e23781SBarry Smith }
769d6e23781SBarry Smith 
770d6e23781SBarry Smith /*MC
771d6e23781SBarry Smith    MatPreallocateSymmetricSetLocalBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be
772d6e23781SBarry Smith        inserted using a local number of the rows and columns
773d6e23781SBarry Smith 
774d6e23781SBarry Smith    Synopsis:
775d6e23781SBarry Smith    #include <petscmat.h>
776d6e23781SBarry Smith    PetscErrorCode MatPreallocateSymmetricSetLocalBlock(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
777bd481603SBarry Smith 
778bd481603SBarry Smith    Not Collective
779bd481603SBarry Smith 
780bd481603SBarry Smith    Input Parameters:
781bd481603SBarry Smith +  map - the mapping between local numbering and global numbering
782bd481603SBarry Smith .  nrows - the number of rows indicated
7831d73ed98SBarry Smith .  rows - the indices of the rows
784bd481603SBarry Smith .  ncols - the number of columns in the matrix
785bd481603SBarry Smith .  cols - the columns indicated
786bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
787bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
788bd481603SBarry Smith 
789bd481603SBarry Smith    Level: intermediate
790bd481603SBarry Smith 
791bd481603SBarry Smith    Notes:
792a7f22e61SSatish Balay     See Users-Manual: ch_performance for more details.
793bd481603SBarry Smith 
794bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
795bd481603SBarry Smith 
796bd481603SBarry Smith   Concepts: preallocation^Matrix
797bd481603SBarry Smith 
798d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateInitialize(),
799d6e23781SBarry Smith           MatPreallocateInitialize(),  MatPreallocateSetLocal()
800bd481603SBarry Smith M*/
801d6e23781SBarry Smith #define MatPreallocateSymmetricSetLocalBlock(map,nrows,rows,ncols,cols,dnz,onz) 0;\
802d3d32019SBarry Smith {\
803c1ac3661SBarry Smith   PetscInt __l;\
804d6e23781SBarry Smith   _4_ierr = ISLocalToGlobalMappingApplyBlock(map,nrows,rows,rows);CHKERRQ(_4_ierr);\
805d6e23781SBarry Smith   _4_ierr = ISLocalToGlobalMappingApplyBlock(map,ncols,cols,cols);CHKERRQ(_4_ierr);\
806d3d32019SBarry Smith   for (__l=0;__l<nrows;__l++) {\
807d6e23781SBarry Smith     _4_ierr = MatPreallocateSymmetricSetBlock((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
808d3d32019SBarry Smith   }\
809d3d32019SBarry Smith }
810bd481603SBarry Smith /*MC
811bd481603SBarry Smith    MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
812bd481603SBarry Smith        inserted using a local number of the rows and columns
813bd481603SBarry Smith 
814bd481603SBarry Smith    Synopsis:
815aaa7dc30SBarry Smith    #include <petscmat.h>
816c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
817bd481603SBarry Smith 
818bd481603SBarry Smith    Not Collective
819bd481603SBarry Smith 
820bd481603SBarry Smith    Input Parameters:
82164075487SBarry Smith +  row - the row
822bd481603SBarry Smith .  ncols - the number of columns in the matrix
823a50f8bf6SHong Zhang -  cols - the columns indicated
824a50f8bf6SHong Zhang 
825a50f8bf6SHong Zhang    Output Parameters:
826a50f8bf6SHong Zhang +  dnz - the array that will be passed to the matrix preallocation routines
827bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
828bd481603SBarry Smith 
829bd481603SBarry Smith    Level: intermediate
830bd481603SBarry Smith 
831bd481603SBarry Smith    Notes:
832a7f22e61SSatish Balay     See Users-Manual: ch_performance for more details.
833bd481603SBarry Smith 
834bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
835bd481603SBarry Smith 
8361620fd73SBarry Smith    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
8371620fd73SBarry Smith 
838bd481603SBarry Smith   Concepts: preallocation^Matrix
839bd481603SBarry Smith 
840d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(),
841d6e23781SBarry Smith           MatPreallocateInitialize(), MatPreallocateSetLocal()
842bd481603SBarry Smith M*/
843c4f061fbSSatish Balay #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\
844c1ac3661SBarry Smith { PetscInt __i; \
845e32f2f54SBarry 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);\
846e32f2f54SBarry 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);\
8477c922b88SBarry Smith   for (__i=0; __i<nc; __i++) {\
84864075487SBarry Smith     if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \
8497cc688d7SBarry Smith     else dnz[row - __rstart]++;\
8507c922b88SBarry Smith   }\
8517c922b88SBarry Smith }
8527c922b88SBarry Smith 
853bd481603SBarry Smith /*MC
854d6e23781SBarry Smith    MatPreallocateSymmetricSetBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be
855bd481603SBarry Smith        inserted using a local number of the rows and columns
856bd481603SBarry Smith 
857bd481603SBarry Smith    Synopsis:
858aaa7dc30SBarry Smith    #include <petscmat.h>
859d6e23781SBarry Smith    PetscErrorCode MatPreallocateSymmetricSetBlock(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
860bd481603SBarry Smith 
861bd481603SBarry Smith    Not Collective
862bd481603SBarry Smith 
863bd481603SBarry Smith    Input Parameters:
864bd481603SBarry Smith +  nrows - the number of rows indicated
8651d73ed98SBarry Smith .  rows - the indices of the rows
866bd481603SBarry Smith .  ncols - the number of columns in the matrix
867bd481603SBarry Smith .  cols - the columns indicated
868bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
869bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
870bd481603SBarry Smith 
871bd481603SBarry Smith    Level: intermediate
872bd481603SBarry Smith 
873bd481603SBarry Smith    Notes:
874a7f22e61SSatish Balay     See Users-Manual: ch_performance for more details.
875bd481603SBarry Smith 
876bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
877bd481603SBarry Smith 
8781620fd73SBarry Smith    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
8791620fd73SBarry Smith 
880bd481603SBarry Smith   Concepts: preallocation^Matrix
881bd481603SBarry Smith 
882d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(),  MatPreallocateInitialize(),
883d6e23781SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocal()
884bd481603SBarry Smith M*/
885d6e23781SBarry Smith #define MatPreallocateSymmetricSetBlock(row,nc,cols,dnz,onz) 0;\
886c1ac3661SBarry Smith { PetscInt __i; \
887d3d32019SBarry Smith   for (__i=0; __i<nc; __i++) {\
888d3d32019SBarry Smith     if (cols[__i] >= __end) onz[row - __rstart]++; \
889d3d32019SBarry Smith     else if (cols[__i] >= row) dnz[row - __rstart]++;\
890d3d32019SBarry Smith   }\
891d3d32019SBarry Smith }
892d3d32019SBarry Smith 
893bd481603SBarry Smith /*MC
89416371a99SBarry Smith    MatPreallocateLocation -  An alternative to MatPreallocationSet() that puts the nonzero locations into the matrix if it exists
89516371a99SBarry Smith 
89616371a99SBarry Smith    Synopsis:
897aaa7dc30SBarry Smith    #include <petscmat.h>
89816371a99SBarry Smith    PetscErrorCode MatPreallocateLocations(Mat A,PetscInt row,PetscInt ncols,PetscInt *cols,PetscInt *dnz,PetscInt *onz)
89916371a99SBarry Smith 
90016371a99SBarry Smith    Not Collective
90116371a99SBarry Smith 
90216371a99SBarry Smith    Input Parameters:
90316371a99SBarry Smith .  A - matrix
90416371a99SBarry Smith .  row - row where values exist (must be local to this process)
90516371a99SBarry Smith .  ncols - number of columns
90616371a99SBarry Smith .  cols - columns with nonzeros
90716371a99SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
90816371a99SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
90916371a99SBarry Smith 
91016371a99SBarry Smith    Level: intermediate
91116371a99SBarry Smith 
91216371a99SBarry Smith    Notes:
913a7f22e61SSatish Balay     See Users-Manual: ch_performance for more details.
91416371a99SBarry Smith 
91516371a99SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
91616371a99SBarry Smith 
91716371a99SBarry Smith    This is a MACRO not a function because it uses a bunch of variables private to the MatPreallocation.... routines.
91816371a99SBarry Smith 
91916371a99SBarry Smith   Concepts: preallocation^Matrix
92016371a99SBarry Smith 
921d6e23781SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(),
922d6e23781SBarry Smith           MatPreallocateSymmetricSetLocalBlock()
92316371a99SBarry Smith M*/
9240298fd71SBarry 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);}
92516371a99SBarry Smith 
92616371a99SBarry Smith 
92716371a99SBarry Smith /*MC
9280d2a7ff2SSatish Balay    MatPreallocateFinalize - Ends the block of code that will count the number of nonzeros per
929bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
930bd481603SBarry Smith 
931bd481603SBarry Smith    Synopsis:
932aaa7dc30SBarry Smith    #include <petscmat.h>
933c1ac3661SBarry Smith    PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz)
934bd481603SBarry Smith 
935bd481603SBarry Smith    Collective on MPI_Comm
936bd481603SBarry Smith 
937bd481603SBarry Smith    Input Parameters:
93816371a99SBarry Smith +  dnz - the array that was be passed to the matrix preallocation routines
939bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
940bd481603SBarry Smith 
941bd481603SBarry Smith    Level: intermediate
942bd481603SBarry Smith 
943bd481603SBarry Smith    Notes:
944a7f22e61SSatish Balay     See Users-Manual: ch_performance for more details.
945bd481603SBarry Smith 
946bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
947bd481603SBarry Smith 
9481620fd73SBarry Smith    This is a MACRO not a function because it closes the { started in MatPreallocateInitialize().
9491620fd73SBarry Smith 
950bd481603SBarry Smith   Concepts: preallocation^Matrix
951bd481603SBarry Smith 
952d6e23781SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(),
953d6e23781SBarry Smith           MatPreallocateSymmetricSetLocalBlock()
954bd481603SBarry Smith M*/
955a89bc211SBarry Smith #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree2(dnz,onz);CHKERRQ(_4_ierr);}
9567c922b88SBarry Smith 
9577b80b807SBarry Smith /* Routines unique to particular data structures */
958014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellGetContext(Mat,void *);
959698d4c6aSKris Buschelman 
960014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInodeAdjustForInodes(Mat,IS*,IS*);
961014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *);
962698d4c6aSKris Buschelman 
963014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetColumnIndices(Mat,PetscInt[]);
964014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetColumnIndices(Mat,PetscInt[]);
965014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
966014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
967014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
968014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJFromTriple(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*,PetscInt,PetscBool);
9697b80b807SBarry Smith 
970a96a251dSBarry Smith #define MAT_SKIP_ALLOCATION -4
971a96a251dSBarry Smith 
972014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
973014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
974014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]);
975273d9f13SBarry Smith 
976014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
977014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
978014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
979014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []);
980014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
981014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]);
982014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
983014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]);
984014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIDenseSetPreallocation(Mat,PetscScalar[]);
985014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqDenseSetPreallocation(Mat,PetscScalar[]);
9869230625dSJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,const PetscInt*[]);
9879230625dSJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,const PetscInt*[]);
988014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAdjCreateNonemptySubcommMat(Mat,Mat*);
989273d9f13SBarry Smith 
9902e1947a5SStefano Zampini PETSC_EXTERN PetscErrorCode MatISSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
991cf0a3239SStefano Zampini PETSC_EXTERN PetscErrorCode MatISSetUpSF(Mat);
992014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqDenseSetLDA(Mat,PetscInt);
993014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDenseGetLocalMatrix(Mat,Mat*);
9941b807ce4Svictorle 
995014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStoreValues(Mat);
996014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRetrieveValues(Mat);
9972e8a6d31SBarry Smith 
998014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDAADSetCtx(Mat,void*);
9993a7fca6bSBarry Smith 
1000014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFindNonzeroRows(Mat,IS*);
10017cfe41e4SPatrick Farrell PETSC_EXTERN PetscErrorCode MatFindZeroRows(Mat,IS*);
10027b80b807SBarry Smith /*
10037b80b807SBarry Smith   These routines are not usually accessed directly, rather solving is
100494b7f48cSBarry Smith   done through the KSP and PC interfaces.
10057b80b807SBarry Smith */
10067b80b807SBarry Smith 
100776bdecfbSBarry Smith /*J
10088f6c3df8SBarry Smith     MatOrderingType - String with the name of a PETSc matrix ordering
1009d9274352SBarry Smith 
1010d9274352SBarry Smith    Level: beginner
1011d9274352SBarry Smith 
1012d9274352SBarry Smith .seealso: MatGetOrdering()
101376bdecfbSBarry Smith J*/
101419fd82e9SBarry Smith typedef const char* MatOrderingType;
10152692d6eeSBarry Smith #define MATORDERINGNATURAL     "natural"
10162692d6eeSBarry Smith #define MATORDERINGND          "nd"
10172692d6eeSBarry Smith #define MATORDERING1WD         "1wd"
10182692d6eeSBarry Smith #define MATORDERINGRCM         "rcm"
10192692d6eeSBarry Smith #define MATORDERINGQMD         "qmd"
10202692d6eeSBarry Smith #define MATORDERINGROWLENGTH   "rowlength"
1021510184a7SJed Brown #define MATORDERINGWBM         "wbm"
1022fc1bef75SJed Brown #define MATORDERINGSPECTRAL    "spectral"
1023312e372aSBarry Smith #define MATORDERINGAMD         "amd"            /* only works if UMFPACK is installed with PETSc */
1024b12f92e5SBarry Smith 
102519fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatGetOrdering(Mat,MatOrderingType,IS*,IS*);
1026140e18c1SBarry Smith PETSC_EXTERN PetscErrorCode MatGetOrderingList(PetscFunctionList*);
1027bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatOrderingRegister(const char[],PetscErrorCode(*)(Mat,MatOrderingType,IS*,IS*));
1028140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatOrderingList;
1029d4fbbf0eSBarry Smith 
1030014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS);
1031fc1bef75SJed Brown PETSC_EXTERN PetscErrorCode MatCreateLaplacian(Mat,PetscReal,PetscBool,Mat*);
1032a2ce50c7SBarry Smith 
1033d91e6319SBarry Smith /*S
1034d90ac83dSHong Zhang     MatFactorShiftType - Numeric Shift.
1035d90ac83dSHong Zhang 
1036d90ac83dSHong Zhang    Level: beginner
1037d90ac83dSHong Zhang 
1038d90ac83dSHong Zhang S*/
1039d90ac83dSHong Zhang typedef enum {MAT_SHIFT_NONE,MAT_SHIFT_NONZERO,MAT_SHIFT_POSITIVE_DEFINITE,MAT_SHIFT_INBLOCKS} MatFactorShiftType;
10406a6fc655SJed Brown PETSC_EXTERN const char *const MatFactorShiftTypes[];
10415e9742b9SJed Brown PETSC_EXTERN const char *const MatFactorShiftTypesDetail[];
1042d90ac83dSHong Zhang 
1043d90ac83dSHong Zhang /*S
10449aa7eafdSHong Zhang     MatFactorError - indicates what type of error in matrix factor
10459aa7eafdSHong Zhang 
10469aa7eafdSHong Zhang     Level: beginner
10479aa7eafdSHong Zhang 
104800e125f8SBarry Smith     Developer Notes: Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
104900e125f8SBarry Smith 
105000e125f8SBarry Smith .seealso: MatGetFactor()
10519aa7eafdSHong Zhang S*/
10525cd7cf9dSHong Zhang typedef enum {MAT_FACTOR_NOERROR,MAT_FACTOR_STRUCT_ZEROPIVOT,MAT_FACTOR_NUMERIC_ZEROPIVOT,MAT_FACTOR_OUTMEMORY,MAT_FACTOR_OTHER} MatFactorError;
10539aa7eafdSHong Zhang 
105400e125f8SBarry Smith PETSC_EXTERN PetscErrorCode MatFactorGetError(Mat,MatFactorError*);
1055b8b68cfdSBarry Smith PETSC_EXTERN PetscErrorCode MatFactorClearError(Mat);
10567b6c816cSBarry Smith PETSC_EXTERN PetscErrorCode MatFactorGetErrorZeroPivot(Mat,PetscReal*,PetscInt*);
105700e125f8SBarry Smith 
10589aa7eafdSHong Zhang /*S
1059422a814eSBarry Smith    MatFactorInfo - Data passed into the matrix factorization routines, and information about the resulting factorization
1060b00f7748SHong Zhang 
106161cad860SBarry Smith    In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE, that is use
106261cad860SBarry Smith $     MatFactorInfo  info(MAT_FACTORINFO_SIZE)
1063b00f7748SHong Zhang 
106415e8a5b3SHong Zhang    Notes: These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC.
1065b00f7748SHong Zhang 
106661cad860SBarry Smith       You can use MatFactorInfoInitialize() to set default values.
106761cad860SBarry Smith 
1068b00f7748SHong Zhang    Level: developer
1069b00f7748SHong Zhang 
1070d7d82daaSBarry Smith .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(),
1071d7d82daaSBarry Smith           MatFactorInfoInitialize()
1072b00f7748SHong Zhang 
1073b00f7748SHong Zhang S*/
1074b00f7748SHong Zhang typedef struct {
107515e8a5b3SHong Zhang   PetscReal     diagonal_fill;  /* force diagonal to fill in if initially not filled */
107685317021SBarry Smith   PetscReal     usedt;
107715e8a5b3SHong Zhang   PetscReal     dt;             /* drop tolerance */
1078b00f7748SHong Zhang   PetscReal     dtcol;          /* tolerance for pivoting */
107915e8a5b3SHong Zhang   PetscReal     dtcount;        /* maximum nonzeros to be allowed per row */
108067e28bfeSBarry Smith   PetscReal     fill;           /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */
1081348344bbSBarry Smith   PetscReal     levels;         /* ICC/ILU(levels) */
1082bcd9e38bSBarry Smith   PetscReal     pivotinblocks;  /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0
1083bcd9e38bSBarry Smith                                    factorization may be faster if do not pivot */
108415e8a5b3SHong Zhang   PetscReal     zeropivot;      /* pivot is called zero if less than this */
1085f4db908eSBarry Smith   PetscReal     shifttype;      /* type of shift added to matrix factor to prevent zero pivots */
1086f4db908eSBarry Smith   PetscReal     shiftamount;     /* how large the shift is */
108715e8a5b3SHong Zhang } MatFactorInfo;
1088ffa6d0a5SLois Curfman McInnes 
1089014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFactorInfoInitialize(MatFactorInfo*);
10909a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatCholeskyFactor(Mat,IS,const MatFactorInfo*);
10919a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatCholeskyFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
10929a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatCholeskyFactorNumeric(Mat,Mat,const MatFactorInfo*);
10939a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatLUFactor(Mat,IS,IS,const MatFactorInfo*);
10949a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatILUFactor(Mat,IS,IS,const MatFactorInfo*);
10959a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatLUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
10969a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatILUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
10979a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatICCFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
10989a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatICCFactor(Mat,IS,const MatFactorInfo*);
10999a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatLUFactorNumeric(Mat,Mat,const MatFactorInfo*);
1100014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*);
1101014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolve(Mat,Vec,Vec);
1102014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatForwardSolve(Mat,Vec,Vec);
1103014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatBackwardSolve(Mat,Vec,Vec);
1104014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveAdd(Mat,Vec,Vec,Vec);
1105014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveTranspose(Mat,Vec,Vec);
1106014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveTransposeAdd(Mat,Vec,Vec,Vec);
1107014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolves(Mat,Vecs,Vecs);
11088e7ba810SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorSetSchurIS(Mat,IS);
11095a05ddb0SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorGetSchurComplement(Mat,Mat*);
11105a05ddb0SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorRestoreSchurComplement(Mat,Mat*);
11115a05ddb0SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorInvertSchurComplement(Mat);
11125a05ddb0SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorCreateSchurComplement(Mat,Mat*);
11135a05ddb0SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorSolveSchurComplement(Mat,Vec,Vec);
11145a05ddb0SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorSolveSchurComplementTranspose(Mat,Vec,Vec);
11156dba178dSStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorFactorizeSchurComplement(Mat);
1116e8ade678SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorSetSchurComplementSolverType(Mat,PetscInt);
1117014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetUnfactored(Mat);
1118bb5a7306SBarry Smith 
1119d91e6319SBarry Smith /*E
1120d91e6319SBarry Smith     MatSORType - What type of (S)SOR to perform
1121bb1eb677SSatish Balay 
1122d91e6319SBarry Smith     Level: beginner
1123d91e6319SBarry Smith 
1124d9274352SBarry Smith    May be bitwise ORd together
1125d9274352SBarry Smith 
1126af0996ceSBarry Smith    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
1127d91e6319SBarry Smith 
11284e7234bfSBarry Smith    MatSORType may be bitwise ORd together, so do not change the numbers
11294e7234bfSBarry Smith 
113041f059aeSBarry Smith .seealso: MatSOR()
1131d91e6319SBarry Smith E*/
1132ee50ffe9SBarry Smith typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3,
1133ee50ffe9SBarry Smith               SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8,
1134ee50ffe9SBarry Smith               SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16,
113584cb2905SBarry Smith               SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType;
1136014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSOR(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
11378ed539a5SBarry Smith 
1138d4fbbf0eSBarry Smith /*
1139639f9d9dSBarry Smith     These routines are for efficiently computing Jacobians via finite differences.
1140639f9d9dSBarry Smith */
1141b12f92e5SBarry Smith 
11427086a01eSPeter Brune /*S
11437086a01eSPeter Brune      MatColoring - Object for managing the coloring of matrices.
11447086a01eSPeter Brune 
11457086a01eSPeter Brune    Level: beginner
11467086a01eSPeter Brune 
11477086a01eSPeter Brune   Concepts: matrix, coloring
11487086a01eSPeter Brune 
11497086a01eSPeter Brune .seealso:  MatFDColoringCreate() ISColoring MatFDColoring
11507086a01eSPeter Brune S*/
11517086a01eSPeter Brune typedef struct _p_MatColoring* MatColoring;
115276bdecfbSBarry Smith /*J
11538f6c3df8SBarry Smith     MatColoringType - String with the name of a PETSc matrix coloring
1154d9274352SBarry Smith 
1155d9274352SBarry Smith    Level: beginner
1156d9274352SBarry Smith 
11577086a01eSPeter Brune .seealso: MatColoringSetType(), MatColoring
115876bdecfbSBarry Smith J*/
11597086a01eSPeter Brune 
116019fd82e9SBarry Smith typedef const  char*           MatColoringType;
11615567d71aSPeter Brune #define MATCOLORINGJP      "jp"
1162b9acb617SPeter Brune #define MATCOLORINGPOWER   "power"
11632692d6eeSBarry Smith #define MATCOLORINGNATURAL "natural"
11642692d6eeSBarry Smith #define MATCOLORINGSL      "sl"
11652692d6eeSBarry Smith #define MATCOLORINGLF      "lf"
11662692d6eeSBarry Smith #define MATCOLORINGID      "id"
11678121bdceSPeter Brune #define MATCOLORINGGREEDY  "greedy"
1168b12f92e5SBarry Smith 
11698ac6da0aSPeter Brune /*E
11708ac6da0aSPeter Brune    MatColoringWeightType - Type of weight scheme
11718ac6da0aSPeter Brune 
11728ac6da0aSPeter Brune     Not Collective
11738ac6da0aSPeter Brune 
11748ac6da0aSPeter Brune +   MAT_COLORING_RANDOM  - Random weights
11758ac6da0aSPeter Brune .   MAT_COLORING_LEXICAL - Lexical weighting based upon global numbering.
11768ac6da0aSPeter Brune -   MAT_COLORING_LF      - Last-first weighting.
11778ac6da0aSPeter Brune 
11788ac6da0aSPeter Brune     Level: intermediate
11798ac6da0aSPeter Brune 
1180af0996ceSBarry Smith    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
11818ac6da0aSPeter Brune 
11828ac6da0aSPeter Brune .seealso: MatCUSPSetFormat(), MatCUSPFormatOperation
11838ac6da0aSPeter Brune E*/
1184301c140bSPeter Brune typedef enum {MAT_COLORING_WEIGHT_RANDOM,MAT_COLORING_WEIGHT_LEXICAL,MAT_COLORING_WEIGHT_LF,MAT_COLORING_WEIGHT_SL} MatColoringWeightType;
11858ac6da0aSPeter Brune 
1186335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringCreate(Mat,MatColoring*);
1187d7f2a307SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringGetDegrees(Mat,PetscInt,PetscInt*);
1188335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringDestroy(MatColoring*);
1189335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringView(MatColoring,PetscViewer);
1190335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetType(MatColoring,MatColoringType);
1191335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetFromOptions(MatColoring);
1192335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetDistance(MatColoring,PetscInt);
1193335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringGetDistance(MatColoring,PetscInt*);
1194335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetMaxColors(MatColoring,PetscInt);
1195335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringGetMaxColors(MatColoring,PetscInt*);
1196335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringApply(MatColoring,ISColoring*);
1197335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringRegister(const char[],PetscErrorCode(*)(MatColoring));
1198014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*);
11998ac6da0aSPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetWeightType(MatColoring,MatColoringWeightType);
1200c29971fcSPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetWeights(MatColoring,PetscReal*,PetscInt*);
12018ac6da0aSPeter Brune PETSC_EXTERN PetscErrorCode MatColoringCreateWeights(MatColoring,PetscReal **,PetscInt **lperm);
1202cb394dc5SBarry Smith PETSC_EXTERN PetscErrorCode MatColoringTestValid(MatColoring,ISColoring);
1203639f9d9dSBarry Smith 
1204d9274352SBarry Smith /*S
1205d9274352SBarry Smith      MatFDColoring - Object for computing a sparse Jacobian via finite differences
1206d9274352SBarry Smith         and coloring
1207639f9d9dSBarry Smith 
1208d9274352SBarry Smith    Level: beginner
1209d9274352SBarry Smith 
1210d9274352SBarry Smith   Concepts: coloring, sparse Jacobian, finite differences
1211d9274352SBarry Smith 
1212d9274352SBarry Smith .seealso:  MatFDColoringCreate()
1213d9274352SBarry Smith S*/
1214e2a1c21fSSatish Balay typedef struct _p_MatFDColoring* MatFDColoring;
1215639f9d9dSBarry Smith 
1216014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringCreate(Mat,ISColoring,MatFDColoring *);
1217014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringDestroy(MatFDColoring*);
1218014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringView(MatFDColoring,PetscViewer);
1219014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*);
1220014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**);
1221014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal);
1222014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetFromOptions(MatFDColoring);
1223d1e9a80fSBarry Smith PETSC_EXTERN PetscErrorCode MatFDColoringApply(Mat,MatFDColoring,Vec,void *);
1224014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetF(MatFDColoring,Vec);
1225014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,PetscInt*[]);
1226f86b9fbaSHong Zhang PETSC_EXTERN PetscErrorCode MatFDColoringSetUp(Mat,ISColoring,MatFDColoring);
1227f86b9fbaSHong Zhang PETSC_EXTERN PetscErrorCode MatFDColoringSetBlockSize(MatFDColoring,PetscInt,PetscInt);
1228f86b9fbaSHong Zhang 
1229b1683b59SHong Zhang 
1230b1683b59SHong Zhang /*S
1231b9af6bddSHong Zhang      MatTransposeColoring - Object for computing a sparse matrix product C=A*B^T via coloring
1232b1683b59SHong Zhang 
1233b1683b59SHong Zhang    Level: beginner
1234b1683b59SHong Zhang 
1235b1683b59SHong Zhang   Concepts: coloring, sparse matrix product
1236b1683b59SHong Zhang 
1237b9af6bddSHong Zhang .seealso:  MatTransposeColoringCreate()
1238b1683b59SHong Zhang S*/
1239b9af6bddSHong Zhang typedef struct _p_MatTransposeColoring* MatTransposeColoring;
1240b1683b59SHong Zhang 
1241014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeColoringCreate(Mat,ISColoring,MatTransposeColoring *);
1242014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransColoringApplySpToDen(MatTransposeColoring,Mat,Mat);
1243014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransColoringApplyDenToSp(MatTransposeColoring,Mat,Mat);
1244014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeColoringDestroy(MatTransposeColoring*);
1245b1683b59SHong Zhang 
1246639f9d9dSBarry Smith /*
12470752156aSBarry Smith     These routines are for partitioning matrices: currently used only
12483eda8832SBarry Smith   for adjacency matrix, MatCreateMPIAdj().
12490752156aSBarry Smith */
1250ca161407SBarry Smith 
1251d9274352SBarry Smith /*S
1252d9274352SBarry Smith      MatPartitioning - Object for managing the partitioning of a matrix or graph
1253d9274352SBarry Smith 
1254d9274352SBarry Smith    Level: beginner
1255d9274352SBarry Smith 
1256d9274352SBarry Smith   Concepts: partitioning
1257d9274352SBarry Smith 
1258743a1026Svictorle .seealso:  MatPartitioningCreate(), MatPartitioningType
1259d9274352SBarry Smith S*/
126091e9ee9fSBarry Smith typedef struct _p_MatPartitioning* MatPartitioning;
1261d9274352SBarry Smith 
126276bdecfbSBarry Smith /*J
12638f6c3df8SBarry Smith     MatPartitioningType - String with the name of a PETSc matrix partitioning
1264d9274352SBarry Smith 
1265d9274352SBarry Smith    Level: beginner
12660d04baf8SBarry Smith dm
1267b547a13bSvictorle .seealso: MatPartitioningCreate(), MatPartitioning
126876bdecfbSBarry Smith J*/
126919fd82e9SBarry Smith typedef const char* MatPartitioningType;
12702692d6eeSBarry Smith #define MATPARTITIONINGCURRENT  "current"
1271aa1e27eaSFande Kong #define MATPARTITIONINGAVERAGE   "average"
12722692d6eeSBarry Smith #define MATPARTITIONINGSQUARE   "square"
12732692d6eeSBarry Smith #define MATPARTITIONINGPARMETIS "parmetis"
12742692d6eeSBarry Smith #define MATPARTITIONINGCHACO    "chaco"
12752692d6eeSBarry Smith #define MATPARTITIONINGPARTY    "party"
127650d91057SBarry Smith #define MATPARTITIONINGPTSCOTCH "ptscotch"
127788d2ac2bSFande Kong #define MATPARTITIONINGHIERARCH  "hierarch"
127871306c60Sjroman 
1279ca161407SBarry Smith 
1280014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningCreate(MPI_Comm,MatPartitioning*);
128119fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningSetType(MatPartitioning,MatPartitioningType);
1282014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetNParts(MatPartitioning,PetscInt);
1283014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetAdjacency(MatPartitioning,Mat);
1284014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]);
1285014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []);
1286014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningApply(MatPartitioning,IS*);
1287014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningDestroy(MatPartitioning*);
12882aabb6bbSBarry Smith 
1289bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningRegister(const char[],PetscErrorCode (*)(MatPartitioning));
129030de9b25SBarry Smith 
1291f1144a30SSatish Balay 
12922bad1931SBarry Smith 
1293014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningView(MatPartitioning,PetscViewer);
1294014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetFromOptions(MatPartitioning);
129519fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningGetType(MatPartitioning,MatPartitioningType*);
1296ca161407SBarry Smith 
129727538973SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetRepartition(MatPartitioning part);
1298014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetCoarseSequential(MatPartitioning);
1299014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningParmetisGetEdgeCut(MatPartitioning, PetscInt *);
13000752156aSBarry Smith 
1301b6956c12SJose Roman typedef enum { MP_CHACO_MULTILEVEL=1,MP_CHACO_SPECTRAL=2,MP_CHACO_LINEAR=4,MP_CHACO_RANDOM=5,MP_CHACO_SCATTERED=6 } MPChacoGlobalType;
13026a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoGlobalTypes[];
1303b6956c12SJose Roman typedef enum { MP_CHACO_KERNIGHAN=1,MP_CHACO_NONE=2 } MPChacoLocalType;
13046a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoLocalTypes[];
1305b6956c12SJose Roman typedef enum { MP_CHACO_LANCZOS=0,MP_CHACO_RQI=1 } MPChacoEigenType;
13066a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoEigenTypes[];
1307b6956c12SJose Roman 
1308014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetGlobal(MatPartitioning,MPChacoGlobalType);
1309014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetGlobal(MatPartitioning,MPChacoGlobalType*);
1310014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetLocal(MatPartitioning,MPChacoLocalType);
1311014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetLocal(MatPartitioning,MPChacoLocalType*);
1312014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal);
1313014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType);
1314014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenSolver(MatPartitioning,MPChacoEigenType*);
1315014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenTol(MatPartitioning,PetscReal);
1316014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenTol(MatPartitioning,PetscReal*);
1317014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenNumber(MatPartitioning,PetscInt);
1318014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenNumber(MatPartitioning,PetscInt*);
131971306c60Sjroman 
132071306c60Sjroman #define MP_PARTY_OPT "opt"
132171306c60Sjroman #define MP_PARTY_LIN "lin"
132271306c60Sjroman #define MP_PARTY_SCA "sca"
132371306c60Sjroman #define MP_PARTY_RAN "ran"
132471306c60Sjroman #define MP_PARTY_GBF "gbf"
132571306c60Sjroman #define MP_PARTY_GCF "gcf"
132671306c60Sjroman #define MP_PARTY_BUB "bub"
132771306c60Sjroman #define MP_PARTY_DEF "def"
1328014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetGlobal(MatPartitioning,const char*);
132971306c60Sjroman #define MP_PARTY_HELPFUL_SETS "hs"
133071306c60Sjroman #define MP_PARTY_KERNIGHAN_LIN "kl"
133171306c60Sjroman #define MP_PARTY_NONE "no"
1332014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetLocal(MatPartitioning,const char*);
1333014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal);
1334014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetBipart(MatPartitioning,PetscBool);
1335014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscBool);
133671306c60Sjroman 
133750d91057SBarry Smith typedef enum { MP_PTSCOTCH_QUALITY,MP_PTSCOTCH_SPEED,MP_PTSCOTCH_BALANCE,MP_PTSCOTCH_SAFETY,MP_PTSCOTCH_SCALABILITY } MPPTScotchStrategyType;
13386a6fc655SJed Brown PETSC_EXTERN const char *const MPPTScotchStrategyTypes[];
1339e0f1cffaSJose Roman 
1340014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetImbalance(MatPartitioning,PetscReal);
1341014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetImbalance(MatPartitioning,PetscReal*);
1342014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetStrategy(MatPartitioning,MPPTScotchStrategyType);
1343014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetStrategy(MatPartitioning,MPPTScotchStrategyType*);
134471306c60Sjroman 
1345b43b03e9SMark F. Adams /*
13466294fa2bSFande Kong  * hierarchical partitioning
13476294fa2bSFande Kong  */
13484e713f55SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalGetFineparts(MatPartitioning,IS*);
13494e713f55SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalGetCoarseparts(MatPartitioning,IS*);
13504e713f55SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalSetNcoarseparts(MatPartitioning,PetscInt);
13514e713f55SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalSetNfineparts(MatPartitioning, PetscInt);
13526294fa2bSFande Kong 
1353014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMeshToVertexGraph(Mat,PetscInt,Mat*);
1354014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMeshToCellGraph(Mat,PetscInt,Mat*);
1355591294e4SBarry Smith 
13560752156aSBarry Smith /*
1357af0996ceSBarry Smith     If you add entries here you must also add them to petsc/finclude/petscmat.h
1358d4fbbf0eSBarry Smith */
13591c1c02c0SLois Curfman McInnes typedef enum { MATOP_SET_VALUES=0,
13601c1c02c0SLois Curfman McInnes                MATOP_GET_ROW=1,
13611c1c02c0SLois Curfman McInnes                MATOP_RESTORE_ROW=2,
13621c1c02c0SLois Curfman McInnes                MATOP_MULT=3,
13631c1c02c0SLois Curfman McInnes                MATOP_MULT_ADD=4,
13647c922b88SBarry Smith                MATOP_MULT_TRANSPOSE=5,
13657c922b88SBarry Smith                MATOP_MULT_TRANSPOSE_ADD=6,
13661c1c02c0SLois Curfman McInnes                MATOP_SOLVE=7,
13671c1c02c0SLois Curfman McInnes                MATOP_SOLVE_ADD=8,
13687c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE=9,
13697c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE_ADD=10,
13701c1c02c0SLois Curfman McInnes                MATOP_LUFACTOR=11,
13711c1c02c0SLois Curfman McInnes                MATOP_CHOLESKYFACTOR=12,
1372a714c06dSBarry Smith                MATOP_SOR=13,
13731c1c02c0SLois Curfman McInnes                MATOP_TRANSPOSE=14,
13741c1c02c0SLois Curfman McInnes                MATOP_GETINFO=15,
13751c1c02c0SLois Curfman McInnes                MATOP_EQUAL=16,
13761c1c02c0SLois Curfman McInnes                MATOP_GET_DIAGONAL=17,
13771c1c02c0SLois Curfman McInnes                MATOP_DIAGONAL_SCALE=18,
13781c1c02c0SLois Curfman McInnes                MATOP_NORM=19,
13791c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_BEGIN=20,
13801c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_END=21,
1381d519adbfSMatthew Knepley                MATOP_SET_OPTION=22,
1382d519adbfSMatthew Knepley                MATOP_ZERO_ENTRIES=23,
1383d519adbfSMatthew Knepley                MATOP_ZERO_ROWS=24,
1384d519adbfSMatthew Knepley                MATOP_LUFACTOR_SYMBOLIC=25,
1385d519adbfSMatthew Knepley                MATOP_LUFACTOR_NUMERIC=26,
1386d519adbfSMatthew Knepley                MATOP_CHOLESKY_FACTOR_SYMBOLIC=27,
1387d519adbfSMatthew Knepley                MATOP_CHOLESKY_FACTOR_NUMERIC=28,
1388d519adbfSMatthew Knepley                MATOP_SETUP_PREALLOCATION=29,
1389d519adbfSMatthew Knepley                MATOP_ILUFACTOR_SYMBOLIC=30,
1390d519adbfSMatthew Knepley                MATOP_ICCFACTOR_SYMBOLIC=31,
1391a5b7ff6bSBarry Smith                MATOP_GET_DIAGONAL_BLOCK=32,
13929d39f69dSJed Brown                /* MATOP_PLACEHOLDER_33=33, */
1393d519adbfSMatthew Knepley                MATOP_DUPLICATE=34,
1394d519adbfSMatthew Knepley                MATOP_FORWARD_SOLVE=35,
1395d519adbfSMatthew Knepley                MATOP_BACKWARD_SOLVE=36,
1396d519adbfSMatthew Knepley                MATOP_ILUFACTOR=37,
1397d519adbfSMatthew Knepley                MATOP_ICCFACTOR=38,
1398d519adbfSMatthew Knepley                MATOP_AXPY=39,
1399d519adbfSMatthew Knepley                MATOP_GET_SUBMATRICES=40,
1400d519adbfSMatthew Knepley                MATOP_INCREASE_OVERLAP=41,
1401d519adbfSMatthew Knepley                MATOP_GET_VALUES=42,
1402d519adbfSMatthew Knepley                MATOP_COPY=43,
1403d519adbfSMatthew Knepley                MATOP_GET_ROW_MAX=44,
1404d519adbfSMatthew Knepley                MATOP_SCALE=45,
1405d519adbfSMatthew Knepley                MATOP_SHIFT=46,
140635153367SBarry Smith                MATOP_DIAGONAL_SET=47,
14079d39f69dSJed Brown                MATOP_ZERO_ROWS_COLUMNS=48,
14089d39f69dSJed Brown                MATOP_SET_RANDOM=49,
1409d519adbfSMatthew Knepley                MATOP_GET_ROW_IJ=50,
1410d519adbfSMatthew Knepley                MATOP_RESTORE_ROW_IJ=51,
1411d519adbfSMatthew Knepley                MATOP_GET_COLUMN_IJ=52,
1412d519adbfSMatthew Knepley                MATOP_RESTORE_COLUMN_IJ=53,
1413d519adbfSMatthew Knepley                MATOP_FDCOLORING_CREATE=54,
1414d519adbfSMatthew Knepley                MATOP_COLORING_PATCH=55,
1415d519adbfSMatthew Knepley                MATOP_SET_UNFACTORED=56,
1416d519adbfSMatthew Knepley                MATOP_PERMUTE=57,
1417d519adbfSMatthew Knepley                MATOP_SET_VALUES_BLOCKED=58,
1418d519adbfSMatthew Knepley                MATOP_GET_SUBMATRIX=59,
1419d519adbfSMatthew Knepley                MATOP_DESTROY=60,
1420d519adbfSMatthew Knepley                MATOP_VIEW=61,
1421d519adbfSMatthew Knepley                MATOP_CONVERT_FROM=62,
14227bab7c10SHong Zhang                MATOP_MATMAT_MULT=63,
14237bab7c10SHong Zhang                MATOP_MATMAT_MULT_SYMBOLIC=64,
14247bab7c10SHong Zhang                MATOP_MATMAT_MULT_NUMERIC=65,
1425d519adbfSMatthew Knepley                MATOP_SET_LOCAL_TO_GLOBAL_MAP=66,
1426d519adbfSMatthew Knepley                MATOP_SET_VALUES_LOCAL=67,
1427d519adbfSMatthew Knepley                MATOP_ZERO_ROWS_LOCAL=68,
1428d519adbfSMatthew Knepley                MATOP_GET_ROW_MAX_ABS=69,
1429d519adbfSMatthew Knepley                MATOP_GET_ROW_MIN_ABS=70,
1430d519adbfSMatthew Knepley                MATOP_CONVERT=71,
1431d519adbfSMatthew Knepley                MATOP_SET_COLORING=72,
14329d39f69dSJed Brown                /* MATOP_PLACEHOLDER_73=73, */
1433d519adbfSMatthew Knepley                MATOP_SET_VALUES_ADIFOR=74,
1434d519adbfSMatthew Knepley                MATOP_FD_COLORING_APPLY=75,
1435d519adbfSMatthew Knepley                MATOP_SET_FROM_OPTIONS=76,
1436bda6c4cbSBarry Smith                MATOP_MULT_CONSTRAINED=77,
1437bda6c4cbSBarry Smith                MATOP_MULT_TRANSPOSE_CONSTRAIN=78,
14389d39f69dSJed Brown                MATOP_FIND_ZERO_DIAGONALS=79,
1439d519adbfSMatthew Knepley                MATOP_MULT_MULTIPLE=80,
1440d519adbfSMatthew Knepley                MATOP_SOLVE_MULTIPLE=81,
1441d519adbfSMatthew Knepley                MATOP_GET_INERTIA=82,
1442d519adbfSMatthew Knepley                MATOP_LOAD=83,
1443d519adbfSMatthew Knepley                MATOP_IS_SYMMETRIC=84,
1444d519adbfSMatthew Knepley                MATOP_IS_HERMITIAN=85,
1445d519adbfSMatthew Knepley                MATOP_IS_STRUCTURALLY_SYMMETRIC=86,
1446820469bcSHong Zhang                MATOP_SET_VALUES_BLOCKEDLOCAL=87,
1447d519adbfSMatthew Knepley                MATOP_GET_VECS=88,
1448d519adbfSMatthew Knepley                MATOP_MAT_MULT=89,
1449d519adbfSMatthew Knepley                MATOP_MAT_MULT_SYMBOLIC=90,
1450d519adbfSMatthew Knepley                MATOP_MAT_MULT_NUMERIC=91,
1451d519adbfSMatthew Knepley                MATOP_PTAP=92,
1452d519adbfSMatthew Knepley                MATOP_PTAP_SYMBOLIC=93,
1453d519adbfSMatthew Knepley                MATOP_PTAP_NUMERIC=94,
1454bda6c4cbSBarry Smith                MATOP_MAT_TRANSPOSE_MULT=95,
1455bda6c4cbSBarry Smith                MATOP_MAT_TRANSPOSE_MULT_SYMBO=96,
1456bda6c4cbSBarry Smith                MATOP_MAT_TRANSPOSE_MULT_NUMER=97,
14579d39f69dSJed Brown                /* MATOP_PLACEHOLDER_98=98, */
14589d39f69dSJed Brown                /* MATOP_PLACEHOLDER_99=99, */
14599d39f69dSJed Brown                /* MATOP_PLACEHOLDER_100=100, */
14609d39f69dSJed Brown                /* MATOP_PLACEHOLDER_101=101, */
1461d519adbfSMatthew Knepley                MATOP_CONJUGATE=102,
14629d39f69dSJed Brown                /* MATOP_PLACEHOLDER_103=103, */
1463d519adbfSMatthew Knepley                MATOP_SET_VALUES_ROW=104,
1464d519adbfSMatthew Knepley                MATOP_REAL_PART=105,
1465bda6c4cbSBarry Smith                MATOP_IMAGINARY_PART=106,
1466bda6c4cbSBarry Smith                MATOP_GET_ROW_UPPER_TRIANGULAR=107,
1467bda6c4cbSBarry Smith                MATOP_RESTORE_ROW_UPPER_TRIANG=108,
1468bda6c4cbSBarry Smith                MATOP_MAT_SOLVE=109,
14690e8d04f7SBarry Smith                MATOP_GET_REDUNDANT_MATRIX=110,
1470d519adbfSMatthew Knepley                MATOP_GET_ROW_MIN=111,
14710e8d04f7SBarry Smith                MATOP_GET_COLUMN_VECTOR=112,
1472d519adbfSMatthew Knepley                MATOP_MISSING_DIAGONAL=113,
14730e8d04f7SBarry Smith                MATOP_GET_SEQ_NONZERO_STRUCTUR=114,
147489c6957cSBarry Smith                MATOP_CREATE=115,
147589c6957cSBarry Smith                MATOP_GET_GHOSTS=116,
14760e8d04f7SBarry Smith                MATOP_GET_LOCAL_SUB_MATRIX=117,
14770e8d04f7SBarry Smith                MATOP_RESTORE_LOCALSUB_MATRIX=118,
1478eeffb40dSHong Zhang                MATOP_MULT_DIAGONAL_BLOCK=119,
14790e8d04f7SBarry Smith                MATOP_HERMITIAN_TRANSPOSE=120,
14800e8d04f7SBarry Smith                MATOP_MULT_HERMITIAN_TRANSPOSE=121,
14810e8d04f7SBarry Smith                MATOP_MULT_HERMITIAN_TRANS_ADD=122,
14820e8d04f7SBarry Smith                MATOP_GET_MULTI_PROC_BLOCK=123,
14839d39f69dSJed Brown                MATOP_FIND_NONZERO_ROWS=124,
14840e8d04f7SBarry Smith                MATOP_GET_COLUMN_NORMS=125,
14859d39f69dSJed Brown                MATOP_INVERT_BLOCK_DIAGONAL=126,
14869d39f69dSJed Brown                /* MATOP_PLACEHOLDER_127=127, */
14870e8d04f7SBarry Smith                MATOP_GET_SUB_MATRICES_PARALLE=128,
14882b8ad9a3SHong Zhang                MATOP_SET_VALUES_BATCH=129,
14890e8d04f7SBarry Smith                MATOP_TRANSPOSE_MAT_MULT=130,
14900e8d04f7SBarry Smith                MATOP_TRANSPOSE_MAT_MULT_SYMBO=131,
1491e9b602ebSSatish Balay                MATOP_TRANSPOSE_MAT_MULT_NUMER=132,
14920e8d04f7SBarry Smith                MATOP_TRANSPOSE_COLORING_CREAT=133,
14930e8d04f7SBarry Smith                MATOP_TRANS_COLORING_APPLY_SPT=134,
14940e8d04f7SBarry Smith                MATOP_TRANS_COLORING_APPLY_DEN=135,
14950e8d04f7SBarry Smith                MATOP_RART=136,
14960e8d04f7SBarry Smith                MATOP_RART_SYMBOLIC=137,
14970e8d04f7SBarry Smith                MATOP_RART_NUMERIC=138,
1498e09a3074SHong Zhang                MATOP_SET_BLOCK_SIZES=139,
1499f9426fe0SMark Adams                MATOP_AYPX=140,
15001919a2e2SJed Brown                MATOP_RESIDUAL=141,
15019c8f2541SHong Zhang                MATOP_FDCOLORING_SETUP=142,
15029c8f2541SHong Zhang                MATOP_MPICONCATENATESEQ=144
1503fae171e0SBarry Smith              } MatOperation;
1504014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatHasOperation(Mat,MatOperation,PetscBool *);
1505014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellSetOperation(Mat,MatOperation,void(*)(void));
1506014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellGetOperation(Mat,MatOperation,void(**)(void));
1507014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellSetContext(Mat,void*);
1508112a2221SBarry Smith 
150990ace30eSBarry Smith /*
151090ace30eSBarry Smith    Codes for matrices stored on disk. By default they are
151190ace30eSBarry Smith    stored in a universal format. By changing the format with
15126a9046bcSBarry Smith    PetscViewerPushFormat(viewer,PETSC_VIEWER_NATIVE); the matrices will
151390ace30eSBarry Smith    be stored in a way natural for the matrix, for example dense matrices
151490ace30eSBarry Smith    would be stored as dense. Matrices stored this way may only be
1515f4403165SShri Abhyankar    read into matrices of the same type.
151690ace30eSBarry Smith */
151790ace30eSBarry Smith #define MATRIX_BINARY_FORMAT_DENSE -1
151890ace30eSBarry Smith 
1519014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetHashTableFactor(Mat,PetscReal);
1520014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatISGetLocalMat(Mat,Mat*);
1521014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatISSetLocalMat(Mat,Mat);
1522b7ce53b6SStefano Zampini PETSC_EXTERN PetscErrorCode MatISGetMPIXAIJ(Mat,MatReuse,Mat*);
15231f4e1ec7SBarry Smith 
1524d9274352SBarry Smith /*S
1525d9274352SBarry Smith      MatNullSpace - Object that removes a null space from a vector, i.e.
1526d9274352SBarry Smith          orthogonalizes the vector to a subsapce
1527d9274352SBarry Smith 
1528f7a9e4ceSBarry Smith    Level: advanced
1529d9274352SBarry Smith 
1530d9274352SBarry Smith   Concepts: matrix; linear operator, null space
1531d9274352SBarry Smith 
15326e1639daSBarry Smith   Users manual sections:
15336e1639daSBarry Smith .   sec_singular
15346e1639daSBarry Smith 
1535d9274352SBarry Smith .seealso:  MatNullSpaceCreate()
1536d9274352SBarry Smith S*/
153774637425SBarry Smith typedef struct _p_MatNullSpace* MatNullSpace;
1538d9274352SBarry Smith 
1539014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceCreate(MPI_Comm,PetscBool ,PetscInt,const Vec[],MatNullSpace*);
1540014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(MatNullSpace,Vec,void*),void*);
1541014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceDestroy(MatNullSpace*);
1542d0195637SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceRemove(MatNullSpace,Vec);
1543014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNullSpace(Mat, MatNullSpace *);
15445fa7ec2dSBarry Smith PETSC_EXTERN PetscErrorCode MatGetTransposeNullSpace(Mat, MatNullSpace *);
15455fa7ec2dSBarry Smith PETSC_EXTERN PetscErrorCode MatSetTransposeNullSpace(Mat,MatNullSpace);
1546014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNullSpace(Mat,MatNullSpace);
1547014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNearNullSpace(Mat,MatNullSpace);
1548014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNearNullSpace(Mat,MatNullSpace*);
1549014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceTest(MatNullSpace,Mat,PetscBool  *);
1550014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceView(MatNullSpace,PetscViewer);
1551014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceGetVecs(MatNullSpace,PetscBool*,PetscInt*,const Vec**);
1552014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceCreateRigidBody(Vec,MatNullSpace*);
155374637425SBarry Smith 
1554014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatReorderingSeqSBAIJ(Mat,IS);
1555014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetHashTableFactor(Mat,PetscReal);
1556014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetColumnIndices(Mat,PetscInt *);
1557014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJInvertBlockDiagonal(Mat);
15583f1d51d7SBarry Smith 
1559014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMAIJ(Mat,PetscInt,Mat*);
1560014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMAIJRedimension(Mat,PetscInt,Mat*);
1561014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMAIJGetAIJ(Mat,Mat*);
1562c4f061fbSSatish Balay 
1563014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatComputeExplicitOperator(Mat,Mat*);
1564b0a32e0cSBarry Smith 
1565014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalScaleLocal(Mat,Vec);
156604f1ad80SBarry Smith 
1567014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*);
1568014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetBase(Mat,Vec,Vec);
1569014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*);
1570014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*));
1571014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec));
1572014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt);
1573014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDResetHHistory(Mat);
1574014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctionError(Mat,PetscReal);
1575014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetPeriod(Mat,PetscInt);
1576014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDGetH(Mat,PetscScalar *);
1577014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetOptionsPrefix(Mat,const char[]);
1578014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*);
1579014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*);
1580e884886fSBarry Smith 
15816370053bSBarry Smith /*S
15826370053bSBarry Smith     MatMFFD - A data structured used to manage the computation of the h differencing parameter for matrix-free
15836370053bSBarry Smith               Jacobian vector products
1584e884886fSBarry Smith 
15856370053bSBarry Smith     Notes: MATMFFD is a specific MatType which uses the MatMFFD data structure
15866370053bSBarry Smith 
15876370053bSBarry Smith            MatMFFD*() methods actually take the Mat as their first argument. Not a MatMFFD data structure
15886370053bSBarry Smith 
15896370053bSBarry Smith     Level: developer
15906370053bSBarry Smith 
15916370053bSBarry Smith .seealso: MATMFFD, MatCreateMFFD(), MatMFFDSetFuction(), MatMFFDSetType(), MatMFFDRegister()
15926370053bSBarry Smith S*/
1593e884886fSBarry Smith typedef struct _p_MatMFFD* MatMFFD;
1594e884886fSBarry Smith 
159576bdecfbSBarry Smith /*J
1596e884886fSBarry Smith     MatMFFDType - algorithm used to compute the h used in computing matrix-vector products via differencing of the function
1597e884886fSBarry Smith 
1598e884886fSBarry Smith    Level: beginner
1599e884886fSBarry Smith 
1600e884886fSBarry Smith .seealso: MatMFFDSetType(), MatMFFDRegister()
160176bdecfbSBarry Smith J*/
160219fd82e9SBarry Smith typedef const char* MatMFFDType;
1603e884886fSBarry Smith #define MATMFFD_DS  "ds"
1604e884886fSBarry Smith #define MATMFFD_WP  "wp"
1605e884886fSBarry Smith 
160619fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDSetType(Mat,MatMFFDType);
1607bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDRegister(const char[],PetscErrorCode (*)(MatMFFD));
1608e884886fSBarry Smith 
1609014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDDSSetUmin(Mat,PetscReal);
1610014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDWPSetComputeNormU(Mat,PetscBool );
1611e884886fSBarry Smith 
161261ab5df0SBarry Smith PETSC_EXTERN PetscErrorCode MatFDColoringSetType(MatFDColoring,MatMFFDType);
161361ab5df0SBarry Smith 
1614014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *);
1615014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *);
16167dbadf16SMatthew Knepley 
161797969023SHong Zhang /*
161897969023SHong Zhang    PETSc interface to MUMPS
161997969023SHong Zhang */
162097969023SHong Zhang #ifdef PETSC_HAVE_MUMPS
1621014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMumpsSetIcntl(Mat,PetscInt,PetscInt);
1622bc6112feSHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetIcntl(Mat,PetscInt,PetscInt*);
1623b9e9ea26SSatish Balay PETSC_EXTERN PetscErrorCode MatMumpsSetCntl(Mat,PetscInt,PetscReal);
1624bc6112feSHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetCntl(Mat,PetscInt,PetscReal*);
1625bc6112feSHong Zhang 
1626ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetInfo(Mat,PetscInt,PetscInt*);
1627ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetInfog(Mat,PetscInt,PetscInt*);
1628ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetRinfo(Mat,PetscInt,PetscReal*);
1629ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetRinfog(Mat,PetscInt,PetscReal*);
163097969023SHong Zhang #endif
163123a5497aSJed Brown 
1632d954db57SHong Zhang /*
1633d615f992SSatish Balay    PETSc interface to Mkl_Pardiso
1634b500a6b7SJose David Bermeol */
1635b500a6b7SJose David Bermeol #ifdef PETSC_HAVE_MKL_PARDISO
1636d615f992SSatish Balay PETSC_EXTERN PetscErrorCode MatMkl_PardisoSetCntl(Mat,PetscInt,PetscInt);
1637b500a6b7SJose David Bermeol #endif
1638b500a6b7SJose David Bermeol 
1639b500a6b7SJose David Bermeol /*
1640d305a81bSVasiliy Kozyrev    PETSc interface to Mkl_CPardiso
1641d305a81bSVasiliy Kozyrev */
1642d305a81bSVasiliy Kozyrev #ifdef PETSC_HAVE_MKL_CPARDISO
1643d305a81bSVasiliy Kozyrev PETSC_EXTERN PetscErrorCode MatMkl_CPardisoSetCntl(Mat,PetscInt,PetscInt);
1644d305a81bSVasiliy Kozyrev #endif
1645d305a81bSVasiliy Kozyrev 
1646d305a81bSVasiliy Kozyrev /*
1647d954db57SHong Zhang    PETSc interface to SUPERLU
1648d954db57SHong Zhang */
1649d954db57SHong Zhang #ifdef PETSC_HAVE_SUPERLU
1650014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSuperluSetILUDropTol(Mat,PetscReal);
1651d954db57SHong Zhang #endif
1652d954db57SHong Zhang 
1653fb785984SHong Zhang /*
1654fb785984SHong Zhang    PETSc interface to SUPERLU_DIST
1655fb785984SHong Zhang */
1656fb785984SHong Zhang #ifdef PETSC_HAVE_SUPERLU_DIST
1657fb785984SHong Zhang PETSC_EXTERN PetscErrorCode MatSuperluDistGetDiagU(Mat,PetscScalar*);
1658fb785984SHong Zhang #endif
1659fb785984SHong Zhang 
1660575704cbSPieter Ghysels /*
1661575704cbSPieter Ghysels    PETSc interface to STRUMPACK
1662575704cbSPieter Ghysels */
1663575704cbSPieter Ghysels #ifdef PETSC_HAVE_STRUMPACK
1664575704cbSPieter Ghysels PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSRelCompTol(Mat,PetscReal);
1665575704cbSPieter Ghysels PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSMinSize(Mat,PetscInt);
1666575704cbSPieter Ghysels PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetColPerm(Mat,PetscBool);
1667575704cbSPieter Ghysels #endif
1668575704cbSPieter Ghysels 
1669575704cbSPieter Ghysels 
1670aa372e3fSPaul Mullowney #ifdef PETSC_HAVE_CUDA
16713f9c0db1SPaul Mullowney /*E
1672e057df02SPaul Mullowney     MatCUSPARSEStorageFormat - indicates the storage format for CUSPARSE (GPU)
16732692e278SPaul Mullowney     matrices.
1674e057df02SPaul Mullowney 
1675e057df02SPaul Mullowney     Not Collective
1676e057df02SPaul Mullowney 
16778468deeeSKarl Rupp +   MAT_CUSPARSE_CSR - Compressed Sparse Row
16782692e278SPaul Mullowney .   MAT_CUSPARSE_ELL - Ellpack (requires CUDA 4.2 or later).
16792692e278SPaul Mullowney -   MAT_CUSPARSE_HYB - Hybrid, a combination of Ellpack and Coordinate format (requires CUDA 4.2 or later).
1680e057df02SPaul Mullowney 
1681e057df02SPaul Mullowney     Level: intermediate
1682e057df02SPaul Mullowney 
1683af0996ceSBarry Smith    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
1684e057df02SPaul Mullowney 
1685e057df02SPaul Mullowney .seealso: MatCUSPARSESetFormat(), MatCUSPARSEFormatOperation
1686e057df02SPaul Mullowney E*/
1687e057df02SPaul Mullowney 
16883f9c0db1SPaul Mullowney typedef enum {MAT_CUSPARSE_CSR, MAT_CUSPARSE_ELL, MAT_CUSPARSE_HYB} MatCUSPARSEStorageFormat;
1689e057df02SPaul Mullowney 
1690e057df02SPaul Mullowney /* these will be strings associated with enumerated type defined above */
1691e057df02SPaul Mullowney PETSC_EXTERN const char *const MatCUSPARSEStorageFormats[];
1692e057df02SPaul Mullowney 
16933f9c0db1SPaul Mullowney /*E
1694e057df02SPaul Mullowney     MatCUSPARSEFormatOperation - indicates the operation of CUSPARSE (GPU)
16952692e278SPaul Mullowney     matrices whose operation should use a particular storage format.
1696e057df02SPaul Mullowney 
1697e057df02SPaul Mullowney     Not Collective
1698e057df02SPaul Mullowney 
16998468deeeSKarl Rupp +   MAT_CUSPARSE_MULT_DIAG - sets the storage format for the diagonal matrix in the parallel MatMult
17008468deeeSKarl Rupp .   MAT_CUSPARSE_MULT_OFFDIAG - sets the storage format for the offdiagonal matrix in the parallel MatMult
17018468deeeSKarl Rupp .   MAT_CUSPARSE_MULT - sets the storage format for the entire matrix in the serial (single GPU) MatMult
17028468deeeSKarl Rupp -   MAT_CUSPARSE_ALL - sets the storage format for all CUSPARSE (GPU) matrices
1703e057df02SPaul Mullowney 
1704e057df02SPaul Mullowney     Level: intermediate
1705e057df02SPaul Mullowney 
1706e057df02SPaul Mullowney .seealso: MatCUSPARSESetFormat(), MatCUSPARSEStorageFormat
1707e057df02SPaul Mullowney E*/
170836d62e41SPaul Mullowney typedef enum {MAT_CUSPARSE_MULT_DIAG, MAT_CUSPARSE_MULT_OFFDIAG, MAT_CUSPARSE_MULT, MAT_CUSPARSE_ALL} MatCUSPARSEFormatOperation;
1709e057df02SPaul Mullowney 
171010e9db80SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
171110e9db80SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCreateAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
1712e057df02SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCUSPARSESetFormat(Mat,MatCUSPARSEFormatOperation,MatCUSPARSEStorageFormat);
17139ae82921SPaul Mullowney #endif
17149ae82921SPaul Mullowney 
171590c53902SBarry Smith #if defined(PETSC_HAVE_CUSP)
1716014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCUSP(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
1717014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateAIJCUSP(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
1718e057df02SPaul Mullowney 
17193f9c0db1SPaul Mullowney /*E
1720e057df02SPaul Mullowney     MatCUSPStorageFormat - indicates the storage format for CUSP (GPU)
172136d62e41SPaul Mullowney     matrices.
1722e057df02SPaul Mullowney 
1723e057df02SPaul Mullowney     Not Collective
1724e057df02SPaul Mullowney 
17258468deeeSKarl Rupp +   MAT_CUSP_CSR - Compressed Sparse Row
17268468deeeSKarl Rupp .   MAT_CUSP_DIA - Diagonal
17278468deeeSKarl Rupp -   MAT_CUSP_ELL - Ellpack
1728e057df02SPaul Mullowney 
1729e057df02SPaul Mullowney     Level: intermediate
1730e057df02SPaul Mullowney 
1731af0996ceSBarry Smith    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
1732e057df02SPaul Mullowney 
1733e057df02SPaul Mullowney .seealso: MatCUSPSetFormat(), MatCUSPFormatOperation
1734e057df02SPaul Mullowney E*/
17353f9c0db1SPaul Mullowney typedef enum {MAT_CUSP_CSR, MAT_CUSP_DIA, MAT_CUSP_ELL} MatCUSPStorageFormat;
1736e057df02SPaul Mullowney 
1737e057df02SPaul Mullowney /* these will be strings associated with enumerated type defined above */
1738e057df02SPaul Mullowney PETSC_EXTERN const char *const MatCUSPStorageFormats[];
1739e057df02SPaul Mullowney 
17403f9c0db1SPaul Mullowney /*E
1741e057df02SPaul Mullowney     MatCUSPFormatOperation - indicates the operation of CUSP (GPU)
174236d62e41SPaul Mullowney     matrices whose operation should use a particular storage format.
1743e057df02SPaul Mullowney 
1744e057df02SPaul Mullowney     Not Collective
1745e057df02SPaul Mullowney 
17468468deeeSKarl Rupp +   MAT_CUSP_MULT_DIAG - sets the storage format for the diagonal matrix in the parallel MatMult
17478468deeeSKarl Rupp .   MAT_CUSP_MULT_OFFDIAG - sets the storage format for the offdiagonal matrix in the parallel MatMult
17488468deeeSKarl Rupp .   MAT_CUSP_MULT - sets the storage format for the entire matrix in the serial (single GPU) MatMult
17498468deeeSKarl Rupp -   MAT_CUSP_ALL - sets the storage format for all CUSP (GPU) matrices
1750e057df02SPaul Mullowney 
1751e057df02SPaul Mullowney     Level: intermediate
1752e057df02SPaul Mullowney 
1753af0996ceSBarry Smith    Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h
1754e057df02SPaul Mullowney 
1755e057df02SPaul Mullowney .seealso: MatCUSPSetFormat(), MatCUSPStorageFormat
1756e057df02SPaul Mullowney E*/
175736d62e41SPaul Mullowney typedef enum {MAT_CUSP_MULT_DIAG, MAT_CUSP_MULT_OFFDIAG, MAT_CUSP_MULT, MAT_CUSP_ALL} MatCUSPFormatOperation;
1758e057df02SPaul Mullowney 
1759e057df02SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCUSPSetFormat(Mat,MatCUSPFormatOperation,MatCUSPStorageFormat);
1760e057df02SPaul Mullowney #endif
176190c53902SBarry Smith 
1762d67ff14aSKarl Rupp #if defined(PETSC_HAVE_VIENNACL)
1763d67ff14aSKarl Rupp PETSC_EXTERN PetscErrorCode MatCreateSeqAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
1764d67ff14aSKarl Rupp PETSC_EXTERN PetscErrorCode MatCreateAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
1765d67ff14aSKarl Rupp #endif
1766d67ff14aSKarl Rupp 
176754efbe56SHong Zhang /*
176854efbe56SHong Zhang    PETSc interface to FFTW
176954efbe56SHong Zhang */
177054efbe56SHong Zhang #if defined(PETSC_HAVE_FFTW)
1771014dd563SJed Brown PETSC_EXTERN PetscErrorCode VecScatterPetscToFFTW(Mat,Vec,Vec);
1772014dd563SJed Brown PETSC_EXTERN PetscErrorCode VecScatterFFTWToPetsc(Mat,Vec,Vec);
17732a7a6963SBarry Smith PETSC_EXTERN PetscErrorCode MatCreateVecsFFTW(Mat,Vec*,Vec*,Vec*);
177454efbe56SHong Zhang #endif
177554efbe56SHong Zhang 
1776014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateNest(MPI_Comm,PetscInt,const IS[],PetscInt,const IS[],const Mat[],Mat*);
1777014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSize(Mat,PetscInt*,PetscInt*);
1778014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetISs(Mat,IS[],IS[]);
1779014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetLocalISs(Mat,IS[],IS[]);
1780014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSubMats(Mat,PetscInt*,PetscInt*,Mat***);
1781014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSubMat(Mat,PetscInt,PetscInt,Mat*);
178219fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatNestSetVecType(Mat,VecType);
1783014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestSetSubMats(Mat,PetscInt,const IS[],PetscInt,const IS[],const Mat[]);
1784014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestSetSubMat(Mat,PetscInt,PetscInt,Mat);
1785d8588912SDave May 
17864325cce7SMatthew G Knepley PETSC_EXTERN PetscErrorCode MatChop(Mat,PetscReal);
1787e0a58c10SMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatComputeBandwidth(Mat,PetscReal,PetscInt*);
17884325cce7SMatthew G Knepley 
1789b541e6a4SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSubdomainsCreateCoalesce(Mat,PetscInt,PetscInt*,IS**);
179053dd7562SDmitry Karpeev 
1791c094ef40SMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatPreallocatorPreallocate(Mat,PetscBool,Mat);
1792c094ef40SMatthew G. Knepley 
1793539c167fSBarry Smith PETSC_INTERN PetscErrorCode MatHeaderMerge(Mat,Mat*);
1794539c167fSBarry Smith PETSC_EXTERN PetscErrorCode MatHeaderReplace(Mat,Mat*);
1795539c167fSBarry Smith 
179623a5497aSJed Brown #endif
1797