xref: /petsc/include/petscmat.h (revision 0d04baf89eca684847eb03c0d4ee6b8dfe3c620a)
12eac72dbSBarry Smith /*
22eac72dbSBarry Smith      Include file for the matrix component of PETSc
32eac72dbSBarry Smith */
40a835dfdSSatish Balay #ifndef __PETSCMAT_H
50a835dfdSSatish Balay #define __PETSCMAT_H
60a835dfdSSatish Balay #include "petscvec.h"
7e9fa29b7SSatish Balay PETSC_EXTERN_CXX_BEGIN
82eac72dbSBarry Smith 
9d9274352SBarry Smith /*S
10d9274352SBarry Smith      Mat - Abstract PETSc matrix object
112eac72dbSBarry Smith 
12d91e6319SBarry Smith    Level: beginner
13d91e6319SBarry Smith 
14d9274352SBarry Smith   Concepts: matrix; linear operator
15d9274352SBarry Smith 
16d9274352SBarry Smith .seealso:  MatCreate(), MatType, MatSetType()
17d9274352SBarry Smith S*/
18d9274352SBarry Smith typedef struct _p_Mat*           Mat;
19d9274352SBarry Smith 
2076bdecfbSBarry Smith /*J
21d9274352SBarry Smith     MatType - String with the name of a PETSc matrix or the creation function
22d9274352SBarry Smith        with an optional dynamic library name, for example
23d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:mymatcreate()
24d9274352SBarry Smith 
25d9274352SBarry Smith    Level: beginner
26d9274352SBarry Smith 
27c7393fdbSBarry Smith .seealso: MatSetType(), Mat, MatSolverPackage
2876bdecfbSBarry Smith J*/
29a313700dSBarry Smith #define MatType char*
30273d9f13SBarry Smith #define MATSAME            "same"
315a11e1b2SBarry Smith #define MATMAIJ            "maij"
32273d9f13SBarry Smith #define MATSEQMAIJ         "seqmaij"
33273d9f13SBarry Smith #define MATMPIMAIJ         "mpimaij"
34273d9f13SBarry Smith #define MATIS              "is"
355a11e1b2SBarry Smith #define MATAIJ             "aij"
36273d9f13SBarry Smith #define MATSEQAIJ          "seqaij"
377d6a0e61SBarry Smith #define MATSEQAIJPTHREAD   "seqaijpthread"
38bf2c1783SBarry Smith #define MATAIJPTHREAD      "aijpthread"
39273d9f13SBarry Smith #define MATMPIAIJ          "mpiaij"
405a11e1b2SBarry Smith #define MATAIJCRL          "aijcrl"
415a11e1b2SBarry Smith #define MATSEQAIJCRL       "seqaijcrl"
425a11e1b2SBarry Smith #define MATMPIAIJCRL       "mpiaijcrl"
438154be41SBarry Smith #define MATAIJCUSP         "aijcusp"
448154be41SBarry Smith #define MATSEQAIJCUSP      "seqaijcusp"
458154be41SBarry Smith #define MATMPIAIJCUSP      "mpiaijcusp"
465a11e1b2SBarry Smith #define MATAIJPERM         "aijperm"
475a11e1b2SBarry Smith #define MATSEQAIJPERM      "seqaijperm"
485a11e1b2SBarry Smith #define MATMPIAIJPERM      "mpiaijperm"
49273d9f13SBarry Smith #define MATSHELL           "shell"
505a11e1b2SBarry Smith #define MATDENSE           "dense"
51209238afSKris Buschelman #define MATSEQDENSE        "seqdense"
52273d9f13SBarry Smith #define MATMPIDENSE        "mpidense"
535a11e1b2SBarry Smith #define MATBAIJ            "baij"
54273d9f13SBarry Smith #define MATSEQBAIJ         "seqbaij"
55273d9f13SBarry Smith #define MATMPIBAIJ         "mpibaij"
56273d9f13SBarry Smith #define MATMPIADJ          "mpiadj"
575a11e1b2SBarry Smith #define MATSBAIJ           "sbaij"
58273d9f13SBarry Smith #define MATSEQSBAIJ        "seqsbaij"
59273d9f13SBarry Smith #define MATMPISBAIJ        "mpisbaij"
60c0cdd4a1SDahai Guo #define MATSEQBSTRM        "seqbstrm"
61c0cdd4a1SDahai Guo #define MATMPIBSTRM        "mpibstrm"
62c0cdd4a1SDahai Guo #define MATBSTRM           "bstrm"
63aa5a9175SDahai Guo #define MATSEQSBSTRM       "seqsbstrm"
64aa5a9175SDahai Guo #define MATMPISBSTRM       "mpisbstrm"
65aa5a9175SDahai Guo #define MATSBSTRM          "sbstrm"
66cebc7f6cSBarry Smith #define MATDAAD            "daad"
67cebc7f6cSBarry Smith #define MATMFFD            "mffd"
68c8a8475eSBarry Smith #define MATNORMAL          "normal"
69ab92ecdeSBarry Smith #define MATLRC             "lrc"
702a6744ebSBarry Smith #define MATSCATTER         "scatter"
71421e10b8SBarry Smith #define MATBLOCKMAT        "blockmat"
72793850ffSBarry Smith #define MATCOMPOSITE       "composite"
731f8c7532SHong Zhang #define MATFFT             "fft"
741f8c7532SHong Zhang #define MATFFTW            "fftw"
75e133240eSMatthew G Knepley #define MATSEQCUFFT        "seqcufft"
76557cca28SSatish Balay #define MATTRANSPOSEMAT    "transpose"
7772ca8751SBarry Smith #define MATSCHURCOMPLEMENT "schurcomplement"
781d6018f0SLisandro Dalcin #define MATPYTHON          "python"
79f91d8e95SBarry Smith #define MATHYPRESTRUCT     "hyprestruct"
80a9e6138eSGlenn Hammond #define MATHYPRESSTRUCT    "hypresstruct"
81ee1cef2cSJed Brown #define MATSUBMATRIX       "submatrix"
822c0dbf93SJed Brown #define MATLOCALREF        "localref"
83d8588912SDave May #define MATNEST            "nest"
8415c1f1baSDmitry Karpeev #define MATIJ              "ij"
85773941d6SBarry Smith 
8676bdecfbSBarry Smith /*J
87c7393fdbSBarry Smith     MatSolverPackage - String with the name of a PETSc matrix solver type.
889989ab13SBarry Smith 
899989ab13SBarry Smith     For example: "petsc" indicates what PETSc provides, "superlu" indicates either
909989ab13SBarry Smith        SuperLU or SuperLU_Dist etc.
919989ab13SBarry Smith 
929989ab13SBarry Smith 
939989ab13SBarry Smith    Level: beginner
949989ab13SBarry Smith 
955c9eb25fSBarry Smith .seealso: MatGetFactor(), Mat, MatSetType(), MatType
9676bdecfbSBarry Smith J*/
97c7393fdbSBarry Smith #define MatSolverPackage char*
982692d6eeSBarry Smith #define MATSOLVERSPOOLES      "spooles"
992692d6eeSBarry Smith #define MATSOLVERSUPERLU      "superlu"
1002692d6eeSBarry Smith #define MATSOLVERSUPERLU_DIST "superlu_dist"
1012692d6eeSBarry Smith #define MATSOLVERUMFPACK      "umfpack"
10220db9a53SJed Brown #define MATSOLVERCHOLMOD      "cholmod"
1032692d6eeSBarry Smith #define MATSOLVERESSL         "essl"
1042692d6eeSBarry Smith #define MATSOLVERLUSOL        "lusol"
1052692d6eeSBarry Smith #define MATSOLVERMUMPS        "mumps"
1062692d6eeSBarry Smith #define MATSOLVERPASTIX       "pastix"
1072692d6eeSBarry Smith #define MATSOLVERMATLAB       "matlab"
1082692d6eeSBarry Smith #define MATSOLVERPETSC        "petsc"
1092692d6eeSBarry Smith #define MATSOLVERPLAPACK      "plapack"
1102692d6eeSBarry Smith #define MATSOLVERBAS          "bas"
111773941d6SBarry Smith 
112aa5a9175SDahai Guo #define MATSOLVERBSTRM        "bstrm"
113aa5a9175SDahai Guo #define MATSOLVERSBSTRM       "sbstrm"
114c0cdd4a1SDahai Guo 
115b24902e0SBarry Smith /*E
116b24902e0SBarry Smith     MatFactorType - indicates what type of factorization is requested
117b24902e0SBarry Smith 
118b24902e0SBarry Smith     Level: beginner
119b24902e0SBarry Smith 
120b24902e0SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
121b24902e0SBarry Smith 
122c7393fdbSBarry Smith .seealso: MatSolverPackage, MatGetFactor()
123b24902e0SBarry Smith E*/
124599ef60dSHong Zhang typedef enum {MAT_FACTOR_NONE, MAT_FACTOR_LU, MAT_FACTOR_CHOLESKY, MAT_FACTOR_ILU, MAT_FACTOR_ICC,MAT_FACTOR_ILUDT} MatFactorType;
12534918c53SJed Brown extern const char *const MatFactorTypes[];
126e92e720dSBarry Smith 
1277087cfbeSBarry Smith extern PetscErrorCode  MatGetFactor(Mat,const MatSolverPackage,MatFactorType,Mat*);
1287087cfbeSBarry Smith extern PetscErrorCode  MatGetFactorAvailable(Mat,const MatSolverPackage,MatFactorType,PetscBool *);
1297087cfbeSBarry Smith extern PetscErrorCode  MatFactorGetSolverPackage(Mat,const MatSolverPackage*);
1307087cfbeSBarry Smith extern PetscErrorCode  MatGetFactorType(Mat,MatFactorType*);
1319989ab13SBarry Smith 
132c06d978dSMatthew Knepley /* Logging support */
1330700a824SBarry Smith #define    MAT_FILE_CLASSID 1211216    /* used to indicate matrices in binary files */
1347087cfbeSBarry Smith extern PetscClassId  MAT_CLASSID;
1357087cfbeSBarry Smith extern PetscClassId  MAT_FDCOLORING_CLASSID;
136b9af6bddSHong Zhang extern PetscClassId  MAT_TRANSPOSECOLORING_CLASSID;
1377087cfbeSBarry Smith extern PetscClassId  MAT_PARTITIONING_CLASSID;
1387087cfbeSBarry Smith extern PetscClassId  MAT_NULLSPACE_CLASSID;
1397087cfbeSBarry Smith extern PetscClassId  MATMFFD_CLASSID;
140c06d978dSMatthew Knepley 
141ceb03754SKris Buschelman /*E
142ceb03754SKris Buschelman     MatReuse - Indicates if matrices obtained from a previous call to MatGetSubMatrices()
143d6eff37eSBarry Smith      or MatGetSubMatrix() are to be reused to store the new matrix values. For MatConvert() is used to indicate
144d6eff37eSBarry Smith      that the input matrix is to be replaced with the converted matrix.
145ceb03754SKris Buschelman 
146ceb03754SKris Buschelman     Level: beginner
147ceb03754SKris Buschelman 
148ceb03754SKris Buschelman    Any additions/changes here MUST also be made in include/finclude/petscmat.h
149ceb03754SKris Buschelman 
1500c451bc4SBarry Smith .seealso: MatGetSubMatrices(), MatGetSubMatrix(), MatDestroyMatrices(), MatConvert()
151ceb03754SKris Buschelman E*/
152dfe085dbSJed Brown typedef enum {MAT_INITIAL_MATRIX,MAT_REUSE_MATRIX,MAT_IGNORE_MATRIX} MatReuse;
153ceb03754SKris Buschelman 
1545494a064SHong Zhang /*E
1555494a064SHong Zhang     MatGetSubMatrixOption - Indicates if matrices obtained from a call to MatGetSubMatrices()
156829201f2SHong Zhang      include the matrix values. Currently it is only used by MatGetSeqNonzerostructure().
1575494a064SHong Zhang 
1585494a064SHong Zhang     Level: beginner
1595494a064SHong Zhang 
160829201f2SHong Zhang .seealso: MatGetSeqNonzerostructure()
1615494a064SHong Zhang E*/
1625494a064SHong Zhang typedef enum {MAT_DO_NOT_GET_VALUES,MAT_GET_VALUES} MatGetSubMatrixOption;
1635494a064SHong Zhang 
1647087cfbeSBarry Smith extern PetscErrorCode  MatInitializePackage(const char[]);
165c06d978dSMatthew Knepley 
1667087cfbeSBarry Smith extern PetscErrorCode  MatCreate(MPI_Comm,Mat*);
167a69a960dSMatthew Knepley PetscPolymorphicFunction(MatCreate,(MPI_Comm comm),(comm,&A),Mat,A)
168a69a960dSMatthew Knepley PetscPolymorphicFunction(MatCreate,(),(PETSC_COMM_WORLD,&A),Mat,A)
1697087cfbeSBarry Smith extern PetscErrorCode  MatSetSizes(Mat,PetscInt,PetscInt,PetscInt,PetscInt);
1707087cfbeSBarry Smith extern PetscErrorCode  MatSetType(Mat,const MatType);
1717087cfbeSBarry Smith extern PetscErrorCode  MatSetFromOptions(Mat);
1727087cfbeSBarry Smith extern PetscErrorCode  MatSetUpPreallocation(Mat);
1737087cfbeSBarry Smith extern PetscErrorCode  MatRegisterAll(const char[]);
1747087cfbeSBarry Smith extern PetscErrorCode  MatRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat));
17501bebe75SBarry Smith extern PetscErrorCode  MatRegisterBaseName(const char[],const char[],const char[]);
1767087cfbeSBarry Smith extern PetscErrorCode  MatSetOptionsPrefix(Mat,const char[]);
1777087cfbeSBarry Smith extern PetscErrorCode  MatAppendOptionsPrefix(Mat,const char[]);
1787087cfbeSBarry Smith extern PetscErrorCode  MatGetOptionsPrefix(Mat,const char*[]);
179f69a0ea3SMatthew Knepley 
18030de9b25SBarry Smith /*MC
18130de9b25SBarry Smith    MatRegisterDynamic - Adds a new matrix type
18230de9b25SBarry Smith 
18330de9b25SBarry Smith    Synopsis:
1841890ba74SBarry Smith    PetscErrorCode MatRegisterDynamic(const char *name,const char *path,const char *name_create,PetscErrorCode (*routine_create)(Mat))
18530de9b25SBarry Smith 
18630de9b25SBarry Smith    Not Collective
18730de9b25SBarry Smith 
18830de9b25SBarry Smith    Input Parameters:
18930de9b25SBarry Smith +  name - name of a new user-defined matrix type
19030de9b25SBarry Smith .  path - path (either absolute or relative) the library containing this solver
19130de9b25SBarry Smith .  name_create - name of routine to create method context
19230de9b25SBarry Smith -  routine_create - routine to create method context
19330de9b25SBarry Smith 
19430de9b25SBarry Smith    Notes:
19530de9b25SBarry Smith    MatRegisterDynamic() may be called multiple times to add several user-defined solvers.
19630de9b25SBarry Smith 
19730de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
19830de9b25SBarry Smith    is ignored.
19930de9b25SBarry Smith 
20030de9b25SBarry Smith    Sample usage:
20130de9b25SBarry Smith .vb
20230de9b25SBarry Smith    MatRegisterDynamic("my_mat",/home/username/my_lib/lib/libO/solaris/mylib.a,
20330de9b25SBarry Smith                "MyMatCreate",MyMatCreate);
20430de9b25SBarry Smith .ve
20530de9b25SBarry Smith 
20630de9b25SBarry Smith    Then, your solver can be chosen with the procedural interface via
20730de9b25SBarry Smith $     MatSetType(Mat,"my_mat")
20830de9b25SBarry Smith    or at runtime via the option
20930de9b25SBarry Smith $     -mat_type my_mat
21030de9b25SBarry Smith 
21130de9b25SBarry Smith    Level: advanced
21230de9b25SBarry Smith 
213ab901514SBarry Smith    Notes: ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values.
21430de9b25SBarry Smith          If your function is not being put into a shared library then use VecRegister() instead
21530de9b25SBarry Smith 
21630de9b25SBarry Smith .keywords: Mat, register
21730de9b25SBarry Smith 
21830de9b25SBarry Smith .seealso: MatRegisterAll(), MatRegisterDestroy()
21930de9b25SBarry Smith 
22030de9b25SBarry Smith M*/
221273d9f13SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
222273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,0)
223273d9f13SBarry Smith #else
224273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,d)
22530de9b25SBarry Smith #endif
22630de9b25SBarry Smith 
227ace3abfcSBarry Smith extern PetscBool  MatRegisterAllCalled;
228b0a32e0cSBarry Smith extern PetscFList MatList;
229b022a5c1SBarry Smith extern PetscFList MatColoringList;
230b022a5c1SBarry Smith extern PetscFList MatPartitioningList;
23128988994SBarry Smith 
2323b224e63SBarry Smith /*E
2333b224e63SBarry Smith     MatStructure - Indicates if the matrix has the same nonzero structure
2343b224e63SBarry Smith 
2353b224e63SBarry Smith     Level: beginner
2363b224e63SBarry Smith 
2373b224e63SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
2383b224e63SBarry Smith 
2393b224e63SBarry Smith .seealso: MatCopy(), KSPSetOperators(), PCSetOperators()
2403b224e63SBarry Smith E*/
241214bc6a2SJed Brown typedef enum {DIFFERENT_NONZERO_PATTERN,SUBSET_NONZERO_PATTERN,SAME_NONZERO_PATTERN,SAME_PRECONDITIONER} MatStructure;
2423b224e63SBarry Smith 
2437087cfbeSBarry Smith extern PetscErrorCode  MatCreateSeqDense(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*);
2447087cfbeSBarry Smith extern PetscErrorCode  MatCreateMPIDense(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*);
2457087cfbeSBarry Smith extern PetscErrorCode  MatCreateSeqAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
246ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,nz,nnz,&A),Mat,A)
247ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,&A),Mat,A)
248ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,0,nnz,&A),Mat,A)
249ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,&A),Mat,A)
250ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,A))
251ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,m,n,0,nnz,A))
252ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,A))
2537087cfbeSBarry Smith extern PetscErrorCode  MatCreateMPIAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
254ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPIAIJ,(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,const PetscInt nnz[],PetscInt onz,const PetscInt onnz[]),(comm,m,n,M,N,nz,nnz,onz,onnz,&A),Mat,A)
255ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPIAIJ,(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,PetscInt nnz),(comm,m,n,M,N,nz,PETSC_NULL,nnz,PETSC_NULL,&A),Mat,A)
256ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPIAIJ,(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt nnz[],const PetscInt onz[]),(comm,m,n,M,N,0,nnz,0,onz,&A),Mat,A)
257ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPIAIJ,(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N),(comm,m,n,M,N,0,PETSC_NULL,0,PETSC_NULL,&A),Mat,A)
258ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPIAIJ,(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,PetscInt nnz,Mat *A),(comm,m,n,M,N,nz,PETSC_NULL,nnz,PETSC_NULL,A))
259ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPIAIJ,(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt nnz[],const PetscInt onz[],Mat *A),(comm,m,n,M,N,0,nnz,0,onz,A))
260ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPIAIJ,(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,Mat *A),(comm,m,n,M,N,0,PETSC_NULL,0,PETSC_NULL,A))
261ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPIAIJ,(PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,const PetscInt nnz[],PetscInt onz,const PetscInt onnz[]),(PETSC_COMM_WORLD,m,n,M,N,nz,nnz,onz,onnz,&A),Mat,A)
262ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPIAIJ,(PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,PetscInt nnz),(PETSC_COMM_WORLD,m,n,M,N,nz,PETSC_NULL,nnz,PETSC_NULL,&A),Mat,A)
263ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPIAIJ,(PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt nnz[],const PetscInt onz[]),(PETSC_COMM_WORLD,m,n,M,N,0,nnz,0,onz,&A),Mat,A)
264ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPIAIJ,(PetscInt m,PetscInt n,PetscInt M,PetscInt N),(PETSC_COMM_WORLD,m,n,M,N,0,PETSC_NULL,0,PETSC_NULL,&A),Mat,A)
265ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPIAIJ,(PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,PetscInt nnz,Mat *A),(PETSC_COMM_WORLD,m,n,M,N,nz,PETSC_NULL,nnz,PETSC_NULL,A))
266ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPIAIJ,(PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt nnz[],const PetscInt onz[],Mat *A),(PETSC_COMM_WORLD,m,n,M,N,0,nnz,0,onz,A))
267ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPIAIJ,(PetscInt m,PetscInt n,PetscInt M,PetscInt N,Mat *A),(PETSC_COMM_WORLD,m,n,M,N,0,PETSC_NULL,0,PETSC_NULL,A))
2687087cfbeSBarry Smith extern PetscErrorCode  MatCreateMPIAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *);
2697087cfbeSBarry Smith extern PetscErrorCode  MatCreateMPIAIJWithSplitArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],PetscInt[],PetscInt[],PetscScalar[],Mat*);
2708d7a6e47SBarry Smith 
2717087cfbeSBarry Smith extern PetscErrorCode  MatCreateSeqBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
272ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,nz,nnz,&A),Mat,A)
273ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A)
274ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A)
275ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A)
276ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A))
277ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A))
278ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A))
2797087cfbeSBarry Smith extern PetscErrorCode  MatCreateMPIBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
280ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPIBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,const PetscInt nnz[],PetscInt onz,const PetscInt onnz[]),(comm,bs,m,n,M,N,nz,nnz,onz,onnz,&A),Mat,A)
281ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPIBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,PetscInt nnz),(comm,bs,m,n,M,N,nz,PETSC_NULL,nnz,PETSC_NULL,&A),Mat,A)
282ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPIBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt nnz[],const PetscInt onz[]),(comm,bs,m,n,M,N,0,nnz,0,onz,&A),Mat,A)
283ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPIBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N),(comm,bs,m,n,M,N,0,PETSC_NULL,0,PETSC_NULL,&A),Mat,A)
284ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPIBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,PetscInt nnz,Mat *A),(comm,bs,m,n,M,N,nz,PETSC_NULL,nnz,PETSC_NULL,A))
285ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPIBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt nnz[],const PetscInt onz[],Mat *A),(comm,bs,m,n,M,N,0,nnz,0,onz,A))
286ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPIBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,Mat *A),(comm,bs,m,n,M,N,0,PETSC_NULL,0,PETSC_NULL,A))
287ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPIBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,const PetscInt nnz[],PetscInt onz,const PetscInt onnz[]),(PETSC_COMM_WORLD,bs,m,n,M,N,nz,nnz,onz,onnz,&A),Mat,A)
288ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPIBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,PetscInt nnz),(PETSC_COMM_WORLD,bs,m,n,M,N,nz,PETSC_NULL,nnz,PETSC_NULL,&A),Mat,A)
289ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPIBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt nnz[],const PetscInt onz[]),(PETSC_COMM_WORLD,bs,m,n,M,N,0,nnz,0,onz,&A),Mat,A)
290ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPIBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N),(PETSC_COMM_WORLD,bs,m,n,M,N,0,PETSC_NULL,0,PETSC_NULL,&A),Mat,A)
291ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPIBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,PetscInt nnz,Mat *A),(PETSC_COMM_WORLD,bs,m,n,M,N,nz,PETSC_NULL,nnz,PETSC_NULL,A))
292ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPIBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt nnz[],const PetscInt onz[],Mat *A),(PETSC_COMM_WORLD,bs,m,n,M,N,0,nnz,0,onz,A))
293ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPIBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,Mat *A),(PETSC_COMM_WORLD,bs,m,n,M,N,0,PETSC_NULL,0,PETSC_NULL,A))
2947087cfbeSBarry Smith extern PetscErrorCode  MatCreateMPIBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat*);
295d21a29f3SJed Brown 
2967087cfbeSBarry Smith extern PetscErrorCode  MatCreateMPIAdj(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscInt[],Mat*);
2977087cfbeSBarry Smith extern PetscErrorCode  MatCreateSeqSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
298ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,nz,nnz,&A),Mat,A)
299ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A)
300ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A)
301ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A)
302ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A))
303ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A))
304ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A))
305d21a29f3SJed Brown 
3067087cfbeSBarry Smith extern PetscErrorCode  MatCreateMPISBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
307ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPISBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,const PetscInt nnz[],PetscInt onz,const PetscInt onnz[]),(comm,bs,m,n,M,N,nz,nnz,onz,onnz,&A),Mat,A)
308ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPISBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,PetscInt nnz),(comm,bs,m,n,M,N,nz,PETSC_NULL,nnz,PETSC_NULL,&A),Mat,A)
309ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPISBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt nnz[],const PetscInt onz[]),(comm,bs,m,n,M,N,0,nnz,0,onz,&A),Mat,A)
310ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPISBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N),(comm,bs,m,n,M,N,0,PETSC_NULL,0,PETSC_NULL,&A),Mat,A)
311ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPISBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,PetscInt nnz,Mat *A),(comm,bs,m,n,M,N,nz,PETSC_NULL,nnz,PETSC_NULL,A))
312ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPISBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt nnz[],const PetscInt onz[],Mat *A),(comm,bs,m,n,M,N,0,nnz,0,onz,A))
313ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPISBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,Mat *A),(comm,bs,m,n,M,N,0,PETSC_NULL,0,PETSC_NULL,A))
314ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPISBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,const PetscInt nnz[],PetscInt onz,const PetscInt onnz[]),(PETSC_COMM_WORLD,bs,m,n,M,N,nz,nnz,onz,onnz,&A),Mat,A)
315ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPISBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,PetscInt nnz),(PETSC_COMM_WORLD,bs,m,n,M,N,nz,PETSC_NULL,nnz,PETSC_NULL,&A),Mat,A)
316ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPISBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt nnz[],const PetscInt onz[]),(PETSC_COMM_WORLD,bs,m,n,M,N,0,nnz,0,onz,&A),Mat,A)
317ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPISBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N),(PETSC_COMM_WORLD,bs,m,n,M,N,0,PETSC_NULL,0,PETSC_NULL,&A),Mat,A)
318ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPISBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,PetscInt nnz,Mat *A),(PETSC_COMM_WORLD,bs,m,n,M,N,nz,PETSC_NULL,nnz,PETSC_NULL,A))
319ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPISBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt nnz[],const PetscInt onz[],Mat *A),(PETSC_COMM_WORLD,bs,m,n,M,N,0,nnz,0,onz,A))
320ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPISBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,Mat *A),(PETSC_COMM_WORLD,bs,m,n,M,N,0,PETSC_NULL,0,PETSC_NULL,A))
3217087cfbeSBarry Smith extern PetscErrorCode  MatCreateMPISBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *);
3227087cfbeSBarry Smith extern PetscErrorCode  MatMPISBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
323dfb205c3SBarry Smith 
3247087cfbeSBarry Smith extern PetscErrorCode  MatCreateShell(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,void *,Mat*);
325ba966639SSatish Balay PetscPolymorphicFunction(MatCreateShell,(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx),(comm,m,n,M,N,ctx,&A),Mat,A)
326ba966639SSatish Balay PetscPolymorphicFunction(MatCreateShell,(PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx),(PETSC_COMM_WORLD,m,n,M,N,ctx,&A),Mat,A)
3277087cfbeSBarry Smith extern PetscErrorCode  MatCreateAdic(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,void (*)(void),Mat*);
3287087cfbeSBarry Smith extern PetscErrorCode  MatCreateNormal(Mat,Mat*);
329ba966639SSatish Balay PetscPolymorphicFunction(MatCreateNormal,(Mat mat),(mat,&A),Mat,A)
3307087cfbeSBarry Smith extern PetscErrorCode  MatCreateLRC(Mat,Mat,Mat,Mat*);
3317087cfbeSBarry Smith extern PetscErrorCode  MatCreateIS(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,ISLocalToGlobalMapping,Mat*);
3327087cfbeSBarry Smith extern PetscErrorCode  MatCreateSeqAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
3337087cfbeSBarry Smith extern PetscErrorCode  MatCreateMPIAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
334c0cdd4a1SDahai Guo 
335c0cdd4a1SDahai Guo extern PetscErrorCode  MatCreateSeqBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
336c0cdd4a1SDahai Guo extern PetscErrorCode  MatCreateMPIBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
337aa5a9175SDahai Guo extern PetscErrorCode  MatCreateSeqSBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
338aa5a9175SDahai Guo extern PetscErrorCode  MatCreateMPISBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
339c0cdd4a1SDahai Guo 
3407087cfbeSBarry Smith extern PetscErrorCode  MatCreateScatter(MPI_Comm,VecScatter,Mat*);
3417087cfbeSBarry Smith extern PetscErrorCode  MatScatterSetVecScatter(Mat,VecScatter);
3427087cfbeSBarry Smith extern PetscErrorCode  MatScatterGetVecScatter(Mat,VecScatter*);
3437087cfbeSBarry Smith extern PetscErrorCode  MatCreateBlockMat(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,Mat*);
3447087cfbeSBarry Smith extern PetscErrorCode  MatCompositeAddMat(Mat,Mat);
3457087cfbeSBarry Smith extern PetscErrorCode  MatCompositeMerge(Mat);
3467087cfbeSBarry Smith extern PetscErrorCode  MatCreateComposite(MPI_Comm,PetscInt,const Mat*,Mat*);
3476d7c1e57SBarry Smith typedef enum {MAT_COMPOSITE_ADDITIVE,MAT_COMPOSITE_MULTIPLICATIVE} MatCompositeType;
3487087cfbeSBarry Smith extern PetscErrorCode  MatCompositeSetType(Mat,MatCompositeType);
3496d7c1e57SBarry Smith 
350dedccee8SHong Zhang extern PetscErrorCode  MatCreateFFT(MPI_Comm,PetscInt,const PetscInt[],const MatType,Mat*);
3517087cfbeSBarry Smith extern PetscErrorCode  MatCreateSeqCUFFT(MPI_Comm,PetscInt,const PetscInt[],Mat*);
352dedccee8SHong Zhang 
3537087cfbeSBarry Smith extern PetscErrorCode  MatCreateTranspose(Mat,Mat*);
3547087cfbeSBarry Smith extern PetscErrorCode  MatCreateSubMatrix(Mat,IS,IS,Mat*);
3557087cfbeSBarry Smith extern PetscErrorCode  MatSubMatrixUpdate(Mat,Mat,IS,IS);
3567087cfbeSBarry Smith extern PetscErrorCode  MatCreateLocalRef(Mat,IS,IS,Mat*);
3571472f72bSBarry Smith 
3587087cfbeSBarry Smith extern PetscErrorCode  MatPythonSetType(Mat,const char[]);
3591d6018f0SLisandro Dalcin 
3607087cfbeSBarry Smith extern PetscErrorCode  MatSetUp(Mat);
361fcfd50ebSBarry Smith extern PetscErrorCode  MatDestroy(Mat*);
36221c89e3eSBarry Smith 
3637087cfbeSBarry Smith extern PetscErrorCode  MatConjugate(Mat);
3647087cfbeSBarry Smith extern PetscErrorCode  MatRealPart(Mat);
3657087cfbeSBarry Smith extern PetscErrorCode  MatImaginaryPart(Mat);
36611bd1e4dSLisandro Dalcin extern PetscErrorCode  MatGetDiagonalBlock(Mat,Mat*);
3677087cfbeSBarry Smith extern PetscErrorCode  MatGetTrace(Mat,PetscScalar*);
368bbead8a2SBarry Smith extern PetscErrorCode  MatInvertBlockDiagonal(Mat,PetscScalar **);
36999cafbc1SBarry Smith 
3708ed539a5SBarry Smith /* ------------------------------------------------------------*/
3717087cfbeSBarry Smith extern PetscErrorCode  MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
3727087cfbeSBarry Smith extern PetscErrorCode  MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
3737087cfbeSBarry Smith extern PetscErrorCode  MatSetValuesRow(Mat,PetscInt,const PetscScalar[]);
3747087cfbeSBarry Smith extern PetscErrorCode  MatSetValuesRowLocal(Mat,PetscInt,const PetscScalar[]);
3752f35a324SMatthew G Knepley extern PetscErrorCode  MatSetValuesBatch(Mat,PetscInt,PetscInt,PetscInt[],const PetscScalar[]);
37684cb2905SBarry Smith 
3772ef4de8bSBarry Smith /*S
3782ef4de8bSBarry Smith      MatStencil - Data structure (C struct) for storing information about a single row or
379be479b30SJed Brown         column of a matrix as indexed on an associated grid.
380be479b30SJed Brown 
381be479b30SJed Brown    Fortran usage is different, see MatSetValuesStencil() for details.
3822ef4de8bSBarry Smith 
3832ef4de8bSBarry Smith    Level: beginner
3842ef4de8bSBarry Smith 
3852ef4de8bSBarry Smith   Concepts: matrix; linear operator
3862ef4de8bSBarry Smith 
387be479b30SJed Brown .seealso:  MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockedStencil()
3882ef4de8bSBarry Smith S*/
389435da068SBarry Smith typedef struct {
390c1ac3661SBarry Smith   PetscInt k,j,i,c;
391435da068SBarry Smith } MatStencil;
3922ef4de8bSBarry Smith 
3937087cfbeSBarry Smith extern PetscErrorCode  MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
3947087cfbeSBarry Smith extern PetscErrorCode  MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
3957087cfbeSBarry Smith extern PetscErrorCode  MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt);
396435da068SBarry Smith 
3977087cfbeSBarry Smith extern PetscErrorCode  MatSetColoring(Mat,ISColoring);
3987087cfbeSBarry Smith extern PetscErrorCode  MatSetValuesAdic(Mat,void*);
3997087cfbeSBarry Smith extern PetscErrorCode  MatSetValuesAdifor(Mat,PetscInt,void*);
4003a7fca6bSBarry Smith 
401d91e6319SBarry Smith /*E
402d91e6319SBarry Smith     MatAssemblyType - Indicates if the matrix is now to be used, or if you plan
403d91e6319SBarry Smith      to continue to add values to it
404d91e6319SBarry Smith 
405d91e6319SBarry Smith     Level: beginner
406d91e6319SBarry Smith 
407d91e6319SBarry Smith .seealso: MatAssemblyBegin(), MatAssemblyEnd()
408d91e6319SBarry Smith E*/
4096d4a8577SBarry Smith typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType;
4107087cfbeSBarry Smith extern PetscErrorCode  MatAssemblyBegin(Mat,MatAssemblyType);
4117087cfbeSBarry Smith extern PetscErrorCode  MatAssemblyEnd(Mat,MatAssemblyType);
4127087cfbeSBarry Smith extern PetscErrorCode  MatAssembled(Mat,PetscBool *);
4134f9c727eSBarry Smith 
4141d73ed98SBarry Smith 
41530de9b25SBarry Smith 
416d91e6319SBarry Smith /*E
417d91e6319SBarry Smith     MatOption - Options that may be set for a matrix and its behavior or storage
418d91e6319SBarry Smith 
419d91e6319SBarry Smith     Level: beginner
420d91e6319SBarry Smith 
4210a835dfdSSatish Balay    Any additions/changes here MUST also be made in include/finclude/petscmat.h
422d91e6319SBarry Smith 
423d91e6319SBarry Smith .seealso: MatSetOption()
424d91e6319SBarry Smith E*/
425512a5fc5SBarry Smith typedef enum {MAT_ROW_ORIENTED,MAT_NEW_NONZERO_LOCATIONS,
4264e0d8c25SBarry Smith               MAT_SYMMETRIC,
4274e0d8c25SBarry Smith               MAT_STRUCTURALLY_SYMMETRIC,
4288047e77bSBarry Smith               MAT_NEW_DIAGONALS,MAT_IGNORE_OFF_PROC_ENTRIES,
4294e0d8c25SBarry Smith               MAT_NEW_NONZERO_LOCATION_ERR,
4304e0d8c25SBarry Smith               MAT_NEW_NONZERO_ALLOCATION_ERR,MAT_USE_HASH_TABLE,
431a9817697SBarry Smith               MAT_KEEP_NONZERO_PATTERN,MAT_IGNORE_ZERO_ENTRIES,
4324e0d8c25SBarry Smith               MAT_USE_INODES,
4334e0d8c25SBarry Smith               MAT_HERMITIAN,
4344e0d8c25SBarry Smith               MAT_SYMMETRY_ETERNAL,
435cd6b891eSBarry Smith               MAT_CHECK_COMPRESSED_ROW,
4364e0d8c25SBarry Smith               MAT_IGNORE_LOWER_TRIANGULAR,MAT_ERROR_LOWER_TRIANGULAR,
43728b2fa4aSMatthew Knepley               MAT_GETROW_UPPERTRIANGULAR,MAT_UNUSED_NONZERO_LOCATION_ERR,
4384cb17eb5SBarry Smith               MAT_SPD,MAT_NO_OFF_PROC_ENTRIES,MAT_NO_OFF_PROC_ZERO_ROWS,
43928b2fa4aSMatthew Knepley               NUM_MAT_OPTIONS} MatOption;
440290bbb0aSBarry Smith extern const char *MatOptions[];
4417087cfbeSBarry Smith extern PetscErrorCode  MatSetOption(Mat,MatOption,PetscBool );
4427087cfbeSBarry Smith extern PetscErrorCode  MatGetType(Mat,const MatType*);
443a313700dSBarry Smith PetscPolymorphicFunction(MatGetType,(Mat mat),(mat,&t),const MatType,t)
44484cb2905SBarry Smith 
4457087cfbeSBarry Smith extern PetscErrorCode  MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]);
4467087cfbeSBarry Smith extern PetscErrorCode  MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
4477087cfbeSBarry Smith extern PetscErrorCode  MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
4487087cfbeSBarry Smith extern PetscErrorCode  MatGetRowUpperTriangular(Mat);
4497087cfbeSBarry Smith extern PetscErrorCode  MatRestoreRowUpperTriangular(Mat);
4507087cfbeSBarry Smith extern PetscErrorCode  MatGetColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
4517087cfbeSBarry Smith extern PetscErrorCode  MatRestoreColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
4527087cfbeSBarry Smith extern PetscErrorCode  MatGetColumnVector(Mat,Vec,PetscInt);
4537087cfbeSBarry Smith extern PetscErrorCode  MatGetArray(Mat,PetscScalar *[]);
454ba966639SSatish Balay PetscPolymorphicFunction(MatGetArray,(Mat mat),(mat,&a),PetscScalar*,a)
4557087cfbeSBarry Smith extern PetscErrorCode  MatRestoreArray(Mat,PetscScalar *[]);
4567087cfbeSBarry Smith extern PetscErrorCode  MatGetBlockSize(Mat,PetscInt *);
457ba966639SSatish Balay PetscPolymorphicFunction(MatGetBlockSize,(Mat mat),(mat,&a),PetscInt,a)
4587087cfbeSBarry Smith extern PetscErrorCode  MatSetBlockSize(Mat,PetscInt);
4597b80b807SBarry Smith 
4601620fd73SBarry Smith 
4617087cfbeSBarry Smith extern PetscErrorCode  MatMult(Mat,Vec,Vec);
4627087cfbeSBarry Smith extern PetscErrorCode  MatMultDiagonalBlock(Mat,Vec,Vec);
4637087cfbeSBarry Smith extern PetscErrorCode  MatMultAdd(Mat,Vec,Vec,Vec);
464ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
4657087cfbeSBarry Smith extern PetscErrorCode  MatMultTranspose(Mat,Vec,Vec);
4667087cfbeSBarry Smith extern PetscErrorCode  MatMultHermitianTranspose(Mat,Vec,Vec);
4677087cfbeSBarry Smith extern PetscErrorCode  MatIsTranspose(Mat,Mat,PetscReal,PetscBool *);
468ace3abfcSBarry Smith PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B,PetscReal tol),(A,B,tol,&t),PetscBool ,t)
469ace3abfcSBarry Smith PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B),(A,B,0,&t),PetscBool ,t)
4707087cfbeSBarry Smith extern PetscErrorCode  MatIsHermitianTranspose(Mat,Mat,PetscReal,PetscBool *);
4717087cfbeSBarry Smith extern PetscErrorCode  MatMultTransposeAdd(Mat,Vec,Vec,Vec);
4727087cfbeSBarry Smith extern PetscErrorCode  MatMultHermitianTransposeAdd(Mat,Vec,Vec,Vec);
473ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultTransposeAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
4747087cfbeSBarry Smith extern PetscErrorCode  MatMultConstrained(Mat,Vec,Vec);
4757087cfbeSBarry Smith extern PetscErrorCode  MatMultTransposeConstrained(Mat,Vec,Vec);
4767087cfbeSBarry Smith extern PetscErrorCode  MatMatSolve(Mat,Mat,Mat);
477f5edf698SHong Zhang 
478d91e6319SBarry Smith /*E
479d91e6319SBarry Smith     MatDuplicateOption - Indicates if a duplicated sparse matrix should have
480d91e6319SBarry Smith   its numerical values copied over or just its nonzero structure.
481d91e6319SBarry Smith 
482d91e6319SBarry Smith     Level: beginner
483d91e6319SBarry Smith 
484d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
485d91e6319SBarry Smith 
48670dcbbb9SBarry Smith $   MAT_SHARE_NONZERO_PATTERN - the i and j arrays in the new matrix will be shared with the original matrix
48770dcbbb9SBarry Smith $                               this also triggers the MAT_DO_NOT_COPY_VALUES option. This is used when you
48870dcbbb9SBarry Smith $                               have several matrices with the same nonzero pattern.
48970dcbbb9SBarry Smith 
490d91e6319SBarry Smith .seealso: MatDuplicate()
491d91e6319SBarry Smith E*/
49270dcbbb9SBarry Smith typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES,MAT_SHARE_NONZERO_PATTERN} MatDuplicateOption;
4932e8a6d31SBarry Smith 
4947087cfbeSBarry Smith extern PetscErrorCode  MatConvert(Mat,const MatType,MatReuse,Mat*);
495a313700dSBarry Smith PetscPolymorphicFunction(MatConvert,(Mat A,const MatType t),(A,t,MAT_INITIAL_MATRIX,&a),Mat,a)
4967087cfbeSBarry Smith extern PetscErrorCode  MatDuplicate(Mat,MatDuplicateOption,Mat*);
497ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A,MatDuplicateOption o),(A,o,&a),Mat,a)
498ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A),(A,MAT_COPY_VALUES,&a),Mat,a)
49994a9d846SBarry Smith 
500cb5b572fSBarry Smith 
5017087cfbeSBarry Smith extern PetscErrorCode  MatCopy(Mat,Mat,MatStructure);
5027087cfbeSBarry Smith extern PetscErrorCode  MatView(Mat,PetscViewer);
5037087cfbeSBarry Smith extern PetscErrorCode  MatIsSymmetric(Mat,PetscReal,PetscBool *);
504ace3abfcSBarry Smith PetscPolymorphicFunction(MatIsSymmetric,(Mat A,PetscReal tol),(A,tol,&t),PetscBool ,t)
505ace3abfcSBarry Smith PetscPolymorphicFunction(MatIsSymmetric,(Mat A),(A,0,&t),PetscBool ,t)
5067087cfbeSBarry Smith extern PetscErrorCode  MatIsStructurallySymmetric(Mat,PetscBool *);
507ace3abfcSBarry Smith PetscPolymorphicFunction(MatIsStructurallySymmetric,(Mat A),(A,&t),PetscBool ,t)
5087087cfbeSBarry Smith extern PetscErrorCode  MatIsHermitian(Mat,PetscReal,PetscBool *);
509ace3abfcSBarry Smith PetscPolymorphicFunction(MatIsHermitian,(Mat A),(A,0,&t),PetscBool ,t)
5107087cfbeSBarry Smith extern PetscErrorCode  MatIsSymmetricKnown(Mat,PetscBool *,PetscBool *);
5117087cfbeSBarry Smith extern PetscErrorCode  MatIsHermitianKnown(Mat,PetscBool *,PetscBool *);
5127087cfbeSBarry Smith extern PetscErrorCode  MatMissingDiagonal(Mat,PetscBool  *,PetscInt *);
5137087cfbeSBarry Smith extern PetscErrorCode  MatLoad(Mat, PetscViewer);
5147b80b807SBarry Smith 
5157087cfbeSBarry Smith extern PetscErrorCode  MatGetRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,PetscInt *[],PetscInt *[],PetscBool  *);
5167087cfbeSBarry Smith extern PetscErrorCode  MatRestoreRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,PetscInt *[],PetscInt *[],PetscBool  *);
5177087cfbeSBarry Smith extern PetscErrorCode  MatGetColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,PetscInt *[],PetscInt *[],PetscBool  *);
5187087cfbeSBarry Smith extern PetscErrorCode  MatRestoreColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,PetscInt *[],PetscInt *[],PetscBool  *);
519d4fbbf0eSBarry Smith 
520d91e6319SBarry Smith /*S
521d91e6319SBarry Smith      MatInfo - Context of matrix information, used with MatGetInfo()
522d91e6319SBarry Smith 
523d91e6319SBarry Smith    In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE
524d91e6319SBarry Smith 
525d91e6319SBarry Smith    Level: intermediate
526d91e6319SBarry Smith 
527d91e6319SBarry Smith   Concepts: matrix^nonzero information
528d91e6319SBarry Smith 
529d9274352SBarry Smith .seealso:  MatGetInfo(), MatInfoType
530d91e6319SBarry Smith S*/
5314e220ebcSLois Curfman McInnes typedef struct {
532b0a32e0cSBarry Smith   PetscLogDouble block_size;                         /* block size */
533b0a32e0cSBarry Smith   PetscLogDouble nz_allocated,nz_used,nz_unneeded;   /* number of nonzeros */
534b0a32e0cSBarry Smith   PetscLogDouble memory;                             /* memory allocated */
535b0a32e0cSBarry Smith   PetscLogDouble assemblies;                         /* number of matrix assemblies called */
536b0a32e0cSBarry Smith   PetscLogDouble mallocs;                            /* number of mallocs during MatSetValues() */
537b0a32e0cSBarry Smith   PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */
538b0a32e0cSBarry Smith   PetscLogDouble factor_mallocs;                     /* number of mallocs during factorization */
5394e220ebcSLois Curfman McInnes } MatInfo;
5404e220ebcSLois Curfman McInnes 
541d9274352SBarry Smith /*E
542d9274352SBarry Smith     MatInfoType - Indicates if you want information about the local part of the matrix,
543d9274352SBarry Smith      the entire parallel matrix or the maximum over all the local parts.
544d9274352SBarry Smith 
545d9274352SBarry Smith     Level: beginner
546d9274352SBarry Smith 
547d9274352SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
548d9274352SBarry Smith 
549d9274352SBarry Smith .seealso: MatGetInfo(), MatInfo
550d9274352SBarry Smith E*/
5517b80b807SBarry Smith typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType;
5527087cfbeSBarry Smith extern PetscErrorCode  MatGetInfo(Mat,MatInfoType,MatInfo*);
5537087cfbeSBarry Smith extern PetscErrorCode  MatGetDiagonal(Mat,Vec);
5547087cfbeSBarry Smith extern PetscErrorCode  MatGetRowMax(Mat,Vec,PetscInt[]);
5557087cfbeSBarry Smith extern PetscErrorCode  MatGetRowMin(Mat,Vec,PetscInt[]);
5567087cfbeSBarry Smith extern PetscErrorCode  MatGetRowMaxAbs(Mat,Vec,PetscInt[]);
5577087cfbeSBarry Smith extern PetscErrorCode  MatGetRowMinAbs(Mat,Vec,PetscInt[]);
5587087cfbeSBarry Smith extern PetscErrorCode  MatGetRowSum(Mat,Vec);
5597087cfbeSBarry Smith extern PetscErrorCode  MatTranspose(Mat,MatReuse,Mat*);
560fc4dec0aSBarry Smith PetscPolymorphicFunction(MatTranspose,(Mat A),(A,MAT_INITIAL_MATRIX,&t),Mat,t)
5617087cfbeSBarry Smith extern PetscErrorCode  MatHermitianTranspose(Mat,MatReuse,Mat*);
5627087cfbeSBarry Smith extern PetscErrorCode  MatPermute(Mat,IS,IS,Mat *);
563ba966639SSatish Balay PetscPolymorphicFunction(MatPermute,(Mat A,IS is1,IS is2),(A,is1,is2,&t),Mat,t)
5647087cfbeSBarry Smith extern PetscErrorCode  MatDiagonalScale(Mat,Vec,Vec);
5657087cfbeSBarry Smith extern PetscErrorCode  MatDiagonalSet(Mat,Vec,InsertMode);
5667087cfbeSBarry Smith extern PetscErrorCode  MatEqual(Mat,Mat,PetscBool *);
567ace3abfcSBarry Smith PetscPolymorphicFunction(MatEqual,(Mat A,Mat B),(A,B,&t),PetscBool ,t)
5687087cfbeSBarry Smith extern PetscErrorCode  MatMultEqual(Mat,Mat,PetscInt,PetscBool *);
5697087cfbeSBarry Smith extern PetscErrorCode  MatMultAddEqual(Mat,Mat,PetscInt,PetscBool *);
5707087cfbeSBarry Smith extern PetscErrorCode  MatMultTransposeEqual(Mat,Mat,PetscInt,PetscBool *);
5717087cfbeSBarry Smith extern PetscErrorCode  MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscBool *);
5727b80b807SBarry Smith 
5737087cfbeSBarry Smith extern PetscErrorCode  MatNorm(Mat,NormType,PetscReal *);
574ba966639SSatish Balay PetscPolymorphicFunction(MatNorm,(Mat A,NormType t),(A,t,&n),PetscReal,n)
57509573ac7SBarry Smith extern PetscErrorCode  MatGetColumnNorms(Mat,NormType,PetscReal *);
5767087cfbeSBarry Smith extern PetscErrorCode  MatZeroEntries(Mat);
5777087cfbeSBarry Smith extern PetscErrorCode  MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
5787087cfbeSBarry Smith extern PetscErrorCode  MatZeroRowsIS(Mat,IS,PetscScalar,Vec,Vec);
5797087cfbeSBarry Smith extern PetscErrorCode  MatZeroRowsStencil(Mat,PetscInt,const MatStencil [],PetscScalar,Vec,Vec);
580ac793be5SBarry Smith extern PetscErrorCode  MatZeroRowsColumnsStencil(Mat,PetscInt,const MatStencil[],PetscScalar,Vec,Vec);
5817087cfbeSBarry Smith extern PetscErrorCode  MatZeroRowsColumns(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
5827087cfbeSBarry Smith extern PetscErrorCode  MatZeroRowsColumnsIS(Mat,IS,PetscScalar,Vec,Vec);
5837b80b807SBarry Smith 
5847087cfbeSBarry Smith extern PetscErrorCode  MatUseScaledForm(Mat,PetscBool );
5857087cfbeSBarry Smith extern PetscErrorCode  MatScaleSystem(Mat,Vec,Vec);
5867087cfbeSBarry Smith extern PetscErrorCode  MatUnScaleSystem(Mat,Vec,Vec);
5875ef9f2a5SBarry Smith 
5887087cfbeSBarry Smith extern PetscErrorCode  MatGetSize(Mat,PetscInt*,PetscInt*);
5897087cfbeSBarry Smith extern PetscErrorCode  MatGetLocalSize(Mat,PetscInt*,PetscInt*);
5907087cfbeSBarry Smith extern PetscErrorCode  MatGetOwnershipRange(Mat,PetscInt*,PetscInt*);
5917087cfbeSBarry Smith extern PetscErrorCode  MatGetOwnershipRanges(Mat,const PetscInt**);
5927087cfbeSBarry Smith extern PetscErrorCode  MatGetOwnershipRangeColumn(Mat,PetscInt*,PetscInt*);
5937087cfbeSBarry Smith extern PetscErrorCode  MatGetOwnershipRangesColumn(Mat,const PetscInt**);
5947b80b807SBarry Smith 
5957087cfbeSBarry Smith extern PetscErrorCode  MatGetSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
5967087cfbeSBarry Smith extern PetscErrorCode  MatGetSubMatricesParallel(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
5977087cfbeSBarry Smith extern PetscErrorCode  MatDestroyMatrices(PetscInt,Mat *[]);
5987087cfbeSBarry Smith extern PetscErrorCode  MatGetSubMatrix(Mat,IS,IS,MatReuse,Mat *);
5997087cfbeSBarry Smith extern PetscErrorCode  MatGetLocalSubMatrix(Mat,IS,IS,Mat*);
6007087cfbeSBarry Smith extern PetscErrorCode  MatRestoreLocalSubMatrix(Mat,IS,IS,Mat*);
6017087cfbeSBarry Smith extern PetscErrorCode  MatGetSeqNonzeroStructure(Mat,Mat*);
6027087cfbeSBarry Smith extern PetscErrorCode  MatDestroySeqNonzeroStructure(Mat*);
6035494a064SHong Zhang 
6047087cfbeSBarry Smith extern PetscErrorCode  MatMerge(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
6059b8102ccSHong Zhang extern PetscErrorCode  MatMergeSymbolic(MPI_Comm,Mat,PetscInt,Mat*);
6069b8102ccSHong Zhang extern PetscErrorCode  MatMergeNumeric(MPI_Comm,Mat,PetscInt,Mat);
6077087cfbeSBarry Smith extern PetscErrorCode  MatMerge_SeqsToMPI(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*);
6087087cfbeSBarry Smith extern PetscErrorCode  MatMerge_SeqsToMPISymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*);
6097087cfbeSBarry Smith extern PetscErrorCode  MatMerge_SeqsToMPINumeric(Mat,Mat);
6104a2b5492SBarry Smith extern PetscErrorCode  MatMPIAIJGetLocalMat(Mat,MatReuse,Mat*);
6114a2b5492SBarry Smith extern PetscErrorCode  MatMPIAIJGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*);
61266bfb163SHong Zhang extern PetscErrorCode  MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,Mat*);
61343eb5e2fSMatthew Knepley #if defined (PETSC_USE_CTABLE)
6147087cfbeSBarry Smith extern PetscErrorCode  MatGetCommunicationStructs(Mat, Vec *, PetscTable *, VecScatter *);
61543eb5e2fSMatthew Knepley #else
6167087cfbeSBarry Smith extern PetscErrorCode  MatGetCommunicationStructs(Mat, Vec *, PetscInt *[], VecScatter *);
61743eb5e2fSMatthew Knepley #endif
6187087cfbeSBarry Smith extern PetscErrorCode  MatGetGhosts(Mat, PetscInt *,const PetscInt *[]);
6198efafbd8SBarry Smith 
6207087cfbeSBarry Smith extern PetscErrorCode  MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt);
6217b80b807SBarry Smith 
6227087cfbeSBarry Smith extern PetscErrorCode  MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*);
6237087cfbeSBarry Smith extern PetscErrorCode  MatMatMultSymbolic(Mat,Mat,PetscReal,Mat*);
6247087cfbeSBarry Smith extern PetscErrorCode  MatMatMultNumeric(Mat,Mat,Mat);
62522440eb1SKris Buschelman 
6267087cfbeSBarry Smith extern PetscErrorCode  MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*);
6277087cfbeSBarry Smith extern PetscErrorCode  MatPtAPSymbolic(Mat,Mat,PetscReal,Mat*);
6287087cfbeSBarry Smith extern PetscErrorCode  MatPtAPNumeric(Mat,Mat,Mat);
629*2b8ad9a3SHong Zhang extern PetscErrorCode  MatRARt(Mat,Mat,MatReuse,PetscReal,Mat*);
630*2b8ad9a3SHong Zhang extern PetscErrorCode  MatRARtSymbolic(Mat,Mat,PetscReal,Mat*);
631*2b8ad9a3SHong Zhang extern PetscErrorCode  MatRARtNumeric(Mat,Mat,Mat);
63222440eb1SKris Buschelman 
63375648e8dSHong Zhang extern PetscErrorCode  MatTransposeMatMult(Mat,Mat,MatReuse,PetscReal,Mat*);
63475648e8dSHong Zhang extern PetscErrorCode  MatTransposetMatMultSymbolic(Mat,Mat,PetscReal,Mat*);
63575648e8dSHong Zhang extern PetscErrorCode  MatTransposetMatMultNumeric(Mat,Mat,Mat);
6366fc122caSHong Zhang extern PetscErrorCode  MatMatTransposeMult(Mat,Mat,MatReuse,PetscReal,Mat*);
6376fc122caSHong Zhang extern PetscErrorCode  MatMatTransposeMultSymbolic(Mat,Mat,PetscReal,Mat*);
6386fc122caSHong Zhang extern PetscErrorCode  MatMatTransposeMultNumeric(Mat,Mat,Mat);
639bc011b1eSHong Zhang 
6407087cfbeSBarry Smith extern PetscErrorCode  MatAXPY(Mat,PetscScalar,Mat,MatStructure);
6417087cfbeSBarry Smith extern PetscErrorCode  MatAYPX(Mat,PetscScalar,Mat,MatStructure);
6421c741599SHong Zhang 
6437087cfbeSBarry Smith extern PetscErrorCode  MatScale(Mat,PetscScalar);
6447087cfbeSBarry Smith extern PetscErrorCode  MatShift(Mat,PetscScalar);
6457b80b807SBarry Smith 
6467087cfbeSBarry Smith extern PetscErrorCode  MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping);
6477087cfbeSBarry Smith extern PetscErrorCode  MatSetLocalToGlobalMappingBlock(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping);
6487087cfbeSBarry Smith extern PetscErrorCode  MatGetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*);
6497087cfbeSBarry Smith extern PetscErrorCode  MatGetLocalToGlobalMappingBlock(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*);
6507087cfbeSBarry Smith extern PetscErrorCode  MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
6517087cfbeSBarry Smith extern PetscErrorCode  MatZeroRowsLocalIS(Mat,IS,PetscScalar,Vec,Vec);
6527087cfbeSBarry Smith extern PetscErrorCode  MatZeroRowsColumnsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
6537087cfbeSBarry Smith extern PetscErrorCode  MatZeroRowsColumnsLocalIS(Mat,IS,PetscScalar,Vec,Vec);
6547087cfbeSBarry Smith extern PetscErrorCode  MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
6557087cfbeSBarry Smith extern PetscErrorCode  MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
656052efed2SBarry Smith 
6577087cfbeSBarry Smith extern PetscErrorCode  MatStashSetInitialSize(Mat,PetscInt,PetscInt);
6587087cfbeSBarry Smith extern PetscErrorCode  MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
65990f02eecSBarry Smith 
6607087cfbeSBarry Smith extern PetscErrorCode  MatInterpolate(Mat,Vec,Vec);
6617087cfbeSBarry Smith extern PetscErrorCode  MatInterpolateAdd(Mat,Vec,Vec,Vec);
662ba966639SSatish Balay PetscPolymorphicSubroutine(MatInterpolateAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
6637087cfbeSBarry Smith extern PetscErrorCode  MatRestrict(Mat,Vec,Vec);
6647087cfbeSBarry Smith extern PetscErrorCode  MatGetVecs(Mat,Vec*,Vec*);
6657087cfbeSBarry Smith extern PetscErrorCode  MatGetRedundantMatrix(Mat,PetscInt,MPI_Comm,PetscInt,MatReuse,Mat*);
6667087cfbeSBarry Smith extern PetscErrorCode  MatGetMultiProcBlock(Mat,MPI_Comm,Mat*);
667fc9bc008SSatish Balay extern PetscErrorCode  MatFindZeroDiagonals(Mat,IS*);
668bd481603SBarry Smith 
669bd481603SBarry Smith /*MC
6701620fd73SBarry Smith    MatSetValue - Set a single entry into a matrix.
6711620fd73SBarry Smith 
6721620fd73SBarry Smith    Not collective
6731620fd73SBarry Smith 
6741620fd73SBarry Smith    Input Parameters:
6751620fd73SBarry Smith +  m - the matrix
6761620fd73SBarry Smith .  row - the row location of the entry
6771620fd73SBarry Smith .  col - the column location of the entry
6781620fd73SBarry Smith .  value - the value to insert
6791620fd73SBarry Smith -  mode - either INSERT_VALUES or ADD_VALUES
6801620fd73SBarry Smith 
6811620fd73SBarry Smith    Notes:
6821620fd73SBarry Smith    For efficiency one should use MatSetValues() and set several or many
6831620fd73SBarry Smith    values simultaneously if possible.
6841620fd73SBarry Smith 
6851620fd73SBarry Smith    Level: beginner
6861620fd73SBarry Smith 
6871620fd73SBarry Smith .seealso: MatSetValues(), MatSetValueLocal()
6881620fd73SBarry Smith M*/
6891620fd73SBarry 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);}
6901620fd73SBarry Smith 
6911620fd73SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatGetValue(Mat v,PetscInt i,PetscInt j,PetscScalar *va) {return MatGetValues(v,1,&i,1,&j,va);}
6921620fd73SBarry Smith 
6931620fd73SBarry 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);}
6941620fd73SBarry Smith 
6951620fd73SBarry Smith /*MC
6960d2a7ff2SSatish Balay    MatPreallocateInitialize - Begins the block of code that will count the number of nonzeros per
697bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
698bd481603SBarry Smith 
699bd481603SBarry Smith    Synopsis:
700c1ac3661SBarry Smith    PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)
701bd481603SBarry Smith 
702bd481603SBarry Smith    Collective on MPI_Comm
703bd481603SBarry Smith 
704bd481603SBarry Smith    Input Parameters:
705bd481603SBarry Smith +  comm - the communicator that will share the eventually allocated matrix
706859fcb39SBarry Smith .  nrows - the number of LOCAL rows in the matrix
707859fcb39SBarry Smith -  ncols - the number of LOCAL columns in the matrix
708bd481603SBarry Smith 
709bd481603SBarry Smith    Output Parameters:
710bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
711bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
712bd481603SBarry Smith 
713bd481603SBarry Smith 
714bd481603SBarry Smith    Level: intermediate
715bd481603SBarry Smith 
716bd481603SBarry Smith    Notes:
7170598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
718bd481603SBarry Smith 
7191d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
720bd481603SBarry Smith 
721bd481603SBarry Smith    Use MatPreallocateInitializeSymmetric() for symmetric matrices (MPISBAIJ matrices)
722bd481603SBarry Smith 
7231620fd73SBarry Smith    This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize().
7241620fd73SBarry Smith 
725bd481603SBarry Smith   Concepts: preallocation^Matrix
726bd481603SBarry Smith 
727bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
728bd481603SBarry Smith           MatPreallocateInitializeSymmetric(), MatPreallocateSymmetricSetLocal()
729bd481603SBarry Smith M*/
730c4f061fbSSatish Balay #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \
7317c922b88SBarry Smith { \
732a89bc211SBarry Smith   PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__start,__end; \
733a89bc211SBarry Smith   _4_ierr = PetscMalloc2(__nrows,PetscInt,&dnz,__nrows,PetscInt,&onz);CHKERRQ(_4_ierr); \
734a89bc211SBarry Smith   _4_ierr = PetscMemzero(dnz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
735a89bc211SBarry Smith   _4_ierr = PetscMemzero(onz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
736c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ctmp;\
737a89bc211SBarry Smith   _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows;
7387c922b88SBarry Smith 
739bd481603SBarry Smith /*MC
7400d2a7ff2SSatish Balay    MatPreallocateSymmetricInitialize - Begins the block of code that will count the number of nonzeros per
741bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
742bd481603SBarry Smith 
743bd481603SBarry Smith    Synopsis:
744c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)
745bd481603SBarry Smith 
746bd481603SBarry Smith    Collective on MPI_Comm
747bd481603SBarry Smith 
748bd481603SBarry Smith    Input Parameters:
749bd481603SBarry Smith +  comm - the communicator that will share the eventually allocated matrix
750859fcb39SBarry Smith .  nrows - the number of LOCAL rows in the matrix
751859fcb39SBarry Smith -  ncols - the number of LOCAL columns in the matrix
752bd481603SBarry Smith 
753bd481603SBarry Smith    Output Parameters:
754bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
755bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
756bd481603SBarry Smith 
757bd481603SBarry Smith 
758bd481603SBarry Smith    Level: intermediate
759bd481603SBarry Smith 
760bd481603SBarry Smith    Notes:
7610598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
762bd481603SBarry Smith 
7631d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
764bd481603SBarry Smith 
7651620fd73SBarry Smith    This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize().
7661620fd73SBarry Smith 
767bd481603SBarry Smith   Concepts: preallocation^Matrix
768bd481603SBarry Smith 
769bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
770bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal()
771bd481603SBarry Smith M*/
772222b16d4SBarry Smith #define MatPreallocateSymmetricInitialize(comm,nrows,ncols,dnz,onz) 0; \
773222b16d4SBarry Smith { \
774a89bc211SBarry Smith   PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__end; \
775a89bc211SBarry Smith   _4_ierr = PetscMalloc2(__nrows,PetscInt,&dnz,__nrows,PetscInt,&onz);CHKERRQ(_4_ierr); \
776a89bc211SBarry Smith   _4_ierr = PetscMemzero(dnz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
777a89bc211SBarry Smith   _4_ierr = PetscMemzero(onz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
778c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr);\
779a89bc211SBarry Smith   _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows;
780222b16d4SBarry Smith 
781bd481603SBarry Smith /*MC
782bd481603SBarry Smith    MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
783bd481603SBarry Smith        inserted using a local number of the rows and columns
784bd481603SBarry Smith 
785bd481603SBarry Smith    Synopsis:
786c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
787bd481603SBarry Smith 
788bd481603SBarry Smith    Not Collective
789bd481603SBarry Smith 
790bd481603SBarry Smith    Input Parameters:
791784ac674SJed Brown +  map - the row mapping from local numbering to global numbering
792bd481603SBarry Smith .  nrows - the number of rows indicated
7931d73ed98SBarry Smith .  rows - the indices of the rows
794784ac674SJed Brown .  cmap - the column mapping from local to global numbering
795bd481603SBarry Smith .  ncols - the number of columns in the matrix
796bd481603SBarry Smith .  cols - the columns indicated
797bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
798bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
799bd481603SBarry Smith 
800bd481603SBarry Smith 
801bd481603SBarry Smith    Level: intermediate
802bd481603SBarry Smith 
803bd481603SBarry Smith    Notes:
8040598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
805bd481603SBarry Smith 
8061d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
807bd481603SBarry Smith 
808bd481603SBarry Smith   Concepts: preallocation^Matrix
809bd481603SBarry Smith 
810bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
811bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal()
812bd481603SBarry Smith M*/
813784ac674SJed Brown #define MatPreallocateSetLocal(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \
814c4f061fbSSatish Balay {\
815c1ac3661SBarry Smith   PetscInt __l;\
816784ac674SJed Brown   _4_ierr = ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\
817784ac674SJed Brown   _4_ierr = ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\
818c4f061fbSSatish Balay   for (__l=0;__l<nrows;__l++) {\
819ef66eb69SBarry Smith     _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
820c4f061fbSSatish Balay   }\
821c4f061fbSSatish Balay }
822c4f061fbSSatish Balay 
823bd481603SBarry Smith /*MC
824bd481603SBarry Smith    MatPreallocateSymmetricSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
825bd481603SBarry Smith        inserted using a local number of the rows and columns
826bd481603SBarry Smith 
827bd481603SBarry Smith    Synopsis:
828c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
829bd481603SBarry Smith 
830bd481603SBarry Smith    Not Collective
831bd481603SBarry Smith 
832bd481603SBarry Smith    Input Parameters:
833bd481603SBarry Smith +  map - the mapping between local numbering and global numbering
834bd481603SBarry Smith .  nrows - the number of rows indicated
8351d73ed98SBarry Smith .  rows - the indices of the rows
836bd481603SBarry Smith .  ncols - the number of columns in the matrix
837bd481603SBarry Smith .  cols - the columns indicated
838bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
839bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
840bd481603SBarry Smith 
841bd481603SBarry Smith 
842bd481603SBarry Smith    Level: intermediate
843bd481603SBarry Smith 
844bd481603SBarry Smith    Notes:
8450598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
846bd481603SBarry Smith 
847bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
848bd481603SBarry Smith 
849bd481603SBarry Smith   Concepts: preallocation^Matrix
850bd481603SBarry Smith 
851bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
852bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
853bd481603SBarry Smith M*/
854d3d32019SBarry Smith #define MatPreallocateSymmetricSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\
855d3d32019SBarry Smith {\
856c1ac3661SBarry Smith   PetscInt __l;\
857d3d32019SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\
858d3d32019SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\
859d3d32019SBarry Smith   for (__l=0;__l<nrows;__l++) {\
860d3d32019SBarry Smith     _4_ierr = MatPreallocateSymmetricSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
861d3d32019SBarry Smith   }\
862d3d32019SBarry Smith }
863d3d32019SBarry Smith 
864bd481603SBarry Smith /*MC
865bd481603SBarry Smith    MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
866bd481603SBarry Smith        inserted using a local number of the rows and columns
867bd481603SBarry Smith 
868bd481603SBarry Smith    Synopsis:
869c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
870bd481603SBarry Smith 
871bd481603SBarry Smith    Not Collective
872bd481603SBarry Smith 
873bd481603SBarry Smith    Input Parameters:
87464075487SBarry Smith +  row - the row
875bd481603SBarry Smith .  ncols - the number of columns in the matrix
876a50f8bf6SHong Zhang -  cols - the columns indicated
877a50f8bf6SHong Zhang 
878a50f8bf6SHong Zhang    Output Parameters:
879a50f8bf6SHong Zhang +  dnz - the array that will be passed to the matrix preallocation routines
880bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
881bd481603SBarry Smith 
882bd481603SBarry Smith 
883bd481603SBarry Smith    Level: intermediate
884bd481603SBarry Smith 
885bd481603SBarry Smith    Notes:
8860598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
887bd481603SBarry Smith 
888bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
889bd481603SBarry Smith 
8901620fd73SBarry Smith    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
8911620fd73SBarry Smith 
892bd481603SBarry Smith   Concepts: preallocation^Matrix
893bd481603SBarry Smith 
894bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
895bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
896bd481603SBarry Smith M*/
897c4f061fbSSatish Balay #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\
898c1ac3661SBarry Smith { PetscInt __i; \
899e32f2f54SBarry 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);\
900e32f2f54SBarry 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);\
9017c922b88SBarry Smith   for (__i=0; __i<nc; __i++) {\
90264075487SBarry Smith     if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \
9037cc688d7SBarry Smith     else dnz[row - __rstart]++;\
9047c922b88SBarry Smith   }\
9057c922b88SBarry Smith }
9067c922b88SBarry Smith 
907bd481603SBarry Smith /*MC
908bd481603SBarry Smith    MatPreallocateSymmetricSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
909bd481603SBarry Smith        inserted using a local number of the rows and columns
910bd481603SBarry Smith 
911bd481603SBarry Smith    Synopsis:
912c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
913bd481603SBarry Smith 
914bd481603SBarry Smith    Not Collective
915bd481603SBarry Smith 
916bd481603SBarry Smith    Input Parameters:
917bd481603SBarry Smith +  nrows - the number of rows indicated
9181d73ed98SBarry Smith .  rows - the indices of the rows
919bd481603SBarry Smith .  ncols - the number of columns in the matrix
920bd481603SBarry Smith .  cols - the columns indicated
921bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
922bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
923bd481603SBarry Smith 
924bd481603SBarry Smith 
925bd481603SBarry Smith    Level: intermediate
926bd481603SBarry Smith 
927bd481603SBarry Smith    Notes:
9280598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
929bd481603SBarry Smith 
930bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
931bd481603SBarry Smith 
9321620fd73SBarry Smith    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
9331620fd73SBarry Smith 
934bd481603SBarry Smith   Concepts: preallocation^Matrix
935bd481603SBarry Smith 
936bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
937bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
938bd481603SBarry Smith M*/
939d3d32019SBarry Smith #define MatPreallocateSymmetricSet(row,nc,cols,dnz,onz) 0;\
940c1ac3661SBarry Smith { PetscInt __i; \
941d3d32019SBarry Smith   for (__i=0; __i<nc; __i++) {\
942d3d32019SBarry Smith     if (cols[__i] >= __end) onz[row - __rstart]++; \
943d3d32019SBarry Smith     else if (cols[__i] >= row) dnz[row - __rstart]++;\
944d3d32019SBarry Smith   }\
945d3d32019SBarry Smith }
946d3d32019SBarry Smith 
947bd481603SBarry Smith /*MC
94816371a99SBarry Smith    MatPreallocateLocation -  An alternative to MatPreallocationSet() that puts the nonzero locations into the matrix if it exists
94916371a99SBarry Smith 
95016371a99SBarry Smith    Synopsis:
95116371a99SBarry Smith    PetscErrorCode MatPreallocateLocations(Mat A,PetscInt row,PetscInt ncols,PetscInt *cols,PetscInt *dnz,PetscInt *onz)
95216371a99SBarry Smith 
95316371a99SBarry Smith    Not Collective
95416371a99SBarry Smith 
95516371a99SBarry Smith    Input Parameters:
95616371a99SBarry Smith .  A - matrix
95716371a99SBarry Smith .  row - row where values exist (must be local to this process)
95816371a99SBarry Smith .  ncols - number of columns
95916371a99SBarry Smith .  cols - columns with nonzeros
96016371a99SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
96116371a99SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
96216371a99SBarry Smith 
96316371a99SBarry Smith 
96416371a99SBarry Smith    Level: intermediate
96516371a99SBarry Smith 
96616371a99SBarry Smith    Notes:
9670598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
96816371a99SBarry Smith 
96916371a99SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
97016371a99SBarry Smith 
97116371a99SBarry Smith    This is a MACRO not a function because it uses a bunch of variables private to the MatPreallocation.... routines.
97216371a99SBarry Smith 
97316371a99SBarry Smith   Concepts: preallocation^Matrix
97416371a99SBarry Smith 
97516371a99SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
97616371a99SBarry Smith           MatPreallocateSymmetricInitialize(), MatPreallocateSymmetricSetLocal()
97716371a99SBarry Smith M*/
97816371a99SBarry Smith #define MatPreallocateLocation(A,row,ncols,cols,dnz,onz) 0;if (A) {ierr = MatSetValues(A,1,&row,ncols,cols,PETSC_NULL,INSERT_VALUES);CHKERRQ(ierr);} else {ierr =  MatPreallocateSet(row,ncols,cols,dnz,onz);CHKERRQ(ierr);}
97916371a99SBarry Smith 
98016371a99SBarry Smith 
98116371a99SBarry Smith /*MC
9820d2a7ff2SSatish Balay    MatPreallocateFinalize - Ends the block of code that will count the number of nonzeros per
983bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
984bd481603SBarry Smith 
985bd481603SBarry Smith    Synopsis:
986c1ac3661SBarry Smith    PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz)
987bd481603SBarry Smith 
988bd481603SBarry Smith    Collective on MPI_Comm
989bd481603SBarry Smith 
990bd481603SBarry Smith    Input Parameters:
99116371a99SBarry Smith +  dnz - the array that was be passed to the matrix preallocation routines
992bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
993bd481603SBarry Smith 
994bd481603SBarry Smith 
995bd481603SBarry Smith    Level: intermediate
996bd481603SBarry Smith 
997bd481603SBarry Smith    Notes:
9980598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
999bd481603SBarry Smith 
1000bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
1001bd481603SBarry Smith 
10021620fd73SBarry Smith    This is a MACRO not a function because it closes the { started in MatPreallocateInitialize().
10031620fd73SBarry Smith 
1004bd481603SBarry Smith   Concepts: preallocation^Matrix
1005bd481603SBarry Smith 
1006bd481603SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
1007bd481603SBarry Smith           MatPreallocateSymmetricInitialize(), MatPreallocateSymmetricSetLocal()
1008bd481603SBarry Smith M*/
1009a89bc211SBarry Smith #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree2(dnz,onz);CHKERRQ(_4_ierr);}
10107c922b88SBarry Smith 
1011bd481603SBarry Smith 
1012bd481603SBarry Smith 
10137b80b807SBarry Smith /* Routines unique to particular data structures */
10148fe8eb27SJed Brown extern PetscErrorCode  MatShellGetContext(Mat,void *);
1015ba966639SSatish Balay PetscPolymorphicFunction(MatShellGetContext,(Mat A),(A,&t),void*,t)
1016698d4c6aSKris Buschelman 
10177087cfbeSBarry Smith extern PetscErrorCode  MatInodeAdjustForInodes(Mat,IS*,IS*);
10187087cfbeSBarry Smith extern PetscErrorCode  MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *);
1019698d4c6aSKris Buschelman 
10207087cfbeSBarry Smith extern PetscErrorCode  MatSeqAIJSetColumnIndices(Mat,PetscInt[]);
10217087cfbeSBarry Smith extern PetscErrorCode  MatSeqBAIJSetColumnIndices(Mat,PetscInt[]);
10227087cfbeSBarry Smith extern PetscErrorCode  MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
10237087cfbeSBarry Smith extern PetscErrorCode  MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
10247087cfbeSBarry Smith extern PetscErrorCode  MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
10251230e6d1SVictor Minden extern PetscErrorCode MatCreateSeqAIJFromTriple(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*,PetscInt,PetscBool);
10267b80b807SBarry Smith 
1027a96a251dSBarry Smith #define MAT_SKIP_ALLOCATION -4
1028a96a251dSBarry Smith 
10297087cfbeSBarry Smith extern PetscErrorCode  MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
1030ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz))
10317087cfbeSBarry Smith extern PetscErrorCode  MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
1032ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqSBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz))
10337087cfbeSBarry Smith extern PetscErrorCode  MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]);
1034ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqAIJSetPreallocation,(Mat A,const PetscInt nnz[]),(A,0,nnz))
1035273d9f13SBarry Smith 
10367087cfbeSBarry Smith extern PetscErrorCode  MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
1037ba966639SSatish Balay PetscPolymorphicSubroutine(MatMPIBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[],const PetscInt onz[]),(A,bs,0,nnz,0,onz))
10387087cfbeSBarry Smith extern PetscErrorCode  MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
10397087cfbeSBarry Smith extern PetscErrorCode  MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
10407087cfbeSBarry Smith extern PetscErrorCode  MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []);
10417087cfbeSBarry Smith extern PetscErrorCode  MatSeqBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
10427087cfbeSBarry Smith extern PetscErrorCode  MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]);
10437087cfbeSBarry Smith extern PetscErrorCode  MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
10447087cfbeSBarry Smith extern PetscErrorCode  MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]);
10457087cfbeSBarry Smith extern PetscErrorCode  MatMPIDenseSetPreallocation(Mat,PetscScalar[]);
10467087cfbeSBarry Smith extern PetscErrorCode  MatSeqDenseSetPreallocation(Mat,PetscScalar[]);
10477087cfbeSBarry Smith extern PetscErrorCode  MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,PetscInt*[]);
10487087cfbeSBarry Smith extern PetscErrorCode  MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,PetscInt*[]);
10497087cfbeSBarry Smith extern PetscErrorCode  MatAdicSetLocalFunction(Mat,void (*)(void));
1050273d9f13SBarry Smith 
10517087cfbeSBarry Smith extern PetscErrorCode  MatSeqDenseSetLDA(Mat,PetscInt);
10527087cfbeSBarry Smith extern PetscErrorCode  MatDenseGetLocalMatrix(Mat,Mat*);
10531b807ce4Svictorle 
10547087cfbeSBarry Smith extern PetscErrorCode  MatStoreValues(Mat);
10557087cfbeSBarry Smith extern PetscErrorCode  MatRetrieveValues(Mat);
10562e8a6d31SBarry Smith 
10577087cfbeSBarry Smith extern PetscErrorCode  MatDAADSetCtx(Mat,void*);
10583a7fca6bSBarry Smith 
1059b3a44c85SBarry Smith extern PetscErrorCode  MatFindNonzeroRows(Mat,IS*);
10607b80b807SBarry Smith /*
10617b80b807SBarry Smith   These routines are not usually accessed directly, rather solving is
106294b7f48cSBarry Smith   done through the KSP and PC interfaces.
10637b80b807SBarry Smith */
10647b80b807SBarry Smith 
106576bdecfbSBarry Smith /*J
1066d9274352SBarry Smith     MatOrderingType - String with the name of a PETSc matrix ordering or the creation function
1067d9274352SBarry Smith        with an optional dynamic library name, for example
1068d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:orderingcreate()
1069d9274352SBarry Smith 
1070d9274352SBarry Smith    Level: beginner
1071d9274352SBarry Smith 
1072e5a9bf91SBarry Smith    Cannot use const because the PC objects manipulate the string
1073e5a9bf91SBarry Smith 
1074d9274352SBarry Smith .seealso: MatGetOrdering()
107576bdecfbSBarry Smith J*/
10763eaad2caSSatish Balay #define MatOrderingType char*
10772692d6eeSBarry Smith #define MATORDERINGNATURAL     "natural"
10782692d6eeSBarry Smith #define MATORDERINGND          "nd"
10792692d6eeSBarry Smith #define MATORDERING1WD         "1wd"
10802692d6eeSBarry Smith #define MATORDERINGRCM         "rcm"
10812692d6eeSBarry Smith #define MATORDERINGQMD         "qmd"
10822692d6eeSBarry Smith #define MATORDERINGROWLENGTH   "rowlength"
1083312e372aSBarry Smith #define MATORDERINGAMD         "amd"            /* only works if UMFPACK is installed with PETSc */
1084b12f92e5SBarry Smith 
10857087cfbeSBarry Smith extern PetscErrorCode  MatGetOrdering(Mat,const MatOrderingType,IS*,IS*);
10867087cfbeSBarry Smith extern PetscErrorCode  MatGetOrderingList(PetscFList *list);
10877087cfbeSBarry Smith extern PetscErrorCode  MatOrderingRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,const MatOrderingType,IS*,IS*));
108830de9b25SBarry Smith 
108930de9b25SBarry Smith /*MC
10901890ba74SBarry Smith    MatOrderingRegisterDynamic - Adds a new sparse matrix ordering to the matrix package.
109130de9b25SBarry Smith 
109230de9b25SBarry Smith    Synopsis:
10931890ba74SBarry Smith    PetscErrorCode MatOrderingRegisterDynamic(const char *name_ordering,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatOrdering))
109430de9b25SBarry Smith 
109530de9b25SBarry Smith    Not Collective
109630de9b25SBarry Smith 
109730de9b25SBarry Smith    Input Parameters:
10982692d6eeSBarry Smith +  sname - name of ordering (for example MATORDERINGND)
109930de9b25SBarry Smith .  path - location of library where creation routine is
110030de9b25SBarry Smith .  name - name of function that creates the ordering type,a string
110130de9b25SBarry Smith -  function - function pointer that creates the ordering
110230de9b25SBarry Smith 
110330de9b25SBarry Smith    Level: developer
110430de9b25SBarry Smith 
110530de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
110630de9b25SBarry Smith    is ignored.
110730de9b25SBarry Smith 
110830de9b25SBarry Smith    Sample usage:
110930de9b25SBarry Smith .vb
111030de9b25SBarry Smith    MatOrderingRegisterDynamic("my_order",/home/username/my_lib/lib/libO/solaris/mylib.a,
111130de9b25SBarry Smith                "MyOrder",MyOrder);
111230de9b25SBarry Smith .ve
111330de9b25SBarry Smith 
111430de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
111530de9b25SBarry Smith $     MatOrderingSetType(part,"my_order)
111630de9b25SBarry Smith    or at runtime via the option
11172401956bSBarry Smith $     -pc_factor_mat_ordering_type my_order
111830de9b25SBarry Smith 
1119ab901514SBarry Smith    ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values.
112030de9b25SBarry Smith 
112130de9b25SBarry Smith .keywords: matrix, ordering, register
112230de9b25SBarry Smith 
112330de9b25SBarry Smith .seealso: MatOrderingRegisterDestroy(), MatOrderingRegisterAll()
112430de9b25SBarry Smith M*/
1125aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1126f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,0)
1127b12f92e5SBarry Smith #else
1128f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,d)
1129b12f92e5SBarry Smith #endif
113030de9b25SBarry Smith 
11317087cfbeSBarry Smith extern PetscErrorCode  MatOrderingRegisterDestroy(void);
11327087cfbeSBarry Smith extern PetscErrorCode  MatOrderingRegisterAll(const char[]);
1133ace3abfcSBarry Smith extern PetscBool  MatOrderingRegisterAllCalled;
1134b0a32e0cSBarry Smith extern PetscFList MatOrderingList;
1135d4fbbf0eSBarry Smith 
11367087cfbeSBarry Smith extern PetscErrorCode  MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS);
1137a2ce50c7SBarry Smith 
1138d91e6319SBarry Smith /*S
1139d90ac83dSHong Zhang     MatFactorShiftType - Numeric Shift.
1140d90ac83dSHong Zhang 
1141d90ac83dSHong Zhang    Level: beginner
1142d90ac83dSHong Zhang 
1143d90ac83dSHong Zhang S*/
1144d90ac83dSHong Zhang typedef enum {MAT_SHIFT_NONE,MAT_SHIFT_NONZERO,MAT_SHIFT_POSITIVE_DEFINITE,MAT_SHIFT_INBLOCKS} MatFactorShiftType;
1145d90ac83dSHong Zhang extern const char *MatFactorShiftTypes[];
1146d90ac83dSHong Zhang 
1147d90ac83dSHong Zhang /*S
11482401956bSBarry Smith    MatFactorInfo - Data passed into the matrix factorization routines
1149b00f7748SHong Zhang 
115061cad860SBarry Smith    In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE, that is use
115161cad860SBarry Smith $     MatFactorInfo  info(MAT_FACTORINFO_SIZE)
1152b00f7748SHong Zhang 
115315e8a5b3SHong Zhang    Notes: These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC.
1154b00f7748SHong Zhang 
115561cad860SBarry Smith       You can use MatFactorInfoInitialize() to set default values.
115661cad860SBarry Smith 
1157b00f7748SHong Zhang    Level: developer
1158b00f7748SHong Zhang 
1159d7d82daaSBarry Smith .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(),
1160d7d82daaSBarry Smith           MatFactorInfoInitialize()
1161b00f7748SHong Zhang 
1162b00f7748SHong Zhang S*/
1163b00f7748SHong Zhang typedef struct {
116415e8a5b3SHong Zhang   PetscReal     diagonal_fill;  /* force diagonal to fill in if initially not filled */
116585317021SBarry Smith   PetscReal     usedt;
116615e8a5b3SHong Zhang   PetscReal     dt;             /* drop tolerance */
1167b00f7748SHong Zhang   PetscReal     dtcol;          /* tolerance for pivoting */
116815e8a5b3SHong Zhang   PetscReal     dtcount;        /* maximum nonzeros to be allowed per row */
116967e28bfeSBarry Smith   PetscReal     fill;           /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */
1170348344bbSBarry Smith   PetscReal     levels;         /* ICC/ILU(levels) */
1171bcd9e38bSBarry Smith   PetscReal     pivotinblocks;  /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0
1172bcd9e38bSBarry Smith                                    factorization may be faster if do not pivot */
117315e8a5b3SHong Zhang   PetscReal     zeropivot;      /* pivot is called zero if less than this */
1174f4db908eSBarry Smith   PetscReal     shifttype;      /* type of shift added to matrix factor to prevent zero pivots */
1175f4db908eSBarry Smith   PetscReal     shiftamount;     /* how large the shift is */
117615e8a5b3SHong Zhang } MatFactorInfo;
1177ffa6d0a5SLois Curfman McInnes 
11787087cfbeSBarry Smith extern PetscErrorCode  MatFactorInfoInitialize(MatFactorInfo*);
11797087cfbeSBarry Smith extern PetscErrorCode  MatCholeskyFactor(Mat,IS,const MatFactorInfo*);
11807087cfbeSBarry Smith extern PetscErrorCode  MatCholeskyFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
11817087cfbeSBarry Smith extern PetscErrorCode  MatCholeskyFactorNumeric(Mat,Mat,const MatFactorInfo*);
11827087cfbeSBarry Smith extern PetscErrorCode  MatLUFactor(Mat,IS,IS,const MatFactorInfo*);
11837087cfbeSBarry Smith extern PetscErrorCode  MatILUFactor(Mat,IS,IS,const MatFactorInfo*);
11847087cfbeSBarry Smith extern PetscErrorCode  MatLUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
11857087cfbeSBarry Smith extern PetscErrorCode  MatILUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
11867087cfbeSBarry Smith extern PetscErrorCode  MatICCFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
11877087cfbeSBarry Smith extern PetscErrorCode  MatICCFactor(Mat,IS,const MatFactorInfo*);
11887087cfbeSBarry Smith extern PetscErrorCode  MatLUFactorNumeric(Mat,Mat,const MatFactorInfo*);
11897087cfbeSBarry Smith extern PetscErrorCode  MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*);
11907087cfbeSBarry Smith extern PetscErrorCode  MatSolve(Mat,Vec,Vec);
11917087cfbeSBarry Smith extern PetscErrorCode  MatForwardSolve(Mat,Vec,Vec);
11927087cfbeSBarry Smith extern PetscErrorCode  MatBackwardSolve(Mat,Vec,Vec);
11937087cfbeSBarry Smith extern PetscErrorCode  MatSolveAdd(Mat,Vec,Vec,Vec);
11947087cfbeSBarry Smith extern PetscErrorCode  MatSolveTranspose(Mat,Vec,Vec);
11957087cfbeSBarry Smith extern PetscErrorCode  MatSolveTransposeAdd(Mat,Vec,Vec,Vec);
11967087cfbeSBarry Smith extern PetscErrorCode  MatSolves(Mat,Vecs,Vecs);
11978ed539a5SBarry Smith 
11987087cfbeSBarry Smith extern PetscErrorCode  MatSetUnfactored(Mat);
1199bb5a7306SBarry Smith 
1200d91e6319SBarry Smith /*E
1201d91e6319SBarry Smith     MatSORType - What type of (S)SOR to perform
1202bb1eb677SSatish Balay 
1203d91e6319SBarry Smith     Level: beginner
1204d91e6319SBarry Smith 
1205d9274352SBarry Smith    May be bitwise ORd together
1206d9274352SBarry Smith 
1207d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
1208d91e6319SBarry Smith 
12094e7234bfSBarry Smith    MatSORType may be bitwise ORd together, so do not change the numbers
12104e7234bfSBarry Smith 
121141f059aeSBarry Smith .seealso: MatSOR()
1212d91e6319SBarry Smith E*/
1213ee50ffe9SBarry Smith typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3,
1214ee50ffe9SBarry Smith               SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8,
1215ee50ffe9SBarry Smith               SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16,
121684cb2905SBarry Smith               SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType;
12177087cfbeSBarry Smith extern PetscErrorCode  MatSOR(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
12188ed539a5SBarry Smith 
1219d4fbbf0eSBarry Smith /*
1220639f9d9dSBarry Smith     These routines are for efficiently computing Jacobians via finite differences.
1221639f9d9dSBarry Smith */
1222b12f92e5SBarry Smith 
122376bdecfbSBarry Smith /*J
1224d9274352SBarry Smith     MatColoringType - String with the name of a PETSc matrix coloring or the creation function
1225d9274352SBarry Smith        with an optional dynamic library name, for example
1226d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:coloringcreate()
1227d9274352SBarry Smith 
1228d9274352SBarry Smith    Level: beginner
1229d9274352SBarry Smith 
1230d9274352SBarry Smith .seealso: MatGetColoring()
123176bdecfbSBarry Smith J*/
1232a313700dSBarry Smith #define MatColoringType char*
12332692d6eeSBarry Smith #define MATCOLORINGNATURAL "natural"
12342692d6eeSBarry Smith #define MATCOLORINGSL      "sl"
12352692d6eeSBarry Smith #define MATCOLORINGLF      "lf"
12362692d6eeSBarry Smith #define MATCOLORINGID      "id"
1237b12f92e5SBarry Smith 
12387087cfbeSBarry Smith extern PetscErrorCode  MatGetColoring(Mat,const MatColoringType,ISColoring*);
12397087cfbeSBarry Smith extern PetscErrorCode  MatColoringRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,MatColoringType,ISColoring *));
124030de9b25SBarry Smith 
124130de9b25SBarry Smith /*MC
124230de9b25SBarry Smith    MatColoringRegisterDynamic - Adds a new sparse matrix coloring to the
124330de9b25SBarry Smith                                matrix package.
124430de9b25SBarry Smith 
124530de9b25SBarry Smith    Synopsis:
12461890ba74SBarry Smith    PetscErrorCode MatColoringRegisterDynamic(const char *name_coloring,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatColoring))
124730de9b25SBarry Smith 
124830de9b25SBarry Smith    Not Collective
124930de9b25SBarry Smith 
125030de9b25SBarry Smith    Input Parameters:
12512692d6eeSBarry Smith +  sname - name of Coloring (for example MATCOLORINGSL)
125230de9b25SBarry Smith .  path - location of library where creation routine is
125330de9b25SBarry Smith .  name - name of function that creates the Coloring type, a string
125430de9b25SBarry Smith -  function - function pointer that creates the coloring
125530de9b25SBarry Smith 
125630de9b25SBarry Smith    Level: developer
125730de9b25SBarry Smith 
125830de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
125930de9b25SBarry Smith    is ignored.
126030de9b25SBarry Smith 
126130de9b25SBarry Smith    Sample usage:
126230de9b25SBarry Smith .vb
126330de9b25SBarry Smith    MatColoringRegisterDynamic("my_color",/home/username/my_lib/lib/libO/solaris/mylib.a,
126430de9b25SBarry Smith                "MyColor",MyColor);
126530de9b25SBarry Smith .ve
126630de9b25SBarry Smith 
126730de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
126830de9b25SBarry Smith $     MatColoringSetType(part,"my_color")
126930de9b25SBarry Smith    or at runtime via the option
127030de9b25SBarry Smith $     -mat_coloring_type my_color
127130de9b25SBarry Smith 
1272ab901514SBarry Smith    $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
127330de9b25SBarry Smith 
127430de9b25SBarry Smith .keywords: matrix, Coloring, register
127530de9b25SBarry Smith 
127630de9b25SBarry Smith .seealso: MatColoringRegisterDestroy(), MatColoringRegisterAll()
127730de9b25SBarry Smith M*/
1278aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1279f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,0)
1280b12f92e5SBarry Smith #else
1281f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,d)
1282b12f92e5SBarry Smith #endif
128330de9b25SBarry Smith 
1284ace3abfcSBarry Smith extern PetscBool  MatColoringRegisterAllCalled;
1285f1144a30SSatish Balay 
12867087cfbeSBarry Smith extern PetscErrorCode  MatColoringRegisterAll(const char[]);
12877087cfbeSBarry Smith extern PetscErrorCode  MatColoringRegisterDestroy(void);
12887087cfbeSBarry Smith extern PetscErrorCode  MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*);
1289639f9d9dSBarry Smith 
1290d9274352SBarry Smith /*S
1291d9274352SBarry Smith      MatFDColoring - Object for computing a sparse Jacobian via finite differences
1292d9274352SBarry Smith         and coloring
1293639f9d9dSBarry Smith 
1294d9274352SBarry Smith    Level: beginner
1295d9274352SBarry Smith 
1296d9274352SBarry Smith   Concepts: coloring, sparse Jacobian, finite differences
1297d9274352SBarry Smith 
1298d9274352SBarry Smith .seealso:  MatFDColoringCreate()
1299d9274352SBarry Smith S*/
1300e2a1c21fSSatish Balay typedef struct _p_MatFDColoring* MatFDColoring;
1301639f9d9dSBarry Smith 
13027087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringCreate(Mat,ISColoring,MatFDColoring *);
1303fcfd50ebSBarry Smith extern PetscErrorCode  MatFDColoringDestroy(MatFDColoring*);
13047087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringView(MatFDColoring,PetscViewer);
13057087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*);
13067087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**);
13077087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal);
13087087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringSetFromOptions(MatFDColoring);
13097087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringApply(Mat,MatFDColoring,Vec,MatStructure*,void *);
13107087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringSetF(MatFDColoring,Vec);
13117087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,PetscInt*[]);
1312b1683b59SHong Zhang 
1313b1683b59SHong Zhang /*S
1314b9af6bddSHong Zhang      MatTransposeColoring - Object for computing a sparse matrix product C=A*B^T via coloring
1315b1683b59SHong Zhang 
1316b1683b59SHong Zhang    Level: beginner
1317b1683b59SHong Zhang 
1318b1683b59SHong Zhang   Concepts: coloring, sparse matrix product
1319b1683b59SHong Zhang 
1320b9af6bddSHong Zhang .seealso:  MatTransposeColoringCreate()
1321b1683b59SHong Zhang S*/
1322b9af6bddSHong Zhang typedef struct _p_MatTransposeColoring* MatTransposeColoring;
1323b1683b59SHong Zhang 
1324b9af6bddSHong Zhang extern PetscErrorCode MatTransposeColoringCreate(Mat,ISColoring,MatTransposeColoring *);
1325b9af6bddSHong Zhang extern PetscErrorCode MatTransColoringApplySpToDen(MatTransposeColoring,Mat,Mat);
1326b9af6bddSHong Zhang extern PetscErrorCode MatTransColoringApplyDenToSp(MatTransposeColoring,Mat,Mat);
1327b9af6bddSHong Zhang extern PetscErrorCode MatTransposeColoringDestroy(MatTransposeColoring*);
1328b1683b59SHong Zhang 
1329639f9d9dSBarry Smith /*
13300752156aSBarry Smith     These routines are for partitioning matrices: currently used only
13313eda8832SBarry Smith   for adjacency matrix, MatCreateMPIAdj().
13320752156aSBarry Smith */
1333ca161407SBarry Smith 
1334d9274352SBarry Smith /*S
1335d9274352SBarry Smith      MatPartitioning - Object for managing the partitioning of a matrix or graph
1336d9274352SBarry Smith 
1337d9274352SBarry Smith    Level: beginner
1338d9274352SBarry Smith 
1339d9274352SBarry Smith   Concepts: partitioning
1340d9274352SBarry Smith 
1341743a1026Svictorle .seealso:  MatPartitioningCreate(), MatPartitioningType
1342d9274352SBarry Smith S*/
134391e9ee9fSBarry Smith typedef struct _p_MatPartitioning* MatPartitioning;
1344d9274352SBarry Smith 
134576bdecfbSBarry Smith /*J
13465bc6e7ccSvictorle     MatPartitioningType - String with the name of a PETSc matrix partitioning or the creation function
1347d9274352SBarry Smith        with an optional dynamic library name, for example
1348d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:partitioningcreate()
1349d9274352SBarry Smith 
1350d9274352SBarry Smith    Level: beginner
13510d04baf8SBarry Smith dm
1352b547a13bSvictorle .seealso: MatPartitioningCreate(), MatPartitioning
135376bdecfbSBarry Smith J*/
1354a313700dSBarry Smith #define MatPartitioningType char*
13552692d6eeSBarry Smith #define MATPARTITIONINGCURRENT  "current"
13562692d6eeSBarry Smith #define MATPARTITIONINGSQUARE   "square"
13572692d6eeSBarry Smith #define MATPARTITIONINGPARMETIS "parmetis"
13582692d6eeSBarry Smith #define MATPARTITIONINGCHACO    "chaco"
13592692d6eeSBarry Smith #define MATPARTITIONINGPARTY    "party"
136050d91057SBarry Smith #define MATPARTITIONINGPTSCOTCH "ptscotch"
136171306c60Sjroman 
1362ca161407SBarry Smith 
13637087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningCreate(MPI_Comm,MatPartitioning*);
13647087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningSetType(MatPartitioning,const MatPartitioningType);
13657087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningSetNParts(MatPartitioning,PetscInt);
13667087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningSetAdjacency(MatPartitioning,Mat);
13677087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]);
13687087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []);
13697087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningApply(MatPartitioning,IS*);
1370fcfd50ebSBarry Smith extern PetscErrorCode  MatPartitioningDestroy(MatPartitioning*);
13712aabb6bbSBarry Smith 
13727087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatPartitioning));
137330de9b25SBarry Smith 
137430de9b25SBarry Smith /*MC
137530de9b25SBarry Smith    MatPartitioningRegisterDynamic - Adds a new sparse matrix partitioning to the
137630de9b25SBarry Smith    matrix package.
137730de9b25SBarry Smith 
137830de9b25SBarry Smith    Synopsis:
13791890ba74SBarry Smith    PetscErrorCode MatPartitioningRegisterDynamic(const char *name_partitioning,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatPartitioning))
138030de9b25SBarry Smith 
138130de9b25SBarry Smith    Not Collective
138230de9b25SBarry Smith 
138330de9b25SBarry Smith    Input Parameters:
13842692d6eeSBarry Smith +  sname - name of partitioning (for example MATPARTITIONINGCURRENT) or parmetis
138530de9b25SBarry Smith .  path - location of library where creation routine is
138630de9b25SBarry Smith .  name - name of function that creates the partitioning type, a string
138730de9b25SBarry Smith -  function - function pointer that creates the partitioning type
138830de9b25SBarry Smith 
138930de9b25SBarry Smith    Level: developer
139030de9b25SBarry Smith 
139130de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
139230de9b25SBarry Smith    is ignored.
139330de9b25SBarry Smith 
139430de9b25SBarry Smith    Sample usage:
139530de9b25SBarry Smith .vb
139630de9b25SBarry Smith    MatPartitioningRegisterDynamic("my_part",/home/username/my_lib/lib/libO/solaris/mylib.a,
139730de9b25SBarry Smith                "MyPartCreate",MyPartCreate);
139830de9b25SBarry Smith .ve
139930de9b25SBarry Smith 
140030de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
140130de9b25SBarry Smith $     MatPartitioningSetType(part,"my_part")
140230de9b25SBarry Smith    or at runtime via the option
140330de9b25SBarry Smith $     -mat_partitioning_type my_part
140430de9b25SBarry Smith 
1405ab901514SBarry Smith    $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
140630de9b25SBarry Smith 
140730de9b25SBarry Smith .keywords: matrix, partitioning, register
140830de9b25SBarry Smith 
140930de9b25SBarry Smith .seealso: MatPartitioningRegisterDestroy(), MatPartitioningRegisterAll()
141030de9b25SBarry Smith M*/
1411aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1412f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,0)
14132aabb6bbSBarry Smith #else
1414f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,d)
14152aabb6bbSBarry Smith #endif
14162aabb6bbSBarry Smith 
1417ace3abfcSBarry Smith extern PetscBool  MatPartitioningRegisterAllCalled;
1418f1144a30SSatish Balay 
14197087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningRegisterAll(const char[]);
14207087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningRegisterDestroy(void);
14212bad1931SBarry Smith 
14227087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningView(MatPartitioning,PetscViewer);
14237087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningSetFromOptions(MatPartitioning);
14247087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningGetType(MatPartitioning,const MatPartitioningType*);
1425ca161407SBarry Smith 
14267087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningParmetisSetCoarseSequential(MatPartitioning);
14277087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningParmetisGetEdgeCut(MatPartitioning, PetscInt *);
14280752156aSBarry Smith 
1429b6956c12SJose Roman typedef enum { MP_CHACO_MULTILEVEL=1,MP_CHACO_SPECTRAL=2,MP_CHACO_LINEAR=4,MP_CHACO_RANDOM=5,MP_CHACO_SCATTERED=6 } MPChacoGlobalType;
1430b6956c12SJose Roman extern const char *MPChacoGlobalTypes[];
1431b6956c12SJose Roman typedef enum { MP_CHACO_KERNIGHAN=1,MP_CHACO_NONE=2 } MPChacoLocalType;
1432b6956c12SJose Roman extern const char *MPChacoLocalTypes[];
1433b6956c12SJose Roman typedef enum { MP_CHACO_LANCZOS=0,MP_CHACO_RQI=1 } MPChacoEigenType;
1434b6956c12SJose Roman extern const char *MPChacoEigenTypes[];
1435b6956c12SJose Roman 
14367087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningChacoSetGlobal(MatPartitioning,MPChacoGlobalType);
1437b6956c12SJose Roman extern PetscErrorCode  MatPartitioningChacoGetGlobal(MatPartitioning,MPChacoGlobalType*);
14387087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningChacoSetLocal(MatPartitioning,MPChacoLocalType);
1439b6956c12SJose Roman extern PetscErrorCode  MatPartitioningChacoGetLocal(MatPartitioning,MPChacoLocalType*);
14407087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal);
14417087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType);
1442b6956c12SJose Roman extern PetscErrorCode  MatPartitioningChacoGetEigenSolver(MatPartitioning,MPChacoEigenType*);
14437087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningChacoSetEigenTol(MatPartitioning,PetscReal);
1444b6956c12SJose Roman extern PetscErrorCode  MatPartitioningChacoGetEigenTol(MatPartitioning,PetscReal*);
14457087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningChacoSetEigenNumber(MatPartitioning,PetscInt);
1446b6956c12SJose Roman extern PetscErrorCode  MatPartitioningChacoGetEigenNumber(MatPartitioning,PetscInt*);
144771306c60Sjroman 
144871306c60Sjroman #define MP_PARTY_OPT "opt"
144971306c60Sjroman #define MP_PARTY_LIN "lin"
145071306c60Sjroman #define MP_PARTY_SCA "sca"
145171306c60Sjroman #define MP_PARTY_RAN "ran"
145271306c60Sjroman #define MP_PARTY_GBF "gbf"
145371306c60Sjroman #define MP_PARTY_GCF "gcf"
145471306c60Sjroman #define MP_PARTY_BUB "bub"
145571306c60Sjroman #define MP_PARTY_DEF "def"
14567087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningPartySetGlobal(MatPartitioning,const char*);
145771306c60Sjroman #define MP_PARTY_HELPFUL_SETS "hs"
145871306c60Sjroman #define MP_PARTY_KERNIGHAN_LIN "kl"
145971306c60Sjroman #define MP_PARTY_NONE "no"
14607087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningPartySetLocal(MatPartitioning,const char*);
14617087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal);
14627087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningPartySetBipart(MatPartitioning,PetscBool);
14637087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscBool);
146471306c60Sjroman 
146550d91057SBarry Smith typedef enum { MP_PTSCOTCH_QUALITY,MP_PTSCOTCH_SPEED,MP_PTSCOTCH_BALANCE,MP_PTSCOTCH_SAFETY,MP_PTSCOTCH_SCALABILITY } MPPTScotchStrategyType;
146650d91057SBarry Smith extern const char *MPPTScotchStrategyTypes[];
1467e0f1cffaSJose Roman 
146850d91057SBarry Smith extern PetscErrorCode MatPartitioningPTScotchSetImbalance(MatPartitioning,PetscReal);
146950d91057SBarry Smith extern PetscErrorCode MatPartitioningPTScotchGetImbalance(MatPartitioning,PetscReal*);
147050d91057SBarry Smith extern PetscErrorCode MatPartitioningPTScotchSetStrategy(MatPartitioning,MPPTScotchStrategyType);
147150d91057SBarry Smith extern PetscErrorCode MatPartitioningPTScotchGetStrategy(MatPartitioning,MPPTScotchStrategyType*);
147271306c60Sjroman 
147309573ac7SBarry Smith extern PetscErrorCode MatMeshToVertexGraph(Mat,PetscInt,Mat*);
147409573ac7SBarry Smith extern PetscErrorCode MatMeshToCellGraph(Mat,PetscInt,Mat*);
1475591294e4SBarry Smith 
14760752156aSBarry Smith /*
14770a835dfdSSatish Balay     If you add entries here you must also add them to finclude/petscmat.h
1478d4fbbf0eSBarry Smith */
14791c1c02c0SLois Curfman McInnes typedef enum { MATOP_SET_VALUES=0,
14801c1c02c0SLois Curfman McInnes                MATOP_GET_ROW=1,
14811c1c02c0SLois Curfman McInnes                MATOP_RESTORE_ROW=2,
14821c1c02c0SLois Curfman McInnes                MATOP_MULT=3,
14831c1c02c0SLois Curfman McInnes                MATOP_MULT_ADD=4,
14847c922b88SBarry Smith                MATOP_MULT_TRANSPOSE=5,
14857c922b88SBarry Smith                MATOP_MULT_TRANSPOSE_ADD=6,
14861c1c02c0SLois Curfman McInnes                MATOP_SOLVE=7,
14871c1c02c0SLois Curfman McInnes                MATOP_SOLVE_ADD=8,
14887c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE=9,
14897c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE_ADD=10,
14901c1c02c0SLois Curfman McInnes                MATOP_LUFACTOR=11,
14911c1c02c0SLois Curfman McInnes                MATOP_CHOLESKYFACTOR=12,
1492a714c06dSBarry Smith                MATOP_SOR=13,
14931c1c02c0SLois Curfman McInnes                MATOP_TRANSPOSE=14,
14941c1c02c0SLois Curfman McInnes                MATOP_GETINFO=15,
14951c1c02c0SLois Curfman McInnes                MATOP_EQUAL=16,
14961c1c02c0SLois Curfman McInnes                MATOP_GET_DIAGONAL=17,
14971c1c02c0SLois Curfman McInnes                MATOP_DIAGONAL_SCALE=18,
14981c1c02c0SLois Curfman McInnes                MATOP_NORM=19,
14991c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_BEGIN=20,
15001c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_END=21,
1501d519adbfSMatthew Knepley                MATOP_SET_OPTION=22,
1502d519adbfSMatthew Knepley                MATOP_ZERO_ENTRIES=23,
1503d519adbfSMatthew Knepley                MATOP_ZERO_ROWS=24,
1504d519adbfSMatthew Knepley                MATOP_LUFACTOR_SYMBOLIC=25,
1505d519adbfSMatthew Knepley                MATOP_LUFACTOR_NUMERIC=26,
1506d519adbfSMatthew Knepley                MATOP_CHOLESKY_FACTOR_SYMBOLIC=27,
1507d519adbfSMatthew Knepley                MATOP_CHOLESKY_FACTOR_NUMERIC=28,
1508d519adbfSMatthew Knepley                MATOP_SETUP_PREALLOCATION=29,
1509d519adbfSMatthew Knepley                MATOP_ILUFACTOR_SYMBOLIC=30,
1510d519adbfSMatthew Knepley                MATOP_ICCFACTOR_SYMBOLIC=31,
1511d519adbfSMatthew Knepley                MATOP_GET_ARRAY=32,
1512d519adbfSMatthew Knepley                MATOP_RESTORE_ARRAY=33,
1513d519adbfSMatthew Knepley                MATOP_DUPLICATE=34,
1514d519adbfSMatthew Knepley                MATOP_FORWARD_SOLVE=35,
1515d519adbfSMatthew Knepley                MATOP_BACKWARD_SOLVE=36,
1516d519adbfSMatthew Knepley                MATOP_ILUFACTOR=37,
1517d519adbfSMatthew Knepley                MATOP_ICCFACTOR=38,
1518d519adbfSMatthew Knepley                MATOP_AXPY=39,
1519d519adbfSMatthew Knepley                MATOP_GET_SUBMATRICES=40,
1520d519adbfSMatthew Knepley                MATOP_INCREASE_OVERLAP=41,
1521d519adbfSMatthew Knepley                MATOP_GET_VALUES=42,
1522d519adbfSMatthew Knepley                MATOP_COPY=43,
1523d519adbfSMatthew Knepley                MATOP_GET_ROW_MAX=44,
1524d519adbfSMatthew Knepley                MATOP_SCALE=45,
1525d519adbfSMatthew Knepley                MATOP_SHIFT=46,
152635153367SBarry Smith                MATOP_DIAGONAL_SET=47,
1527d519adbfSMatthew Knepley                MATOP_ILUDT_FACTOR=48,
1528d519adbfSMatthew Knepley                MATOP_SET_BLOCK_SIZE=49,
1529d519adbfSMatthew Knepley                MATOP_GET_ROW_IJ=50,
1530d519adbfSMatthew Knepley                MATOP_RESTORE_ROW_IJ=51,
1531d519adbfSMatthew Knepley                MATOP_GET_COLUMN_IJ=52,
1532d519adbfSMatthew Knepley                MATOP_RESTORE_COLUMN_IJ=53,
1533d519adbfSMatthew Knepley                MATOP_FDCOLORING_CREATE=54,
1534d519adbfSMatthew Knepley                MATOP_COLORING_PATCH=55,
1535d519adbfSMatthew Knepley                MATOP_SET_UNFACTORED=56,
1536d519adbfSMatthew Knepley                MATOP_PERMUTE=57,
1537d519adbfSMatthew Knepley                MATOP_SET_VALUES_BLOCKED=58,
1538d519adbfSMatthew Knepley                MATOP_GET_SUBMATRIX=59,
1539d519adbfSMatthew Knepley                MATOP_DESTROY=60,
1540d519adbfSMatthew Knepley                MATOP_VIEW=61,
1541d519adbfSMatthew Knepley                MATOP_CONVERT_FROM=62,
1542d519adbfSMatthew Knepley                MATOP_USE_SCALED_FORM=63,
1543d519adbfSMatthew Knepley                MATOP_SCALE_SYSTEM=64,
1544d519adbfSMatthew Knepley                MATOP_UNSCALE_SYSTEM=65,
1545d519adbfSMatthew Knepley                MATOP_SET_LOCAL_TO_GLOBAL_MAP=66,
1546d519adbfSMatthew Knepley                MATOP_SET_VALUES_LOCAL=67,
1547d519adbfSMatthew Knepley                MATOP_ZERO_ROWS_LOCAL=68,
1548d519adbfSMatthew Knepley                MATOP_GET_ROW_MAX_ABS=69,
1549d519adbfSMatthew Knepley                MATOP_GET_ROW_MIN_ABS=70,
1550d519adbfSMatthew Knepley                MATOP_CONVERT=71,
1551d519adbfSMatthew Knepley                MATOP_SET_COLORING=72,
1552d519adbfSMatthew Knepley                MATOP_SET_VALUES_ADIC=73,
1553d519adbfSMatthew Knepley                MATOP_SET_VALUES_ADIFOR=74,
1554d519adbfSMatthew Knepley                MATOP_FD_COLORING_APPLY=75,
1555d519adbfSMatthew Knepley                MATOP_SET_FROM_OPTIONS=76,
1556d519adbfSMatthew Knepley                MATOP_MULT_CON=77,
1557d519adbfSMatthew Knepley                MATOP_MULT_TRANSPOSE_CON=78,
1558d519adbfSMatthew Knepley                MATOP_PERMUTE_SPARSIFY=79,
1559d519adbfSMatthew Knepley                MATOP_MULT_MULTIPLE=80,
1560d519adbfSMatthew Knepley                MATOP_SOLVE_MULTIPLE=81,
1561d519adbfSMatthew Knepley                MATOP_GET_INERTIA=82,
1562d519adbfSMatthew Knepley                MATOP_LOAD=83,
1563d519adbfSMatthew Knepley                MATOP_IS_SYMMETRIC=84,
1564d519adbfSMatthew Knepley                MATOP_IS_HERMITIAN=85,
1565d519adbfSMatthew Knepley                MATOP_IS_STRUCTURALLY_SYMMETRIC=86,
156641f059aeSBarry Smith                MATOP_DUMMY=87,
1567d519adbfSMatthew Knepley                MATOP_GET_VECS=88,
1568d519adbfSMatthew Knepley                MATOP_MAT_MULT=89,
1569d519adbfSMatthew Knepley                MATOP_MAT_MULT_SYMBOLIC=90,
1570d519adbfSMatthew Knepley                MATOP_MAT_MULT_NUMERIC=91,
1571d519adbfSMatthew Knepley                MATOP_PTAP=92,
1572d519adbfSMatthew Knepley                MATOP_PTAP_SYMBOLIC=93,
1573d519adbfSMatthew Knepley                MATOP_PTAP_NUMERIC=94,
1574d519adbfSMatthew Knepley                MATOP_MAT_MULTTRANSPOSE=95,
15751763ac4eSSatish Balay                MATOP_MAT_MULTTRANSPOSE_SYM=96,
15761763ac4eSSatish Balay                MATOP_MAT_MULTTRANSPOSE_NUM=97,
1577d519adbfSMatthew Knepley                MATOP_PTAP_SYMBOLIC_SEQAIJ=98,
1578d519adbfSMatthew Knepley                MATOP_PTAP_NUMERIC_SEQAIJ=99,
1579d519adbfSMatthew Knepley                MATOP_PTAP_SYMBOLIC_MPIAIJ=100,
1580d519adbfSMatthew Knepley                MATOP_PTAP_NUMERIC_MPIAIJ=101,
1581d519adbfSMatthew Knepley                MATOP_CONJUGATE=102,
1582d519adbfSMatthew Knepley                MATOP_SET_SIZES=103,
1583d519adbfSMatthew Knepley                MATOP_SET_VALUES_ROW=104,
1584d519adbfSMatthew Knepley                MATOP_REAL_PART=105,
1585d519adbfSMatthew Knepley                MATOP_IMAG_PART=106,
1586d519adbfSMatthew Knepley                MATOP_GET_ROW_UTRIANGULAR=107,
1587d519adbfSMatthew Knepley                MATOP_RESTORE_ROW_UTRIANGULAR=108,
1588d519adbfSMatthew Knepley                MATOP_MATSOLVE=109,
1589d519adbfSMatthew Knepley                MATOP_GET_REDUNDANTMATRIX=110,
1590d519adbfSMatthew Knepley                MATOP_GET_ROW_MIN=111,
1591d519adbfSMatthew Knepley                MATOP_GET_COLUMN_VEC=112,
1592d519adbfSMatthew Knepley                MATOP_MISSING_DIAGONAL=113,
1593d519adbfSMatthew Knepley                MATOP_MATGETSEQNONZEROSTRUCTURE=114,
159489c6957cSBarry Smith                MATOP_CREATE=115,
159589c6957cSBarry Smith                MATOP_GET_GHOSTS=116,
1596ea38565cSJed Brown                MATOP_GET_LOCALSUBMATRIX=117,
1597ea38565cSJed Brown                MATOP_RESTORE_LOCALSUBMATRIX=118,
1598eeffb40dSHong Zhang                MATOP_MULT_DIAGONAL_BLOCK=119,
1599547795f9SHong Zhang                MATOP_HERMITIANTRANSPOSE=120,
1600547795f9SHong Zhang                MATOP_MULTHERMITIANTRANSPOSE=121,
1601fbdbba38SShri Abhyankar                MATOP_MULTHERMITIANTRANSPOSEADD=122,
1602b9614d88SDmitry Karpeev                MATOP_GETMULTIPROCBLOCK=123,
16030716a85fSBarry Smith                MATOP_GETCOLUMNNORMS=125,
160437868618SMatthew G Knepley 	       MATOP_GET_SUBMATRICES_PARALLEL=128,
1605*2b8ad9a3SHong Zhang                MATOP_SET_VALUES_BATCH=129,
1606*2b8ad9a3SHong Zhang                MATOP_TRANSPOSEMATMULT=130,
1607*2b8ad9a3SHong Zhang                MATOP_TRANSPOSEMATMULT_SYMBOLIC=131,
1608*2b8ad9a3SHong Zhang                MATOP_TRANSPOSEMATMULT_NUMERIC=132,
1609*2b8ad9a3SHong Zhang                MATOP_TRANSPOSECOLORING_CREATE=133,
1610*2b8ad9a3SHong Zhang                MATOP_TRANSCOLORING_APPLY_SPTODEN=134,
1611*2b8ad9a3SHong Zhang                MATOP_TRANSCOLORING_APPLY_DENTOSP=135,
1612*2b8ad9a3SHong Zhang                MATOP_RARt=136,
1613*2b8ad9a3SHong Zhang                MATOP_RARt_SYMBOLIC=137,
1614*2b8ad9a3SHong Zhang                MATOP_RARt_NUMERIC=138
1615fae171e0SBarry Smith              } MatOperation;
16167087cfbeSBarry Smith extern PetscErrorCode  MatHasOperation(Mat,MatOperation,PetscBool *);
16177087cfbeSBarry Smith extern PetscErrorCode  MatShellSetOperation(Mat,MatOperation,void(*)(void));
16187087cfbeSBarry Smith extern PetscErrorCode  MatShellGetOperation(Mat,MatOperation,void(**)(void));
16197087cfbeSBarry Smith extern PetscErrorCode  MatShellSetContext(Mat,void*);
1620112a2221SBarry Smith 
162190ace30eSBarry Smith /*
162290ace30eSBarry Smith    Codes for matrices stored on disk. By default they are
162390ace30eSBarry Smith    stored in a universal format. By changing the format with
16247973ad04SBarry Smith    PetscViewerSetFormat(viewer,PETSC_VIEWER_NATIVE); the matrices will
162590ace30eSBarry Smith    be stored in a way natural for the matrix, for example dense matrices
162690ace30eSBarry Smith    would be stored as dense. Matrices stored this way may only be
1627f4403165SShri Abhyankar    read into matrices of the same type.
162890ace30eSBarry Smith */
162990ace30eSBarry Smith #define MATRIX_BINARY_FORMAT_DENSE -1
163090ace30eSBarry Smith 
16317087cfbeSBarry Smith extern PetscErrorCode  MatMPIBAIJSetHashTableFactor(Mat,PetscReal);
16327087cfbeSBarry Smith extern PetscErrorCode  MatISGetLocalMat(Mat,Mat*);
16331f4e1ec7SBarry Smith 
1634d9274352SBarry Smith /*S
1635d9274352SBarry Smith      MatNullSpace - Object that removes a null space from a vector, i.e.
1636d9274352SBarry Smith          orthogonalizes the vector to a subsapce
1637d9274352SBarry Smith 
1638f7a9e4ceSBarry Smith    Level: advanced
1639d9274352SBarry Smith 
1640d9274352SBarry Smith   Concepts: matrix; linear operator, null space
1641d9274352SBarry Smith 
16426e1639daSBarry Smith   Users manual sections:
16436e1639daSBarry Smith .   sec_singular
16446e1639daSBarry Smith 
1645d9274352SBarry Smith .seealso:  MatNullSpaceCreate()
1646d9274352SBarry Smith S*/
164774637425SBarry Smith typedef struct _p_MatNullSpace* MatNullSpace;
1648d9274352SBarry Smith 
16497087cfbeSBarry Smith extern PetscErrorCode  MatNullSpaceCreate(MPI_Comm,PetscBool ,PetscInt,const Vec[],MatNullSpace*);
16507087cfbeSBarry Smith extern PetscErrorCode  MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(MatNullSpace,Vec,void*),void*);
1651d34fcf5fSBarry Smith extern PetscErrorCode  MatNullSpaceDestroy(MatNullSpace*);
16527087cfbeSBarry Smith extern PetscErrorCode  MatNullSpaceRemove(MatNullSpace,Vec,Vec*);
16539c6675edSMatthew G Knepley extern PetscErrorCode  MatGetNullSpace(Mat, MatNullSpace *);
16541e633543SBarry Smith extern PetscErrorCode  MatSetNullSpace(Mat,MatNullSpace);
165510077b9fSMatthew G Knepley extern PetscErrorCode  MatSetNearNullSpace(Mat,MatNullSpace);
16567087cfbeSBarry Smith extern PetscErrorCode  MatNullSpaceTest(MatNullSpace,Mat,PetscBool  *);
1657b717e993SJed Brown extern PetscErrorCode  MatNullSpaceView(MatNullSpace,PetscViewer);
165874637425SBarry Smith 
16597087cfbeSBarry Smith extern PetscErrorCode  MatReorderingSeqSBAIJ(Mat,IS);
16607087cfbeSBarry Smith extern PetscErrorCode  MatMPISBAIJSetHashTableFactor(Mat,PetscReal);
16617087cfbeSBarry Smith extern PetscErrorCode  MatSeqSBAIJSetColumnIndices(Mat,PetscInt *);
16627087cfbeSBarry Smith extern PetscErrorCode  MatSeqBAIJInvertBlockDiagonal(Mat);
16633f1d51d7SBarry Smith 
16647087cfbeSBarry Smith extern PetscErrorCode  MatCreateMAIJ(Mat,PetscInt,Mat*);
16657087cfbeSBarry Smith extern PetscErrorCode  MatMAIJRedimension(Mat,PetscInt,Mat*);
16667087cfbeSBarry Smith extern PetscErrorCode  MatMAIJGetAIJ(Mat,Mat*);
1667c4f061fbSSatish Balay 
16687087cfbeSBarry Smith extern PetscErrorCode  MatComputeExplicitOperator(Mat,Mat*);
1669b0a32e0cSBarry Smith 
16707087cfbeSBarry Smith extern PetscErrorCode  MatDiagonalScaleLocal(Mat,Vec);
167104f1ad80SBarry Smith 
16727087cfbeSBarry Smith extern PetscErrorCode  MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*);
16737087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetBase(Mat,Vec,Vec);
16747087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*);
16757087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*));
16767087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec));
16777087cfbeSBarry Smith extern PetscErrorCode  MatMFFDAddNullSpace(Mat,MatNullSpace);
16787087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt);
16797087cfbeSBarry Smith extern PetscErrorCode  MatMFFDResetHHistory(Mat);
16807087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetFunctionError(Mat,PetscReal);
16817087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetPeriod(Mat,PetscInt);
16827087cfbeSBarry Smith extern PetscErrorCode  MatMFFDGetH(Mat,PetscScalar *);
16837087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetOptionsPrefix(Mat,const char[]);
16847087cfbeSBarry Smith extern PetscErrorCode  MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*);
16857087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*);
1686e884886fSBarry Smith 
16876370053bSBarry Smith /*S
16886370053bSBarry Smith     MatMFFD - A data structured used to manage the computation of the h differencing parameter for matrix-free
16896370053bSBarry Smith               Jacobian vector products
1690e884886fSBarry Smith 
16916370053bSBarry Smith     Notes: MATMFFD is a specific MatType which uses the MatMFFD data structure
16926370053bSBarry Smith 
16936370053bSBarry Smith            MatMFFD*() methods actually take the Mat as their first argument. Not a MatMFFD data structure
16946370053bSBarry Smith 
16956370053bSBarry Smith     Level: developer
16966370053bSBarry Smith 
16976370053bSBarry Smith .seealso: MATMFFD, MatCreateMFFD(), MatMFFDSetFuction(), MatMFFDSetType(), MatMFFDRegister()
16986370053bSBarry Smith S*/
1699e884886fSBarry Smith typedef struct _p_MatMFFD* MatMFFD;
1700e884886fSBarry Smith 
170176bdecfbSBarry Smith /*J
1702e884886fSBarry Smith     MatMFFDType - algorithm used to compute the h used in computing matrix-vector products via differencing of the function
1703e884886fSBarry Smith 
1704e884886fSBarry Smith    Level: beginner
1705e884886fSBarry Smith 
1706e884886fSBarry Smith .seealso: MatMFFDSetType(), MatMFFDRegister()
170776bdecfbSBarry Smith J*/
1708a313700dSBarry Smith #define MatMFFDType char*
1709e884886fSBarry Smith #define MATMFFD_DS  "ds"
1710e884886fSBarry Smith #define MATMFFD_WP  "wp"
1711e884886fSBarry Smith 
17127087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetType(Mat,const MatMFFDType);
17137087cfbeSBarry Smith extern PetscErrorCode  MatMFFDRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatMFFD));
1714e884886fSBarry Smith 
1715e884886fSBarry Smith /*MC
1716e884886fSBarry Smith    MatMFFDRegisterDynamic - Adds a method to the MatMFFD registry.
1717e884886fSBarry Smith 
1718e884886fSBarry Smith    Synopsis:
17191890ba74SBarry Smith    PetscErrorCode MatMFFDRegisterDynamic(const char *name_solver,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatMFFD))
1720e884886fSBarry Smith 
1721e884886fSBarry Smith    Not Collective
1722e884886fSBarry Smith 
1723e884886fSBarry Smith    Input Parameters:
1724e884886fSBarry Smith +  name_solver - name of a new user-defined compute-h module
1725e884886fSBarry Smith .  path - path (either absolute or relative) the library containing this solver
1726e884886fSBarry Smith .  name_create - name of routine to create method context
1727e884886fSBarry Smith -  routine_create - routine to create method context
1728e884886fSBarry Smith 
1729e884886fSBarry Smith    Level: developer
1730e884886fSBarry Smith 
1731e884886fSBarry Smith    Notes:
1732e884886fSBarry Smith    MatMFFDRegisterDynamic() may be called multiple times to add several user-defined solvers.
1733e884886fSBarry Smith 
1734e884886fSBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
1735e884886fSBarry Smith    is ignored.
1736e884886fSBarry Smith 
1737e884886fSBarry Smith    Sample usage:
1738e884886fSBarry Smith .vb
1739e884886fSBarry Smith    MatMFFDRegisterDynamic("my_h",/home/username/my_lib/lib/libO/solaris/mylib.a,
1740e884886fSBarry Smith                "MyHCreate",MyHCreate);
1741e884886fSBarry Smith .ve
1742e884886fSBarry Smith 
1743e884886fSBarry Smith    Then, your solver can be chosen with the procedural interface via
1744e884886fSBarry Smith $     MatMFFDSetType(mfctx,"my_h")
1745e884886fSBarry Smith    or at runtime via the option
1746e884886fSBarry Smith $     -snes_mf_type my_h
1747e884886fSBarry Smith 
1748e884886fSBarry Smith .keywords: MatMFFD, register
1749e884886fSBarry Smith 
1750e884886fSBarry Smith .seealso: MatMFFDRegisterAll(), MatMFFDRegisterDestroy()
1751e884886fSBarry Smith M*/
1752e884886fSBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1753e884886fSBarry Smith #define MatMFFDRegisterDynamic(a,b,c,d) MatMFFDRegister(a,b,c,0)
1754e884886fSBarry Smith #else
1755e884886fSBarry Smith #define MatMFFDRegisterDynamic(a,b,c,d) MatMFFDRegister(a,b,c,d)
1756e884886fSBarry Smith #endif
1757e884886fSBarry Smith 
17587087cfbeSBarry Smith extern PetscErrorCode  MatMFFDRegisterAll(const char[]);
17597087cfbeSBarry Smith extern PetscErrorCode  MatMFFDRegisterDestroy(void);
17607087cfbeSBarry Smith extern PetscErrorCode  MatMFFDDSSetUmin(Mat,PetscReal);
17617087cfbeSBarry Smith extern PetscErrorCode  MatMFFDWPSetComputeNormU(Mat,PetscBool );
1762e884886fSBarry Smith 
1763e884886fSBarry Smith 
17647087cfbeSBarry Smith extern PetscErrorCode  PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *);
17657087cfbeSBarry Smith extern PetscErrorCode  PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *);
17667dbadf16SMatthew Knepley 
176797969023SHong Zhang /*
176897969023SHong Zhang    PETSc interface to MUMPS
176997969023SHong Zhang */
177097969023SHong Zhang #ifdef PETSC_HAVE_MUMPS
1771f250808bSBarry Smith extern PetscErrorCode  MatMumpsSetIcntl(Mat,PetscInt,PetscInt);
177297969023SHong Zhang #endif
177323a5497aSJed Brown 
1774d954db57SHong Zhang /*
1775d954db57SHong Zhang    PETSc interface to SUPERLU
1776d954db57SHong Zhang */
1777d954db57SHong Zhang #ifdef PETSC_HAVE_SUPERLU
1778f250808bSBarry Smith extern PetscErrorCode  MatSuperluSetILUDropTol(Mat,PetscReal);
1779d954db57SHong Zhang #endif
1780d954db57SHong Zhang 
178190c53902SBarry Smith #if defined(PETSC_HAVE_CUSP)
178290c53902SBarry Smith extern PetscErrorCode  MatCreateSeqAIJCUSP(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
178390c53902SBarry Smith extern PetscErrorCode  MatCreateMPIAIJCUSP(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
178490c53902SBarry Smith #endif
178590c53902SBarry Smith 
178654efbe56SHong Zhang /*
178754efbe56SHong Zhang    PETSc interface to FFTW
178854efbe56SHong Zhang */
178954efbe56SHong Zhang #if defined(PETSC_HAVE_FFTW)
179074a26884SAmlan Barua extern PetscErrorCode VecScatterPetscToFFTW(Mat,Vec,Vec);
179174a26884SAmlan Barua extern PetscErrorCode VecScatterFFTWToPetsc(Mat,Vec,Vec);
17924be45526SHong Zhang extern PetscErrorCode MatGetVecsFFTW(Mat,Vec*,Vec*,Vec*);
179354efbe56SHong Zhang #endif
179454efbe56SHong Zhang 
1795c8883902SJed Brown extern PetscErrorCode MatCreateNest(MPI_Comm,PetscInt,const IS[],PetscInt,const IS[],const Mat[],Mat*);
1796c8883902SJed Brown extern PetscErrorCode MatNestGetSize(Mat,PetscInt*,PetscInt*);
1797c8883902SJed Brown extern PetscErrorCode MatNestGetSubMats(Mat,PetscInt*,PetscInt*,Mat***);
1798c8883902SJed Brown extern PetscErrorCode MatNestGetSubMat(Mat,PetscInt,PetscInt,Mat*);
17997087cfbeSBarry Smith extern PetscErrorCode MatNestSetVecType(Mat,const VecType);
1800c8883902SJed Brown extern PetscErrorCode MatNestSetSubMats(Mat,PetscInt,const IS[],PetscInt,const IS[],const Mat[]);
18010782ca92SJed Brown extern PetscErrorCode MatNestSetSubMat(Mat,PetscInt,PetscInt,Mat);
1802d8588912SDave May 
180315c1f1baSDmitry Karpeev /*
180415c1f1baSDmitry Karpeev  MatIJ:
180515c1f1baSDmitry Karpeev  An unweighted directed pseudograph
180615c1f1baSDmitry Karpeev  An interpretation of this matrix as a (pseudo)graph allows us to define additional operations on it:
180715c1f1baSDmitry Karpeev  A MatIJ can act on sparse arrays: arrays of indices, or index arrays of integers, scalars, or integer-scalar pairs
180815c1f1baSDmitry Karpeev  by mapping the indices to the indices connected to them by the (pseudo)graph ed
180915c1f1baSDmitry Karpeev  */
181015c1f1baSDmitry Karpeev typedef enum {MATIJ_LOCAL, MATIJ_GLOBAL} MatIJIndexType;
181115c1f1baSDmitry Karpeev extern  PetscErrorCode MatIJSetMultivalued(Mat, PetscBool);
181215c1f1baSDmitry Karpeev extern  PetscErrorCode MatIJGetMultivalued(Mat, PetscBool*);
181315c1f1baSDmitry Karpeev extern  PetscErrorCode MatIJSetEdges(Mat, PetscInt, const PetscInt*, const PetscInt*);
181415c1f1baSDmitry Karpeev extern  PetscErrorCode MatIJGetEdges(Mat, PetscInt *, PetscInt **, PetscInt **);
181515c1f1baSDmitry Karpeev extern  PetscErrorCode MatIJSetEdgesIS(Mat, IS, IS);
181615c1f1baSDmitry Karpeev extern  PetscErrorCode MatIJGetEdgesIS(Mat, IS*, IS*);
181715c1f1baSDmitry Karpeev extern  PetscErrorCode MatIJGetRowSizes(Mat, MatIJIndexType, PetscInt, const PetscInt *, PetscInt **);
181815c1f1baSDmitry Karpeev extern  PetscErrorCode MatIJGetMinRowSize(Mat, PetscInt *);
181915c1f1baSDmitry Karpeev extern  PetscErrorCode MatIJGetMaxRowSize(Mat, PetscInt *);
182015c1f1baSDmitry Karpeev extern  PetscErrorCode MatIJGetSupport(Mat,  PetscInt *, PetscInt **);
182115c1f1baSDmitry Karpeev extern  PetscErrorCode MatIJGetSupportIS(Mat, IS *);
182215c1f1baSDmitry Karpeev extern  PetscErrorCode MatIJGetImage(Mat, PetscInt*, PetscInt**);
182315c1f1baSDmitry Karpeev extern  PetscErrorCode MatIJGetImageIS(Mat, IS *);
182415c1f1baSDmitry Karpeev extern  PetscErrorCode MatIJGetSupportSize(Mat, PetscInt *);
182515c1f1baSDmitry Karpeev extern  PetscErrorCode MatIJGetImageSize(Mat, PetscInt *);
182615c1f1baSDmitry Karpeev 
182715c1f1baSDmitry Karpeev extern  PetscErrorCode MatIJBinRenumber(Mat, Mat*);
182815c1f1baSDmitry Karpeev 
18299b67c4b3SDmitry Karpeev extern  PetscErrorCode MatIJMap(Mat, MatIJIndexType, PetscInt,const PetscInt*,const PetscInt*,const PetscScalar*, MatIJIndexType,PetscInt*,PetscInt**,PetscInt**,PetscScalar**,PetscInt**);
18309b67c4b3SDmitry Karpeev extern  PetscErrorCode MatIJBin(Mat, MatIJIndexType, PetscInt,const PetscInt*,const PetscInt*,const PetscScalar*,PetscInt*,PetscInt**,PetscInt**,PetscScalar**,PetscInt**);
18319b67c4b3SDmitry Karpeev extern  PetscErrorCode MatIJBinMap(Mat,Mat, MatIJIndexType,PetscInt,const PetscInt*,const PetscInt*,const PetscScalar*,MatIJIndexType,PetscInt*,PetscInt**,PetscInt**,PetscScalar**,PetscInt**);
183215c1f1baSDmitry Karpeev 
183323a5497aSJed Brown PETSC_EXTERN_CXX_END
183423a5497aSJed Brown #endif
1835