xref: /petsc/include/petscmat.h (revision 76bdecfb13ff6dcb09bbf73e6c929f447c427e71)
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 
20*76bdecfbSBarry Smith /*J
21d9274352SBarry Smith     MatType - String with the name of a PETSc matrix or the creation function
22d9274352SBarry Smith        with an optional dynamic library name, for example
23d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:mymatcreate()
24d9274352SBarry Smith 
25d9274352SBarry Smith    Level: beginner
26d9274352SBarry Smith 
27c7393fdbSBarry Smith .seealso: MatSetType(), Mat, MatSolverPackage
28*76bdecfbSBarry Smith J*/
29a313700dSBarry Smith #define MatType char*
30273d9f13SBarry Smith #define MATSAME            "same"
315a11e1b2SBarry Smith #define MATMAIJ            "maij"
32273d9f13SBarry Smith #define MATSEQMAIJ         "seqmaij"
33273d9f13SBarry Smith #define MATMPIMAIJ         "mpimaij"
34273d9f13SBarry Smith #define MATIS              "is"
355a11e1b2SBarry Smith #define MATAIJ             "aij"
36273d9f13SBarry Smith #define MATSEQAIJ          "seqaij"
377d6a0e61SBarry Smith #define MATSEQAIJPTHREAD   "seqaijpthread"
38bf2c1783SBarry Smith #define MATAIJPTHREAD      "aijpthread"
39273d9f13SBarry Smith #define MATMPIAIJ          "mpiaij"
405a11e1b2SBarry Smith #define MATAIJCRL          "aijcrl"
415a11e1b2SBarry Smith #define MATSEQAIJCRL       "seqaijcrl"
425a11e1b2SBarry Smith #define MATMPIAIJCRL       "mpiaijcrl"
438154be41SBarry Smith #define MATAIJCUSP         "aijcusp"
448154be41SBarry Smith #define MATSEQAIJCUSP      "seqaijcusp"
458154be41SBarry Smith #define MATMPIAIJCUSP      "mpiaijcusp"
465a11e1b2SBarry Smith #define MATAIJPERM         "aijperm"
475a11e1b2SBarry Smith #define MATSEQAIJPERM      "seqaijperm"
485a11e1b2SBarry Smith #define MATMPIAIJPERM      "mpiaijperm"
49273d9f13SBarry Smith #define MATSHELL           "shell"
505a11e1b2SBarry Smith #define MATDENSE           "dense"
51209238afSKris Buschelman #define MATSEQDENSE        "seqdense"
52273d9f13SBarry Smith #define MATMPIDENSE        "mpidense"
535a11e1b2SBarry Smith #define MATBAIJ            "baij"
54273d9f13SBarry Smith #define MATSEQBAIJ         "seqbaij"
55273d9f13SBarry Smith #define MATMPIBAIJ         "mpibaij"
56273d9f13SBarry Smith #define MATMPIADJ          "mpiadj"
575a11e1b2SBarry Smith #define MATSBAIJ           "sbaij"
58273d9f13SBarry Smith #define MATSEQSBAIJ        "seqsbaij"
59273d9f13SBarry Smith #define MATMPISBAIJ        "mpisbaij"
60c0cdd4a1SDahai Guo #define MATSEQBSTRM        "seqbstrm"
61c0cdd4a1SDahai Guo #define MATMPIBSTRM        "mpibstrm"
62c0cdd4a1SDahai Guo #define MATBSTRM           "bstrm"
63aa5a9175SDahai Guo #define MATSEQSBSTRM       "seqsbstrm"
64aa5a9175SDahai Guo #define MATMPISBSTRM       "mpisbstrm"
65aa5a9175SDahai Guo #define MATSBSTRM          "sbstrm"
66cebc7f6cSBarry Smith #define MATDAAD            "daad"
67cebc7f6cSBarry Smith #define MATMFFD            "mffd"
68c8a8475eSBarry Smith #define MATNORMAL          "normal"
69ab92ecdeSBarry Smith #define MATLRC             "lrc"
702a6744ebSBarry Smith #define MATSCATTER         "scatter"
71421e10b8SBarry Smith #define MATBLOCKMAT        "blockmat"
72793850ffSBarry Smith #define MATCOMPOSITE       "composite"
731f8c7532SHong Zhang #define MATFFT             "fft"
741f8c7532SHong Zhang #define MATFFTW            "fftw"
75e133240eSMatthew G Knepley #define MATSEQCUFFT        "seqcufft"
76557cca28SSatish Balay #define MATTRANSPOSEMAT    "transpose"
7772ca8751SBarry Smith #define MATSCHURCOMPLEMENT "schurcomplement"
781d6018f0SLisandro Dalcin #define MATPYTHON          "python"
79f91d8e95SBarry Smith #define MATHYPRESTRUCT     "hyprestruct"
80a9e6138eSGlenn Hammond #define MATHYPRESSTRUCT    "hypresstruct"
81ee1cef2cSJed Brown #define MATSUBMATRIX       "submatrix"
822c0dbf93SJed Brown #define MATLOCALREF        "localref"
83d8588912SDave May #define MATNEST            "nest"
84773941d6SBarry Smith 
85*76bdecfbSBarry Smith /*J
86c7393fdbSBarry Smith     MatSolverPackage - String with the name of a PETSc matrix solver type.
879989ab13SBarry Smith 
889989ab13SBarry Smith     For example: "petsc" indicates what PETSc provides, "superlu" indicates either
899989ab13SBarry Smith        SuperLU or SuperLU_Dist etc.
909989ab13SBarry Smith 
919989ab13SBarry Smith 
929989ab13SBarry Smith    Level: beginner
939989ab13SBarry Smith 
945c9eb25fSBarry Smith .seealso: MatGetFactor(), Mat, MatSetType(), MatType
95*76bdecfbSBarry Smith J*/
96c7393fdbSBarry Smith #define MatSolverPackage char*
972692d6eeSBarry Smith #define MATSOLVERSPOOLES      "spooles"
982692d6eeSBarry Smith #define MATSOLVERSUPERLU      "superlu"
992692d6eeSBarry Smith #define MATSOLVERSUPERLU_DIST "superlu_dist"
1002692d6eeSBarry Smith #define MATSOLVERUMFPACK      "umfpack"
10120db9a53SJed Brown #define MATSOLVERCHOLMOD      "cholmod"
1022692d6eeSBarry Smith #define MATSOLVERESSL         "essl"
1032692d6eeSBarry Smith #define MATSOLVERLUSOL        "lusol"
1042692d6eeSBarry Smith #define MATSOLVERMUMPS        "mumps"
1052692d6eeSBarry Smith #define MATSOLVERPASTIX       "pastix"
1062692d6eeSBarry Smith #define MATSOLVERMATLAB       "matlab"
1072692d6eeSBarry Smith #define MATSOLVERPETSC        "petsc"
1082692d6eeSBarry Smith #define MATSOLVERPLAPACK      "plapack"
1092692d6eeSBarry Smith #define MATSOLVERBAS          "bas"
110773941d6SBarry Smith 
111aa5a9175SDahai Guo #define MATSOLVERBSTRM        "bstrm"
112aa5a9175SDahai Guo #define MATSOLVERSBSTRM       "sbstrm"
113c0cdd4a1SDahai Guo 
114b24902e0SBarry Smith /*E
115b24902e0SBarry Smith     MatFactorType - indicates what type of factorization is requested
116b24902e0SBarry Smith 
117b24902e0SBarry Smith     Level: beginner
118b24902e0SBarry Smith 
119b24902e0SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
120b24902e0SBarry Smith 
121c7393fdbSBarry Smith .seealso: MatSolverPackage, MatGetFactor()
122b24902e0SBarry Smith E*/
123599ef60dSHong Zhang typedef enum {MAT_FACTOR_NONE, MAT_FACTOR_LU, MAT_FACTOR_CHOLESKY, MAT_FACTOR_ILU, MAT_FACTOR_ICC,MAT_FACTOR_ILUDT} MatFactorType;
12434918c53SJed Brown extern const char *const MatFactorTypes[];
125e92e720dSBarry Smith 
1267087cfbeSBarry Smith extern PetscErrorCode  MatGetFactor(Mat,const MatSolverPackage,MatFactorType,Mat*);
1277087cfbeSBarry Smith extern PetscErrorCode  MatGetFactorAvailable(Mat,const MatSolverPackage,MatFactorType,PetscBool *);
1287087cfbeSBarry Smith extern PetscErrorCode  MatFactorGetSolverPackage(Mat,const MatSolverPackage*);
1297087cfbeSBarry Smith extern PetscErrorCode  MatGetFactorType(Mat,MatFactorType*);
1309989ab13SBarry Smith 
131c06d978dSMatthew Knepley /* Logging support */
1320700a824SBarry Smith #define    MAT_FILE_CLASSID 1211216    /* used to indicate matrices in binary files */
1337087cfbeSBarry Smith extern PetscClassId  MAT_CLASSID;
1347087cfbeSBarry Smith extern PetscClassId  MAT_FDCOLORING_CLASSID;
1357087cfbeSBarry Smith extern PetscClassId  MAT_PARTITIONING_CLASSID;
1367087cfbeSBarry Smith extern PetscClassId  MAT_NULLSPACE_CLASSID;
1377087cfbeSBarry Smith extern PetscClassId  MATMFFD_CLASSID;
138c06d978dSMatthew Knepley 
139ceb03754SKris Buschelman /*E
140ceb03754SKris Buschelman     MatReuse - Indicates if matrices obtained from a previous call to MatGetSubMatrices()
141d6eff37eSBarry Smith      or MatGetSubMatrix() are to be reused to store the new matrix values. For MatConvert() is used to indicate
142d6eff37eSBarry Smith      that the input matrix is to be replaced with the converted matrix.
143ceb03754SKris Buschelman 
144ceb03754SKris Buschelman     Level: beginner
145ceb03754SKris Buschelman 
146ceb03754SKris Buschelman    Any additions/changes here MUST also be made in include/finclude/petscmat.h
147ceb03754SKris Buschelman 
1480c451bc4SBarry Smith .seealso: MatGetSubMatrices(), MatGetSubMatrix(), MatDestroyMatrices(), MatConvert()
149ceb03754SKris Buschelman E*/
150dfe085dbSJed Brown typedef enum {MAT_INITIAL_MATRIX,MAT_REUSE_MATRIX,MAT_IGNORE_MATRIX} MatReuse;
151ceb03754SKris Buschelman 
1525494a064SHong Zhang /*E
1535494a064SHong Zhang     MatGetSubMatrixOption - Indicates if matrices obtained from a call to MatGetSubMatrices()
154829201f2SHong Zhang      include the matrix values. Currently it is only used by MatGetSeqNonzerostructure().
1555494a064SHong Zhang 
1565494a064SHong Zhang     Level: beginner
1575494a064SHong Zhang 
158829201f2SHong Zhang .seealso: MatGetSeqNonzerostructure()
1595494a064SHong Zhang E*/
1605494a064SHong Zhang typedef enum {MAT_DO_NOT_GET_VALUES,MAT_GET_VALUES} MatGetSubMatrixOption;
1615494a064SHong Zhang 
1627087cfbeSBarry Smith extern PetscErrorCode  MatInitializePackage(const char[]);
163c06d978dSMatthew Knepley 
1647087cfbeSBarry Smith extern PetscErrorCode  MatCreate(MPI_Comm,Mat*);
165a69a960dSMatthew Knepley PetscPolymorphicFunction(MatCreate,(MPI_Comm comm),(comm,&A),Mat,A)
166a69a960dSMatthew Knepley PetscPolymorphicFunction(MatCreate,(),(PETSC_COMM_WORLD,&A),Mat,A)
1677087cfbeSBarry Smith extern PetscErrorCode  MatSetSizes(Mat,PetscInt,PetscInt,PetscInt,PetscInt);
1687087cfbeSBarry Smith extern PetscErrorCode  MatSetType(Mat,const MatType);
1697087cfbeSBarry Smith extern PetscErrorCode  MatSetFromOptions(Mat);
1707087cfbeSBarry Smith extern PetscErrorCode  MatSetUpPreallocation(Mat);
1717087cfbeSBarry Smith extern PetscErrorCode  MatRegisterAll(const char[]);
1727087cfbeSBarry Smith extern PetscErrorCode  MatRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat));
17301bebe75SBarry Smith extern PetscErrorCode  MatRegisterBaseName(const char[],const char[],const char[]);
1747087cfbeSBarry Smith extern PetscErrorCode  MatSetOptionsPrefix(Mat,const char[]);
1757087cfbeSBarry Smith extern PetscErrorCode  MatAppendOptionsPrefix(Mat,const char[]);
1767087cfbeSBarry Smith extern PetscErrorCode  MatGetOptionsPrefix(Mat,const char*[]);
177f69a0ea3SMatthew Knepley 
17830de9b25SBarry Smith /*MC
17930de9b25SBarry Smith    MatRegisterDynamic - Adds a new matrix type
18030de9b25SBarry Smith 
18130de9b25SBarry Smith    Synopsis:
1821890ba74SBarry Smith    PetscErrorCode MatRegisterDynamic(const char *name,const char *path,const char *name_create,PetscErrorCode (*routine_create)(Mat))
18330de9b25SBarry Smith 
18430de9b25SBarry Smith    Not Collective
18530de9b25SBarry Smith 
18630de9b25SBarry Smith    Input Parameters:
18730de9b25SBarry Smith +  name - name of a new user-defined matrix type
18830de9b25SBarry Smith .  path - path (either absolute or relative) the library containing this solver
18930de9b25SBarry Smith .  name_create - name of routine to create method context
19030de9b25SBarry Smith -  routine_create - routine to create method context
19130de9b25SBarry Smith 
19230de9b25SBarry Smith    Notes:
19330de9b25SBarry Smith    MatRegisterDynamic() may be called multiple times to add several user-defined solvers.
19430de9b25SBarry Smith 
19530de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
19630de9b25SBarry Smith    is ignored.
19730de9b25SBarry Smith 
19830de9b25SBarry Smith    Sample usage:
19930de9b25SBarry Smith .vb
20030de9b25SBarry Smith    MatRegisterDynamic("my_mat",/home/username/my_lib/lib/libO/solaris/mylib.a,
20130de9b25SBarry Smith                "MyMatCreate",MyMatCreate);
20230de9b25SBarry Smith .ve
20330de9b25SBarry Smith 
20430de9b25SBarry Smith    Then, your solver can be chosen with the procedural interface via
20530de9b25SBarry Smith $     MatSetType(Mat,"my_mat")
20630de9b25SBarry Smith    or at runtime via the option
20730de9b25SBarry Smith $     -mat_type my_mat
20830de9b25SBarry Smith 
20930de9b25SBarry Smith    Level: advanced
21030de9b25SBarry Smith 
211ab901514SBarry Smith    Notes: ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values.
21230de9b25SBarry Smith          If your function is not being put into a shared library then use VecRegister() instead
21330de9b25SBarry Smith 
21430de9b25SBarry Smith .keywords: Mat, register
21530de9b25SBarry Smith 
21630de9b25SBarry Smith .seealso: MatRegisterAll(), MatRegisterDestroy()
21730de9b25SBarry Smith 
21830de9b25SBarry Smith M*/
219273d9f13SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
220273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,0)
221273d9f13SBarry Smith #else
222273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,d)
22330de9b25SBarry Smith #endif
22430de9b25SBarry Smith 
225ace3abfcSBarry Smith extern PetscBool  MatRegisterAllCalled;
226b0a32e0cSBarry Smith extern PetscFList MatList;
227b022a5c1SBarry Smith extern PetscFList MatColoringList;
228b022a5c1SBarry Smith extern PetscFList MatPartitioningList;
22928988994SBarry Smith 
2303b224e63SBarry Smith /*E
2313b224e63SBarry Smith     MatStructure - Indicates if the matrix has the same nonzero structure
2323b224e63SBarry Smith 
2333b224e63SBarry Smith     Level: beginner
2343b224e63SBarry Smith 
2353b224e63SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
2363b224e63SBarry Smith 
2373b224e63SBarry Smith .seealso: MatCopy(), KSPSetOperators(), PCSetOperators()
2383b224e63SBarry Smith E*/
239214bc6a2SJed Brown typedef enum {DIFFERENT_NONZERO_PATTERN,SUBSET_NONZERO_PATTERN,SAME_NONZERO_PATTERN,SAME_PRECONDITIONER} MatStructure;
2403b224e63SBarry Smith 
2417087cfbeSBarry Smith extern PetscErrorCode  MatCreateSeqDense(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*);
2427087cfbeSBarry Smith extern PetscErrorCode  MatCreateMPIDense(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*);
2437087cfbeSBarry Smith extern PetscErrorCode  MatCreateSeqAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
244ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,nz,nnz,&A),Mat,A)
245ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,&A),Mat,A)
246ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,0,nnz,&A),Mat,A)
247ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,&A),Mat,A)
248ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,A))
249ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,m,n,0,nnz,A))
250ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,A))
2517087cfbeSBarry Smith extern PetscErrorCode  MatCreateMPIAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
252ba966639SSatish 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)
253ba966639SSatish 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)
254ba966639SSatish 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)
255ba966639SSatish 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)
256ba966639SSatish 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))
257ba966639SSatish 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))
258ba966639SSatish 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))
259ba966639SSatish 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)
260ba966639SSatish 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)
261ba966639SSatish 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)
262ba966639SSatish 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)
263ba966639SSatish 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))
264ba966639SSatish 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))
265ba966639SSatish 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))
2667087cfbeSBarry Smith extern PetscErrorCode  MatCreateMPIAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *);
2677087cfbeSBarry Smith extern PetscErrorCode  MatCreateMPIAIJWithSplitArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],PetscInt[],PetscInt[],PetscScalar[],Mat*);
2688d7a6e47SBarry Smith 
2697087cfbeSBarry Smith extern PetscErrorCode  MatCreateSeqBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
270ba966639SSatish 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)
271ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A)
272ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A)
273ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A)
274ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A))
275ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A))
276ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A))
2777087cfbeSBarry Smith extern PetscErrorCode  MatCreateMPIBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
278ba966639SSatish 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)
279ba966639SSatish 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)
280ba966639SSatish 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)
281ba966639SSatish 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)
282ba966639SSatish 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))
283ba966639SSatish 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))
284ba966639SSatish 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))
285ba966639SSatish 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)
286ba966639SSatish 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)
287ba966639SSatish 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)
288ba966639SSatish 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)
289ba966639SSatish 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))
290ba966639SSatish 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))
291ba966639SSatish 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))
2927087cfbeSBarry Smith extern PetscErrorCode  MatCreateMPIBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat*);
293d21a29f3SJed Brown 
2947087cfbeSBarry Smith extern PetscErrorCode  MatCreateMPIAdj(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscInt[],Mat*);
2957087cfbeSBarry Smith extern PetscErrorCode  MatCreateSeqSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
296ba966639SSatish 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)
297ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A)
298ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A)
299ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A)
300ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A))
301ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A))
302ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A))
303d21a29f3SJed Brown 
3047087cfbeSBarry Smith extern PetscErrorCode  MatCreateMPISBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
305ba966639SSatish 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)
306ba966639SSatish 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)
307ba966639SSatish 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)
308ba966639SSatish 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)
309ba966639SSatish 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))
310ba966639SSatish 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))
311ba966639SSatish 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))
312ba966639SSatish 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)
313ba966639SSatish 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)
314ba966639SSatish 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)
315ba966639SSatish 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)
316ba966639SSatish 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))
317ba966639SSatish 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))
318ba966639SSatish 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))
3197087cfbeSBarry Smith extern PetscErrorCode  MatCreateMPISBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *);
3207087cfbeSBarry Smith extern PetscErrorCode  MatMPISBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
321dfb205c3SBarry Smith 
3227087cfbeSBarry Smith extern PetscErrorCode  MatCreateShell(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,void *,Mat*);
323ba966639SSatish 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)
324ba966639SSatish Balay PetscPolymorphicFunction(MatCreateShell,(PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx),(PETSC_COMM_WORLD,m,n,M,N,ctx,&A),Mat,A)
3257087cfbeSBarry Smith extern PetscErrorCode  MatCreateAdic(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,void (*)(void),Mat*);
3267087cfbeSBarry Smith extern PetscErrorCode  MatCreateNormal(Mat,Mat*);
327ba966639SSatish Balay PetscPolymorphicFunction(MatCreateNormal,(Mat mat),(mat,&A),Mat,A)
3287087cfbeSBarry Smith extern PetscErrorCode  MatCreateLRC(Mat,Mat,Mat,Mat*);
3297087cfbeSBarry Smith extern PetscErrorCode  MatCreateIS(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,ISLocalToGlobalMapping,Mat*);
3307087cfbeSBarry Smith extern PetscErrorCode  MatCreateSeqAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
3317087cfbeSBarry Smith extern PetscErrorCode  MatCreateMPIAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
332c0cdd4a1SDahai Guo 
333c0cdd4a1SDahai Guo extern PetscErrorCode  MatCreateSeqBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
334c0cdd4a1SDahai Guo extern PetscErrorCode  MatCreateMPIBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
335aa5a9175SDahai Guo extern PetscErrorCode  MatCreateSeqSBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
336aa5a9175SDahai Guo extern PetscErrorCode  MatCreateMPISBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
337c0cdd4a1SDahai Guo 
3387087cfbeSBarry Smith extern PetscErrorCode  MatCreateScatter(MPI_Comm,VecScatter,Mat*);
3397087cfbeSBarry Smith extern PetscErrorCode  MatScatterSetVecScatter(Mat,VecScatter);
3407087cfbeSBarry Smith extern PetscErrorCode  MatScatterGetVecScatter(Mat,VecScatter*);
3417087cfbeSBarry Smith extern PetscErrorCode  MatCreateBlockMat(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,Mat*);
3427087cfbeSBarry Smith extern PetscErrorCode  MatCompositeAddMat(Mat,Mat);
3437087cfbeSBarry Smith extern PetscErrorCode  MatCompositeMerge(Mat);
3447087cfbeSBarry Smith extern PetscErrorCode  MatCreateComposite(MPI_Comm,PetscInt,const Mat*,Mat*);
3456d7c1e57SBarry Smith typedef enum {MAT_COMPOSITE_ADDITIVE,MAT_COMPOSITE_MULTIPLICATIVE} MatCompositeType;
3467087cfbeSBarry Smith extern PetscErrorCode  MatCompositeSetType(Mat,MatCompositeType);
3476d7c1e57SBarry Smith 
348dedccee8SHong Zhang extern PetscErrorCode  MatCreateFFT(MPI_Comm,PetscInt,const PetscInt[],const MatType,Mat*);
3497087cfbeSBarry Smith extern PetscErrorCode  MatCreateSeqCUFFT(MPI_Comm,PetscInt,const PetscInt[],Mat*);
350dedccee8SHong Zhang 
3517087cfbeSBarry Smith extern PetscErrorCode  MatCreateTranspose(Mat,Mat*);
3527087cfbeSBarry Smith extern PetscErrorCode  MatCreateSubMatrix(Mat,IS,IS,Mat*);
3537087cfbeSBarry Smith extern PetscErrorCode  MatSubMatrixUpdate(Mat,Mat,IS,IS);
3547087cfbeSBarry Smith extern PetscErrorCode  MatCreateLocalRef(Mat,IS,IS,Mat*);
3551472f72bSBarry Smith 
3567087cfbeSBarry Smith extern PetscErrorCode  MatPythonSetType(Mat,const char[]);
3571d6018f0SLisandro Dalcin 
3587087cfbeSBarry Smith extern PetscErrorCode  MatSetUp(Mat);
359fcfd50ebSBarry Smith extern PetscErrorCode  MatDestroy(Mat*);
36021c89e3eSBarry Smith 
3617087cfbeSBarry Smith extern PetscErrorCode  MatConjugate(Mat);
3627087cfbeSBarry Smith extern PetscErrorCode  MatRealPart(Mat);
3637087cfbeSBarry Smith extern PetscErrorCode  MatImaginaryPart(Mat);
36411bd1e4dSLisandro Dalcin extern PetscErrorCode  MatGetDiagonalBlock(Mat,Mat*);
3657087cfbeSBarry Smith extern PetscErrorCode  MatGetTrace(Mat,PetscScalar*);
366bbead8a2SBarry Smith extern PetscErrorCode  MatInvertBlockDiagonal(Mat,PetscScalar **);
36799cafbc1SBarry Smith 
3688ed539a5SBarry Smith /* ------------------------------------------------------------*/
3697087cfbeSBarry Smith extern PetscErrorCode  MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
3707087cfbeSBarry Smith extern PetscErrorCode  MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
3717087cfbeSBarry Smith extern PetscErrorCode  MatSetValuesRow(Mat,PetscInt,const PetscScalar[]);
3727087cfbeSBarry Smith extern PetscErrorCode  MatSetValuesRowLocal(Mat,PetscInt,const PetscScalar[]);
3732f35a324SMatthew G Knepley extern PetscErrorCode  MatSetValuesBatch(Mat,PetscInt,PetscInt,PetscInt[],const PetscScalar[]);
37484cb2905SBarry Smith 
3752ef4de8bSBarry Smith /*S
3762ef4de8bSBarry Smith      MatStencil - Data structure (C struct) for storing information about a single row or
3772ef4de8bSBarry Smith         column of a matrix as index on an associated grid.
3782ef4de8bSBarry Smith 
3792ef4de8bSBarry Smith    Level: beginner
3802ef4de8bSBarry Smith 
3812ef4de8bSBarry Smith   Concepts: matrix; linear operator
3822ef4de8bSBarry Smith 
383d7bd93baSBarry Smith .seealso:  MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockStencil()
3842ef4de8bSBarry Smith S*/
385435da068SBarry Smith typedef struct {
386c1ac3661SBarry Smith   PetscInt k,j,i,c;
387435da068SBarry Smith } MatStencil;
3882ef4de8bSBarry Smith 
3897087cfbeSBarry Smith extern PetscErrorCode  MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
3907087cfbeSBarry Smith extern PetscErrorCode  MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
3917087cfbeSBarry Smith extern PetscErrorCode  MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt);
392435da068SBarry Smith 
3937087cfbeSBarry Smith extern PetscErrorCode  MatSetColoring(Mat,ISColoring);
3947087cfbeSBarry Smith extern PetscErrorCode  MatSetValuesAdic(Mat,void*);
3957087cfbeSBarry Smith extern PetscErrorCode  MatSetValuesAdifor(Mat,PetscInt,void*);
3963a7fca6bSBarry Smith 
397d91e6319SBarry Smith /*E
398d91e6319SBarry Smith     MatAssemblyType - Indicates if the matrix is now to be used, or if you plan
399d91e6319SBarry Smith      to continue to add values to it
400d91e6319SBarry Smith 
401d91e6319SBarry Smith     Level: beginner
402d91e6319SBarry Smith 
403d91e6319SBarry Smith .seealso: MatAssemblyBegin(), MatAssemblyEnd()
404d91e6319SBarry Smith E*/
4056d4a8577SBarry Smith typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType;
4067087cfbeSBarry Smith extern PetscErrorCode  MatAssemblyBegin(Mat,MatAssemblyType);
4077087cfbeSBarry Smith extern PetscErrorCode  MatAssemblyEnd(Mat,MatAssemblyType);
4087087cfbeSBarry Smith extern PetscErrorCode  MatAssembled(Mat,PetscBool *);
4094f9c727eSBarry Smith 
4101d73ed98SBarry Smith 
41130de9b25SBarry Smith 
412d91e6319SBarry Smith /*E
413d91e6319SBarry Smith     MatOption - Options that may be set for a matrix and its behavior or storage
414d91e6319SBarry Smith 
415d91e6319SBarry Smith     Level: beginner
416d91e6319SBarry Smith 
4170a835dfdSSatish Balay    Any additions/changes here MUST also be made in include/finclude/petscmat.h
418d91e6319SBarry Smith 
419d91e6319SBarry Smith .seealso: MatSetOption()
420d91e6319SBarry Smith E*/
421512a5fc5SBarry Smith typedef enum {MAT_ROW_ORIENTED,MAT_NEW_NONZERO_LOCATIONS,
4224e0d8c25SBarry Smith               MAT_SYMMETRIC,
4234e0d8c25SBarry Smith               MAT_STRUCTURALLY_SYMMETRIC,
4248047e77bSBarry Smith               MAT_NEW_DIAGONALS,MAT_IGNORE_OFF_PROC_ENTRIES,
4254e0d8c25SBarry Smith               MAT_NEW_NONZERO_LOCATION_ERR,
4264e0d8c25SBarry Smith               MAT_NEW_NONZERO_ALLOCATION_ERR,MAT_USE_HASH_TABLE,
427a9817697SBarry Smith               MAT_KEEP_NONZERO_PATTERN,MAT_IGNORE_ZERO_ENTRIES,
4284e0d8c25SBarry Smith               MAT_USE_INODES,
4294e0d8c25SBarry Smith               MAT_HERMITIAN,
4304e0d8c25SBarry Smith               MAT_SYMMETRY_ETERNAL,
431cd6b891eSBarry Smith               MAT_CHECK_COMPRESSED_ROW,
4324e0d8c25SBarry Smith               MAT_IGNORE_LOWER_TRIANGULAR,MAT_ERROR_LOWER_TRIANGULAR,
43328b2fa4aSMatthew Knepley               MAT_GETROW_UPPERTRIANGULAR,MAT_UNUSED_NONZERO_LOCATION_ERR,
4344cb17eb5SBarry Smith               MAT_SPD,MAT_NO_OFF_PROC_ENTRIES,MAT_NO_OFF_PROC_ZERO_ROWS,
43528b2fa4aSMatthew Knepley               NUM_MAT_OPTIONS} MatOption;
436290bbb0aSBarry Smith extern const char *MatOptions[];
4377087cfbeSBarry Smith extern PetscErrorCode  MatSetOption(Mat,MatOption,PetscBool );
4387087cfbeSBarry Smith extern PetscErrorCode  MatGetType(Mat,const MatType*);
439a313700dSBarry Smith PetscPolymorphicFunction(MatGetType,(Mat mat),(mat,&t),const MatType,t)
44084cb2905SBarry Smith 
4417087cfbeSBarry Smith extern PetscErrorCode  MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]);
4427087cfbeSBarry Smith extern PetscErrorCode  MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
4437087cfbeSBarry Smith extern PetscErrorCode  MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
4447087cfbeSBarry Smith extern PetscErrorCode  MatGetRowUpperTriangular(Mat);
4457087cfbeSBarry Smith extern PetscErrorCode  MatRestoreRowUpperTriangular(Mat);
4467087cfbeSBarry Smith extern PetscErrorCode  MatGetColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
4477087cfbeSBarry Smith extern PetscErrorCode  MatRestoreColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
4487087cfbeSBarry Smith extern PetscErrorCode  MatGetColumnVector(Mat,Vec,PetscInt);
4497087cfbeSBarry Smith extern PetscErrorCode  MatGetArray(Mat,PetscScalar *[]);
450ba966639SSatish Balay PetscPolymorphicFunction(MatGetArray,(Mat mat),(mat,&a),PetscScalar*,a)
4517087cfbeSBarry Smith extern PetscErrorCode  MatRestoreArray(Mat,PetscScalar *[]);
4527087cfbeSBarry Smith extern PetscErrorCode  MatGetBlockSize(Mat,PetscInt *);
453ba966639SSatish Balay PetscPolymorphicFunction(MatGetBlockSize,(Mat mat),(mat,&a),PetscInt,a)
4547087cfbeSBarry Smith extern PetscErrorCode  MatSetBlockSize(Mat,PetscInt);
4557b80b807SBarry Smith 
4561620fd73SBarry Smith 
4577087cfbeSBarry Smith extern PetscErrorCode  MatMult(Mat,Vec,Vec);
4587087cfbeSBarry Smith extern PetscErrorCode  MatMultDiagonalBlock(Mat,Vec,Vec);
4597087cfbeSBarry Smith extern PetscErrorCode  MatMultAdd(Mat,Vec,Vec,Vec);
460ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
4617087cfbeSBarry Smith extern PetscErrorCode  MatMultTranspose(Mat,Vec,Vec);
4627087cfbeSBarry Smith extern PetscErrorCode  MatMultHermitianTranspose(Mat,Vec,Vec);
4637087cfbeSBarry Smith extern PetscErrorCode  MatIsTranspose(Mat,Mat,PetscReal,PetscBool *);
464ace3abfcSBarry Smith PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B,PetscReal tol),(A,B,tol,&t),PetscBool ,t)
465ace3abfcSBarry Smith PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B),(A,B,0,&t),PetscBool ,t)
4667087cfbeSBarry Smith extern PetscErrorCode  MatIsHermitianTranspose(Mat,Mat,PetscReal,PetscBool *);
4677087cfbeSBarry Smith extern PetscErrorCode  MatMultTransposeAdd(Mat,Vec,Vec,Vec);
4687087cfbeSBarry Smith extern PetscErrorCode  MatMultHermitianTransposeAdd(Mat,Vec,Vec,Vec);
469ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultTransposeAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
4707087cfbeSBarry Smith extern PetscErrorCode  MatMultConstrained(Mat,Vec,Vec);
4717087cfbeSBarry Smith extern PetscErrorCode  MatMultTransposeConstrained(Mat,Vec,Vec);
4727087cfbeSBarry Smith extern PetscErrorCode  MatMatSolve(Mat,Mat,Mat);
473f5edf698SHong Zhang 
474d91e6319SBarry Smith /*E
475d91e6319SBarry Smith     MatDuplicateOption - Indicates if a duplicated sparse matrix should have
476d91e6319SBarry Smith   its numerical values copied over or just its nonzero structure.
477d91e6319SBarry Smith 
478d91e6319SBarry Smith     Level: beginner
479d91e6319SBarry Smith 
480d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
481d91e6319SBarry Smith 
48270dcbbb9SBarry Smith $   MAT_SHARE_NONZERO_PATTERN - the i and j arrays in the new matrix will be shared with the original matrix
48370dcbbb9SBarry Smith $                               this also triggers the MAT_DO_NOT_COPY_VALUES option. This is used when you
48470dcbbb9SBarry Smith $                               have several matrices with the same nonzero pattern.
48570dcbbb9SBarry Smith 
486d91e6319SBarry Smith .seealso: MatDuplicate()
487d91e6319SBarry Smith E*/
48870dcbbb9SBarry Smith typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES,MAT_SHARE_NONZERO_PATTERN} MatDuplicateOption;
4892e8a6d31SBarry Smith 
4907087cfbeSBarry Smith extern PetscErrorCode  MatConvert(Mat,const MatType,MatReuse,Mat*);
491a313700dSBarry Smith PetscPolymorphicFunction(MatConvert,(Mat A,const MatType t),(A,t,MAT_INITIAL_MATRIX,&a),Mat,a)
4927087cfbeSBarry Smith extern PetscErrorCode  MatDuplicate(Mat,MatDuplicateOption,Mat*);
493ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A,MatDuplicateOption o),(A,o,&a),Mat,a)
494ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A),(A,MAT_COPY_VALUES,&a),Mat,a)
49594a9d846SBarry Smith 
496cb5b572fSBarry Smith 
4977087cfbeSBarry Smith extern PetscErrorCode  MatCopy(Mat,Mat,MatStructure);
4987087cfbeSBarry Smith extern PetscErrorCode  MatView(Mat,PetscViewer);
4997087cfbeSBarry Smith extern PetscErrorCode  MatIsSymmetric(Mat,PetscReal,PetscBool *);
500ace3abfcSBarry Smith PetscPolymorphicFunction(MatIsSymmetric,(Mat A,PetscReal tol),(A,tol,&t),PetscBool ,t)
501ace3abfcSBarry Smith PetscPolymorphicFunction(MatIsSymmetric,(Mat A),(A,0,&t),PetscBool ,t)
5027087cfbeSBarry Smith extern PetscErrorCode  MatIsStructurallySymmetric(Mat,PetscBool *);
503ace3abfcSBarry Smith PetscPolymorphicFunction(MatIsStructurallySymmetric,(Mat A),(A,&t),PetscBool ,t)
5047087cfbeSBarry Smith extern PetscErrorCode  MatIsHermitian(Mat,PetscReal,PetscBool *);
505ace3abfcSBarry Smith PetscPolymorphicFunction(MatIsHermitian,(Mat A),(A,0,&t),PetscBool ,t)
5067087cfbeSBarry Smith extern PetscErrorCode  MatIsSymmetricKnown(Mat,PetscBool *,PetscBool *);
5077087cfbeSBarry Smith extern PetscErrorCode  MatIsHermitianKnown(Mat,PetscBool *,PetscBool *);
5087087cfbeSBarry Smith extern PetscErrorCode  MatMissingDiagonal(Mat,PetscBool  *,PetscInt *);
5097087cfbeSBarry Smith extern PetscErrorCode  MatLoad(Mat, PetscViewer);
5107b80b807SBarry Smith 
5117087cfbeSBarry Smith extern PetscErrorCode  MatGetRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,PetscInt *[],PetscInt *[],PetscBool  *);
5127087cfbeSBarry Smith extern PetscErrorCode  MatRestoreRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,PetscInt *[],PetscInt *[],PetscBool  *);
5137087cfbeSBarry Smith extern PetscErrorCode  MatGetColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,PetscInt *[],PetscInt *[],PetscBool  *);
5147087cfbeSBarry Smith extern PetscErrorCode  MatRestoreColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,PetscInt *[],PetscInt *[],PetscBool  *);
515d4fbbf0eSBarry Smith 
516d91e6319SBarry Smith /*S
517d91e6319SBarry Smith      MatInfo - Context of matrix information, used with MatGetInfo()
518d91e6319SBarry Smith 
519d91e6319SBarry Smith    In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE
520d91e6319SBarry Smith 
521d91e6319SBarry Smith    Level: intermediate
522d91e6319SBarry Smith 
523d91e6319SBarry Smith   Concepts: matrix^nonzero information
524d91e6319SBarry Smith 
525d9274352SBarry Smith .seealso:  MatGetInfo(), MatInfoType
526d91e6319SBarry Smith S*/
5274e220ebcSLois Curfman McInnes typedef struct {
528b0a32e0cSBarry Smith   PetscLogDouble block_size;                         /* block size */
529b0a32e0cSBarry Smith   PetscLogDouble nz_allocated,nz_used,nz_unneeded;   /* number of nonzeros */
530b0a32e0cSBarry Smith   PetscLogDouble memory;                             /* memory allocated */
531b0a32e0cSBarry Smith   PetscLogDouble assemblies;                         /* number of matrix assemblies called */
532b0a32e0cSBarry Smith   PetscLogDouble mallocs;                            /* number of mallocs during MatSetValues() */
533b0a32e0cSBarry Smith   PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */
534b0a32e0cSBarry Smith   PetscLogDouble factor_mallocs;                     /* number of mallocs during factorization */
5354e220ebcSLois Curfman McInnes } MatInfo;
5364e220ebcSLois Curfman McInnes 
537d9274352SBarry Smith /*E
538d9274352SBarry Smith     MatInfoType - Indicates if you want information about the local part of the matrix,
539d9274352SBarry Smith      the entire parallel matrix or the maximum over all the local parts.
540d9274352SBarry Smith 
541d9274352SBarry Smith     Level: beginner
542d9274352SBarry Smith 
543d9274352SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
544d9274352SBarry Smith 
545d9274352SBarry Smith .seealso: MatGetInfo(), MatInfo
546d9274352SBarry Smith E*/
5477b80b807SBarry Smith typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType;
5487087cfbeSBarry Smith extern PetscErrorCode  MatGetInfo(Mat,MatInfoType,MatInfo*);
5497087cfbeSBarry Smith extern PetscErrorCode  MatGetDiagonal(Mat,Vec);
5507087cfbeSBarry Smith extern PetscErrorCode  MatGetRowMax(Mat,Vec,PetscInt[]);
5517087cfbeSBarry Smith extern PetscErrorCode  MatGetRowMin(Mat,Vec,PetscInt[]);
5527087cfbeSBarry Smith extern PetscErrorCode  MatGetRowMaxAbs(Mat,Vec,PetscInt[]);
5537087cfbeSBarry Smith extern PetscErrorCode  MatGetRowMinAbs(Mat,Vec,PetscInt[]);
5547087cfbeSBarry Smith extern PetscErrorCode  MatGetRowSum(Mat,Vec);
5557087cfbeSBarry Smith extern PetscErrorCode  MatTranspose(Mat,MatReuse,Mat*);
556fc4dec0aSBarry Smith PetscPolymorphicFunction(MatTranspose,(Mat A),(A,MAT_INITIAL_MATRIX,&t),Mat,t)
5577087cfbeSBarry Smith extern PetscErrorCode  MatHermitianTranspose(Mat,MatReuse,Mat*);
5587087cfbeSBarry Smith extern PetscErrorCode  MatPermute(Mat,IS,IS,Mat *);
559ba966639SSatish Balay PetscPolymorphicFunction(MatPermute,(Mat A,IS is1,IS is2),(A,is1,is2,&t),Mat,t)
5607087cfbeSBarry Smith extern PetscErrorCode  MatDiagonalScale(Mat,Vec,Vec);
5617087cfbeSBarry Smith extern PetscErrorCode  MatDiagonalSet(Mat,Vec,InsertMode);
5627087cfbeSBarry Smith extern PetscErrorCode  MatEqual(Mat,Mat,PetscBool *);
563ace3abfcSBarry Smith PetscPolymorphicFunction(MatEqual,(Mat A,Mat B),(A,B,&t),PetscBool ,t)
5647087cfbeSBarry Smith extern PetscErrorCode  MatMultEqual(Mat,Mat,PetscInt,PetscBool *);
5657087cfbeSBarry Smith extern PetscErrorCode  MatMultAddEqual(Mat,Mat,PetscInt,PetscBool *);
5667087cfbeSBarry Smith extern PetscErrorCode  MatMultTransposeEqual(Mat,Mat,PetscInt,PetscBool *);
5677087cfbeSBarry Smith extern PetscErrorCode  MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscBool *);
5687b80b807SBarry Smith 
5697087cfbeSBarry Smith extern PetscErrorCode  MatNorm(Mat,NormType,PetscReal *);
570ba966639SSatish Balay PetscPolymorphicFunction(MatNorm,(Mat A,NormType t),(A,t,&n),PetscReal,n)
57109573ac7SBarry Smith extern PetscErrorCode  MatGetColumnNorms(Mat,NormType,PetscReal *);
5727087cfbeSBarry Smith extern PetscErrorCode  MatZeroEntries(Mat);
5737087cfbeSBarry Smith extern PetscErrorCode  MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
5747087cfbeSBarry Smith extern PetscErrorCode  MatZeroRowsIS(Mat,IS,PetscScalar,Vec,Vec);
5757087cfbeSBarry Smith extern PetscErrorCode  MatZeroRowsStencil(Mat,PetscInt,const MatStencil [],PetscScalar,Vec,Vec);
5767087cfbeSBarry Smith extern PetscErrorCode  MatZeroRowsColumns(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
5777087cfbeSBarry Smith extern PetscErrorCode  MatZeroRowsColumnsIS(Mat,IS,PetscScalar,Vec,Vec);
5787b80b807SBarry Smith 
5797087cfbeSBarry Smith extern PetscErrorCode  MatUseScaledForm(Mat,PetscBool );
5807087cfbeSBarry Smith extern PetscErrorCode  MatScaleSystem(Mat,Vec,Vec);
5817087cfbeSBarry Smith extern PetscErrorCode  MatUnScaleSystem(Mat,Vec,Vec);
5825ef9f2a5SBarry Smith 
5837087cfbeSBarry Smith extern PetscErrorCode  MatGetSize(Mat,PetscInt*,PetscInt*);
5847087cfbeSBarry Smith extern PetscErrorCode  MatGetLocalSize(Mat,PetscInt*,PetscInt*);
5857087cfbeSBarry Smith extern PetscErrorCode  MatGetOwnershipRange(Mat,PetscInt*,PetscInt*);
5867087cfbeSBarry Smith extern PetscErrorCode  MatGetOwnershipRanges(Mat,const PetscInt**);
5877087cfbeSBarry Smith extern PetscErrorCode  MatGetOwnershipRangeColumn(Mat,PetscInt*,PetscInt*);
5887087cfbeSBarry Smith extern PetscErrorCode  MatGetOwnershipRangesColumn(Mat,const PetscInt**);
5897b80b807SBarry Smith 
5907087cfbeSBarry Smith extern PetscErrorCode  MatGetSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
5917087cfbeSBarry Smith extern PetscErrorCode  MatGetSubMatricesParallel(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
5927087cfbeSBarry Smith extern PetscErrorCode  MatDestroyMatrices(PetscInt,Mat *[]);
5937087cfbeSBarry Smith extern PetscErrorCode  MatGetSubMatrix(Mat,IS,IS,MatReuse,Mat *);
5947087cfbeSBarry Smith extern PetscErrorCode  MatGetLocalSubMatrix(Mat,IS,IS,Mat*);
5957087cfbeSBarry Smith extern PetscErrorCode  MatRestoreLocalSubMatrix(Mat,IS,IS,Mat*);
5967087cfbeSBarry Smith extern PetscErrorCode  MatGetSeqNonzeroStructure(Mat,Mat*);
5977087cfbeSBarry Smith extern PetscErrorCode  MatDestroySeqNonzeroStructure(Mat*);
5985494a064SHong Zhang 
5997087cfbeSBarry Smith extern PetscErrorCode  MatMerge(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
6007087cfbeSBarry Smith extern PetscErrorCode  MatMerge_SeqsToMPI(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*);
6017087cfbeSBarry Smith extern PetscErrorCode  MatMerge_SeqsToMPISymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*);
6027087cfbeSBarry Smith extern PetscErrorCode  MatMerge_SeqsToMPINumeric(Mat,Mat);
6034a2b5492SBarry Smith extern PetscErrorCode  MatMPIAIJGetLocalMat(Mat,MatReuse,Mat*);
6044a2b5492SBarry Smith extern PetscErrorCode  MatMPIAIJGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*);
6057087cfbeSBarry Smith extern PetscErrorCode  MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,PetscInt*,Mat*);
6067087cfbeSBarry Smith extern PetscErrorCode  MatGetBrowsOfAoCols(Mat,Mat,MatReuse,PetscInt**,PetscInt**,MatScalar**,Mat*);
60743eb5e2fSMatthew Knepley #if defined (PETSC_USE_CTABLE)
608cd116e26SSatish Balay #include "petscctable.h"
6097087cfbeSBarry Smith extern PetscErrorCode  MatGetCommunicationStructs(Mat, Vec *, PetscTable *, VecScatter *);
61043eb5e2fSMatthew Knepley #else
6117087cfbeSBarry Smith extern PetscErrorCode  MatGetCommunicationStructs(Mat, Vec *, PetscInt *[], VecScatter *);
61243eb5e2fSMatthew Knepley #endif
6137087cfbeSBarry Smith extern PetscErrorCode  MatGetGhosts(Mat, PetscInt *,const PetscInt *[]);
6148efafbd8SBarry Smith 
6157087cfbeSBarry Smith extern PetscErrorCode  MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt);
6167b80b807SBarry Smith 
6177087cfbeSBarry Smith extern PetscErrorCode  MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*);
6187087cfbeSBarry Smith extern PetscErrorCode  MatMatMultSymbolic(Mat,Mat,PetscReal,Mat*);
6197087cfbeSBarry Smith extern PetscErrorCode  MatMatMultNumeric(Mat,Mat,Mat);
62022440eb1SKris Buschelman 
6217087cfbeSBarry Smith extern PetscErrorCode  MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*);
6227087cfbeSBarry Smith extern PetscErrorCode  MatPtAPSymbolic(Mat,Mat,PetscReal,Mat*);
6237087cfbeSBarry Smith extern PetscErrorCode  MatPtAPNumeric(Mat,Mat,Mat);
62422440eb1SKris Buschelman 
6257087cfbeSBarry Smith extern PetscErrorCode  MatMatMultTranspose(Mat,Mat,MatReuse,PetscReal,Mat*);
6267087cfbeSBarry Smith extern PetscErrorCode  MatMatMultTransposeSymbolic(Mat,Mat,PetscReal,Mat*);
6277087cfbeSBarry Smith extern PetscErrorCode  MatMatMultTransposeNumeric(Mat,Mat,Mat);
628bc011b1eSHong Zhang 
6297087cfbeSBarry Smith extern PetscErrorCode  MatAXPY(Mat,PetscScalar,Mat,MatStructure);
6307087cfbeSBarry Smith extern PetscErrorCode  MatAYPX(Mat,PetscScalar,Mat,MatStructure);
6311c741599SHong Zhang 
6327087cfbeSBarry Smith extern PetscErrorCode  MatScale(Mat,PetscScalar);
6337087cfbeSBarry Smith extern PetscErrorCode  MatShift(Mat,PetscScalar);
6347b80b807SBarry Smith 
6357087cfbeSBarry Smith extern PetscErrorCode  MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping);
6367087cfbeSBarry Smith extern PetscErrorCode  MatSetLocalToGlobalMappingBlock(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping);
6377087cfbeSBarry Smith extern PetscErrorCode  MatGetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*);
6387087cfbeSBarry Smith extern PetscErrorCode  MatGetLocalToGlobalMappingBlock(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*);
6397087cfbeSBarry Smith extern PetscErrorCode  MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
6407087cfbeSBarry Smith extern PetscErrorCode  MatZeroRowsLocalIS(Mat,IS,PetscScalar,Vec,Vec);
6417087cfbeSBarry Smith extern PetscErrorCode  MatZeroRowsColumnsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
6427087cfbeSBarry Smith extern PetscErrorCode  MatZeroRowsColumnsLocalIS(Mat,IS,PetscScalar,Vec,Vec);
6437087cfbeSBarry Smith extern PetscErrorCode  MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
6447087cfbeSBarry Smith extern PetscErrorCode  MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
645052efed2SBarry Smith 
6467087cfbeSBarry Smith extern PetscErrorCode  MatStashSetInitialSize(Mat,PetscInt,PetscInt);
6477087cfbeSBarry Smith extern PetscErrorCode  MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
64890f02eecSBarry Smith 
6497087cfbeSBarry Smith extern PetscErrorCode  MatInterpolate(Mat,Vec,Vec);
6507087cfbeSBarry Smith extern PetscErrorCode  MatInterpolateAdd(Mat,Vec,Vec,Vec);
651ba966639SSatish Balay PetscPolymorphicSubroutine(MatInterpolateAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
6527087cfbeSBarry Smith extern PetscErrorCode  MatRestrict(Mat,Vec,Vec);
6537087cfbeSBarry Smith extern PetscErrorCode  MatGetVecs(Mat,Vec*,Vec*);
6547087cfbeSBarry Smith extern PetscErrorCode  MatGetRedundantMatrix(Mat,PetscInt,MPI_Comm,PetscInt,MatReuse,Mat*);
6557087cfbeSBarry Smith extern PetscErrorCode  MatGetMultiProcBlock(Mat,MPI_Comm,Mat*);
656fc9bc008SSatish Balay extern PetscErrorCode  MatFindZeroDiagonals(Mat,IS*);
657bd481603SBarry Smith 
658bd481603SBarry Smith /*MC
6591620fd73SBarry Smith    MatSetValue - Set a single entry into a matrix.
6601620fd73SBarry Smith 
6611620fd73SBarry Smith    Not collective
6621620fd73SBarry Smith 
6631620fd73SBarry Smith    Input Parameters:
6641620fd73SBarry Smith +  m - the matrix
6651620fd73SBarry Smith .  row - the row location of the entry
6661620fd73SBarry Smith .  col - the column location of the entry
6671620fd73SBarry Smith .  value - the value to insert
6681620fd73SBarry Smith -  mode - either INSERT_VALUES or ADD_VALUES
6691620fd73SBarry Smith 
6701620fd73SBarry Smith    Notes:
6711620fd73SBarry Smith    For efficiency one should use MatSetValues() and set several or many
6721620fd73SBarry Smith    values simultaneously if possible.
6731620fd73SBarry Smith 
6741620fd73SBarry Smith    Level: beginner
6751620fd73SBarry Smith 
6761620fd73SBarry Smith .seealso: MatSetValues(), MatSetValueLocal()
6771620fd73SBarry Smith M*/
6781620fd73SBarry 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);}
6791620fd73SBarry Smith 
6801620fd73SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatGetValue(Mat v,PetscInt i,PetscInt j,PetscScalar *va) {return MatGetValues(v,1,&i,1,&j,va);}
6811620fd73SBarry Smith 
6821620fd73SBarry 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);}
6831620fd73SBarry Smith 
6841620fd73SBarry Smith /*MC
6850d2a7ff2SSatish Balay    MatPreallocateInitialize - Begins the block of code that will count the number of nonzeros per
686bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
687bd481603SBarry Smith 
688bd481603SBarry Smith    Synopsis:
689c1ac3661SBarry Smith    PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)
690bd481603SBarry Smith 
691bd481603SBarry Smith    Collective on MPI_Comm
692bd481603SBarry Smith 
693bd481603SBarry Smith    Input Parameters:
694bd481603SBarry Smith +  comm - the communicator that will share the eventually allocated matrix
695859fcb39SBarry Smith .  nrows - the number of LOCAL rows in the matrix
696859fcb39SBarry Smith -  ncols - the number of LOCAL columns in the matrix
697bd481603SBarry Smith 
698bd481603SBarry Smith    Output Parameters:
699bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
700bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
701bd481603SBarry Smith 
702bd481603SBarry Smith 
703bd481603SBarry Smith    Level: intermediate
704bd481603SBarry Smith 
705bd481603SBarry Smith    Notes:
7060598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
707bd481603SBarry Smith 
7081d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
709bd481603SBarry Smith 
710bd481603SBarry Smith    Use MatPreallocateInitializeSymmetric() for symmetric matrices (MPISBAIJ matrices)
711bd481603SBarry Smith 
7121620fd73SBarry Smith    This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize().
7131620fd73SBarry Smith 
714bd481603SBarry Smith   Concepts: preallocation^Matrix
715bd481603SBarry Smith 
716bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
717bd481603SBarry Smith           MatPreallocateInitializeSymmetric(), MatPreallocateSymmetricSetLocal()
718bd481603SBarry Smith M*/
719c4f061fbSSatish Balay #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \
7207c922b88SBarry Smith { \
721a89bc211SBarry Smith   PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__start,__end; \
722a89bc211SBarry Smith   _4_ierr = PetscMalloc2(__nrows,PetscInt,&dnz,__nrows,PetscInt,&onz);CHKERRQ(_4_ierr); \
723a89bc211SBarry Smith   _4_ierr = PetscMemzero(dnz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
724a89bc211SBarry Smith   _4_ierr = PetscMemzero(onz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
725c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ctmp;\
726a89bc211SBarry Smith   _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows;
7277c922b88SBarry Smith 
728bd481603SBarry Smith /*MC
7290d2a7ff2SSatish Balay    MatPreallocateSymmetricInitialize - Begins the block of code that will count the number of nonzeros per
730bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
731bd481603SBarry Smith 
732bd481603SBarry Smith    Synopsis:
733c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)
734bd481603SBarry Smith 
735bd481603SBarry Smith    Collective on MPI_Comm
736bd481603SBarry Smith 
737bd481603SBarry Smith    Input Parameters:
738bd481603SBarry Smith +  comm - the communicator that will share the eventually allocated matrix
739859fcb39SBarry Smith .  nrows - the number of LOCAL rows in the matrix
740859fcb39SBarry Smith -  ncols - the number of LOCAL columns in the matrix
741bd481603SBarry Smith 
742bd481603SBarry Smith    Output Parameters:
743bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
744bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
745bd481603SBarry Smith 
746bd481603SBarry Smith 
747bd481603SBarry Smith    Level: intermediate
748bd481603SBarry Smith 
749bd481603SBarry Smith    Notes:
7500598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
751bd481603SBarry Smith 
7521d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
753bd481603SBarry Smith 
7541620fd73SBarry Smith    This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize().
7551620fd73SBarry Smith 
756bd481603SBarry Smith   Concepts: preallocation^Matrix
757bd481603SBarry Smith 
758bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
759bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal()
760bd481603SBarry Smith M*/
761222b16d4SBarry Smith #define MatPreallocateSymmetricInitialize(comm,nrows,ncols,dnz,onz) 0; \
762222b16d4SBarry Smith { \
763a89bc211SBarry Smith   PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__end; \
764a89bc211SBarry Smith   _4_ierr = PetscMalloc2(__nrows,PetscInt,&dnz,__nrows,PetscInt,&onz);CHKERRQ(_4_ierr); \
765a89bc211SBarry Smith   _4_ierr = PetscMemzero(dnz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
766a89bc211SBarry Smith   _4_ierr = PetscMemzero(onz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
767c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr);\
768a89bc211SBarry Smith   _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows;
769222b16d4SBarry Smith 
770bd481603SBarry Smith /*MC
771bd481603SBarry Smith    MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
772bd481603SBarry Smith        inserted using a local number of the rows and columns
773bd481603SBarry Smith 
774bd481603SBarry Smith    Synopsis:
775c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
776bd481603SBarry Smith 
777bd481603SBarry Smith    Not Collective
778bd481603SBarry Smith 
779bd481603SBarry Smith    Input Parameters:
780784ac674SJed Brown +  map - the row mapping from local numbering to global numbering
781bd481603SBarry Smith .  nrows - the number of rows indicated
7821d73ed98SBarry Smith .  rows - the indices of the rows
783784ac674SJed Brown .  cmap - the column mapping from local to global numbering
784bd481603SBarry Smith .  ncols - the number of columns in the matrix
785bd481603SBarry Smith .  cols - the columns indicated
786bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
787bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
788bd481603SBarry Smith 
789bd481603SBarry Smith 
790bd481603SBarry Smith    Level: intermediate
791bd481603SBarry Smith 
792bd481603SBarry Smith    Notes:
7930598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
794bd481603SBarry Smith 
7951d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
796bd481603SBarry Smith 
797bd481603SBarry Smith   Concepts: preallocation^Matrix
798bd481603SBarry Smith 
799bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
800bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal()
801bd481603SBarry Smith M*/
802784ac674SJed Brown #define MatPreallocateSetLocal(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \
803c4f061fbSSatish Balay {\
804c1ac3661SBarry Smith   PetscInt __l;\
805784ac674SJed Brown   _4_ierr = ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\
806784ac674SJed Brown   _4_ierr = ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\
807c4f061fbSSatish Balay   for (__l=0;__l<nrows;__l++) {\
808ef66eb69SBarry Smith     _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
809c4f061fbSSatish Balay   }\
810c4f061fbSSatish Balay }
811c4f061fbSSatish Balay 
812bd481603SBarry Smith /*MC
813bd481603SBarry Smith    MatPreallocateSymmetricSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
814bd481603SBarry Smith        inserted using a local number of the rows and columns
815bd481603SBarry Smith 
816bd481603SBarry Smith    Synopsis:
817c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
818bd481603SBarry Smith 
819bd481603SBarry Smith    Not Collective
820bd481603SBarry Smith 
821bd481603SBarry Smith    Input Parameters:
822bd481603SBarry Smith +  map - the mapping between local numbering and global numbering
823bd481603SBarry Smith .  nrows - the number of rows indicated
8241d73ed98SBarry Smith .  rows - the indices of the rows
825bd481603SBarry Smith .  ncols - the number of columns in the matrix
826bd481603SBarry Smith .  cols - the columns indicated
827bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
828bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
829bd481603SBarry Smith 
830bd481603SBarry Smith 
831bd481603SBarry Smith    Level: intermediate
832bd481603SBarry Smith 
833bd481603SBarry Smith    Notes:
8340598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
835bd481603SBarry Smith 
836bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
837bd481603SBarry Smith 
838bd481603SBarry Smith   Concepts: preallocation^Matrix
839bd481603SBarry Smith 
840bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
841bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
842bd481603SBarry Smith M*/
843d3d32019SBarry Smith #define MatPreallocateSymmetricSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\
844d3d32019SBarry Smith {\
845c1ac3661SBarry Smith   PetscInt __l;\
846d3d32019SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\
847d3d32019SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\
848d3d32019SBarry Smith   for (__l=0;__l<nrows;__l++) {\
849d3d32019SBarry Smith     _4_ierr = MatPreallocateSymmetricSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
850d3d32019SBarry Smith   }\
851d3d32019SBarry Smith }
852d3d32019SBarry Smith 
853bd481603SBarry Smith /*MC
854bd481603SBarry Smith    MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
855bd481603SBarry Smith        inserted using a local number of the rows and columns
856bd481603SBarry Smith 
857bd481603SBarry Smith    Synopsis:
858c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
859bd481603SBarry Smith 
860bd481603SBarry Smith    Not Collective
861bd481603SBarry Smith 
862bd481603SBarry Smith    Input Parameters:
86364075487SBarry Smith +  row - the row
864bd481603SBarry Smith .  ncols - the number of columns in the matrix
865a50f8bf6SHong Zhang -  cols - the columns indicated
866a50f8bf6SHong Zhang 
867a50f8bf6SHong Zhang    Output Parameters:
868a50f8bf6SHong Zhang +  dnz - the array that will be passed to the matrix preallocation routines
869bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
870bd481603SBarry Smith 
871bd481603SBarry Smith 
872bd481603SBarry Smith    Level: intermediate
873bd481603SBarry Smith 
874bd481603SBarry Smith    Notes:
8750598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
876bd481603SBarry Smith 
877bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
878bd481603SBarry Smith 
8791620fd73SBarry Smith    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
8801620fd73SBarry Smith 
881bd481603SBarry Smith   Concepts: preallocation^Matrix
882bd481603SBarry Smith 
883bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
884bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
885bd481603SBarry Smith M*/
886c4f061fbSSatish Balay #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\
887c1ac3661SBarry Smith { PetscInt __i; \
888e32f2f54SBarry 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);\
889e32f2f54SBarry 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);\
8907c922b88SBarry Smith   for (__i=0; __i<nc; __i++) {\
89164075487SBarry Smith     if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \
8927cc688d7SBarry Smith     else dnz[row - __rstart]++;\
8937c922b88SBarry Smith   }\
8947c922b88SBarry Smith }
8957c922b88SBarry Smith 
896bd481603SBarry Smith /*MC
897bd481603SBarry Smith    MatPreallocateSymmetricSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
898bd481603SBarry Smith        inserted using a local number of the rows and columns
899bd481603SBarry Smith 
900bd481603SBarry Smith    Synopsis:
901c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
902bd481603SBarry Smith 
903bd481603SBarry Smith    Not Collective
904bd481603SBarry Smith 
905bd481603SBarry Smith    Input Parameters:
906bd481603SBarry Smith +  nrows - the number of rows indicated
9071d73ed98SBarry Smith .  rows - the indices of the rows
908bd481603SBarry Smith .  ncols - the number of columns in the matrix
909bd481603SBarry Smith .  cols - the columns indicated
910bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
911bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
912bd481603SBarry Smith 
913bd481603SBarry Smith 
914bd481603SBarry Smith    Level: intermediate
915bd481603SBarry Smith 
916bd481603SBarry Smith    Notes:
9170598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
918bd481603SBarry Smith 
919bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
920bd481603SBarry Smith 
9211620fd73SBarry Smith    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
9221620fd73SBarry Smith 
923bd481603SBarry Smith   Concepts: preallocation^Matrix
924bd481603SBarry Smith 
925bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
926bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
927bd481603SBarry Smith M*/
928d3d32019SBarry Smith #define MatPreallocateSymmetricSet(row,nc,cols,dnz,onz) 0;\
929c1ac3661SBarry Smith { PetscInt __i; \
930d3d32019SBarry Smith   for (__i=0; __i<nc; __i++) {\
931d3d32019SBarry Smith     if (cols[__i] >= __end) onz[row - __rstart]++; \
932d3d32019SBarry Smith     else if (cols[__i] >= row) dnz[row - __rstart]++;\
933d3d32019SBarry Smith   }\
934d3d32019SBarry Smith }
935d3d32019SBarry Smith 
936bd481603SBarry Smith /*MC
93716371a99SBarry Smith    MatPreallocateLocation -  An alternative to MatPreallocationSet() that puts the nonzero locations into the matrix if it exists
93816371a99SBarry Smith 
93916371a99SBarry Smith    Synopsis:
94016371a99SBarry Smith    PetscErrorCode MatPreallocateLocations(Mat A,PetscInt row,PetscInt ncols,PetscInt *cols,PetscInt *dnz,PetscInt *onz)
94116371a99SBarry Smith 
94216371a99SBarry Smith    Not Collective
94316371a99SBarry Smith 
94416371a99SBarry Smith    Input Parameters:
94516371a99SBarry Smith .  A - matrix
94616371a99SBarry Smith .  row - row where values exist (must be local to this process)
94716371a99SBarry Smith .  ncols - number of columns
94816371a99SBarry Smith .  cols - columns with nonzeros
94916371a99SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
95016371a99SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
95116371a99SBarry Smith 
95216371a99SBarry Smith 
95316371a99SBarry Smith    Level: intermediate
95416371a99SBarry Smith 
95516371a99SBarry Smith    Notes:
9560598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
95716371a99SBarry Smith 
95816371a99SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
95916371a99SBarry Smith 
96016371a99SBarry Smith    This is a MACRO not a function because it uses a bunch of variables private to the MatPreallocation.... routines.
96116371a99SBarry Smith 
96216371a99SBarry Smith   Concepts: preallocation^Matrix
96316371a99SBarry Smith 
96416371a99SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
96516371a99SBarry Smith           MatPreallocateSymmetricInitialize(), MatPreallocateSymmetricSetLocal()
96616371a99SBarry Smith M*/
96716371a99SBarry 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);}
96816371a99SBarry Smith 
96916371a99SBarry Smith 
97016371a99SBarry Smith /*MC
9710d2a7ff2SSatish Balay    MatPreallocateFinalize - Ends the block of code that will count the number of nonzeros per
972bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
973bd481603SBarry Smith 
974bd481603SBarry Smith    Synopsis:
975c1ac3661SBarry Smith    PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz)
976bd481603SBarry Smith 
977bd481603SBarry Smith    Collective on MPI_Comm
978bd481603SBarry Smith 
979bd481603SBarry Smith    Input Parameters:
98016371a99SBarry Smith +  dnz - the array that was be passed to the matrix preallocation routines
981bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
982bd481603SBarry Smith 
983bd481603SBarry Smith 
984bd481603SBarry Smith    Level: intermediate
985bd481603SBarry Smith 
986bd481603SBarry Smith    Notes:
9870598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
988bd481603SBarry Smith 
989bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
990bd481603SBarry Smith 
9911620fd73SBarry Smith    This is a MACRO not a function because it closes the { started in MatPreallocateInitialize().
9921620fd73SBarry Smith 
993bd481603SBarry Smith   Concepts: preallocation^Matrix
994bd481603SBarry Smith 
995bd481603SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
996bd481603SBarry Smith           MatPreallocateSymmetricInitialize(), MatPreallocateSymmetricSetLocal()
997bd481603SBarry Smith M*/
998a89bc211SBarry Smith #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree2(dnz,onz);CHKERRQ(_4_ierr);}
9997c922b88SBarry Smith 
1000bd481603SBarry Smith 
1001bd481603SBarry Smith 
10027b80b807SBarry Smith /* Routines unique to particular data structures */
10038fe8eb27SJed Brown extern PetscErrorCode  MatShellGetContext(Mat,void *);
1004ba966639SSatish Balay PetscPolymorphicFunction(MatShellGetContext,(Mat A),(A,&t),void*,t)
1005698d4c6aSKris Buschelman 
10067087cfbeSBarry Smith extern PetscErrorCode  MatInodeAdjustForInodes(Mat,IS*,IS*);
10077087cfbeSBarry Smith extern PetscErrorCode  MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *);
1008698d4c6aSKris Buschelman 
10097087cfbeSBarry Smith extern PetscErrorCode  MatSeqAIJSetColumnIndices(Mat,PetscInt[]);
10107087cfbeSBarry Smith extern PetscErrorCode  MatSeqBAIJSetColumnIndices(Mat,PetscInt[]);
10117087cfbeSBarry Smith extern PetscErrorCode  MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
10127087cfbeSBarry Smith extern PetscErrorCode  MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
10137087cfbeSBarry Smith extern PetscErrorCode  MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
10147b80b807SBarry Smith 
1015a96a251dSBarry Smith #define MAT_SKIP_ALLOCATION -4
1016a96a251dSBarry Smith 
10177087cfbeSBarry Smith extern PetscErrorCode  MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
1018ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz))
10197087cfbeSBarry Smith extern PetscErrorCode  MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
1020ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqSBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz))
10217087cfbeSBarry Smith extern PetscErrorCode  MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]);
1022ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqAIJSetPreallocation,(Mat A,const PetscInt nnz[]),(A,0,nnz))
1023273d9f13SBarry Smith 
10247087cfbeSBarry Smith extern PetscErrorCode  MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
1025ba966639SSatish Balay PetscPolymorphicSubroutine(MatMPIBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[],const PetscInt onz[]),(A,bs,0,nnz,0,onz))
10267087cfbeSBarry Smith extern PetscErrorCode  MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
10277087cfbeSBarry Smith extern PetscErrorCode  MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
10287087cfbeSBarry Smith extern PetscErrorCode  MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []);
10297087cfbeSBarry Smith extern PetscErrorCode  MatSeqBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
10307087cfbeSBarry Smith extern PetscErrorCode  MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]);
10317087cfbeSBarry Smith extern PetscErrorCode  MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
10327087cfbeSBarry Smith extern PetscErrorCode  MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]);
10337087cfbeSBarry Smith extern PetscErrorCode  MatMPIDenseSetPreallocation(Mat,PetscScalar[]);
10347087cfbeSBarry Smith extern PetscErrorCode  MatSeqDenseSetPreallocation(Mat,PetscScalar[]);
10357087cfbeSBarry Smith extern PetscErrorCode  MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,PetscInt*[]);
10367087cfbeSBarry Smith extern PetscErrorCode  MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,PetscInt*[]);
10377087cfbeSBarry Smith extern PetscErrorCode  MatAdicSetLocalFunction(Mat,void (*)(void));
1038273d9f13SBarry Smith 
10397087cfbeSBarry Smith extern PetscErrorCode  MatSeqDenseSetLDA(Mat,PetscInt);
10407087cfbeSBarry Smith extern PetscErrorCode  MatDenseGetLocalMatrix(Mat,Mat*);
10411b807ce4Svictorle 
10427087cfbeSBarry Smith extern PetscErrorCode  MatStoreValues(Mat);
10437087cfbeSBarry Smith extern PetscErrorCode  MatRetrieveValues(Mat);
10442e8a6d31SBarry Smith 
10457087cfbeSBarry Smith extern PetscErrorCode  MatDAADSetCtx(Mat,void*);
10463a7fca6bSBarry Smith 
1047b3a44c85SBarry Smith extern PetscErrorCode  MatFindNonzeroRows(Mat,IS*);
10487b80b807SBarry Smith /*
10497b80b807SBarry Smith   These routines are not usually accessed directly, rather solving is
105094b7f48cSBarry Smith   done through the KSP and PC interfaces.
10517b80b807SBarry Smith */
10527b80b807SBarry Smith 
1053*76bdecfbSBarry Smith /*J
1054d9274352SBarry Smith     MatOrderingType - String with the name of a PETSc matrix ordering or the creation function
1055d9274352SBarry Smith        with an optional dynamic library name, for example
1056d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:orderingcreate()
1057d9274352SBarry Smith 
1058d9274352SBarry Smith    Level: beginner
1059d9274352SBarry Smith 
1060e5a9bf91SBarry Smith    Cannot use const because the PC objects manipulate the string
1061e5a9bf91SBarry Smith 
1062d9274352SBarry Smith .seealso: MatGetOrdering()
1063*76bdecfbSBarry Smith J*/
10643eaad2caSSatish Balay #define MatOrderingType char*
10652692d6eeSBarry Smith #define MATORDERINGNATURAL     "natural"
10662692d6eeSBarry Smith #define MATORDERINGND          "nd"
10672692d6eeSBarry Smith #define MATORDERING1WD         "1wd"
10682692d6eeSBarry Smith #define MATORDERINGRCM         "rcm"
10692692d6eeSBarry Smith #define MATORDERINGQMD         "qmd"
10702692d6eeSBarry Smith #define MATORDERINGROWLENGTH   "rowlength"
1071312e372aSBarry Smith #define MATORDERINGAMD         "amd"            /* only works if UMFPACK is installed with PETSc */
1072b12f92e5SBarry Smith 
10737087cfbeSBarry Smith extern PetscErrorCode  MatGetOrdering(Mat,const MatOrderingType,IS*,IS*);
10747087cfbeSBarry Smith extern PetscErrorCode  MatGetOrderingList(PetscFList *list);
10757087cfbeSBarry Smith extern PetscErrorCode  MatOrderingRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,const MatOrderingType,IS*,IS*));
107630de9b25SBarry Smith 
107730de9b25SBarry Smith /*MC
10781890ba74SBarry Smith    MatOrderingRegisterDynamic - Adds a new sparse matrix ordering to the matrix package.
107930de9b25SBarry Smith 
108030de9b25SBarry Smith    Synopsis:
10811890ba74SBarry Smith    PetscErrorCode MatOrderingRegisterDynamic(const char *name_ordering,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatOrdering))
108230de9b25SBarry Smith 
108330de9b25SBarry Smith    Not Collective
108430de9b25SBarry Smith 
108530de9b25SBarry Smith    Input Parameters:
10862692d6eeSBarry Smith +  sname - name of ordering (for example MATORDERINGND)
108730de9b25SBarry Smith .  path - location of library where creation routine is
108830de9b25SBarry Smith .  name - name of function that creates the ordering type,a string
108930de9b25SBarry Smith -  function - function pointer that creates the ordering
109030de9b25SBarry Smith 
109130de9b25SBarry Smith    Level: developer
109230de9b25SBarry Smith 
109330de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
109430de9b25SBarry Smith    is ignored.
109530de9b25SBarry Smith 
109630de9b25SBarry Smith    Sample usage:
109730de9b25SBarry Smith .vb
109830de9b25SBarry Smith    MatOrderingRegisterDynamic("my_order",/home/username/my_lib/lib/libO/solaris/mylib.a,
109930de9b25SBarry Smith                "MyOrder",MyOrder);
110030de9b25SBarry Smith .ve
110130de9b25SBarry Smith 
110230de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
110330de9b25SBarry Smith $     MatOrderingSetType(part,"my_order)
110430de9b25SBarry Smith    or at runtime via the option
11052401956bSBarry Smith $     -pc_factor_mat_ordering_type my_order
110630de9b25SBarry Smith 
1107ab901514SBarry Smith    ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values.
110830de9b25SBarry Smith 
110930de9b25SBarry Smith .keywords: matrix, ordering, register
111030de9b25SBarry Smith 
111130de9b25SBarry Smith .seealso: MatOrderingRegisterDestroy(), MatOrderingRegisterAll()
111230de9b25SBarry Smith M*/
1113aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1114f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,0)
1115b12f92e5SBarry Smith #else
1116f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,d)
1117b12f92e5SBarry Smith #endif
111830de9b25SBarry Smith 
11197087cfbeSBarry Smith extern PetscErrorCode  MatOrderingRegisterDestroy(void);
11207087cfbeSBarry Smith extern PetscErrorCode  MatOrderingRegisterAll(const char[]);
1121ace3abfcSBarry Smith extern PetscBool  MatOrderingRegisterAllCalled;
1122b0a32e0cSBarry Smith extern PetscFList MatOrderingList;
1123d4fbbf0eSBarry Smith 
11247087cfbeSBarry Smith extern PetscErrorCode  MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS);
1125a2ce50c7SBarry Smith 
1126d91e6319SBarry Smith /*S
1127d90ac83dSHong Zhang     MatFactorShiftType - Numeric Shift.
1128d90ac83dSHong Zhang 
1129d90ac83dSHong Zhang    Level: beginner
1130d90ac83dSHong Zhang 
1131d90ac83dSHong Zhang S*/
1132d90ac83dSHong Zhang typedef enum {MAT_SHIFT_NONE,MAT_SHIFT_NONZERO,MAT_SHIFT_POSITIVE_DEFINITE,MAT_SHIFT_INBLOCKS} MatFactorShiftType;
1133d90ac83dSHong Zhang extern const char *MatFactorShiftTypes[];
1134d90ac83dSHong Zhang 
1135d90ac83dSHong Zhang /*S
11362401956bSBarry Smith    MatFactorInfo - Data passed into the matrix factorization routines
1137b00f7748SHong Zhang 
113861cad860SBarry Smith    In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE, that is use
113961cad860SBarry Smith $     MatFactorInfo  info(MAT_FACTORINFO_SIZE)
1140b00f7748SHong Zhang 
114115e8a5b3SHong Zhang    Notes: These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC.
1142b00f7748SHong Zhang 
114361cad860SBarry Smith       You can use MatFactorInfoInitialize() to set default values.
114461cad860SBarry Smith 
1145b00f7748SHong Zhang    Level: developer
1146b00f7748SHong Zhang 
1147d7d82daaSBarry Smith .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(),
1148d7d82daaSBarry Smith           MatFactorInfoInitialize()
1149b00f7748SHong Zhang 
1150b00f7748SHong Zhang S*/
1151b00f7748SHong Zhang typedef struct {
115215e8a5b3SHong Zhang   PetscReal     diagonal_fill;  /* force diagonal to fill in if initially not filled */
115385317021SBarry Smith   PetscReal     usedt;
115415e8a5b3SHong Zhang   PetscReal     dt;             /* drop tolerance */
1155b00f7748SHong Zhang   PetscReal     dtcol;          /* tolerance for pivoting */
115615e8a5b3SHong Zhang   PetscReal     dtcount;        /* maximum nonzeros to be allowed per row */
115767e28bfeSBarry Smith   PetscReal     fill;           /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */
1158348344bbSBarry Smith   PetscReal     levels;         /* ICC/ILU(levels) */
1159bcd9e38bSBarry Smith   PetscReal     pivotinblocks;  /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0
1160bcd9e38bSBarry Smith                                    factorization may be faster if do not pivot */
116115e8a5b3SHong Zhang   PetscReal     zeropivot;      /* pivot is called zero if less than this */
1162f4db908eSBarry Smith   PetscReal     shifttype;      /* type of shift added to matrix factor to prevent zero pivots */
1163f4db908eSBarry Smith   PetscReal     shiftamount;     /* how large the shift is */
116415e8a5b3SHong Zhang } MatFactorInfo;
1165ffa6d0a5SLois Curfman McInnes 
11667087cfbeSBarry Smith extern PetscErrorCode  MatFactorInfoInitialize(MatFactorInfo*);
11677087cfbeSBarry Smith extern PetscErrorCode  MatCholeskyFactor(Mat,IS,const MatFactorInfo*);
11687087cfbeSBarry Smith extern PetscErrorCode  MatCholeskyFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
11697087cfbeSBarry Smith extern PetscErrorCode  MatCholeskyFactorNumeric(Mat,Mat,const MatFactorInfo*);
11707087cfbeSBarry Smith extern PetscErrorCode  MatLUFactor(Mat,IS,IS,const MatFactorInfo*);
11717087cfbeSBarry Smith extern PetscErrorCode  MatILUFactor(Mat,IS,IS,const MatFactorInfo*);
11727087cfbeSBarry Smith extern PetscErrorCode  MatLUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
11737087cfbeSBarry Smith extern PetscErrorCode  MatILUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
11747087cfbeSBarry Smith extern PetscErrorCode  MatICCFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
11757087cfbeSBarry Smith extern PetscErrorCode  MatICCFactor(Mat,IS,const MatFactorInfo*);
11767087cfbeSBarry Smith extern PetscErrorCode  MatLUFactorNumeric(Mat,Mat,const MatFactorInfo*);
11777087cfbeSBarry Smith extern PetscErrorCode  MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*);
11787087cfbeSBarry Smith extern PetscErrorCode  MatSolve(Mat,Vec,Vec);
11797087cfbeSBarry Smith extern PetscErrorCode  MatForwardSolve(Mat,Vec,Vec);
11807087cfbeSBarry Smith extern PetscErrorCode  MatBackwardSolve(Mat,Vec,Vec);
11817087cfbeSBarry Smith extern PetscErrorCode  MatSolveAdd(Mat,Vec,Vec,Vec);
11827087cfbeSBarry Smith extern PetscErrorCode  MatSolveTranspose(Mat,Vec,Vec);
11837087cfbeSBarry Smith extern PetscErrorCode  MatSolveTransposeAdd(Mat,Vec,Vec,Vec);
11847087cfbeSBarry Smith extern PetscErrorCode  MatSolves(Mat,Vecs,Vecs);
11858ed539a5SBarry Smith 
11867087cfbeSBarry Smith extern PetscErrorCode  MatSetUnfactored(Mat);
1187bb5a7306SBarry Smith 
1188d91e6319SBarry Smith /*E
1189d91e6319SBarry Smith     MatSORType - What type of (S)SOR to perform
1190bb1eb677SSatish Balay 
1191d91e6319SBarry Smith     Level: beginner
1192d91e6319SBarry Smith 
1193d9274352SBarry Smith    May be bitwise ORd together
1194d9274352SBarry Smith 
1195d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
1196d91e6319SBarry Smith 
11974e7234bfSBarry Smith    MatSORType may be bitwise ORd together, so do not change the numbers
11984e7234bfSBarry Smith 
119941f059aeSBarry Smith .seealso: MatSOR()
1200d91e6319SBarry Smith E*/
1201ee50ffe9SBarry Smith typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3,
1202ee50ffe9SBarry Smith               SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8,
1203ee50ffe9SBarry Smith               SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16,
120484cb2905SBarry Smith               SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType;
12057087cfbeSBarry Smith extern PetscErrorCode  MatSOR(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
12068ed539a5SBarry Smith 
1207d4fbbf0eSBarry Smith /*
1208639f9d9dSBarry Smith     These routines are for efficiently computing Jacobians via finite differences.
1209639f9d9dSBarry Smith */
1210b12f92e5SBarry Smith 
1211*76bdecfbSBarry Smith /*J
1212d9274352SBarry Smith     MatColoringType - String with the name of a PETSc matrix coloring or the creation function
1213d9274352SBarry Smith        with an optional dynamic library name, for example
1214d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:coloringcreate()
1215d9274352SBarry Smith 
1216d9274352SBarry Smith    Level: beginner
1217d9274352SBarry Smith 
1218d9274352SBarry Smith .seealso: MatGetColoring()
1219*76bdecfbSBarry Smith J*/
1220a313700dSBarry Smith #define MatColoringType char*
12212692d6eeSBarry Smith #define MATCOLORINGNATURAL "natural"
12222692d6eeSBarry Smith #define MATCOLORINGSL      "sl"
12232692d6eeSBarry Smith #define MATCOLORINGLF      "lf"
12242692d6eeSBarry Smith #define MATCOLORINGID      "id"
1225b12f92e5SBarry Smith 
12267087cfbeSBarry Smith extern PetscErrorCode  MatGetColoring(Mat,const MatColoringType,ISColoring*);
12277087cfbeSBarry Smith extern PetscErrorCode  MatColoringRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,MatColoringType,ISColoring *));
122830de9b25SBarry Smith 
122930de9b25SBarry Smith /*MC
123030de9b25SBarry Smith    MatColoringRegisterDynamic - Adds a new sparse matrix coloring to the
123130de9b25SBarry Smith                                matrix package.
123230de9b25SBarry Smith 
123330de9b25SBarry Smith    Synopsis:
12341890ba74SBarry Smith    PetscErrorCode MatColoringRegisterDynamic(const char *name_coloring,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatColoring))
123530de9b25SBarry Smith 
123630de9b25SBarry Smith    Not Collective
123730de9b25SBarry Smith 
123830de9b25SBarry Smith    Input Parameters:
12392692d6eeSBarry Smith +  sname - name of Coloring (for example MATCOLORINGSL)
124030de9b25SBarry Smith .  path - location of library where creation routine is
124130de9b25SBarry Smith .  name - name of function that creates the Coloring type, a string
124230de9b25SBarry Smith -  function - function pointer that creates the coloring
124330de9b25SBarry Smith 
124430de9b25SBarry Smith    Level: developer
124530de9b25SBarry Smith 
124630de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
124730de9b25SBarry Smith    is ignored.
124830de9b25SBarry Smith 
124930de9b25SBarry Smith    Sample usage:
125030de9b25SBarry Smith .vb
125130de9b25SBarry Smith    MatColoringRegisterDynamic("my_color",/home/username/my_lib/lib/libO/solaris/mylib.a,
125230de9b25SBarry Smith                "MyColor",MyColor);
125330de9b25SBarry Smith .ve
125430de9b25SBarry Smith 
125530de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
125630de9b25SBarry Smith $     MatColoringSetType(part,"my_color")
125730de9b25SBarry Smith    or at runtime via the option
125830de9b25SBarry Smith $     -mat_coloring_type my_color
125930de9b25SBarry Smith 
1260ab901514SBarry Smith    $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
126130de9b25SBarry Smith 
126230de9b25SBarry Smith .keywords: matrix, Coloring, register
126330de9b25SBarry Smith 
126430de9b25SBarry Smith .seealso: MatColoringRegisterDestroy(), MatColoringRegisterAll()
126530de9b25SBarry Smith M*/
1266aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1267f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,0)
1268b12f92e5SBarry Smith #else
1269f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,d)
1270b12f92e5SBarry Smith #endif
127130de9b25SBarry Smith 
1272ace3abfcSBarry Smith extern PetscBool  MatColoringRegisterAllCalled;
1273f1144a30SSatish Balay 
12747087cfbeSBarry Smith extern PetscErrorCode  MatColoringRegisterAll(const char[]);
12757087cfbeSBarry Smith extern PetscErrorCode  MatColoringRegisterDestroy(void);
12767087cfbeSBarry Smith extern PetscErrorCode  MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*);
1277639f9d9dSBarry Smith 
1278d9274352SBarry Smith /*S
1279d9274352SBarry Smith      MatFDColoring - Object for computing a sparse Jacobian via finite differences
1280d9274352SBarry Smith         and coloring
1281639f9d9dSBarry Smith 
1282d9274352SBarry Smith    Level: beginner
1283d9274352SBarry Smith 
1284d9274352SBarry Smith   Concepts: coloring, sparse Jacobian, finite differences
1285d9274352SBarry Smith 
1286d9274352SBarry Smith .seealso:  MatFDColoringCreate()
1287d9274352SBarry Smith S*/
1288e2a1c21fSSatish Balay typedef struct _p_MatFDColoring* MatFDColoring;
1289639f9d9dSBarry Smith 
12907087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringCreate(Mat,ISColoring,MatFDColoring *);
1291fcfd50ebSBarry Smith extern PetscErrorCode  MatFDColoringDestroy(MatFDColoring*);
12927087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringView(MatFDColoring,PetscViewer);
12937087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*);
12947087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**);
12957087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal);
12967087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringSetFromOptions(MatFDColoring);
12977087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringApply(Mat,MatFDColoring,Vec,MatStructure*,void *);
12987087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringApplyTS(Mat,MatFDColoring,PetscReal,Vec,MatStructure*,void *);
12997087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringSetF(MatFDColoring,Vec);
13007087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,PetscInt*[]);
1301639f9d9dSBarry Smith /*
13020752156aSBarry Smith     These routines are for partitioning matrices: currently used only
13033eda8832SBarry Smith   for adjacency matrix, MatCreateMPIAdj().
13040752156aSBarry Smith */
1305ca161407SBarry Smith 
1306d9274352SBarry Smith /*S
1307d9274352SBarry Smith      MatPartitioning - Object for managing the partitioning of a matrix or graph
1308d9274352SBarry Smith 
1309d9274352SBarry Smith    Level: beginner
1310d9274352SBarry Smith 
1311d9274352SBarry Smith   Concepts: partitioning
1312d9274352SBarry Smith 
1313743a1026Svictorle .seealso:  MatPartitioningCreate(), MatPartitioningType
1314d9274352SBarry Smith S*/
131591e9ee9fSBarry Smith typedef struct _p_MatPartitioning* MatPartitioning;
1316d9274352SBarry Smith 
1317*76bdecfbSBarry Smith /*J
13185bc6e7ccSvictorle     MatPartitioningType - String with the name of a PETSc matrix partitioning or the creation function
1319d9274352SBarry Smith        with an optional dynamic library name, for example
1320d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:partitioningcreate()
1321d9274352SBarry Smith 
1322d9274352SBarry Smith    Level: beginner
1323d9274352SBarry Smith 
1324b547a13bSvictorle .seealso: MatPartitioningCreate(), MatPartitioning
1325*76bdecfbSBarry Smith J*/
1326a313700dSBarry Smith #define MatPartitioningType char*
13272692d6eeSBarry Smith #define MATPARTITIONINGCURRENT  "current"
13282692d6eeSBarry Smith #define MATPARTITIONINGSQUARE   "square"
13292692d6eeSBarry Smith #define MATPARTITIONINGPARMETIS "parmetis"
13302692d6eeSBarry Smith #define MATPARTITIONINGCHACO    "chaco"
13312692d6eeSBarry Smith #define MATPARTITIONINGPARTY    "party"
133250d91057SBarry Smith #define MATPARTITIONINGPTSCOTCH "ptscotch"
133371306c60Sjroman 
1334ca161407SBarry Smith 
13357087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningCreate(MPI_Comm,MatPartitioning*);
13367087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningSetType(MatPartitioning,const MatPartitioningType);
13377087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningSetNParts(MatPartitioning,PetscInt);
13387087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningSetAdjacency(MatPartitioning,Mat);
13397087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]);
13407087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []);
13417087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningApply(MatPartitioning,IS*);
1342fcfd50ebSBarry Smith extern PetscErrorCode  MatPartitioningDestroy(MatPartitioning*);
13432aabb6bbSBarry Smith 
13447087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatPartitioning));
134530de9b25SBarry Smith 
134630de9b25SBarry Smith /*MC
134730de9b25SBarry Smith    MatPartitioningRegisterDynamic - Adds a new sparse matrix partitioning to the
134830de9b25SBarry Smith    matrix package.
134930de9b25SBarry Smith 
135030de9b25SBarry Smith    Synopsis:
13511890ba74SBarry Smith    PetscErrorCode MatPartitioningRegisterDynamic(const char *name_partitioning,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatPartitioning))
135230de9b25SBarry Smith 
135330de9b25SBarry Smith    Not Collective
135430de9b25SBarry Smith 
135530de9b25SBarry Smith    Input Parameters:
13562692d6eeSBarry Smith +  sname - name of partitioning (for example MATPARTITIONINGCURRENT) or parmetis
135730de9b25SBarry Smith .  path - location of library where creation routine is
135830de9b25SBarry Smith .  name - name of function that creates the partitioning type, a string
135930de9b25SBarry Smith -  function - function pointer that creates the partitioning type
136030de9b25SBarry Smith 
136130de9b25SBarry Smith    Level: developer
136230de9b25SBarry Smith 
136330de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
136430de9b25SBarry Smith    is ignored.
136530de9b25SBarry Smith 
136630de9b25SBarry Smith    Sample usage:
136730de9b25SBarry Smith .vb
136830de9b25SBarry Smith    MatPartitioningRegisterDynamic("my_part",/home/username/my_lib/lib/libO/solaris/mylib.a,
136930de9b25SBarry Smith                "MyPartCreate",MyPartCreate);
137030de9b25SBarry Smith .ve
137130de9b25SBarry Smith 
137230de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
137330de9b25SBarry Smith $     MatPartitioningSetType(part,"my_part")
137430de9b25SBarry Smith    or at runtime via the option
137530de9b25SBarry Smith $     -mat_partitioning_type my_part
137630de9b25SBarry Smith 
1377ab901514SBarry Smith    $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
137830de9b25SBarry Smith 
137930de9b25SBarry Smith .keywords: matrix, partitioning, register
138030de9b25SBarry Smith 
138130de9b25SBarry Smith .seealso: MatPartitioningRegisterDestroy(), MatPartitioningRegisterAll()
138230de9b25SBarry Smith M*/
1383aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1384f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,0)
13852aabb6bbSBarry Smith #else
1386f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,d)
13872aabb6bbSBarry Smith #endif
13882aabb6bbSBarry Smith 
1389ace3abfcSBarry Smith extern PetscBool  MatPartitioningRegisterAllCalled;
1390f1144a30SSatish Balay 
13917087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningRegisterAll(const char[]);
13927087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningRegisterDestroy(void);
13932bad1931SBarry Smith 
13947087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningView(MatPartitioning,PetscViewer);
13957087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningSetFromOptions(MatPartitioning);
13967087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningGetType(MatPartitioning,const MatPartitioningType*);
1397ca161407SBarry Smith 
13987087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningParmetisSetCoarseSequential(MatPartitioning);
13997087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningParmetisGetEdgeCut(MatPartitioning, PetscInt *);
14000752156aSBarry Smith 
1401b6956c12SJose Roman typedef enum { MP_CHACO_MULTILEVEL=1,MP_CHACO_SPECTRAL=2,MP_CHACO_LINEAR=4,MP_CHACO_RANDOM=5,MP_CHACO_SCATTERED=6 } MPChacoGlobalType;
1402b6956c12SJose Roman extern const char *MPChacoGlobalTypes[];
1403b6956c12SJose Roman typedef enum { MP_CHACO_KERNIGHAN=1,MP_CHACO_NONE=2 } MPChacoLocalType;
1404b6956c12SJose Roman extern const char *MPChacoLocalTypes[];
1405b6956c12SJose Roman typedef enum { MP_CHACO_LANCZOS=0,MP_CHACO_RQI=1 } MPChacoEigenType;
1406b6956c12SJose Roman extern const char *MPChacoEigenTypes[];
1407b6956c12SJose Roman 
14087087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningChacoSetGlobal(MatPartitioning,MPChacoGlobalType);
1409b6956c12SJose Roman extern PetscErrorCode  MatPartitioningChacoGetGlobal(MatPartitioning,MPChacoGlobalType*);
14107087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningChacoSetLocal(MatPartitioning,MPChacoLocalType);
1411b6956c12SJose Roman extern PetscErrorCode  MatPartitioningChacoGetLocal(MatPartitioning,MPChacoLocalType*);
14127087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal);
14137087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType);
1414b6956c12SJose Roman extern PetscErrorCode  MatPartitioningChacoGetEigenSolver(MatPartitioning,MPChacoEigenType*);
14157087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningChacoSetEigenTol(MatPartitioning,PetscReal);
1416b6956c12SJose Roman extern PetscErrorCode  MatPartitioningChacoGetEigenTol(MatPartitioning,PetscReal*);
14177087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningChacoSetEigenNumber(MatPartitioning,PetscInt);
1418b6956c12SJose Roman extern PetscErrorCode  MatPartitioningChacoGetEigenNumber(MatPartitioning,PetscInt*);
141971306c60Sjroman 
142071306c60Sjroman #define MP_PARTY_OPT "opt"
142171306c60Sjroman #define MP_PARTY_LIN "lin"
142271306c60Sjroman #define MP_PARTY_SCA "sca"
142371306c60Sjroman #define MP_PARTY_RAN "ran"
142471306c60Sjroman #define MP_PARTY_GBF "gbf"
142571306c60Sjroman #define MP_PARTY_GCF "gcf"
142671306c60Sjroman #define MP_PARTY_BUB "bub"
142771306c60Sjroman #define MP_PARTY_DEF "def"
14287087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningPartySetGlobal(MatPartitioning,const char*);
142971306c60Sjroman #define MP_PARTY_HELPFUL_SETS "hs"
143071306c60Sjroman #define MP_PARTY_KERNIGHAN_LIN "kl"
143171306c60Sjroman #define MP_PARTY_NONE "no"
14327087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningPartySetLocal(MatPartitioning,const char*);
14337087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal);
14347087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningPartySetBipart(MatPartitioning,PetscBool);
14357087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscBool);
143671306c60Sjroman 
143750d91057SBarry Smith typedef enum { MP_PTSCOTCH_QUALITY,MP_PTSCOTCH_SPEED,MP_PTSCOTCH_BALANCE,MP_PTSCOTCH_SAFETY,MP_PTSCOTCH_SCALABILITY } MPPTScotchStrategyType;
143850d91057SBarry Smith extern const char *MPPTScotchStrategyTypes[];
1439e0f1cffaSJose Roman 
144050d91057SBarry Smith extern PetscErrorCode MatPartitioningPTScotchSetImbalance(MatPartitioning,PetscReal);
144150d91057SBarry Smith extern PetscErrorCode MatPartitioningPTScotchGetImbalance(MatPartitioning,PetscReal*);
144250d91057SBarry Smith extern PetscErrorCode MatPartitioningPTScotchSetStrategy(MatPartitioning,MPPTScotchStrategyType);
144350d91057SBarry Smith extern PetscErrorCode MatPartitioningPTScotchGetStrategy(MatPartitioning,MPPTScotchStrategyType*);
144471306c60Sjroman 
144509573ac7SBarry Smith extern PetscErrorCode MatMeshToVertexGraph(Mat,PetscInt,Mat*);
144609573ac7SBarry Smith extern PetscErrorCode MatMeshToCellGraph(Mat,PetscInt,Mat*);
1447591294e4SBarry Smith 
14480752156aSBarry Smith /*
14490a835dfdSSatish Balay     If you add entries here you must also add them to finclude/petscmat.h
1450d4fbbf0eSBarry Smith */
14511c1c02c0SLois Curfman McInnes typedef enum { MATOP_SET_VALUES=0,
14521c1c02c0SLois Curfman McInnes                MATOP_GET_ROW=1,
14531c1c02c0SLois Curfman McInnes                MATOP_RESTORE_ROW=2,
14541c1c02c0SLois Curfman McInnes                MATOP_MULT=3,
14551c1c02c0SLois Curfman McInnes                MATOP_MULT_ADD=4,
14567c922b88SBarry Smith                MATOP_MULT_TRANSPOSE=5,
14577c922b88SBarry Smith                MATOP_MULT_TRANSPOSE_ADD=6,
14581c1c02c0SLois Curfman McInnes                MATOP_SOLVE=7,
14591c1c02c0SLois Curfman McInnes                MATOP_SOLVE_ADD=8,
14607c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE=9,
14617c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE_ADD=10,
14621c1c02c0SLois Curfman McInnes                MATOP_LUFACTOR=11,
14631c1c02c0SLois Curfman McInnes                MATOP_CHOLESKYFACTOR=12,
1464a714c06dSBarry Smith                MATOP_SOR=13,
14651c1c02c0SLois Curfman McInnes                MATOP_TRANSPOSE=14,
14661c1c02c0SLois Curfman McInnes                MATOP_GETINFO=15,
14671c1c02c0SLois Curfman McInnes                MATOP_EQUAL=16,
14681c1c02c0SLois Curfman McInnes                MATOP_GET_DIAGONAL=17,
14691c1c02c0SLois Curfman McInnes                MATOP_DIAGONAL_SCALE=18,
14701c1c02c0SLois Curfman McInnes                MATOP_NORM=19,
14711c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_BEGIN=20,
14721c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_END=21,
1473d519adbfSMatthew Knepley                MATOP_SET_OPTION=22,
1474d519adbfSMatthew Knepley                MATOP_ZERO_ENTRIES=23,
1475d519adbfSMatthew Knepley                MATOP_ZERO_ROWS=24,
1476d519adbfSMatthew Knepley                MATOP_LUFACTOR_SYMBOLIC=25,
1477d519adbfSMatthew Knepley                MATOP_LUFACTOR_NUMERIC=26,
1478d519adbfSMatthew Knepley                MATOP_CHOLESKY_FACTOR_SYMBOLIC=27,
1479d519adbfSMatthew Knepley                MATOP_CHOLESKY_FACTOR_NUMERIC=28,
1480d519adbfSMatthew Knepley                MATOP_SETUP_PREALLOCATION=29,
1481d519adbfSMatthew Knepley                MATOP_ILUFACTOR_SYMBOLIC=30,
1482d519adbfSMatthew Knepley                MATOP_ICCFACTOR_SYMBOLIC=31,
1483d519adbfSMatthew Knepley                MATOP_GET_ARRAY=32,
1484d519adbfSMatthew Knepley                MATOP_RESTORE_ARRAY=33,
1485d519adbfSMatthew Knepley                MATOP_DUPLICATE=34,
1486d519adbfSMatthew Knepley                MATOP_FORWARD_SOLVE=35,
1487d519adbfSMatthew Knepley                MATOP_BACKWARD_SOLVE=36,
1488d519adbfSMatthew Knepley                MATOP_ILUFACTOR=37,
1489d519adbfSMatthew Knepley                MATOP_ICCFACTOR=38,
1490d519adbfSMatthew Knepley                MATOP_AXPY=39,
1491d519adbfSMatthew Knepley                MATOP_GET_SUBMATRICES=40,
1492d519adbfSMatthew Knepley                MATOP_INCREASE_OVERLAP=41,
1493d519adbfSMatthew Knepley                MATOP_GET_VALUES=42,
1494d519adbfSMatthew Knepley                MATOP_COPY=43,
1495d519adbfSMatthew Knepley                MATOP_GET_ROW_MAX=44,
1496d519adbfSMatthew Knepley                MATOP_SCALE=45,
1497d519adbfSMatthew Knepley                MATOP_SHIFT=46,
149835153367SBarry Smith                MATOP_DIAGONAL_SET=47,
1499d519adbfSMatthew Knepley                MATOP_ILUDT_FACTOR=48,
1500d519adbfSMatthew Knepley                MATOP_SET_BLOCK_SIZE=49,
1501d519adbfSMatthew Knepley                MATOP_GET_ROW_IJ=50,
1502d519adbfSMatthew Knepley                MATOP_RESTORE_ROW_IJ=51,
1503d519adbfSMatthew Knepley                MATOP_GET_COLUMN_IJ=52,
1504d519adbfSMatthew Knepley                MATOP_RESTORE_COLUMN_IJ=53,
1505d519adbfSMatthew Knepley                MATOP_FDCOLORING_CREATE=54,
1506d519adbfSMatthew Knepley                MATOP_COLORING_PATCH=55,
1507d519adbfSMatthew Knepley                MATOP_SET_UNFACTORED=56,
1508d519adbfSMatthew Knepley                MATOP_PERMUTE=57,
1509d519adbfSMatthew Knepley                MATOP_SET_VALUES_BLOCKED=58,
1510d519adbfSMatthew Knepley                MATOP_GET_SUBMATRIX=59,
1511d519adbfSMatthew Knepley                MATOP_DESTROY=60,
1512d519adbfSMatthew Knepley                MATOP_VIEW=61,
1513d519adbfSMatthew Knepley                MATOP_CONVERT_FROM=62,
1514d519adbfSMatthew Knepley                MATOP_USE_SCALED_FORM=63,
1515d519adbfSMatthew Knepley                MATOP_SCALE_SYSTEM=64,
1516d519adbfSMatthew Knepley                MATOP_UNSCALE_SYSTEM=65,
1517d519adbfSMatthew Knepley                MATOP_SET_LOCAL_TO_GLOBAL_MAP=66,
1518d519adbfSMatthew Knepley                MATOP_SET_VALUES_LOCAL=67,
1519d519adbfSMatthew Knepley                MATOP_ZERO_ROWS_LOCAL=68,
1520d519adbfSMatthew Knepley                MATOP_GET_ROW_MAX_ABS=69,
1521d519adbfSMatthew Knepley                MATOP_GET_ROW_MIN_ABS=70,
1522d519adbfSMatthew Knepley                MATOP_CONVERT=71,
1523d519adbfSMatthew Knepley                MATOP_SET_COLORING=72,
1524d519adbfSMatthew Knepley                MATOP_SET_VALUES_ADIC=73,
1525d519adbfSMatthew Knepley                MATOP_SET_VALUES_ADIFOR=74,
1526d519adbfSMatthew Knepley                MATOP_FD_COLORING_APPLY=75,
1527d519adbfSMatthew Knepley                MATOP_SET_FROM_OPTIONS=76,
1528d519adbfSMatthew Knepley                MATOP_MULT_CON=77,
1529d519adbfSMatthew Knepley                MATOP_MULT_TRANSPOSE_CON=78,
1530d519adbfSMatthew Knepley                MATOP_PERMUTE_SPARSIFY=79,
1531d519adbfSMatthew Knepley                MATOP_MULT_MULTIPLE=80,
1532d519adbfSMatthew Knepley                MATOP_SOLVE_MULTIPLE=81,
1533d519adbfSMatthew Knepley                MATOP_GET_INERTIA=82,
1534d519adbfSMatthew Knepley                MATOP_LOAD=83,
1535d519adbfSMatthew Knepley                MATOP_IS_SYMMETRIC=84,
1536d519adbfSMatthew Knepley                MATOP_IS_HERMITIAN=85,
1537d519adbfSMatthew Knepley                MATOP_IS_STRUCTURALLY_SYMMETRIC=86,
153841f059aeSBarry Smith                MATOP_DUMMY=87,
1539d519adbfSMatthew Knepley                MATOP_GET_VECS=88,
1540d519adbfSMatthew Knepley                MATOP_MAT_MULT=89,
1541d519adbfSMatthew Knepley                MATOP_MAT_MULT_SYMBOLIC=90,
1542d519adbfSMatthew Knepley                MATOP_MAT_MULT_NUMERIC=91,
1543d519adbfSMatthew Knepley                MATOP_PTAP=92,
1544d519adbfSMatthew Knepley                MATOP_PTAP_SYMBOLIC=93,
1545d519adbfSMatthew Knepley                MATOP_PTAP_NUMERIC=94,
1546d519adbfSMatthew Knepley                MATOP_MAT_MULTTRANSPOSE=95,
15471763ac4eSSatish Balay                MATOP_MAT_MULTTRANSPOSE_SYM=96,
15481763ac4eSSatish Balay                MATOP_MAT_MULTTRANSPOSE_NUM=97,
1549d519adbfSMatthew Knepley                MATOP_PTAP_SYMBOLIC_SEQAIJ=98,
1550d519adbfSMatthew Knepley                MATOP_PTAP_NUMERIC_SEQAIJ=99,
1551d519adbfSMatthew Knepley                MATOP_PTAP_SYMBOLIC_MPIAIJ=100,
1552d519adbfSMatthew Knepley                MATOP_PTAP_NUMERIC_MPIAIJ=101,
1553d519adbfSMatthew Knepley                MATOP_CONJUGATE=102,
1554d519adbfSMatthew Knepley                MATOP_SET_SIZES=103,
1555d519adbfSMatthew Knepley                MATOP_SET_VALUES_ROW=104,
1556d519adbfSMatthew Knepley                MATOP_REAL_PART=105,
1557d519adbfSMatthew Knepley                MATOP_IMAG_PART=106,
1558d519adbfSMatthew Knepley                MATOP_GET_ROW_UTRIANGULAR=107,
1559d519adbfSMatthew Knepley                MATOP_RESTORE_ROW_UTRIANGULAR=108,
1560d519adbfSMatthew Knepley                MATOP_MATSOLVE=109,
1561d519adbfSMatthew Knepley                MATOP_GET_REDUNDANTMATRIX=110,
1562d519adbfSMatthew Knepley                MATOP_GET_ROW_MIN=111,
1563d519adbfSMatthew Knepley                MATOP_GET_COLUMN_VEC=112,
1564d519adbfSMatthew Knepley                MATOP_MISSING_DIAGONAL=113,
1565d519adbfSMatthew Knepley                MATOP_MATGETSEQNONZEROSTRUCTURE=114,
156689c6957cSBarry Smith                MATOP_CREATE=115,
156789c6957cSBarry Smith                MATOP_GET_GHOSTS=116,
1568ea38565cSJed Brown                MATOP_GET_LOCALSUBMATRIX=117,
1569ea38565cSJed Brown                MATOP_RESTORE_LOCALSUBMATRIX=118,
1570eeffb40dSHong Zhang                MATOP_MULT_DIAGONAL_BLOCK=119,
1571547795f9SHong Zhang                MATOP_HERMITIANTRANSPOSE=120,
1572547795f9SHong Zhang                MATOP_MULTHERMITIANTRANSPOSE=121,
1573fbdbba38SShri Abhyankar                MATOP_MULTHERMITIANTRANSPOSEADD=122,
1574b9614d88SDmitry Karpeev                MATOP_GETMULTIPROCBLOCK=123,
15750716a85fSBarry Smith                MATOP_GETCOLUMNNORMS=125,
157637868618SMatthew G Knepley 	       MATOP_GET_SUBMATRICES_PARALLEL=128,
157737868618SMatthew G Knepley                MATOP_SET_VALUES_BATCH=129
1578fae171e0SBarry Smith              } MatOperation;
15797087cfbeSBarry Smith extern PetscErrorCode  MatHasOperation(Mat,MatOperation,PetscBool *);
15807087cfbeSBarry Smith extern PetscErrorCode  MatShellSetOperation(Mat,MatOperation,void(*)(void));
15817087cfbeSBarry Smith extern PetscErrorCode  MatShellGetOperation(Mat,MatOperation,void(**)(void));
15827087cfbeSBarry Smith extern PetscErrorCode  MatShellSetContext(Mat,void*);
1583112a2221SBarry Smith 
158490ace30eSBarry Smith /*
158590ace30eSBarry Smith    Codes for matrices stored on disk. By default they are
158690ace30eSBarry Smith    stored in a universal format. By changing the format with
15877973ad04SBarry Smith    PetscViewerSetFormat(viewer,PETSC_VIEWER_NATIVE); the matrices will
158890ace30eSBarry Smith    be stored in a way natural for the matrix, for example dense matrices
158990ace30eSBarry Smith    would be stored as dense. Matrices stored this way may only be
1590f4403165SShri Abhyankar    read into matrices of the same type.
159190ace30eSBarry Smith */
159290ace30eSBarry Smith #define MATRIX_BINARY_FORMAT_DENSE -1
159390ace30eSBarry Smith 
15947087cfbeSBarry Smith extern PetscErrorCode  MatMPIBAIJSetHashTableFactor(Mat,PetscReal);
15957087cfbeSBarry Smith extern PetscErrorCode  MatISGetLocalMat(Mat,Mat*);
15961f4e1ec7SBarry Smith 
1597d9274352SBarry Smith /*S
1598d9274352SBarry Smith      MatNullSpace - Object that removes a null space from a vector, i.e.
1599d9274352SBarry Smith          orthogonalizes the vector to a subsapce
1600d9274352SBarry Smith 
1601f7a9e4ceSBarry Smith    Level: advanced
1602d9274352SBarry Smith 
1603d9274352SBarry Smith   Concepts: matrix; linear operator, null space
1604d9274352SBarry Smith 
16056e1639daSBarry Smith   Users manual sections:
16066e1639daSBarry Smith .   sec_singular
16076e1639daSBarry Smith 
1608d9274352SBarry Smith .seealso:  MatNullSpaceCreate()
1609d9274352SBarry Smith S*/
161074637425SBarry Smith typedef struct _p_MatNullSpace* MatNullSpace;
1611d9274352SBarry Smith 
16127087cfbeSBarry Smith extern PetscErrorCode  MatNullSpaceCreate(MPI_Comm,PetscBool ,PetscInt,const Vec[],MatNullSpace*);
16137087cfbeSBarry Smith extern PetscErrorCode  MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(MatNullSpace,Vec,void*),void*);
1614d34fcf5fSBarry Smith extern PetscErrorCode  MatNullSpaceDestroy(MatNullSpace*);
16157087cfbeSBarry Smith extern PetscErrorCode  MatNullSpaceRemove(MatNullSpace,Vec,Vec*);
16167087cfbeSBarry Smith extern PetscErrorCode  MatNullSpaceAttach(Mat,MatNullSpace);
16177087cfbeSBarry Smith extern PetscErrorCode  MatNullSpaceTest(MatNullSpace,Mat,PetscBool  *);
1618b717e993SJed Brown extern PetscErrorCode  MatNullSpaceView(MatNullSpace,PetscViewer);
161974637425SBarry Smith 
16207087cfbeSBarry Smith extern PetscErrorCode  MatReorderingSeqSBAIJ(Mat,IS);
16217087cfbeSBarry Smith extern PetscErrorCode  MatMPISBAIJSetHashTableFactor(Mat,PetscReal);
16227087cfbeSBarry Smith extern PetscErrorCode  MatSeqSBAIJSetColumnIndices(Mat,PetscInt *);
16237087cfbeSBarry Smith extern PetscErrorCode  MatSeqBAIJInvertBlockDiagonal(Mat);
16243f1d51d7SBarry Smith 
16257087cfbeSBarry Smith extern PetscErrorCode  MatCreateMAIJ(Mat,PetscInt,Mat*);
16267087cfbeSBarry Smith extern PetscErrorCode  MatMAIJRedimension(Mat,PetscInt,Mat*);
16277087cfbeSBarry Smith extern PetscErrorCode  MatMAIJGetAIJ(Mat,Mat*);
1628c4f061fbSSatish Balay 
16297087cfbeSBarry Smith extern PetscErrorCode  MatComputeExplicitOperator(Mat,Mat*);
1630b0a32e0cSBarry Smith 
16317087cfbeSBarry Smith extern PetscErrorCode  MatDiagonalScaleLocal(Mat,Vec);
163204f1ad80SBarry Smith 
16337087cfbeSBarry Smith extern PetscErrorCode  MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*);
16347087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetBase(Mat,Vec,Vec);
16357087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*);
16367087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*));
16377087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec));
16387087cfbeSBarry Smith extern PetscErrorCode  MatMFFDAddNullSpace(Mat,MatNullSpace);
16397087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt);
16407087cfbeSBarry Smith extern PetscErrorCode  MatMFFDResetHHistory(Mat);
16417087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetFunctionError(Mat,PetscReal);
16427087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetPeriod(Mat,PetscInt);
16437087cfbeSBarry Smith extern PetscErrorCode  MatMFFDGetH(Mat,PetscScalar *);
16447087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetOptionsPrefix(Mat,const char[]);
16457087cfbeSBarry Smith extern PetscErrorCode  MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*);
16467087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*);
1647e884886fSBarry Smith 
16486370053bSBarry Smith /*S
16496370053bSBarry Smith     MatMFFD - A data structured used to manage the computation of the h differencing parameter for matrix-free
16506370053bSBarry Smith               Jacobian vector products
1651e884886fSBarry Smith 
16526370053bSBarry Smith     Notes: MATMFFD is a specific MatType which uses the MatMFFD data structure
16536370053bSBarry Smith 
16546370053bSBarry Smith            MatMFFD*() methods actually take the Mat as their first argument. Not a MatMFFD data structure
16556370053bSBarry Smith 
16566370053bSBarry Smith     Level: developer
16576370053bSBarry Smith 
16586370053bSBarry Smith .seealso: MATMFFD, MatCreateMFFD(), MatMFFDSetFuction(), MatMFFDSetType(), MatMFFDRegister()
16596370053bSBarry Smith S*/
1660e884886fSBarry Smith typedef struct _p_MatMFFD* MatMFFD;
1661e884886fSBarry Smith 
1662*76bdecfbSBarry Smith /*J
1663e884886fSBarry Smith     MatMFFDType - algorithm used to compute the h used in computing matrix-vector products via differencing of the function
1664e884886fSBarry Smith 
1665e884886fSBarry Smith    Level: beginner
1666e884886fSBarry Smith 
1667e884886fSBarry Smith .seealso: MatMFFDSetType(), MatMFFDRegister()
1668*76bdecfbSBarry Smith J*/
1669a313700dSBarry Smith #define MatMFFDType char*
1670e884886fSBarry Smith #define MATMFFD_DS  "ds"
1671e884886fSBarry Smith #define MATMFFD_WP  "wp"
1672e884886fSBarry Smith 
16737087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetType(Mat,const MatMFFDType);
16747087cfbeSBarry Smith extern PetscErrorCode  MatMFFDRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatMFFD));
1675e884886fSBarry Smith 
1676e884886fSBarry Smith /*MC
1677e884886fSBarry Smith    MatMFFDRegisterDynamic - Adds a method to the MatMFFD registry.
1678e884886fSBarry Smith 
1679e884886fSBarry Smith    Synopsis:
16801890ba74SBarry Smith    PetscErrorCode MatMFFDRegisterDynamic(const char *name_solver,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatMFFD))
1681e884886fSBarry Smith 
1682e884886fSBarry Smith    Not Collective
1683e884886fSBarry Smith 
1684e884886fSBarry Smith    Input Parameters:
1685e884886fSBarry Smith +  name_solver - name of a new user-defined compute-h module
1686e884886fSBarry Smith .  path - path (either absolute or relative) the library containing this solver
1687e884886fSBarry Smith .  name_create - name of routine to create method context
1688e884886fSBarry Smith -  routine_create - routine to create method context
1689e884886fSBarry Smith 
1690e884886fSBarry Smith    Level: developer
1691e884886fSBarry Smith 
1692e884886fSBarry Smith    Notes:
1693e884886fSBarry Smith    MatMFFDRegisterDynamic() may be called multiple times to add several user-defined solvers.
1694e884886fSBarry Smith 
1695e884886fSBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
1696e884886fSBarry Smith    is ignored.
1697e884886fSBarry Smith 
1698e884886fSBarry Smith    Sample usage:
1699e884886fSBarry Smith .vb
1700e884886fSBarry Smith    MatMFFDRegisterDynamic("my_h",/home/username/my_lib/lib/libO/solaris/mylib.a,
1701e884886fSBarry Smith                "MyHCreate",MyHCreate);
1702e884886fSBarry Smith .ve
1703e884886fSBarry Smith 
1704e884886fSBarry Smith    Then, your solver can be chosen with the procedural interface via
1705e884886fSBarry Smith $     MatMFFDSetType(mfctx,"my_h")
1706e884886fSBarry Smith    or at runtime via the option
1707e884886fSBarry Smith $     -snes_mf_type my_h
1708e884886fSBarry Smith 
1709e884886fSBarry Smith .keywords: MatMFFD, register
1710e884886fSBarry Smith 
1711e884886fSBarry Smith .seealso: MatMFFDRegisterAll(), MatMFFDRegisterDestroy()
1712e884886fSBarry Smith M*/
1713e884886fSBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1714e884886fSBarry Smith #define MatMFFDRegisterDynamic(a,b,c,d) MatMFFDRegister(a,b,c,0)
1715e884886fSBarry Smith #else
1716e884886fSBarry Smith #define MatMFFDRegisterDynamic(a,b,c,d) MatMFFDRegister(a,b,c,d)
1717e884886fSBarry Smith #endif
1718e884886fSBarry Smith 
17197087cfbeSBarry Smith extern PetscErrorCode  MatMFFDRegisterAll(const char[]);
17207087cfbeSBarry Smith extern PetscErrorCode  MatMFFDRegisterDestroy(void);
17217087cfbeSBarry Smith extern PetscErrorCode  MatMFFDDSSetUmin(Mat,PetscReal);
17227087cfbeSBarry Smith extern PetscErrorCode  MatMFFDWPSetComputeNormU(Mat,PetscBool );
1723e884886fSBarry Smith 
1724e884886fSBarry Smith 
17257087cfbeSBarry Smith extern PetscErrorCode  PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *);
17267087cfbeSBarry Smith extern PetscErrorCode  PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *);
17277dbadf16SMatthew Knepley 
172897969023SHong Zhang /*
172997969023SHong Zhang    PETSc interface to MUMPS
173097969023SHong Zhang */
173197969023SHong Zhang #ifdef PETSC_HAVE_MUMPS
1732f250808bSBarry Smith extern PetscErrorCode  MatMumpsSetIcntl(Mat,PetscInt,PetscInt);
173397969023SHong Zhang #endif
173423a5497aSJed Brown 
1735d954db57SHong Zhang /*
1736d954db57SHong Zhang    PETSc interface to SUPERLU
1737d954db57SHong Zhang */
1738d954db57SHong Zhang #ifdef PETSC_HAVE_SUPERLU
1739f250808bSBarry Smith extern PetscErrorCode  MatSuperluSetILUDropTol(Mat,PetscReal);
1740d954db57SHong Zhang #endif
1741d954db57SHong Zhang 
174290c53902SBarry Smith #if defined(PETSC_HAVE_CUSP)
174390c53902SBarry Smith extern PetscErrorCode  MatCreateSeqAIJCUSP(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
174490c53902SBarry Smith extern PetscErrorCode  MatCreateMPIAIJCUSP(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
174590c53902SBarry Smith #endif
174690c53902SBarry Smith 
174754efbe56SHong Zhang /*
174854efbe56SHong Zhang    PETSc interface to FFTW
174954efbe56SHong Zhang */
175054efbe56SHong Zhang #if defined(PETSC_HAVE_FFTW)
175174a26884SAmlan Barua extern PetscErrorCode VecScatterPetscToFFTW(Mat,Vec,Vec);
175274a26884SAmlan Barua extern PetscErrorCode VecScatterFFTWToPetsc(Mat,Vec,Vec);
17534be45526SHong Zhang extern PetscErrorCode MatGetVecsFFTW(Mat,Vec*,Vec*,Vec*);
175454efbe56SHong Zhang #endif
175554efbe56SHong Zhang 
1756c8883902SJed Brown extern PetscErrorCode MatCreateNest(MPI_Comm,PetscInt,const IS[],PetscInt,const IS[],const Mat[],Mat*);
1757c8883902SJed Brown extern PetscErrorCode MatNestGetSize(Mat,PetscInt*,PetscInt*);
1758c8883902SJed Brown extern PetscErrorCode MatNestGetSubMats(Mat,PetscInt*,PetscInt*,Mat***);
1759c8883902SJed Brown extern PetscErrorCode MatNestGetSubMat(Mat,PetscInt,PetscInt,Mat*);
17607087cfbeSBarry Smith extern PetscErrorCode MatNestSetVecType(Mat,const VecType);
1761c8883902SJed Brown extern PetscErrorCode MatNestSetSubMats(Mat,PetscInt,const IS[],PetscInt,const IS[],const Mat[]);
17620782ca92SJed Brown extern PetscErrorCode MatNestSetSubMat(Mat,PetscInt,PetscInt,Mat);
1763d8588912SDave May 
176423a5497aSJed Brown PETSC_EXTERN_CXX_END
176523a5497aSJed Brown #endif
1766