xref: /petsc/include/petscmat.h (revision d615f992f49c992098f6d6dcf1af9a293ec6e320)
12eac72dbSBarry Smith /*
22eac72dbSBarry Smith      Include file for the matrix component of PETSc
32eac72dbSBarry Smith */
40a835dfdSSatish Balay #ifndef __PETSCMAT_H
50a835dfdSSatish Balay #define __PETSCMAT_H
62c8e378dSBarry Smith #include <petscvec.h>
72eac72dbSBarry Smith 
8d9274352SBarry Smith /*S
98f6c3df8SBarry Smith      Mat - Abstract PETSc matrix object used to manage all linear operators in PETSc, even those without
108f6c3df8SBarry Smith            an explicit sparse representation (such as matrix-free operators)
112eac72dbSBarry Smith 
12d91e6319SBarry Smith    Level: beginner
13d91e6319SBarry Smith 
14d9274352SBarry Smith   Concepts: matrix; linear operator
15d9274352SBarry Smith 
168f6c3df8SBarry Smith .seealso:  MatCreate(), MatType, MatSetType(), MatDestroy()
17d9274352SBarry Smith S*/
18d9274352SBarry Smith typedef struct _p_Mat*           Mat;
19d9274352SBarry Smith 
2076bdecfbSBarry Smith /*J
218f6c3df8SBarry Smith     MatType - String with the name of a PETSc matrix type
22d9274352SBarry Smith 
23d9274352SBarry Smith    Level: beginner
24d9274352SBarry Smith 
258f6c3df8SBarry Smith .seealso: MatSetType(), Mat, MatSolverPackage, MatRegister()
2676bdecfbSBarry Smith J*/
2719fd82e9SBarry Smith typedef const char* MatType;
28273d9f13SBarry Smith #define MATSAME            "same"
295a11e1b2SBarry Smith #define MATMAIJ            "maij"
30273d9f13SBarry Smith #define MATSEQMAIJ         "seqmaij"
31273d9f13SBarry Smith #define MATMPIMAIJ         "mpimaij"
32273d9f13SBarry Smith #define MATIS              "is"
335a11e1b2SBarry Smith #define MATAIJ             "aij"
34273d9f13SBarry Smith #define MATSEQAIJ          "seqaij"
357d6a0e61SBarry Smith #define MATSEQAIJPTHREAD   "seqaijpthread"
36bf2c1783SBarry Smith #define MATAIJPTHREAD      "aijpthread"
37273d9f13SBarry Smith #define MATMPIAIJ          "mpiaij"
385a11e1b2SBarry Smith #define MATAIJCRL          "aijcrl"
395a11e1b2SBarry Smith #define MATSEQAIJCRL       "seqaijcrl"
405a11e1b2SBarry Smith #define MATMPIAIJCRL       "mpiaijcrl"
418154be41SBarry Smith #define MATAIJCUSP         "aijcusp"
428154be41SBarry Smith #define MATSEQAIJCUSP      "seqaijcusp"
438154be41SBarry Smith #define MATMPIAIJCUSP      "mpiaijcusp"
449ae82921SPaul Mullowney #define MATAIJCUSPARSE     "aijcusparse"
459ae82921SPaul Mullowney #define MATSEQAIJCUSPARSE  "seqaijcusparse"
469ae82921SPaul Mullowney #define MATMPIAIJCUSPARSE  "mpiaijcusparse"
47d67ff14aSKarl Rupp #define MATAIJVIENNACL     "aijviennacl"
48d67ff14aSKarl Rupp #define MATSEQAIJVIENNACL  "seqaijviennacl"
49d67ff14aSKarl Rupp #define MATMPIAIJVIENNACL  "mpiaijviennacl"
505a11e1b2SBarry Smith #define MATAIJPERM         "aijperm"
515a11e1b2SBarry Smith #define MATSEQAIJPERM      "seqaijperm"
525a11e1b2SBarry Smith #define MATMPIAIJPERM      "mpiaijperm"
53273d9f13SBarry Smith #define MATSHELL           "shell"
545a11e1b2SBarry Smith #define MATDENSE           "dense"
55209238afSKris Buschelman #define MATSEQDENSE        "seqdense"
56273d9f13SBarry Smith #define MATMPIDENSE        "mpidense"
57db31f6deSJed Brown #define MATELEMENTAL       "elemental"
585a11e1b2SBarry Smith #define MATBAIJ            "baij"
59273d9f13SBarry Smith #define MATSEQBAIJ         "seqbaij"
60273d9f13SBarry Smith #define MATMPIBAIJ         "mpibaij"
61273d9f13SBarry Smith #define MATMPIADJ          "mpiadj"
625a11e1b2SBarry Smith #define MATSBAIJ           "sbaij"
63273d9f13SBarry Smith #define MATSEQSBAIJ        "seqsbaij"
64273d9f13SBarry Smith #define MATMPISBAIJ        "mpisbaij"
65c0cdd4a1SDahai Guo #define MATSEQBSTRM        "seqbstrm"
66c0cdd4a1SDahai Guo #define MATMPIBSTRM        "mpibstrm"
67c0cdd4a1SDahai Guo #define MATBSTRM           "bstrm"
68aa5a9175SDahai Guo #define MATSEQSBSTRM       "seqsbstrm"
69aa5a9175SDahai Guo #define MATMPISBSTRM       "mpisbstrm"
70aa5a9175SDahai Guo #define MATSBSTRM          "sbstrm"
71cebc7f6cSBarry Smith #define MATDAAD            "daad"
72cebc7f6cSBarry Smith #define MATMFFD            "mffd"
73c8a8475eSBarry Smith #define MATNORMAL          "normal"
74ab92ecdeSBarry Smith #define MATLRC             "lrc"
752a6744ebSBarry Smith #define MATSCATTER         "scatter"
76421e10b8SBarry Smith #define MATBLOCKMAT        "blockmat"
77793850ffSBarry Smith #define MATCOMPOSITE       "composite"
781f8c7532SHong Zhang #define MATFFT             "fft"
791f8c7532SHong Zhang #define MATFFTW            "fftw"
80e133240eSMatthew G Knepley #define MATSEQCUFFT        "seqcufft"
81557cca28SSatish Balay #define MATTRANSPOSEMAT    "transpose"
8272ca8751SBarry Smith #define MATSCHURCOMPLEMENT "schurcomplement"
831d6018f0SLisandro Dalcin #define MATPYTHON          "python"
84f91d8e95SBarry Smith #define MATHYPRESTRUCT     "hyprestruct"
85a9e6138eSGlenn Hammond #define MATHYPRESSTRUCT    "hypresstruct"
86ee1cef2cSJed Brown #define MATSUBMATRIX       "submatrix"
872c0dbf93SJed Brown #define MATLOCALREF        "localref"
88d8588912SDave May #define MATNEST            "nest"
89773941d6SBarry Smith 
9076bdecfbSBarry Smith /*J
91c7393fdbSBarry Smith     MatSolverPackage - String with the name of a PETSc matrix solver type.
929989ab13SBarry Smith 
939989ab13SBarry Smith     For example: "petsc" indicates what PETSc provides, "superlu" indicates either
949989ab13SBarry Smith        SuperLU or SuperLU_Dist etc.
959989ab13SBarry Smith 
969989ab13SBarry Smith 
979989ab13SBarry Smith    Level: beginner
989989ab13SBarry Smith 
995c9eb25fSBarry Smith .seealso: MatGetFactor(), Mat, MatSetType(), MatType
10076bdecfbSBarry Smith J*/
101c7393fdbSBarry Smith #define MatSolverPackage char*
1022692d6eeSBarry Smith #define MATSOLVERSUPERLU      "superlu"
1032692d6eeSBarry Smith #define MATSOLVERSUPERLU_DIST "superlu_dist"
1042692d6eeSBarry Smith #define MATSOLVERUMFPACK      "umfpack"
10520db9a53SJed Brown #define MATSOLVERCHOLMOD      "cholmod"
1062692d6eeSBarry Smith #define MATSOLVERESSL         "essl"
1072692d6eeSBarry Smith #define MATSOLVERLUSOL        "lusol"
1082692d6eeSBarry Smith #define MATSOLVERMUMPS        "mumps"
109*d615f992SSatish Balay #define MATSOLVERMKL_PARDISO  "mkl_pardiso"
1102692d6eeSBarry Smith #define MATSOLVERPASTIX       "pastix"
1112692d6eeSBarry Smith #define MATSOLVERMATLAB       "matlab"
1122692d6eeSBarry Smith #define MATSOLVERPETSC        "petsc"
1132692d6eeSBarry Smith #define MATSOLVERBAS          "bas"
1149ae82921SPaul Mullowney #define MATSOLVERCUSPARSE     "cusparse"
115aa5a9175SDahai Guo #define MATSOLVERBSTRM        "bstrm"
116aa5a9175SDahai Guo #define MATSOLVERSBSTRM       "sbstrm"
11715767789SHong Zhang #define MATSOLVERELEMENTAL    "elemental"
11817f1a0eaSHong Zhang #define MATSOLVERCLIQUE       "clique"
11969e15a41SShri Abhyankar #define MATSOLVERKLU          "klu"
120c0cdd4a1SDahai Guo 
121b24902e0SBarry Smith /*E
122b24902e0SBarry Smith     MatFactorType - indicates what type of factorization is requested
123b24902e0SBarry Smith 
124b24902e0SBarry Smith     Level: beginner
125b24902e0SBarry Smith 
126b24902e0SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
127b24902e0SBarry Smith 
128c7393fdbSBarry Smith .seealso: MatSolverPackage, MatGetFactor()
129b24902e0SBarry Smith E*/
130599ef60dSHong Zhang typedef enum {MAT_FACTOR_NONE, MAT_FACTOR_LU, MAT_FACTOR_CHOLESKY, MAT_FACTOR_ILU, MAT_FACTOR_ICC,MAT_FACTOR_ILUDT} MatFactorType;
131014dd563SJed Brown PETSC_EXTERN const char *const MatFactorTypes[];
132e92e720dSBarry Smith 
133014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactor(Mat,const MatSolverPackage,MatFactorType,Mat*);
134014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactorAvailable(Mat,const MatSolverPackage,MatFactorType,PetscBool *);
135014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFactorGetSolverPackage(Mat,const MatSolverPackage*);
136014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactorType(Mat,MatFactorType*);
1379989ab13SBarry Smith 
138c06d978dSMatthew Knepley /* Logging support */
1390700a824SBarry Smith #define    MAT_FILE_CLASSID 1211216    /* used to indicate matrices in binary files */
140014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_CLASSID;
141335efc43SPeter Brune PETSC_EXTERN PetscClassId MAT_COLORING_CLASSID;
142014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_FDCOLORING_CLASSID;
143014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_TRANSPOSECOLORING_CLASSID;
144014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_PARTITIONING_CLASSID;
145014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_COARSEN_CLASSID;
146014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_NULLSPACE_CLASSID;
147014dd563SJed Brown PETSC_EXTERN PetscClassId MATMFFD_CLASSID;
148c06d978dSMatthew Knepley 
149ceb03754SKris Buschelman /*E
150ceb03754SKris Buschelman     MatReuse - Indicates if matrices obtained from a previous call to MatGetSubMatrices()
151d6eff37eSBarry Smith      or MatGetSubMatrix() are to be reused to store the new matrix values. For MatConvert() is used to indicate
152d6eff37eSBarry Smith      that the input matrix is to be replaced with the converted matrix.
153ceb03754SKris Buschelman 
154ceb03754SKris Buschelman     Level: beginner
155ceb03754SKris Buschelman 
156ceb03754SKris Buschelman    Any additions/changes here MUST also be made in include/finclude/petscmat.h
157ceb03754SKris Buschelman 
1580c451bc4SBarry Smith .seealso: MatGetSubMatrices(), MatGetSubMatrix(), MatDestroyMatrices(), MatConvert()
159ceb03754SKris Buschelman E*/
160dfe085dbSJed Brown typedef enum {MAT_INITIAL_MATRIX,MAT_REUSE_MATRIX,MAT_IGNORE_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);
178ce1779c8SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatViewFromOptions(Mat A,const char prefix[],const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,prefix,name);}
179607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatRegisterAll(void);
180bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatRegister(const char[],PetscErrorCode(*)(Mat));
181014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRegisterBaseName(const char[],const char[],const char[]);
182014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetOptionsPrefix(Mat,const char[]);
183014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAppendOptionsPrefix(Mat,const char[]);
184014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOptionsPrefix(Mat,const char*[]);
185f69a0ea3SMatthew Knepley 
186014dd563SJed Brown PETSC_EXTERN PetscBool         MatRegisterAllCalled;
187140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatList;
188140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatColoringList;
189140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatPartitioningList;
190140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatCoarsenList;
19128988994SBarry Smith 
1923b224e63SBarry Smith /*E
193aa6c7ce3SBarry Smith     MatStructure - Indicates if two matrices have the same nonzero structure
1943b224e63SBarry Smith 
1953b224e63SBarry Smith     Level: beginner
1963b224e63SBarry Smith 
1973b224e63SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
1983b224e63SBarry Smith 
199aa6c7ce3SBarry Smith .seealso: MatCopy(), MatAXPY()
2003b224e63SBarry Smith E*/
201aa6c7ce3SBarry Smith typedef enum {DIFFERENT_NONZERO_PATTERN,SUBSET_NONZERO_PATTERN,SAME_NONZERO_PATTERN} MatStructure;
2023b224e63SBarry Smith 
203014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqDense(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*);
204014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateDense(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*);
205014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
206014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
207014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *);
208014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithSplitArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],PetscInt[],PetscInt[],PetscScalar[],Mat*);
2098d7a6e47SBarry Smith 
210014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
211014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
212014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat*);
213d21a29f3SJed Brown 
214014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAdj(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscInt[],Mat*);
215014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
216d21a29f3SJed Brown 
217014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
218014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPISBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *);
21938f409ebSLisandro Dalcin PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
220014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
2214cce697cSJed Brown PETSC_EXTERN PetscErrorCode MatXAIJSetPreallocation(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscInt[],const PetscInt[]);
222dfb205c3SBarry Smith 
223014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateShell(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,void *,Mat*);
224014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateNormal(Mat,Mat*);
225014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateLRC(Mat,Mat,Mat,Mat*);
226014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateIS(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,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 MatCreateSeqBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
231014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
232014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqSBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
233014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPISBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
234c0cdd4a1SDahai Guo 
235014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateScatter(MPI_Comm,VecScatter,Mat*);
236014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScatterSetVecScatter(Mat,VecScatter);
237014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScatterGetVecScatter(Mat,VecScatter*);
238014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateBlockMat(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,Mat*);
239014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeAddMat(Mat,Mat);
240014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeMerge(Mat);
241014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateComposite(MPI_Comm,PetscInt,const Mat*,Mat*);
2426d7c1e57SBarry Smith typedef enum {MAT_COMPOSITE_ADDITIVE,MAT_COMPOSITE_MULTIPLICATIVE} MatCompositeType;
243014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeSetType(Mat,MatCompositeType);
2446d7c1e57SBarry Smith 
24519fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatCreateFFT(MPI_Comm,PetscInt,const PetscInt[],MatType,Mat*);
246014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqCUFFT(MPI_Comm,PetscInt,const PetscInt[],Mat*);
247dedccee8SHong Zhang 
248014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateTranspose(Mat,Mat*);
249014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSubMatrix(Mat,IS,IS,Mat*);
250014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSubMatrixUpdate(Mat,Mat,IS,IS);
251014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateLocalRef(Mat,IS,IS,Mat*);
2521472f72bSBarry Smith 
253014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPythonSetType(Mat,const char[]);
2541d6018f0SLisandro Dalcin 
255014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetUp(Mat);
256014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroy(Mat*);
257e56f5c9eSBarry Smith PETSC_EXTERN PetscErrorCode MatGetNonzeroState(Mat,PetscObjectState*);
25821c89e3eSBarry Smith 
259014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatConjugate(Mat);
260014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRealPart(Mat);
261014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatImaginaryPart(Mat);
262014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetDiagonalBlock(Mat,Mat*);
263014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetTrace(Mat,PetscScalar*);
264713ccfa9SJed Brown PETSC_EXTERN PetscErrorCode MatInvertBlockDiagonal(Mat,const PetscScalar **);
26599cafbc1SBarry Smith 
2668ed539a5SBarry Smith /* ------------------------------------------------------------*/
267014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
268014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
269014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesRow(Mat,PetscInt,const PetscScalar[]);
270014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesRowLocal(Mat,PetscInt,const PetscScalar[]);
271014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBatch(Mat,PetscInt,PetscInt,PetscInt[],const PetscScalar[]);
27273a71a0fSBarry Smith PETSC_EXTERN PetscErrorCode MatSetRandom(Mat,PetscRandom);
27384cb2905SBarry Smith 
2742ef4de8bSBarry Smith /*S
2752ef4de8bSBarry Smith      MatStencil - Data structure (C struct) for storing information about a single row or
276be479b30SJed Brown         column of a matrix as indexed on an associated grid.
277be479b30SJed Brown 
278be479b30SJed Brown    Fortran usage is different, see MatSetValuesStencil() for details.
2792ef4de8bSBarry Smith 
2802ef4de8bSBarry Smith    Level: beginner
2812ef4de8bSBarry Smith 
2822ef4de8bSBarry Smith   Concepts: matrix; linear operator
2832ef4de8bSBarry Smith 
284be479b30SJed Brown .seealso:  MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockedStencil()
2852ef4de8bSBarry Smith S*/
286435da068SBarry Smith typedef struct {
287c1ac3661SBarry Smith   PetscInt k,j,i,c;
288435da068SBarry Smith } MatStencil;
2892ef4de8bSBarry Smith 
290014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
291014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
292014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt);
293435da068SBarry Smith 
294014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetColoring(Mat,ISColoring);
295014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesAdifor(Mat,PetscInt,void*);
2963a7fca6bSBarry Smith 
297d91e6319SBarry Smith /*E
298d91e6319SBarry Smith     MatAssemblyType - Indicates if the matrix is now to be used, or if you plan
299d91e6319SBarry Smith      to continue to add values to it
300d91e6319SBarry Smith 
301d91e6319SBarry Smith     Level: beginner
302d91e6319SBarry Smith 
303d91e6319SBarry Smith .seealso: MatAssemblyBegin(), MatAssemblyEnd()
304d91e6319SBarry Smith E*/
3056d4a8577SBarry Smith typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType;
306014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssemblyBegin(Mat,MatAssemblyType);
307014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssemblyEnd(Mat,MatAssemblyType);
308014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssembled(Mat,PetscBool *);
3094f9c727eSBarry Smith 
3101d73ed98SBarry Smith 
31130de9b25SBarry Smith 
312d91e6319SBarry Smith /*E
313d91e6319SBarry Smith     MatOption - Options that may be set for a matrix and its behavior or storage
314d91e6319SBarry Smith 
315d91e6319SBarry Smith     Level: beginner
316d91e6319SBarry Smith 
3170a835dfdSSatish Balay    Any additions/changes here MUST also be made in include/finclude/petscmat.h
318d91e6319SBarry Smith 
319382164b0SBarry Smith    Developer Notes: Entries that are negative need not be called collectively by all processes.
320382164b0SBarry Smith 
321d91e6319SBarry Smith .seealso: MatSetOption()
322d91e6319SBarry Smith E*/
32392d4d306SBarry Smith typedef enum {MAT_OPTION_MIN = -5,
32492d4d306SBarry Smith               MAT_NEW_NONZERO_LOCATION_ERR = -4,
32592d4d306SBarry Smith               MAT_UNUSED_NONZERO_LOCATION_ERR = -3,
32692d4d306SBarry Smith               MAT_NEW_NONZERO_ALLOCATION_ERR = -2,
32792d4d306SBarry Smith               MAT_ROW_ORIENTED = -1,
328382164b0SBarry Smith               MAT_SYMMETRIC = 1,
329382164b0SBarry Smith               MAT_STRUCTURALLY_SYMMETRIC = 2,
330382164b0SBarry Smith               MAT_NEW_DIAGONALS = 3,
331382164b0SBarry Smith               MAT_IGNORE_OFF_PROC_ENTRIES = 4,
332382164b0SBarry Smith               MAT_USE_HASH_TABLE = 5,
333382164b0SBarry Smith               MAT_KEEP_NONZERO_PATTERN = 6,
334382164b0SBarry Smith               MAT_IGNORE_ZERO_ENTRIES = 7,
335382164b0SBarry Smith               MAT_USE_INODES = 8,
336382164b0SBarry Smith               MAT_HERMITIAN = 9,
337382164b0SBarry Smith               MAT_SYMMETRY_ETERNAL = 10,
33811e456e1SBarry Smith               MAT_DUMMY = 11,
339382164b0SBarry Smith               MAT_IGNORE_LOWER_TRIANGULAR = 12,
340382164b0SBarry Smith               MAT_ERROR_LOWER_TRIANGULAR = 13,
341382164b0SBarry Smith               MAT_GETROW_UPPERTRIANGULAR = 14,
342382164b0SBarry Smith               MAT_SPD = 15,
34392d4d306SBarry Smith               MAT_NO_OFF_PROC_ZERO_ROWS = 16,
34492d4d306SBarry Smith               MAT_NO_OFF_PROC_ENTRIES = 17,
34592d4d306SBarry Smith               MAT_NEW_NONZERO_LOCATIONS = 18,
34692d4d306SBarry Smith               MAT_OPTION_MAX = 19} MatOption;
347e057df02SPaul Mullowney 
348014dd563SJed Brown PETSC_EXTERN const char *MatOptions[];
349014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetOption(Mat,MatOption,PetscBool );
35019fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatGetType(Mat,MatType*);
35184cb2905SBarry Smith 
352014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]);
353014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
354014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
355014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowUpperTriangular(Mat);
356014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRowUpperTriangular(Mat);
357014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
358014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
359014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnVector(Mat,Vec,PetscInt);
3608c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJGetArray(Mat,PetscScalar *[]);
3618c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJRestoreArray(Mat,PetscScalar *[]);
362bd04181cSBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJSetValuesLocalFast(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
3638c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatDenseGetArray(Mat,PetscScalar *[]);
3648c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatDenseRestoreArray(Mat,PetscScalar *[]);
365014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBlockSize(Mat,PetscInt *);
366014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSize(Mat,PetscInt);
367014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBlockSizes(Mat,PetscInt *,PetscInt *);
368014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSizes(Mat,PetscInt,PetscInt);
36933d57670SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSizesFromMats(Mat,Mat,Mat);
370014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNThreads(Mat,PetscInt);
371014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNThreads(Mat,PetscInt*);
3721620fd73SBarry Smith 
373014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMult(Mat,Vec,Vec);
374014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultDiagonalBlock(Mat,Vec,Vec);
375014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultAdd(Mat,Vec,Vec,Vec);
376014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTranspose(Mat,Vec,Vec);
377014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultHermitianTranspose(Mat,Vec,Vec);
378014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsTranspose(Mat,Mat,PetscReal,PetscBool *);
379014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitianTranspose(Mat,Mat,PetscReal,PetscBool *);
380014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeAdd(Mat,Vec,Vec,Vec);
381014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultHermitianTransposeAdd(Mat,Vec,Vec,Vec);
382014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultConstrained(Mat,Vec,Vec);
383014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeConstrained(Mat,Vec,Vec);
384014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatSolve(Mat,Mat,Mat);
385f9426fe0SMark Adams PETSC_EXTERN PetscErrorCode MatResidual(Mat,Vec,Vec,Vec);
386f5edf698SHong Zhang 
387d91e6319SBarry Smith /*E
388d91e6319SBarry Smith     MatDuplicateOption - Indicates if a duplicated sparse matrix should have
389d91e6319SBarry Smith   its numerical values copied over or just its nonzero structure.
390d91e6319SBarry Smith 
391d91e6319SBarry Smith     Level: beginner
392d91e6319SBarry Smith 
393d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
394d91e6319SBarry Smith 
39570dcbbb9SBarry Smith $   MAT_SHARE_NONZERO_PATTERN - the i and j arrays in the new matrix will be shared with the original matrix
39670dcbbb9SBarry Smith $                               this also triggers the MAT_DO_NOT_COPY_VALUES option. This is used when you
39770dcbbb9SBarry Smith $                               have several matrices with the same nonzero pattern.
39870dcbbb9SBarry Smith 
399d91e6319SBarry Smith .seealso: MatDuplicate()
400d91e6319SBarry Smith E*/
40170dcbbb9SBarry Smith typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES,MAT_SHARE_NONZERO_PATTERN} MatDuplicateOption;
4022e8a6d31SBarry Smith 
40319fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatConvert(Mat,MatType,MatReuse,Mat*);
404014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDuplicate(Mat,MatDuplicateOption,Mat*);
40594a9d846SBarry Smith 
406cb5b572fSBarry Smith 
407014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCopy(Mat,Mat,MatStructure);
408014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatView(Mat,PetscViewer);
409014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsSymmetric(Mat,PetscReal,PetscBool *);
410014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsStructurallySymmetric(Mat,PetscBool *);
411014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitian(Mat,PetscReal,PetscBool *);
412014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsSymmetricKnown(Mat,PetscBool *,PetscBool *);
413014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitianKnown(Mat,PetscBool *,PetscBool *);
414014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMissingDiagonal(Mat,PetscBool  *,PetscInt *);
415014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLoad(Mat, PetscViewer);
4167b80b807SBarry Smith 
4171a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool  *);
4181a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool  *);
4191a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool  *);
4201a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool  *);
421d4fbbf0eSBarry Smith 
422d91e6319SBarry Smith /*S
423d91e6319SBarry Smith      MatInfo - Context of matrix information, used with MatGetInfo()
424d91e6319SBarry Smith 
425d91e6319SBarry Smith    In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE
426d91e6319SBarry Smith 
427d91e6319SBarry Smith    Level: intermediate
428d91e6319SBarry Smith 
429d91e6319SBarry Smith   Concepts: matrix^nonzero information
430d91e6319SBarry Smith 
431d9274352SBarry Smith .seealso:  MatGetInfo(), MatInfoType
432d91e6319SBarry Smith S*/
4334e220ebcSLois Curfman McInnes typedef struct {
434b0a32e0cSBarry Smith   PetscLogDouble block_size;                         /* block size */
435b0a32e0cSBarry Smith   PetscLogDouble nz_allocated,nz_used,nz_unneeded;   /* number of nonzeros */
436b0a32e0cSBarry Smith   PetscLogDouble memory;                             /* memory allocated */
437b0a32e0cSBarry Smith   PetscLogDouble assemblies;                         /* number of matrix assemblies called */
438b0a32e0cSBarry Smith   PetscLogDouble mallocs;                            /* number of mallocs during MatSetValues() */
439b0a32e0cSBarry Smith   PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */
440b0a32e0cSBarry Smith   PetscLogDouble factor_mallocs;                     /* number of mallocs during factorization */
4414e220ebcSLois Curfman McInnes } MatInfo;
4424e220ebcSLois Curfman McInnes 
443d9274352SBarry Smith /*E
444d9274352SBarry Smith     MatInfoType - Indicates if you want information about the local part of the matrix,
445d9274352SBarry Smith      the entire parallel matrix or the maximum over all the local parts.
446d9274352SBarry Smith 
447d9274352SBarry Smith     Level: beginner
448d9274352SBarry Smith 
449d9274352SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
450d9274352SBarry Smith 
451d9274352SBarry Smith .seealso: MatGetInfo(), MatInfo
452d9274352SBarry Smith E*/
4537b80b807SBarry Smith typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType;
454014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetInfo(Mat,MatInfoType,MatInfo*);
455014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetDiagonal(Mat,Vec);
456014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMax(Mat,Vec,PetscInt[]);
457014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMin(Mat,Vec,PetscInt[]);
458014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMaxAbs(Mat,Vec,PetscInt[]);
459014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMinAbs(Mat,Vec,PetscInt[]);
460014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowSum(Mat,Vec);
461014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTranspose(Mat,MatReuse,Mat*);
462014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatHermitianTranspose(Mat,MatReuse,Mat*);
463014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPermute(Mat,IS,IS,Mat *);
464014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalScale(Mat,Vec,Vec);
465014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalSet(Mat,Vec,InsertMode);
466014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatEqual(Mat,Mat,PetscBool *);
467014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultEqual(Mat,Mat,PetscInt,PetscBool *);
468014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultAddEqual(Mat,Mat,PetscInt,PetscBool *);
469014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeEqual(Mat,Mat,PetscInt,PetscBool *);
470014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscBool *);
4717b80b807SBarry Smith 
472014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNorm(Mat,NormType,PetscReal *);
473014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnNorms(Mat,NormType,PetscReal *);
474014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroEntries(Mat);
475014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
476014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsIS(Mat,IS,PetscScalar,Vec,Vec);
477014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsStencil(Mat,PetscInt,const MatStencil [],PetscScalar,Vec,Vec);
478014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsStencil(Mat,PetscInt,const MatStencil[],PetscScalar,Vec,Vec);
479014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumns(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
480014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsIS(Mat,IS,PetscScalar,Vec,Vec);
4817b80b807SBarry Smith 
482014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSize(Mat,PetscInt*,PetscInt*);
483014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalSize(Mat,PetscInt*,PetscInt*);
484014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRange(Mat,PetscInt*,PetscInt*);
485014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRanges(Mat,const PetscInt**);
486014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRangeColumn(Mat,PetscInt*,PetscInt*);
487014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRangesColumn(Mat,const PetscInt**);
488566876eaSJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipIS(Mat,IS*,IS*);
4895ef9f2a5SBarry Smith 
490014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
491014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSubMatricesParallel(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
492014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroyMatrices(PetscInt,Mat *[]);
493014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSubMatrix(Mat,IS,IS,MatReuse,Mat *);
494014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalSubMatrix(Mat,IS,IS,Mat*);
495014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreLocalSubMatrix(Mat,IS,IS,Mat*);
496014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSeqNonzeroStructure(Mat,Mat*);
497014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroySeqNonzeroStructure(Mat*);
4987b80b807SBarry Smith 
499014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJConcatenateSeqAIJ(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
500014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJConcatenateSeqAIJSymbolic(MPI_Comm,Mat,PetscInt,Mat*);
501014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJConcatenateSeqAIJNumeric(MPI_Comm,Mat,PetscInt,Mat);
502014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJ(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*);
503014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJSymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*);
504014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJNumeric(Mat,Mat);
505014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMat(Mat,MatReuse,Mat*);
506014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*);
507014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,Mat*);
508014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetGhosts(Mat, PetscInt *,const PetscInt *[]);
5098efafbd8SBarry Smith 
510014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt);
5117b80b807SBarry Smith 
512014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*);
513014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMultSymbolic(Mat,Mat,PetscReal,Mat*);
514014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMultNumeric(Mat,Mat,Mat);
51522440eb1SKris Buschelman 
5167bab7c10SHong Zhang PETSC_EXTERN PetscErrorCode MatMatMatMult(Mat,Mat,Mat,MatReuse,PetscReal,Mat*);
5177bab7c10SHong Zhang PETSC_EXTERN PetscErrorCode MatMatMatMultSymbolic(Mat,Mat,Mat,PetscReal,Mat*);
5187bab7c10SHong Zhang PETSC_EXTERN PetscErrorCode MatMatMatMultNumeric(Mat,Mat,Mat,Mat);
5197bab7c10SHong Zhang 
520014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*);
521014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAPSymbolic(Mat,Mat,PetscReal,Mat*);
522014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAPNumeric(Mat,Mat,Mat);
523014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARt(Mat,Mat,MatReuse,PetscReal,Mat*);
524014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARtSymbolic(Mat,Mat,PetscReal,Mat*);
525014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARtNumeric(Mat,Mat,Mat);
52622440eb1SKris Buschelman 
527014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeMatMult(Mat,Mat,MatReuse,PetscReal,Mat*);
528014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposetMatMultSymbolic(Mat,Mat,PetscReal,Mat*);
529014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposetMatMultNumeric(Mat,Mat,Mat);
530014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMult(Mat,Mat,MatReuse,PetscReal,Mat*);
531014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMultSymbolic(Mat,Mat,PetscReal,Mat*);
532014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMultNumeric(Mat,Mat,Mat);
533bc011b1eSHong Zhang 
534014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAXPY(Mat,PetscScalar,Mat,MatStructure);
535014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAYPX(Mat,PetscScalar,Mat,MatStructure);
5361c741599SHong Zhang 
537014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScale(Mat,PetscScalar);
538014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShift(Mat,PetscScalar);
5397b80b807SBarry Smith 
540014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping);
541014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*);
542a93ff8c4SPeter Brune PETSC_EXTERN PetscErrorCode MatGetLayouts(Mat,PetscLayout*,PetscLayout*);
543014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
544014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsLocalIS(Mat,IS,PetscScalar,Vec,Vec);
545014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
546014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocalIS(Mat,IS,PetscScalar,Vec,Vec);
547014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
548014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
549052efed2SBarry Smith 
550014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStashSetInitialSize(Mat,PetscInt,PetscInt);
551014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
55290f02eecSBarry Smith 
553014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInterpolate(Mat,Vec,Vec);
554014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInterpolateAdd(Mat,Vec,Vec,Vec);
555014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestrict(Mat,Vec,Vec);
556014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetVecs(Mat,Vec*,Vec*);
557b2bf6370SHong Zhang PETSC_EXTERN PetscErrorCode MatGetRedundantMatrix(Mat,PetscInt,MPI_Comm,MatReuse,Mat*);
558014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetMultiProcBlock(Mat,MPI_Comm,MatReuse,Mat*);
559014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFindZeroDiagonals(Mat,IS*);
5603a062f41SBarry Smith PETSC_EXTERN PetscErrorCode MatFindOffBlockDiagonalEntries(Mat,IS*);
561bd481603SBarry Smith 
562bd481603SBarry Smith /*MC
5631620fd73SBarry Smith    MatSetValue - Set a single entry into a matrix.
5641620fd73SBarry Smith 
5651620fd73SBarry Smith    Not collective
5661620fd73SBarry Smith 
5671620fd73SBarry Smith    Input Parameters:
5681620fd73SBarry Smith +  m - the matrix
5691620fd73SBarry Smith .  row - the row location of the entry
5701620fd73SBarry Smith .  col - the column location of the entry
5711620fd73SBarry Smith .  value - the value to insert
5721620fd73SBarry Smith -  mode - either INSERT_VALUES or ADD_VALUES
5731620fd73SBarry Smith 
5741620fd73SBarry Smith    Notes:
5751620fd73SBarry Smith    For efficiency one should use MatSetValues() and set several or many
5761620fd73SBarry Smith    values simultaneously if possible.
5771620fd73SBarry Smith 
5781620fd73SBarry Smith    Level: beginner
5791620fd73SBarry Smith 
5801620fd73SBarry Smith .seealso: MatSetValues(), MatSetValueLocal()
5811620fd73SBarry Smith M*/
5821620fd73SBarry 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);}
5831620fd73SBarry Smith 
5841620fd73SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatGetValue(Mat v,PetscInt i,PetscInt j,PetscScalar *va) {return MatGetValues(v,1,&i,1,&j,va);}
5851620fd73SBarry Smith 
5861620fd73SBarry 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);}
5871620fd73SBarry Smith 
5881620fd73SBarry Smith /*MC
5890d2a7ff2SSatish Balay    MatPreallocateInitialize - Begins the block of code that will count the number of nonzeros per
590bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
591bd481603SBarry Smith 
592bd481603SBarry Smith    Synopsis:
593aaa7dc30SBarry Smith    #include <petscmat.h>
594c1ac3661SBarry Smith    PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)
595bd481603SBarry Smith 
596bd481603SBarry Smith    Collective on MPI_Comm
597bd481603SBarry Smith 
598bd481603SBarry Smith    Input Parameters:
599bd481603SBarry Smith +  comm - the communicator that will share the eventually allocated matrix
600859fcb39SBarry Smith .  nrows - the number of LOCAL rows in the matrix
601859fcb39SBarry Smith -  ncols - the number of LOCAL columns in the matrix
602bd481603SBarry Smith 
603bd481603SBarry Smith    Output Parameters:
604bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
605bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
606bd481603SBarry Smith 
607bd481603SBarry Smith    Level: intermediate
608bd481603SBarry Smith 
609bd481603SBarry Smith    Notes:
610a7f22e61SSatish Balay     See Users-Manual: ch_performance for more details.
611bd481603SBarry Smith 
6121d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
613bd481603SBarry Smith 
614bd481603SBarry Smith    Use MatPreallocateInitializeSymmetric() for symmetric matrices (MPISBAIJ matrices)
615bd481603SBarry Smith 
6161620fd73SBarry Smith    This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize().
6171620fd73SBarry Smith 
618bd481603SBarry Smith   Concepts: preallocation^Matrix
619bd481603SBarry Smith 
620d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(),
621d6e23781SBarry Smith           MatPreallocateInitializeSymmetric(), MatPreallocateSymmetricSetLocalBlock()
622bd481603SBarry Smith M*/
623c4f061fbSSatish Balay #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \
6247c922b88SBarry Smith { \
625a89bc211SBarry Smith   PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__start,__end; \
6261795a4d1SJed Brown   _4_ierr = PetscCalloc2(__nrows,&dnz,__nrows,&onz);CHKERRQ(_4_ierr); \
6271795a4d1SJed Brown   __start = 0; __end = __start;                                         \
628c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ctmp;\
629a89bc211SBarry Smith   _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows;
6307c922b88SBarry Smith 
631bd481603SBarry Smith /*MC
632bd481603SBarry Smith    MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
633bd481603SBarry Smith        inserted using a local number of the rows and columns
634bd481603SBarry Smith 
635bd481603SBarry Smith    Synopsis:
636aaa7dc30SBarry Smith    #include <petscmat.h>
637c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
638bd481603SBarry Smith 
639bd481603SBarry Smith    Not Collective
640bd481603SBarry Smith 
641bd481603SBarry Smith    Input Parameters:
642784ac674SJed Brown +  map - the row mapping from local numbering to global numbering
643bd481603SBarry Smith .  nrows - the number of rows indicated
6441d73ed98SBarry Smith .  rows - the indices of the rows
645784ac674SJed Brown .  cmap - the column mapping from local to global numbering
646bd481603SBarry Smith .  ncols - the number of columns in the matrix
647bd481603SBarry Smith .  cols - the columns indicated
648bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
649bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
650bd481603SBarry Smith 
651bd481603SBarry Smith    Level: intermediate
652bd481603SBarry Smith 
653bd481603SBarry Smith    Notes:
654a7f22e61SSatish Balay     See Users-Manual: ch_performance for more details.
655bd481603SBarry Smith 
6561d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
657bd481603SBarry Smith 
658bd481603SBarry Smith   Concepts: preallocation^Matrix
659bd481603SBarry Smith 
660d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(),
661d6e23781SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock()
662bd481603SBarry Smith M*/
663784ac674SJed Brown #define MatPreallocateSetLocal(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \
664c4f061fbSSatish Balay {\
665c1ac3661SBarry Smith   PetscInt __l;\
666784ac674SJed Brown   _4_ierr = ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\
667784ac674SJed Brown   _4_ierr = ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\
668c4f061fbSSatish Balay   for (__l=0;__l<nrows;__l++) {\
669ef66eb69SBarry Smith     _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
670c4f061fbSSatish Balay   }\
671c4f061fbSSatish Balay }
672c4f061fbSSatish Balay 
673bd481603SBarry Smith /*MC
674d6e23781SBarry Smith    MatPreallocateSetLocalBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be
675bd481603SBarry Smith        inserted using a local number of the rows and columns
676bd481603SBarry Smith 
677bd481603SBarry Smith    Synopsis:
678aaa7dc30SBarry Smith    #include <petscmat.h>
679d6e23781SBarry Smith    PetscErrorCode MatPreallocateSetLocalBlock(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
680d6e23781SBarry Smith 
681d6e23781SBarry Smith    Not Collective
682d6e23781SBarry Smith 
683d6e23781SBarry Smith    Input Parameters:
684d6e23781SBarry Smith +  map - the row mapping from local numbering to global numbering
685d6e23781SBarry Smith .  nrows - the number of rows indicated
686d6e23781SBarry Smith .  rows - the indices of the rows
687d6e23781SBarry Smith .  cmap - the column mapping from local to global numbering
688d6e23781SBarry Smith .  ncols - the number of columns in the matrix
689d6e23781SBarry Smith .  cols - the columns indicated
690d6e23781SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
691d6e23781SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
692d6e23781SBarry Smith 
693d6e23781SBarry Smith    Level: intermediate
694d6e23781SBarry Smith 
695d6e23781SBarry Smith    Notes:
696d6e23781SBarry Smith     See Users-Manual: ch_performance for more details.
697d6e23781SBarry Smith 
698d6e23781SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
699d6e23781SBarry Smith 
700d6e23781SBarry Smith   Concepts: preallocation^Matrix
701d6e23781SBarry Smith 
702d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(),
703d6e23781SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock()
704d6e23781SBarry Smith M*/
705d6e23781SBarry Smith #define MatPreallocateSetLocalBlock(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \
706d6e23781SBarry Smith {\
707d6e23781SBarry Smith   PetscInt __l;\
708d6e23781SBarry Smith   _4_ierr = ISLocalToGlobalMappingApplyBlock(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\
709d6e23781SBarry Smith   _4_ierr = ISLocalToGlobalMappingApplyBlock(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\
710d6e23781SBarry Smith   for (__l=0;__l<nrows;__l++) {\
711d6e23781SBarry Smith     _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
712d6e23781SBarry Smith   }\
713d6e23781SBarry Smith }
714d6e23781SBarry Smith 
715d6e23781SBarry Smith /*MC
716d6e23781SBarry Smith    MatPreallocateSymmetricSetLocalBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be
717d6e23781SBarry Smith        inserted using a local number of the rows and columns
718d6e23781SBarry Smith 
719d6e23781SBarry Smith    Synopsis:
720d6e23781SBarry Smith    #include <petscmat.h>
721d6e23781SBarry Smith    PetscErrorCode MatPreallocateSymmetricSetLocalBlock(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
722bd481603SBarry Smith 
723bd481603SBarry Smith    Not Collective
724bd481603SBarry Smith 
725bd481603SBarry Smith    Input Parameters:
726bd481603SBarry Smith +  map - the mapping between local numbering and global numbering
727bd481603SBarry Smith .  nrows - the number of rows indicated
7281d73ed98SBarry Smith .  rows - the indices of the rows
729bd481603SBarry Smith .  ncols - the number of columns in the matrix
730bd481603SBarry Smith .  cols - the columns indicated
731bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
732bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
733bd481603SBarry Smith 
734bd481603SBarry Smith    Level: intermediate
735bd481603SBarry Smith 
736bd481603SBarry Smith    Notes:
737a7f22e61SSatish Balay     See Users-Manual: ch_performance for more details.
738bd481603SBarry Smith 
739bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
740bd481603SBarry Smith 
741bd481603SBarry Smith   Concepts: preallocation^Matrix
742bd481603SBarry Smith 
743d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateInitialize(),
744d6e23781SBarry Smith           MatPreallocateInitialize(),  MatPreallocateSetLocal()
745bd481603SBarry Smith M*/
746d6e23781SBarry Smith #define MatPreallocateSymmetricSetLocalBlock(map,nrows,rows,ncols,cols,dnz,onz) 0;\
747d3d32019SBarry Smith {\
748c1ac3661SBarry Smith   PetscInt __l;\
749d6e23781SBarry Smith   _4_ierr = ISLocalToGlobalMappingApplyBlock(map,nrows,rows,rows);CHKERRQ(_4_ierr);\
750d6e23781SBarry Smith   _4_ierr = ISLocalToGlobalMappingApplyBlock(map,ncols,cols,cols);CHKERRQ(_4_ierr);\
751d3d32019SBarry Smith   for (__l=0;__l<nrows;__l++) {\
752d6e23781SBarry Smith     _4_ierr = MatPreallocateSymmetricSetBlock((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
753d3d32019SBarry Smith   }\
754d3d32019SBarry Smith }
755bd481603SBarry Smith /*MC
756bd481603SBarry Smith    MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
757bd481603SBarry Smith        inserted using a local number of the rows and columns
758bd481603SBarry Smith 
759bd481603SBarry Smith    Synopsis:
760aaa7dc30SBarry Smith    #include <petscmat.h>
761c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
762bd481603SBarry Smith 
763bd481603SBarry Smith    Not Collective
764bd481603SBarry Smith 
765bd481603SBarry Smith    Input Parameters:
76664075487SBarry Smith +  row - the row
767bd481603SBarry Smith .  ncols - the number of columns in the matrix
768a50f8bf6SHong Zhang -  cols - the columns indicated
769a50f8bf6SHong Zhang 
770a50f8bf6SHong Zhang    Output Parameters:
771a50f8bf6SHong Zhang +  dnz - the array that will be passed to the matrix preallocation routines
772bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
773bd481603SBarry Smith 
774bd481603SBarry Smith    Level: intermediate
775bd481603SBarry Smith 
776bd481603SBarry Smith    Notes:
777a7f22e61SSatish Balay     See Users-Manual: ch_performance for more details.
778bd481603SBarry Smith 
779bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
780bd481603SBarry Smith 
7811620fd73SBarry Smith    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
7821620fd73SBarry Smith 
783bd481603SBarry Smith   Concepts: preallocation^Matrix
784bd481603SBarry Smith 
785d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(),
786d6e23781SBarry Smith           MatPreallocateInitialize(), MatPreallocateSetLocal()
787bd481603SBarry Smith M*/
788c4f061fbSSatish Balay #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\
789c1ac3661SBarry Smith { PetscInt __i; \
790e32f2f54SBarry 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);\
791e32f2f54SBarry 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);\
7927c922b88SBarry Smith   for (__i=0; __i<nc; __i++) {\
79364075487SBarry Smith     if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \
7947cc688d7SBarry Smith     else dnz[row - __rstart]++;\
7957c922b88SBarry Smith   }\
7967c922b88SBarry Smith }
7977c922b88SBarry Smith 
798bd481603SBarry Smith /*MC
799d6e23781SBarry Smith    MatPreallocateSymmetricSetBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be
800bd481603SBarry Smith        inserted using a local number of the rows and columns
801bd481603SBarry Smith 
802bd481603SBarry Smith    Synopsis:
803aaa7dc30SBarry Smith    #include <petscmat.h>
804d6e23781SBarry Smith    PetscErrorCode MatPreallocateSymmetricSetBlock(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
805bd481603SBarry Smith 
806bd481603SBarry Smith    Not Collective
807bd481603SBarry Smith 
808bd481603SBarry Smith    Input Parameters:
809bd481603SBarry Smith +  nrows - the number of rows indicated
8101d73ed98SBarry Smith .  rows - the indices of the rows
811bd481603SBarry Smith .  ncols - the number of columns in the matrix
812bd481603SBarry Smith .  cols - the columns indicated
813bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
814bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
815bd481603SBarry Smith 
816bd481603SBarry Smith    Level: intermediate
817bd481603SBarry Smith 
818bd481603SBarry Smith    Notes:
819a7f22e61SSatish Balay     See Users-Manual: ch_performance for more details.
820bd481603SBarry Smith 
821bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
822bd481603SBarry Smith 
8231620fd73SBarry Smith    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
8241620fd73SBarry Smith 
825bd481603SBarry Smith   Concepts: preallocation^Matrix
826bd481603SBarry Smith 
827d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(),  MatPreallocateInitialize(),
828d6e23781SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocal()
829bd481603SBarry Smith M*/
830d6e23781SBarry Smith #define MatPreallocateSymmetricSetBlock(row,nc,cols,dnz,onz) 0;\
831c1ac3661SBarry Smith { PetscInt __i; \
832d3d32019SBarry Smith   for (__i=0; __i<nc; __i++) {\
833d3d32019SBarry Smith     if (cols[__i] >= __end) onz[row - __rstart]++; \
834d3d32019SBarry Smith     else if (cols[__i] >= row) dnz[row - __rstart]++;\
835d3d32019SBarry Smith   }\
836d3d32019SBarry Smith }
837d3d32019SBarry Smith 
838bd481603SBarry Smith /*MC
83916371a99SBarry Smith    MatPreallocateLocation -  An alternative to MatPreallocationSet() that puts the nonzero locations into the matrix if it exists
84016371a99SBarry Smith 
84116371a99SBarry Smith    Synopsis:
842aaa7dc30SBarry Smith    #include <petscmat.h>
84316371a99SBarry Smith    PetscErrorCode MatPreallocateLocations(Mat A,PetscInt row,PetscInt ncols,PetscInt *cols,PetscInt *dnz,PetscInt *onz)
84416371a99SBarry Smith 
84516371a99SBarry Smith    Not Collective
84616371a99SBarry Smith 
84716371a99SBarry Smith    Input Parameters:
84816371a99SBarry Smith .  A - matrix
84916371a99SBarry Smith .  row - row where values exist (must be local to this process)
85016371a99SBarry Smith .  ncols - number of columns
85116371a99SBarry Smith .  cols - columns with nonzeros
85216371a99SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
85316371a99SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
85416371a99SBarry Smith 
85516371a99SBarry Smith    Level: intermediate
85616371a99SBarry Smith 
85716371a99SBarry Smith    Notes:
858a7f22e61SSatish Balay     See Users-Manual: ch_performance for more details.
85916371a99SBarry Smith 
86016371a99SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
86116371a99SBarry Smith 
86216371a99SBarry Smith    This is a MACRO not a function because it uses a bunch of variables private to the MatPreallocation.... routines.
86316371a99SBarry Smith 
86416371a99SBarry Smith   Concepts: preallocation^Matrix
86516371a99SBarry Smith 
866d6e23781SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(),
867d6e23781SBarry Smith           MatPreallocateSymmetricSetLocalBlock()
86816371a99SBarry Smith M*/
8690298fd71SBarry 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);}
87016371a99SBarry Smith 
87116371a99SBarry Smith 
87216371a99SBarry Smith /*MC
8730d2a7ff2SSatish Balay    MatPreallocateFinalize - Ends the block of code that will count the number of nonzeros per
874bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
875bd481603SBarry Smith 
876bd481603SBarry Smith    Synopsis:
877aaa7dc30SBarry Smith    #include <petscmat.h>
878c1ac3661SBarry Smith    PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz)
879bd481603SBarry Smith 
880bd481603SBarry Smith    Collective on MPI_Comm
881bd481603SBarry Smith 
882bd481603SBarry Smith    Input Parameters:
88316371a99SBarry Smith +  dnz - the array that was be passed to the matrix preallocation routines
884bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
885bd481603SBarry Smith 
886bd481603SBarry Smith    Level: intermediate
887bd481603SBarry Smith 
888bd481603SBarry Smith    Notes:
889a7f22e61SSatish Balay     See Users-Manual: ch_performance for more details.
890bd481603SBarry Smith 
891bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
892bd481603SBarry Smith 
8931620fd73SBarry Smith    This is a MACRO not a function because it closes the { started in MatPreallocateInitialize().
8941620fd73SBarry Smith 
895bd481603SBarry Smith   Concepts: preallocation^Matrix
896bd481603SBarry Smith 
897d6e23781SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(),
898d6e23781SBarry Smith           MatPreallocateSymmetricSetLocalBlock()
899bd481603SBarry Smith M*/
900a89bc211SBarry Smith #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree2(dnz,onz);CHKERRQ(_4_ierr);}
9017c922b88SBarry Smith 
9027b80b807SBarry Smith /* Routines unique to particular data structures */
903014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellGetContext(Mat,void *);
904698d4c6aSKris Buschelman 
905014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInodeAdjustForInodes(Mat,IS*,IS*);
906014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *);
907698d4c6aSKris Buschelman 
908014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetColumnIndices(Mat,PetscInt[]);
909014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetColumnIndices(Mat,PetscInt[]);
910014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
911014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
912014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
913014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJFromTriple(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*,PetscInt,PetscBool);
9147b80b807SBarry Smith 
915a96a251dSBarry Smith #define MAT_SKIP_ALLOCATION -4
916a96a251dSBarry Smith 
917014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
918014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
919014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]);
920273d9f13SBarry Smith 
921014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
922014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
923014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
924014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []);
925014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
926014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]);
927014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
928014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]);
929014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIDenseSetPreallocation(Mat,PetscScalar[]);
930014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqDenseSetPreallocation(Mat,PetscScalar[]);
9319230625dSJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,const PetscInt*[]);
9329230625dSJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,const PetscInt*[]);
933014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAdjCreateNonemptySubcommMat(Mat,Mat*);
934273d9f13SBarry Smith 
935014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqDenseSetLDA(Mat,PetscInt);
936014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDenseGetLocalMatrix(Mat,Mat*);
9371b807ce4Svictorle 
938014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStoreValues(Mat);
939014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRetrieveValues(Mat);
9402e8a6d31SBarry Smith 
941014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDAADSetCtx(Mat,void*);
9423a7fca6bSBarry Smith 
943014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFindNonzeroRows(Mat,IS*);
9447b80b807SBarry Smith /*
9457b80b807SBarry Smith   These routines are not usually accessed directly, rather solving is
94694b7f48cSBarry Smith   done through the KSP and PC interfaces.
9477b80b807SBarry Smith */
9487b80b807SBarry Smith 
94976bdecfbSBarry Smith /*J
9508f6c3df8SBarry Smith     MatOrderingType - String with the name of a PETSc matrix ordering
951d9274352SBarry Smith 
952d9274352SBarry Smith    Level: beginner
953d9274352SBarry Smith 
954d9274352SBarry Smith .seealso: MatGetOrdering()
95576bdecfbSBarry Smith J*/
95619fd82e9SBarry Smith typedef const char* MatOrderingType;
9572692d6eeSBarry Smith #define MATORDERINGNATURAL     "natural"
9582692d6eeSBarry Smith #define MATORDERINGND          "nd"
9592692d6eeSBarry Smith #define MATORDERING1WD         "1wd"
9602692d6eeSBarry Smith #define MATORDERINGRCM         "rcm"
9612692d6eeSBarry Smith #define MATORDERINGQMD         "qmd"
9622692d6eeSBarry Smith #define MATORDERINGROWLENGTH   "rowlength"
963510184a7SJed Brown #define MATORDERINGWBM         "wbm"
964fc1bef75SJed Brown #define MATORDERINGSPECTRAL    "spectral"
965312e372aSBarry Smith #define MATORDERINGAMD         "amd"            /* only works if UMFPACK is installed with PETSc */
966b12f92e5SBarry Smith 
96719fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatGetOrdering(Mat,MatOrderingType,IS*,IS*);
968140e18c1SBarry Smith PETSC_EXTERN PetscErrorCode MatGetOrderingList(PetscFunctionList*);
969bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatOrderingRegister(const char[],PetscErrorCode(*)(Mat,MatOrderingType,IS*,IS*));
970607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatOrderingRegisterAll(void);
971014dd563SJed Brown PETSC_EXTERN PetscBool         MatOrderingRegisterAllCalled;
972140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatOrderingList;
973d4fbbf0eSBarry Smith 
974014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS);
975fc1bef75SJed Brown PETSC_EXTERN PetscErrorCode MatCreateLaplacian(Mat,PetscReal,PetscBool,Mat*);
976a2ce50c7SBarry Smith 
977d91e6319SBarry Smith /*S
978d90ac83dSHong Zhang     MatFactorShiftType - Numeric Shift.
979d90ac83dSHong Zhang 
980d90ac83dSHong Zhang    Level: beginner
981d90ac83dSHong Zhang 
982d90ac83dSHong Zhang S*/
983d90ac83dSHong Zhang typedef enum {MAT_SHIFT_NONE,MAT_SHIFT_NONZERO,MAT_SHIFT_POSITIVE_DEFINITE,MAT_SHIFT_INBLOCKS} MatFactorShiftType;
9846a6fc655SJed Brown PETSC_EXTERN const char *const MatFactorShiftTypes[];
9855e9742b9SJed Brown PETSC_EXTERN const char *const MatFactorShiftTypesDetail[];
986d90ac83dSHong Zhang 
987d90ac83dSHong Zhang /*S
9882401956bSBarry Smith    MatFactorInfo - Data passed into the matrix factorization routines
989b00f7748SHong Zhang 
99061cad860SBarry Smith    In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE, that is use
99161cad860SBarry Smith $     MatFactorInfo  info(MAT_FACTORINFO_SIZE)
992b00f7748SHong Zhang 
99315e8a5b3SHong Zhang    Notes: These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC.
994b00f7748SHong Zhang 
99561cad860SBarry Smith       You can use MatFactorInfoInitialize() to set default values.
99661cad860SBarry Smith 
997b00f7748SHong Zhang    Level: developer
998b00f7748SHong Zhang 
999d7d82daaSBarry Smith .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(),
1000d7d82daaSBarry Smith           MatFactorInfoInitialize()
1001b00f7748SHong Zhang 
1002b00f7748SHong Zhang S*/
1003b00f7748SHong Zhang typedef struct {
100415e8a5b3SHong Zhang   PetscReal     diagonal_fill;  /* force diagonal to fill in if initially not filled */
100585317021SBarry Smith   PetscReal     usedt;
100615e8a5b3SHong Zhang   PetscReal     dt;             /* drop tolerance */
1007b00f7748SHong Zhang   PetscReal     dtcol;          /* tolerance for pivoting */
100815e8a5b3SHong Zhang   PetscReal     dtcount;        /* maximum nonzeros to be allowed per row */
100967e28bfeSBarry Smith   PetscReal     fill;           /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */
1010348344bbSBarry Smith   PetscReal     levels;         /* ICC/ILU(levels) */
1011bcd9e38bSBarry Smith   PetscReal     pivotinblocks;  /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0
1012bcd9e38bSBarry Smith                                    factorization may be faster if do not pivot */
101315e8a5b3SHong Zhang   PetscReal     zeropivot;      /* pivot is called zero if less than this */
1014f4db908eSBarry Smith   PetscReal     shifttype;      /* type of shift added to matrix factor to prevent zero pivots */
1015f4db908eSBarry Smith   PetscReal     shiftamount;     /* how large the shift is */
101615e8a5b3SHong Zhang } MatFactorInfo;
1017ffa6d0a5SLois Curfman McInnes 
1018014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFactorInfoInitialize(MatFactorInfo*);
1019014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCholeskyFactor(Mat,IS,const MatFactorInfo*);
1020014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCholeskyFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
1021014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCholeskyFactorNumeric(Mat,Mat,const MatFactorInfo*);
1022014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLUFactor(Mat,IS,IS,const MatFactorInfo*);
1023014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatILUFactor(Mat,IS,IS,const MatFactorInfo*);
1024014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
1025014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatILUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
1026014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatICCFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
1027014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatICCFactor(Mat,IS,const MatFactorInfo*);
1028014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLUFactorNumeric(Mat,Mat,const MatFactorInfo*);
1029014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*);
1030014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolve(Mat,Vec,Vec);
1031014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatForwardSolve(Mat,Vec,Vec);
1032014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatBackwardSolve(Mat,Vec,Vec);
1033014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveAdd(Mat,Vec,Vec,Vec);
1034014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveTranspose(Mat,Vec,Vec);
1035014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveTransposeAdd(Mat,Vec,Vec,Vec);
1036014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolves(Mat,Vecs,Vecs);
10378ed539a5SBarry Smith 
1038014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetUnfactored(Mat);
1039bb5a7306SBarry Smith 
1040d91e6319SBarry Smith /*E
1041d91e6319SBarry Smith     MatSORType - What type of (S)SOR to perform
1042bb1eb677SSatish Balay 
1043d91e6319SBarry Smith     Level: beginner
1044d91e6319SBarry Smith 
1045d9274352SBarry Smith    May be bitwise ORd together
1046d9274352SBarry Smith 
1047d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
1048d91e6319SBarry Smith 
10494e7234bfSBarry Smith    MatSORType may be bitwise ORd together, so do not change the numbers
10504e7234bfSBarry Smith 
105141f059aeSBarry Smith .seealso: MatSOR()
1052d91e6319SBarry Smith E*/
1053ee50ffe9SBarry Smith typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3,
1054ee50ffe9SBarry Smith               SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8,
1055ee50ffe9SBarry Smith               SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16,
105684cb2905SBarry Smith               SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType;
1057014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSOR(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
10588ed539a5SBarry Smith 
1059d4fbbf0eSBarry Smith /*
1060639f9d9dSBarry Smith     These routines are for efficiently computing Jacobians via finite differences.
1061639f9d9dSBarry Smith */
1062b12f92e5SBarry Smith 
10637086a01eSPeter Brune /*S
10647086a01eSPeter Brune      MatColoring - Object for managing the coloring of matrices.
10657086a01eSPeter Brune 
10667086a01eSPeter Brune    Level: beginner
10677086a01eSPeter Brune 
10687086a01eSPeter Brune   Concepts: matrix, coloring
10697086a01eSPeter Brune 
10707086a01eSPeter Brune .seealso:  MatFDColoringCreate() ISColoring MatFDColoring
10717086a01eSPeter Brune S*/
10727086a01eSPeter Brune typedef struct _p_MatColoring* MatColoring;
107376bdecfbSBarry Smith /*J
10748f6c3df8SBarry Smith     MatColoringType - String with the name of a PETSc matrix coloring
1075d9274352SBarry Smith 
1076d9274352SBarry Smith    Level: beginner
1077d9274352SBarry Smith 
10787086a01eSPeter Brune .seealso: MatColoringSetType(), MatColoring
107976bdecfbSBarry Smith J*/
10807086a01eSPeter Brune 
108119fd82e9SBarry Smith typedef const  char*           MatColoringType;
10825567d71aSPeter Brune #define MATCOLORINGJP      "jp"
1083b9acb617SPeter Brune #define MATCOLORINGPOWER   "power"
10842692d6eeSBarry Smith #define MATCOLORINGNATURAL "natural"
10852692d6eeSBarry Smith #define MATCOLORINGSL      "sl"
10862692d6eeSBarry Smith #define MATCOLORINGLF      "lf"
10872692d6eeSBarry Smith #define MATCOLORINGID      "id"
10888121bdceSPeter Brune #define MATCOLORINGGREEDY  "greedy"
1089b12f92e5SBarry Smith 
10908ac6da0aSPeter Brune /*E
10918ac6da0aSPeter Brune    MatColoringWeightType - Type of weight scheme
10928ac6da0aSPeter Brune 
10938ac6da0aSPeter Brune     Not Collective
10948ac6da0aSPeter Brune 
10958ac6da0aSPeter Brune +   MAT_COLORING_RANDOM  - Random weights
10968ac6da0aSPeter Brune .   MAT_COLORING_LEXICAL - Lexical weighting based upon global numbering.
10978ac6da0aSPeter Brune -   MAT_COLORING_LF      - Last-first weighting.
10988ac6da0aSPeter Brune 
10998ac6da0aSPeter Brune     Level: intermediate
11008ac6da0aSPeter Brune 
11018ac6da0aSPeter Brune    Any additions/changes here MUST also be made in include/finclude/petscmat.h
11028ac6da0aSPeter Brune 
11038ac6da0aSPeter Brune .seealso: MatCUSPSetFormat(), MatCUSPFormatOperation
11048ac6da0aSPeter Brune E*/
1105301c140bSPeter Brune typedef enum {MAT_COLORING_WEIGHT_RANDOM,MAT_COLORING_WEIGHT_LEXICAL,MAT_COLORING_WEIGHT_LF,MAT_COLORING_WEIGHT_SL} MatColoringWeightType;
11068ac6da0aSPeter Brune 
1107335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringCreate(Mat,MatColoring*);
1108d7f2a307SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringGetDegrees(Mat,PetscInt,PetscInt*);
1109335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringDestroy(MatColoring*);
1110335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringView(MatColoring,PetscViewer);
1111335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetType(MatColoring,MatColoringType);
1112335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetFromOptions(MatColoring);
1113335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetDistance(MatColoring,PetscInt);
1114335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringGetDistance(MatColoring,PetscInt*);
1115335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetMaxColors(MatColoring,PetscInt);
1116335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringGetMaxColors(MatColoring,PetscInt*);
1117335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringApply(MatColoring,ISColoring*);
1118607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatColoringRegisterAll(void);
1119335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringRegister(const char[],PetscErrorCode(*)(MatColoring));
1120335efc43SPeter Brune PETSC_EXTERN PetscBool MatColoringRegisterAllCalled;
1121014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*);
11228ac6da0aSPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetWeightType(MatColoring,MatColoringWeightType);
1123c29971fcSPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetWeights(MatColoring,PetscReal*,PetscInt*);
11248ac6da0aSPeter Brune PETSC_EXTERN PetscErrorCode MatColoringCreateWeights(MatColoring,PetscReal **,PetscInt **lperm);
1125639f9d9dSBarry Smith 
1126d9274352SBarry Smith /*S
1127d9274352SBarry Smith      MatFDColoring - Object for computing a sparse Jacobian via finite differences
1128d9274352SBarry Smith         and coloring
1129639f9d9dSBarry Smith 
1130d9274352SBarry Smith    Level: beginner
1131d9274352SBarry Smith 
1132d9274352SBarry Smith   Concepts: coloring, sparse Jacobian, finite differences
1133d9274352SBarry Smith 
1134d9274352SBarry Smith .seealso:  MatFDColoringCreate()
1135d9274352SBarry Smith S*/
1136e2a1c21fSSatish Balay typedef struct _p_MatFDColoring* MatFDColoring;
1137639f9d9dSBarry Smith 
1138014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringCreate(Mat,ISColoring,MatFDColoring *);
1139014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringDestroy(MatFDColoring*);
1140014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringView(MatFDColoring,PetscViewer);
1141014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*);
1142014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**);
1143014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal);
1144014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetFromOptions(MatFDColoring);
1145d1e9a80fSBarry Smith PETSC_EXTERN PetscErrorCode MatFDColoringApply(Mat,MatFDColoring,Vec,void *);
1146014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetF(MatFDColoring,Vec);
1147014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,PetscInt*[]);
1148f86b9fbaSHong Zhang PETSC_EXTERN PetscErrorCode MatFDColoringSetUp(Mat,ISColoring,MatFDColoring);
1149f86b9fbaSHong Zhang PETSC_EXTERN PetscErrorCode MatFDColoringSetBlockSize(MatFDColoring,PetscInt,PetscInt);
1150f86b9fbaSHong Zhang 
1151b1683b59SHong Zhang 
1152b1683b59SHong Zhang /*S
1153b9af6bddSHong Zhang      MatTransposeColoring - Object for computing a sparse matrix product C=A*B^T via coloring
1154b1683b59SHong Zhang 
1155b1683b59SHong Zhang    Level: beginner
1156b1683b59SHong Zhang 
1157b1683b59SHong Zhang   Concepts: coloring, sparse matrix product
1158b1683b59SHong Zhang 
1159b9af6bddSHong Zhang .seealso:  MatTransposeColoringCreate()
1160b1683b59SHong Zhang S*/
1161b9af6bddSHong Zhang typedef struct _p_MatTransposeColoring* MatTransposeColoring;
1162b1683b59SHong Zhang 
1163014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeColoringCreate(Mat,ISColoring,MatTransposeColoring *);
1164014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransColoringApplySpToDen(MatTransposeColoring,Mat,Mat);
1165014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransColoringApplyDenToSp(MatTransposeColoring,Mat,Mat);
1166014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeColoringDestroy(MatTransposeColoring*);
1167b1683b59SHong Zhang 
1168639f9d9dSBarry Smith /*
11690752156aSBarry Smith     These routines are for partitioning matrices: currently used only
11703eda8832SBarry Smith   for adjacency matrix, MatCreateMPIAdj().
11710752156aSBarry Smith */
1172ca161407SBarry Smith 
1173d9274352SBarry Smith /*S
1174d9274352SBarry Smith      MatPartitioning - Object for managing the partitioning of a matrix or graph
1175d9274352SBarry Smith 
1176d9274352SBarry Smith    Level: beginner
1177d9274352SBarry Smith 
1178d9274352SBarry Smith   Concepts: partitioning
1179d9274352SBarry Smith 
1180743a1026Svictorle .seealso:  MatPartitioningCreate(), MatPartitioningType
1181d9274352SBarry Smith S*/
118291e9ee9fSBarry Smith typedef struct _p_MatPartitioning* MatPartitioning;
1183d9274352SBarry Smith 
118476bdecfbSBarry Smith /*J
11858f6c3df8SBarry Smith     MatPartitioningType - String with the name of a PETSc matrix partitioning
1186d9274352SBarry Smith 
1187d9274352SBarry Smith    Level: beginner
11880d04baf8SBarry Smith dm
1189b547a13bSvictorle .seealso: MatPartitioningCreate(), MatPartitioning
119076bdecfbSBarry Smith J*/
119119fd82e9SBarry Smith typedef const char* MatPartitioningType;
11922692d6eeSBarry Smith #define MATPARTITIONINGCURRENT  "current"
11932692d6eeSBarry Smith #define MATPARTITIONINGSQUARE   "square"
11942692d6eeSBarry Smith #define MATPARTITIONINGPARMETIS "parmetis"
11952692d6eeSBarry Smith #define MATPARTITIONINGCHACO    "chaco"
11962692d6eeSBarry Smith #define MATPARTITIONINGPARTY    "party"
119750d91057SBarry Smith #define MATPARTITIONINGPTSCOTCH "ptscotch"
119871306c60Sjroman 
1199ca161407SBarry Smith 
1200014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningCreate(MPI_Comm,MatPartitioning*);
120119fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningSetType(MatPartitioning,MatPartitioningType);
1202014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetNParts(MatPartitioning,PetscInt);
1203014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetAdjacency(MatPartitioning,Mat);
1204014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]);
1205014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []);
1206014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningApply(MatPartitioning,IS*);
1207014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningDestroy(MatPartitioning*);
12082aabb6bbSBarry Smith 
1209bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningRegister(const char[],PetscErrorCode (*)(MatPartitioning));
121030de9b25SBarry Smith 
1211014dd563SJed Brown PETSC_EXTERN PetscBool MatPartitioningRegisterAllCalled;
1212f1144a30SSatish Balay 
1213607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningRegisterAll(void);
12142bad1931SBarry Smith 
1215014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningView(MatPartitioning,PetscViewer);
1216014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetFromOptions(MatPartitioning);
121719fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningGetType(MatPartitioning,MatPartitioningType*);
1218ca161407SBarry Smith 
1219014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetCoarseSequential(MatPartitioning);
1220014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningParmetisGetEdgeCut(MatPartitioning, PetscInt *);
12210752156aSBarry Smith 
1222b6956c12SJose Roman typedef enum { MP_CHACO_MULTILEVEL=1,MP_CHACO_SPECTRAL=2,MP_CHACO_LINEAR=4,MP_CHACO_RANDOM=5,MP_CHACO_SCATTERED=6 } MPChacoGlobalType;
12236a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoGlobalTypes[];
1224b6956c12SJose Roman typedef enum { MP_CHACO_KERNIGHAN=1,MP_CHACO_NONE=2 } MPChacoLocalType;
12256a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoLocalTypes[];
1226b6956c12SJose Roman typedef enum { MP_CHACO_LANCZOS=0,MP_CHACO_RQI=1 } MPChacoEigenType;
12276a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoEigenTypes[];
1228b6956c12SJose Roman 
1229014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetGlobal(MatPartitioning,MPChacoGlobalType);
1230014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetGlobal(MatPartitioning,MPChacoGlobalType*);
1231014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetLocal(MatPartitioning,MPChacoLocalType);
1232014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetLocal(MatPartitioning,MPChacoLocalType*);
1233014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal);
1234014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType);
1235014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenSolver(MatPartitioning,MPChacoEigenType*);
1236014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenTol(MatPartitioning,PetscReal);
1237014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenTol(MatPartitioning,PetscReal*);
1238014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenNumber(MatPartitioning,PetscInt);
1239014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenNumber(MatPartitioning,PetscInt*);
124071306c60Sjroman 
124171306c60Sjroman #define MP_PARTY_OPT "opt"
124271306c60Sjroman #define MP_PARTY_LIN "lin"
124371306c60Sjroman #define MP_PARTY_SCA "sca"
124471306c60Sjroman #define MP_PARTY_RAN "ran"
124571306c60Sjroman #define MP_PARTY_GBF "gbf"
124671306c60Sjroman #define MP_PARTY_GCF "gcf"
124771306c60Sjroman #define MP_PARTY_BUB "bub"
124871306c60Sjroman #define MP_PARTY_DEF "def"
1249014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetGlobal(MatPartitioning,const char*);
125071306c60Sjroman #define MP_PARTY_HELPFUL_SETS "hs"
125171306c60Sjroman #define MP_PARTY_KERNIGHAN_LIN "kl"
125271306c60Sjroman #define MP_PARTY_NONE "no"
1253014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetLocal(MatPartitioning,const char*);
1254014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal);
1255014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetBipart(MatPartitioning,PetscBool);
1256014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscBool);
125771306c60Sjroman 
125850d91057SBarry Smith typedef enum { MP_PTSCOTCH_QUALITY,MP_PTSCOTCH_SPEED,MP_PTSCOTCH_BALANCE,MP_PTSCOTCH_SAFETY,MP_PTSCOTCH_SCALABILITY } MPPTScotchStrategyType;
12596a6fc655SJed Brown PETSC_EXTERN const char *const MPPTScotchStrategyTypes[];
1260e0f1cffaSJose Roman 
1261014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetImbalance(MatPartitioning,PetscReal);
1262014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetImbalance(MatPartitioning,PetscReal*);
1263014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetStrategy(MatPartitioning,MPPTScotchStrategyType);
1264014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetStrategy(MatPartitioning,MPPTScotchStrategyType*);
126571306c60Sjroman 
1266b43b03e9SMark F. Adams /*
1267b43b03e9SMark F. Adams     These routines are for coarsening matrices:
1268b43b03e9SMark F. Adams */
1269b43b03e9SMark F. Adams 
1270b43b03e9SMark F. Adams /*S
1271b43b03e9SMark F. Adams      MatCoarsen - Object for managing the coarsening of a graph (symmetric matrix)
1272b43b03e9SMark F. Adams 
1273b43b03e9SMark F. Adams    Level: beginner
1274b43b03e9SMark F. Adams 
1275b43b03e9SMark F. Adams   Concepts: coarsen
1276b43b03e9SMark F. Adams 
1277b43b03e9SMark F. Adams .seealso:  MatCoarsenCreate), MatCoarsenType
1278b43b03e9SMark F. Adams S*/
1279b43b03e9SMark F. Adams typedef struct _p_MatCoarsen* MatCoarsen;
1280b43b03e9SMark F. Adams 
1281b43b03e9SMark F. Adams /*J
12828f6c3df8SBarry Smith     MatCoarsenType - String with the name of a PETSc matrix coarsen
1283b43b03e9SMark F. Adams 
1284b43b03e9SMark F. Adams    Level: beginner
12858f6c3df8SBarry Smith 
1286b43b03e9SMark F. Adams .seealso: MatCoarsenCreate(), MatCoarsen
1287b43b03e9SMark F. Adams J*/
128819fd82e9SBarry Smith typedef const char* MatCoarsenType;
1289b43b03e9SMark F. Adams #define MATCOARSENMIS  "mis"
12909057884aSMark F. Adams #define MATCOARSENHEM  "hem"
1291b43b03e9SMark F. Adams 
12920cbbd2e1SMark F. Adams /* linked list for aggregates */
129341b27cdeSMark F. Adams typedef struct _PetscCDIntNd{
129441b27cdeSMark F. Adams   struct _PetscCDIntNd *next;
12950cbbd2e1SMark F. Adams   PetscInt             gid;
129641b27cdeSMark F. Adams }PetscCDIntNd;
12970cbbd2e1SMark F. Adams 
12980cbbd2e1SMark F. Adams /* only used by node pool */
129941b27cdeSMark F. Adams typedef struct _PetscCDArrNd{
130041b27cdeSMark F. Adams   struct _PetscCDArrNd *next;
130141b27cdeSMark F. Adams   struct _PetscCDIntNd *array;
130241b27cdeSMark F. Adams }PetscCDArrNd;
13030cbbd2e1SMark F. Adams 
13040cbbd2e1SMark F. Adams typedef struct _PetscCoarsenData{
130582c86c8fSBarry Smith   PetscCDArrNd pool_list;  /* node pool */
130641b27cdeSMark F. Adams   PetscCDIntNd *new_node;
13070cbbd2e1SMark F. Adams   PetscInt     new_left;
13080cbbd2e1SMark F. Adams   PetscInt     chk_sz;
130941b27cdeSMark F. Adams   PetscCDIntNd *extra_nodes;
131082c86c8fSBarry Smith   PetscCDIntNd **array;  /* Array of lists */
13110cbbd2e1SMark F. Adams   PetscInt     size;
131282c86c8fSBarry Smith   Mat          mat;  /* cache a Mat for communication data */
13130cbbd2e1SMark F. Adams }PetscCoarsenData;
13140cbbd2e1SMark F. Adams 
1315014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenCreate(MPI_Comm,MatCoarsen*);
131619fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenSetType(MatCoarsen,MatCoarsenType);
1317014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetAdjacency(MatCoarsen,Mat);
1318014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetGreedyOrdering(MatCoarsen,const IS);
1319014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetStrictAggs(MatCoarsen,PetscBool);
1320014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetVerbose(MatCoarsen,PetscInt);
1321014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenGetData( MatCoarsen, PetscCoarsenData ** );
1322014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenApply(MatCoarsen);
1323014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenDestroy(MatCoarsen*);
1324b43b03e9SMark F. Adams 
1325bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenRegister(const char[],PetscErrorCode (*)(MatCoarsen));
1326b43b03e9SMark F. Adams 
1327014dd563SJed Brown PETSC_EXTERN PetscBool MatCoarsenRegisterAllCalled;
1328b43b03e9SMark F. Adams 
1329607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenRegisterAll(void);
1330b43b03e9SMark F. Adams 
1331014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenView(MatCoarsen,PetscViewer);
1332014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetFromOptions(MatCoarsen);
133319fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenGetType(MatCoarsen,MatCoarsenType*);
1334ce1779c8SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatCoarsenViewFromOptions(MatCoarsen A,const char prefix[],const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,prefix,name);}
1335b43b03e9SMark F. Adams 
1336014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMeshToVertexGraph(Mat,PetscInt,Mat*);
1337014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMeshToCellGraph(Mat,PetscInt,Mat*);
1338591294e4SBarry Smith 
13390752156aSBarry Smith /*
13400a835dfdSSatish Balay     If you add entries here you must also add them to finclude/petscmat.h
1341d4fbbf0eSBarry Smith */
13421c1c02c0SLois Curfman McInnes typedef enum { MATOP_SET_VALUES=0,
13431c1c02c0SLois Curfman McInnes                MATOP_GET_ROW=1,
13441c1c02c0SLois Curfman McInnes                MATOP_RESTORE_ROW=2,
13451c1c02c0SLois Curfman McInnes                MATOP_MULT=3,
13461c1c02c0SLois Curfman McInnes                MATOP_MULT_ADD=4,
13477c922b88SBarry Smith                MATOP_MULT_TRANSPOSE=5,
13487c922b88SBarry Smith                MATOP_MULT_TRANSPOSE_ADD=6,
13491c1c02c0SLois Curfman McInnes                MATOP_SOLVE=7,
13501c1c02c0SLois Curfman McInnes                MATOP_SOLVE_ADD=8,
13517c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE=9,
13527c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE_ADD=10,
13531c1c02c0SLois Curfman McInnes                MATOP_LUFACTOR=11,
13541c1c02c0SLois Curfman McInnes                MATOP_CHOLESKYFACTOR=12,
1355a714c06dSBarry Smith                MATOP_SOR=13,
13561c1c02c0SLois Curfman McInnes                MATOP_TRANSPOSE=14,
13571c1c02c0SLois Curfman McInnes                MATOP_GETINFO=15,
13581c1c02c0SLois Curfman McInnes                MATOP_EQUAL=16,
13591c1c02c0SLois Curfman McInnes                MATOP_GET_DIAGONAL=17,
13601c1c02c0SLois Curfman McInnes                MATOP_DIAGONAL_SCALE=18,
13611c1c02c0SLois Curfman McInnes                MATOP_NORM=19,
13621c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_BEGIN=20,
13631c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_END=21,
1364d519adbfSMatthew Knepley                MATOP_SET_OPTION=22,
1365d519adbfSMatthew Knepley                MATOP_ZERO_ENTRIES=23,
1366d519adbfSMatthew Knepley                MATOP_ZERO_ROWS=24,
1367d519adbfSMatthew Knepley                MATOP_LUFACTOR_SYMBOLIC=25,
1368d519adbfSMatthew Knepley                MATOP_LUFACTOR_NUMERIC=26,
1369d519adbfSMatthew Knepley                MATOP_CHOLESKY_FACTOR_SYMBOLIC=27,
1370d519adbfSMatthew Knepley                MATOP_CHOLESKY_FACTOR_NUMERIC=28,
1371d519adbfSMatthew Knepley                MATOP_SETUP_PREALLOCATION=29,
1372d519adbfSMatthew Knepley                MATOP_ILUFACTOR_SYMBOLIC=30,
1373d519adbfSMatthew Knepley                MATOP_ICCFACTOR_SYMBOLIC=31,
13749d39f69dSJed Brown                /* MATOP_PLACEHOLDER_32=32, */
13759d39f69dSJed Brown                /* MATOP_PLACEHOLDER_33=33, */
1376d519adbfSMatthew Knepley                MATOP_DUPLICATE=34,
1377d519adbfSMatthew Knepley                MATOP_FORWARD_SOLVE=35,
1378d519adbfSMatthew Knepley                MATOP_BACKWARD_SOLVE=36,
1379d519adbfSMatthew Knepley                MATOP_ILUFACTOR=37,
1380d519adbfSMatthew Knepley                MATOP_ICCFACTOR=38,
1381d519adbfSMatthew Knepley                MATOP_AXPY=39,
1382d519adbfSMatthew Knepley                MATOP_GET_SUBMATRICES=40,
1383d519adbfSMatthew Knepley                MATOP_INCREASE_OVERLAP=41,
1384d519adbfSMatthew Knepley                MATOP_GET_VALUES=42,
1385d519adbfSMatthew Knepley                MATOP_COPY=43,
1386d519adbfSMatthew Knepley                MATOP_GET_ROW_MAX=44,
1387d519adbfSMatthew Knepley                MATOP_SCALE=45,
1388d519adbfSMatthew Knepley                MATOP_SHIFT=46,
138935153367SBarry Smith                MATOP_DIAGONAL_SET=47,
13909d39f69dSJed Brown                MATOP_ZERO_ROWS_COLUMNS=48,
13919d39f69dSJed Brown                MATOP_SET_RANDOM=49,
1392d519adbfSMatthew Knepley                MATOP_GET_ROW_IJ=50,
1393d519adbfSMatthew Knepley                MATOP_RESTORE_ROW_IJ=51,
1394d519adbfSMatthew Knepley                MATOP_GET_COLUMN_IJ=52,
1395d519adbfSMatthew Knepley                MATOP_RESTORE_COLUMN_IJ=53,
1396d519adbfSMatthew Knepley                MATOP_FDCOLORING_CREATE=54,
1397d519adbfSMatthew Knepley                MATOP_COLORING_PATCH=55,
1398d519adbfSMatthew Knepley                MATOP_SET_UNFACTORED=56,
1399d519adbfSMatthew Knepley                MATOP_PERMUTE=57,
1400d519adbfSMatthew Knepley                MATOP_SET_VALUES_BLOCKED=58,
1401d519adbfSMatthew Knepley                MATOP_GET_SUBMATRIX=59,
1402d519adbfSMatthew Knepley                MATOP_DESTROY=60,
1403d519adbfSMatthew Knepley                MATOP_VIEW=61,
1404d519adbfSMatthew Knepley                MATOP_CONVERT_FROM=62,
14057bab7c10SHong Zhang                MATOP_MATMAT_MULT=63,
14067bab7c10SHong Zhang                MATOP_MATMAT_MULT_SYMBOLIC=64,
14077bab7c10SHong Zhang                MATOP_MATMAT_MULT_NUMERIC=65,
1408d519adbfSMatthew Knepley                MATOP_SET_LOCAL_TO_GLOBAL_MAP=66,
1409d519adbfSMatthew Knepley                MATOP_SET_VALUES_LOCAL=67,
1410d519adbfSMatthew Knepley                MATOP_ZERO_ROWS_LOCAL=68,
1411d519adbfSMatthew Knepley                MATOP_GET_ROW_MAX_ABS=69,
1412d519adbfSMatthew Knepley                MATOP_GET_ROW_MIN_ABS=70,
1413d519adbfSMatthew Knepley                MATOP_CONVERT=71,
1414d519adbfSMatthew Knepley                MATOP_SET_COLORING=72,
14159d39f69dSJed Brown                /* MATOP_PLACEHOLDER_73=73, */
1416d519adbfSMatthew Knepley                MATOP_SET_VALUES_ADIFOR=74,
1417d519adbfSMatthew Knepley                MATOP_FD_COLORING_APPLY=75,
1418d519adbfSMatthew Knepley                MATOP_SET_FROM_OPTIONS=76,
1419bda6c4cbSBarry Smith                MATOP_MULT_CONSTRAINED=77,
1420bda6c4cbSBarry Smith                MATOP_MULT_TRANSPOSE_CONSTRAIN=78,
14219d39f69dSJed Brown                MATOP_FIND_ZERO_DIAGONALS=79,
1422d519adbfSMatthew Knepley                MATOP_MULT_MULTIPLE=80,
1423d519adbfSMatthew Knepley                MATOP_SOLVE_MULTIPLE=81,
1424d519adbfSMatthew Knepley                MATOP_GET_INERTIA=82,
1425d519adbfSMatthew Knepley                MATOP_LOAD=83,
1426d519adbfSMatthew Knepley                MATOP_IS_SYMMETRIC=84,
1427d519adbfSMatthew Knepley                MATOP_IS_HERMITIAN=85,
1428d519adbfSMatthew Knepley                MATOP_IS_STRUCTURALLY_SYMMETRIC=86,
1429820469bcSHong Zhang                MATOP_SET_VALUES_BLOCKEDLOCAL=87,
1430d519adbfSMatthew Knepley                MATOP_GET_VECS=88,
1431d519adbfSMatthew Knepley                MATOP_MAT_MULT=89,
1432d519adbfSMatthew Knepley                MATOP_MAT_MULT_SYMBOLIC=90,
1433d519adbfSMatthew Knepley                MATOP_MAT_MULT_NUMERIC=91,
1434d519adbfSMatthew Knepley                MATOP_PTAP=92,
1435d519adbfSMatthew Knepley                MATOP_PTAP_SYMBOLIC=93,
1436d519adbfSMatthew Knepley                MATOP_PTAP_NUMERIC=94,
1437bda6c4cbSBarry Smith                MATOP_MAT_TRANSPOSE_MULT=95,
1438bda6c4cbSBarry Smith                MATOP_MAT_TRANSPOSE_MULT_SYMBO=96,
1439bda6c4cbSBarry Smith                MATOP_MAT_TRANSPOSE_MULT_NUMER=97,
14409d39f69dSJed Brown                /* MATOP_PLACEHOLDER_98=98, */
14419d39f69dSJed Brown                /* MATOP_PLACEHOLDER_99=99, */
14429d39f69dSJed Brown                /* MATOP_PLACEHOLDER_100=100, */
14439d39f69dSJed Brown                /* MATOP_PLACEHOLDER_101=101, */
1444d519adbfSMatthew Knepley                MATOP_CONJUGATE=102,
14459d39f69dSJed Brown                /* MATOP_PLACEHOLDER_103=103, */
1446d519adbfSMatthew Knepley                MATOP_SET_VALUES_ROW=104,
1447d519adbfSMatthew Knepley                MATOP_REAL_PART=105,
1448bda6c4cbSBarry Smith                MATOP_IMAGINARY_PART=106,
1449bda6c4cbSBarry Smith                MATOP_GET_ROW_UPPER_TRIANGULAR=107,
1450bda6c4cbSBarry Smith                MATOP_RESTORE_ROW_UPPER_TRIANG=108,
1451bda6c4cbSBarry Smith                MATOP_MAT_SOLVE=109,
14520e8d04f7SBarry Smith                MATOP_GET_REDUNDANT_MATRIX=110,
1453d519adbfSMatthew Knepley                MATOP_GET_ROW_MIN=111,
14540e8d04f7SBarry Smith                MATOP_GET_COLUMN_VECTOR=112,
1455d519adbfSMatthew Knepley                MATOP_MISSING_DIAGONAL=113,
14560e8d04f7SBarry Smith                MATOP_GET_SEQ_NONZERO_STRUCTUR=114,
145789c6957cSBarry Smith                MATOP_CREATE=115,
145889c6957cSBarry Smith                MATOP_GET_GHOSTS=116,
14590e8d04f7SBarry Smith                MATOP_GET_LOCAL_SUB_MATRIX=117,
14600e8d04f7SBarry Smith                MATOP_RESTORE_LOCALSUB_MATRIX=118,
1461eeffb40dSHong Zhang                MATOP_MULT_DIAGONAL_BLOCK=119,
14620e8d04f7SBarry Smith                MATOP_HERMITIAN_TRANSPOSE=120,
14630e8d04f7SBarry Smith                MATOP_MULT_HERMITIAN_TRANSPOSE=121,
14640e8d04f7SBarry Smith                MATOP_MULT_HERMITIAN_TRANS_ADD=122,
14650e8d04f7SBarry Smith                MATOP_GET_MULTI_PROC_BLOCK=123,
14669d39f69dSJed Brown                MATOP_FIND_NONZERO_ROWS=124,
14670e8d04f7SBarry Smith                MATOP_GET_COLUMN_NORMS=125,
14689d39f69dSJed Brown                MATOP_INVERT_BLOCK_DIAGONAL=126,
14699d39f69dSJed Brown                /* MATOP_PLACEHOLDER_127=127, */
14700e8d04f7SBarry Smith                MATOP_GET_SUB_MATRICES_PARALLE=128,
14712b8ad9a3SHong Zhang                MATOP_SET_VALUES_BATCH=129,
14720e8d04f7SBarry Smith                MATOP_TRANSPOSE_MAT_MULT=130,
14730e8d04f7SBarry Smith                MATOP_TRANSPOSE_MAT_MULT_SYMBO=131,
14740e8d04f7SBarry Smith                MATOP_TRANSPOSE_MAT_MULT_NUMER=132,
14750e8d04f7SBarry Smith                MATOP_TRANSPOSE_COLORING_CREAT=133,
14760e8d04f7SBarry Smith                MATOP_TRANS_COLORING_APPLY_SPT=134,
14770e8d04f7SBarry Smith                MATOP_TRANS_COLORING_APPLY_DEN=135,
14780e8d04f7SBarry Smith                MATOP_RART=136,
14790e8d04f7SBarry Smith                MATOP_RART_SYMBOLIC=137,
14800e8d04f7SBarry Smith                MATOP_RART_NUMERIC=138,
1481e09a3074SHong Zhang                MATOP_SET_BLOCK_SIZES=139,
1482f9426fe0SMark Adams                MATOP_AYPX=140,
14831919a2e2SJed Brown                MATOP_RESIDUAL=141,
14841919a2e2SJed Brown                MATOP_FDCOLORING_SETUP= 142
1485fae171e0SBarry Smith              } MatOperation;
1486014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatHasOperation(Mat,MatOperation,PetscBool *);
1487014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellSetOperation(Mat,MatOperation,void(*)(void));
1488014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellGetOperation(Mat,MatOperation,void(**)(void));
1489014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellSetContext(Mat,void*);
1490112a2221SBarry Smith 
149190ace30eSBarry Smith /*
149290ace30eSBarry Smith    Codes for matrices stored on disk. By default they are
149390ace30eSBarry Smith    stored in a universal format. By changing the format with
14947973ad04SBarry Smith    PetscViewerSetFormat(viewer,PETSC_VIEWER_NATIVE); the matrices will
149590ace30eSBarry Smith    be stored in a way natural for the matrix, for example dense matrices
149690ace30eSBarry Smith    would be stored as dense. Matrices stored this way may only be
1497f4403165SShri Abhyankar    read into matrices of the same type.
149890ace30eSBarry Smith */
149990ace30eSBarry Smith #define MATRIX_BINARY_FORMAT_DENSE -1
150090ace30eSBarry Smith 
1501014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetHashTableFactor(Mat,PetscReal);
1502014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatISGetLocalMat(Mat,Mat*);
1503014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatISSetLocalMat(Mat,Mat);
15041f4e1ec7SBarry Smith 
1505d9274352SBarry Smith /*S
1506d9274352SBarry Smith      MatNullSpace - Object that removes a null space from a vector, i.e.
1507d9274352SBarry Smith          orthogonalizes the vector to a subsapce
1508d9274352SBarry Smith 
1509f7a9e4ceSBarry Smith    Level: advanced
1510d9274352SBarry Smith 
1511d9274352SBarry Smith   Concepts: matrix; linear operator, null space
1512d9274352SBarry Smith 
15136e1639daSBarry Smith   Users manual sections:
15146e1639daSBarry Smith .   sec_singular
15156e1639daSBarry Smith 
1516d9274352SBarry Smith .seealso:  MatNullSpaceCreate()
1517d9274352SBarry Smith S*/
151874637425SBarry Smith typedef struct _p_MatNullSpace* MatNullSpace;
1519d9274352SBarry Smith 
1520014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceCreate(MPI_Comm,PetscBool ,PetscInt,const Vec[],MatNullSpace*);
1521014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(MatNullSpace,Vec,void*),void*);
1522014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceDestroy(MatNullSpace*);
1523d0195637SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceRemove(MatNullSpace,Vec);
1524014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNullSpace(Mat, MatNullSpace *);
1525014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNullSpace(Mat,MatNullSpace);
1526014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNearNullSpace(Mat,MatNullSpace);
1527014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNearNullSpace(Mat,MatNullSpace*);
1528014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceTest(MatNullSpace,Mat,PetscBool  *);
1529014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceView(MatNullSpace,PetscViewer);
1530014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceGetVecs(MatNullSpace,PetscBool*,PetscInt*,const Vec**);
1531014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceCreateRigidBody(Vec,MatNullSpace*);
153274637425SBarry Smith 
1533014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatReorderingSeqSBAIJ(Mat,IS);
1534014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetHashTableFactor(Mat,PetscReal);
1535014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetColumnIndices(Mat,PetscInt *);
1536014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJInvertBlockDiagonal(Mat);
15373f1d51d7SBarry Smith 
1538014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMAIJ(Mat,PetscInt,Mat*);
1539014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMAIJRedimension(Mat,PetscInt,Mat*);
1540014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMAIJGetAIJ(Mat,Mat*);
1541c4f061fbSSatish Balay 
1542014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatComputeExplicitOperator(Mat,Mat*);
1543b0a32e0cSBarry Smith 
1544014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalScaleLocal(Mat,Vec);
154504f1ad80SBarry Smith 
1546014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*);
1547014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetBase(Mat,Vec,Vec);
1548014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*);
1549014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*));
1550014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec));
1551014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDAddNullSpace(Mat,MatNullSpace);
1552014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt);
1553014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDResetHHistory(Mat);
1554014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctionError(Mat,PetscReal);
1555014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetPeriod(Mat,PetscInt);
1556014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDGetH(Mat,PetscScalar *);
1557014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetOptionsPrefix(Mat,const char[]);
1558014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*);
1559014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*);
1560e884886fSBarry Smith 
15616370053bSBarry Smith /*S
15626370053bSBarry Smith     MatMFFD - A data structured used to manage the computation of the h differencing parameter for matrix-free
15636370053bSBarry Smith               Jacobian vector products
1564e884886fSBarry Smith 
15656370053bSBarry Smith     Notes: MATMFFD is a specific MatType which uses the MatMFFD data structure
15666370053bSBarry Smith 
15676370053bSBarry Smith            MatMFFD*() methods actually take the Mat as their first argument. Not a MatMFFD data structure
15686370053bSBarry Smith 
15696370053bSBarry Smith     Level: developer
15706370053bSBarry Smith 
15716370053bSBarry Smith .seealso: MATMFFD, MatCreateMFFD(), MatMFFDSetFuction(), MatMFFDSetType(), MatMFFDRegister()
15726370053bSBarry Smith S*/
1573e884886fSBarry Smith typedef struct _p_MatMFFD* MatMFFD;
1574e884886fSBarry Smith 
157576bdecfbSBarry Smith /*J
1576e884886fSBarry Smith     MatMFFDType - algorithm used to compute the h used in computing matrix-vector products via differencing of the function
1577e884886fSBarry Smith 
1578e884886fSBarry Smith    Level: beginner
1579e884886fSBarry Smith 
1580e884886fSBarry Smith .seealso: MatMFFDSetType(), MatMFFDRegister()
158176bdecfbSBarry Smith J*/
158219fd82e9SBarry Smith typedef const char* MatMFFDType;
1583e884886fSBarry Smith #define MATMFFD_DS  "ds"
1584e884886fSBarry Smith #define MATMFFD_WP  "wp"
1585e884886fSBarry Smith 
158619fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDSetType(Mat,MatMFFDType);
1587bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDRegister(const char[],PetscErrorCode (*)(MatMFFD));
1588e884886fSBarry Smith 
1589607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDRegisterAll(void);
1590014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDDSSetUmin(Mat,PetscReal);
1591014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDWPSetComputeNormU(Mat,PetscBool );
1592e884886fSBarry Smith 
1593014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *);
1594014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *);
15957dbadf16SMatthew Knepley 
159697969023SHong Zhang /*
159797969023SHong Zhang    PETSc interface to MUMPS
159897969023SHong Zhang */
159997969023SHong Zhang #ifdef PETSC_HAVE_MUMPS
1600014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMumpsSetIcntl(Mat,PetscInt,PetscInt);
1601bc6112feSHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetIcntl(Mat,PetscInt,PetscInt*);
1602b9e9ea26SSatish Balay PETSC_EXTERN PetscErrorCode MatMumpsSetCntl(Mat,PetscInt,PetscReal);
1603bc6112feSHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetCntl(Mat,PetscInt,PetscReal*);
1604bc6112feSHong Zhang 
1605ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetInfo(Mat,PetscInt,PetscInt*);
1606ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetInfog(Mat,PetscInt,PetscInt*);
1607ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetRinfo(Mat,PetscInt,PetscReal*);
1608ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetRinfog(Mat,PetscInt,PetscReal*);
160997969023SHong Zhang #endif
161023a5497aSJed Brown 
1611d954db57SHong Zhang /*
1612*d615f992SSatish Balay    PETSc interface to Mkl_Pardiso
1613b500a6b7SJose David Bermeol */
1614b500a6b7SJose David Bermeol #ifdef PETSC_HAVE_MKL_PARDISO
1615*d615f992SSatish Balay PETSC_EXTERN PetscErrorCode MatMkl_PardisoSetCntl(Mat,PetscInt,PetscInt);
1616b500a6b7SJose David Bermeol #endif
1617b500a6b7SJose David Bermeol 
1618b500a6b7SJose David Bermeol /*
1619d954db57SHong Zhang    PETSc interface to SUPERLU
1620d954db57SHong Zhang */
1621d954db57SHong Zhang #ifdef PETSC_HAVE_SUPERLU
1622014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSuperluSetILUDropTol(Mat,PetscReal);
1623d954db57SHong Zhang #endif
1624d954db57SHong Zhang 
1625aa372e3fSPaul Mullowney #ifdef PETSC_HAVE_CUDA
16263f9c0db1SPaul Mullowney /*E
1627e057df02SPaul Mullowney     MatCUSPARSEStorageFormat - indicates the storage format for CUSPARSE (GPU)
16282692e278SPaul Mullowney     matrices.
1629e057df02SPaul Mullowney 
1630e057df02SPaul Mullowney     Not Collective
1631e057df02SPaul Mullowney 
16328468deeeSKarl Rupp +   MAT_CUSPARSE_CSR - Compressed Sparse Row
16332692e278SPaul Mullowney .   MAT_CUSPARSE_ELL - Ellpack (requires CUDA 4.2 or later).
16342692e278SPaul Mullowney -   MAT_CUSPARSE_HYB - Hybrid, a combination of Ellpack and Coordinate format (requires CUDA 4.2 or later).
1635e057df02SPaul Mullowney 
1636e057df02SPaul Mullowney     Level: intermediate
1637e057df02SPaul Mullowney 
1638e057df02SPaul Mullowney    Any additions/changes here MUST also be made in include/finclude/petscmat.h
1639e057df02SPaul Mullowney 
1640e057df02SPaul Mullowney .seealso: MatCUSPARSESetFormat(), MatCUSPARSEFormatOperation
1641e057df02SPaul Mullowney E*/
1642e057df02SPaul Mullowney 
16433f9c0db1SPaul Mullowney typedef enum {MAT_CUSPARSE_CSR, MAT_CUSPARSE_ELL, MAT_CUSPARSE_HYB} MatCUSPARSEStorageFormat;
1644e057df02SPaul Mullowney 
1645e057df02SPaul Mullowney /* these will be strings associated with enumerated type defined above */
1646e057df02SPaul Mullowney PETSC_EXTERN const char *const MatCUSPARSEStorageFormats[];
1647e057df02SPaul Mullowney 
16483f9c0db1SPaul Mullowney /*E
1649e057df02SPaul Mullowney     MatCUSPARSEFormatOperation - indicates the operation of CUSPARSE (GPU)
16502692e278SPaul Mullowney     matrices whose operation should use a particular storage format.
1651e057df02SPaul Mullowney 
1652e057df02SPaul Mullowney     Not Collective
1653e057df02SPaul Mullowney 
16548468deeeSKarl Rupp +   MAT_CUSPARSE_MULT_DIAG - sets the storage format for the diagonal matrix in the parallel MatMult
16558468deeeSKarl Rupp .   MAT_CUSPARSE_MULT_OFFDIAG - sets the storage format for the offdiagonal matrix in the parallel MatMult
16568468deeeSKarl Rupp .   MAT_CUSPARSE_MULT - sets the storage format for the entire matrix in the serial (single GPU) MatMult
16578468deeeSKarl Rupp -   MAT_CUSPARSE_ALL - sets the storage format for all CUSPARSE (GPU) matrices
1658e057df02SPaul Mullowney 
1659e057df02SPaul Mullowney     Level: intermediate
1660e057df02SPaul Mullowney 
1661e057df02SPaul Mullowney .seealso: MatCUSPARSESetFormat(), MatCUSPARSEStorageFormat
1662e057df02SPaul Mullowney E*/
166336d62e41SPaul Mullowney typedef enum {MAT_CUSPARSE_MULT_DIAG, MAT_CUSPARSE_MULT_OFFDIAG, MAT_CUSPARSE_MULT, MAT_CUSPARSE_ALL} MatCUSPARSEFormatOperation;
1664e057df02SPaul Mullowney 
166510e9db80SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
166610e9db80SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCreateAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
1667e057df02SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCUSPARSESetFormat(Mat,MatCUSPARSEFormatOperation,MatCUSPARSEStorageFormat);
16689ae82921SPaul Mullowney #endif
16699ae82921SPaul Mullowney 
167090c53902SBarry Smith #if defined(PETSC_HAVE_CUSP)
1671014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCUSP(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
1672014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateAIJCUSP(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
1673e057df02SPaul Mullowney 
16743f9c0db1SPaul Mullowney /*E
1675e057df02SPaul Mullowney     MatCUSPStorageFormat - indicates the storage format for CUSP (GPU)
167636d62e41SPaul Mullowney     matrices.
1677e057df02SPaul Mullowney 
1678e057df02SPaul Mullowney     Not Collective
1679e057df02SPaul Mullowney 
16808468deeeSKarl Rupp +   MAT_CUSP_CSR - Compressed Sparse Row
16818468deeeSKarl Rupp .   MAT_CUSP_DIA - Diagonal
16828468deeeSKarl Rupp -   MAT_CUSP_ELL - Ellpack
1683e057df02SPaul Mullowney 
1684e057df02SPaul Mullowney     Level: intermediate
1685e057df02SPaul Mullowney 
1686e057df02SPaul Mullowney    Any additions/changes here MUST also be made in include/finclude/petscmat.h
1687e057df02SPaul Mullowney 
1688e057df02SPaul Mullowney .seealso: MatCUSPSetFormat(), MatCUSPFormatOperation
1689e057df02SPaul Mullowney E*/
16903f9c0db1SPaul Mullowney typedef enum {MAT_CUSP_CSR, MAT_CUSP_DIA, MAT_CUSP_ELL} MatCUSPStorageFormat;
1691e057df02SPaul Mullowney 
1692e057df02SPaul Mullowney /* these will be strings associated with enumerated type defined above */
1693e057df02SPaul Mullowney PETSC_EXTERN const char *const MatCUSPStorageFormats[];
1694e057df02SPaul Mullowney 
16953f9c0db1SPaul Mullowney /*E
1696e057df02SPaul Mullowney     MatCUSPFormatOperation - indicates the operation of CUSP (GPU)
169736d62e41SPaul Mullowney     matrices whose operation should use a particular storage format.
1698e057df02SPaul Mullowney 
1699e057df02SPaul Mullowney     Not Collective
1700e057df02SPaul Mullowney 
17018468deeeSKarl Rupp +   MAT_CUSP_MULT_DIAG - sets the storage format for the diagonal matrix in the parallel MatMult
17028468deeeSKarl Rupp .   MAT_CUSP_MULT_OFFDIAG - sets the storage format for the offdiagonal matrix in the parallel MatMult
17038468deeeSKarl Rupp .   MAT_CUSP_MULT - sets the storage format for the entire matrix in the serial (single GPU) MatMult
17048468deeeSKarl Rupp -   MAT_CUSP_ALL - sets the storage format for all CUSP (GPU) matrices
1705e057df02SPaul Mullowney 
1706e057df02SPaul Mullowney     Level: intermediate
1707e057df02SPaul Mullowney 
1708e057df02SPaul Mullowney    Any additions/changes here MUST also be made in include/finclude/petscmat.h
1709e057df02SPaul Mullowney 
1710e057df02SPaul Mullowney .seealso: MatCUSPSetFormat(), MatCUSPStorageFormat
1711e057df02SPaul Mullowney E*/
171236d62e41SPaul Mullowney typedef enum {MAT_CUSP_MULT_DIAG, MAT_CUSP_MULT_OFFDIAG, MAT_CUSP_MULT, MAT_CUSP_ALL} MatCUSPFormatOperation;
1713e057df02SPaul Mullowney 
1714e057df02SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCUSPSetFormat(Mat,MatCUSPFormatOperation,MatCUSPStorageFormat);
1715e057df02SPaul Mullowney #endif
171690c53902SBarry Smith 
1717d67ff14aSKarl Rupp #if defined(PETSC_HAVE_VIENNACL)
1718d67ff14aSKarl Rupp PETSC_EXTERN PetscErrorCode MatCreateSeqAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
1719d67ff14aSKarl Rupp PETSC_EXTERN PetscErrorCode MatCreateAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
1720d67ff14aSKarl Rupp #endif
1721d67ff14aSKarl Rupp 
172254efbe56SHong Zhang /*
172354efbe56SHong Zhang    PETSc interface to FFTW
172454efbe56SHong Zhang */
172554efbe56SHong Zhang #if defined(PETSC_HAVE_FFTW)
1726014dd563SJed Brown PETSC_EXTERN PetscErrorCode VecScatterPetscToFFTW(Mat,Vec,Vec);
1727014dd563SJed Brown PETSC_EXTERN PetscErrorCode VecScatterFFTWToPetsc(Mat,Vec,Vec);
1728014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetVecsFFTW(Mat,Vec*,Vec*,Vec*);
172954efbe56SHong Zhang #endif
173054efbe56SHong Zhang 
1731db31f6deSJed Brown #if defined(PETSC_HAVE_ELEMENTAL)
1732607a6623SBarry Smith PETSC_EXTERN PetscErrorCode PetscElementalInitializePackage(void);
1733db31f6deSJed Brown PETSC_EXTERN PetscErrorCode PetscElementalFinalizePackage(void);
1734db31f6deSJed Brown #endif
1735db31f6deSJed Brown 
1736014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateNest(MPI_Comm,PetscInt,const IS[],PetscInt,const IS[],const Mat[],Mat*);
1737014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSize(Mat,PetscInt*,PetscInt*);
1738014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetISs(Mat,IS[],IS[]);
1739014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetLocalISs(Mat,IS[],IS[]);
1740014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSubMats(Mat,PetscInt*,PetscInt*,Mat***);
1741014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSubMat(Mat,PetscInt,PetscInt,Mat*);
174219fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatNestSetVecType(Mat,VecType);
1743014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestSetSubMats(Mat,PetscInt,const IS[],PetscInt,const IS[],const Mat[]);
1744014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestSetSubMat(Mat,PetscInt,PetscInt,Mat);
1745d8588912SDave May 
17464325cce7SMatthew G Knepley PETSC_EXTERN PetscErrorCode MatChop(Mat,PetscReal);
1747e0a58c10SMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatComputeBandwidth(Mat,PetscReal,PetscInt*);
17484325cce7SMatthew G Knepley 
174923a5497aSJed Brown #endif
1750