xref: /petsc/include/petscmat.h (revision 1e633543be291cfcfe1cde705e1b44f9f92cb5cc)
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 
20d9274352SBarry Smith /*E
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
28d91e6319SBarry Smith E*/
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 
61c0cdd4a1SDahai Guo #define MATSEQBSTRM        "seqbstrm"
62c0cdd4a1SDahai Guo #define MATMPIBSTRM        "mpibstrm"
63c0cdd4a1SDahai Guo #define MATBSTRM           "bstrm"
64aa5a9175SDahai Guo #define MATSEQSBSTRM       "seqsbstrm"
65aa5a9175SDahai Guo #define MATMPISBSTRM       "mpisbstrm"
66aa5a9175SDahai Guo #define MATSBSTRM          "sbstrm"
67c0cdd4a1SDahai Guo 
68cebc7f6cSBarry Smith #define MATDAAD            "daad"
69cebc7f6cSBarry Smith #define MATMFFD            "mffd"
70c8a8475eSBarry Smith #define MATNORMAL          "normal"
71ab92ecdeSBarry Smith #define MATLRC             "lrc"
722a6744ebSBarry Smith #define MATSCATTER         "scatter"
73421e10b8SBarry Smith #define MATBLOCKMAT        "blockmat"
74793850ffSBarry Smith #define MATCOMPOSITE       "composite"
751f8c7532SHong Zhang #define MATFFT             "fft"
761f8c7532SHong Zhang #define MATFFTW            "fftw"
77e133240eSMatthew G Knepley #define MATSEQCUFFT        "seqcufft"
78557cca28SSatish Balay #define MATTRANSPOSEMAT    "transpose"
7972ca8751SBarry Smith #define MATSCHURCOMPLEMENT "schurcomplement"
801d6018f0SLisandro Dalcin #define MATPYTHON          "python"
81f91d8e95SBarry Smith #define MATHYPRESTRUCT     "hyprestruct"
82a9e6138eSGlenn Hammond #define MATHYPRESSTRUCT    "hypresstruct"
83ee1cef2cSJed Brown #define MATSUBMATRIX       "submatrix"
842c0dbf93SJed Brown #define MATLOCALREF        "localref"
85d8588912SDave May #define MATNEST            "nest"
86773941d6SBarry Smith 
879989ab13SBarry Smith /*E
88c7393fdbSBarry Smith     MatSolverPackage - String with the name of a PETSc matrix solver type.
899989ab13SBarry Smith 
909989ab13SBarry Smith     For example: "petsc" indicates what PETSc provides, "superlu" indicates either
919989ab13SBarry Smith        SuperLU or SuperLU_Dist etc.
929989ab13SBarry Smith 
939989ab13SBarry Smith 
949989ab13SBarry Smith    Level: beginner
959989ab13SBarry Smith 
965c9eb25fSBarry Smith .seealso: MatGetFactor(), Mat, MatSetType(), MatType
979989ab13SBarry Smith E*/
98c7393fdbSBarry Smith #define MatSolverPackage char*
992692d6eeSBarry Smith #define MATSOLVERSPOOLES      "spooles"
1002692d6eeSBarry Smith #define MATSOLVERSUPERLU      "superlu"
1012692d6eeSBarry Smith #define MATSOLVERSUPERLU_DIST "superlu_dist"
1022692d6eeSBarry Smith #define MATSOLVERUMFPACK      "umfpack"
10320db9a53SJed Brown #define MATSOLVERCHOLMOD      "cholmod"
1042692d6eeSBarry Smith #define MATSOLVERESSL         "essl"
1052692d6eeSBarry Smith #define MATSOLVERLUSOL        "lusol"
1062692d6eeSBarry Smith #define MATSOLVERMUMPS        "mumps"
1072692d6eeSBarry Smith #define MATSOLVERPASTIX       "pastix"
1082692d6eeSBarry Smith #define MATSOLVERMATLAB       "matlab"
1092692d6eeSBarry Smith #define MATSOLVERPETSC        "petsc"
1102692d6eeSBarry Smith #define MATSOLVERPLAPACK      "plapack"
1112692d6eeSBarry Smith #define MATSOLVERBAS          "bas"
112773941d6SBarry Smith 
113aa5a9175SDahai Guo #define MATSOLVERBSTRM        "bstrm"
114aa5a9175SDahai Guo #define MATSOLVERSBSTRM       "sbstrm"
115c0cdd4a1SDahai Guo 
116b24902e0SBarry Smith /*E
117b24902e0SBarry Smith     MatFactorType - indicates what type of factorization is requested
118b24902e0SBarry Smith 
119b24902e0SBarry Smith     Level: beginner
120b24902e0SBarry Smith 
121b24902e0SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
122b24902e0SBarry Smith 
123c7393fdbSBarry Smith .seealso: MatSolverPackage, MatGetFactor()
124b24902e0SBarry Smith E*/
125599ef60dSHong Zhang typedef enum {MAT_FACTOR_NONE, MAT_FACTOR_LU, MAT_FACTOR_CHOLESKY, MAT_FACTOR_ILU, MAT_FACTOR_ICC,MAT_FACTOR_ILUDT} MatFactorType;
12634918c53SJed Brown extern const char *const MatFactorTypes[];
127e92e720dSBarry Smith 
1287087cfbeSBarry Smith extern PetscErrorCode  MatGetFactor(Mat,const MatSolverPackage,MatFactorType,Mat*);
1297087cfbeSBarry Smith extern PetscErrorCode  MatGetFactorAvailable(Mat,const MatSolverPackage,MatFactorType,PetscBool *);
1307087cfbeSBarry Smith extern PetscErrorCode  MatFactorGetSolverPackage(Mat,const MatSolverPackage*);
1317087cfbeSBarry Smith extern PetscErrorCode  MatGetFactorType(Mat,MatFactorType*);
1329989ab13SBarry Smith 
133c06d978dSMatthew Knepley /* Logging support */
1340700a824SBarry Smith #define    MAT_FILE_CLASSID 1211216    /* used to indicate matrices in binary files */
1357087cfbeSBarry Smith extern PetscClassId  MAT_CLASSID;
1367087cfbeSBarry Smith extern PetscClassId  MAT_FDCOLORING_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
3792ef4de8bSBarry Smith         column of a matrix as index on an associated grid.
3802ef4de8bSBarry Smith 
3812ef4de8bSBarry Smith    Level: beginner
3822ef4de8bSBarry Smith 
3832ef4de8bSBarry Smith   Concepts: matrix; linear operator
3842ef4de8bSBarry Smith 
385d7bd93baSBarry Smith .seealso:  MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockStencil()
3862ef4de8bSBarry Smith S*/
387435da068SBarry Smith typedef struct {
388c1ac3661SBarry Smith   PetscInt k,j,i,c;
389435da068SBarry Smith } MatStencil;
3902ef4de8bSBarry Smith 
3917087cfbeSBarry Smith extern PetscErrorCode  MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
3927087cfbeSBarry Smith extern PetscErrorCode  MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
3937087cfbeSBarry Smith extern PetscErrorCode  MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt);
394435da068SBarry Smith 
3957087cfbeSBarry Smith extern PetscErrorCode  MatSetColoring(Mat,ISColoring);
3967087cfbeSBarry Smith extern PetscErrorCode  MatSetValuesAdic(Mat,void*);
3977087cfbeSBarry Smith extern PetscErrorCode  MatSetValuesAdifor(Mat,PetscInt,void*);
3983a7fca6bSBarry Smith 
399d91e6319SBarry Smith /*E
400d91e6319SBarry Smith     MatAssemblyType - Indicates if the matrix is now to be used, or if you plan
401d91e6319SBarry Smith      to continue to add values to it
402d91e6319SBarry Smith 
403d91e6319SBarry Smith     Level: beginner
404d91e6319SBarry Smith 
405d91e6319SBarry Smith .seealso: MatAssemblyBegin(), MatAssemblyEnd()
406d91e6319SBarry Smith E*/
4076d4a8577SBarry Smith typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType;
4087087cfbeSBarry Smith extern PetscErrorCode  MatAssemblyBegin(Mat,MatAssemblyType);
4097087cfbeSBarry Smith extern PetscErrorCode  MatAssemblyEnd(Mat,MatAssemblyType);
4107087cfbeSBarry Smith extern PetscErrorCode  MatAssembled(Mat,PetscBool *);
4114f9c727eSBarry Smith 
4121d73ed98SBarry Smith 
41330de9b25SBarry Smith 
414d91e6319SBarry Smith /*E
415d91e6319SBarry Smith     MatOption - Options that may be set for a matrix and its behavior or storage
416d91e6319SBarry Smith 
417d91e6319SBarry Smith     Level: beginner
418d91e6319SBarry Smith 
4190a835dfdSSatish Balay    Any additions/changes here MUST also be made in include/finclude/petscmat.h
420d91e6319SBarry Smith 
421d91e6319SBarry Smith .seealso: MatSetOption()
422d91e6319SBarry Smith E*/
423512a5fc5SBarry Smith typedef enum {MAT_ROW_ORIENTED,MAT_NEW_NONZERO_LOCATIONS,
4244e0d8c25SBarry Smith               MAT_SYMMETRIC,
4254e0d8c25SBarry Smith               MAT_STRUCTURALLY_SYMMETRIC,
4268047e77bSBarry Smith               MAT_NEW_DIAGONALS,MAT_IGNORE_OFF_PROC_ENTRIES,
4274e0d8c25SBarry Smith               MAT_NEW_NONZERO_LOCATION_ERR,
4284e0d8c25SBarry Smith               MAT_NEW_NONZERO_ALLOCATION_ERR,MAT_USE_HASH_TABLE,
429a9817697SBarry Smith               MAT_KEEP_NONZERO_PATTERN,MAT_IGNORE_ZERO_ENTRIES,
4304e0d8c25SBarry Smith               MAT_USE_INODES,
4314e0d8c25SBarry Smith               MAT_HERMITIAN,
4324e0d8c25SBarry Smith               MAT_SYMMETRY_ETERNAL,
433cd6b891eSBarry Smith               MAT_CHECK_COMPRESSED_ROW,
4344e0d8c25SBarry Smith               MAT_IGNORE_LOWER_TRIANGULAR,MAT_ERROR_LOWER_TRIANGULAR,
43528b2fa4aSMatthew Knepley               MAT_GETROW_UPPERTRIANGULAR,MAT_UNUSED_NONZERO_LOCATION_ERR,
4364cb17eb5SBarry Smith               MAT_SPD,MAT_NO_OFF_PROC_ENTRIES,MAT_NO_OFF_PROC_ZERO_ROWS,
43728b2fa4aSMatthew Knepley               NUM_MAT_OPTIONS} MatOption;
438290bbb0aSBarry Smith extern const char *MatOptions[];
4397087cfbeSBarry Smith extern PetscErrorCode  MatSetOption(Mat,MatOption,PetscBool );
4407087cfbeSBarry Smith extern PetscErrorCode  MatGetType(Mat,const MatType*);
441a313700dSBarry Smith PetscPolymorphicFunction(MatGetType,(Mat mat),(mat,&t),const MatType,t)
44284cb2905SBarry Smith 
4437087cfbeSBarry Smith extern PetscErrorCode  MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]);
4447087cfbeSBarry Smith extern PetscErrorCode  MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
4457087cfbeSBarry Smith extern PetscErrorCode  MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
4467087cfbeSBarry Smith extern PetscErrorCode  MatGetRowUpperTriangular(Mat);
4477087cfbeSBarry Smith extern PetscErrorCode  MatRestoreRowUpperTriangular(Mat);
4487087cfbeSBarry Smith extern PetscErrorCode  MatGetColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
4497087cfbeSBarry Smith extern PetscErrorCode  MatRestoreColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
4507087cfbeSBarry Smith extern PetscErrorCode  MatGetColumnVector(Mat,Vec,PetscInt);
4517087cfbeSBarry Smith extern PetscErrorCode  MatGetArray(Mat,PetscScalar *[]);
452ba966639SSatish Balay PetscPolymorphicFunction(MatGetArray,(Mat mat),(mat,&a),PetscScalar*,a)
4537087cfbeSBarry Smith extern PetscErrorCode  MatRestoreArray(Mat,PetscScalar *[]);
4547087cfbeSBarry Smith extern PetscErrorCode  MatGetBlockSize(Mat,PetscInt *);
455ba966639SSatish Balay PetscPolymorphicFunction(MatGetBlockSize,(Mat mat),(mat,&a),PetscInt,a)
4567087cfbeSBarry Smith extern PetscErrorCode  MatSetBlockSize(Mat,PetscInt);
4577b80b807SBarry Smith 
4581620fd73SBarry Smith 
4597087cfbeSBarry Smith extern PetscErrorCode  MatMult(Mat,Vec,Vec);
4607087cfbeSBarry Smith extern PetscErrorCode  MatMultDiagonalBlock(Mat,Vec,Vec);
4617087cfbeSBarry Smith extern PetscErrorCode  MatMultAdd(Mat,Vec,Vec,Vec);
462ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
4637087cfbeSBarry Smith extern PetscErrorCode  MatMultTranspose(Mat,Vec,Vec);
4647087cfbeSBarry Smith extern PetscErrorCode  MatMultHermitianTranspose(Mat,Vec,Vec);
4657087cfbeSBarry Smith extern PetscErrorCode  MatIsTranspose(Mat,Mat,PetscReal,PetscBool *);
466ace3abfcSBarry Smith PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B,PetscReal tol),(A,B,tol,&t),PetscBool ,t)
467ace3abfcSBarry Smith PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B),(A,B,0,&t),PetscBool ,t)
4687087cfbeSBarry Smith extern PetscErrorCode  MatIsHermitianTranspose(Mat,Mat,PetscReal,PetscBool *);
4697087cfbeSBarry Smith extern PetscErrorCode  MatMultTransposeAdd(Mat,Vec,Vec,Vec);
4707087cfbeSBarry Smith extern PetscErrorCode  MatMultHermitianTransposeAdd(Mat,Vec,Vec,Vec);
471ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultTransposeAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
4727087cfbeSBarry Smith extern PetscErrorCode  MatMultConstrained(Mat,Vec,Vec);
4737087cfbeSBarry Smith extern PetscErrorCode  MatMultTransposeConstrained(Mat,Vec,Vec);
4747087cfbeSBarry Smith extern PetscErrorCode  MatMatSolve(Mat,Mat,Mat);
475f5edf698SHong Zhang 
476d91e6319SBarry Smith /*E
477d91e6319SBarry Smith     MatDuplicateOption - Indicates if a duplicated sparse matrix should have
478d91e6319SBarry Smith   its numerical values copied over or just its nonzero structure.
479d91e6319SBarry Smith 
480d91e6319SBarry Smith     Level: beginner
481d91e6319SBarry Smith 
482d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
483d91e6319SBarry Smith 
48470dcbbb9SBarry Smith $   MAT_SHARE_NONZERO_PATTERN - the i and j arrays in the new matrix will be shared with the original matrix
48570dcbbb9SBarry Smith $                               this also triggers the MAT_DO_NOT_COPY_VALUES option. This is used when you
48670dcbbb9SBarry Smith $                               have several matrices with the same nonzero pattern.
48770dcbbb9SBarry Smith 
488d91e6319SBarry Smith .seealso: MatDuplicate()
489d91e6319SBarry Smith E*/
49070dcbbb9SBarry Smith typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES,MAT_SHARE_NONZERO_PATTERN} MatDuplicateOption;
4912e8a6d31SBarry Smith 
4927087cfbeSBarry Smith extern PetscErrorCode  MatConvert(Mat,const MatType,MatReuse,Mat*);
493a313700dSBarry Smith PetscPolymorphicFunction(MatConvert,(Mat A,const MatType t),(A,t,MAT_INITIAL_MATRIX,&a),Mat,a)
4947087cfbeSBarry Smith extern PetscErrorCode  MatDuplicate(Mat,MatDuplicateOption,Mat*);
495ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A,MatDuplicateOption o),(A,o,&a),Mat,a)
496ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A),(A,MAT_COPY_VALUES,&a),Mat,a)
49794a9d846SBarry Smith 
498cb5b572fSBarry Smith 
4997087cfbeSBarry Smith extern PetscErrorCode  MatCopy(Mat,Mat,MatStructure);
5007087cfbeSBarry Smith extern PetscErrorCode  MatView(Mat,PetscViewer);
5017087cfbeSBarry Smith extern PetscErrorCode  MatIsSymmetric(Mat,PetscReal,PetscBool *);
502ace3abfcSBarry Smith PetscPolymorphicFunction(MatIsSymmetric,(Mat A,PetscReal tol),(A,tol,&t),PetscBool ,t)
503ace3abfcSBarry Smith PetscPolymorphicFunction(MatIsSymmetric,(Mat A),(A,0,&t),PetscBool ,t)
5047087cfbeSBarry Smith extern PetscErrorCode  MatIsStructurallySymmetric(Mat,PetscBool *);
505ace3abfcSBarry Smith PetscPolymorphicFunction(MatIsStructurallySymmetric,(Mat A),(A,&t),PetscBool ,t)
5067087cfbeSBarry Smith extern PetscErrorCode  MatIsHermitian(Mat,PetscReal,PetscBool *);
507ace3abfcSBarry Smith PetscPolymorphicFunction(MatIsHermitian,(Mat A),(A,0,&t),PetscBool ,t)
5087087cfbeSBarry Smith extern PetscErrorCode  MatIsSymmetricKnown(Mat,PetscBool *,PetscBool *);
5097087cfbeSBarry Smith extern PetscErrorCode  MatIsHermitianKnown(Mat,PetscBool *,PetscBool *);
5107087cfbeSBarry Smith extern PetscErrorCode  MatMissingDiagonal(Mat,PetscBool  *,PetscInt *);
5117087cfbeSBarry Smith extern PetscErrorCode  MatLoad(Mat, PetscViewer);
5127b80b807SBarry Smith 
5137087cfbeSBarry Smith extern PetscErrorCode  MatGetRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,PetscInt *[],PetscInt *[],PetscBool  *);
5147087cfbeSBarry Smith extern PetscErrorCode  MatRestoreRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,PetscInt *[],PetscInt *[],PetscBool  *);
5157087cfbeSBarry Smith extern PetscErrorCode  MatGetColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,PetscInt *[],PetscInt *[],PetscBool  *);
5167087cfbeSBarry Smith extern PetscErrorCode  MatRestoreColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,PetscInt *[],PetscInt *[],PetscBool  *);
517d4fbbf0eSBarry Smith 
518d91e6319SBarry Smith /*S
519d91e6319SBarry Smith      MatInfo - Context of matrix information, used with MatGetInfo()
520d91e6319SBarry Smith 
521d91e6319SBarry Smith    In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE
522d91e6319SBarry Smith 
523d91e6319SBarry Smith    Level: intermediate
524d91e6319SBarry Smith 
525d91e6319SBarry Smith   Concepts: matrix^nonzero information
526d91e6319SBarry Smith 
527d9274352SBarry Smith .seealso:  MatGetInfo(), MatInfoType
528d91e6319SBarry Smith S*/
5294e220ebcSLois Curfman McInnes typedef struct {
530b0a32e0cSBarry Smith   PetscLogDouble block_size;                         /* block size */
531b0a32e0cSBarry Smith   PetscLogDouble nz_allocated,nz_used,nz_unneeded;   /* number of nonzeros */
532b0a32e0cSBarry Smith   PetscLogDouble memory;                             /* memory allocated */
533b0a32e0cSBarry Smith   PetscLogDouble assemblies;                         /* number of matrix assemblies called */
534b0a32e0cSBarry Smith   PetscLogDouble mallocs;                            /* number of mallocs during MatSetValues() */
535b0a32e0cSBarry Smith   PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */
536b0a32e0cSBarry Smith   PetscLogDouble factor_mallocs;                     /* number of mallocs during factorization */
5374e220ebcSLois Curfman McInnes } MatInfo;
5384e220ebcSLois Curfman McInnes 
539d9274352SBarry Smith /*E
540d9274352SBarry Smith     MatInfoType - Indicates if you want information about the local part of the matrix,
541d9274352SBarry Smith      the entire parallel matrix or the maximum over all the local parts.
542d9274352SBarry Smith 
543d9274352SBarry Smith     Level: beginner
544d9274352SBarry Smith 
545d9274352SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
546d9274352SBarry Smith 
547d9274352SBarry Smith .seealso: MatGetInfo(), MatInfo
548d9274352SBarry Smith E*/
5497b80b807SBarry Smith typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType;
5507087cfbeSBarry Smith extern PetscErrorCode  MatGetInfo(Mat,MatInfoType,MatInfo*);
5517087cfbeSBarry Smith extern PetscErrorCode  MatGetDiagonal(Mat,Vec);
5527087cfbeSBarry Smith extern PetscErrorCode  MatGetRowMax(Mat,Vec,PetscInt[]);
5537087cfbeSBarry Smith extern PetscErrorCode  MatGetRowMin(Mat,Vec,PetscInt[]);
5547087cfbeSBarry Smith extern PetscErrorCode  MatGetRowMaxAbs(Mat,Vec,PetscInt[]);
5557087cfbeSBarry Smith extern PetscErrorCode  MatGetRowMinAbs(Mat,Vec,PetscInt[]);
5567087cfbeSBarry Smith extern PetscErrorCode  MatGetRowSum(Mat,Vec);
5577087cfbeSBarry Smith extern PetscErrorCode  MatTranspose(Mat,MatReuse,Mat*);
558fc4dec0aSBarry Smith PetscPolymorphicFunction(MatTranspose,(Mat A),(A,MAT_INITIAL_MATRIX,&t),Mat,t)
5597087cfbeSBarry Smith extern PetscErrorCode  MatHermitianTranspose(Mat,MatReuse,Mat*);
5607087cfbeSBarry Smith extern PetscErrorCode  MatPermute(Mat,IS,IS,Mat *);
561ba966639SSatish Balay PetscPolymorphicFunction(MatPermute,(Mat A,IS is1,IS is2),(A,is1,is2,&t),Mat,t)
5627087cfbeSBarry Smith extern PetscErrorCode  MatDiagonalScale(Mat,Vec,Vec);
5637087cfbeSBarry Smith extern PetscErrorCode  MatDiagonalSet(Mat,Vec,InsertMode);
5647087cfbeSBarry Smith extern PetscErrorCode  MatEqual(Mat,Mat,PetscBool *);
565ace3abfcSBarry Smith PetscPolymorphicFunction(MatEqual,(Mat A,Mat B),(A,B,&t),PetscBool ,t)
5667087cfbeSBarry Smith extern PetscErrorCode  MatMultEqual(Mat,Mat,PetscInt,PetscBool *);
5677087cfbeSBarry Smith extern PetscErrorCode  MatMultAddEqual(Mat,Mat,PetscInt,PetscBool *);
5687087cfbeSBarry Smith extern PetscErrorCode  MatMultTransposeEqual(Mat,Mat,PetscInt,PetscBool *);
5697087cfbeSBarry Smith extern PetscErrorCode  MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscBool *);
5707b80b807SBarry Smith 
5717087cfbeSBarry Smith extern PetscErrorCode  MatNorm(Mat,NormType,PetscReal *);
572ba966639SSatish Balay PetscPolymorphicFunction(MatNorm,(Mat A,NormType t),(A,t,&n),PetscReal,n)
57309573ac7SBarry Smith extern PetscErrorCode  MatGetColumnNorms(Mat,NormType,PetscReal *);
5747087cfbeSBarry Smith extern PetscErrorCode  MatZeroEntries(Mat);
5757087cfbeSBarry Smith extern PetscErrorCode  MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
5767087cfbeSBarry Smith extern PetscErrorCode  MatZeroRowsIS(Mat,IS,PetscScalar,Vec,Vec);
5777087cfbeSBarry Smith extern PetscErrorCode  MatZeroRowsStencil(Mat,PetscInt,const MatStencil [],PetscScalar,Vec,Vec);
5787087cfbeSBarry Smith extern PetscErrorCode  MatZeroRowsColumns(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
5797087cfbeSBarry Smith extern PetscErrorCode  MatZeroRowsColumnsIS(Mat,IS,PetscScalar,Vec,Vec);
5807b80b807SBarry Smith 
5817087cfbeSBarry Smith extern PetscErrorCode  MatUseScaledForm(Mat,PetscBool );
5827087cfbeSBarry Smith extern PetscErrorCode  MatScaleSystem(Mat,Vec,Vec);
5837087cfbeSBarry Smith extern PetscErrorCode  MatUnScaleSystem(Mat,Vec,Vec);
5845ef9f2a5SBarry Smith 
5857087cfbeSBarry Smith extern PetscErrorCode  MatGetSize(Mat,PetscInt*,PetscInt*);
5867087cfbeSBarry Smith extern PetscErrorCode  MatGetLocalSize(Mat,PetscInt*,PetscInt*);
5877087cfbeSBarry Smith extern PetscErrorCode  MatGetOwnershipRange(Mat,PetscInt*,PetscInt*);
5887087cfbeSBarry Smith extern PetscErrorCode  MatGetOwnershipRanges(Mat,const PetscInt**);
5897087cfbeSBarry Smith extern PetscErrorCode  MatGetOwnershipRangeColumn(Mat,PetscInt*,PetscInt*);
5907087cfbeSBarry Smith extern PetscErrorCode  MatGetOwnershipRangesColumn(Mat,const PetscInt**);
5917b80b807SBarry Smith 
5927087cfbeSBarry Smith extern PetscErrorCode  MatGetSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
5937087cfbeSBarry Smith extern PetscErrorCode  MatGetSubMatricesParallel(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
5947087cfbeSBarry Smith extern PetscErrorCode  MatDestroyMatrices(PetscInt,Mat *[]);
5957087cfbeSBarry Smith extern PetscErrorCode  MatGetSubMatrix(Mat,IS,IS,MatReuse,Mat *);
5967087cfbeSBarry Smith extern PetscErrorCode  MatGetLocalSubMatrix(Mat,IS,IS,Mat*);
5977087cfbeSBarry Smith extern PetscErrorCode  MatRestoreLocalSubMatrix(Mat,IS,IS,Mat*);
5987087cfbeSBarry Smith extern PetscErrorCode  MatGetSeqNonzeroStructure(Mat,Mat*);
5997087cfbeSBarry Smith extern PetscErrorCode  MatDestroySeqNonzeroStructure(Mat*);
6005494a064SHong Zhang 
6017087cfbeSBarry Smith extern PetscErrorCode  MatMerge(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
6027087cfbeSBarry Smith extern PetscErrorCode  MatMerge_SeqsToMPI(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*);
6037087cfbeSBarry Smith extern PetscErrorCode  MatMerge_SeqsToMPISymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*);
6047087cfbeSBarry Smith extern PetscErrorCode  MatMerge_SeqsToMPINumeric(Mat,Mat);
6054a2b5492SBarry Smith extern PetscErrorCode  MatMPIAIJGetLocalMat(Mat,MatReuse,Mat*);
6064a2b5492SBarry Smith extern PetscErrorCode  MatMPIAIJGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*);
6077087cfbeSBarry Smith extern PetscErrorCode  MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,PetscInt*,Mat*);
6087087cfbeSBarry Smith extern PetscErrorCode  MatGetBrowsOfAoCols(Mat,Mat,MatReuse,PetscInt**,PetscInt**,MatScalar**,Mat*);
60943eb5e2fSMatthew Knepley #if defined (PETSC_USE_CTABLE)
610cd116e26SSatish Balay #include "petscctable.h"
6117087cfbeSBarry Smith extern PetscErrorCode  MatGetCommunicationStructs(Mat, Vec *, PetscTable *, VecScatter *);
61243eb5e2fSMatthew Knepley #else
6137087cfbeSBarry Smith extern PetscErrorCode  MatGetCommunicationStructs(Mat, Vec *, PetscInt *[], VecScatter *);
61443eb5e2fSMatthew Knepley #endif
6157087cfbeSBarry Smith extern PetscErrorCode  MatGetGhosts(Mat, PetscInt *,const PetscInt *[]);
6168efafbd8SBarry Smith 
6177087cfbeSBarry Smith extern PetscErrorCode  MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt);
6187b80b807SBarry Smith 
6197087cfbeSBarry Smith extern PetscErrorCode  MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*);
6207087cfbeSBarry Smith extern PetscErrorCode  MatMatMultSymbolic(Mat,Mat,PetscReal,Mat*);
6217087cfbeSBarry Smith extern PetscErrorCode  MatMatMultNumeric(Mat,Mat,Mat);
62222440eb1SKris Buschelman 
6237087cfbeSBarry Smith extern PetscErrorCode  MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*);
6247087cfbeSBarry Smith extern PetscErrorCode  MatPtAPSymbolic(Mat,Mat,PetscReal,Mat*);
6257087cfbeSBarry Smith extern PetscErrorCode  MatPtAPNumeric(Mat,Mat,Mat);
62622440eb1SKris Buschelman 
6277087cfbeSBarry Smith extern PetscErrorCode  MatMatMultTranspose(Mat,Mat,MatReuse,PetscReal,Mat*);
6287087cfbeSBarry Smith extern PetscErrorCode  MatMatMultTransposeSymbolic(Mat,Mat,PetscReal,Mat*);
6297087cfbeSBarry Smith extern PetscErrorCode  MatMatMultTransposeNumeric(Mat,Mat,Mat);
630bc011b1eSHong Zhang 
6317087cfbeSBarry Smith extern PetscErrorCode  MatAXPY(Mat,PetscScalar,Mat,MatStructure);
6327087cfbeSBarry Smith extern PetscErrorCode  MatAYPX(Mat,PetscScalar,Mat,MatStructure);
6331c741599SHong Zhang 
6347087cfbeSBarry Smith extern PetscErrorCode  MatScale(Mat,PetscScalar);
6357087cfbeSBarry Smith extern PetscErrorCode  MatShift(Mat,PetscScalar);
6367b80b807SBarry Smith 
6377087cfbeSBarry Smith extern PetscErrorCode  MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping);
6387087cfbeSBarry Smith extern PetscErrorCode  MatSetLocalToGlobalMappingBlock(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping);
6397087cfbeSBarry Smith extern PetscErrorCode  MatGetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*);
6407087cfbeSBarry Smith extern PetscErrorCode  MatGetLocalToGlobalMappingBlock(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*);
6417087cfbeSBarry Smith extern PetscErrorCode  MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
6427087cfbeSBarry Smith extern PetscErrorCode  MatZeroRowsLocalIS(Mat,IS,PetscScalar,Vec,Vec);
6437087cfbeSBarry Smith extern PetscErrorCode  MatZeroRowsColumnsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
6447087cfbeSBarry Smith extern PetscErrorCode  MatZeroRowsColumnsLocalIS(Mat,IS,PetscScalar,Vec,Vec);
6457087cfbeSBarry Smith extern PetscErrorCode  MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
6467087cfbeSBarry Smith extern PetscErrorCode  MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
647052efed2SBarry Smith 
6487087cfbeSBarry Smith extern PetscErrorCode  MatStashSetInitialSize(Mat,PetscInt,PetscInt);
6497087cfbeSBarry Smith extern PetscErrorCode  MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
65090f02eecSBarry Smith 
6517087cfbeSBarry Smith extern PetscErrorCode  MatInterpolate(Mat,Vec,Vec);
6527087cfbeSBarry Smith extern PetscErrorCode  MatInterpolateAdd(Mat,Vec,Vec,Vec);
653ba966639SSatish Balay PetscPolymorphicSubroutine(MatInterpolateAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
6547087cfbeSBarry Smith extern PetscErrorCode  MatRestrict(Mat,Vec,Vec);
6557087cfbeSBarry Smith extern PetscErrorCode  MatGetVecs(Mat,Vec*,Vec*);
6567087cfbeSBarry Smith extern PetscErrorCode  MatGetRedundantMatrix(Mat,PetscInt,MPI_Comm,PetscInt,MatReuse,Mat*);
6577087cfbeSBarry Smith extern PetscErrorCode  MatGetMultiProcBlock(Mat,MPI_Comm,Mat*);
658fc9bc008SSatish Balay extern PetscErrorCode  MatFindZeroDiagonals(Mat,IS*);
659bd481603SBarry Smith 
660bd481603SBarry Smith /*MC
6611620fd73SBarry Smith    MatSetValue - Set a single entry into a matrix.
6621620fd73SBarry Smith 
6631620fd73SBarry Smith    Not collective
6641620fd73SBarry Smith 
6651620fd73SBarry Smith    Input Parameters:
6661620fd73SBarry Smith +  m - the matrix
6671620fd73SBarry Smith .  row - the row location of the entry
6681620fd73SBarry Smith .  col - the column location of the entry
6691620fd73SBarry Smith .  value - the value to insert
6701620fd73SBarry Smith -  mode - either INSERT_VALUES or ADD_VALUES
6711620fd73SBarry Smith 
6721620fd73SBarry Smith    Notes:
6731620fd73SBarry Smith    For efficiency one should use MatSetValues() and set several or many
6741620fd73SBarry Smith    values simultaneously if possible.
6751620fd73SBarry Smith 
6761620fd73SBarry Smith    Level: beginner
6771620fd73SBarry Smith 
6781620fd73SBarry Smith .seealso: MatSetValues(), MatSetValueLocal()
6791620fd73SBarry Smith M*/
6801620fd73SBarry 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);}
6811620fd73SBarry Smith 
6821620fd73SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatGetValue(Mat v,PetscInt i,PetscInt j,PetscScalar *va) {return MatGetValues(v,1,&i,1,&j,va);}
6831620fd73SBarry Smith 
6841620fd73SBarry 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);}
6851620fd73SBarry Smith 
6861620fd73SBarry Smith /*MC
6870d2a7ff2SSatish Balay    MatPreallocateInitialize - Begins the block of code that will count the number of nonzeros per
688bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
689bd481603SBarry Smith 
690bd481603SBarry Smith    Synopsis:
691c1ac3661SBarry Smith    PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)
692bd481603SBarry Smith 
693bd481603SBarry Smith    Collective on MPI_Comm
694bd481603SBarry Smith 
695bd481603SBarry Smith    Input Parameters:
696bd481603SBarry Smith +  comm - the communicator that will share the eventually allocated matrix
697859fcb39SBarry Smith .  nrows - the number of LOCAL rows in the matrix
698859fcb39SBarry Smith -  ncols - the number of LOCAL columns in the matrix
699bd481603SBarry Smith 
700bd481603SBarry Smith    Output Parameters:
701bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
702bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
703bd481603SBarry Smith 
704bd481603SBarry Smith 
705bd481603SBarry Smith    Level: intermediate
706bd481603SBarry Smith 
707bd481603SBarry Smith    Notes:
7080598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
709bd481603SBarry Smith 
7101d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
711bd481603SBarry Smith 
712bd481603SBarry Smith    Use MatPreallocateInitializeSymmetric() for symmetric matrices (MPISBAIJ matrices)
713bd481603SBarry Smith 
7141620fd73SBarry Smith    This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize().
7151620fd73SBarry Smith 
716bd481603SBarry Smith   Concepts: preallocation^Matrix
717bd481603SBarry Smith 
718bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
719bd481603SBarry Smith           MatPreallocateInitializeSymmetric(), MatPreallocateSymmetricSetLocal()
720bd481603SBarry Smith M*/
721c4f061fbSSatish Balay #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \
7227c922b88SBarry Smith { \
723a89bc211SBarry Smith   PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__start,__end; \
724a89bc211SBarry Smith   _4_ierr = PetscMalloc2(__nrows,PetscInt,&dnz,__nrows,PetscInt,&onz);CHKERRQ(_4_ierr); \
725a89bc211SBarry Smith   _4_ierr = PetscMemzero(dnz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
726a89bc211SBarry Smith   _4_ierr = PetscMemzero(onz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
727c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ctmp;\
728a89bc211SBarry Smith   _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows;
7297c922b88SBarry Smith 
730bd481603SBarry Smith /*MC
7310d2a7ff2SSatish Balay    MatPreallocateSymmetricInitialize - Begins the block of code that will count the number of nonzeros per
732bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
733bd481603SBarry Smith 
734bd481603SBarry Smith    Synopsis:
735c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)
736bd481603SBarry Smith 
737bd481603SBarry Smith    Collective on MPI_Comm
738bd481603SBarry Smith 
739bd481603SBarry Smith    Input Parameters:
740bd481603SBarry Smith +  comm - the communicator that will share the eventually allocated matrix
741859fcb39SBarry Smith .  nrows - the number of LOCAL rows in the matrix
742859fcb39SBarry Smith -  ncols - the number of LOCAL columns in the matrix
743bd481603SBarry Smith 
744bd481603SBarry Smith    Output Parameters:
745bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
746bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
747bd481603SBarry Smith 
748bd481603SBarry Smith 
749bd481603SBarry Smith    Level: intermediate
750bd481603SBarry Smith 
751bd481603SBarry Smith    Notes:
7520598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
753bd481603SBarry Smith 
7541d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
755bd481603SBarry Smith 
7561620fd73SBarry Smith    This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize().
7571620fd73SBarry Smith 
758bd481603SBarry Smith   Concepts: preallocation^Matrix
759bd481603SBarry Smith 
760bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
761bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal()
762bd481603SBarry Smith M*/
763222b16d4SBarry Smith #define MatPreallocateSymmetricInitialize(comm,nrows,ncols,dnz,onz) 0; \
764222b16d4SBarry Smith { \
765a89bc211SBarry Smith   PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__end; \
766a89bc211SBarry Smith   _4_ierr = PetscMalloc2(__nrows,PetscInt,&dnz,__nrows,PetscInt,&onz);CHKERRQ(_4_ierr); \
767a89bc211SBarry Smith   _4_ierr = PetscMemzero(dnz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
768a89bc211SBarry Smith   _4_ierr = PetscMemzero(onz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
769c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr);\
770a89bc211SBarry Smith   _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows;
771222b16d4SBarry Smith 
772bd481603SBarry Smith /*MC
773bd481603SBarry Smith    MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
774bd481603SBarry Smith        inserted using a local number of the rows and columns
775bd481603SBarry Smith 
776bd481603SBarry Smith    Synopsis:
777c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
778bd481603SBarry Smith 
779bd481603SBarry Smith    Not Collective
780bd481603SBarry Smith 
781bd481603SBarry Smith    Input Parameters:
782784ac674SJed Brown +  map - the row mapping from local numbering to global numbering
783bd481603SBarry Smith .  nrows - the number of rows indicated
7841d73ed98SBarry Smith .  rows - the indices of the rows
785784ac674SJed Brown .  cmap - the column mapping from local to global numbering
786bd481603SBarry Smith .  ncols - the number of columns in the matrix
787bd481603SBarry Smith .  cols - the columns indicated
788bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
789bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
790bd481603SBarry Smith 
791bd481603SBarry Smith 
792bd481603SBarry Smith    Level: intermediate
793bd481603SBarry Smith 
794bd481603SBarry Smith    Notes:
7950598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
796bd481603SBarry Smith 
7971d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
798bd481603SBarry Smith 
799bd481603SBarry Smith   Concepts: preallocation^Matrix
800bd481603SBarry Smith 
801bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
802bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal()
803bd481603SBarry Smith M*/
804784ac674SJed Brown #define MatPreallocateSetLocal(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \
805c4f061fbSSatish Balay {\
806c1ac3661SBarry Smith   PetscInt __l;\
807784ac674SJed Brown   _4_ierr = ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\
808784ac674SJed Brown   _4_ierr = ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\
809c4f061fbSSatish Balay   for (__l=0;__l<nrows;__l++) {\
810ef66eb69SBarry Smith     _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
811c4f061fbSSatish Balay   }\
812c4f061fbSSatish Balay }
813c4f061fbSSatish Balay 
814bd481603SBarry Smith /*MC
815bd481603SBarry Smith    MatPreallocateSymmetricSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
816bd481603SBarry Smith        inserted using a local number of the rows and columns
817bd481603SBarry Smith 
818bd481603SBarry Smith    Synopsis:
819c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
820bd481603SBarry Smith 
821bd481603SBarry Smith    Not Collective
822bd481603SBarry Smith 
823bd481603SBarry Smith    Input Parameters:
824bd481603SBarry Smith +  map - the mapping between local numbering and global numbering
825bd481603SBarry Smith .  nrows - the number of rows indicated
8261d73ed98SBarry Smith .  rows - the indices of the rows
827bd481603SBarry Smith .  ncols - the number of columns in the matrix
828bd481603SBarry Smith .  cols - the columns indicated
829bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
830bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
831bd481603SBarry Smith 
832bd481603SBarry Smith 
833bd481603SBarry Smith    Level: intermediate
834bd481603SBarry Smith 
835bd481603SBarry Smith    Notes:
8360598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
837bd481603SBarry Smith 
838bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
839bd481603SBarry Smith 
840bd481603SBarry Smith   Concepts: preallocation^Matrix
841bd481603SBarry Smith 
842bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
843bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
844bd481603SBarry Smith M*/
845d3d32019SBarry Smith #define MatPreallocateSymmetricSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\
846d3d32019SBarry Smith {\
847c1ac3661SBarry Smith   PetscInt __l;\
848d3d32019SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\
849d3d32019SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\
850d3d32019SBarry Smith   for (__l=0;__l<nrows;__l++) {\
851d3d32019SBarry Smith     _4_ierr = MatPreallocateSymmetricSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
852d3d32019SBarry Smith   }\
853d3d32019SBarry Smith }
854d3d32019SBarry Smith 
855bd481603SBarry Smith /*MC
856bd481603SBarry Smith    MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
857bd481603SBarry Smith        inserted using a local number of the rows and columns
858bd481603SBarry Smith 
859bd481603SBarry Smith    Synopsis:
860c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
861bd481603SBarry Smith 
862bd481603SBarry Smith    Not Collective
863bd481603SBarry Smith 
864bd481603SBarry Smith    Input Parameters:
86564075487SBarry Smith +  row - the row
866bd481603SBarry Smith .  ncols - the number of columns in the matrix
867a50f8bf6SHong Zhang -  cols - the columns indicated
868a50f8bf6SHong Zhang 
869a50f8bf6SHong Zhang    Output Parameters:
870a50f8bf6SHong Zhang +  dnz - the array that will be passed to the matrix preallocation routines
871bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
872bd481603SBarry Smith 
873bd481603SBarry Smith 
874bd481603SBarry Smith    Level: intermediate
875bd481603SBarry Smith 
876bd481603SBarry Smith    Notes:
8770598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
878bd481603SBarry Smith 
879bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
880bd481603SBarry Smith 
8811620fd73SBarry Smith    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
8821620fd73SBarry Smith 
883bd481603SBarry Smith   Concepts: preallocation^Matrix
884bd481603SBarry Smith 
885bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
886bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
887bd481603SBarry Smith M*/
888c4f061fbSSatish Balay #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\
889c1ac3661SBarry Smith { PetscInt __i; \
890e32f2f54SBarry 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);\
891e32f2f54SBarry 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);\
8927c922b88SBarry Smith   for (__i=0; __i<nc; __i++) {\
89364075487SBarry Smith     if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \
8947cc688d7SBarry Smith     else dnz[row - __rstart]++;\
8957c922b88SBarry Smith   }\
8967c922b88SBarry Smith }
8977c922b88SBarry Smith 
898bd481603SBarry Smith /*MC
899bd481603SBarry Smith    MatPreallocateSymmetricSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
900bd481603SBarry Smith        inserted using a local number of the rows and columns
901bd481603SBarry Smith 
902bd481603SBarry Smith    Synopsis:
903c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
904bd481603SBarry Smith 
905bd481603SBarry Smith    Not Collective
906bd481603SBarry Smith 
907bd481603SBarry Smith    Input Parameters:
908bd481603SBarry Smith +  nrows - the number of rows indicated
9091d73ed98SBarry Smith .  rows - the indices of the rows
910bd481603SBarry Smith .  ncols - the number of columns in the matrix
911bd481603SBarry Smith .  cols - the columns indicated
912bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
913bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
914bd481603SBarry Smith 
915bd481603SBarry Smith 
916bd481603SBarry Smith    Level: intermediate
917bd481603SBarry Smith 
918bd481603SBarry Smith    Notes:
9190598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
920bd481603SBarry Smith 
921bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
922bd481603SBarry Smith 
9231620fd73SBarry Smith    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
9241620fd73SBarry Smith 
925bd481603SBarry Smith   Concepts: preallocation^Matrix
926bd481603SBarry Smith 
927bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
928bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
929bd481603SBarry Smith M*/
930d3d32019SBarry Smith #define MatPreallocateSymmetricSet(row,nc,cols,dnz,onz) 0;\
931c1ac3661SBarry Smith { PetscInt __i; \
932d3d32019SBarry Smith   for (__i=0; __i<nc; __i++) {\
933d3d32019SBarry Smith     if (cols[__i] >= __end) onz[row - __rstart]++; \
934d3d32019SBarry Smith     else if (cols[__i] >= row) dnz[row - __rstart]++;\
935d3d32019SBarry Smith   }\
936d3d32019SBarry Smith }
937d3d32019SBarry Smith 
938bd481603SBarry Smith /*MC
93916371a99SBarry Smith    MatPreallocateLocation -  An alternative to MatPreallocationSet() that puts the nonzero locations into the matrix if it exists
94016371a99SBarry Smith 
94116371a99SBarry Smith    Synopsis:
94216371a99SBarry Smith    PetscErrorCode MatPreallocateLocations(Mat A,PetscInt row,PetscInt ncols,PetscInt *cols,PetscInt *dnz,PetscInt *onz)
94316371a99SBarry Smith 
94416371a99SBarry Smith    Not Collective
94516371a99SBarry Smith 
94616371a99SBarry Smith    Input Parameters:
94716371a99SBarry Smith .  A - matrix
94816371a99SBarry Smith .  row - row where values exist (must be local to this process)
94916371a99SBarry Smith .  ncols - number of columns
95016371a99SBarry Smith .  cols - columns with nonzeros
95116371a99SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
95216371a99SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
95316371a99SBarry Smith 
95416371a99SBarry Smith 
95516371a99SBarry Smith    Level: intermediate
95616371a99SBarry Smith 
95716371a99SBarry Smith    Notes:
9580598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
95916371a99SBarry Smith 
96016371a99SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
96116371a99SBarry Smith 
96216371a99SBarry Smith    This is a MACRO not a function because it uses a bunch of variables private to the MatPreallocation.... routines.
96316371a99SBarry Smith 
96416371a99SBarry Smith   Concepts: preallocation^Matrix
96516371a99SBarry Smith 
96616371a99SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
96716371a99SBarry Smith           MatPreallocateSymmetricInitialize(), MatPreallocateSymmetricSetLocal()
96816371a99SBarry Smith M*/
96916371a99SBarry 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);}
97016371a99SBarry Smith 
97116371a99SBarry Smith 
97216371a99SBarry Smith /*MC
9730d2a7ff2SSatish Balay    MatPreallocateFinalize - Ends the block of code that will count the number of nonzeros per
974bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
975bd481603SBarry Smith 
976bd481603SBarry Smith    Synopsis:
977c1ac3661SBarry Smith    PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz)
978bd481603SBarry Smith 
979bd481603SBarry Smith    Collective on MPI_Comm
980bd481603SBarry Smith 
981bd481603SBarry Smith    Input Parameters:
98216371a99SBarry Smith +  dnz - the array that was be passed to the matrix preallocation routines
983bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
984bd481603SBarry Smith 
985bd481603SBarry Smith 
986bd481603SBarry Smith    Level: intermediate
987bd481603SBarry Smith 
988bd481603SBarry Smith    Notes:
9890598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
990bd481603SBarry Smith 
991bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
992bd481603SBarry Smith 
9931620fd73SBarry Smith    This is a MACRO not a function because it closes the { started in MatPreallocateInitialize().
9941620fd73SBarry Smith 
995bd481603SBarry Smith   Concepts: preallocation^Matrix
996bd481603SBarry Smith 
997bd481603SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
998bd481603SBarry Smith           MatPreallocateSymmetricInitialize(), MatPreallocateSymmetricSetLocal()
999bd481603SBarry Smith M*/
1000a89bc211SBarry Smith #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree2(dnz,onz);CHKERRQ(_4_ierr);}
10017c922b88SBarry Smith 
1002bd481603SBarry Smith 
1003bd481603SBarry Smith 
10047b80b807SBarry Smith /* Routines unique to particular data structures */
10058fe8eb27SJed Brown extern PetscErrorCode  MatShellGetContext(Mat,void *);
1006ba966639SSatish Balay PetscPolymorphicFunction(MatShellGetContext,(Mat A),(A,&t),void*,t)
1007698d4c6aSKris Buschelman 
10087087cfbeSBarry Smith extern PetscErrorCode  MatInodeAdjustForInodes(Mat,IS*,IS*);
10097087cfbeSBarry Smith extern PetscErrorCode  MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *);
1010698d4c6aSKris Buschelman 
10117087cfbeSBarry Smith extern PetscErrorCode  MatSeqAIJSetColumnIndices(Mat,PetscInt[]);
10127087cfbeSBarry Smith extern PetscErrorCode  MatSeqBAIJSetColumnIndices(Mat,PetscInt[]);
10137087cfbeSBarry Smith extern PetscErrorCode  MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
10147087cfbeSBarry Smith extern PetscErrorCode  MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
10157087cfbeSBarry Smith extern PetscErrorCode  MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
10167b80b807SBarry Smith 
1017a96a251dSBarry Smith #define MAT_SKIP_ALLOCATION -4
1018a96a251dSBarry Smith 
10197087cfbeSBarry Smith extern PetscErrorCode  MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
1020ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz))
10217087cfbeSBarry Smith extern PetscErrorCode  MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
1022ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqSBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz))
10237087cfbeSBarry Smith extern PetscErrorCode  MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]);
1024ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqAIJSetPreallocation,(Mat A,const PetscInt nnz[]),(A,0,nnz))
1025273d9f13SBarry Smith 
10267087cfbeSBarry Smith extern PetscErrorCode  MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
1027ba966639SSatish Balay PetscPolymorphicSubroutine(MatMPIBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[],const PetscInt onz[]),(A,bs,0,nnz,0,onz))
10287087cfbeSBarry Smith extern PetscErrorCode  MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
10297087cfbeSBarry Smith extern PetscErrorCode  MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
10307087cfbeSBarry Smith extern PetscErrorCode  MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []);
10317087cfbeSBarry Smith extern PetscErrorCode  MatSeqBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
10327087cfbeSBarry Smith extern PetscErrorCode  MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]);
10337087cfbeSBarry Smith extern PetscErrorCode  MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
10347087cfbeSBarry Smith extern PetscErrorCode  MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]);
10357087cfbeSBarry Smith extern PetscErrorCode  MatMPIDenseSetPreallocation(Mat,PetscScalar[]);
10367087cfbeSBarry Smith extern PetscErrorCode  MatSeqDenseSetPreallocation(Mat,PetscScalar[]);
10377087cfbeSBarry Smith extern PetscErrorCode  MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,PetscInt*[]);
10387087cfbeSBarry Smith extern PetscErrorCode  MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,PetscInt*[]);
10397087cfbeSBarry Smith extern PetscErrorCode  MatAdicSetLocalFunction(Mat,void (*)(void));
1040273d9f13SBarry Smith 
10417087cfbeSBarry Smith extern PetscErrorCode  MatSeqDenseSetLDA(Mat,PetscInt);
10427087cfbeSBarry Smith extern PetscErrorCode  MatDenseGetLocalMatrix(Mat,Mat*);
10431b807ce4Svictorle 
10447087cfbeSBarry Smith extern PetscErrorCode  MatStoreValues(Mat);
10457087cfbeSBarry Smith extern PetscErrorCode  MatRetrieveValues(Mat);
10462e8a6d31SBarry Smith 
10477087cfbeSBarry Smith extern PetscErrorCode  MatDAADSetCtx(Mat,void*);
10483a7fca6bSBarry Smith 
1049b3a44c85SBarry Smith extern PetscErrorCode  MatFindNonzeroRows(Mat,IS*);
10507b80b807SBarry Smith /*
10517b80b807SBarry Smith   These routines are not usually accessed directly, rather solving is
105294b7f48cSBarry Smith   done through the KSP and PC interfaces.
10537b80b807SBarry Smith */
10547b80b807SBarry Smith 
1055d9274352SBarry Smith /*E
1056d9274352SBarry Smith     MatOrderingType - String with the name of a PETSc matrix ordering or the creation function
1057d9274352SBarry Smith        with an optional dynamic library name, for example
1058d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:orderingcreate()
1059d9274352SBarry Smith 
1060d9274352SBarry Smith    Level: beginner
1061d9274352SBarry Smith 
1062e5a9bf91SBarry Smith    Cannot use const because the PC objects manipulate the string
1063e5a9bf91SBarry Smith 
1064d9274352SBarry Smith .seealso: MatGetOrdering()
1065d9274352SBarry Smith E*/
10663eaad2caSSatish Balay #define MatOrderingType char*
10672692d6eeSBarry Smith #define MATORDERINGNATURAL     "natural"
10682692d6eeSBarry Smith #define MATORDERINGND          "nd"
10692692d6eeSBarry Smith #define MATORDERING1WD         "1wd"
10702692d6eeSBarry Smith #define MATORDERINGRCM         "rcm"
10712692d6eeSBarry Smith #define MATORDERINGQMD         "qmd"
10722692d6eeSBarry Smith #define MATORDERINGROWLENGTH   "rowlength"
1073312e372aSBarry Smith #define MATORDERINGAMD         "amd"            /* only works if UMFPACK is installed with PETSc */
1074b12f92e5SBarry Smith 
10757087cfbeSBarry Smith extern PetscErrorCode  MatGetOrdering(Mat,const MatOrderingType,IS*,IS*);
10767087cfbeSBarry Smith extern PetscErrorCode  MatGetOrderingList(PetscFList *list);
10777087cfbeSBarry Smith extern PetscErrorCode  MatOrderingRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,const MatOrderingType,IS*,IS*));
107830de9b25SBarry Smith 
107930de9b25SBarry Smith /*MC
10801890ba74SBarry Smith    MatOrderingRegisterDynamic - Adds a new sparse matrix ordering to the matrix package.
108130de9b25SBarry Smith 
108230de9b25SBarry Smith    Synopsis:
10831890ba74SBarry Smith    PetscErrorCode MatOrderingRegisterDynamic(const char *name_ordering,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatOrdering))
108430de9b25SBarry Smith 
108530de9b25SBarry Smith    Not Collective
108630de9b25SBarry Smith 
108730de9b25SBarry Smith    Input Parameters:
10882692d6eeSBarry Smith +  sname - name of ordering (for example MATORDERINGND)
108930de9b25SBarry Smith .  path - location of library where creation routine is
109030de9b25SBarry Smith .  name - name of function that creates the ordering type,a string
109130de9b25SBarry Smith -  function - function pointer that creates the ordering
109230de9b25SBarry Smith 
109330de9b25SBarry Smith    Level: developer
109430de9b25SBarry Smith 
109530de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
109630de9b25SBarry Smith    is ignored.
109730de9b25SBarry Smith 
109830de9b25SBarry Smith    Sample usage:
109930de9b25SBarry Smith .vb
110030de9b25SBarry Smith    MatOrderingRegisterDynamic("my_order",/home/username/my_lib/lib/libO/solaris/mylib.a,
110130de9b25SBarry Smith                "MyOrder",MyOrder);
110230de9b25SBarry Smith .ve
110330de9b25SBarry Smith 
110430de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
110530de9b25SBarry Smith $     MatOrderingSetType(part,"my_order)
110630de9b25SBarry Smith    or at runtime via the option
11072401956bSBarry Smith $     -pc_factor_mat_ordering_type my_order
110830de9b25SBarry Smith 
1109ab901514SBarry Smith    ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values.
111030de9b25SBarry Smith 
111130de9b25SBarry Smith .keywords: matrix, ordering, register
111230de9b25SBarry Smith 
111330de9b25SBarry Smith .seealso: MatOrderingRegisterDestroy(), MatOrderingRegisterAll()
111430de9b25SBarry Smith M*/
1115aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1116f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,0)
1117b12f92e5SBarry Smith #else
1118f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,d)
1119b12f92e5SBarry Smith #endif
112030de9b25SBarry Smith 
11217087cfbeSBarry Smith extern PetscErrorCode  MatOrderingRegisterDestroy(void);
11227087cfbeSBarry Smith extern PetscErrorCode  MatOrderingRegisterAll(const char[]);
1123ace3abfcSBarry Smith extern PetscBool  MatOrderingRegisterAllCalled;
1124b0a32e0cSBarry Smith extern PetscFList MatOrderingList;
1125d4fbbf0eSBarry Smith 
11267087cfbeSBarry Smith extern PetscErrorCode  MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS);
1127a2ce50c7SBarry Smith 
1128d91e6319SBarry Smith /*S
1129d90ac83dSHong Zhang     MatFactorShiftType - Numeric Shift.
1130d90ac83dSHong Zhang 
1131d90ac83dSHong Zhang    Level: beginner
1132d90ac83dSHong Zhang 
1133d90ac83dSHong Zhang S*/
1134d90ac83dSHong Zhang typedef enum {MAT_SHIFT_NONE,MAT_SHIFT_NONZERO,MAT_SHIFT_POSITIVE_DEFINITE,MAT_SHIFT_INBLOCKS} MatFactorShiftType;
1135d90ac83dSHong Zhang extern const char *MatFactorShiftTypes[];
1136d90ac83dSHong Zhang 
1137d90ac83dSHong Zhang /*S
11382401956bSBarry Smith    MatFactorInfo - Data passed into the matrix factorization routines
1139b00f7748SHong Zhang 
114061cad860SBarry Smith    In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE, that is use
114161cad860SBarry Smith $     MatFactorInfo  info(MAT_FACTORINFO_SIZE)
1142b00f7748SHong Zhang 
114315e8a5b3SHong Zhang    Notes: These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC.
1144b00f7748SHong Zhang 
114561cad860SBarry Smith       You can use MatFactorInfoInitialize() to set default values.
114661cad860SBarry Smith 
1147b00f7748SHong Zhang    Level: developer
1148b00f7748SHong Zhang 
1149d7d82daaSBarry Smith .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(),
1150d7d82daaSBarry Smith           MatFactorInfoInitialize()
1151b00f7748SHong Zhang 
1152b00f7748SHong Zhang S*/
1153b00f7748SHong Zhang typedef struct {
115415e8a5b3SHong Zhang   PetscReal     diagonal_fill;  /* force diagonal to fill in if initially not filled */
115585317021SBarry Smith   PetscReal     usedt;
115615e8a5b3SHong Zhang   PetscReal     dt;             /* drop tolerance */
1157b00f7748SHong Zhang   PetscReal     dtcol;          /* tolerance for pivoting */
115815e8a5b3SHong Zhang   PetscReal     dtcount;        /* maximum nonzeros to be allowed per row */
115967e28bfeSBarry Smith   PetscReal     fill;           /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */
1160348344bbSBarry Smith   PetscReal     levels;         /* ICC/ILU(levels) */
1161bcd9e38bSBarry Smith   PetscReal     pivotinblocks;  /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0
1162bcd9e38bSBarry Smith                                    factorization may be faster if do not pivot */
116315e8a5b3SHong Zhang   PetscReal     zeropivot;      /* pivot is called zero if less than this */
1164f4db908eSBarry Smith   PetscReal     shifttype;      /* type of shift added to matrix factor to prevent zero pivots */
1165f4db908eSBarry Smith   PetscReal     shiftamount;     /* how large the shift is */
116615e8a5b3SHong Zhang } MatFactorInfo;
1167ffa6d0a5SLois Curfman McInnes 
11687087cfbeSBarry Smith extern PetscErrorCode  MatFactorInfoInitialize(MatFactorInfo*);
11697087cfbeSBarry Smith extern PetscErrorCode  MatCholeskyFactor(Mat,IS,const MatFactorInfo*);
11707087cfbeSBarry Smith extern PetscErrorCode  MatCholeskyFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
11717087cfbeSBarry Smith extern PetscErrorCode  MatCholeskyFactorNumeric(Mat,Mat,const MatFactorInfo*);
11727087cfbeSBarry Smith extern PetscErrorCode  MatLUFactor(Mat,IS,IS,const MatFactorInfo*);
11737087cfbeSBarry Smith extern PetscErrorCode  MatILUFactor(Mat,IS,IS,const MatFactorInfo*);
11747087cfbeSBarry Smith extern PetscErrorCode  MatLUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
11757087cfbeSBarry Smith extern PetscErrorCode  MatILUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
11767087cfbeSBarry Smith extern PetscErrorCode  MatICCFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
11777087cfbeSBarry Smith extern PetscErrorCode  MatICCFactor(Mat,IS,const MatFactorInfo*);
11787087cfbeSBarry Smith extern PetscErrorCode  MatLUFactorNumeric(Mat,Mat,const MatFactorInfo*);
11797087cfbeSBarry Smith extern PetscErrorCode  MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*);
11807087cfbeSBarry Smith extern PetscErrorCode  MatSolve(Mat,Vec,Vec);
11817087cfbeSBarry Smith extern PetscErrorCode  MatForwardSolve(Mat,Vec,Vec);
11827087cfbeSBarry Smith extern PetscErrorCode  MatBackwardSolve(Mat,Vec,Vec);
11837087cfbeSBarry Smith extern PetscErrorCode  MatSolveAdd(Mat,Vec,Vec,Vec);
11847087cfbeSBarry Smith extern PetscErrorCode  MatSolveTranspose(Mat,Vec,Vec);
11857087cfbeSBarry Smith extern PetscErrorCode  MatSolveTransposeAdd(Mat,Vec,Vec,Vec);
11867087cfbeSBarry Smith extern PetscErrorCode  MatSolves(Mat,Vecs,Vecs);
11878ed539a5SBarry Smith 
11887087cfbeSBarry Smith extern PetscErrorCode  MatSetUnfactored(Mat);
1189bb5a7306SBarry Smith 
1190d91e6319SBarry Smith /*E
1191d91e6319SBarry Smith     MatSORType - What type of (S)SOR to perform
1192bb1eb677SSatish Balay 
1193d91e6319SBarry Smith     Level: beginner
1194d91e6319SBarry Smith 
1195d9274352SBarry Smith    May be bitwise ORd together
1196d9274352SBarry Smith 
1197d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
1198d91e6319SBarry Smith 
11994e7234bfSBarry Smith    MatSORType may be bitwise ORd together, so do not change the numbers
12004e7234bfSBarry Smith 
120141f059aeSBarry Smith .seealso: MatSOR()
1202d91e6319SBarry Smith E*/
1203ee50ffe9SBarry Smith typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3,
1204ee50ffe9SBarry Smith               SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8,
1205ee50ffe9SBarry Smith               SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16,
120684cb2905SBarry Smith               SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType;
12077087cfbeSBarry Smith extern PetscErrorCode  MatSOR(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
12088ed539a5SBarry Smith 
1209d4fbbf0eSBarry Smith /*
1210639f9d9dSBarry Smith     These routines are for efficiently computing Jacobians via finite differences.
1211639f9d9dSBarry Smith */
1212b12f92e5SBarry Smith 
1213d9274352SBarry Smith /*E
1214d9274352SBarry Smith     MatColoringType - String with the name of a PETSc matrix coloring or the creation function
1215d9274352SBarry Smith        with an optional dynamic library name, for example
1216d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:coloringcreate()
1217d9274352SBarry Smith 
1218d9274352SBarry Smith    Level: beginner
1219d9274352SBarry Smith 
1220d9274352SBarry Smith .seealso: MatGetColoring()
1221d9274352SBarry Smith E*/
1222a313700dSBarry Smith #define MatColoringType char*
12232692d6eeSBarry Smith #define MATCOLORINGNATURAL "natural"
12242692d6eeSBarry Smith #define MATCOLORINGSL      "sl"
12252692d6eeSBarry Smith #define MATCOLORINGLF      "lf"
12262692d6eeSBarry Smith #define MATCOLORINGID      "id"
1227b12f92e5SBarry Smith 
12287087cfbeSBarry Smith extern PetscErrorCode  MatGetColoring(Mat,const MatColoringType,ISColoring*);
12297087cfbeSBarry Smith extern PetscErrorCode  MatColoringRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,MatColoringType,ISColoring *));
123030de9b25SBarry Smith 
123130de9b25SBarry Smith /*MC
123230de9b25SBarry Smith    MatColoringRegisterDynamic - Adds a new sparse matrix coloring to the
123330de9b25SBarry Smith                                matrix package.
123430de9b25SBarry Smith 
123530de9b25SBarry Smith    Synopsis:
12361890ba74SBarry Smith    PetscErrorCode MatColoringRegisterDynamic(const char *name_coloring,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatColoring))
123730de9b25SBarry Smith 
123830de9b25SBarry Smith    Not Collective
123930de9b25SBarry Smith 
124030de9b25SBarry Smith    Input Parameters:
12412692d6eeSBarry Smith +  sname - name of Coloring (for example MATCOLORINGSL)
124230de9b25SBarry Smith .  path - location of library where creation routine is
124330de9b25SBarry Smith .  name - name of function that creates the Coloring type, a string
124430de9b25SBarry Smith -  function - function pointer that creates the coloring
124530de9b25SBarry Smith 
124630de9b25SBarry Smith    Level: developer
124730de9b25SBarry Smith 
124830de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
124930de9b25SBarry Smith    is ignored.
125030de9b25SBarry Smith 
125130de9b25SBarry Smith    Sample usage:
125230de9b25SBarry Smith .vb
125330de9b25SBarry Smith    MatColoringRegisterDynamic("my_color",/home/username/my_lib/lib/libO/solaris/mylib.a,
125430de9b25SBarry Smith                "MyColor",MyColor);
125530de9b25SBarry Smith .ve
125630de9b25SBarry Smith 
125730de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
125830de9b25SBarry Smith $     MatColoringSetType(part,"my_color")
125930de9b25SBarry Smith    or at runtime via the option
126030de9b25SBarry Smith $     -mat_coloring_type my_color
126130de9b25SBarry Smith 
1262ab901514SBarry Smith    $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
126330de9b25SBarry Smith 
126430de9b25SBarry Smith .keywords: matrix, Coloring, register
126530de9b25SBarry Smith 
126630de9b25SBarry Smith .seealso: MatColoringRegisterDestroy(), MatColoringRegisterAll()
126730de9b25SBarry Smith M*/
1268aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1269f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,0)
1270b12f92e5SBarry Smith #else
1271f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,d)
1272b12f92e5SBarry Smith #endif
127330de9b25SBarry Smith 
1274ace3abfcSBarry Smith extern PetscBool  MatColoringRegisterAllCalled;
1275f1144a30SSatish Balay 
12767087cfbeSBarry Smith extern PetscErrorCode  MatColoringRegisterAll(const char[]);
12777087cfbeSBarry Smith extern PetscErrorCode  MatColoringRegisterDestroy(void);
12787087cfbeSBarry Smith extern PetscErrorCode  MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*);
1279639f9d9dSBarry Smith 
1280d9274352SBarry Smith /*S
1281d9274352SBarry Smith      MatFDColoring - Object for computing a sparse Jacobian via finite differences
1282d9274352SBarry Smith         and coloring
1283639f9d9dSBarry Smith 
1284d9274352SBarry Smith    Level: beginner
1285d9274352SBarry Smith 
1286d9274352SBarry Smith   Concepts: coloring, sparse Jacobian, finite differences
1287d9274352SBarry Smith 
1288d9274352SBarry Smith .seealso:  MatFDColoringCreate()
1289d9274352SBarry Smith S*/
1290e2a1c21fSSatish Balay typedef struct _p_MatFDColoring* MatFDColoring;
1291639f9d9dSBarry Smith 
12927087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringCreate(Mat,ISColoring,MatFDColoring *);
1293fcfd50ebSBarry Smith extern PetscErrorCode  MatFDColoringDestroy(MatFDColoring*);
12947087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringView(MatFDColoring,PetscViewer);
12957087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*);
12967087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**);
12977087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal);
12987087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringSetFromOptions(MatFDColoring);
12997087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringApply(Mat,MatFDColoring,Vec,MatStructure*,void *);
13007087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringApplyTS(Mat,MatFDColoring,PetscReal,Vec,MatStructure*,void *);
13017087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringSetF(MatFDColoring,Vec);
13027087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,PetscInt*[]);
1303639f9d9dSBarry Smith /*
13040752156aSBarry Smith     These routines are for partitioning matrices: currently used only
13053eda8832SBarry Smith   for adjacency matrix, MatCreateMPIAdj().
13060752156aSBarry Smith */
1307ca161407SBarry Smith 
1308d9274352SBarry Smith /*S
1309d9274352SBarry Smith      MatPartitioning - Object for managing the partitioning of a matrix or graph
1310d9274352SBarry Smith 
1311d9274352SBarry Smith    Level: beginner
1312d9274352SBarry Smith 
1313d9274352SBarry Smith   Concepts: partitioning
1314d9274352SBarry Smith 
1315743a1026Svictorle .seealso:  MatPartitioningCreate(), MatPartitioningType
1316d9274352SBarry Smith S*/
131791e9ee9fSBarry Smith typedef struct _p_MatPartitioning* MatPartitioning;
1318d9274352SBarry Smith 
1319d9274352SBarry Smith /*E
13205bc6e7ccSvictorle     MatPartitioningType - String with the name of a PETSc matrix partitioning or the creation function
1321d9274352SBarry Smith        with an optional dynamic library name, for example
1322d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:partitioningcreate()
1323d9274352SBarry Smith 
1324d9274352SBarry Smith    Level: beginner
1325d9274352SBarry Smith 
1326b547a13bSvictorle .seealso: MatPartitioningCreate(), MatPartitioning
1327d9274352SBarry Smith E*/
1328a313700dSBarry Smith #define MatPartitioningType char*
13292692d6eeSBarry Smith #define MATPARTITIONINGCURRENT  "current"
13302692d6eeSBarry Smith #define MATPARTITIONINGSQUARE   "square"
13312692d6eeSBarry Smith #define MATPARTITIONINGPARMETIS "parmetis"
13322692d6eeSBarry Smith #define MATPARTITIONINGCHACO    "chaco"
13332692d6eeSBarry Smith #define MATPARTITIONINGPARTY    "party"
133450d91057SBarry Smith #define MATPARTITIONINGPTSCOTCH "ptscotch"
133571306c60Sjroman 
1336ca161407SBarry Smith 
13377087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningCreate(MPI_Comm,MatPartitioning*);
13387087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningSetType(MatPartitioning,const MatPartitioningType);
13397087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningSetNParts(MatPartitioning,PetscInt);
13407087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningSetAdjacency(MatPartitioning,Mat);
13417087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]);
13427087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []);
13437087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningApply(MatPartitioning,IS*);
1344fcfd50ebSBarry Smith extern PetscErrorCode  MatPartitioningDestroy(MatPartitioning*);
13452aabb6bbSBarry Smith 
13467087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatPartitioning));
134730de9b25SBarry Smith 
134830de9b25SBarry Smith /*MC
134930de9b25SBarry Smith    MatPartitioningRegisterDynamic - Adds a new sparse matrix partitioning to the
135030de9b25SBarry Smith    matrix package.
135130de9b25SBarry Smith 
135230de9b25SBarry Smith    Synopsis:
13531890ba74SBarry Smith    PetscErrorCode MatPartitioningRegisterDynamic(const char *name_partitioning,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatPartitioning))
135430de9b25SBarry Smith 
135530de9b25SBarry Smith    Not Collective
135630de9b25SBarry Smith 
135730de9b25SBarry Smith    Input Parameters:
13582692d6eeSBarry Smith +  sname - name of partitioning (for example MATPARTITIONINGCURRENT) or parmetis
135930de9b25SBarry Smith .  path - location of library where creation routine is
136030de9b25SBarry Smith .  name - name of function that creates the partitioning type, a string
136130de9b25SBarry Smith -  function - function pointer that creates the partitioning type
136230de9b25SBarry Smith 
136330de9b25SBarry Smith    Level: developer
136430de9b25SBarry Smith 
136530de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
136630de9b25SBarry Smith    is ignored.
136730de9b25SBarry Smith 
136830de9b25SBarry Smith    Sample usage:
136930de9b25SBarry Smith .vb
137030de9b25SBarry Smith    MatPartitioningRegisterDynamic("my_part",/home/username/my_lib/lib/libO/solaris/mylib.a,
137130de9b25SBarry Smith                "MyPartCreate",MyPartCreate);
137230de9b25SBarry Smith .ve
137330de9b25SBarry Smith 
137430de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
137530de9b25SBarry Smith $     MatPartitioningSetType(part,"my_part")
137630de9b25SBarry Smith    or at runtime via the option
137730de9b25SBarry Smith $     -mat_partitioning_type my_part
137830de9b25SBarry Smith 
1379ab901514SBarry Smith    $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
138030de9b25SBarry Smith 
138130de9b25SBarry Smith .keywords: matrix, partitioning, register
138230de9b25SBarry Smith 
138330de9b25SBarry Smith .seealso: MatPartitioningRegisterDestroy(), MatPartitioningRegisterAll()
138430de9b25SBarry Smith M*/
1385aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1386f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,0)
13872aabb6bbSBarry Smith #else
1388f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,d)
13892aabb6bbSBarry Smith #endif
13902aabb6bbSBarry Smith 
1391ace3abfcSBarry Smith extern PetscBool  MatPartitioningRegisterAllCalled;
1392f1144a30SSatish Balay 
13937087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningRegisterAll(const char[]);
13947087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningRegisterDestroy(void);
13952bad1931SBarry Smith 
13967087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningView(MatPartitioning,PetscViewer);
13977087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningSetFromOptions(MatPartitioning);
13987087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningGetType(MatPartitioning,const MatPartitioningType*);
1399ca161407SBarry Smith 
14007087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningParmetisSetCoarseSequential(MatPartitioning);
14017087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningParmetisGetEdgeCut(MatPartitioning, PetscInt *);
14020752156aSBarry Smith 
1403b6956c12SJose Roman typedef enum { MP_CHACO_MULTILEVEL=1,MP_CHACO_SPECTRAL=2,MP_CHACO_LINEAR=4,MP_CHACO_RANDOM=5,MP_CHACO_SCATTERED=6 } MPChacoGlobalType;
1404b6956c12SJose Roman extern const char *MPChacoGlobalTypes[];
1405b6956c12SJose Roman typedef enum { MP_CHACO_KERNIGHAN=1,MP_CHACO_NONE=2 } MPChacoLocalType;
1406b6956c12SJose Roman extern const char *MPChacoLocalTypes[];
1407b6956c12SJose Roman typedef enum { MP_CHACO_LANCZOS=0,MP_CHACO_RQI=1 } MPChacoEigenType;
1408b6956c12SJose Roman extern const char *MPChacoEigenTypes[];
1409b6956c12SJose Roman 
14107087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningChacoSetGlobal(MatPartitioning,MPChacoGlobalType);
1411b6956c12SJose Roman extern PetscErrorCode  MatPartitioningChacoGetGlobal(MatPartitioning,MPChacoGlobalType*);
14127087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningChacoSetLocal(MatPartitioning,MPChacoLocalType);
1413b6956c12SJose Roman extern PetscErrorCode  MatPartitioningChacoGetLocal(MatPartitioning,MPChacoLocalType*);
14147087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal);
14157087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType);
1416b6956c12SJose Roman extern PetscErrorCode  MatPartitioningChacoGetEigenSolver(MatPartitioning,MPChacoEigenType*);
14177087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningChacoSetEigenTol(MatPartitioning,PetscReal);
1418b6956c12SJose Roman extern PetscErrorCode  MatPartitioningChacoGetEigenTol(MatPartitioning,PetscReal*);
14197087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningChacoSetEigenNumber(MatPartitioning,PetscInt);
1420b6956c12SJose Roman extern PetscErrorCode  MatPartitioningChacoGetEigenNumber(MatPartitioning,PetscInt*);
142171306c60Sjroman 
142271306c60Sjroman #define MP_PARTY_OPT "opt"
142371306c60Sjroman #define MP_PARTY_LIN "lin"
142471306c60Sjroman #define MP_PARTY_SCA "sca"
142571306c60Sjroman #define MP_PARTY_RAN "ran"
142671306c60Sjroman #define MP_PARTY_GBF "gbf"
142771306c60Sjroman #define MP_PARTY_GCF "gcf"
142871306c60Sjroman #define MP_PARTY_BUB "bub"
142971306c60Sjroman #define MP_PARTY_DEF "def"
14307087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningPartySetGlobal(MatPartitioning,const char*);
143171306c60Sjroman #define MP_PARTY_HELPFUL_SETS "hs"
143271306c60Sjroman #define MP_PARTY_KERNIGHAN_LIN "kl"
143371306c60Sjroman #define MP_PARTY_NONE "no"
14347087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningPartySetLocal(MatPartitioning,const char*);
14357087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal);
14367087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningPartySetBipart(MatPartitioning,PetscBool);
14377087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscBool);
143871306c60Sjroman 
143950d91057SBarry Smith typedef enum { MP_PTSCOTCH_QUALITY,MP_PTSCOTCH_SPEED,MP_PTSCOTCH_BALANCE,MP_PTSCOTCH_SAFETY,MP_PTSCOTCH_SCALABILITY } MPPTScotchStrategyType;
144050d91057SBarry Smith extern const char *MPPTScotchStrategyTypes[];
1441e0f1cffaSJose Roman 
144250d91057SBarry Smith extern PetscErrorCode MatPartitioningPTScotchSetImbalance(MatPartitioning,PetscReal);
144350d91057SBarry Smith extern PetscErrorCode MatPartitioningPTScotchGetImbalance(MatPartitioning,PetscReal*);
144450d91057SBarry Smith extern PetscErrorCode MatPartitioningPTScotchSetStrategy(MatPartitioning,MPPTScotchStrategyType);
144550d91057SBarry Smith extern PetscErrorCode MatPartitioningPTScotchGetStrategy(MatPartitioning,MPPTScotchStrategyType*);
144671306c60Sjroman 
144709573ac7SBarry Smith extern PetscErrorCode MatMeshToVertexGraph(Mat,PetscInt,Mat*);
144809573ac7SBarry Smith extern PetscErrorCode MatMeshToCellGraph(Mat,PetscInt,Mat*);
1449591294e4SBarry Smith 
14500752156aSBarry Smith /*
14510a835dfdSSatish Balay     If you add entries here you must also add them to finclude/petscmat.h
1452d4fbbf0eSBarry Smith */
14531c1c02c0SLois Curfman McInnes typedef enum { MATOP_SET_VALUES=0,
14541c1c02c0SLois Curfman McInnes                MATOP_GET_ROW=1,
14551c1c02c0SLois Curfman McInnes                MATOP_RESTORE_ROW=2,
14561c1c02c0SLois Curfman McInnes                MATOP_MULT=3,
14571c1c02c0SLois Curfman McInnes                MATOP_MULT_ADD=4,
14587c922b88SBarry Smith                MATOP_MULT_TRANSPOSE=5,
14597c922b88SBarry Smith                MATOP_MULT_TRANSPOSE_ADD=6,
14601c1c02c0SLois Curfman McInnes                MATOP_SOLVE=7,
14611c1c02c0SLois Curfman McInnes                MATOP_SOLVE_ADD=8,
14627c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE=9,
14637c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE_ADD=10,
14641c1c02c0SLois Curfman McInnes                MATOP_LUFACTOR=11,
14651c1c02c0SLois Curfman McInnes                MATOP_CHOLESKYFACTOR=12,
1466a714c06dSBarry Smith                MATOP_SOR=13,
14671c1c02c0SLois Curfman McInnes                MATOP_TRANSPOSE=14,
14681c1c02c0SLois Curfman McInnes                MATOP_GETINFO=15,
14691c1c02c0SLois Curfman McInnes                MATOP_EQUAL=16,
14701c1c02c0SLois Curfman McInnes                MATOP_GET_DIAGONAL=17,
14711c1c02c0SLois Curfman McInnes                MATOP_DIAGONAL_SCALE=18,
14721c1c02c0SLois Curfman McInnes                MATOP_NORM=19,
14731c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_BEGIN=20,
14741c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_END=21,
1475d519adbfSMatthew Knepley                MATOP_SET_OPTION=22,
1476d519adbfSMatthew Knepley                MATOP_ZERO_ENTRIES=23,
1477d519adbfSMatthew Knepley                MATOP_ZERO_ROWS=24,
1478d519adbfSMatthew Knepley                MATOP_LUFACTOR_SYMBOLIC=25,
1479d519adbfSMatthew Knepley                MATOP_LUFACTOR_NUMERIC=26,
1480d519adbfSMatthew Knepley                MATOP_CHOLESKY_FACTOR_SYMBOLIC=27,
1481d519adbfSMatthew Knepley                MATOP_CHOLESKY_FACTOR_NUMERIC=28,
1482d519adbfSMatthew Knepley                MATOP_SETUP_PREALLOCATION=29,
1483d519adbfSMatthew Knepley                MATOP_ILUFACTOR_SYMBOLIC=30,
1484d519adbfSMatthew Knepley                MATOP_ICCFACTOR_SYMBOLIC=31,
1485d519adbfSMatthew Knepley                MATOP_GET_ARRAY=32,
1486d519adbfSMatthew Knepley                MATOP_RESTORE_ARRAY=33,
1487d519adbfSMatthew Knepley                MATOP_DUPLICATE=34,
1488d519adbfSMatthew Knepley                MATOP_FORWARD_SOLVE=35,
1489d519adbfSMatthew Knepley                MATOP_BACKWARD_SOLVE=36,
1490d519adbfSMatthew Knepley                MATOP_ILUFACTOR=37,
1491d519adbfSMatthew Knepley                MATOP_ICCFACTOR=38,
1492d519adbfSMatthew Knepley                MATOP_AXPY=39,
1493d519adbfSMatthew Knepley                MATOP_GET_SUBMATRICES=40,
1494d519adbfSMatthew Knepley                MATOP_INCREASE_OVERLAP=41,
1495d519adbfSMatthew Knepley                MATOP_GET_VALUES=42,
1496d519adbfSMatthew Knepley                MATOP_COPY=43,
1497d519adbfSMatthew Knepley                MATOP_GET_ROW_MAX=44,
1498d519adbfSMatthew Knepley                MATOP_SCALE=45,
1499d519adbfSMatthew Knepley                MATOP_SHIFT=46,
150035153367SBarry Smith                MATOP_DIAGONAL_SET=47,
1501d519adbfSMatthew Knepley                MATOP_ILUDT_FACTOR=48,
1502d519adbfSMatthew Knepley                MATOP_SET_BLOCK_SIZE=49,
1503d519adbfSMatthew Knepley                MATOP_GET_ROW_IJ=50,
1504d519adbfSMatthew Knepley                MATOP_RESTORE_ROW_IJ=51,
1505d519adbfSMatthew Knepley                MATOP_GET_COLUMN_IJ=52,
1506d519adbfSMatthew Knepley                MATOP_RESTORE_COLUMN_IJ=53,
1507d519adbfSMatthew Knepley                MATOP_FDCOLORING_CREATE=54,
1508d519adbfSMatthew Knepley                MATOP_COLORING_PATCH=55,
1509d519adbfSMatthew Knepley                MATOP_SET_UNFACTORED=56,
1510d519adbfSMatthew Knepley                MATOP_PERMUTE=57,
1511d519adbfSMatthew Knepley                MATOP_SET_VALUES_BLOCKED=58,
1512d519adbfSMatthew Knepley                MATOP_GET_SUBMATRIX=59,
1513d519adbfSMatthew Knepley                MATOP_DESTROY=60,
1514d519adbfSMatthew Knepley                MATOP_VIEW=61,
1515d519adbfSMatthew Knepley                MATOP_CONVERT_FROM=62,
1516d519adbfSMatthew Knepley                MATOP_USE_SCALED_FORM=63,
1517d519adbfSMatthew Knepley                MATOP_SCALE_SYSTEM=64,
1518d519adbfSMatthew Knepley                MATOP_UNSCALE_SYSTEM=65,
1519d519adbfSMatthew Knepley                MATOP_SET_LOCAL_TO_GLOBAL_MAP=66,
1520d519adbfSMatthew Knepley                MATOP_SET_VALUES_LOCAL=67,
1521d519adbfSMatthew Knepley                MATOP_ZERO_ROWS_LOCAL=68,
1522d519adbfSMatthew Knepley                MATOP_GET_ROW_MAX_ABS=69,
1523d519adbfSMatthew Knepley                MATOP_GET_ROW_MIN_ABS=70,
1524d519adbfSMatthew Knepley                MATOP_CONVERT=71,
1525d519adbfSMatthew Knepley                MATOP_SET_COLORING=72,
1526d519adbfSMatthew Knepley                MATOP_SET_VALUES_ADIC=73,
1527d519adbfSMatthew Knepley                MATOP_SET_VALUES_ADIFOR=74,
1528d519adbfSMatthew Knepley                MATOP_FD_COLORING_APPLY=75,
1529d519adbfSMatthew Knepley                MATOP_SET_FROM_OPTIONS=76,
1530d519adbfSMatthew Knepley                MATOP_MULT_CON=77,
1531d519adbfSMatthew Knepley                MATOP_MULT_TRANSPOSE_CON=78,
1532d519adbfSMatthew Knepley                MATOP_PERMUTE_SPARSIFY=79,
1533d519adbfSMatthew Knepley                MATOP_MULT_MULTIPLE=80,
1534d519adbfSMatthew Knepley                MATOP_SOLVE_MULTIPLE=81,
1535d519adbfSMatthew Knepley                MATOP_GET_INERTIA=82,
1536d519adbfSMatthew Knepley                MATOP_LOAD=83,
1537d519adbfSMatthew Knepley                MATOP_IS_SYMMETRIC=84,
1538d519adbfSMatthew Knepley                MATOP_IS_HERMITIAN=85,
1539d519adbfSMatthew Knepley                MATOP_IS_STRUCTURALLY_SYMMETRIC=86,
154041f059aeSBarry Smith                MATOP_DUMMY=87,
1541d519adbfSMatthew Knepley                MATOP_GET_VECS=88,
1542d519adbfSMatthew Knepley                MATOP_MAT_MULT=89,
1543d519adbfSMatthew Knepley                MATOP_MAT_MULT_SYMBOLIC=90,
1544d519adbfSMatthew Knepley                MATOP_MAT_MULT_NUMERIC=91,
1545d519adbfSMatthew Knepley                MATOP_PTAP=92,
1546d519adbfSMatthew Knepley                MATOP_PTAP_SYMBOLIC=93,
1547d519adbfSMatthew Knepley                MATOP_PTAP_NUMERIC=94,
1548d519adbfSMatthew Knepley                MATOP_MAT_MULTTRANSPOSE=95,
15491763ac4eSSatish Balay                MATOP_MAT_MULTTRANSPOSE_SYM=96,
15501763ac4eSSatish Balay                MATOP_MAT_MULTTRANSPOSE_NUM=97,
1551d519adbfSMatthew Knepley                MATOP_PTAP_SYMBOLIC_SEQAIJ=98,
1552d519adbfSMatthew Knepley                MATOP_PTAP_NUMERIC_SEQAIJ=99,
1553d519adbfSMatthew Knepley                MATOP_PTAP_SYMBOLIC_MPIAIJ=100,
1554d519adbfSMatthew Knepley                MATOP_PTAP_NUMERIC_MPIAIJ=101,
1555d519adbfSMatthew Knepley                MATOP_CONJUGATE=102,
1556d519adbfSMatthew Knepley                MATOP_SET_SIZES=103,
1557d519adbfSMatthew Knepley                MATOP_SET_VALUES_ROW=104,
1558d519adbfSMatthew Knepley                MATOP_REAL_PART=105,
1559d519adbfSMatthew Knepley                MATOP_IMAG_PART=106,
1560d519adbfSMatthew Knepley                MATOP_GET_ROW_UTRIANGULAR=107,
1561d519adbfSMatthew Knepley                MATOP_RESTORE_ROW_UTRIANGULAR=108,
1562d519adbfSMatthew Knepley                MATOP_MATSOLVE=109,
1563d519adbfSMatthew Knepley                MATOP_GET_REDUNDANTMATRIX=110,
1564d519adbfSMatthew Knepley                MATOP_GET_ROW_MIN=111,
1565d519adbfSMatthew Knepley                MATOP_GET_COLUMN_VEC=112,
1566d519adbfSMatthew Knepley                MATOP_MISSING_DIAGONAL=113,
1567d519adbfSMatthew Knepley                MATOP_MATGETSEQNONZEROSTRUCTURE=114,
156889c6957cSBarry Smith                MATOP_CREATE=115,
156989c6957cSBarry Smith                MATOP_GET_GHOSTS=116,
1570ea38565cSJed Brown                MATOP_GET_LOCALSUBMATRIX=117,
1571ea38565cSJed Brown                MATOP_RESTORE_LOCALSUBMATRIX=118,
1572eeffb40dSHong Zhang                MATOP_MULT_DIAGONAL_BLOCK=119,
1573547795f9SHong Zhang                MATOP_HERMITIANTRANSPOSE=120,
1574547795f9SHong Zhang                MATOP_MULTHERMITIANTRANSPOSE=121,
1575fbdbba38SShri Abhyankar                MATOP_MULTHERMITIANTRANSPOSEADD=122,
1576b9614d88SDmitry Karpeev                MATOP_GETMULTIPROCBLOCK=123,
15770716a85fSBarry Smith                MATOP_GETCOLUMNNORMS=125,
157837868618SMatthew G Knepley 	       MATOP_GET_SUBMATRICES_PARALLEL=128,
157937868618SMatthew G Knepley                MATOP_SET_VALUES_BATCH=129
1580fae171e0SBarry Smith              } MatOperation;
15817087cfbeSBarry Smith extern PetscErrorCode  MatHasOperation(Mat,MatOperation,PetscBool *);
15827087cfbeSBarry Smith extern PetscErrorCode  MatShellSetOperation(Mat,MatOperation,void(*)(void));
15837087cfbeSBarry Smith extern PetscErrorCode  MatShellGetOperation(Mat,MatOperation,void(**)(void));
15847087cfbeSBarry Smith extern PetscErrorCode  MatShellSetContext(Mat,void*);
1585112a2221SBarry Smith 
158690ace30eSBarry Smith /*
158790ace30eSBarry Smith    Codes for matrices stored on disk. By default they are
158890ace30eSBarry Smith    stored in a universal format. By changing the format with
15897973ad04SBarry Smith    PetscViewerSetFormat(viewer,PETSC_VIEWER_NATIVE); the matrices will
159090ace30eSBarry Smith    be stored in a way natural for the matrix, for example dense matrices
159190ace30eSBarry Smith    would be stored as dense. Matrices stored this way may only be
1592f4403165SShri Abhyankar    read into matrices of the same type.
159390ace30eSBarry Smith */
159490ace30eSBarry Smith #define MATRIX_BINARY_FORMAT_DENSE -1
159590ace30eSBarry Smith 
15967087cfbeSBarry Smith extern PetscErrorCode  MatMPIBAIJSetHashTableFactor(Mat,PetscReal);
15977087cfbeSBarry Smith extern PetscErrorCode  MatISGetLocalMat(Mat,Mat*);
15981f4e1ec7SBarry Smith 
1599d9274352SBarry Smith /*S
1600d9274352SBarry Smith      MatNullSpace - Object that removes a null space from a vector, i.e.
1601d9274352SBarry Smith          orthogonalizes the vector to a subsapce
1602d9274352SBarry Smith 
1603f7a9e4ceSBarry Smith    Level: advanced
1604d9274352SBarry Smith 
1605d9274352SBarry Smith   Concepts: matrix; linear operator, null space
1606d9274352SBarry Smith 
16076e1639daSBarry Smith   Users manual sections:
16086e1639daSBarry Smith .   sec_singular
16096e1639daSBarry Smith 
1610d9274352SBarry Smith .seealso:  MatNullSpaceCreate()
1611d9274352SBarry Smith S*/
161274637425SBarry Smith typedef struct _p_MatNullSpace* MatNullSpace;
1613d9274352SBarry Smith 
16147087cfbeSBarry Smith extern PetscErrorCode  MatNullSpaceCreate(MPI_Comm,PetscBool ,PetscInt,const Vec[],MatNullSpace*);
16157087cfbeSBarry Smith extern PetscErrorCode  MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(MatNullSpace,Vec,void*),void*);
1616d34fcf5fSBarry Smith extern PetscErrorCode  MatNullSpaceDestroy(MatNullSpace*);
16177087cfbeSBarry Smith extern PetscErrorCode  MatNullSpaceRemove(MatNullSpace,Vec,Vec*);
1618*1e633543SBarry Smith extern PetscErrorCode  MatSetNullSpace(Mat,MatNullSpace);
16197087cfbeSBarry Smith extern PetscErrorCode  MatNullSpaceTest(MatNullSpace,Mat,PetscBool  *);
1620b717e993SJed Brown extern PetscErrorCode  MatNullSpaceView(MatNullSpace,PetscViewer);
162174637425SBarry Smith 
16227087cfbeSBarry Smith extern PetscErrorCode  MatReorderingSeqSBAIJ(Mat,IS);
16237087cfbeSBarry Smith extern PetscErrorCode  MatMPISBAIJSetHashTableFactor(Mat,PetscReal);
16247087cfbeSBarry Smith extern PetscErrorCode  MatSeqSBAIJSetColumnIndices(Mat,PetscInt *);
16257087cfbeSBarry Smith extern PetscErrorCode  MatSeqBAIJInvertBlockDiagonal(Mat);
16263f1d51d7SBarry Smith 
16277087cfbeSBarry Smith extern PetscErrorCode  MatCreateMAIJ(Mat,PetscInt,Mat*);
16287087cfbeSBarry Smith extern PetscErrorCode  MatMAIJRedimension(Mat,PetscInt,Mat*);
16297087cfbeSBarry Smith extern PetscErrorCode  MatMAIJGetAIJ(Mat,Mat*);
1630c4f061fbSSatish Balay 
16317087cfbeSBarry Smith extern PetscErrorCode  MatComputeExplicitOperator(Mat,Mat*);
1632b0a32e0cSBarry Smith 
16337087cfbeSBarry Smith extern PetscErrorCode  MatDiagonalScaleLocal(Mat,Vec);
163404f1ad80SBarry Smith 
16357087cfbeSBarry Smith extern PetscErrorCode  MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*);
16367087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetBase(Mat,Vec,Vec);
16377087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*);
16387087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*));
16397087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec));
16407087cfbeSBarry Smith extern PetscErrorCode  MatMFFDAddNullSpace(Mat,MatNullSpace);
16417087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt);
16427087cfbeSBarry Smith extern PetscErrorCode  MatMFFDResetHHistory(Mat);
16437087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetFunctionError(Mat,PetscReal);
16447087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetPeriod(Mat,PetscInt);
16457087cfbeSBarry Smith extern PetscErrorCode  MatMFFDGetH(Mat,PetscScalar *);
16467087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetOptionsPrefix(Mat,const char[]);
16477087cfbeSBarry Smith extern PetscErrorCode  MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*);
16487087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*);
1649e884886fSBarry Smith 
16506370053bSBarry Smith /*S
16516370053bSBarry Smith     MatMFFD - A data structured used to manage the computation of the h differencing parameter for matrix-free
16526370053bSBarry Smith               Jacobian vector products
1653e884886fSBarry Smith 
16546370053bSBarry Smith     Notes: MATMFFD is a specific MatType which uses the MatMFFD data structure
16556370053bSBarry Smith 
16566370053bSBarry Smith            MatMFFD*() methods actually take the Mat as their first argument. Not a MatMFFD data structure
16576370053bSBarry Smith 
16586370053bSBarry Smith     Level: developer
16596370053bSBarry Smith 
16606370053bSBarry Smith .seealso: MATMFFD, MatCreateMFFD(), MatMFFDSetFuction(), MatMFFDSetType(), MatMFFDRegister()
16616370053bSBarry Smith S*/
1662e884886fSBarry Smith typedef struct _p_MatMFFD* MatMFFD;
1663e884886fSBarry Smith 
1664e884886fSBarry Smith /*E
1665e884886fSBarry Smith     MatMFFDType - algorithm used to compute the h used in computing matrix-vector products via differencing of the function
1666e884886fSBarry Smith 
1667e884886fSBarry Smith    Level: beginner
1668e884886fSBarry Smith 
1669e884886fSBarry Smith .seealso: MatMFFDSetType(), MatMFFDRegister()
1670e884886fSBarry Smith E*/
1671a313700dSBarry Smith #define MatMFFDType char*
1672e884886fSBarry Smith #define MATMFFD_DS  "ds"
1673e884886fSBarry Smith #define MATMFFD_WP  "wp"
1674e884886fSBarry Smith 
16757087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetType(Mat,const MatMFFDType);
16767087cfbeSBarry Smith extern PetscErrorCode  MatMFFDRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatMFFD));
1677e884886fSBarry Smith 
1678e884886fSBarry Smith /*MC
1679e884886fSBarry Smith    MatMFFDRegisterDynamic - Adds a method to the MatMFFD registry.
1680e884886fSBarry Smith 
1681e884886fSBarry Smith    Synopsis:
16821890ba74SBarry Smith    PetscErrorCode MatMFFDRegisterDynamic(const char *name_solver,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatMFFD))
1683e884886fSBarry Smith 
1684e884886fSBarry Smith    Not Collective
1685e884886fSBarry Smith 
1686e884886fSBarry Smith    Input Parameters:
1687e884886fSBarry Smith +  name_solver - name of a new user-defined compute-h module
1688e884886fSBarry Smith .  path - path (either absolute or relative) the library containing this solver
1689e884886fSBarry Smith .  name_create - name of routine to create method context
1690e884886fSBarry Smith -  routine_create - routine to create method context
1691e884886fSBarry Smith 
1692e884886fSBarry Smith    Level: developer
1693e884886fSBarry Smith 
1694e884886fSBarry Smith    Notes:
1695e884886fSBarry Smith    MatMFFDRegisterDynamic() may be called multiple times to add several user-defined solvers.
1696e884886fSBarry Smith 
1697e884886fSBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
1698e884886fSBarry Smith    is ignored.
1699e884886fSBarry Smith 
1700e884886fSBarry Smith    Sample usage:
1701e884886fSBarry Smith .vb
1702e884886fSBarry Smith    MatMFFDRegisterDynamic("my_h",/home/username/my_lib/lib/libO/solaris/mylib.a,
1703e884886fSBarry Smith                "MyHCreate",MyHCreate);
1704e884886fSBarry Smith .ve
1705e884886fSBarry Smith 
1706e884886fSBarry Smith    Then, your solver can be chosen with the procedural interface via
1707e884886fSBarry Smith $     MatMFFDSetType(mfctx,"my_h")
1708e884886fSBarry Smith    or at runtime via the option
1709e884886fSBarry Smith $     -snes_mf_type my_h
1710e884886fSBarry Smith 
1711e884886fSBarry Smith .keywords: MatMFFD, register
1712e884886fSBarry Smith 
1713e884886fSBarry Smith .seealso: MatMFFDRegisterAll(), MatMFFDRegisterDestroy()
1714e884886fSBarry Smith M*/
1715e884886fSBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1716e884886fSBarry Smith #define MatMFFDRegisterDynamic(a,b,c,d) MatMFFDRegister(a,b,c,0)
1717e884886fSBarry Smith #else
1718e884886fSBarry Smith #define MatMFFDRegisterDynamic(a,b,c,d) MatMFFDRegister(a,b,c,d)
1719e884886fSBarry Smith #endif
1720e884886fSBarry Smith 
17217087cfbeSBarry Smith extern PetscErrorCode  MatMFFDRegisterAll(const char[]);
17227087cfbeSBarry Smith extern PetscErrorCode  MatMFFDRegisterDestroy(void);
17237087cfbeSBarry Smith extern PetscErrorCode  MatMFFDDSSetUmin(Mat,PetscReal);
17247087cfbeSBarry Smith extern PetscErrorCode  MatMFFDWPSetComputeNormU(Mat,PetscBool );
1725e884886fSBarry Smith 
1726e884886fSBarry Smith 
17277087cfbeSBarry Smith extern PetscErrorCode  PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *);
17287087cfbeSBarry Smith extern PetscErrorCode  PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *);
17297dbadf16SMatthew Knepley 
173097969023SHong Zhang /*
173197969023SHong Zhang    PETSc interface to MUMPS
173297969023SHong Zhang */
173397969023SHong Zhang #ifdef PETSC_HAVE_MUMPS
1734f250808bSBarry Smith extern PetscErrorCode  MatMumpsSetIcntl(Mat,PetscInt,PetscInt);
173597969023SHong Zhang #endif
173623a5497aSJed Brown 
1737d954db57SHong Zhang /*
1738d954db57SHong Zhang    PETSc interface to SUPERLU
1739d954db57SHong Zhang */
1740d954db57SHong Zhang #ifdef PETSC_HAVE_SUPERLU
1741f250808bSBarry Smith extern PetscErrorCode  MatSuperluSetILUDropTol(Mat,PetscReal);
1742d954db57SHong Zhang #endif
1743d954db57SHong Zhang 
174490c53902SBarry Smith #if defined(PETSC_HAVE_CUSP)
174590c53902SBarry Smith extern PetscErrorCode  MatCreateSeqAIJCUSP(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
174690c53902SBarry Smith extern PetscErrorCode  MatCreateMPIAIJCUSP(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
174790c53902SBarry Smith #endif
174890c53902SBarry Smith 
174954efbe56SHong Zhang /*
175054efbe56SHong Zhang    PETSc interface to FFTW
175154efbe56SHong Zhang */
175254efbe56SHong Zhang #if defined(PETSC_HAVE_FFTW)
175374a26884SAmlan Barua extern PetscErrorCode VecScatterPetscToFFTW(Mat,Vec,Vec);
175474a26884SAmlan Barua extern PetscErrorCode VecScatterFFTWToPetsc(Mat,Vec,Vec);
17554be45526SHong Zhang extern PetscErrorCode MatGetVecsFFTW(Mat,Vec*,Vec*,Vec*);
175654efbe56SHong Zhang #endif
175754efbe56SHong Zhang 
1758c8883902SJed Brown extern PetscErrorCode MatCreateNest(MPI_Comm,PetscInt,const IS[],PetscInt,const IS[],const Mat[],Mat*);
1759c8883902SJed Brown extern PetscErrorCode MatNestGetSize(Mat,PetscInt*,PetscInt*);
1760c8883902SJed Brown extern PetscErrorCode MatNestGetSubMats(Mat,PetscInt*,PetscInt*,Mat***);
1761c8883902SJed Brown extern PetscErrorCode MatNestGetSubMat(Mat,PetscInt,PetscInt,Mat*);
17627087cfbeSBarry Smith extern PetscErrorCode MatNestSetVecType(Mat,const VecType);
1763c8883902SJed Brown extern PetscErrorCode MatNestSetSubMats(Mat,PetscInt,const IS[],PetscInt,const IS[],const Mat[]);
17640782ca92SJed Brown extern PetscErrorCode MatNestSetSubMat(Mat,PetscInt,PetscInt,Mat);
1765d8588912SDave May 
176623a5497aSJed Brown PETSC_EXTERN_CXX_END
176723a5497aSJed Brown #endif
1768