xref: /petsc/include/petscmat.h (revision c0cdd4a1e0711c7ed0175bf9f41df2ec9c8de89c)
12eac72dbSBarry Smith /*
22eac72dbSBarry Smith      Include file for the matrix component of PETSc
32eac72dbSBarry Smith */
40a835dfdSSatish Balay #ifndef __PETSCMAT_H
50a835dfdSSatish Balay #define __PETSCMAT_H
60a835dfdSSatish Balay #include "petscvec.h"
7e9fa29b7SSatish Balay PETSC_EXTERN_CXX_BEGIN
82eac72dbSBarry Smith 
9d9274352SBarry Smith /*S
10d9274352SBarry Smith      Mat - Abstract PETSc matrix object
112eac72dbSBarry Smith 
12d91e6319SBarry Smith    Level: beginner
13d91e6319SBarry Smith 
14d9274352SBarry Smith   Concepts: matrix; linear operator
15d9274352SBarry Smith 
16d9274352SBarry Smith .seealso:  MatCreate(), MatType, MatSetType()
17d9274352SBarry Smith S*/
18d9274352SBarry Smith typedef struct _p_Mat*           Mat;
19d9274352SBarry Smith 
20d9274352SBarry Smith /*E
21d9274352SBarry Smith     MatType - String with the name of a PETSc matrix or the creation function
22d9274352SBarry Smith        with an optional dynamic library name, for example
23d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:mymatcreate()
24d9274352SBarry Smith 
25d9274352SBarry Smith    Level: beginner
26d9274352SBarry Smith 
27c7393fdbSBarry Smith .seealso: MatSetType(), Mat, MatSolverPackage
28d91e6319SBarry Smith E*/
29a313700dSBarry Smith #define MatType char*
30273d9f13SBarry Smith #define MATSAME            "same"
315a11e1b2SBarry Smith #define MATMAIJ            "maij"
32273d9f13SBarry Smith #define MATSEQMAIJ           "seqmaij"
33273d9f13SBarry Smith #define MATMPIMAIJ           "mpimaij"
34273d9f13SBarry Smith #define MATIS              "is"
355a11e1b2SBarry Smith #define MATAIJ             "aij"
36273d9f13SBarry Smith #define MATSEQAIJ            "seqaij"
37273d9f13SBarry Smith #define MATMPIAIJ            "mpiaij"
385a11e1b2SBarry Smith #define MATAIJCRL              "aijcrl"
395a11e1b2SBarry Smith #define MATSEQAIJCRL             "seqaijcrl"
405a11e1b2SBarry Smith #define MATMPIAIJCRL             "mpiaijcrl"
418154be41SBarry Smith #define MATAIJCUSP             "aijcusp"
428154be41SBarry Smith #define MATSEQAIJCUSP            "seqaijcusp"
438154be41SBarry Smith #define MATMPIAIJCUSP            "mpiaijcusp"
445a11e1b2SBarry Smith #define MATAIJPERM             "aijperm"
455a11e1b2SBarry Smith #define MATSEQAIJPERM            "seqaijperm"
465a11e1b2SBarry Smith #define MATMPIAIJPERM            "mpiaijperm"
47273d9f13SBarry Smith #define MATSHELL           "shell"
485a11e1b2SBarry Smith #define MATDENSE           "dense"
49209238afSKris Buschelman #define MATSEQDENSE          "seqdense"
50273d9f13SBarry Smith #define MATMPIDENSE          "mpidense"
515a11e1b2SBarry Smith #define MATBAIJ            "baij"
52273d9f13SBarry Smith #define MATSEQBAIJ           "seqbaij"
53273d9f13SBarry Smith #define MATMPIBAIJ           "mpibaij"
54273d9f13SBarry Smith #define MATMPIADJ          "mpiadj"
555a11e1b2SBarry Smith #define MATSBAIJ           "sbaij"
56273d9f13SBarry Smith #define MATSEQSBAIJ          "seqsbaij"
57273d9f13SBarry Smith #define MATMPISBAIJ          "mpisbaij"
58*c0cdd4a1SDahai Guo 
59*c0cdd4a1SDahai Guo #define MATSEQBSTRM        "seqbstrm"
60*c0cdd4a1SDahai Guo #define MATMPIBSTRM        "mpibstrm"
61*c0cdd4a1SDahai Guo #define MATBSTRM           "bstrm"
62*c0cdd4a1SDahai Guo 
63cebc7f6cSBarry Smith #define MATDAAD            "daad"
64cebc7f6cSBarry Smith #define MATMFFD            "mffd"
65c8a8475eSBarry Smith #define MATNORMAL          "normal"
66ab92ecdeSBarry Smith #define MATLRC             "lrc"
672a6744ebSBarry Smith #define MATSCATTER         "scatter"
68421e10b8SBarry Smith #define MATBLOCKMAT        "blockmat"
69793850ffSBarry Smith #define MATCOMPOSITE       "composite"
701f8c7532SHong Zhang #define MATFFT             "fft"
711f8c7532SHong Zhang #define MATFFTW              "fftw"
72e133240eSMatthew G Knepley #define MATSEQCUFFT          "seqcufft"
73557cca28SSatish Balay #define MATTRANSPOSEMAT    "transpose"
7472ca8751SBarry Smith #define MATSCHURCOMPLEMENT "schurcomplement"
751d6018f0SLisandro Dalcin #define MATPYTHON          "python"
76f91d8e95SBarry Smith #define MATHYPRESTRUCT     "hyprestruct"
77a9e6138eSGlenn Hammond #define MATHYPRESSTRUCT    "hypresstruct"
78ee1cef2cSJed Brown #define MATSUBMATRIX       "submatrix"
792c0dbf93SJed Brown #define MATLOCALREF        "localref"
80d8588912SDave May #define MATNEST            "nest"
81773941d6SBarry Smith 
829989ab13SBarry Smith /*E
83c7393fdbSBarry Smith     MatSolverPackage - String with the name of a PETSc matrix solver type.
849989ab13SBarry Smith 
859989ab13SBarry Smith     For example: "petsc" indicates what PETSc provides, "superlu" indicates either
869989ab13SBarry Smith        SuperLU or SuperLU_Dist etc.
879989ab13SBarry Smith 
889989ab13SBarry Smith 
899989ab13SBarry Smith    Level: beginner
909989ab13SBarry Smith 
915c9eb25fSBarry Smith .seealso: MatGetFactor(), Mat, MatSetType(), MatType
929989ab13SBarry Smith E*/
93c7393fdbSBarry Smith #define MatSolverPackage char*
942692d6eeSBarry Smith #define MATSOLVERSPOOLES      "spooles"
952692d6eeSBarry Smith #define MATSOLVERSUPERLU      "superlu"
962692d6eeSBarry Smith #define MATSOLVERSUPERLU_DIST "superlu_dist"
972692d6eeSBarry Smith #define MATSOLVERUMFPACK      "umfpack"
9820db9a53SJed Brown #define MATSOLVERCHOLMOD      "cholmod"
992692d6eeSBarry Smith #define MATSOLVERESSL         "essl"
1002692d6eeSBarry Smith #define MATSOLVERLUSOL        "lusol"
1012692d6eeSBarry Smith #define MATSOLVERMUMPS        "mumps"
1022692d6eeSBarry Smith #define MATSOLVERPASTIX       "pastix"
1032692d6eeSBarry Smith #define MATSOLVERDSCPACK      "dscpack"
1042692d6eeSBarry Smith #define MATSOLVERMATLAB       "matlab"
1052692d6eeSBarry Smith #define MATSOLVERPETSC        "petsc"
1062692d6eeSBarry Smith #define MATSOLVERPLAPACK      "plapack"
1072692d6eeSBarry Smith #define MATSOLVERBAS          "bas"
108773941d6SBarry Smith 
109*c0cdd4a1SDahai Guo #define MAT_SOLVER_BSTRM       "bstrm"
110*c0cdd4a1SDahai Guo 
111b24902e0SBarry Smith /*E
112b24902e0SBarry Smith     MatFactorType - indicates what type of factorization is requested
113b24902e0SBarry Smith 
114b24902e0SBarry Smith     Level: beginner
115b24902e0SBarry Smith 
116b24902e0SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
117b24902e0SBarry Smith 
118c7393fdbSBarry Smith .seealso: MatSolverPackage, MatGetFactor()
119b24902e0SBarry Smith E*/
120599ef60dSHong Zhang typedef enum {MAT_FACTOR_NONE, MAT_FACTOR_LU, MAT_FACTOR_CHOLESKY, MAT_FACTOR_ILU, MAT_FACTOR_ICC,MAT_FACTOR_ILUDT} MatFactorType;
12134918c53SJed Brown extern const char *const MatFactorTypes[];
122e92e720dSBarry Smith 
1237087cfbeSBarry Smith extern PetscErrorCode  MatGetFactor(Mat,const MatSolverPackage,MatFactorType,Mat*);
1247087cfbeSBarry Smith extern PetscErrorCode  MatGetFactorAvailable(Mat,const MatSolverPackage,MatFactorType,PetscBool *);
1257087cfbeSBarry Smith extern PetscErrorCode  MatFactorGetSolverPackage(Mat,const MatSolverPackage*);
1267087cfbeSBarry Smith extern PetscErrorCode  MatGetFactorType(Mat,MatFactorType*);
1279989ab13SBarry Smith 
128c06d978dSMatthew Knepley /* Logging support */
1290700a824SBarry Smith #define    MAT_FILE_CLASSID 1211216    /* used to indicate matrices in binary files */
1307087cfbeSBarry Smith extern PetscClassId  MAT_CLASSID;
1317087cfbeSBarry Smith extern PetscClassId  MAT_FDCOLORING_CLASSID;
1327087cfbeSBarry Smith extern PetscClassId  MAT_PARTITIONING_CLASSID;
1337087cfbeSBarry Smith extern PetscClassId  MAT_NULLSPACE_CLASSID;
1347087cfbeSBarry Smith extern PetscClassId  MATMFFD_CLASSID;
135c06d978dSMatthew Knepley 
136ceb03754SKris Buschelman /*E
137ceb03754SKris Buschelman     MatReuse - Indicates if matrices obtained from a previous call to MatGetSubMatrices()
138d6eff37eSBarry Smith      or MatGetSubMatrix() are to be reused to store the new matrix values. For MatConvert() is used to indicate
139d6eff37eSBarry Smith      that the input matrix is to be replaced with the converted matrix.
140ceb03754SKris Buschelman 
141ceb03754SKris Buschelman     Level: beginner
142ceb03754SKris Buschelman 
143ceb03754SKris Buschelman    Any additions/changes here MUST also be made in include/finclude/petscmat.h
144ceb03754SKris Buschelman 
1450c451bc4SBarry Smith .seealso: MatGetSubMatrices(), MatGetSubMatrix(), MatDestroyMatrices(), MatConvert()
146ceb03754SKris Buschelman E*/
147dfe085dbSJed Brown typedef enum {MAT_INITIAL_MATRIX,MAT_REUSE_MATRIX,MAT_IGNORE_MATRIX} MatReuse;
148ceb03754SKris Buschelman 
1495494a064SHong Zhang /*E
1505494a064SHong Zhang     MatGetSubMatrixOption - Indicates if matrices obtained from a call to MatGetSubMatrices()
151829201f2SHong Zhang      include the matrix values. Currently it is only used by MatGetSeqNonzerostructure().
1525494a064SHong Zhang 
1535494a064SHong Zhang     Level: beginner
1545494a064SHong Zhang 
155829201f2SHong Zhang .seealso: MatGetSeqNonzerostructure()
1565494a064SHong Zhang E*/
1575494a064SHong Zhang typedef enum {MAT_DO_NOT_GET_VALUES,MAT_GET_VALUES} MatGetSubMatrixOption;
1585494a064SHong Zhang 
1597087cfbeSBarry Smith extern PetscErrorCode  MatInitializePackage(const char[]);
160c06d978dSMatthew Knepley 
1617087cfbeSBarry Smith extern PetscErrorCode  MatCreate(MPI_Comm,Mat*);
162a69a960dSMatthew Knepley PetscPolymorphicFunction(MatCreate,(MPI_Comm comm),(comm,&A),Mat,A)
163a69a960dSMatthew Knepley PetscPolymorphicFunction(MatCreate,(),(PETSC_COMM_WORLD,&A),Mat,A)
1647087cfbeSBarry Smith extern PetscErrorCode  MatSetSizes(Mat,PetscInt,PetscInt,PetscInt,PetscInt);
1657087cfbeSBarry Smith extern PetscErrorCode  MatSetType(Mat,const MatType);
1667087cfbeSBarry Smith extern PetscErrorCode  MatSetFromOptions(Mat);
1677087cfbeSBarry Smith extern PetscErrorCode  MatSetUpPreallocation(Mat);
1687087cfbeSBarry Smith extern PetscErrorCode  MatRegisterAll(const char[]);
1697087cfbeSBarry Smith extern PetscErrorCode  MatRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat));
17030de9b25SBarry Smith 
1717087cfbeSBarry Smith extern PetscErrorCode  MatSetOptionsPrefix(Mat,const char[]);
1727087cfbeSBarry Smith extern PetscErrorCode  MatAppendOptionsPrefix(Mat,const char[]);
1737087cfbeSBarry Smith extern PetscErrorCode  MatGetOptionsPrefix(Mat,const char*[]);
174f69a0ea3SMatthew Knepley 
17530de9b25SBarry Smith /*MC
17630de9b25SBarry Smith    MatRegisterDynamic - Adds a new matrix type
17730de9b25SBarry Smith 
17830de9b25SBarry Smith    Synopsis:
1791890ba74SBarry Smith    PetscErrorCode MatRegisterDynamic(const char *name,const char *path,const char *name_create,PetscErrorCode (*routine_create)(Mat))
18030de9b25SBarry Smith 
18130de9b25SBarry Smith    Not Collective
18230de9b25SBarry Smith 
18330de9b25SBarry Smith    Input Parameters:
18430de9b25SBarry Smith +  name - name of a new user-defined matrix type
18530de9b25SBarry Smith .  path - path (either absolute or relative) the library containing this solver
18630de9b25SBarry Smith .  name_create - name of routine to create method context
18730de9b25SBarry Smith -  routine_create - routine to create method context
18830de9b25SBarry Smith 
18930de9b25SBarry Smith    Notes:
19030de9b25SBarry Smith    MatRegisterDynamic() may be called multiple times to add several user-defined solvers.
19130de9b25SBarry Smith 
19230de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
19330de9b25SBarry Smith    is ignored.
19430de9b25SBarry Smith 
19530de9b25SBarry Smith    Sample usage:
19630de9b25SBarry Smith .vb
19730de9b25SBarry Smith    MatRegisterDynamic("my_mat",/home/username/my_lib/lib/libO/solaris/mylib.a,
19830de9b25SBarry Smith                "MyMatCreate",MyMatCreate);
19930de9b25SBarry Smith .ve
20030de9b25SBarry Smith 
20130de9b25SBarry Smith    Then, your solver can be chosen with the procedural interface via
20230de9b25SBarry Smith $     MatSetType(Mat,"my_mat")
20330de9b25SBarry Smith    or at runtime via the option
20430de9b25SBarry Smith $     -mat_type my_mat
20530de9b25SBarry Smith 
20630de9b25SBarry Smith    Level: advanced
20730de9b25SBarry Smith 
208ab901514SBarry Smith    Notes: ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values.
20930de9b25SBarry Smith          If your function is not being put into a shared library then use VecRegister() instead
21030de9b25SBarry Smith 
21130de9b25SBarry Smith .keywords: Mat, register
21230de9b25SBarry Smith 
21330de9b25SBarry Smith .seealso: MatRegisterAll(), MatRegisterDestroy()
21430de9b25SBarry Smith 
21530de9b25SBarry Smith M*/
216273d9f13SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
217273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,0)
218273d9f13SBarry Smith #else
219273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,d)
22030de9b25SBarry Smith #endif
22130de9b25SBarry Smith 
222ace3abfcSBarry Smith extern PetscBool  MatRegisterAllCalled;
223b0a32e0cSBarry Smith extern PetscFList MatList;
224b022a5c1SBarry Smith extern PetscFList MatColoringList;
225b022a5c1SBarry Smith extern PetscFList MatPartitioningList;
22628988994SBarry Smith 
2273b224e63SBarry Smith /*E
2283b224e63SBarry Smith     MatStructure - Indicates if the matrix has the same nonzero structure
2293b224e63SBarry Smith 
2303b224e63SBarry Smith     Level: beginner
2313b224e63SBarry Smith 
2323b224e63SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
2333b224e63SBarry Smith 
2343b224e63SBarry Smith .seealso: MatCopy(), KSPSetOperators(), PCSetOperators()
2353b224e63SBarry Smith E*/
2363b224e63SBarry Smith typedef enum {SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER,SUBSET_NONZERO_PATTERN} MatStructure;
2373b224e63SBarry Smith 
2387087cfbeSBarry Smith extern PetscErrorCode  MatCreateSeqDense(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*);
2397087cfbeSBarry Smith extern PetscErrorCode  MatCreateMPIDense(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*);
2407087cfbeSBarry Smith extern PetscErrorCode  MatCreateSeqAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
241ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,nz,nnz,&A),Mat,A)
242ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,&A),Mat,A)
243ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,0,nnz,&A),Mat,A)
244ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,&A),Mat,A)
245ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,A))
246ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,m,n,0,nnz,A))
247ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,A))
2487087cfbeSBarry Smith extern PetscErrorCode  MatCreateMPIAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
249ba966639SSatish 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)
250ba966639SSatish 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)
251ba966639SSatish 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)
252ba966639SSatish 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)
253ba966639SSatish 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))
254ba966639SSatish 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))
255ba966639SSatish 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))
256ba966639SSatish 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)
257ba966639SSatish 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)
258ba966639SSatish 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)
259ba966639SSatish 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)
260ba966639SSatish 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))
261ba966639SSatish 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))
262ba966639SSatish 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))
2637087cfbeSBarry Smith extern PetscErrorCode  MatCreateMPIAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *);
2647087cfbeSBarry Smith extern PetscErrorCode  MatCreateMPIAIJWithSplitArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],PetscInt[],PetscInt[],PetscScalar[],Mat*);
2658d7a6e47SBarry Smith 
2667087cfbeSBarry Smith extern PetscErrorCode  MatCreateSeqBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
267ba966639SSatish 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)
268ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A)
269ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A)
270ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A)
271ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A))
272ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A))
273ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A))
2747087cfbeSBarry Smith extern PetscErrorCode  MatCreateMPIBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
275ba966639SSatish 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)
276ba966639SSatish 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)
277ba966639SSatish 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)
278ba966639SSatish 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)
279ba966639SSatish 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))
280ba966639SSatish 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))
281ba966639SSatish 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))
282ba966639SSatish 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)
283ba966639SSatish 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)
284ba966639SSatish 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)
285ba966639SSatish 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)
286ba966639SSatish 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))
287ba966639SSatish 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))
288ba966639SSatish 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))
2897087cfbeSBarry Smith extern PetscErrorCode  MatCreateMPIBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat*);
290d21a29f3SJed Brown 
2917087cfbeSBarry Smith extern PetscErrorCode  MatCreateMPIAdj(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscInt[],Mat*);
2927087cfbeSBarry Smith extern PetscErrorCode  MatCreateSeqSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
293ba966639SSatish 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)
294ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A)
295ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A)
296ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A)
297ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A))
298ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A))
299ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A))
300d21a29f3SJed Brown 
3017087cfbeSBarry Smith extern PetscErrorCode  MatCreateMPISBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
302ba966639SSatish 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)
303ba966639SSatish 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)
304ba966639SSatish 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)
305ba966639SSatish 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)
306ba966639SSatish 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))
307ba966639SSatish 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))
308ba966639SSatish 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))
309ba966639SSatish 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)
310ba966639SSatish 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)
311ba966639SSatish 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)
312ba966639SSatish 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)
313ba966639SSatish 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))
314ba966639SSatish 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))
315ba966639SSatish 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))
3167087cfbeSBarry Smith extern PetscErrorCode  MatCreateMPISBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *);
3177087cfbeSBarry Smith extern PetscErrorCode  MatMPISBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
318dfb205c3SBarry Smith 
3197087cfbeSBarry Smith extern PetscErrorCode  MatCreateShell(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,void *,Mat*);
320ba966639SSatish 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)
321ba966639SSatish Balay PetscPolymorphicFunction(MatCreateShell,(PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx),(PETSC_COMM_WORLD,m,n,M,N,ctx,&A),Mat,A)
3227087cfbeSBarry Smith extern PetscErrorCode  MatCreateAdic(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,void (*)(void),Mat*);
3237087cfbeSBarry Smith extern PetscErrorCode  MatCreateNormal(Mat,Mat*);
324ba966639SSatish Balay PetscPolymorphicFunction(MatCreateNormal,(Mat mat),(mat,&A),Mat,A)
3257087cfbeSBarry Smith extern PetscErrorCode  MatCreateLRC(Mat,Mat,Mat,Mat*);
3267087cfbeSBarry Smith extern PetscErrorCode  MatCreateIS(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,ISLocalToGlobalMapping,Mat*);
3277087cfbeSBarry Smith extern PetscErrorCode  MatCreateSeqAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
3287087cfbeSBarry Smith extern PetscErrorCode  MatCreateMPIAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
329*c0cdd4a1SDahai Guo 
330*c0cdd4a1SDahai Guo extern PetscErrorCode  MatCreateSeqBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
331*c0cdd4a1SDahai Guo extern PetscErrorCode  MatCreateMPIBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
332*c0cdd4a1SDahai Guo 
3337087cfbeSBarry Smith extern PetscErrorCode  MatCreateScatter(MPI_Comm,VecScatter,Mat*);
3347087cfbeSBarry Smith extern PetscErrorCode  MatScatterSetVecScatter(Mat,VecScatter);
3357087cfbeSBarry Smith extern PetscErrorCode  MatScatterGetVecScatter(Mat,VecScatter*);
3367087cfbeSBarry Smith extern PetscErrorCode  MatCreateBlockMat(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,Mat*);
3377087cfbeSBarry Smith extern PetscErrorCode  MatCompositeAddMat(Mat,Mat);
3387087cfbeSBarry Smith extern PetscErrorCode  MatCompositeMerge(Mat);
3397087cfbeSBarry Smith extern PetscErrorCode  MatCreateComposite(MPI_Comm,PetscInt,const Mat*,Mat*);
3406d7c1e57SBarry Smith typedef enum {MAT_COMPOSITE_ADDITIVE,MAT_COMPOSITE_MULTIPLICATIVE} MatCompositeType;
3417087cfbeSBarry Smith extern PetscErrorCode  MatCompositeSetType(Mat,MatCompositeType);
3426d7c1e57SBarry Smith 
343dedccee8SHong Zhang extern PetscErrorCode  MatCreateFFT(MPI_Comm,PetscInt,const PetscInt[],const MatType,Mat*);
3447087cfbeSBarry Smith extern PetscErrorCode  MatCreateSeqCUFFT(MPI_Comm,PetscInt,const PetscInt[],Mat*);
345dedccee8SHong Zhang 
3467087cfbeSBarry Smith extern PetscErrorCode  MatCreateTranspose(Mat,Mat*);
3477087cfbeSBarry Smith extern PetscErrorCode  MatCreateSubMatrix(Mat,IS,IS,Mat*);
3487087cfbeSBarry Smith extern PetscErrorCode  MatSubMatrixUpdate(Mat,Mat,IS,IS);
3497087cfbeSBarry Smith extern PetscErrorCode  MatCreateLocalRef(Mat,IS,IS,Mat*);
3501472f72bSBarry Smith 
3517087cfbeSBarry Smith extern PetscErrorCode  MatPythonSetType(Mat,const char[]);
3521d6018f0SLisandro Dalcin 
3537087cfbeSBarry Smith extern PetscErrorCode  MatSetUp(Mat);
354fcfd50ebSBarry Smith extern PetscErrorCode  MatDestroy(Mat*);
35521c89e3eSBarry Smith 
3567087cfbeSBarry Smith extern PetscErrorCode  MatConjugate(Mat);
3577087cfbeSBarry Smith extern PetscErrorCode  MatRealPart(Mat);
3587087cfbeSBarry Smith extern PetscErrorCode  MatImaginaryPart(Mat);
3597087cfbeSBarry Smith extern PetscErrorCode  MatGetDiagonalBlock(Mat,PetscBool *,MatReuse,Mat*);
3607087cfbeSBarry Smith extern PetscErrorCode  MatGetTrace(Mat,PetscScalar*);
36199cafbc1SBarry Smith 
3628ed539a5SBarry Smith /* ------------------------------------------------------------*/
3637087cfbeSBarry Smith extern PetscErrorCode  MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
3647087cfbeSBarry Smith extern PetscErrorCode  MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
3657087cfbeSBarry Smith extern PetscErrorCode  MatSetValuesRow(Mat,PetscInt,const PetscScalar[]);
3667087cfbeSBarry Smith extern PetscErrorCode  MatSetValuesRowLocal(Mat,PetscInt,const PetscScalar[]);
36784cb2905SBarry Smith 
3682ef4de8bSBarry Smith /*S
3692ef4de8bSBarry Smith      MatStencil - Data structure (C struct) for storing information about a single row or
3702ef4de8bSBarry Smith         column of a matrix as index on an associated grid.
3712ef4de8bSBarry Smith 
3722ef4de8bSBarry Smith    Level: beginner
3732ef4de8bSBarry Smith 
3742ef4de8bSBarry Smith   Concepts: matrix; linear operator
3752ef4de8bSBarry Smith 
376d7bd93baSBarry Smith .seealso:  MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockStencil()
3772ef4de8bSBarry Smith S*/
378435da068SBarry Smith typedef struct {
379c1ac3661SBarry Smith   PetscInt k,j,i,c;
380435da068SBarry Smith } MatStencil;
3812ef4de8bSBarry Smith 
3827087cfbeSBarry Smith extern PetscErrorCode  MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
3837087cfbeSBarry Smith extern PetscErrorCode  MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
3847087cfbeSBarry Smith extern PetscErrorCode  MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt);
385435da068SBarry Smith 
3867087cfbeSBarry Smith extern PetscErrorCode  MatSetColoring(Mat,ISColoring);
3877087cfbeSBarry Smith extern PetscErrorCode  MatSetValuesAdic(Mat,void*);
3887087cfbeSBarry Smith extern PetscErrorCode  MatSetValuesAdifor(Mat,PetscInt,void*);
3893a7fca6bSBarry Smith 
390d91e6319SBarry Smith /*E
391d91e6319SBarry Smith     MatAssemblyType - Indicates if the matrix is now to be used, or if you plan
392d91e6319SBarry Smith      to continue to add values to it
393d91e6319SBarry Smith 
394d91e6319SBarry Smith     Level: beginner
395d91e6319SBarry Smith 
396d91e6319SBarry Smith .seealso: MatAssemblyBegin(), MatAssemblyEnd()
397d91e6319SBarry Smith E*/
3986d4a8577SBarry Smith typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType;
3997087cfbeSBarry Smith extern PetscErrorCode  MatAssemblyBegin(Mat,MatAssemblyType);
4007087cfbeSBarry Smith extern PetscErrorCode  MatAssemblyEnd(Mat,MatAssemblyType);
4017087cfbeSBarry Smith extern PetscErrorCode  MatAssembled(Mat,PetscBool *);
4024f9c727eSBarry Smith 
4031d73ed98SBarry Smith 
40430de9b25SBarry Smith 
405d91e6319SBarry Smith /*E
406d91e6319SBarry Smith     MatOption - Options that may be set for a matrix and its behavior or storage
407d91e6319SBarry Smith 
408d91e6319SBarry Smith     Level: beginner
409d91e6319SBarry Smith 
4100a835dfdSSatish Balay    Any additions/changes here MUST also be made in include/finclude/petscmat.h
411d91e6319SBarry Smith 
412d91e6319SBarry Smith .seealso: MatSetOption()
413d91e6319SBarry Smith E*/
414512a5fc5SBarry Smith typedef enum {MAT_ROW_ORIENTED,MAT_NEW_NONZERO_LOCATIONS,
4154e0d8c25SBarry Smith               MAT_SYMMETRIC,
4164e0d8c25SBarry Smith               MAT_STRUCTURALLY_SYMMETRIC,
4178047e77bSBarry Smith               MAT_NEW_DIAGONALS,MAT_IGNORE_OFF_PROC_ENTRIES,
4184e0d8c25SBarry Smith               MAT_NEW_NONZERO_LOCATION_ERR,
4194e0d8c25SBarry Smith               MAT_NEW_NONZERO_ALLOCATION_ERR,MAT_USE_HASH_TABLE,
420a9817697SBarry Smith               MAT_KEEP_NONZERO_PATTERN,MAT_IGNORE_ZERO_ENTRIES,
4214e0d8c25SBarry Smith               MAT_USE_INODES,
4224e0d8c25SBarry Smith               MAT_HERMITIAN,
4234e0d8c25SBarry Smith               MAT_SYMMETRY_ETERNAL,
424cd6b891eSBarry Smith               MAT_CHECK_COMPRESSED_ROW,
4254e0d8c25SBarry Smith               MAT_IGNORE_LOWER_TRIANGULAR,MAT_ERROR_LOWER_TRIANGULAR,
42628b2fa4aSMatthew Knepley               MAT_GETROW_UPPERTRIANGULAR,MAT_UNUSED_NONZERO_LOCATION_ERR,
4274cb17eb5SBarry Smith               MAT_SPD,MAT_NO_OFF_PROC_ENTRIES,MAT_NO_OFF_PROC_ZERO_ROWS,
42828b2fa4aSMatthew Knepley               NUM_MAT_OPTIONS} MatOption;
429290bbb0aSBarry Smith extern const char *MatOptions[];
4307087cfbeSBarry Smith extern PetscErrorCode  MatSetOption(Mat,MatOption,PetscBool );
4317087cfbeSBarry Smith extern PetscErrorCode  MatGetType(Mat,const MatType*);
432a313700dSBarry Smith PetscPolymorphicFunction(MatGetType,(Mat mat),(mat,&t),const MatType,t)
43384cb2905SBarry Smith 
4347087cfbeSBarry Smith extern PetscErrorCode  MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]);
4357087cfbeSBarry Smith extern PetscErrorCode  MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
4367087cfbeSBarry Smith extern PetscErrorCode  MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
4377087cfbeSBarry Smith extern PetscErrorCode  MatGetRowUpperTriangular(Mat);
4387087cfbeSBarry Smith extern PetscErrorCode  MatRestoreRowUpperTriangular(Mat);
4397087cfbeSBarry Smith extern PetscErrorCode  MatGetColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
4407087cfbeSBarry Smith extern PetscErrorCode  MatRestoreColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
4417087cfbeSBarry Smith extern PetscErrorCode  MatGetColumnVector(Mat,Vec,PetscInt);
4427087cfbeSBarry Smith extern PetscErrorCode  MatGetArray(Mat,PetscScalar *[]);
443ba966639SSatish Balay PetscPolymorphicFunction(MatGetArray,(Mat mat),(mat,&a),PetscScalar*,a)
4447087cfbeSBarry Smith extern PetscErrorCode  MatRestoreArray(Mat,PetscScalar *[]);
4457087cfbeSBarry Smith extern PetscErrorCode  MatGetBlockSize(Mat,PetscInt *);
446ba966639SSatish Balay PetscPolymorphicFunction(MatGetBlockSize,(Mat mat),(mat,&a),PetscInt,a)
4477087cfbeSBarry Smith extern PetscErrorCode  MatSetBlockSize(Mat,PetscInt);
4487b80b807SBarry Smith 
4491620fd73SBarry Smith 
4507087cfbeSBarry Smith extern PetscErrorCode  MatMult(Mat,Vec,Vec);
4517087cfbeSBarry Smith extern PetscErrorCode  MatMultDiagonalBlock(Mat,Vec,Vec);
4527087cfbeSBarry Smith extern PetscErrorCode  MatMultAdd(Mat,Vec,Vec,Vec);
453ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
4547087cfbeSBarry Smith extern PetscErrorCode  MatMultTranspose(Mat,Vec,Vec);
4557087cfbeSBarry Smith extern PetscErrorCode  MatMultHermitianTranspose(Mat,Vec,Vec);
4567087cfbeSBarry Smith extern PetscErrorCode  MatIsTranspose(Mat,Mat,PetscReal,PetscBool *);
457ace3abfcSBarry Smith PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B,PetscReal tol),(A,B,tol,&t),PetscBool ,t)
458ace3abfcSBarry Smith PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B),(A,B,0,&t),PetscBool ,t)
4597087cfbeSBarry Smith extern PetscErrorCode  MatIsHermitianTranspose(Mat,Mat,PetscReal,PetscBool *);
4607087cfbeSBarry Smith extern PetscErrorCode  MatMultTransposeAdd(Mat,Vec,Vec,Vec);
4617087cfbeSBarry Smith extern PetscErrorCode  MatMultHermitianTransposeAdd(Mat,Vec,Vec,Vec);
462ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultTransposeAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
4637087cfbeSBarry Smith extern PetscErrorCode  MatMultConstrained(Mat,Vec,Vec);
4647087cfbeSBarry Smith extern PetscErrorCode  MatMultTransposeConstrained(Mat,Vec,Vec);
4657087cfbeSBarry Smith extern PetscErrorCode  MatMatSolve(Mat,Mat,Mat);
466f5edf698SHong Zhang 
467d91e6319SBarry Smith /*E
468d91e6319SBarry Smith     MatDuplicateOption - Indicates if a duplicated sparse matrix should have
469d91e6319SBarry Smith   its numerical values copied over or just its nonzero structure.
470d91e6319SBarry Smith 
471d91e6319SBarry Smith     Level: beginner
472d91e6319SBarry Smith 
473d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
474d91e6319SBarry Smith 
47570dcbbb9SBarry Smith $   MAT_SHARE_NONZERO_PATTERN - the i and j arrays in the new matrix will be shared with the original matrix
47670dcbbb9SBarry Smith $                               this also triggers the MAT_DO_NOT_COPY_VALUES option. This is used when you
47770dcbbb9SBarry Smith $                               have several matrices with the same nonzero pattern.
47870dcbbb9SBarry Smith 
479d91e6319SBarry Smith .seealso: MatDuplicate()
480d91e6319SBarry Smith E*/
48170dcbbb9SBarry Smith typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES,MAT_SHARE_NONZERO_PATTERN} MatDuplicateOption;
4822e8a6d31SBarry Smith 
4837087cfbeSBarry Smith extern PetscErrorCode  MatConvert(Mat,const MatType,MatReuse,Mat*);
484a313700dSBarry Smith PetscPolymorphicFunction(MatConvert,(Mat A,const MatType t),(A,t,MAT_INITIAL_MATRIX,&a),Mat,a)
4857087cfbeSBarry Smith extern PetscErrorCode  MatDuplicate(Mat,MatDuplicateOption,Mat*);
486ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A,MatDuplicateOption o),(A,o,&a),Mat,a)
487ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A),(A,MAT_COPY_VALUES,&a),Mat,a)
48894a9d846SBarry Smith 
489cb5b572fSBarry Smith 
4907087cfbeSBarry Smith extern PetscErrorCode  MatCopy(Mat,Mat,MatStructure);
4917087cfbeSBarry Smith extern PetscErrorCode  MatView(Mat,PetscViewer);
4927087cfbeSBarry Smith extern PetscErrorCode  MatIsSymmetric(Mat,PetscReal,PetscBool *);
493ace3abfcSBarry Smith PetscPolymorphicFunction(MatIsSymmetric,(Mat A,PetscReal tol),(A,tol,&t),PetscBool ,t)
494ace3abfcSBarry Smith PetscPolymorphicFunction(MatIsSymmetric,(Mat A),(A,0,&t),PetscBool ,t)
4957087cfbeSBarry Smith extern PetscErrorCode  MatIsStructurallySymmetric(Mat,PetscBool *);
496ace3abfcSBarry Smith PetscPolymorphicFunction(MatIsStructurallySymmetric,(Mat A),(A,&t),PetscBool ,t)
4977087cfbeSBarry Smith extern PetscErrorCode  MatIsHermitian(Mat,PetscReal,PetscBool *);
498ace3abfcSBarry Smith PetscPolymorphicFunction(MatIsHermitian,(Mat A),(A,0,&t),PetscBool ,t)
4997087cfbeSBarry Smith extern PetscErrorCode  MatIsSymmetricKnown(Mat,PetscBool *,PetscBool *);
5007087cfbeSBarry Smith extern PetscErrorCode  MatIsHermitianKnown(Mat,PetscBool *,PetscBool *);
5017087cfbeSBarry Smith extern PetscErrorCode  MatMissingDiagonal(Mat,PetscBool  *,PetscInt *);
5027087cfbeSBarry Smith extern PetscErrorCode  MatLoad(Mat, PetscViewer);
5037b80b807SBarry Smith 
5047087cfbeSBarry Smith extern PetscErrorCode  MatGetRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,PetscInt *[],PetscInt *[],PetscBool  *);
5057087cfbeSBarry Smith extern PetscErrorCode  MatRestoreRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,PetscInt *[],PetscInt *[],PetscBool  *);
5067087cfbeSBarry Smith extern PetscErrorCode  MatGetColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,PetscInt *[],PetscInt *[],PetscBool  *);
5077087cfbeSBarry Smith extern PetscErrorCode  MatRestoreColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,PetscInt *[],PetscInt *[],PetscBool  *);
508d4fbbf0eSBarry Smith 
509d91e6319SBarry Smith /*S
510d91e6319SBarry Smith      MatInfo - Context of matrix information, used with MatGetInfo()
511d91e6319SBarry Smith 
512d91e6319SBarry Smith    In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE
513d91e6319SBarry Smith 
514d91e6319SBarry Smith    Level: intermediate
515d91e6319SBarry Smith 
516d91e6319SBarry Smith   Concepts: matrix^nonzero information
517d91e6319SBarry Smith 
518d9274352SBarry Smith .seealso:  MatGetInfo(), MatInfoType
519d91e6319SBarry Smith S*/
5204e220ebcSLois Curfman McInnes typedef struct {
521b0a32e0cSBarry Smith   PetscLogDouble block_size;                         /* block size */
522b0a32e0cSBarry Smith   PetscLogDouble nz_allocated,nz_used,nz_unneeded;   /* number of nonzeros */
523b0a32e0cSBarry Smith   PetscLogDouble memory;                             /* memory allocated */
524b0a32e0cSBarry Smith   PetscLogDouble assemblies;                         /* number of matrix assemblies called */
525b0a32e0cSBarry Smith   PetscLogDouble mallocs;                            /* number of mallocs during MatSetValues() */
526b0a32e0cSBarry Smith   PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */
527b0a32e0cSBarry Smith   PetscLogDouble factor_mallocs;                     /* number of mallocs during factorization */
5284e220ebcSLois Curfman McInnes } MatInfo;
5294e220ebcSLois Curfman McInnes 
530d9274352SBarry Smith /*E
531d9274352SBarry Smith     MatInfoType - Indicates if you want information about the local part of the matrix,
532d9274352SBarry Smith      the entire parallel matrix or the maximum over all the local parts.
533d9274352SBarry Smith 
534d9274352SBarry Smith     Level: beginner
535d9274352SBarry Smith 
536d9274352SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
537d9274352SBarry Smith 
538d9274352SBarry Smith .seealso: MatGetInfo(), MatInfo
539d9274352SBarry Smith E*/
5407b80b807SBarry Smith typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType;
5417087cfbeSBarry Smith extern PetscErrorCode  MatGetInfo(Mat,MatInfoType,MatInfo*);
5427087cfbeSBarry Smith extern PetscErrorCode  MatGetDiagonal(Mat,Vec);
5437087cfbeSBarry Smith extern PetscErrorCode  MatGetRowMax(Mat,Vec,PetscInt[]);
5447087cfbeSBarry Smith extern PetscErrorCode  MatGetRowMin(Mat,Vec,PetscInt[]);
5457087cfbeSBarry Smith extern PetscErrorCode  MatGetRowMaxAbs(Mat,Vec,PetscInt[]);
5467087cfbeSBarry Smith extern PetscErrorCode  MatGetRowMinAbs(Mat,Vec,PetscInt[]);
5477087cfbeSBarry Smith extern PetscErrorCode  MatGetRowSum(Mat,Vec);
5487087cfbeSBarry Smith extern PetscErrorCode  MatTranspose(Mat,MatReuse,Mat*);
549fc4dec0aSBarry Smith PetscPolymorphicFunction(MatTranspose,(Mat A),(A,MAT_INITIAL_MATRIX,&t),Mat,t)
5507087cfbeSBarry Smith extern PetscErrorCode  MatHermitianTranspose(Mat,MatReuse,Mat*);
5517087cfbeSBarry Smith extern PetscErrorCode  MatPermute(Mat,IS,IS,Mat *);
552ba966639SSatish Balay PetscPolymorphicFunction(MatPermute,(Mat A,IS is1,IS is2),(A,is1,is2,&t),Mat,t)
5537087cfbeSBarry Smith extern PetscErrorCode  MatDiagonalScale(Mat,Vec,Vec);
5547087cfbeSBarry Smith extern PetscErrorCode  MatDiagonalSet(Mat,Vec,InsertMode);
5557087cfbeSBarry Smith extern PetscErrorCode  MatEqual(Mat,Mat,PetscBool *);
556ace3abfcSBarry Smith PetscPolymorphicFunction(MatEqual,(Mat A,Mat B),(A,B,&t),PetscBool ,t)
5577087cfbeSBarry Smith extern PetscErrorCode  MatMultEqual(Mat,Mat,PetscInt,PetscBool *);
5587087cfbeSBarry Smith extern PetscErrorCode  MatMultAddEqual(Mat,Mat,PetscInt,PetscBool *);
5597087cfbeSBarry Smith extern PetscErrorCode  MatMultTransposeEqual(Mat,Mat,PetscInt,PetscBool *);
5607087cfbeSBarry Smith extern PetscErrorCode  MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscBool *);
5617b80b807SBarry Smith 
5627087cfbeSBarry Smith extern PetscErrorCode  MatNorm(Mat,NormType,PetscReal *);
563ba966639SSatish Balay PetscPolymorphicFunction(MatNorm,(Mat A,NormType t),(A,t,&n),PetscReal,n)
56409573ac7SBarry Smith extern PetscErrorCode MatGetColumnNorms(Mat,NormType,PetscReal *);
5657087cfbeSBarry Smith extern PetscErrorCode  MatZeroEntries(Mat);
5667087cfbeSBarry Smith extern PetscErrorCode  MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
5677087cfbeSBarry Smith extern PetscErrorCode  MatZeroRowsIS(Mat,IS,PetscScalar,Vec,Vec);
5687087cfbeSBarry Smith extern PetscErrorCode  MatZeroRowsStencil(Mat,PetscInt,const MatStencil [],PetscScalar,Vec,Vec);
5697087cfbeSBarry Smith extern PetscErrorCode  MatZeroRowsColumns(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
5707087cfbeSBarry Smith extern PetscErrorCode  MatZeroRowsColumnsIS(Mat,IS,PetscScalar,Vec,Vec);
5717b80b807SBarry Smith 
5727087cfbeSBarry Smith extern PetscErrorCode  MatUseScaledForm(Mat,PetscBool );
5737087cfbeSBarry Smith extern PetscErrorCode  MatScaleSystem(Mat,Vec,Vec);
5747087cfbeSBarry Smith extern PetscErrorCode  MatUnScaleSystem(Mat,Vec,Vec);
5755ef9f2a5SBarry Smith 
5767087cfbeSBarry Smith extern PetscErrorCode  MatGetSize(Mat,PetscInt*,PetscInt*);
5777087cfbeSBarry Smith extern PetscErrorCode  MatGetLocalSize(Mat,PetscInt*,PetscInt*);
5787087cfbeSBarry Smith extern PetscErrorCode  MatGetOwnershipRange(Mat,PetscInt*,PetscInt*);
5797087cfbeSBarry Smith extern PetscErrorCode  MatGetOwnershipRanges(Mat,const PetscInt**);
5807087cfbeSBarry Smith extern PetscErrorCode  MatGetOwnershipRangeColumn(Mat,PetscInt*,PetscInt*);
5817087cfbeSBarry Smith extern PetscErrorCode  MatGetOwnershipRangesColumn(Mat,const PetscInt**);
5827b80b807SBarry Smith 
5837087cfbeSBarry Smith extern PetscErrorCode  MatGetSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
5847087cfbeSBarry Smith extern PetscErrorCode  MatGetSubMatricesParallel(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
5857087cfbeSBarry Smith extern PetscErrorCode  MatDestroyMatrices(PetscInt,Mat *[]);
5867087cfbeSBarry Smith extern PetscErrorCode  MatGetSubMatrix(Mat,IS,IS,MatReuse,Mat *);
5877087cfbeSBarry Smith extern PetscErrorCode  MatGetLocalSubMatrix(Mat,IS,IS,Mat*);
5887087cfbeSBarry Smith extern PetscErrorCode  MatRestoreLocalSubMatrix(Mat,IS,IS,Mat*);
5897087cfbeSBarry Smith extern PetscErrorCode  MatGetSeqNonzeroStructure(Mat,Mat*);
5907087cfbeSBarry Smith extern PetscErrorCode  MatDestroySeqNonzeroStructure(Mat*);
5915494a064SHong Zhang 
5927087cfbeSBarry Smith extern PetscErrorCode  MatMerge(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
5937087cfbeSBarry Smith extern PetscErrorCode  MatMerge_SeqsToMPI(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*);
5947087cfbeSBarry Smith extern PetscErrorCode  MatMerge_SeqsToMPISymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*);
5957087cfbeSBarry Smith extern PetscErrorCode  MatMerge_SeqsToMPINumeric(Mat,Mat);
5964a2b5492SBarry Smith extern PetscErrorCode  MatMPIAIJGetLocalMat(Mat,MatReuse,Mat*);
5974a2b5492SBarry Smith extern PetscErrorCode  MatMPIAIJGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*);
5987087cfbeSBarry Smith extern PetscErrorCode  MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,PetscInt*,Mat*);
5997087cfbeSBarry Smith extern PetscErrorCode  MatGetBrowsOfAoCols(Mat,Mat,MatReuse,PetscInt**,PetscInt**,MatScalar**,Mat*);
60043eb5e2fSMatthew Knepley #if defined (PETSC_USE_CTABLE)
601cd116e26SSatish Balay #include "petscctable.h"
6027087cfbeSBarry Smith extern PetscErrorCode  MatGetCommunicationStructs(Mat, Vec *, PetscTable *, VecScatter *);
60343eb5e2fSMatthew Knepley #else
6047087cfbeSBarry Smith extern PetscErrorCode  MatGetCommunicationStructs(Mat, Vec *, PetscInt *[], VecScatter *);
60543eb5e2fSMatthew Knepley #endif
6067087cfbeSBarry Smith extern PetscErrorCode  MatGetGhosts(Mat, PetscInt *,const PetscInt *[]);
6078efafbd8SBarry Smith 
6087087cfbeSBarry Smith extern PetscErrorCode  MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt);
6097b80b807SBarry Smith 
6107087cfbeSBarry Smith extern PetscErrorCode  MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*);
6117087cfbeSBarry Smith extern PetscErrorCode  MatMatMultSymbolic(Mat,Mat,PetscReal,Mat*);
6127087cfbeSBarry Smith extern PetscErrorCode  MatMatMultNumeric(Mat,Mat,Mat);
61322440eb1SKris Buschelman 
6147087cfbeSBarry Smith extern PetscErrorCode  MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*);
6157087cfbeSBarry Smith extern PetscErrorCode  MatPtAPSymbolic(Mat,Mat,PetscReal,Mat*);
6167087cfbeSBarry Smith extern PetscErrorCode  MatPtAPNumeric(Mat,Mat,Mat);
61722440eb1SKris Buschelman 
6187087cfbeSBarry Smith extern PetscErrorCode  MatMatMultTranspose(Mat,Mat,MatReuse,PetscReal,Mat*);
6197087cfbeSBarry Smith extern PetscErrorCode  MatMatMultTransposeSymbolic(Mat,Mat,PetscReal,Mat*);
6207087cfbeSBarry Smith extern PetscErrorCode  MatMatMultTransposeNumeric(Mat,Mat,Mat);
621bc011b1eSHong Zhang 
6227087cfbeSBarry Smith extern PetscErrorCode  MatAXPY(Mat,PetscScalar,Mat,MatStructure);
6237087cfbeSBarry Smith extern PetscErrorCode  MatAYPX(Mat,PetscScalar,Mat,MatStructure);
6241c741599SHong Zhang 
6257087cfbeSBarry Smith extern PetscErrorCode  MatScale(Mat,PetscScalar);
6267087cfbeSBarry Smith extern PetscErrorCode  MatShift(Mat,PetscScalar);
6277b80b807SBarry Smith 
6287087cfbeSBarry Smith extern PetscErrorCode  MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping);
6297087cfbeSBarry Smith extern PetscErrorCode  MatSetLocalToGlobalMappingBlock(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping);
6307087cfbeSBarry Smith extern PetscErrorCode  MatGetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*);
6317087cfbeSBarry Smith extern PetscErrorCode  MatGetLocalToGlobalMappingBlock(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*);
6327087cfbeSBarry Smith extern PetscErrorCode  MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
6337087cfbeSBarry Smith extern PetscErrorCode  MatZeroRowsLocalIS(Mat,IS,PetscScalar,Vec,Vec);
6347087cfbeSBarry Smith extern PetscErrorCode  MatZeroRowsColumnsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
6357087cfbeSBarry Smith extern PetscErrorCode  MatZeroRowsColumnsLocalIS(Mat,IS,PetscScalar,Vec,Vec);
6367087cfbeSBarry Smith extern PetscErrorCode  MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
6377087cfbeSBarry Smith extern PetscErrorCode  MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
638052efed2SBarry Smith 
6397087cfbeSBarry Smith extern PetscErrorCode  MatStashSetInitialSize(Mat,PetscInt,PetscInt);
6407087cfbeSBarry Smith extern PetscErrorCode  MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
64190f02eecSBarry Smith 
6427087cfbeSBarry Smith extern PetscErrorCode  MatInterpolate(Mat,Vec,Vec);
6437087cfbeSBarry Smith extern PetscErrorCode  MatInterpolateAdd(Mat,Vec,Vec,Vec);
644ba966639SSatish Balay PetscPolymorphicSubroutine(MatInterpolateAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
6457087cfbeSBarry Smith extern PetscErrorCode  MatRestrict(Mat,Vec,Vec);
6467087cfbeSBarry Smith extern PetscErrorCode  MatGetVecs(Mat,Vec*,Vec*);
6477087cfbeSBarry Smith extern PetscErrorCode  MatGetRedundantMatrix(Mat,PetscInt,MPI_Comm,PetscInt,MatReuse,Mat*);
6487087cfbeSBarry Smith extern PetscErrorCode  MatGetMultiProcBlock(Mat,MPI_Comm,Mat*);
649fc9bc008SSatish Balay extern PetscErrorCode  MatFindZeroDiagonals(Mat,IS*);
650bd481603SBarry Smith 
651bd481603SBarry Smith /*MC
6521620fd73SBarry Smith    MatSetValue - Set a single entry into a matrix.
6531620fd73SBarry Smith 
6541620fd73SBarry Smith    Not collective
6551620fd73SBarry Smith 
6561620fd73SBarry Smith    Input Parameters:
6571620fd73SBarry Smith +  m - the matrix
6581620fd73SBarry Smith .  row - the row location of the entry
6591620fd73SBarry Smith .  col - the column location of the entry
6601620fd73SBarry Smith .  value - the value to insert
6611620fd73SBarry Smith -  mode - either INSERT_VALUES or ADD_VALUES
6621620fd73SBarry Smith 
6631620fd73SBarry Smith    Notes:
6641620fd73SBarry Smith    For efficiency one should use MatSetValues() and set several or many
6651620fd73SBarry Smith    values simultaneously if possible.
6661620fd73SBarry Smith 
6671620fd73SBarry Smith    Level: beginner
6681620fd73SBarry Smith 
6691620fd73SBarry Smith .seealso: MatSetValues(), MatSetValueLocal()
6701620fd73SBarry Smith M*/
6711620fd73SBarry 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);}
6721620fd73SBarry Smith 
6731620fd73SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatGetValue(Mat v,PetscInt i,PetscInt j,PetscScalar *va) {return MatGetValues(v,1,&i,1,&j,va);}
6741620fd73SBarry Smith 
6751620fd73SBarry 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);}
6761620fd73SBarry Smith 
6771620fd73SBarry Smith /*MC
6780d2a7ff2SSatish Balay    MatPreallocateInitialize - Begins the block of code that will count the number of nonzeros per
679bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
680bd481603SBarry Smith 
681bd481603SBarry Smith    Synopsis:
682c1ac3661SBarry Smith    PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)
683bd481603SBarry Smith 
684bd481603SBarry Smith    Collective on MPI_Comm
685bd481603SBarry Smith 
686bd481603SBarry Smith    Input Parameters:
687bd481603SBarry Smith +  comm - the communicator that will share the eventually allocated matrix
688859fcb39SBarry Smith .  nrows - the number of LOCAL rows in the matrix
689859fcb39SBarry Smith -  ncols - the number of LOCAL columns in the matrix
690bd481603SBarry Smith 
691bd481603SBarry Smith    Output Parameters:
692bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
693bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
694bd481603SBarry Smith 
695bd481603SBarry Smith 
696bd481603SBarry Smith    Level: intermediate
697bd481603SBarry Smith 
698bd481603SBarry Smith    Notes:
6990598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
700bd481603SBarry Smith 
7011d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
702bd481603SBarry Smith 
703bd481603SBarry Smith    Use MatPreallocateInitializeSymmetric() for symmetric matrices (MPISBAIJ matrices)
704bd481603SBarry Smith 
7051620fd73SBarry Smith    This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize().
7061620fd73SBarry Smith 
707bd481603SBarry Smith   Concepts: preallocation^Matrix
708bd481603SBarry Smith 
709bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
710bd481603SBarry Smith           MatPreallocateInitializeSymmetric(), MatPreallocateSymmetricSetLocal()
711bd481603SBarry Smith M*/
712c4f061fbSSatish Balay #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \
7137c922b88SBarry Smith { \
714a89bc211SBarry Smith   PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__start,__end; \
715a89bc211SBarry Smith   _4_ierr = PetscMalloc2(__nrows,PetscInt,&dnz,__nrows,PetscInt,&onz);CHKERRQ(_4_ierr); \
716a89bc211SBarry Smith   _4_ierr = PetscMemzero(dnz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
717a89bc211SBarry Smith   _4_ierr = PetscMemzero(onz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
718c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ctmp;\
719a89bc211SBarry Smith   _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows;
7207c922b88SBarry Smith 
721bd481603SBarry Smith /*MC
7220d2a7ff2SSatish Balay    MatPreallocateSymmetricInitialize - Begins the block of code that will count the number of nonzeros per
723bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
724bd481603SBarry Smith 
725bd481603SBarry Smith    Synopsis:
726c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)
727bd481603SBarry Smith 
728bd481603SBarry Smith    Collective on MPI_Comm
729bd481603SBarry Smith 
730bd481603SBarry Smith    Input Parameters:
731bd481603SBarry Smith +  comm - the communicator that will share the eventually allocated matrix
732859fcb39SBarry Smith .  nrows - the number of LOCAL rows in the matrix
733859fcb39SBarry Smith -  ncols - the number of LOCAL columns in the matrix
734bd481603SBarry Smith 
735bd481603SBarry Smith    Output Parameters:
736bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
737bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
738bd481603SBarry Smith 
739bd481603SBarry Smith 
740bd481603SBarry Smith    Level: intermediate
741bd481603SBarry Smith 
742bd481603SBarry Smith    Notes:
7430598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
744bd481603SBarry Smith 
7451d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
746bd481603SBarry Smith 
7471620fd73SBarry Smith    This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize().
7481620fd73SBarry Smith 
749bd481603SBarry Smith   Concepts: preallocation^Matrix
750bd481603SBarry Smith 
751bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
752bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal()
753bd481603SBarry Smith M*/
754222b16d4SBarry Smith #define MatPreallocateSymmetricInitialize(comm,nrows,ncols,dnz,onz) 0; \
755222b16d4SBarry Smith { \
756a89bc211SBarry Smith   PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__end; \
757a89bc211SBarry Smith   _4_ierr = PetscMalloc2(__nrows,PetscInt,&dnz,__nrows,PetscInt,&onz);CHKERRQ(_4_ierr); \
758a89bc211SBarry Smith   _4_ierr = PetscMemzero(dnz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
759a89bc211SBarry Smith   _4_ierr = PetscMemzero(onz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
760c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr);\
761a89bc211SBarry Smith   _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows;
762222b16d4SBarry Smith 
763bd481603SBarry Smith /*MC
764bd481603SBarry Smith    MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
765bd481603SBarry Smith        inserted using a local number of the rows and columns
766bd481603SBarry Smith 
767bd481603SBarry Smith    Synopsis:
768c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
769bd481603SBarry Smith 
770bd481603SBarry Smith    Not Collective
771bd481603SBarry Smith 
772bd481603SBarry Smith    Input Parameters:
773784ac674SJed Brown +  map - the row mapping from local numbering to global numbering
774bd481603SBarry Smith .  nrows - the number of rows indicated
7751d73ed98SBarry Smith .  rows - the indices of the rows
776784ac674SJed Brown .  cmap - the column mapping from local to global numbering
777bd481603SBarry Smith .  ncols - the number of columns in the matrix
778bd481603SBarry Smith .  cols - the columns indicated
779bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
780bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
781bd481603SBarry Smith 
782bd481603SBarry Smith 
783bd481603SBarry Smith    Level: intermediate
784bd481603SBarry Smith 
785bd481603SBarry Smith    Notes:
7860598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
787bd481603SBarry Smith 
7881d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
789bd481603SBarry Smith 
790bd481603SBarry Smith   Concepts: preallocation^Matrix
791bd481603SBarry Smith 
792bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
793bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal()
794bd481603SBarry Smith M*/
795784ac674SJed Brown #define MatPreallocateSetLocal(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \
796c4f061fbSSatish Balay {\
797c1ac3661SBarry Smith   PetscInt __l;\
798784ac674SJed Brown   _4_ierr = ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\
799784ac674SJed Brown   _4_ierr = ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\
800c4f061fbSSatish Balay   for (__l=0;__l<nrows;__l++) {\
801ef66eb69SBarry Smith     _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
802c4f061fbSSatish Balay   }\
803c4f061fbSSatish Balay }
804c4f061fbSSatish Balay 
805bd481603SBarry Smith /*MC
806bd481603SBarry Smith    MatPreallocateSymmetricSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
807bd481603SBarry Smith        inserted using a local number of the rows and columns
808bd481603SBarry Smith 
809bd481603SBarry Smith    Synopsis:
810c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
811bd481603SBarry Smith 
812bd481603SBarry Smith    Not Collective
813bd481603SBarry Smith 
814bd481603SBarry Smith    Input Parameters:
815bd481603SBarry Smith +  map - the mapping between local numbering and global numbering
816bd481603SBarry Smith .  nrows - the number of rows indicated
8171d73ed98SBarry Smith .  rows - the indices of the rows
818bd481603SBarry Smith .  ncols - the number of columns in the matrix
819bd481603SBarry Smith .  cols - the columns indicated
820bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
821bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
822bd481603SBarry Smith 
823bd481603SBarry Smith 
824bd481603SBarry Smith    Level: intermediate
825bd481603SBarry Smith 
826bd481603SBarry Smith    Notes:
8270598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
828bd481603SBarry Smith 
829bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
830bd481603SBarry Smith 
831bd481603SBarry Smith   Concepts: preallocation^Matrix
832bd481603SBarry Smith 
833bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
834bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
835bd481603SBarry Smith M*/
836d3d32019SBarry Smith #define MatPreallocateSymmetricSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\
837d3d32019SBarry Smith {\
838c1ac3661SBarry Smith   PetscInt __l;\
839d3d32019SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\
840d3d32019SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\
841d3d32019SBarry Smith   for (__l=0;__l<nrows;__l++) {\
842d3d32019SBarry Smith     _4_ierr = MatPreallocateSymmetricSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
843d3d32019SBarry Smith   }\
844d3d32019SBarry Smith }
845d3d32019SBarry Smith 
846bd481603SBarry Smith /*MC
847bd481603SBarry Smith    MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
848bd481603SBarry Smith        inserted using a local number of the rows and columns
849bd481603SBarry Smith 
850bd481603SBarry Smith    Synopsis:
851c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
852bd481603SBarry Smith 
853bd481603SBarry Smith    Not Collective
854bd481603SBarry Smith 
855bd481603SBarry Smith    Input Parameters:
85664075487SBarry Smith +  row - the row
857bd481603SBarry Smith .  ncols - the number of columns in the matrix
858a50f8bf6SHong Zhang -  cols - the columns indicated
859a50f8bf6SHong Zhang 
860a50f8bf6SHong Zhang    Output Parameters:
861a50f8bf6SHong Zhang +  dnz - the array that will be passed to the matrix preallocation routines
862bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
863bd481603SBarry Smith 
864bd481603SBarry Smith 
865bd481603SBarry Smith    Level: intermediate
866bd481603SBarry Smith 
867bd481603SBarry Smith    Notes:
8680598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
869bd481603SBarry Smith 
870bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
871bd481603SBarry Smith 
8721620fd73SBarry Smith    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
8731620fd73SBarry Smith 
874bd481603SBarry Smith   Concepts: preallocation^Matrix
875bd481603SBarry Smith 
876bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
877bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
878bd481603SBarry Smith M*/
879c4f061fbSSatish Balay #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\
880c1ac3661SBarry Smith { PetscInt __i; \
881e32f2f54SBarry 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);\
882e32f2f54SBarry 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);\
8837c922b88SBarry Smith   for (__i=0; __i<nc; __i++) {\
88464075487SBarry Smith     if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \
8857cc688d7SBarry Smith     else dnz[row - __rstart]++;\
8867c922b88SBarry Smith   }\
8877c922b88SBarry Smith }
8887c922b88SBarry Smith 
889bd481603SBarry Smith /*MC
890bd481603SBarry Smith    MatPreallocateSymmetricSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
891bd481603SBarry Smith        inserted using a local number of the rows and columns
892bd481603SBarry Smith 
893bd481603SBarry Smith    Synopsis:
894c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
895bd481603SBarry Smith 
896bd481603SBarry Smith    Not Collective
897bd481603SBarry Smith 
898bd481603SBarry Smith    Input Parameters:
899bd481603SBarry Smith +  nrows - the number of rows indicated
9001d73ed98SBarry Smith .  rows - the indices of the rows
901bd481603SBarry Smith .  ncols - the number of columns in the matrix
902bd481603SBarry Smith .  cols - the columns indicated
903bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
904bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
905bd481603SBarry Smith 
906bd481603SBarry Smith 
907bd481603SBarry Smith    Level: intermediate
908bd481603SBarry Smith 
909bd481603SBarry Smith    Notes:
9100598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
911bd481603SBarry Smith 
912bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
913bd481603SBarry Smith 
9141620fd73SBarry Smith    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
9151620fd73SBarry Smith 
916bd481603SBarry Smith   Concepts: preallocation^Matrix
917bd481603SBarry Smith 
918bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
919bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
920bd481603SBarry Smith M*/
921d3d32019SBarry Smith #define MatPreallocateSymmetricSet(row,nc,cols,dnz,onz) 0;\
922c1ac3661SBarry Smith { PetscInt __i; \
923d3d32019SBarry Smith   for (__i=0; __i<nc; __i++) {\
924d3d32019SBarry Smith     if (cols[__i] >= __end) onz[row - __rstart]++; \
925d3d32019SBarry Smith     else if (cols[__i] >= row) dnz[row - __rstart]++;\
926d3d32019SBarry Smith   }\
927d3d32019SBarry Smith }
928d3d32019SBarry Smith 
929bd481603SBarry Smith /*MC
93016371a99SBarry Smith    MatPreallocateLocation -  An alternative to MatPreallocationSet() that puts the nonzero locations into the matrix if it exists
93116371a99SBarry Smith 
93216371a99SBarry Smith    Synopsis:
93316371a99SBarry Smith    PetscErrorCode MatPreallocateLocations(Mat A,PetscInt row,PetscInt ncols,PetscInt *cols,PetscInt *dnz,PetscInt *onz)
93416371a99SBarry Smith 
93516371a99SBarry Smith    Not Collective
93616371a99SBarry Smith 
93716371a99SBarry Smith    Input Parameters:
93816371a99SBarry Smith .  A - matrix
93916371a99SBarry Smith .  row - row where values exist (must be local to this process)
94016371a99SBarry Smith .  ncols - number of columns
94116371a99SBarry Smith .  cols - columns with nonzeros
94216371a99SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
94316371a99SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
94416371a99SBarry Smith 
94516371a99SBarry Smith 
94616371a99SBarry Smith    Level: intermediate
94716371a99SBarry Smith 
94816371a99SBarry Smith    Notes:
9490598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
95016371a99SBarry Smith 
95116371a99SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
95216371a99SBarry Smith 
95316371a99SBarry Smith    This is a MACRO not a function because it uses a bunch of variables private to the MatPreallocation.... routines.
95416371a99SBarry Smith 
95516371a99SBarry Smith   Concepts: preallocation^Matrix
95616371a99SBarry Smith 
95716371a99SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
95816371a99SBarry Smith           MatPreallocateSymmetricInitialize(), MatPreallocateSymmetricSetLocal()
95916371a99SBarry Smith M*/
96016371a99SBarry 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);}
96116371a99SBarry Smith 
96216371a99SBarry Smith 
96316371a99SBarry Smith /*MC
9640d2a7ff2SSatish Balay    MatPreallocateFinalize - Ends the block of code that will count the number of nonzeros per
965bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
966bd481603SBarry Smith 
967bd481603SBarry Smith    Synopsis:
968c1ac3661SBarry Smith    PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz)
969bd481603SBarry Smith 
970bd481603SBarry Smith    Collective on MPI_Comm
971bd481603SBarry Smith 
972bd481603SBarry Smith    Input Parameters:
97316371a99SBarry Smith +  dnz - the array that was be passed to the matrix preallocation routines
974bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
975bd481603SBarry Smith 
976bd481603SBarry Smith 
977bd481603SBarry Smith    Level: intermediate
978bd481603SBarry Smith 
979bd481603SBarry Smith    Notes:
9800598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
981bd481603SBarry Smith 
982bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
983bd481603SBarry Smith 
9841620fd73SBarry Smith    This is a MACRO not a function because it closes the { started in MatPreallocateInitialize().
9851620fd73SBarry Smith 
986bd481603SBarry Smith   Concepts: preallocation^Matrix
987bd481603SBarry Smith 
988bd481603SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
989bd481603SBarry Smith           MatPreallocateSymmetricInitialize(), MatPreallocateSymmetricSetLocal()
990bd481603SBarry Smith M*/
991a89bc211SBarry Smith #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree2(dnz,onz);CHKERRQ(_4_ierr);}
9927c922b88SBarry Smith 
993bd481603SBarry Smith 
994bd481603SBarry Smith 
9957b80b807SBarry Smith /* Routines unique to particular data structures */
9967087cfbeSBarry Smith extern PetscErrorCode  MatShellGetContext(Mat,void **);
997ba966639SSatish Balay PetscPolymorphicFunction(MatShellGetContext,(Mat A),(A,&t),void*,t)
998698d4c6aSKris Buschelman 
9997087cfbeSBarry Smith extern PetscErrorCode  MatInodeAdjustForInodes(Mat,IS*,IS*);
10007087cfbeSBarry Smith extern PetscErrorCode  MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *);
1001698d4c6aSKris Buschelman 
10027087cfbeSBarry Smith extern PetscErrorCode  MatSeqAIJSetColumnIndices(Mat,PetscInt[]);
10037087cfbeSBarry Smith extern PetscErrorCode  MatSeqBAIJSetColumnIndices(Mat,PetscInt[]);
10047087cfbeSBarry Smith extern PetscErrorCode  MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
10057087cfbeSBarry Smith extern PetscErrorCode  MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
10067087cfbeSBarry Smith extern PetscErrorCode  MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
10077b80b807SBarry Smith 
1008a96a251dSBarry Smith #define MAT_SKIP_ALLOCATION -4
1009a96a251dSBarry Smith 
10107087cfbeSBarry Smith extern PetscErrorCode  MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
1011ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz))
10127087cfbeSBarry Smith extern PetscErrorCode  MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
1013ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqSBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz))
10147087cfbeSBarry Smith extern PetscErrorCode  MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]);
1015ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqAIJSetPreallocation,(Mat A,const PetscInt nnz[]),(A,0,nnz))
1016273d9f13SBarry Smith 
10177087cfbeSBarry Smith extern PetscErrorCode  MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
1018ba966639SSatish Balay PetscPolymorphicSubroutine(MatMPIBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[],const PetscInt onz[]),(A,bs,0,nnz,0,onz))
10197087cfbeSBarry Smith extern PetscErrorCode  MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
10207087cfbeSBarry Smith extern PetscErrorCode  MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
10217087cfbeSBarry Smith extern PetscErrorCode  MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []);
10227087cfbeSBarry Smith extern PetscErrorCode  MatSeqBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
10237087cfbeSBarry Smith extern PetscErrorCode  MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]);
10247087cfbeSBarry Smith extern PetscErrorCode  MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
10257087cfbeSBarry Smith extern PetscErrorCode  MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]);
10267087cfbeSBarry Smith extern PetscErrorCode  MatMPIDenseSetPreallocation(Mat,PetscScalar[]);
10277087cfbeSBarry Smith extern PetscErrorCode  MatSeqDenseSetPreallocation(Mat,PetscScalar[]);
10287087cfbeSBarry Smith extern PetscErrorCode  MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,PetscInt*[]);
10297087cfbeSBarry Smith extern PetscErrorCode  MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,PetscInt*[]);
10307087cfbeSBarry Smith extern PetscErrorCode  MatAdicSetLocalFunction(Mat,void (*)(void));
1031273d9f13SBarry Smith 
10327087cfbeSBarry Smith extern PetscErrorCode  MatSeqDenseSetLDA(Mat,PetscInt);
10337087cfbeSBarry Smith extern PetscErrorCode  MatDenseGetLocalMatrix(Mat,Mat*);
10341b807ce4Svictorle 
10357087cfbeSBarry Smith extern PetscErrorCode  MatStoreValues(Mat);
10367087cfbeSBarry Smith extern PetscErrorCode  MatRetrieveValues(Mat);
10372e8a6d31SBarry Smith 
10387087cfbeSBarry Smith extern PetscErrorCode  MatDAADSetCtx(Mat,void*);
10393a7fca6bSBarry Smith 
1040b3a44c85SBarry Smith extern PetscErrorCode  MatFindNonzeroRows(Mat,IS*);
10417b80b807SBarry Smith /*
10427b80b807SBarry Smith   These routines are not usually accessed directly, rather solving is
104394b7f48cSBarry Smith   done through the KSP and PC interfaces.
10447b80b807SBarry Smith */
10457b80b807SBarry Smith 
1046d9274352SBarry Smith /*E
1047d9274352SBarry Smith     MatOrderingType - String with the name of a PETSc matrix ordering or the creation function
1048d9274352SBarry Smith        with an optional dynamic library name, for example
1049d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:orderingcreate()
1050d9274352SBarry Smith 
1051d9274352SBarry Smith    Level: beginner
1052d9274352SBarry Smith 
1053e5a9bf91SBarry Smith    Cannot use const because the PC objects manipulate the string
1054e5a9bf91SBarry Smith 
1055d9274352SBarry Smith .seealso: MatGetOrdering()
1056d9274352SBarry Smith E*/
10573eaad2caSSatish Balay #define MatOrderingType char*
10582692d6eeSBarry Smith #define MATORDERINGNATURAL     "natural"
10592692d6eeSBarry Smith #define MATORDERINGND          "nd"
10602692d6eeSBarry Smith #define MATORDERING1WD         "1wd"
10612692d6eeSBarry Smith #define MATORDERINGRCM         "rcm"
10622692d6eeSBarry Smith #define MATORDERINGQMD         "qmd"
10632692d6eeSBarry Smith #define MATORDERINGROWLENGTH   "rowlength"
10642692d6eeSBarry Smith #define MATORDERINGDSC_ND      "dsc_nd"         /* these three are only for DSCPACK, see its documentation for details */
10652692d6eeSBarry Smith #define MATORDERINGDSC_MMD     "dsc_mmd"
10662692d6eeSBarry Smith #define MATORDERINGDSC_MDF     "dsc_mdf"
1067312e372aSBarry Smith #define MATORDERINGAMD         "amd"            /* only works if UMFPACK is installed with PETSc */
1068b12f92e5SBarry Smith 
10697087cfbeSBarry Smith extern PetscErrorCode  MatGetOrdering(Mat,const MatOrderingType,IS*,IS*);
10707087cfbeSBarry Smith extern PetscErrorCode  MatGetOrderingList(PetscFList *list);
10717087cfbeSBarry Smith extern PetscErrorCode  MatOrderingRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,const MatOrderingType,IS*,IS*));
107230de9b25SBarry Smith 
107330de9b25SBarry Smith /*MC
10741890ba74SBarry Smith    MatOrderingRegisterDynamic - Adds a new sparse matrix ordering to the matrix package.
107530de9b25SBarry Smith 
107630de9b25SBarry Smith    Synopsis:
10771890ba74SBarry Smith    PetscErrorCode MatOrderingRegisterDynamic(const char *name_ordering,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatOrdering))
107830de9b25SBarry Smith 
107930de9b25SBarry Smith    Not Collective
108030de9b25SBarry Smith 
108130de9b25SBarry Smith    Input Parameters:
10822692d6eeSBarry Smith +  sname - name of ordering (for example MATORDERINGND)
108330de9b25SBarry Smith .  path - location of library where creation routine is
108430de9b25SBarry Smith .  name - name of function that creates the ordering type,a string
108530de9b25SBarry Smith -  function - function pointer that creates the ordering
108630de9b25SBarry Smith 
108730de9b25SBarry Smith    Level: developer
108830de9b25SBarry Smith 
108930de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
109030de9b25SBarry Smith    is ignored.
109130de9b25SBarry Smith 
109230de9b25SBarry Smith    Sample usage:
109330de9b25SBarry Smith .vb
109430de9b25SBarry Smith    MatOrderingRegisterDynamic("my_order",/home/username/my_lib/lib/libO/solaris/mylib.a,
109530de9b25SBarry Smith                "MyOrder",MyOrder);
109630de9b25SBarry Smith .ve
109730de9b25SBarry Smith 
109830de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
109930de9b25SBarry Smith $     MatOrderingSetType(part,"my_order)
110030de9b25SBarry Smith    or at runtime via the option
11012401956bSBarry Smith $     -pc_factor_mat_ordering_type my_order
110230de9b25SBarry Smith 
1103ab901514SBarry Smith    ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values.
110430de9b25SBarry Smith 
110530de9b25SBarry Smith .keywords: matrix, ordering, register
110630de9b25SBarry Smith 
110730de9b25SBarry Smith .seealso: MatOrderingRegisterDestroy(), MatOrderingRegisterAll()
110830de9b25SBarry Smith M*/
1109aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1110f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,0)
1111b12f92e5SBarry Smith #else
1112f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,d)
1113b12f92e5SBarry Smith #endif
111430de9b25SBarry Smith 
11157087cfbeSBarry Smith extern PetscErrorCode  MatOrderingRegisterDestroy(void);
11167087cfbeSBarry Smith extern PetscErrorCode  MatOrderingRegisterAll(const char[]);
1117ace3abfcSBarry Smith extern PetscBool  MatOrderingRegisterAllCalled;
1118b0a32e0cSBarry Smith extern PetscFList MatOrderingList;
1119d4fbbf0eSBarry Smith 
11207087cfbeSBarry Smith extern PetscErrorCode  MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS);
1121a2ce50c7SBarry Smith 
1122d91e6319SBarry Smith /*S
1123d90ac83dSHong Zhang     MatFactorShiftType - Numeric Shift.
1124d90ac83dSHong Zhang 
1125d90ac83dSHong Zhang    Level: beginner
1126d90ac83dSHong Zhang 
1127d90ac83dSHong Zhang S*/
1128d90ac83dSHong Zhang typedef enum {MAT_SHIFT_NONE,MAT_SHIFT_NONZERO,MAT_SHIFT_POSITIVE_DEFINITE,MAT_SHIFT_INBLOCKS} MatFactorShiftType;
1129d90ac83dSHong Zhang extern const char *MatFactorShiftTypes[];
1130d90ac83dSHong Zhang 
1131d90ac83dSHong Zhang /*S
11322401956bSBarry Smith    MatFactorInfo - Data passed into the matrix factorization routines
1133b00f7748SHong Zhang 
113461cad860SBarry Smith    In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE, that is use
113561cad860SBarry Smith $     MatFactorInfo  info(MAT_FACTORINFO_SIZE)
1136b00f7748SHong Zhang 
113715e8a5b3SHong Zhang    Notes: These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC.
1138b00f7748SHong Zhang 
113961cad860SBarry Smith       You can use MatFactorInfoInitialize() to set default values.
114061cad860SBarry Smith 
1141b00f7748SHong Zhang    Level: developer
1142b00f7748SHong Zhang 
1143d7d82daaSBarry Smith .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(),
1144d7d82daaSBarry Smith           MatFactorInfoInitialize()
1145b00f7748SHong Zhang 
1146b00f7748SHong Zhang S*/
1147b00f7748SHong Zhang typedef struct {
114815e8a5b3SHong Zhang   PetscReal     diagonal_fill;  /* force diagonal to fill in if initially not filled */
114985317021SBarry Smith   PetscReal     usedt;
115015e8a5b3SHong Zhang   PetscReal     dt;             /* drop tolerance */
1151b00f7748SHong Zhang   PetscReal     dtcol;          /* tolerance for pivoting */
115215e8a5b3SHong Zhang   PetscReal     dtcount;        /* maximum nonzeros to be allowed per row */
115367e28bfeSBarry Smith   PetscReal     fill;           /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */
1154348344bbSBarry Smith   PetscReal     levels;         /* ICC/ILU(levels) */
1155bcd9e38bSBarry Smith   PetscReal     pivotinblocks;  /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0
1156bcd9e38bSBarry Smith                                    factorization may be faster if do not pivot */
115715e8a5b3SHong Zhang   PetscReal     zeropivot;      /* pivot is called zero if less than this */
1158f4db908eSBarry Smith   PetscReal     shifttype;      /* type of shift added to matrix factor to prevent zero pivots */
1159f4db908eSBarry Smith   PetscReal     shiftamount;     /* how large the shift is */
116015e8a5b3SHong Zhang } MatFactorInfo;
1161ffa6d0a5SLois Curfman McInnes 
11627087cfbeSBarry Smith extern PetscErrorCode  MatFactorInfoInitialize(MatFactorInfo*);
11637087cfbeSBarry Smith extern PetscErrorCode  MatCholeskyFactor(Mat,IS,const MatFactorInfo*);
11647087cfbeSBarry Smith extern PetscErrorCode  MatCholeskyFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
11657087cfbeSBarry Smith extern PetscErrorCode  MatCholeskyFactorNumeric(Mat,Mat,const MatFactorInfo*);
11667087cfbeSBarry Smith extern PetscErrorCode  MatLUFactor(Mat,IS,IS,const MatFactorInfo*);
11677087cfbeSBarry Smith extern PetscErrorCode  MatILUFactor(Mat,IS,IS,const MatFactorInfo*);
11687087cfbeSBarry Smith extern PetscErrorCode  MatLUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
11697087cfbeSBarry Smith extern PetscErrorCode  MatILUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
11707087cfbeSBarry Smith extern PetscErrorCode  MatICCFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
11717087cfbeSBarry Smith extern PetscErrorCode  MatICCFactor(Mat,IS,const MatFactorInfo*);
11727087cfbeSBarry Smith extern PetscErrorCode  MatLUFactorNumeric(Mat,Mat,const MatFactorInfo*);
11737087cfbeSBarry Smith extern PetscErrorCode  MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*);
11747087cfbeSBarry Smith extern PetscErrorCode  MatSolve(Mat,Vec,Vec);
11757087cfbeSBarry Smith extern PetscErrorCode  MatForwardSolve(Mat,Vec,Vec);
11767087cfbeSBarry Smith extern PetscErrorCode  MatBackwardSolve(Mat,Vec,Vec);
11777087cfbeSBarry Smith extern PetscErrorCode  MatSolveAdd(Mat,Vec,Vec,Vec);
11787087cfbeSBarry Smith extern PetscErrorCode  MatSolveTranspose(Mat,Vec,Vec);
11797087cfbeSBarry Smith extern PetscErrorCode  MatSolveTransposeAdd(Mat,Vec,Vec,Vec);
11807087cfbeSBarry Smith extern PetscErrorCode  MatSolves(Mat,Vecs,Vecs);
11818ed539a5SBarry Smith 
11827087cfbeSBarry Smith extern PetscErrorCode  MatSetUnfactored(Mat);
1183bb5a7306SBarry Smith 
1184d91e6319SBarry Smith /*E
1185d91e6319SBarry Smith     MatSORType - What type of (S)SOR to perform
1186bb1eb677SSatish Balay 
1187d91e6319SBarry Smith     Level: beginner
1188d91e6319SBarry Smith 
1189d9274352SBarry Smith    May be bitwise ORd together
1190d9274352SBarry Smith 
1191d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
1192d91e6319SBarry Smith 
11934e7234bfSBarry Smith    MatSORType may be bitwise ORd together, so do not change the numbers
11944e7234bfSBarry Smith 
119541f059aeSBarry Smith .seealso: MatSOR()
1196d91e6319SBarry Smith E*/
1197ee50ffe9SBarry Smith typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3,
1198ee50ffe9SBarry Smith               SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8,
1199ee50ffe9SBarry Smith               SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16,
120084cb2905SBarry Smith               SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType;
12017087cfbeSBarry Smith extern PetscErrorCode  MatSOR(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
12028ed539a5SBarry Smith 
1203d4fbbf0eSBarry Smith /*
1204639f9d9dSBarry Smith     These routines are for efficiently computing Jacobians via finite differences.
1205639f9d9dSBarry Smith */
1206b12f92e5SBarry Smith 
1207d9274352SBarry Smith /*E
1208d9274352SBarry Smith     MatColoringType - String with the name of a PETSc matrix coloring or the creation function
1209d9274352SBarry Smith        with an optional dynamic library name, for example
1210d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:coloringcreate()
1211d9274352SBarry Smith 
1212d9274352SBarry Smith    Level: beginner
1213d9274352SBarry Smith 
1214d9274352SBarry Smith .seealso: MatGetColoring()
1215d9274352SBarry Smith E*/
1216a313700dSBarry Smith #define MatColoringType char*
12172692d6eeSBarry Smith #define MATCOLORINGNATURAL "natural"
12182692d6eeSBarry Smith #define MATCOLORINGSL      "sl"
12192692d6eeSBarry Smith #define MATCOLORINGLF      "lf"
12202692d6eeSBarry Smith #define MATCOLORINGID      "id"
1221b12f92e5SBarry Smith 
12227087cfbeSBarry Smith extern PetscErrorCode  MatGetColoring(Mat,const MatColoringType,ISColoring*);
12237087cfbeSBarry Smith extern PetscErrorCode  MatColoringRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,MatColoringType,ISColoring *));
122430de9b25SBarry Smith 
122530de9b25SBarry Smith /*MC
122630de9b25SBarry Smith    MatColoringRegisterDynamic - Adds a new sparse matrix coloring to the
122730de9b25SBarry Smith                                matrix package.
122830de9b25SBarry Smith 
122930de9b25SBarry Smith    Synopsis:
12301890ba74SBarry Smith    PetscErrorCode MatColoringRegisterDynamic(const char *name_coloring,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatColoring))
123130de9b25SBarry Smith 
123230de9b25SBarry Smith    Not Collective
123330de9b25SBarry Smith 
123430de9b25SBarry Smith    Input Parameters:
12352692d6eeSBarry Smith +  sname - name of Coloring (for example MATCOLORINGSL)
123630de9b25SBarry Smith .  path - location of library where creation routine is
123730de9b25SBarry Smith .  name - name of function that creates the Coloring type, a string
123830de9b25SBarry Smith -  function - function pointer that creates the coloring
123930de9b25SBarry Smith 
124030de9b25SBarry Smith    Level: developer
124130de9b25SBarry Smith 
124230de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
124330de9b25SBarry Smith    is ignored.
124430de9b25SBarry Smith 
124530de9b25SBarry Smith    Sample usage:
124630de9b25SBarry Smith .vb
124730de9b25SBarry Smith    MatColoringRegisterDynamic("my_color",/home/username/my_lib/lib/libO/solaris/mylib.a,
124830de9b25SBarry Smith                "MyColor",MyColor);
124930de9b25SBarry Smith .ve
125030de9b25SBarry Smith 
125130de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
125230de9b25SBarry Smith $     MatColoringSetType(part,"my_color")
125330de9b25SBarry Smith    or at runtime via the option
125430de9b25SBarry Smith $     -mat_coloring_type my_color
125530de9b25SBarry Smith 
1256ab901514SBarry Smith    $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
125730de9b25SBarry Smith 
125830de9b25SBarry Smith .keywords: matrix, Coloring, register
125930de9b25SBarry Smith 
126030de9b25SBarry Smith .seealso: MatColoringRegisterDestroy(), MatColoringRegisterAll()
126130de9b25SBarry Smith M*/
1262aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1263f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,0)
1264b12f92e5SBarry Smith #else
1265f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,d)
1266b12f92e5SBarry Smith #endif
126730de9b25SBarry Smith 
1268ace3abfcSBarry Smith extern PetscBool  MatColoringRegisterAllCalled;
1269f1144a30SSatish Balay 
12707087cfbeSBarry Smith extern PetscErrorCode  MatColoringRegisterAll(const char[]);
12717087cfbeSBarry Smith extern PetscErrorCode  MatColoringRegisterDestroy(void);
12727087cfbeSBarry Smith extern PetscErrorCode  MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*);
1273639f9d9dSBarry Smith 
1274d9274352SBarry Smith /*S
1275d9274352SBarry Smith      MatFDColoring - Object for computing a sparse Jacobian via finite differences
1276d9274352SBarry Smith         and coloring
1277639f9d9dSBarry Smith 
1278d9274352SBarry Smith    Level: beginner
1279d9274352SBarry Smith 
1280d9274352SBarry Smith   Concepts: coloring, sparse Jacobian, finite differences
1281d9274352SBarry Smith 
1282d9274352SBarry Smith .seealso:  MatFDColoringCreate()
1283d9274352SBarry Smith S*/
1284e2a1c21fSSatish Balay typedef struct _p_MatFDColoring* MatFDColoring;
1285639f9d9dSBarry Smith 
12867087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringCreate(Mat,ISColoring,MatFDColoring *);
1287fcfd50ebSBarry Smith extern PetscErrorCode  MatFDColoringDestroy(MatFDColoring*);
12887087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringView(MatFDColoring,PetscViewer);
12897087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*);
12907087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**);
12917087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal);
12927087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringSetFromOptions(MatFDColoring);
12937087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringApply(Mat,MatFDColoring,Vec,MatStructure*,void *);
12947087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringApplyTS(Mat,MatFDColoring,PetscReal,Vec,MatStructure*,void *);
12957087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringSetF(MatFDColoring,Vec);
12967087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,PetscInt*[]);
1297639f9d9dSBarry Smith /*
12980752156aSBarry Smith     These routines are for partitioning matrices: currently used only
12993eda8832SBarry Smith   for adjacency matrix, MatCreateMPIAdj().
13000752156aSBarry Smith */
1301ca161407SBarry Smith 
1302d9274352SBarry Smith /*S
1303d9274352SBarry Smith      MatPartitioning - Object for managing the partitioning of a matrix or graph
1304d9274352SBarry Smith 
1305d9274352SBarry Smith    Level: beginner
1306d9274352SBarry Smith 
1307d9274352SBarry Smith   Concepts: partitioning
1308d9274352SBarry Smith 
1309743a1026Svictorle .seealso:  MatPartitioningCreate(), MatPartitioningType
1310d9274352SBarry Smith S*/
131191e9ee9fSBarry Smith typedef struct _p_MatPartitioning* MatPartitioning;
1312d9274352SBarry Smith 
1313d9274352SBarry Smith /*E
13145bc6e7ccSvictorle     MatPartitioningType - String with the name of a PETSc matrix partitioning or the creation function
1315d9274352SBarry Smith        with an optional dynamic library name, for example
1316d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:partitioningcreate()
1317d9274352SBarry Smith 
1318d9274352SBarry Smith    Level: beginner
1319d9274352SBarry Smith 
1320b547a13bSvictorle .seealso: MatPartitioningCreate(), MatPartitioning
1321d9274352SBarry Smith E*/
1322a313700dSBarry Smith #define MatPartitioningType char*
13232692d6eeSBarry Smith #define MATPARTITIONINGCURRENT  "current"
13242692d6eeSBarry Smith #define MATPARTITIONINGSQUARE   "square"
13252692d6eeSBarry Smith #define MATPARTITIONINGPARMETIS "parmetis"
13262692d6eeSBarry Smith #define MATPARTITIONINGCHACO    "chaco"
13272692d6eeSBarry Smith #define MATPARTITIONINGJOSTLE   "jostle"
13282692d6eeSBarry Smith #define MATPARTITIONINGPARTY    "party"
13292692d6eeSBarry Smith #define MATPARTITIONINGSCOTCH   "scotch"
133071306c60Sjroman 
1331ca161407SBarry Smith 
13327087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningCreate(MPI_Comm,MatPartitioning*);
13337087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningSetType(MatPartitioning,const MatPartitioningType);
13347087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningSetNParts(MatPartitioning,PetscInt);
13357087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningSetAdjacency(MatPartitioning,Mat);
13367087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]);
13377087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []);
13387087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningApply(MatPartitioning,IS*);
1339fcfd50ebSBarry Smith extern PetscErrorCode  MatPartitioningDestroy(MatPartitioning*);
13402aabb6bbSBarry Smith 
13417087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatPartitioning));
134230de9b25SBarry Smith 
134330de9b25SBarry Smith /*MC
134430de9b25SBarry Smith    MatPartitioningRegisterDynamic - Adds a new sparse matrix partitioning to the
134530de9b25SBarry Smith    matrix package.
134630de9b25SBarry Smith 
134730de9b25SBarry Smith    Synopsis:
13481890ba74SBarry Smith    PetscErrorCode MatPartitioningRegisterDynamic(const char *name_partitioning,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatPartitioning))
134930de9b25SBarry Smith 
135030de9b25SBarry Smith    Not Collective
135130de9b25SBarry Smith 
135230de9b25SBarry Smith    Input Parameters:
13532692d6eeSBarry Smith +  sname - name of partitioning (for example MATPARTITIONINGCURRENT) or parmetis
135430de9b25SBarry Smith .  path - location of library where creation routine is
135530de9b25SBarry Smith .  name - name of function that creates the partitioning type, a string
135630de9b25SBarry Smith -  function - function pointer that creates the partitioning type
135730de9b25SBarry Smith 
135830de9b25SBarry Smith    Level: developer
135930de9b25SBarry Smith 
136030de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
136130de9b25SBarry Smith    is ignored.
136230de9b25SBarry Smith 
136330de9b25SBarry Smith    Sample usage:
136430de9b25SBarry Smith .vb
136530de9b25SBarry Smith    MatPartitioningRegisterDynamic("my_part",/home/username/my_lib/lib/libO/solaris/mylib.a,
136630de9b25SBarry Smith                "MyPartCreate",MyPartCreate);
136730de9b25SBarry Smith .ve
136830de9b25SBarry Smith 
136930de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
137030de9b25SBarry Smith $     MatPartitioningSetType(part,"my_part")
137130de9b25SBarry Smith    or at runtime via the option
137230de9b25SBarry Smith $     -mat_partitioning_type my_part
137330de9b25SBarry Smith 
1374ab901514SBarry Smith    $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
137530de9b25SBarry Smith 
137630de9b25SBarry Smith .keywords: matrix, partitioning, register
137730de9b25SBarry Smith 
137830de9b25SBarry Smith .seealso: MatPartitioningRegisterDestroy(), MatPartitioningRegisterAll()
137930de9b25SBarry Smith M*/
1380aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1381f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,0)
13822aabb6bbSBarry Smith #else
1383f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,d)
13842aabb6bbSBarry Smith #endif
13852aabb6bbSBarry Smith 
1386ace3abfcSBarry Smith extern PetscBool  MatPartitioningRegisterAllCalled;
1387f1144a30SSatish Balay 
13887087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningRegisterAll(const char[]);
13897087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningRegisterDestroy(void);
13902bad1931SBarry Smith 
13917087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningView(MatPartitioning,PetscViewer);
13927087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningSetFromOptions(MatPartitioning);
13937087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningGetType(MatPartitioning,const MatPartitioningType*);
1394ca161407SBarry Smith 
13957087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningParmetisSetCoarseSequential(MatPartitioning);
13967087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningParmetisGetEdgeCut(MatPartitioning, PetscInt *);
13970752156aSBarry Smith 
13987087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningJostleSetCoarseLevel(MatPartitioning,PetscReal);
13997087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningJostleSetCoarseSequential(MatPartitioning);
140071306c60Sjroman 
1401954b3ec6SBarry Smith typedef enum {MP_CHACO_MULTILEVEL_KL,MP_CHACO_SPECTRAL,MP_CHACO_LINEAR,MP_CHACO_RANDOM, MP_CHACO_SCATTERED} MPChacoGlobalType;
14027087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningChacoSetGlobal(MatPartitioning, MPChacoGlobalType);
140371306c60Sjroman typedef enum { MP_CHACO_KERNIGHAN_LIN, MP_CHACO_NONE } MPChacoLocalType;
14047087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningChacoSetLocal(MatPartitioning, MPChacoLocalType);
14057087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal);
140671306c60Sjroman typedef enum { MP_CHACO_LANCZOS, MP_CHACO_RQI_SYMMLQ } MPChacoEigenType;
14077087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType);
14087087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningChacoSetEigenTol(MatPartitioning, PetscReal);
14097087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningChacoSetEigenNumber(MatPartitioning, PetscInt);
141071306c60Sjroman 
141171306c60Sjroman #define MP_PARTY_OPT "opt"
141271306c60Sjroman #define MP_PARTY_LIN "lin"
141371306c60Sjroman #define MP_PARTY_SCA "sca"
141471306c60Sjroman #define MP_PARTY_RAN "ran"
141571306c60Sjroman #define MP_PARTY_GBF "gbf"
141671306c60Sjroman #define MP_PARTY_GCF "gcf"
141771306c60Sjroman #define MP_PARTY_BUB "bub"
141871306c60Sjroman #define MP_PARTY_DEF "def"
14197087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningPartySetGlobal(MatPartitioning, const char*);
142071306c60Sjroman #define MP_PARTY_HELPFUL_SETS "hs"
142171306c60Sjroman #define MP_PARTY_KERNIGHAN_LIN "kl"
142271306c60Sjroman #define MP_PARTY_NONE "no"
14237087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningPartySetLocal(MatPartitioning, const char*);
14247087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal);
14257087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningPartySetBipart(MatPartitioning,PetscBool );
14267087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscBool );
142771306c60Sjroman 
142871306c60Sjroman typedef enum { MP_SCOTCH_GREEDY, MP_SCOTCH_GPS, MP_SCOTCH_GR_GPS } MPScotchGlobalType;
14297087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningScotchSetArch(MatPartitioning,const char*);
14307087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningScotchSetMultilevel(MatPartitioning);
14317087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningScotchSetGlobal(MatPartitioning,MPScotchGlobalType);
14327087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningScotchSetCoarseLevel(MatPartitioning,PetscReal);
14337087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningScotchSetHostList(MatPartitioning,const char*);
143471306c60Sjroman typedef enum { MP_SCOTCH_KERNIGHAN_LIN, MP_SCOTCH_NONE } MPScotchLocalType;
14357087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningScotchSetLocal(MatPartitioning,MPScotchLocalType);
14367087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningScotchSetMapping(MatPartitioning);
14377087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningScotchSetStrategy(MatPartitioning,char*);
143871306c60Sjroman 
143909573ac7SBarry Smith extern PetscErrorCode MatMeshToVertexGraph(Mat,PetscInt,Mat*);
144009573ac7SBarry Smith extern PetscErrorCode MatMeshToCellGraph(Mat,PetscInt,Mat*);
1441591294e4SBarry Smith 
14420752156aSBarry Smith /*
14430a835dfdSSatish Balay     If you add entries here you must also add them to finclude/petscmat.h
1444d4fbbf0eSBarry Smith */
14451c1c02c0SLois Curfman McInnes typedef enum { MATOP_SET_VALUES=0,
14461c1c02c0SLois Curfman McInnes                MATOP_GET_ROW=1,
14471c1c02c0SLois Curfman McInnes                MATOP_RESTORE_ROW=2,
14481c1c02c0SLois Curfman McInnes                MATOP_MULT=3,
14491c1c02c0SLois Curfman McInnes                MATOP_MULT_ADD=4,
14507c922b88SBarry Smith                MATOP_MULT_TRANSPOSE=5,
14517c922b88SBarry Smith                MATOP_MULT_TRANSPOSE_ADD=6,
14521c1c02c0SLois Curfman McInnes                MATOP_SOLVE=7,
14531c1c02c0SLois Curfman McInnes                MATOP_SOLVE_ADD=8,
14547c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE=9,
14557c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE_ADD=10,
14561c1c02c0SLois Curfman McInnes                MATOP_LUFACTOR=11,
14571c1c02c0SLois Curfman McInnes                MATOP_CHOLESKYFACTOR=12,
1458a714c06dSBarry Smith                MATOP_SOR=13,
14591c1c02c0SLois Curfman McInnes                MATOP_TRANSPOSE=14,
14601c1c02c0SLois Curfman McInnes                MATOP_GETINFO=15,
14611c1c02c0SLois Curfman McInnes                MATOP_EQUAL=16,
14621c1c02c0SLois Curfman McInnes                MATOP_GET_DIAGONAL=17,
14631c1c02c0SLois Curfman McInnes                MATOP_DIAGONAL_SCALE=18,
14641c1c02c0SLois Curfman McInnes                MATOP_NORM=19,
14651c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_BEGIN=20,
14661c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_END=21,
1467d519adbfSMatthew Knepley                MATOP_SET_OPTION=22,
1468d519adbfSMatthew Knepley                MATOP_ZERO_ENTRIES=23,
1469d519adbfSMatthew Knepley                MATOP_ZERO_ROWS=24,
1470d519adbfSMatthew Knepley                MATOP_LUFACTOR_SYMBOLIC=25,
1471d519adbfSMatthew Knepley                MATOP_LUFACTOR_NUMERIC=26,
1472d519adbfSMatthew Knepley                MATOP_CHOLESKY_FACTOR_SYMBOLIC=27,
1473d519adbfSMatthew Knepley                MATOP_CHOLESKY_FACTOR_NUMERIC=28,
1474d519adbfSMatthew Knepley                MATOP_SETUP_PREALLOCATION=29,
1475d519adbfSMatthew Knepley                MATOP_ILUFACTOR_SYMBOLIC=30,
1476d519adbfSMatthew Knepley                MATOP_ICCFACTOR_SYMBOLIC=31,
1477d519adbfSMatthew Knepley                MATOP_GET_ARRAY=32,
1478d519adbfSMatthew Knepley                MATOP_RESTORE_ARRAY=33,
1479d519adbfSMatthew Knepley                MATOP_DUPLICATE=34,
1480d519adbfSMatthew Knepley                MATOP_FORWARD_SOLVE=35,
1481d519adbfSMatthew Knepley                MATOP_BACKWARD_SOLVE=36,
1482d519adbfSMatthew Knepley                MATOP_ILUFACTOR=37,
1483d519adbfSMatthew Knepley                MATOP_ICCFACTOR=38,
1484d519adbfSMatthew Knepley                MATOP_AXPY=39,
1485d519adbfSMatthew Knepley                MATOP_GET_SUBMATRICES=40,
1486d519adbfSMatthew Knepley                MATOP_INCREASE_OVERLAP=41,
1487d519adbfSMatthew Knepley                MATOP_GET_VALUES=42,
1488d519adbfSMatthew Knepley                MATOP_COPY=43,
1489d519adbfSMatthew Knepley                MATOP_GET_ROW_MAX=44,
1490d519adbfSMatthew Knepley                MATOP_SCALE=45,
1491d519adbfSMatthew Knepley                MATOP_SHIFT=46,
149235153367SBarry Smith                MATOP_DIAGONAL_SET=47,
1493d519adbfSMatthew Knepley                MATOP_ILUDT_FACTOR=48,
1494d519adbfSMatthew Knepley                MATOP_SET_BLOCK_SIZE=49,
1495d519adbfSMatthew Knepley                MATOP_GET_ROW_IJ=50,
1496d519adbfSMatthew Knepley                MATOP_RESTORE_ROW_IJ=51,
1497d519adbfSMatthew Knepley                MATOP_GET_COLUMN_IJ=52,
1498d519adbfSMatthew Knepley                MATOP_RESTORE_COLUMN_IJ=53,
1499d519adbfSMatthew Knepley                MATOP_FDCOLORING_CREATE=54,
1500d519adbfSMatthew Knepley                MATOP_COLORING_PATCH=55,
1501d519adbfSMatthew Knepley                MATOP_SET_UNFACTORED=56,
1502d519adbfSMatthew Knepley                MATOP_PERMUTE=57,
1503d519adbfSMatthew Knepley                MATOP_SET_VALUES_BLOCKED=58,
1504d519adbfSMatthew Knepley                MATOP_GET_SUBMATRIX=59,
1505d519adbfSMatthew Knepley                MATOP_DESTROY=60,
1506d519adbfSMatthew Knepley                MATOP_VIEW=61,
1507d519adbfSMatthew Knepley                MATOP_CONVERT_FROM=62,
1508d519adbfSMatthew Knepley                MATOP_USE_SCALED_FORM=63,
1509d519adbfSMatthew Knepley                MATOP_SCALE_SYSTEM=64,
1510d519adbfSMatthew Knepley                MATOP_UNSCALE_SYSTEM=65,
1511d519adbfSMatthew Knepley                MATOP_SET_LOCAL_TO_GLOBAL_MAP=66,
1512d519adbfSMatthew Knepley                MATOP_SET_VALUES_LOCAL=67,
1513d519adbfSMatthew Knepley                MATOP_ZERO_ROWS_LOCAL=68,
1514d519adbfSMatthew Knepley                MATOP_GET_ROW_MAX_ABS=69,
1515d519adbfSMatthew Knepley                MATOP_GET_ROW_MIN_ABS=70,
1516d519adbfSMatthew Knepley                MATOP_CONVERT=71,
1517d519adbfSMatthew Knepley                MATOP_SET_COLORING=72,
1518d519adbfSMatthew Knepley                MATOP_SET_VALUES_ADIC=73,
1519d519adbfSMatthew Knepley                MATOP_SET_VALUES_ADIFOR=74,
1520d519adbfSMatthew Knepley                MATOP_FD_COLORING_APPLY=75,
1521d519adbfSMatthew Knepley                MATOP_SET_FROM_OPTIONS=76,
1522d519adbfSMatthew Knepley                MATOP_MULT_CON=77,
1523d519adbfSMatthew Knepley                MATOP_MULT_TRANSPOSE_CON=78,
1524d519adbfSMatthew Knepley                MATOP_PERMUTE_SPARSIFY=79,
1525d519adbfSMatthew Knepley                MATOP_MULT_MULTIPLE=80,
1526d519adbfSMatthew Knepley                MATOP_SOLVE_MULTIPLE=81,
1527d519adbfSMatthew Knepley                MATOP_GET_INERTIA=82,
1528d519adbfSMatthew Knepley                MATOP_LOAD=83,
1529d519adbfSMatthew Knepley                MATOP_IS_SYMMETRIC=84,
1530d519adbfSMatthew Knepley                MATOP_IS_HERMITIAN=85,
1531d519adbfSMatthew Knepley                MATOP_IS_STRUCTURALLY_SYMMETRIC=86,
153241f059aeSBarry Smith                MATOP_DUMMY=87,
1533d519adbfSMatthew Knepley                MATOP_GET_VECS=88,
1534d519adbfSMatthew Knepley                MATOP_MAT_MULT=89,
1535d519adbfSMatthew Knepley                MATOP_MAT_MULT_SYMBOLIC=90,
1536d519adbfSMatthew Knepley                MATOP_MAT_MULT_NUMERIC=91,
1537d519adbfSMatthew Knepley                MATOP_PTAP=92,
1538d519adbfSMatthew Knepley                MATOP_PTAP_SYMBOLIC=93,
1539d519adbfSMatthew Knepley                MATOP_PTAP_NUMERIC=94,
1540d519adbfSMatthew Knepley                MATOP_MAT_MULTTRANSPOSE=95,
15411763ac4eSSatish Balay                MATOP_MAT_MULTTRANSPOSE_SYM=96,
15421763ac4eSSatish Balay                MATOP_MAT_MULTTRANSPOSE_NUM=97,
1543d519adbfSMatthew Knepley                MATOP_PTAP_SYMBOLIC_SEQAIJ=98,
1544d519adbfSMatthew Knepley                MATOP_PTAP_NUMERIC_SEQAIJ=99,
1545d519adbfSMatthew Knepley                MATOP_PTAP_SYMBOLIC_MPIAIJ=100,
1546d519adbfSMatthew Knepley                MATOP_PTAP_NUMERIC_MPIAIJ=101,
1547d519adbfSMatthew Knepley                MATOP_CONJUGATE=102,
1548d519adbfSMatthew Knepley                MATOP_SET_SIZES=103,
1549d519adbfSMatthew Knepley                MATOP_SET_VALUES_ROW=104,
1550d519adbfSMatthew Knepley                MATOP_REAL_PART=105,
1551d519adbfSMatthew Knepley                MATOP_IMAG_PART=106,
1552d519adbfSMatthew Knepley                MATOP_GET_ROW_UTRIANGULAR=107,
1553d519adbfSMatthew Knepley                MATOP_RESTORE_ROW_UTRIANGULAR=108,
1554d519adbfSMatthew Knepley                MATOP_MATSOLVE=109,
1555d519adbfSMatthew Knepley                MATOP_GET_REDUNDANTMATRIX=110,
1556d519adbfSMatthew Knepley                MATOP_GET_ROW_MIN=111,
1557d519adbfSMatthew Knepley                MATOP_GET_COLUMN_VEC=112,
1558d519adbfSMatthew Knepley                MATOP_MISSING_DIAGONAL=113,
1559d519adbfSMatthew Knepley                MATOP_MATGETSEQNONZEROSTRUCTURE=114,
156089c6957cSBarry Smith                MATOP_CREATE=115,
156189c6957cSBarry Smith                MATOP_GET_GHOSTS=116,
1562ea38565cSJed Brown                MATOP_GET_LOCALSUBMATRIX=117,
1563ea38565cSJed Brown                MATOP_RESTORE_LOCALSUBMATRIX=118,
1564eeffb40dSHong Zhang                MATOP_MULT_DIAGONAL_BLOCK=119,
1565547795f9SHong Zhang                MATOP_HERMITIANTRANSPOSE=120,
1566547795f9SHong Zhang                MATOP_MULTHERMITIANTRANSPOSE=121,
1567fbdbba38SShri Abhyankar                MATOP_MULTHERMITIANTRANSPOSEADD=122,
1568b9614d88SDmitry Karpeev                MATOP_GETMULTIPROCBLOCK=123,
1569b9614d88SDmitry Karpeev 	       MATOP_GET_SUBMATRICES_PARALLEL=128
1570fae171e0SBarry Smith              } MatOperation;
15717087cfbeSBarry Smith extern PetscErrorCode  MatHasOperation(Mat,MatOperation,PetscBool *);
15727087cfbeSBarry Smith extern PetscErrorCode  MatShellSetOperation(Mat,MatOperation,void(*)(void));
15737087cfbeSBarry Smith extern PetscErrorCode  MatShellGetOperation(Mat,MatOperation,void(**)(void));
15747087cfbeSBarry Smith extern PetscErrorCode  MatShellSetContext(Mat,void*);
1575112a2221SBarry Smith 
157690ace30eSBarry Smith /*
157790ace30eSBarry Smith    Codes for matrices stored on disk. By default they are
157890ace30eSBarry Smith    stored in a universal format. By changing the format with
15797973ad04SBarry Smith    PetscViewerSetFormat(viewer,PETSC_VIEWER_NATIVE); the matrices will
158090ace30eSBarry Smith    be stored in a way natural for the matrix, for example dense matrices
158190ace30eSBarry Smith    would be stored as dense. Matrices stored this way may only be
1582f4403165SShri Abhyankar    read into matrices of the same type.
158390ace30eSBarry Smith */
158490ace30eSBarry Smith #define MATRIX_BINARY_FORMAT_DENSE -1
158590ace30eSBarry Smith 
15867087cfbeSBarry Smith extern PetscErrorCode  MatMPIBAIJSetHashTableFactor(Mat,PetscReal);
15877087cfbeSBarry Smith extern PetscErrorCode  MatISGetLocalMat(Mat,Mat*);
15881f4e1ec7SBarry Smith 
1589d9274352SBarry Smith /*S
1590d9274352SBarry Smith      MatNullSpace - Object that removes a null space from a vector, i.e.
1591d9274352SBarry Smith          orthogonalizes the vector to a subsapce
1592d9274352SBarry Smith 
1593f7a9e4ceSBarry Smith    Level: advanced
1594d9274352SBarry Smith 
1595d9274352SBarry Smith   Concepts: matrix; linear operator, null space
1596d9274352SBarry Smith 
15976e1639daSBarry Smith   Users manual sections:
15986e1639daSBarry Smith .   sec_singular
15996e1639daSBarry Smith 
1600d9274352SBarry Smith .seealso:  MatNullSpaceCreate()
1601d9274352SBarry Smith S*/
160274637425SBarry Smith typedef struct _p_MatNullSpace* MatNullSpace;
1603d9274352SBarry Smith 
16047087cfbeSBarry Smith extern PetscErrorCode  MatNullSpaceCreate(MPI_Comm,PetscBool ,PetscInt,const Vec[],MatNullSpace*);
16057087cfbeSBarry Smith extern PetscErrorCode  MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(MatNullSpace,Vec,void*),void*);
1606d34fcf5fSBarry Smith extern PetscErrorCode  MatNullSpaceDestroy(MatNullSpace*);
16077087cfbeSBarry Smith extern PetscErrorCode  MatNullSpaceRemove(MatNullSpace,Vec,Vec*);
16087087cfbeSBarry Smith extern PetscErrorCode  MatNullSpaceAttach(Mat,MatNullSpace);
16097087cfbeSBarry Smith extern PetscErrorCode  MatNullSpaceTest(MatNullSpace,Mat,PetscBool  *);
1610b717e993SJed Brown extern PetscErrorCode  MatNullSpaceView(MatNullSpace,PetscViewer);
161174637425SBarry Smith 
16127087cfbeSBarry Smith extern PetscErrorCode  MatReorderingSeqSBAIJ(Mat,IS);
16137087cfbeSBarry Smith extern PetscErrorCode  MatMPISBAIJSetHashTableFactor(Mat,PetscReal);
16147087cfbeSBarry Smith extern PetscErrorCode  MatSeqSBAIJSetColumnIndices(Mat,PetscInt *);
16157087cfbeSBarry Smith extern PetscErrorCode  MatSeqBAIJInvertBlockDiagonal(Mat);
16163f1d51d7SBarry Smith 
16177087cfbeSBarry Smith extern PetscErrorCode  MatCreateMAIJ(Mat,PetscInt,Mat*);
16187087cfbeSBarry Smith extern PetscErrorCode  MatMAIJRedimension(Mat,PetscInt,Mat*);
16197087cfbeSBarry Smith extern PetscErrorCode  MatMAIJGetAIJ(Mat,Mat*);
1620c4f061fbSSatish Balay 
16217087cfbeSBarry Smith extern PetscErrorCode  MatComputeExplicitOperator(Mat,Mat*);
1622b0a32e0cSBarry Smith 
16237087cfbeSBarry Smith extern PetscErrorCode  MatDiagonalScaleLocal(Mat,Vec);
162404f1ad80SBarry Smith 
16257087cfbeSBarry Smith extern PetscErrorCode  MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*);
16267087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetBase(Mat,Vec,Vec);
16277087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*);
16287087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*));
16297087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec));
16307087cfbeSBarry Smith extern PetscErrorCode  MatMFFDAddNullSpace(Mat,MatNullSpace);
16317087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt);
16327087cfbeSBarry Smith extern PetscErrorCode  MatMFFDResetHHistory(Mat);
16337087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetFunctionError(Mat,PetscReal);
16347087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetPeriod(Mat,PetscInt);
16357087cfbeSBarry Smith extern PetscErrorCode  MatMFFDGetH(Mat,PetscScalar *);
16367087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetOptionsPrefix(Mat,const char[]);
16377087cfbeSBarry Smith extern PetscErrorCode  MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*);
16387087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*);
1639e884886fSBarry Smith 
16406370053bSBarry Smith /*S
16416370053bSBarry Smith     MatMFFD - A data structured used to manage the computation of the h differencing parameter for matrix-free
16426370053bSBarry Smith               Jacobian vector products
1643e884886fSBarry Smith 
16446370053bSBarry Smith     Notes: MATMFFD is a specific MatType which uses the MatMFFD data structure
16456370053bSBarry Smith 
16466370053bSBarry Smith            MatMFFD*() methods actually take the Mat as their first argument. Not a MatMFFD data structure
16476370053bSBarry Smith 
16486370053bSBarry Smith     Level: developer
16496370053bSBarry Smith 
16506370053bSBarry Smith .seealso: MATMFFD, MatCreateMFFD(), MatMFFDSetFuction(), MatMFFDSetType(), MatMFFDRegister()
16516370053bSBarry Smith S*/
1652e884886fSBarry Smith typedef struct _p_MatMFFD* MatMFFD;
1653e884886fSBarry Smith 
1654e884886fSBarry Smith /*E
1655e884886fSBarry Smith     MatMFFDType - algorithm used to compute the h used in computing matrix-vector products via differencing of the function
1656e884886fSBarry Smith 
1657e884886fSBarry Smith    Level: beginner
1658e884886fSBarry Smith 
1659e884886fSBarry Smith .seealso: MatMFFDSetType(), MatMFFDRegister()
1660e884886fSBarry Smith E*/
1661a313700dSBarry Smith #define MatMFFDType char*
1662e884886fSBarry Smith #define MATMFFD_DS  "ds"
1663e884886fSBarry Smith #define MATMFFD_WP  "wp"
1664e884886fSBarry Smith 
16657087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetType(Mat,const MatMFFDType);
16667087cfbeSBarry Smith extern PetscErrorCode  MatMFFDRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatMFFD));
1667e884886fSBarry Smith 
1668e884886fSBarry Smith /*MC
1669e884886fSBarry Smith    MatMFFDRegisterDynamic - Adds a method to the MatMFFD registry.
1670e884886fSBarry Smith 
1671e884886fSBarry Smith    Synopsis:
16721890ba74SBarry Smith    PetscErrorCode MatMFFDRegisterDynamic(const char *name_solver,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatMFFD))
1673e884886fSBarry Smith 
1674e884886fSBarry Smith    Not Collective
1675e884886fSBarry Smith 
1676e884886fSBarry Smith    Input Parameters:
1677e884886fSBarry Smith +  name_solver - name of a new user-defined compute-h module
1678e884886fSBarry Smith .  path - path (either absolute or relative) the library containing this solver
1679e884886fSBarry Smith .  name_create - name of routine to create method context
1680e884886fSBarry Smith -  routine_create - routine to create method context
1681e884886fSBarry Smith 
1682e884886fSBarry Smith    Level: developer
1683e884886fSBarry Smith 
1684e884886fSBarry Smith    Notes:
1685e884886fSBarry Smith    MatMFFDRegisterDynamic() may be called multiple times to add several user-defined solvers.
1686e884886fSBarry Smith 
1687e884886fSBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
1688e884886fSBarry Smith    is ignored.
1689e884886fSBarry Smith 
1690e884886fSBarry Smith    Sample usage:
1691e884886fSBarry Smith .vb
1692e884886fSBarry Smith    MatMFFDRegisterDynamic("my_h",/home/username/my_lib/lib/libO/solaris/mylib.a,
1693e884886fSBarry Smith                "MyHCreate",MyHCreate);
1694e884886fSBarry Smith .ve
1695e884886fSBarry Smith 
1696e884886fSBarry Smith    Then, your solver can be chosen with the procedural interface via
1697e884886fSBarry Smith $     MatMFFDSetType(mfctx,"my_h")
1698e884886fSBarry Smith    or at runtime via the option
1699e884886fSBarry Smith $     -snes_mf_type my_h
1700e884886fSBarry Smith 
1701e884886fSBarry Smith .keywords: MatMFFD, register
1702e884886fSBarry Smith 
1703e884886fSBarry Smith .seealso: MatMFFDRegisterAll(), MatMFFDRegisterDestroy()
1704e884886fSBarry Smith M*/
1705e884886fSBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1706e884886fSBarry Smith #define MatMFFDRegisterDynamic(a,b,c,d) MatMFFDRegister(a,b,c,0)
1707e884886fSBarry Smith #else
1708e884886fSBarry Smith #define MatMFFDRegisterDynamic(a,b,c,d) MatMFFDRegister(a,b,c,d)
1709e884886fSBarry Smith #endif
1710e884886fSBarry Smith 
17117087cfbeSBarry Smith extern PetscErrorCode  MatMFFDRegisterAll(const char[]);
17127087cfbeSBarry Smith extern PetscErrorCode  MatMFFDRegisterDestroy(void);
17137087cfbeSBarry Smith extern PetscErrorCode  MatMFFDDSSetUmin(Mat,PetscReal);
17147087cfbeSBarry Smith extern PetscErrorCode  MatMFFDWPSetComputeNormU(Mat,PetscBool );
1715e884886fSBarry Smith 
1716e884886fSBarry Smith 
17177087cfbeSBarry Smith extern PetscErrorCode  PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *);
17187087cfbeSBarry Smith extern PetscErrorCode  PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *);
17197dbadf16SMatthew Knepley 
172097969023SHong Zhang /*
172197969023SHong Zhang    PETSc interface to MUMPS
172297969023SHong Zhang */
172397969023SHong Zhang #ifdef PETSC_HAVE_MUMPS
1724f250808bSBarry Smith extern PetscErrorCode  MatMumpsSetIcntl(Mat,PetscInt,PetscInt);
172597969023SHong Zhang #endif
172623a5497aSJed Brown 
1727d954db57SHong Zhang /*
1728d954db57SHong Zhang    PETSc interface to SUPERLU
1729d954db57SHong Zhang */
1730d954db57SHong Zhang #ifdef PETSC_HAVE_SUPERLU
1731f250808bSBarry Smith extern PetscErrorCode  MatSuperluSetILUDropTol(Mat,PetscReal);
1732d954db57SHong Zhang #endif
1733d954db57SHong Zhang 
1734c8883902SJed Brown extern PetscErrorCode MatCreateNest(MPI_Comm,PetscInt,const IS[],PetscInt,const IS[],const Mat[],Mat*);
1735c8883902SJed Brown extern PetscErrorCode MatNestGetSize(Mat,PetscInt*,PetscInt*);
1736c8883902SJed Brown extern PetscErrorCode MatNestGetSubMats(Mat,PetscInt*,PetscInt*,Mat***);
1737c8883902SJed Brown extern PetscErrorCode MatNestGetSubMat(Mat,PetscInt,PetscInt,Mat*);
17387087cfbeSBarry Smith extern PetscErrorCode MatNestSetVecType(Mat,const VecType);
1739c8883902SJed Brown extern PetscErrorCode MatNestSetSubMats(Mat,PetscInt,const IS[],PetscInt,const IS[],const Mat[]);
1740d8588912SDave May 
174123a5497aSJed Brown PETSC_EXTERN_CXX_END
174223a5497aSJed Brown #endif
1743