xref: /petsc/include/petscmat.h (revision dedccee8f225f91da5577b2f35bba0c3cd8d7504)
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"
58cebc7f6cSBarry Smith #define MATDAAD            "daad"
59cebc7f6cSBarry Smith #define MATMFFD            "mffd"
60c8a8475eSBarry Smith #define MATNORMAL          "normal"
61ab92ecdeSBarry Smith #define MATLRC             "lrc"
622a6744ebSBarry Smith #define MATSCATTER         "scatter"
63421e10b8SBarry Smith #define MATBLOCKMAT        "blockmat"
64793850ffSBarry Smith #define MATCOMPOSITE       "composite"
651f8c7532SHong Zhang #define MATFFT             "fft"
661f8c7532SHong Zhang #define MATFFTW              "fftw"
67e133240eSMatthew G Knepley #define MATSEQCUFFT          "seqcufft"
68557cca28SSatish Balay #define MATTRANSPOSEMAT    "transpose"
6972ca8751SBarry Smith #define MATSCHURCOMPLEMENT "schurcomplement"
701d6018f0SLisandro Dalcin #define MATPYTHON          "python"
71f91d8e95SBarry Smith #define MATHYPRESTRUCT     "hyprestruct"
72a9e6138eSGlenn Hammond #define MATHYPRESSTRUCT    "hypresstruct"
73ee1cef2cSJed Brown #define MATSUBMATRIX       "submatrix"
742c0dbf93SJed Brown #define MATLOCALREF        "localref"
75d8588912SDave May #define MATNEST            "nest"
76773941d6SBarry Smith 
779989ab13SBarry Smith /*E
78c7393fdbSBarry Smith     MatSolverPackage - String with the name of a PETSc matrix solver type.
799989ab13SBarry Smith 
809989ab13SBarry Smith     For example: "petsc" indicates what PETSc provides, "superlu" indicates either
819989ab13SBarry Smith        SuperLU or SuperLU_Dist etc.
829989ab13SBarry Smith 
839989ab13SBarry Smith 
849989ab13SBarry Smith    Level: beginner
859989ab13SBarry Smith 
865c9eb25fSBarry Smith .seealso: MatGetFactor(), Mat, MatSetType(), MatType
879989ab13SBarry Smith E*/
88c7393fdbSBarry Smith #define MatSolverPackage char*
892692d6eeSBarry Smith #define MATSOLVERSPOOLES      "spooles"
902692d6eeSBarry Smith #define MATSOLVERSUPERLU      "superlu"
912692d6eeSBarry Smith #define MATSOLVERSUPERLU_DIST "superlu_dist"
922692d6eeSBarry Smith #define MATSOLVERUMFPACK      "umfpack"
9320db9a53SJed Brown #define MATSOLVERCHOLMOD      "cholmod"
942692d6eeSBarry Smith #define MATSOLVERESSL         "essl"
952692d6eeSBarry Smith #define MATSOLVERLUSOL        "lusol"
962692d6eeSBarry Smith #define MATSOLVERMUMPS        "mumps"
972692d6eeSBarry Smith #define MATSOLVERPASTIX       "pastix"
982692d6eeSBarry Smith #define MATSOLVERDSCPACK      "dscpack"
992692d6eeSBarry Smith #define MATSOLVERMATLAB       "matlab"
1002692d6eeSBarry Smith #define MATSOLVERPETSC        "petsc"
1012692d6eeSBarry Smith #define MATSOLVERPLAPACK      "plapack"
1022692d6eeSBarry Smith #define MATSOLVERBAS          "bas"
103773941d6SBarry Smith 
104b24902e0SBarry Smith /*E
105b24902e0SBarry Smith     MatFactorType - indicates what type of factorization is requested
106b24902e0SBarry Smith 
107b24902e0SBarry Smith     Level: beginner
108b24902e0SBarry Smith 
109b24902e0SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
110b24902e0SBarry Smith 
111c7393fdbSBarry Smith .seealso: MatSolverPackage, MatGetFactor()
112b24902e0SBarry Smith E*/
113599ef60dSHong Zhang typedef enum {MAT_FACTOR_NONE, MAT_FACTOR_LU, MAT_FACTOR_CHOLESKY, MAT_FACTOR_ILU, MAT_FACTOR_ICC,MAT_FACTOR_ILUDT} MatFactorType;
11434918c53SJed Brown extern const char *const MatFactorTypes[];
115e92e720dSBarry Smith 
1167087cfbeSBarry Smith extern PetscErrorCode  MatGetFactor(Mat,const MatSolverPackage,MatFactorType,Mat*);
1177087cfbeSBarry Smith extern PetscErrorCode  MatGetFactorAvailable(Mat,const MatSolverPackage,MatFactorType,PetscBool *);
1187087cfbeSBarry Smith extern PetscErrorCode  MatFactorGetSolverPackage(Mat,const MatSolverPackage*);
1197087cfbeSBarry Smith extern PetscErrorCode  MatGetFactorType(Mat,MatFactorType*);
1209989ab13SBarry Smith 
121c06d978dSMatthew Knepley /* Logging support */
1220700a824SBarry Smith #define    MAT_FILE_CLASSID 1211216    /* used to indicate matrices in binary files */
1237087cfbeSBarry Smith extern PetscClassId  MAT_CLASSID;
1247087cfbeSBarry Smith extern PetscClassId  MAT_FDCOLORING_CLASSID;
1257087cfbeSBarry Smith extern PetscClassId  MAT_PARTITIONING_CLASSID;
1267087cfbeSBarry Smith extern PetscClassId  MAT_NULLSPACE_CLASSID;
1277087cfbeSBarry Smith extern PetscClassId  MATMFFD_CLASSID;
128c06d978dSMatthew Knepley 
129ceb03754SKris Buschelman /*E
130ceb03754SKris Buschelman     MatReuse - Indicates if matrices obtained from a previous call to MatGetSubMatrices()
131d6eff37eSBarry Smith      or MatGetSubMatrix() are to be reused to store the new matrix values. For MatConvert() is used to indicate
132d6eff37eSBarry Smith      that the input matrix is to be replaced with the converted matrix.
133ceb03754SKris Buschelman 
134ceb03754SKris Buschelman     Level: beginner
135ceb03754SKris Buschelman 
136ceb03754SKris Buschelman    Any additions/changes here MUST also be made in include/finclude/petscmat.h
137ceb03754SKris Buschelman 
1380c451bc4SBarry Smith .seealso: MatGetSubMatrices(), MatGetSubMatrix(), MatDestroyMatrices(), MatConvert()
139ceb03754SKris Buschelman E*/
140dfe085dbSJed Brown typedef enum {MAT_INITIAL_MATRIX,MAT_REUSE_MATRIX,MAT_IGNORE_MATRIX} MatReuse;
141ceb03754SKris Buschelman 
1425494a064SHong Zhang /*E
1435494a064SHong Zhang     MatGetSubMatrixOption - Indicates if matrices obtained from a call to MatGetSubMatrices()
144829201f2SHong Zhang      include the matrix values. Currently it is only used by MatGetSeqNonzerostructure().
1455494a064SHong Zhang 
1465494a064SHong Zhang     Level: beginner
1475494a064SHong Zhang 
148829201f2SHong Zhang .seealso: MatGetSeqNonzerostructure()
1495494a064SHong Zhang E*/
1505494a064SHong Zhang typedef enum {MAT_DO_NOT_GET_VALUES,MAT_GET_VALUES} MatGetSubMatrixOption;
1515494a064SHong Zhang 
1527087cfbeSBarry Smith extern PetscErrorCode  MatInitializePackage(const char[]);
153c06d978dSMatthew Knepley 
1547087cfbeSBarry Smith extern PetscErrorCode  MatCreate(MPI_Comm,Mat*);
155a69a960dSMatthew Knepley PetscPolymorphicFunction(MatCreate,(MPI_Comm comm),(comm,&A),Mat,A)
156a69a960dSMatthew Knepley PetscPolymorphicFunction(MatCreate,(),(PETSC_COMM_WORLD,&A),Mat,A)
1577087cfbeSBarry Smith extern PetscErrorCode  MatSetSizes(Mat,PetscInt,PetscInt,PetscInt,PetscInt);
1587087cfbeSBarry Smith extern PetscErrorCode  MatSetType(Mat,const MatType);
1597087cfbeSBarry Smith extern PetscErrorCode  MatSetFromOptions(Mat);
1607087cfbeSBarry Smith extern PetscErrorCode  MatSetUpPreallocation(Mat);
1617087cfbeSBarry Smith extern PetscErrorCode  MatRegisterAll(const char[]);
1627087cfbeSBarry Smith extern PetscErrorCode  MatRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat));
16330de9b25SBarry Smith 
1647087cfbeSBarry Smith extern PetscErrorCode  MatSetOptionsPrefix(Mat,const char[]);
1657087cfbeSBarry Smith extern PetscErrorCode  MatAppendOptionsPrefix(Mat,const char[]);
1667087cfbeSBarry Smith extern PetscErrorCode  MatGetOptionsPrefix(Mat,const char*[]);
167f69a0ea3SMatthew Knepley 
16830de9b25SBarry Smith /*MC
16930de9b25SBarry Smith    MatRegisterDynamic - Adds a new matrix type
17030de9b25SBarry Smith 
17130de9b25SBarry Smith    Synopsis:
1721890ba74SBarry Smith    PetscErrorCode MatRegisterDynamic(const char *name,const char *path,const char *name_create,PetscErrorCode (*routine_create)(Mat))
17330de9b25SBarry Smith 
17430de9b25SBarry Smith    Not Collective
17530de9b25SBarry Smith 
17630de9b25SBarry Smith    Input Parameters:
17730de9b25SBarry Smith +  name - name of a new user-defined matrix type
17830de9b25SBarry Smith .  path - path (either absolute or relative) the library containing this solver
17930de9b25SBarry Smith .  name_create - name of routine to create method context
18030de9b25SBarry Smith -  routine_create - routine to create method context
18130de9b25SBarry Smith 
18230de9b25SBarry Smith    Notes:
18330de9b25SBarry Smith    MatRegisterDynamic() may be called multiple times to add several user-defined solvers.
18430de9b25SBarry Smith 
18530de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
18630de9b25SBarry Smith    is ignored.
18730de9b25SBarry Smith 
18830de9b25SBarry Smith    Sample usage:
18930de9b25SBarry Smith .vb
19030de9b25SBarry Smith    MatRegisterDynamic("my_mat",/home/username/my_lib/lib/libO/solaris/mylib.a,
19130de9b25SBarry Smith                "MyMatCreate",MyMatCreate);
19230de9b25SBarry Smith .ve
19330de9b25SBarry Smith 
19430de9b25SBarry Smith    Then, your solver can be chosen with the procedural interface via
19530de9b25SBarry Smith $     MatSetType(Mat,"my_mat")
19630de9b25SBarry Smith    or at runtime via the option
19730de9b25SBarry Smith $     -mat_type my_mat
19830de9b25SBarry Smith 
19930de9b25SBarry Smith    Level: advanced
20030de9b25SBarry Smith 
201ab901514SBarry Smith    Notes: ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values.
20230de9b25SBarry Smith          If your function is not being put into a shared library then use VecRegister() instead
20330de9b25SBarry Smith 
20430de9b25SBarry Smith .keywords: Mat, register
20530de9b25SBarry Smith 
20630de9b25SBarry Smith .seealso: MatRegisterAll(), MatRegisterDestroy()
20730de9b25SBarry Smith 
20830de9b25SBarry Smith M*/
209273d9f13SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
210273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,0)
211273d9f13SBarry Smith #else
212273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,d)
21330de9b25SBarry Smith #endif
21430de9b25SBarry Smith 
215ace3abfcSBarry Smith extern PetscBool  MatRegisterAllCalled;
216b0a32e0cSBarry Smith extern PetscFList MatList;
217b022a5c1SBarry Smith extern PetscFList MatColoringList;
218b022a5c1SBarry Smith extern PetscFList MatPartitioningList;
21928988994SBarry Smith 
2203b224e63SBarry Smith /*E
2213b224e63SBarry Smith     MatStructure - Indicates if the matrix has the same nonzero structure
2223b224e63SBarry Smith 
2233b224e63SBarry Smith     Level: beginner
2243b224e63SBarry Smith 
2253b224e63SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
2263b224e63SBarry Smith 
2273b224e63SBarry Smith .seealso: MatCopy(), KSPSetOperators(), PCSetOperators()
2283b224e63SBarry Smith E*/
2293b224e63SBarry Smith typedef enum {SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER,SUBSET_NONZERO_PATTERN} MatStructure;
2303b224e63SBarry Smith 
2317087cfbeSBarry Smith extern PetscErrorCode  MatCreateSeqDense(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*);
2327087cfbeSBarry Smith extern PetscErrorCode  MatCreateMPIDense(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*);
2337087cfbeSBarry Smith extern PetscErrorCode  MatCreateSeqAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
234ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,nz,nnz,&A),Mat,A)
235ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,&A),Mat,A)
236ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,0,nnz,&A),Mat,A)
237ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,&A),Mat,A)
238ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,A))
239ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,m,n,0,nnz,A))
240ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,A))
2417087cfbeSBarry Smith extern PetscErrorCode  MatCreateMPIAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
242ba966639SSatish 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)
243ba966639SSatish 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)
244ba966639SSatish 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)
245ba966639SSatish 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)
246ba966639SSatish 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))
247ba966639SSatish 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))
248ba966639SSatish 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))
249ba966639SSatish 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)
250ba966639SSatish 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)
251ba966639SSatish 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)
252ba966639SSatish 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)
253ba966639SSatish 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))
254ba966639SSatish 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))
255ba966639SSatish 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))
2567087cfbeSBarry Smith extern PetscErrorCode  MatCreateMPIAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *);
2577087cfbeSBarry Smith extern PetscErrorCode  MatCreateMPIAIJWithSplitArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],PetscInt[],PetscInt[],PetscScalar[],Mat*);
2588d7a6e47SBarry Smith 
2597087cfbeSBarry Smith extern PetscErrorCode  MatCreateSeqBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
260ba966639SSatish 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)
261ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A)
262ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A)
263ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A)
264ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A))
265ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A))
266ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A))
2677087cfbeSBarry Smith extern PetscErrorCode  MatCreateMPIBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
268ba966639SSatish 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)
269ba966639SSatish 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)
270ba966639SSatish 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)
271ba966639SSatish 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)
272ba966639SSatish 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))
273ba966639SSatish 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))
274ba966639SSatish 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))
275ba966639SSatish 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)
276ba966639SSatish 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)
277ba966639SSatish 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)
278ba966639SSatish 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)
279ba966639SSatish 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))
280ba966639SSatish 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))
281ba966639SSatish 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))
2827087cfbeSBarry Smith extern PetscErrorCode  MatCreateMPIBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat*);
283d21a29f3SJed Brown 
2847087cfbeSBarry Smith extern PetscErrorCode  MatCreateMPIAdj(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscInt[],Mat*);
2857087cfbeSBarry Smith extern PetscErrorCode  MatCreateSeqSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
286ba966639SSatish 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)
287ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A)
288ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A)
289ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A)
290ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A))
291ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A))
292ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A))
293d21a29f3SJed Brown 
2947087cfbeSBarry Smith extern PetscErrorCode  MatCreateMPISBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
295ba966639SSatish 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)
296ba966639SSatish 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)
297ba966639SSatish 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)
298ba966639SSatish 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)
299ba966639SSatish 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))
300ba966639SSatish 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))
301ba966639SSatish 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))
302ba966639SSatish 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)
303ba966639SSatish 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)
304ba966639SSatish 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)
305ba966639SSatish 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)
306ba966639SSatish 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))
307ba966639SSatish 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))
308ba966639SSatish 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))
3097087cfbeSBarry Smith extern PetscErrorCode  MatCreateMPISBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *);
3107087cfbeSBarry Smith extern PetscErrorCode  MatMPISBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
311dfb205c3SBarry Smith 
3127087cfbeSBarry Smith extern PetscErrorCode  MatCreateShell(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,void *,Mat*);
313ba966639SSatish 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)
314ba966639SSatish Balay PetscPolymorphicFunction(MatCreateShell,(PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx),(PETSC_COMM_WORLD,m,n,M,N,ctx,&A),Mat,A)
3157087cfbeSBarry Smith extern PetscErrorCode  MatCreateAdic(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,void (*)(void),Mat*);
3167087cfbeSBarry Smith extern PetscErrorCode  MatCreateNormal(Mat,Mat*);
317ba966639SSatish Balay PetscPolymorphicFunction(MatCreateNormal,(Mat mat),(mat,&A),Mat,A)
3187087cfbeSBarry Smith extern PetscErrorCode  MatCreateLRC(Mat,Mat,Mat,Mat*);
3197087cfbeSBarry Smith extern PetscErrorCode  MatCreateIS(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,ISLocalToGlobalMapping,Mat*);
3207087cfbeSBarry Smith extern PetscErrorCode  MatCreateSeqAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
3217087cfbeSBarry Smith extern PetscErrorCode  MatCreateMPIAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
3227087cfbeSBarry Smith extern PetscErrorCode  MatCreateScatter(MPI_Comm,VecScatter,Mat*);
3237087cfbeSBarry Smith extern PetscErrorCode  MatScatterSetVecScatter(Mat,VecScatter);
3247087cfbeSBarry Smith extern PetscErrorCode  MatScatterGetVecScatter(Mat,VecScatter*);
3257087cfbeSBarry Smith extern PetscErrorCode  MatCreateBlockMat(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,Mat*);
3267087cfbeSBarry Smith extern PetscErrorCode  MatCompositeAddMat(Mat,Mat);
3277087cfbeSBarry Smith extern PetscErrorCode  MatCompositeMerge(Mat);
3287087cfbeSBarry Smith extern PetscErrorCode  MatCreateComposite(MPI_Comm,PetscInt,const Mat*,Mat*);
3296d7c1e57SBarry Smith typedef enum {MAT_COMPOSITE_ADDITIVE,MAT_COMPOSITE_MULTIPLICATIVE} MatCompositeType;
3307087cfbeSBarry Smith extern PetscErrorCode  MatCompositeSetType(Mat,MatCompositeType);
3316d7c1e57SBarry Smith 
332*dedccee8SHong Zhang extern PetscErrorCode  MatCreateFFT(MPI_Comm,PetscInt,const PetscInt[],const MatType,Mat*);
3337087cfbeSBarry Smith extern PetscErrorCode  MatCreateSeqCUFFT(MPI_Comm,PetscInt,const PetscInt[],Mat*);
334*dedccee8SHong Zhang 
3357087cfbeSBarry Smith extern PetscErrorCode  MatCreateTranspose(Mat,Mat*);
3367087cfbeSBarry Smith extern PetscErrorCode  MatCreateSubMatrix(Mat,IS,IS,Mat*);
3377087cfbeSBarry Smith extern PetscErrorCode  MatSubMatrixUpdate(Mat,Mat,IS,IS);
3387087cfbeSBarry Smith extern PetscErrorCode  MatCreateLocalRef(Mat,IS,IS,Mat*);
3391472f72bSBarry Smith 
3407087cfbeSBarry Smith extern PetscErrorCode  MatCreatePython(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const char[],Mat*);
3417087cfbeSBarry Smith extern PetscErrorCode  MatPythonSetType(Mat,const char[]);
3421d6018f0SLisandro Dalcin 
3437087cfbeSBarry Smith extern PetscErrorCode  MatSetUp(Mat);
3447087cfbeSBarry Smith extern PetscErrorCode  MatDestroy(Mat);
34521c89e3eSBarry Smith 
3467087cfbeSBarry Smith extern PetscErrorCode  MatConjugate(Mat);
3477087cfbeSBarry Smith extern PetscErrorCode  MatRealPart(Mat);
3487087cfbeSBarry Smith extern PetscErrorCode  MatImaginaryPart(Mat);
3497087cfbeSBarry Smith extern PetscErrorCode  MatGetDiagonalBlock(Mat,PetscBool *,MatReuse,Mat*);
3507087cfbeSBarry Smith extern PetscErrorCode  MatGetTrace(Mat,PetscScalar*);
35199cafbc1SBarry Smith 
3528ed539a5SBarry Smith /* ------------------------------------------------------------*/
3537087cfbeSBarry Smith extern PetscErrorCode  MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
3547087cfbeSBarry Smith extern PetscErrorCode  MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
3557087cfbeSBarry Smith extern PetscErrorCode  MatSetValuesRow(Mat,PetscInt,const PetscScalar[]);
3567087cfbeSBarry Smith extern PetscErrorCode  MatSetValuesRowLocal(Mat,PetscInt,const PetscScalar[]);
35784cb2905SBarry Smith 
3582ef4de8bSBarry Smith /*S
3592ef4de8bSBarry Smith      MatStencil - Data structure (C struct) for storing information about a single row or
3602ef4de8bSBarry Smith         column of a matrix as index on an associated grid.
3612ef4de8bSBarry Smith 
3622ef4de8bSBarry Smith    Level: beginner
3632ef4de8bSBarry Smith 
3642ef4de8bSBarry Smith   Concepts: matrix; linear operator
3652ef4de8bSBarry Smith 
366d7bd93baSBarry Smith .seealso:  MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockStencil()
3672ef4de8bSBarry Smith S*/
368435da068SBarry Smith typedef struct {
369c1ac3661SBarry Smith   PetscInt k,j,i,c;
370435da068SBarry Smith } MatStencil;
3712ef4de8bSBarry Smith 
3727087cfbeSBarry Smith extern PetscErrorCode  MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
3737087cfbeSBarry Smith extern PetscErrorCode  MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
3747087cfbeSBarry Smith extern PetscErrorCode  MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt);
375435da068SBarry Smith 
3767087cfbeSBarry Smith extern PetscErrorCode  MatSetColoring(Mat,ISColoring);
3777087cfbeSBarry Smith extern PetscErrorCode  MatSetValuesAdic(Mat,void*);
3787087cfbeSBarry Smith extern PetscErrorCode  MatSetValuesAdifor(Mat,PetscInt,void*);
3793a7fca6bSBarry Smith 
380d91e6319SBarry Smith /*E
381d91e6319SBarry Smith     MatAssemblyType - Indicates if the matrix is now to be used, or if you plan
382d91e6319SBarry Smith      to continue to add values to it
383d91e6319SBarry Smith 
384d91e6319SBarry Smith     Level: beginner
385d91e6319SBarry Smith 
386d91e6319SBarry Smith .seealso: MatAssemblyBegin(), MatAssemblyEnd()
387d91e6319SBarry Smith E*/
3886d4a8577SBarry Smith typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType;
3897087cfbeSBarry Smith extern PetscErrorCode  MatAssemblyBegin(Mat,MatAssemblyType);
3907087cfbeSBarry Smith extern PetscErrorCode  MatAssemblyEnd(Mat,MatAssemblyType);
3917087cfbeSBarry Smith extern PetscErrorCode  MatAssembled(Mat,PetscBool *);
3924f9c727eSBarry Smith 
3931d73ed98SBarry Smith 
39430de9b25SBarry Smith 
395d91e6319SBarry Smith /*E
396d91e6319SBarry Smith     MatOption - Options that may be set for a matrix and its behavior or storage
397d91e6319SBarry Smith 
398d91e6319SBarry Smith     Level: beginner
399d91e6319SBarry Smith 
4000a835dfdSSatish Balay    Any additions/changes here MUST also be made in include/finclude/petscmat.h
401d91e6319SBarry Smith 
402d91e6319SBarry Smith .seealso: MatSetOption()
403d91e6319SBarry Smith E*/
404512a5fc5SBarry Smith typedef enum {MAT_ROW_ORIENTED,MAT_NEW_NONZERO_LOCATIONS,
4054e0d8c25SBarry Smith               MAT_SYMMETRIC,
4064e0d8c25SBarry Smith               MAT_STRUCTURALLY_SYMMETRIC,
4078047e77bSBarry Smith               MAT_NEW_DIAGONALS,MAT_IGNORE_OFF_PROC_ENTRIES,
4084e0d8c25SBarry Smith               MAT_NEW_NONZERO_LOCATION_ERR,
4094e0d8c25SBarry Smith               MAT_NEW_NONZERO_ALLOCATION_ERR,MAT_USE_HASH_TABLE,
410a9817697SBarry Smith               MAT_KEEP_NONZERO_PATTERN,MAT_IGNORE_ZERO_ENTRIES,
4114e0d8c25SBarry Smith               MAT_USE_INODES,
4124e0d8c25SBarry Smith               MAT_HERMITIAN,
4134e0d8c25SBarry Smith               MAT_SYMMETRY_ETERNAL,
414cd6b891eSBarry Smith               MAT_CHECK_COMPRESSED_ROW,
4154e0d8c25SBarry Smith               MAT_IGNORE_LOWER_TRIANGULAR,MAT_ERROR_LOWER_TRIANGULAR,
41628b2fa4aSMatthew Knepley               MAT_GETROW_UPPERTRIANGULAR,MAT_UNUSED_NONZERO_LOCATION_ERR,
4174cb17eb5SBarry Smith               MAT_SPD,MAT_NO_OFF_PROC_ENTRIES,MAT_NO_OFF_PROC_ZERO_ROWS,
41828b2fa4aSMatthew Knepley               NUM_MAT_OPTIONS} MatOption;
419290bbb0aSBarry Smith extern const char *MatOptions[];
4207087cfbeSBarry Smith extern PetscErrorCode  MatSetOption(Mat,MatOption,PetscBool );
4217087cfbeSBarry Smith extern PetscErrorCode  MatGetType(Mat,const MatType*);
422a313700dSBarry Smith PetscPolymorphicFunction(MatGetType,(Mat mat),(mat,&t),const MatType,t)
42384cb2905SBarry Smith 
4247087cfbeSBarry Smith extern PetscErrorCode  MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]);
4257087cfbeSBarry Smith extern PetscErrorCode  MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
4267087cfbeSBarry Smith extern PetscErrorCode  MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
4277087cfbeSBarry Smith extern PetscErrorCode  MatGetRowUpperTriangular(Mat);
4287087cfbeSBarry Smith extern PetscErrorCode  MatRestoreRowUpperTriangular(Mat);
4297087cfbeSBarry Smith extern PetscErrorCode  MatGetColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
4307087cfbeSBarry Smith extern PetscErrorCode  MatRestoreColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
4317087cfbeSBarry Smith extern PetscErrorCode  MatGetColumnVector(Mat,Vec,PetscInt);
4327087cfbeSBarry Smith extern PetscErrorCode  MatGetArray(Mat,PetscScalar *[]);
433ba966639SSatish Balay PetscPolymorphicFunction(MatGetArray,(Mat mat),(mat,&a),PetscScalar*,a)
4347087cfbeSBarry Smith extern PetscErrorCode  MatRestoreArray(Mat,PetscScalar *[]);
4357087cfbeSBarry Smith extern PetscErrorCode  MatGetBlockSize(Mat,PetscInt *);
436ba966639SSatish Balay PetscPolymorphicFunction(MatGetBlockSize,(Mat mat),(mat,&a),PetscInt,a)
4377087cfbeSBarry Smith extern PetscErrorCode  MatSetBlockSize(Mat,PetscInt);
4387b80b807SBarry Smith 
4391620fd73SBarry Smith 
4407087cfbeSBarry Smith extern PetscErrorCode  MatMult(Mat,Vec,Vec);
4417087cfbeSBarry Smith extern PetscErrorCode  MatMultDiagonalBlock(Mat,Vec,Vec);
4427087cfbeSBarry Smith extern PetscErrorCode  MatMultAdd(Mat,Vec,Vec,Vec);
443ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
4447087cfbeSBarry Smith extern PetscErrorCode  MatMultTranspose(Mat,Vec,Vec);
4457087cfbeSBarry Smith extern PetscErrorCode  MatMultHermitianTranspose(Mat,Vec,Vec);
4467087cfbeSBarry Smith extern PetscErrorCode  MatIsTranspose(Mat,Mat,PetscReal,PetscBool *);
447ace3abfcSBarry Smith PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B,PetscReal tol),(A,B,tol,&t),PetscBool ,t)
448ace3abfcSBarry Smith PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B),(A,B,0,&t),PetscBool ,t)
4497087cfbeSBarry Smith extern PetscErrorCode  MatIsHermitianTranspose(Mat,Mat,PetscReal,PetscBool *);
4507087cfbeSBarry Smith extern PetscErrorCode  MatMultTransposeAdd(Mat,Vec,Vec,Vec);
4517087cfbeSBarry Smith extern PetscErrorCode  MatMultHermitianTransposeAdd(Mat,Vec,Vec,Vec);
452ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultTransposeAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
4537087cfbeSBarry Smith extern PetscErrorCode  MatMultConstrained(Mat,Vec,Vec);
4547087cfbeSBarry Smith extern PetscErrorCode  MatMultTransposeConstrained(Mat,Vec,Vec);
4557087cfbeSBarry Smith extern PetscErrorCode  MatMatSolve(Mat,Mat,Mat);
456f5edf698SHong Zhang 
457d91e6319SBarry Smith /*E
458d91e6319SBarry Smith     MatDuplicateOption - Indicates if a duplicated sparse matrix should have
459d91e6319SBarry Smith   its numerical values copied over or just its nonzero structure.
460d91e6319SBarry Smith 
461d91e6319SBarry Smith     Level: beginner
462d91e6319SBarry Smith 
463d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
464d91e6319SBarry Smith 
46570dcbbb9SBarry Smith $   MAT_SHARE_NONZERO_PATTERN - the i and j arrays in the new matrix will be shared with the original matrix
46670dcbbb9SBarry Smith $                               this also triggers the MAT_DO_NOT_COPY_VALUES option. This is used when you
46770dcbbb9SBarry Smith $                               have several matrices with the same nonzero pattern.
46870dcbbb9SBarry Smith 
469d91e6319SBarry Smith .seealso: MatDuplicate()
470d91e6319SBarry Smith E*/
47170dcbbb9SBarry Smith typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES,MAT_SHARE_NONZERO_PATTERN} MatDuplicateOption;
4722e8a6d31SBarry Smith 
4737087cfbeSBarry Smith extern PetscErrorCode  MatConvert(Mat,const MatType,MatReuse,Mat*);
474a313700dSBarry Smith PetscPolymorphicFunction(MatConvert,(Mat A,const MatType t),(A,t,MAT_INITIAL_MATRIX,&a),Mat,a)
4757087cfbeSBarry Smith extern PetscErrorCode  MatDuplicate(Mat,MatDuplicateOption,Mat*);
476ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A,MatDuplicateOption o),(A,o,&a),Mat,a)
477ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A),(A,MAT_COPY_VALUES,&a),Mat,a)
47894a9d846SBarry Smith 
479cb5b572fSBarry Smith 
4807087cfbeSBarry Smith extern PetscErrorCode  MatCopy(Mat,Mat,MatStructure);
4817087cfbeSBarry Smith extern PetscErrorCode  MatView(Mat,PetscViewer);
4827087cfbeSBarry Smith extern PetscErrorCode  MatIsSymmetric(Mat,PetscReal,PetscBool *);
483ace3abfcSBarry Smith PetscPolymorphicFunction(MatIsSymmetric,(Mat A,PetscReal tol),(A,tol,&t),PetscBool ,t)
484ace3abfcSBarry Smith PetscPolymorphicFunction(MatIsSymmetric,(Mat A),(A,0,&t),PetscBool ,t)
4857087cfbeSBarry Smith extern PetscErrorCode  MatIsStructurallySymmetric(Mat,PetscBool *);
486ace3abfcSBarry Smith PetscPolymorphicFunction(MatIsStructurallySymmetric,(Mat A),(A,&t),PetscBool ,t)
4877087cfbeSBarry Smith extern PetscErrorCode  MatIsHermitian(Mat,PetscReal,PetscBool *);
488ace3abfcSBarry Smith PetscPolymorphicFunction(MatIsHermitian,(Mat A),(A,0,&t),PetscBool ,t)
4897087cfbeSBarry Smith extern PetscErrorCode  MatIsSymmetricKnown(Mat,PetscBool *,PetscBool *);
4907087cfbeSBarry Smith extern PetscErrorCode  MatIsHermitianKnown(Mat,PetscBool *,PetscBool *);
4917087cfbeSBarry Smith extern PetscErrorCode  MatMissingDiagonal(Mat,PetscBool  *,PetscInt *);
4927087cfbeSBarry Smith extern PetscErrorCode  MatLoad(Mat, PetscViewer);
4937b80b807SBarry Smith 
4947087cfbeSBarry Smith extern PetscErrorCode  MatGetRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,PetscInt *[],PetscInt *[],PetscBool  *);
4957087cfbeSBarry Smith extern PetscErrorCode  MatRestoreRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,PetscInt *[],PetscInt *[],PetscBool  *);
4967087cfbeSBarry Smith extern PetscErrorCode  MatGetColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,PetscInt *[],PetscInt *[],PetscBool  *);
4977087cfbeSBarry Smith extern PetscErrorCode  MatRestoreColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,PetscInt *[],PetscInt *[],PetscBool  *);
498d4fbbf0eSBarry Smith 
499d91e6319SBarry Smith /*S
500d91e6319SBarry Smith      MatInfo - Context of matrix information, used with MatGetInfo()
501d91e6319SBarry Smith 
502d91e6319SBarry Smith    In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE
503d91e6319SBarry Smith 
504d91e6319SBarry Smith    Level: intermediate
505d91e6319SBarry Smith 
506d91e6319SBarry Smith   Concepts: matrix^nonzero information
507d91e6319SBarry Smith 
508d9274352SBarry Smith .seealso:  MatGetInfo(), MatInfoType
509d91e6319SBarry Smith S*/
5104e220ebcSLois Curfman McInnes typedef struct {
511b0a32e0cSBarry Smith   PetscLogDouble block_size;                         /* block size */
512b0a32e0cSBarry Smith   PetscLogDouble nz_allocated,nz_used,nz_unneeded;   /* number of nonzeros */
513b0a32e0cSBarry Smith   PetscLogDouble memory;                             /* memory allocated */
514b0a32e0cSBarry Smith   PetscLogDouble assemblies;                         /* number of matrix assemblies called */
515b0a32e0cSBarry Smith   PetscLogDouble mallocs;                            /* number of mallocs during MatSetValues() */
516b0a32e0cSBarry Smith   PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */
517b0a32e0cSBarry Smith   PetscLogDouble factor_mallocs;                     /* number of mallocs during factorization */
5184e220ebcSLois Curfman McInnes } MatInfo;
5194e220ebcSLois Curfman McInnes 
520d9274352SBarry Smith /*E
521d9274352SBarry Smith     MatInfoType - Indicates if you want information about the local part of the matrix,
522d9274352SBarry Smith      the entire parallel matrix or the maximum over all the local parts.
523d9274352SBarry Smith 
524d9274352SBarry Smith     Level: beginner
525d9274352SBarry Smith 
526d9274352SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
527d9274352SBarry Smith 
528d9274352SBarry Smith .seealso: MatGetInfo(), MatInfo
529d9274352SBarry Smith E*/
5307b80b807SBarry Smith typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType;
5317087cfbeSBarry Smith extern PetscErrorCode  MatGetInfo(Mat,MatInfoType,MatInfo*);
5327087cfbeSBarry Smith extern PetscErrorCode  MatGetDiagonal(Mat,Vec);
5337087cfbeSBarry Smith extern PetscErrorCode  MatGetRowMax(Mat,Vec,PetscInt[]);
5347087cfbeSBarry Smith extern PetscErrorCode  MatGetRowMin(Mat,Vec,PetscInt[]);
5357087cfbeSBarry Smith extern PetscErrorCode  MatGetRowMaxAbs(Mat,Vec,PetscInt[]);
5367087cfbeSBarry Smith extern PetscErrorCode  MatGetRowMinAbs(Mat,Vec,PetscInt[]);
5377087cfbeSBarry Smith extern PetscErrorCode  MatGetRowSum(Mat,Vec);
5387087cfbeSBarry Smith extern PetscErrorCode  MatTranspose(Mat,MatReuse,Mat*);
539fc4dec0aSBarry Smith PetscPolymorphicFunction(MatTranspose,(Mat A),(A,MAT_INITIAL_MATRIX,&t),Mat,t)
5407087cfbeSBarry Smith extern PetscErrorCode  MatHermitianTranspose(Mat,MatReuse,Mat*);
5417087cfbeSBarry Smith extern PetscErrorCode  MatPermute(Mat,IS,IS,Mat *);
542ba966639SSatish Balay PetscPolymorphicFunction(MatPermute,(Mat A,IS is1,IS is2),(A,is1,is2,&t),Mat,t)
5437087cfbeSBarry Smith extern PetscErrorCode  MatDiagonalScale(Mat,Vec,Vec);
5447087cfbeSBarry Smith extern PetscErrorCode  MatDiagonalSet(Mat,Vec,InsertMode);
5457087cfbeSBarry Smith extern PetscErrorCode  MatEqual(Mat,Mat,PetscBool *);
546ace3abfcSBarry Smith PetscPolymorphicFunction(MatEqual,(Mat A,Mat B),(A,B,&t),PetscBool ,t)
5477087cfbeSBarry Smith extern PetscErrorCode  MatMultEqual(Mat,Mat,PetscInt,PetscBool *);
5487087cfbeSBarry Smith extern PetscErrorCode  MatMultAddEqual(Mat,Mat,PetscInt,PetscBool *);
5497087cfbeSBarry Smith extern PetscErrorCode  MatMultTransposeEqual(Mat,Mat,PetscInt,PetscBool *);
5507087cfbeSBarry Smith extern PetscErrorCode  MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscBool *);
5517b80b807SBarry Smith 
5527087cfbeSBarry Smith extern PetscErrorCode  MatNorm(Mat,NormType,PetscReal *);
553ba966639SSatish Balay PetscPolymorphicFunction(MatNorm,(Mat A,NormType t),(A,t,&n),PetscReal,n)
55409573ac7SBarry Smith extern PetscErrorCode MatGetColumnNorms(Mat,NormType,PetscReal *);
5557087cfbeSBarry Smith extern PetscErrorCode  MatZeroEntries(Mat);
5567087cfbeSBarry Smith extern PetscErrorCode  MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
5577087cfbeSBarry Smith extern PetscErrorCode  MatZeroRowsIS(Mat,IS,PetscScalar,Vec,Vec);
5587087cfbeSBarry Smith extern PetscErrorCode  MatZeroRowsStencil(Mat,PetscInt,const MatStencil [],PetscScalar,Vec,Vec);
5597087cfbeSBarry Smith extern PetscErrorCode  MatZeroRowsColumns(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
5607087cfbeSBarry Smith extern PetscErrorCode  MatZeroRowsColumnsIS(Mat,IS,PetscScalar,Vec,Vec);
5617b80b807SBarry Smith 
5627087cfbeSBarry Smith extern PetscErrorCode  MatUseScaledForm(Mat,PetscBool );
5637087cfbeSBarry Smith extern PetscErrorCode  MatScaleSystem(Mat,Vec,Vec);
5647087cfbeSBarry Smith extern PetscErrorCode  MatUnScaleSystem(Mat,Vec,Vec);
5655ef9f2a5SBarry Smith 
5667087cfbeSBarry Smith extern PetscErrorCode  MatGetSize(Mat,PetscInt*,PetscInt*);
5677087cfbeSBarry Smith extern PetscErrorCode  MatGetLocalSize(Mat,PetscInt*,PetscInt*);
5687087cfbeSBarry Smith extern PetscErrorCode  MatGetOwnershipRange(Mat,PetscInt*,PetscInt*);
5697087cfbeSBarry Smith extern PetscErrorCode  MatGetOwnershipRanges(Mat,const PetscInt**);
5707087cfbeSBarry Smith extern PetscErrorCode  MatGetOwnershipRangeColumn(Mat,PetscInt*,PetscInt*);
5717087cfbeSBarry Smith extern PetscErrorCode  MatGetOwnershipRangesColumn(Mat,const PetscInt**);
5727b80b807SBarry Smith 
5737087cfbeSBarry Smith extern PetscErrorCode  MatGetSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
5747087cfbeSBarry Smith extern PetscErrorCode  MatGetSubMatricesParallel(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
5757087cfbeSBarry Smith extern PetscErrorCode  MatDestroyMatrices(PetscInt,Mat *[]);
5767087cfbeSBarry Smith extern PetscErrorCode  MatGetSubMatrix(Mat,IS,IS,MatReuse,Mat *);
5777087cfbeSBarry Smith extern PetscErrorCode  MatGetLocalSubMatrix(Mat,IS,IS,Mat*);
5787087cfbeSBarry Smith extern PetscErrorCode  MatRestoreLocalSubMatrix(Mat,IS,IS,Mat*);
5797087cfbeSBarry Smith extern PetscErrorCode  MatGetSeqNonzeroStructure(Mat,Mat*);
5807087cfbeSBarry Smith extern PetscErrorCode  MatDestroySeqNonzeroStructure(Mat*);
5815494a064SHong Zhang 
5827087cfbeSBarry Smith extern PetscErrorCode  MatMerge(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
5837087cfbeSBarry Smith extern PetscErrorCode  MatMerge_SeqsToMPI(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*);
5847087cfbeSBarry Smith extern PetscErrorCode  MatMerge_SeqsToMPISymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*);
5857087cfbeSBarry Smith extern PetscErrorCode  MatMerge_SeqsToMPINumeric(Mat,Mat);
5867087cfbeSBarry Smith extern PetscErrorCode  MatDestroy_MPIAIJ_SeqsToMPI(Mat);
5877087cfbeSBarry Smith extern PetscErrorCode  MatGetLocalMat(Mat,MatReuse,Mat*);
5887087cfbeSBarry Smith extern PetscErrorCode  MatGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*);
5897087cfbeSBarry Smith extern PetscErrorCode  MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,PetscInt*,Mat*);
5907087cfbeSBarry Smith extern PetscErrorCode  MatGetBrowsOfAoCols(Mat,Mat,MatReuse,PetscInt**,PetscInt**,MatScalar**,Mat*);
59143eb5e2fSMatthew Knepley #if defined (PETSC_USE_CTABLE)
592cd116e26SSatish Balay #include "petscctable.h"
5937087cfbeSBarry Smith extern PetscErrorCode  MatGetCommunicationStructs(Mat, Vec *, PetscTable *, VecScatter *);
59443eb5e2fSMatthew Knepley #else
5957087cfbeSBarry Smith extern PetscErrorCode  MatGetCommunicationStructs(Mat, Vec *, PetscInt *[], VecScatter *);
59643eb5e2fSMatthew Knepley #endif
5977087cfbeSBarry Smith extern PetscErrorCode  MatGetGhosts(Mat, PetscInt *,const PetscInt *[]);
5988efafbd8SBarry Smith 
5997087cfbeSBarry Smith extern PetscErrorCode  MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt);
6007b80b807SBarry Smith 
6017087cfbeSBarry Smith extern PetscErrorCode  MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*);
6027087cfbeSBarry Smith extern PetscErrorCode  MatMatMultSymbolic(Mat,Mat,PetscReal,Mat*);
6037087cfbeSBarry Smith extern PetscErrorCode  MatMatMultNumeric(Mat,Mat,Mat);
60422440eb1SKris Buschelman 
6057087cfbeSBarry Smith extern PetscErrorCode  MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*);
6067087cfbeSBarry Smith extern PetscErrorCode  MatPtAPSymbolic(Mat,Mat,PetscReal,Mat*);
6077087cfbeSBarry Smith extern PetscErrorCode  MatPtAPNumeric(Mat,Mat,Mat);
60822440eb1SKris Buschelman 
6097087cfbeSBarry Smith extern PetscErrorCode  MatMatMultTranspose(Mat,Mat,MatReuse,PetscReal,Mat*);
6107087cfbeSBarry Smith extern PetscErrorCode  MatMatMultTransposeSymbolic(Mat,Mat,PetscReal,Mat*);
6117087cfbeSBarry Smith extern PetscErrorCode  MatMatMultTransposeNumeric(Mat,Mat,Mat);
612bc011b1eSHong Zhang 
6137087cfbeSBarry Smith extern PetscErrorCode  MatAXPY(Mat,PetscScalar,Mat,MatStructure);
6147087cfbeSBarry Smith extern PetscErrorCode  MatAYPX(Mat,PetscScalar,Mat,MatStructure);
6151c741599SHong Zhang 
6167087cfbeSBarry Smith extern PetscErrorCode  MatScale(Mat,PetscScalar);
6177087cfbeSBarry Smith extern PetscErrorCode  MatShift(Mat,PetscScalar);
6187b80b807SBarry Smith 
6197087cfbeSBarry Smith extern PetscErrorCode  MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping);
6207087cfbeSBarry Smith extern PetscErrorCode  MatSetLocalToGlobalMappingBlock(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping);
6217087cfbeSBarry Smith extern PetscErrorCode  MatGetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*);
6227087cfbeSBarry Smith extern PetscErrorCode  MatGetLocalToGlobalMappingBlock(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*);
6237087cfbeSBarry Smith extern PetscErrorCode  MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
6247087cfbeSBarry Smith extern PetscErrorCode  MatZeroRowsLocalIS(Mat,IS,PetscScalar,Vec,Vec);
6257087cfbeSBarry Smith extern PetscErrorCode  MatZeroRowsColumnsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
6267087cfbeSBarry Smith extern PetscErrorCode  MatZeroRowsColumnsLocalIS(Mat,IS,PetscScalar,Vec,Vec);
6277087cfbeSBarry Smith extern PetscErrorCode  MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
6287087cfbeSBarry Smith extern PetscErrorCode  MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
629052efed2SBarry Smith 
6307087cfbeSBarry Smith extern PetscErrorCode  MatStashSetInitialSize(Mat,PetscInt,PetscInt);
6317087cfbeSBarry Smith extern PetscErrorCode  MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
63290f02eecSBarry Smith 
6337087cfbeSBarry Smith extern PetscErrorCode  MatInterpolate(Mat,Vec,Vec);
6347087cfbeSBarry Smith extern PetscErrorCode  MatInterpolateAdd(Mat,Vec,Vec,Vec);
635ba966639SSatish Balay PetscPolymorphicSubroutine(MatInterpolateAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
6367087cfbeSBarry Smith extern PetscErrorCode  MatRestrict(Mat,Vec,Vec);
6377087cfbeSBarry Smith extern PetscErrorCode  MatGetVecs(Mat,Vec*,Vec*);
6387087cfbeSBarry Smith extern PetscErrorCode  MatGetRedundantMatrix(Mat,PetscInt,MPI_Comm,PetscInt,MatReuse,Mat*);
6397087cfbeSBarry Smith extern PetscErrorCode  MatGetMultiProcBlock(Mat,MPI_Comm,Mat*);
640fc9bc008SSatish Balay extern PetscErrorCode  MatFindZeroDiagonals(Mat,IS*);
641bd481603SBarry Smith 
642bd481603SBarry Smith /*MC
6431620fd73SBarry Smith    MatSetValue - Set a single entry into a matrix.
6441620fd73SBarry Smith 
6451620fd73SBarry Smith    Not collective
6461620fd73SBarry Smith 
6471620fd73SBarry Smith    Input Parameters:
6481620fd73SBarry Smith +  m - the matrix
6491620fd73SBarry Smith .  row - the row location of the entry
6501620fd73SBarry Smith .  col - the column location of the entry
6511620fd73SBarry Smith .  value - the value to insert
6521620fd73SBarry Smith -  mode - either INSERT_VALUES or ADD_VALUES
6531620fd73SBarry Smith 
6541620fd73SBarry Smith    Notes:
6551620fd73SBarry Smith    For efficiency one should use MatSetValues() and set several or many
6561620fd73SBarry Smith    values simultaneously if possible.
6571620fd73SBarry Smith 
6581620fd73SBarry Smith    Level: beginner
6591620fd73SBarry Smith 
6601620fd73SBarry Smith .seealso: MatSetValues(), MatSetValueLocal()
6611620fd73SBarry Smith M*/
6621620fd73SBarry 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);}
6631620fd73SBarry Smith 
6641620fd73SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatGetValue(Mat v,PetscInt i,PetscInt j,PetscScalar *va) {return MatGetValues(v,1,&i,1,&j,va);}
6651620fd73SBarry Smith 
6661620fd73SBarry 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);}
6671620fd73SBarry Smith 
6681620fd73SBarry Smith /*MC
6690d2a7ff2SSatish Balay    MatPreallocateInitialize - Begins the block of code that will count the number of nonzeros per
670bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
671bd481603SBarry Smith 
672bd481603SBarry Smith    Synopsis:
673c1ac3661SBarry Smith    PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)
674bd481603SBarry Smith 
675bd481603SBarry Smith    Collective on MPI_Comm
676bd481603SBarry Smith 
677bd481603SBarry Smith    Input Parameters:
678bd481603SBarry Smith +  comm - the communicator that will share the eventually allocated matrix
679859fcb39SBarry Smith .  nrows - the number of LOCAL rows in the matrix
680859fcb39SBarry Smith -  ncols - the number of LOCAL columns in the matrix
681bd481603SBarry Smith 
682bd481603SBarry Smith    Output Parameters:
683bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
684bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
685bd481603SBarry Smith 
686bd481603SBarry Smith 
687bd481603SBarry Smith    Level: intermediate
688bd481603SBarry Smith 
689bd481603SBarry Smith    Notes:
6900598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
691bd481603SBarry Smith 
6921d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
693bd481603SBarry Smith 
694bd481603SBarry Smith    Use MatPreallocateInitializeSymmetric() for symmetric matrices (MPISBAIJ matrices)
695bd481603SBarry Smith 
6961620fd73SBarry Smith    This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize().
6971620fd73SBarry Smith 
698bd481603SBarry Smith   Concepts: preallocation^Matrix
699bd481603SBarry Smith 
700bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
701bd481603SBarry Smith           MatPreallocateInitializeSymmetric(), MatPreallocateSymmetricSetLocal()
702bd481603SBarry Smith M*/
703c4f061fbSSatish Balay #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \
7047c922b88SBarry Smith { \
705a89bc211SBarry Smith   PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__start,__end; \
706a89bc211SBarry Smith   _4_ierr = PetscMalloc2(__nrows,PetscInt,&dnz,__nrows,PetscInt,&onz);CHKERRQ(_4_ierr); \
707a89bc211SBarry Smith   _4_ierr = PetscMemzero(dnz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
708a89bc211SBarry Smith   _4_ierr = PetscMemzero(onz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
709c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ctmp;\
710a89bc211SBarry Smith   _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows;
7117c922b88SBarry Smith 
712bd481603SBarry Smith /*MC
7130d2a7ff2SSatish Balay    MatPreallocateSymmetricInitialize - Begins the block of code that will count the number of nonzeros per
714bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
715bd481603SBarry Smith 
716bd481603SBarry Smith    Synopsis:
717c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)
718bd481603SBarry Smith 
719bd481603SBarry Smith    Collective on MPI_Comm
720bd481603SBarry Smith 
721bd481603SBarry Smith    Input Parameters:
722bd481603SBarry Smith +  comm - the communicator that will share the eventually allocated matrix
723859fcb39SBarry Smith .  nrows - the number of LOCAL rows in the matrix
724859fcb39SBarry Smith -  ncols - the number of LOCAL columns in the matrix
725bd481603SBarry Smith 
726bd481603SBarry Smith    Output Parameters:
727bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
728bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
729bd481603SBarry Smith 
730bd481603SBarry Smith 
731bd481603SBarry Smith    Level: intermediate
732bd481603SBarry Smith 
733bd481603SBarry Smith    Notes:
7340598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
735bd481603SBarry Smith 
7361d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
737bd481603SBarry Smith 
7381620fd73SBarry Smith    This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize().
7391620fd73SBarry Smith 
740bd481603SBarry Smith   Concepts: preallocation^Matrix
741bd481603SBarry Smith 
742bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
743bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal()
744bd481603SBarry Smith M*/
745222b16d4SBarry Smith #define MatPreallocateSymmetricInitialize(comm,nrows,ncols,dnz,onz) 0; \
746222b16d4SBarry Smith { \
747a89bc211SBarry Smith   PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__end; \
748a89bc211SBarry Smith   _4_ierr = PetscMalloc2(__nrows,PetscInt,&dnz,__nrows,PetscInt,&onz);CHKERRQ(_4_ierr); \
749a89bc211SBarry Smith   _4_ierr = PetscMemzero(dnz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
750a89bc211SBarry Smith   _4_ierr = PetscMemzero(onz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
751c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr);\
752a89bc211SBarry Smith   _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows;
753222b16d4SBarry Smith 
754bd481603SBarry Smith /*MC
755bd481603SBarry Smith    MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
756bd481603SBarry Smith        inserted using a local number of the rows and columns
757bd481603SBarry Smith 
758bd481603SBarry Smith    Synopsis:
759c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
760bd481603SBarry Smith 
761bd481603SBarry Smith    Not Collective
762bd481603SBarry Smith 
763bd481603SBarry Smith    Input Parameters:
764784ac674SJed Brown +  map - the row mapping from local numbering to global numbering
765bd481603SBarry Smith .  nrows - the number of rows indicated
7661d73ed98SBarry Smith .  rows - the indices of the rows
767784ac674SJed Brown .  cmap - the column mapping from local to global numbering
768bd481603SBarry Smith .  ncols - the number of columns in the matrix
769bd481603SBarry Smith .  cols - the columns indicated
770bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
771bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
772bd481603SBarry Smith 
773bd481603SBarry Smith 
774bd481603SBarry Smith    Level: intermediate
775bd481603SBarry Smith 
776bd481603SBarry Smith    Notes:
7770598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
778bd481603SBarry Smith 
7791d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
780bd481603SBarry Smith 
781bd481603SBarry Smith   Concepts: preallocation^Matrix
782bd481603SBarry Smith 
783bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
784bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal()
785bd481603SBarry Smith M*/
786784ac674SJed Brown #define MatPreallocateSetLocal(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \
787c4f061fbSSatish Balay {\
788c1ac3661SBarry Smith   PetscInt __l;\
789784ac674SJed Brown   _4_ierr = ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\
790784ac674SJed Brown   _4_ierr = ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\
791c4f061fbSSatish Balay   for (__l=0;__l<nrows;__l++) {\
792ef66eb69SBarry Smith     _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
793c4f061fbSSatish Balay   }\
794c4f061fbSSatish Balay }
795c4f061fbSSatish Balay 
796bd481603SBarry Smith /*MC
797bd481603SBarry Smith    MatPreallocateSymmetricSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
798bd481603SBarry Smith        inserted using a local number of the rows and columns
799bd481603SBarry Smith 
800bd481603SBarry Smith    Synopsis:
801c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
802bd481603SBarry Smith 
803bd481603SBarry Smith    Not Collective
804bd481603SBarry Smith 
805bd481603SBarry Smith    Input Parameters:
806bd481603SBarry Smith +  map - the mapping between local numbering and global numbering
807bd481603SBarry Smith .  nrows - the number of rows indicated
8081d73ed98SBarry Smith .  rows - the indices of the rows
809bd481603SBarry Smith .  ncols - the number of columns in the matrix
810bd481603SBarry Smith .  cols - the columns indicated
811bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
812bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
813bd481603SBarry Smith 
814bd481603SBarry Smith 
815bd481603SBarry Smith    Level: intermediate
816bd481603SBarry Smith 
817bd481603SBarry Smith    Notes:
8180598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
819bd481603SBarry Smith 
820bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
821bd481603SBarry Smith 
822bd481603SBarry Smith   Concepts: preallocation^Matrix
823bd481603SBarry Smith 
824bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
825bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
826bd481603SBarry Smith M*/
827d3d32019SBarry Smith #define MatPreallocateSymmetricSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\
828d3d32019SBarry Smith {\
829c1ac3661SBarry Smith   PetscInt __l;\
830d3d32019SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\
831d3d32019SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\
832d3d32019SBarry Smith   for (__l=0;__l<nrows;__l++) {\
833d3d32019SBarry Smith     _4_ierr = MatPreallocateSymmetricSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
834d3d32019SBarry Smith   }\
835d3d32019SBarry Smith }
836d3d32019SBarry Smith 
837bd481603SBarry Smith /*MC
838bd481603SBarry Smith    MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
839bd481603SBarry Smith        inserted using a local number of the rows and columns
840bd481603SBarry Smith 
841bd481603SBarry Smith    Synopsis:
842c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
843bd481603SBarry Smith 
844bd481603SBarry Smith    Not Collective
845bd481603SBarry Smith 
846bd481603SBarry Smith    Input Parameters:
84764075487SBarry Smith +  row - the row
848bd481603SBarry Smith .  ncols - the number of columns in the matrix
849a50f8bf6SHong Zhang -  cols - the columns indicated
850a50f8bf6SHong Zhang 
851a50f8bf6SHong Zhang    Output Parameters:
852a50f8bf6SHong Zhang +  dnz - the array that will be passed to the matrix preallocation routines
853bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
854bd481603SBarry Smith 
855bd481603SBarry Smith 
856bd481603SBarry Smith    Level: intermediate
857bd481603SBarry Smith 
858bd481603SBarry Smith    Notes:
8590598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
860bd481603SBarry Smith 
861bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
862bd481603SBarry Smith 
8631620fd73SBarry Smith    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
8641620fd73SBarry Smith 
865bd481603SBarry Smith   Concepts: preallocation^Matrix
866bd481603SBarry Smith 
867bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
868bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
869bd481603SBarry Smith M*/
870c4f061fbSSatish Balay #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\
871c1ac3661SBarry Smith { PetscInt __i; \
872e32f2f54SBarry 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);\
873e32f2f54SBarry 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);\
8747c922b88SBarry Smith   for (__i=0; __i<nc; __i++) {\
87564075487SBarry Smith     if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \
8767cc688d7SBarry Smith     else dnz[row - __rstart]++;\
8777c922b88SBarry Smith   }\
8787c922b88SBarry Smith }
8797c922b88SBarry Smith 
880bd481603SBarry Smith /*MC
881bd481603SBarry Smith    MatPreallocateSymmetricSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
882bd481603SBarry Smith        inserted using a local number of the rows and columns
883bd481603SBarry Smith 
884bd481603SBarry Smith    Synopsis:
885c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
886bd481603SBarry Smith 
887bd481603SBarry Smith    Not Collective
888bd481603SBarry Smith 
889bd481603SBarry Smith    Input Parameters:
890bd481603SBarry Smith +  nrows - the number of rows indicated
8911d73ed98SBarry Smith .  rows - the indices of the rows
892bd481603SBarry Smith .  ncols - the number of columns in the matrix
893bd481603SBarry Smith .  cols - the columns indicated
894bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
895bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
896bd481603SBarry Smith 
897bd481603SBarry Smith 
898bd481603SBarry Smith    Level: intermediate
899bd481603SBarry Smith 
900bd481603SBarry Smith    Notes:
9010598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
902bd481603SBarry Smith 
903bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
904bd481603SBarry Smith 
9051620fd73SBarry Smith    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
9061620fd73SBarry Smith 
907bd481603SBarry Smith   Concepts: preallocation^Matrix
908bd481603SBarry Smith 
909bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
910bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
911bd481603SBarry Smith M*/
912d3d32019SBarry Smith #define MatPreallocateSymmetricSet(row,nc,cols,dnz,onz) 0;\
913c1ac3661SBarry Smith { PetscInt __i; \
914d3d32019SBarry Smith   for (__i=0; __i<nc; __i++) {\
915d3d32019SBarry Smith     if (cols[__i] >= __end) onz[row - __rstart]++; \
916d3d32019SBarry Smith     else if (cols[__i] >= row) dnz[row - __rstart]++;\
917d3d32019SBarry Smith   }\
918d3d32019SBarry Smith }
919d3d32019SBarry Smith 
920bd481603SBarry Smith /*MC
92116371a99SBarry Smith    MatPreallocateLocation -  An alternative to MatPreallocationSet() that puts the nonzero locations into the matrix if it exists
92216371a99SBarry Smith 
92316371a99SBarry Smith    Synopsis:
92416371a99SBarry Smith    PetscErrorCode MatPreallocateLocations(Mat A,PetscInt row,PetscInt ncols,PetscInt *cols,PetscInt *dnz,PetscInt *onz)
92516371a99SBarry Smith 
92616371a99SBarry Smith    Not Collective
92716371a99SBarry Smith 
92816371a99SBarry Smith    Input Parameters:
92916371a99SBarry Smith .  A - matrix
93016371a99SBarry Smith .  row - row where values exist (must be local to this process)
93116371a99SBarry Smith .  ncols - number of columns
93216371a99SBarry Smith .  cols - columns with nonzeros
93316371a99SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
93416371a99SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
93516371a99SBarry Smith 
93616371a99SBarry Smith 
93716371a99SBarry Smith    Level: intermediate
93816371a99SBarry Smith 
93916371a99SBarry Smith    Notes:
9400598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
94116371a99SBarry Smith 
94216371a99SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
94316371a99SBarry Smith 
94416371a99SBarry Smith    This is a MACRO not a function because it uses a bunch of variables private to the MatPreallocation.... routines.
94516371a99SBarry Smith 
94616371a99SBarry Smith   Concepts: preallocation^Matrix
94716371a99SBarry Smith 
94816371a99SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
94916371a99SBarry Smith           MatPreallocateSymmetricInitialize(), MatPreallocateSymmetricSetLocal()
95016371a99SBarry Smith M*/
95116371a99SBarry 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);}
95216371a99SBarry Smith 
95316371a99SBarry Smith 
95416371a99SBarry Smith /*MC
9550d2a7ff2SSatish Balay    MatPreallocateFinalize - Ends the block of code that will count the number of nonzeros per
956bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
957bd481603SBarry Smith 
958bd481603SBarry Smith    Synopsis:
959c1ac3661SBarry Smith    PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz)
960bd481603SBarry Smith 
961bd481603SBarry Smith    Collective on MPI_Comm
962bd481603SBarry Smith 
963bd481603SBarry Smith    Input Parameters:
96416371a99SBarry Smith +  dnz - the array that was be passed to the matrix preallocation routines
965bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
966bd481603SBarry Smith 
967bd481603SBarry Smith 
968bd481603SBarry Smith    Level: intermediate
969bd481603SBarry Smith 
970bd481603SBarry Smith    Notes:
9710598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
972bd481603SBarry Smith 
973bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
974bd481603SBarry Smith 
9751620fd73SBarry Smith    This is a MACRO not a function because it closes the { started in MatPreallocateInitialize().
9761620fd73SBarry Smith 
977bd481603SBarry Smith   Concepts: preallocation^Matrix
978bd481603SBarry Smith 
979bd481603SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
980bd481603SBarry Smith           MatPreallocateSymmetricInitialize(), MatPreallocateSymmetricSetLocal()
981bd481603SBarry Smith M*/
982a89bc211SBarry Smith #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree2(dnz,onz);CHKERRQ(_4_ierr);}
9837c922b88SBarry Smith 
984bd481603SBarry Smith 
985bd481603SBarry Smith 
9867b80b807SBarry Smith /* Routines unique to particular data structures */
9877087cfbeSBarry Smith extern PetscErrorCode  MatShellGetContext(Mat,void **);
988ba966639SSatish Balay PetscPolymorphicFunction(MatShellGetContext,(Mat A),(A,&t),void*,t)
989698d4c6aSKris Buschelman 
9907087cfbeSBarry Smith extern PetscErrorCode  MatInodeAdjustForInodes(Mat,IS*,IS*);
9917087cfbeSBarry Smith extern PetscErrorCode  MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *);
992698d4c6aSKris Buschelman 
9937087cfbeSBarry Smith extern PetscErrorCode  MatSeqAIJSetColumnIndices(Mat,PetscInt[]);
9947087cfbeSBarry Smith extern PetscErrorCode  MatSeqBAIJSetColumnIndices(Mat,PetscInt[]);
9957087cfbeSBarry Smith extern PetscErrorCode  MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
9967087cfbeSBarry Smith extern PetscErrorCode  MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
9977087cfbeSBarry Smith extern PetscErrorCode  MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
9987b80b807SBarry Smith 
999a96a251dSBarry Smith #define MAT_SKIP_ALLOCATION -4
1000a96a251dSBarry Smith 
10017087cfbeSBarry Smith extern PetscErrorCode  MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
1002ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz))
10037087cfbeSBarry Smith extern PetscErrorCode  MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
1004ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqSBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz))
10057087cfbeSBarry Smith extern PetscErrorCode  MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]);
1006ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqAIJSetPreallocation,(Mat A,const PetscInt nnz[]),(A,0,nnz))
1007273d9f13SBarry Smith 
10087087cfbeSBarry Smith extern PetscErrorCode  MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
1009ba966639SSatish Balay PetscPolymorphicSubroutine(MatMPIBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[],const PetscInt onz[]),(A,bs,0,nnz,0,onz))
10107087cfbeSBarry Smith extern PetscErrorCode  MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
10117087cfbeSBarry Smith extern PetscErrorCode  MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
10127087cfbeSBarry Smith extern PetscErrorCode  MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []);
10137087cfbeSBarry Smith extern PetscErrorCode  MatSeqBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
10147087cfbeSBarry Smith extern PetscErrorCode  MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]);
10157087cfbeSBarry Smith extern PetscErrorCode  MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
10167087cfbeSBarry Smith extern PetscErrorCode  MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]);
10177087cfbeSBarry Smith extern PetscErrorCode  MatMPIDenseSetPreallocation(Mat,PetscScalar[]);
10187087cfbeSBarry Smith extern PetscErrorCode  MatSeqDenseSetPreallocation(Mat,PetscScalar[]);
10197087cfbeSBarry Smith extern PetscErrorCode  MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,PetscInt*[]);
10207087cfbeSBarry Smith extern PetscErrorCode  MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,PetscInt*[]);
10217087cfbeSBarry Smith extern PetscErrorCode  MatAdicSetLocalFunction(Mat,void (*)(void));
1022273d9f13SBarry Smith 
10237087cfbeSBarry Smith extern PetscErrorCode  MatSeqDenseSetLDA(Mat,PetscInt);
10247087cfbeSBarry Smith extern PetscErrorCode  MatDenseGetLocalMatrix(Mat,Mat*);
10251b807ce4Svictorle 
10267087cfbeSBarry Smith extern PetscErrorCode  MatStoreValues(Mat);
10277087cfbeSBarry Smith extern PetscErrorCode  MatRetrieveValues(Mat);
10282e8a6d31SBarry Smith 
10297087cfbeSBarry Smith extern PetscErrorCode  MatDAADSetCtx(Mat,void*);
10303a7fca6bSBarry Smith 
1031b3a44c85SBarry Smith extern PetscErrorCode  MatFindNonzeroRows(Mat,IS*);
10327b80b807SBarry Smith /*
10337b80b807SBarry Smith   These routines are not usually accessed directly, rather solving is
103494b7f48cSBarry Smith   done through the KSP and PC interfaces.
10357b80b807SBarry Smith */
10367b80b807SBarry Smith 
1037d9274352SBarry Smith /*E
1038d9274352SBarry Smith     MatOrderingType - String with the name of a PETSc matrix ordering or the creation function
1039d9274352SBarry Smith        with an optional dynamic library name, for example
1040d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:orderingcreate()
1041d9274352SBarry Smith 
1042d9274352SBarry Smith    Level: beginner
1043d9274352SBarry Smith 
1044e5a9bf91SBarry Smith    Cannot use const because the PC objects manipulate the string
1045e5a9bf91SBarry Smith 
1046d9274352SBarry Smith .seealso: MatGetOrdering()
1047d9274352SBarry Smith E*/
10483eaad2caSSatish Balay #define MatOrderingType char*
10492692d6eeSBarry Smith #define MATORDERINGNATURAL     "natural"
10502692d6eeSBarry Smith #define MATORDERINGND          "nd"
10512692d6eeSBarry Smith #define MATORDERING1WD         "1wd"
10522692d6eeSBarry Smith #define MATORDERINGRCM         "rcm"
10532692d6eeSBarry Smith #define MATORDERINGQMD         "qmd"
10542692d6eeSBarry Smith #define MATORDERINGROWLENGTH   "rowlength"
10552692d6eeSBarry Smith #define MATORDERINGDSC_ND      "dsc_nd"         /* these three are only for DSCPACK, see its documentation for details */
10562692d6eeSBarry Smith #define MATORDERINGDSC_MMD     "dsc_mmd"
10572692d6eeSBarry Smith #define MATORDERINGDSC_MDF     "dsc_mdf"
1058312e372aSBarry Smith #define MATORDERINGAMD         "amd"            /* only works if UMFPACK is installed with PETSc */
1059b12f92e5SBarry Smith 
10607087cfbeSBarry Smith extern PetscErrorCode  MatGetOrdering(Mat,const MatOrderingType,IS*,IS*);
10617087cfbeSBarry Smith extern PetscErrorCode  MatGetOrderingList(PetscFList *list);
10627087cfbeSBarry Smith extern PetscErrorCode  MatOrderingRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,const MatOrderingType,IS*,IS*));
106330de9b25SBarry Smith 
106430de9b25SBarry Smith /*MC
10651890ba74SBarry Smith    MatOrderingRegisterDynamic - Adds a new sparse matrix ordering to the matrix package.
106630de9b25SBarry Smith 
106730de9b25SBarry Smith    Synopsis:
10681890ba74SBarry Smith    PetscErrorCode MatOrderingRegisterDynamic(const char *name_ordering,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatOrdering))
106930de9b25SBarry Smith 
107030de9b25SBarry Smith    Not Collective
107130de9b25SBarry Smith 
107230de9b25SBarry Smith    Input Parameters:
10732692d6eeSBarry Smith +  sname - name of ordering (for example MATORDERINGND)
107430de9b25SBarry Smith .  path - location of library where creation routine is
107530de9b25SBarry Smith .  name - name of function that creates the ordering type,a string
107630de9b25SBarry Smith -  function - function pointer that creates the ordering
107730de9b25SBarry Smith 
107830de9b25SBarry Smith    Level: developer
107930de9b25SBarry Smith 
108030de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
108130de9b25SBarry Smith    is ignored.
108230de9b25SBarry Smith 
108330de9b25SBarry Smith    Sample usage:
108430de9b25SBarry Smith .vb
108530de9b25SBarry Smith    MatOrderingRegisterDynamic("my_order",/home/username/my_lib/lib/libO/solaris/mylib.a,
108630de9b25SBarry Smith                "MyOrder",MyOrder);
108730de9b25SBarry Smith .ve
108830de9b25SBarry Smith 
108930de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
109030de9b25SBarry Smith $     MatOrderingSetType(part,"my_order)
109130de9b25SBarry Smith    or at runtime via the option
10922401956bSBarry Smith $     -pc_factor_mat_ordering_type my_order
109330de9b25SBarry Smith 
1094ab901514SBarry Smith    ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values.
109530de9b25SBarry Smith 
109630de9b25SBarry Smith .keywords: matrix, ordering, register
109730de9b25SBarry Smith 
109830de9b25SBarry Smith .seealso: MatOrderingRegisterDestroy(), MatOrderingRegisterAll()
109930de9b25SBarry Smith M*/
1100aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1101f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,0)
1102b12f92e5SBarry Smith #else
1103f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,d)
1104b12f92e5SBarry Smith #endif
110530de9b25SBarry Smith 
11067087cfbeSBarry Smith extern PetscErrorCode  MatOrderingRegisterDestroy(void);
11077087cfbeSBarry Smith extern PetscErrorCode  MatOrderingRegisterAll(const char[]);
1108ace3abfcSBarry Smith extern PetscBool  MatOrderingRegisterAllCalled;
1109b0a32e0cSBarry Smith extern PetscFList MatOrderingList;
1110d4fbbf0eSBarry Smith 
11117087cfbeSBarry Smith extern PetscErrorCode  MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS);
1112a2ce50c7SBarry Smith 
1113d91e6319SBarry Smith /*S
1114d90ac83dSHong Zhang     MatFactorShiftType - Numeric Shift.
1115d90ac83dSHong Zhang 
1116d90ac83dSHong Zhang    Level: beginner
1117d90ac83dSHong Zhang 
1118d90ac83dSHong Zhang S*/
1119d90ac83dSHong Zhang typedef enum {MAT_SHIFT_NONE,MAT_SHIFT_NONZERO,MAT_SHIFT_POSITIVE_DEFINITE,MAT_SHIFT_INBLOCKS} MatFactorShiftType;
1120d90ac83dSHong Zhang extern const char *MatFactorShiftTypes[];
1121d90ac83dSHong Zhang 
1122d90ac83dSHong Zhang /*S
11232401956bSBarry Smith    MatFactorInfo - Data passed into the matrix factorization routines
1124b00f7748SHong Zhang 
112561cad860SBarry Smith    In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE, that is use
112661cad860SBarry Smith $     MatFactorInfo  info(MAT_FACTORINFO_SIZE)
1127b00f7748SHong Zhang 
112815e8a5b3SHong Zhang    Notes: These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC.
1129b00f7748SHong Zhang 
113061cad860SBarry Smith       You can use MatFactorInfoInitialize() to set default values.
113161cad860SBarry Smith 
1132b00f7748SHong Zhang    Level: developer
1133b00f7748SHong Zhang 
1134d7d82daaSBarry Smith .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(),
1135d7d82daaSBarry Smith           MatFactorInfoInitialize()
1136b00f7748SHong Zhang 
1137b00f7748SHong Zhang S*/
1138b00f7748SHong Zhang typedef struct {
113915e8a5b3SHong Zhang   PetscReal     diagonal_fill;  /* force diagonal to fill in if initially not filled */
114085317021SBarry Smith   PetscReal     usedt;
114115e8a5b3SHong Zhang   PetscReal     dt;             /* drop tolerance */
1142b00f7748SHong Zhang   PetscReal     dtcol;          /* tolerance for pivoting */
114315e8a5b3SHong Zhang   PetscReal     dtcount;        /* maximum nonzeros to be allowed per row */
114467e28bfeSBarry Smith   PetscReal     fill;           /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */
1145348344bbSBarry Smith   PetscReal     levels;         /* ICC/ILU(levels) */
1146bcd9e38bSBarry Smith   PetscReal     pivotinblocks;  /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0
1147bcd9e38bSBarry Smith                                    factorization may be faster if do not pivot */
114815e8a5b3SHong Zhang   PetscReal     zeropivot;      /* pivot is called zero if less than this */
1149f4db908eSBarry Smith   PetscReal     shifttype;      /* type of shift added to matrix factor to prevent zero pivots */
1150f4db908eSBarry Smith   PetscReal     shiftamount;     /* how large the shift is */
115115e8a5b3SHong Zhang } MatFactorInfo;
1152ffa6d0a5SLois Curfman McInnes 
11537087cfbeSBarry Smith extern PetscErrorCode  MatFactorInfoInitialize(MatFactorInfo*);
11547087cfbeSBarry Smith extern PetscErrorCode  MatCholeskyFactor(Mat,IS,const MatFactorInfo*);
11557087cfbeSBarry Smith extern PetscErrorCode  MatCholeskyFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
11567087cfbeSBarry Smith extern PetscErrorCode  MatCholeskyFactorNumeric(Mat,Mat,const MatFactorInfo*);
11577087cfbeSBarry Smith extern PetscErrorCode  MatLUFactor(Mat,IS,IS,const MatFactorInfo*);
11587087cfbeSBarry Smith extern PetscErrorCode  MatILUFactor(Mat,IS,IS,const MatFactorInfo*);
11597087cfbeSBarry Smith extern PetscErrorCode  MatLUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
11607087cfbeSBarry Smith extern PetscErrorCode  MatILUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
11617087cfbeSBarry Smith extern PetscErrorCode  MatICCFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
11627087cfbeSBarry Smith extern PetscErrorCode  MatICCFactor(Mat,IS,const MatFactorInfo*);
11637087cfbeSBarry Smith extern PetscErrorCode  MatLUFactorNumeric(Mat,Mat,const MatFactorInfo*);
11647087cfbeSBarry Smith extern PetscErrorCode  MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*);
11657087cfbeSBarry Smith extern PetscErrorCode  MatSolve(Mat,Vec,Vec);
11667087cfbeSBarry Smith extern PetscErrorCode  MatForwardSolve(Mat,Vec,Vec);
11677087cfbeSBarry Smith extern PetscErrorCode  MatBackwardSolve(Mat,Vec,Vec);
11687087cfbeSBarry Smith extern PetscErrorCode  MatSolveAdd(Mat,Vec,Vec,Vec);
11697087cfbeSBarry Smith extern PetscErrorCode  MatSolveTranspose(Mat,Vec,Vec);
11707087cfbeSBarry Smith extern PetscErrorCode  MatSolveTransposeAdd(Mat,Vec,Vec,Vec);
11717087cfbeSBarry Smith extern PetscErrorCode  MatSolves(Mat,Vecs,Vecs);
11728ed539a5SBarry Smith 
11737087cfbeSBarry Smith extern PetscErrorCode  MatSetUnfactored(Mat);
1174bb5a7306SBarry Smith 
1175d91e6319SBarry Smith /*E
1176d91e6319SBarry Smith     MatSORType - What type of (S)SOR to perform
1177bb1eb677SSatish Balay 
1178d91e6319SBarry Smith     Level: beginner
1179d91e6319SBarry Smith 
1180d9274352SBarry Smith    May be bitwise ORd together
1181d9274352SBarry Smith 
1182d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
1183d91e6319SBarry Smith 
11844e7234bfSBarry Smith    MatSORType may be bitwise ORd together, so do not change the numbers
11854e7234bfSBarry Smith 
118641f059aeSBarry Smith .seealso: MatSOR()
1187d91e6319SBarry Smith E*/
1188ee50ffe9SBarry Smith typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3,
1189ee50ffe9SBarry Smith               SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8,
1190ee50ffe9SBarry Smith               SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16,
119184cb2905SBarry Smith               SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType;
11927087cfbeSBarry Smith extern PetscErrorCode  MatSOR(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
11938ed539a5SBarry Smith 
1194d4fbbf0eSBarry Smith /*
1195639f9d9dSBarry Smith     These routines are for efficiently computing Jacobians via finite differences.
1196639f9d9dSBarry Smith */
1197b12f92e5SBarry Smith 
1198d9274352SBarry Smith /*E
1199d9274352SBarry Smith     MatColoringType - String with the name of a PETSc matrix coloring or the creation function
1200d9274352SBarry Smith        with an optional dynamic library name, for example
1201d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:coloringcreate()
1202d9274352SBarry Smith 
1203d9274352SBarry Smith    Level: beginner
1204d9274352SBarry Smith 
1205d9274352SBarry Smith .seealso: MatGetColoring()
1206d9274352SBarry Smith E*/
1207a313700dSBarry Smith #define MatColoringType char*
12082692d6eeSBarry Smith #define MATCOLORINGNATURAL "natural"
12092692d6eeSBarry Smith #define MATCOLORINGSL      "sl"
12102692d6eeSBarry Smith #define MATCOLORINGLF      "lf"
12112692d6eeSBarry Smith #define MATCOLORINGID      "id"
1212b12f92e5SBarry Smith 
12137087cfbeSBarry Smith extern PetscErrorCode  MatGetColoring(Mat,const MatColoringType,ISColoring*);
12147087cfbeSBarry Smith extern PetscErrorCode  MatColoringRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,MatColoringType,ISColoring *));
121530de9b25SBarry Smith 
121630de9b25SBarry Smith /*MC
121730de9b25SBarry Smith    MatColoringRegisterDynamic - Adds a new sparse matrix coloring to the
121830de9b25SBarry Smith                                matrix package.
121930de9b25SBarry Smith 
122030de9b25SBarry Smith    Synopsis:
12211890ba74SBarry Smith    PetscErrorCode MatColoringRegisterDynamic(const char *name_coloring,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatColoring))
122230de9b25SBarry Smith 
122330de9b25SBarry Smith    Not Collective
122430de9b25SBarry Smith 
122530de9b25SBarry Smith    Input Parameters:
12262692d6eeSBarry Smith +  sname - name of Coloring (for example MATCOLORINGSL)
122730de9b25SBarry Smith .  path - location of library where creation routine is
122830de9b25SBarry Smith .  name - name of function that creates the Coloring type, a string
122930de9b25SBarry Smith -  function - function pointer that creates the coloring
123030de9b25SBarry Smith 
123130de9b25SBarry Smith    Level: developer
123230de9b25SBarry Smith 
123330de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
123430de9b25SBarry Smith    is ignored.
123530de9b25SBarry Smith 
123630de9b25SBarry Smith    Sample usage:
123730de9b25SBarry Smith .vb
123830de9b25SBarry Smith    MatColoringRegisterDynamic("my_color",/home/username/my_lib/lib/libO/solaris/mylib.a,
123930de9b25SBarry Smith                "MyColor",MyColor);
124030de9b25SBarry Smith .ve
124130de9b25SBarry Smith 
124230de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
124330de9b25SBarry Smith $     MatColoringSetType(part,"my_color")
124430de9b25SBarry Smith    or at runtime via the option
124530de9b25SBarry Smith $     -mat_coloring_type my_color
124630de9b25SBarry Smith 
1247ab901514SBarry Smith    $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
124830de9b25SBarry Smith 
124930de9b25SBarry Smith .keywords: matrix, Coloring, register
125030de9b25SBarry Smith 
125130de9b25SBarry Smith .seealso: MatColoringRegisterDestroy(), MatColoringRegisterAll()
125230de9b25SBarry Smith M*/
1253aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1254f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,0)
1255b12f92e5SBarry Smith #else
1256f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,d)
1257b12f92e5SBarry Smith #endif
125830de9b25SBarry Smith 
1259ace3abfcSBarry Smith extern PetscBool  MatColoringRegisterAllCalled;
1260f1144a30SSatish Balay 
12617087cfbeSBarry Smith extern PetscErrorCode  MatColoringRegisterAll(const char[]);
12627087cfbeSBarry Smith extern PetscErrorCode  MatColoringRegisterDestroy(void);
12637087cfbeSBarry Smith extern PetscErrorCode  MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*);
1264639f9d9dSBarry Smith 
1265d9274352SBarry Smith /*S
1266d9274352SBarry Smith      MatFDColoring - Object for computing a sparse Jacobian via finite differences
1267d9274352SBarry Smith         and coloring
1268639f9d9dSBarry Smith 
1269d9274352SBarry Smith    Level: beginner
1270d9274352SBarry Smith 
1271d9274352SBarry Smith   Concepts: coloring, sparse Jacobian, finite differences
1272d9274352SBarry Smith 
1273d9274352SBarry Smith .seealso:  MatFDColoringCreate()
1274d9274352SBarry Smith S*/
1275e2a1c21fSSatish Balay typedef struct _p_MatFDColoring* MatFDColoring;
1276639f9d9dSBarry Smith 
12777087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringCreate(Mat,ISColoring,MatFDColoring *);
12787087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringDestroy(MatFDColoring);
12797087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringView(MatFDColoring,PetscViewer);
12807087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*);
12817087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**);
12827087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal);
12837087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringSetFromOptions(MatFDColoring);
12847087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringApply(Mat,MatFDColoring,Vec,MatStructure*,void *);
12857087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringApplyTS(Mat,MatFDColoring,PetscReal,Vec,MatStructure*,void *);
12867087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringSetF(MatFDColoring,Vec);
12877087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,PetscInt*[]);
1288639f9d9dSBarry Smith /*
12890752156aSBarry Smith     These routines are for partitioning matrices: currently used only
12903eda8832SBarry Smith   for adjacency matrix, MatCreateMPIAdj().
12910752156aSBarry Smith */
1292ca161407SBarry Smith 
1293d9274352SBarry Smith /*S
1294d9274352SBarry Smith      MatPartitioning - Object for managing the partitioning of a matrix or graph
1295d9274352SBarry Smith 
1296d9274352SBarry Smith    Level: beginner
1297d9274352SBarry Smith 
1298d9274352SBarry Smith   Concepts: partitioning
1299d9274352SBarry Smith 
1300743a1026Svictorle .seealso:  MatPartitioningCreate(), MatPartitioningType
1301d9274352SBarry Smith S*/
130291e9ee9fSBarry Smith typedef struct _p_MatPartitioning* MatPartitioning;
1303d9274352SBarry Smith 
1304d9274352SBarry Smith /*E
13055bc6e7ccSvictorle     MatPartitioningType - String with the name of a PETSc matrix partitioning or the creation function
1306d9274352SBarry Smith        with an optional dynamic library name, for example
1307d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:partitioningcreate()
1308d9274352SBarry Smith 
1309d9274352SBarry Smith    Level: beginner
1310d9274352SBarry Smith 
1311b547a13bSvictorle .seealso: MatPartitioningCreate(), MatPartitioning
1312d9274352SBarry Smith E*/
1313a313700dSBarry Smith #define MatPartitioningType char*
13142692d6eeSBarry Smith #define MATPARTITIONINGCURRENT  "current"
13152692d6eeSBarry Smith #define MATPARTITIONINGSQUARE   "square"
13162692d6eeSBarry Smith #define MATPARTITIONINGPARMETIS "parmetis"
13172692d6eeSBarry Smith #define MATPARTITIONINGCHACO    "chaco"
13182692d6eeSBarry Smith #define MATPARTITIONINGJOSTLE   "jostle"
13192692d6eeSBarry Smith #define MATPARTITIONINGPARTY    "party"
13202692d6eeSBarry Smith #define MATPARTITIONINGSCOTCH   "scotch"
132171306c60Sjroman 
1322ca161407SBarry Smith 
13237087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningCreate(MPI_Comm,MatPartitioning*);
13247087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningSetType(MatPartitioning,const MatPartitioningType);
13257087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningSetNParts(MatPartitioning,PetscInt);
13267087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningSetAdjacency(MatPartitioning,Mat);
13277087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]);
13287087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []);
13297087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningApply(MatPartitioning,IS*);
13307087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningDestroy(MatPartitioning);
13312aabb6bbSBarry Smith 
13327087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatPartitioning));
133330de9b25SBarry Smith 
133430de9b25SBarry Smith /*MC
133530de9b25SBarry Smith    MatPartitioningRegisterDynamic - Adds a new sparse matrix partitioning to the
133630de9b25SBarry Smith    matrix package.
133730de9b25SBarry Smith 
133830de9b25SBarry Smith    Synopsis:
13391890ba74SBarry Smith    PetscErrorCode MatPartitioningRegisterDynamic(const char *name_partitioning,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatPartitioning))
134030de9b25SBarry Smith 
134130de9b25SBarry Smith    Not Collective
134230de9b25SBarry Smith 
134330de9b25SBarry Smith    Input Parameters:
13442692d6eeSBarry Smith +  sname - name of partitioning (for example MATPARTITIONINGCURRENT) or parmetis
134530de9b25SBarry Smith .  path - location of library where creation routine is
134630de9b25SBarry Smith .  name - name of function that creates the partitioning type, a string
134730de9b25SBarry Smith -  function - function pointer that creates the partitioning type
134830de9b25SBarry Smith 
134930de9b25SBarry Smith    Level: developer
135030de9b25SBarry Smith 
135130de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
135230de9b25SBarry Smith    is ignored.
135330de9b25SBarry Smith 
135430de9b25SBarry Smith    Sample usage:
135530de9b25SBarry Smith .vb
135630de9b25SBarry Smith    MatPartitioningRegisterDynamic("my_part",/home/username/my_lib/lib/libO/solaris/mylib.a,
135730de9b25SBarry Smith                "MyPartCreate",MyPartCreate);
135830de9b25SBarry Smith .ve
135930de9b25SBarry Smith 
136030de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
136130de9b25SBarry Smith $     MatPartitioningSetType(part,"my_part")
136230de9b25SBarry Smith    or at runtime via the option
136330de9b25SBarry Smith $     -mat_partitioning_type my_part
136430de9b25SBarry Smith 
1365ab901514SBarry Smith    $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
136630de9b25SBarry Smith 
136730de9b25SBarry Smith .keywords: matrix, partitioning, register
136830de9b25SBarry Smith 
136930de9b25SBarry Smith .seealso: MatPartitioningRegisterDestroy(), MatPartitioningRegisterAll()
137030de9b25SBarry Smith M*/
1371aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1372f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,0)
13732aabb6bbSBarry Smith #else
1374f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,d)
13752aabb6bbSBarry Smith #endif
13762aabb6bbSBarry Smith 
1377ace3abfcSBarry Smith extern PetscBool  MatPartitioningRegisterAllCalled;
1378f1144a30SSatish Balay 
13797087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningRegisterAll(const char[]);
13807087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningRegisterDestroy(void);
13812bad1931SBarry Smith 
13827087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningView(MatPartitioning,PetscViewer);
13837087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningSetFromOptions(MatPartitioning);
13847087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningGetType(MatPartitioning,const MatPartitioningType*);
1385ca161407SBarry Smith 
13867087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningParmetisSetCoarseSequential(MatPartitioning);
13877087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningParmetisGetEdgeCut(MatPartitioning, PetscInt *);
13880752156aSBarry Smith 
13897087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningJostleSetCoarseLevel(MatPartitioning,PetscReal);
13907087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningJostleSetCoarseSequential(MatPartitioning);
139171306c60Sjroman 
1392954b3ec6SBarry Smith typedef enum {MP_CHACO_MULTILEVEL_KL,MP_CHACO_SPECTRAL,MP_CHACO_LINEAR,MP_CHACO_RANDOM, MP_CHACO_SCATTERED} MPChacoGlobalType;
13937087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningChacoSetGlobal(MatPartitioning, MPChacoGlobalType);
139471306c60Sjroman typedef enum { MP_CHACO_KERNIGHAN_LIN, MP_CHACO_NONE } MPChacoLocalType;
13957087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningChacoSetLocal(MatPartitioning, MPChacoLocalType);
13967087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal);
139771306c60Sjroman typedef enum { MP_CHACO_LANCZOS, MP_CHACO_RQI_SYMMLQ } MPChacoEigenType;
13987087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType);
13997087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningChacoSetEigenTol(MatPartitioning, PetscReal);
14007087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningChacoSetEigenNumber(MatPartitioning, PetscInt);
140171306c60Sjroman 
140271306c60Sjroman #define MP_PARTY_OPT "opt"
140371306c60Sjroman #define MP_PARTY_LIN "lin"
140471306c60Sjroman #define MP_PARTY_SCA "sca"
140571306c60Sjroman #define MP_PARTY_RAN "ran"
140671306c60Sjroman #define MP_PARTY_GBF "gbf"
140771306c60Sjroman #define MP_PARTY_GCF "gcf"
140871306c60Sjroman #define MP_PARTY_BUB "bub"
140971306c60Sjroman #define MP_PARTY_DEF "def"
14107087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningPartySetGlobal(MatPartitioning, const char*);
141171306c60Sjroman #define MP_PARTY_HELPFUL_SETS "hs"
141271306c60Sjroman #define MP_PARTY_KERNIGHAN_LIN "kl"
141371306c60Sjroman #define MP_PARTY_NONE "no"
14147087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningPartySetLocal(MatPartitioning, const char*);
14157087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal);
14167087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningPartySetBipart(MatPartitioning,PetscBool );
14177087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscBool );
141871306c60Sjroman 
141971306c60Sjroman typedef enum { MP_SCOTCH_GREEDY, MP_SCOTCH_GPS, MP_SCOTCH_GR_GPS } MPScotchGlobalType;
14207087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningScotchSetArch(MatPartitioning,const char*);
14217087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningScotchSetMultilevel(MatPartitioning);
14227087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningScotchSetGlobal(MatPartitioning,MPScotchGlobalType);
14237087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningScotchSetCoarseLevel(MatPartitioning,PetscReal);
14247087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningScotchSetHostList(MatPartitioning,const char*);
142571306c60Sjroman typedef enum { MP_SCOTCH_KERNIGHAN_LIN, MP_SCOTCH_NONE } MPScotchLocalType;
14267087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningScotchSetLocal(MatPartitioning,MPScotchLocalType);
14277087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningScotchSetMapping(MatPartitioning);
14287087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningScotchSetStrategy(MatPartitioning,char*);
142971306c60Sjroman 
143009573ac7SBarry Smith extern PetscErrorCode MatMeshToVertexGraph(Mat,PetscInt,Mat*);
143109573ac7SBarry Smith extern PetscErrorCode MatMeshToCellGraph(Mat,PetscInt,Mat*);
1432591294e4SBarry Smith 
14330752156aSBarry Smith /*
14340a835dfdSSatish Balay     If you add entries here you must also add them to finclude/petscmat.h
1435d4fbbf0eSBarry Smith */
14361c1c02c0SLois Curfman McInnes typedef enum { MATOP_SET_VALUES=0,
14371c1c02c0SLois Curfman McInnes                MATOP_GET_ROW=1,
14381c1c02c0SLois Curfman McInnes                MATOP_RESTORE_ROW=2,
14391c1c02c0SLois Curfman McInnes                MATOP_MULT=3,
14401c1c02c0SLois Curfman McInnes                MATOP_MULT_ADD=4,
14417c922b88SBarry Smith                MATOP_MULT_TRANSPOSE=5,
14427c922b88SBarry Smith                MATOP_MULT_TRANSPOSE_ADD=6,
14431c1c02c0SLois Curfman McInnes                MATOP_SOLVE=7,
14441c1c02c0SLois Curfman McInnes                MATOP_SOLVE_ADD=8,
14457c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE=9,
14467c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE_ADD=10,
14471c1c02c0SLois Curfman McInnes                MATOP_LUFACTOR=11,
14481c1c02c0SLois Curfman McInnes                MATOP_CHOLESKYFACTOR=12,
1449a714c06dSBarry Smith                MATOP_SOR=13,
14501c1c02c0SLois Curfman McInnes                MATOP_TRANSPOSE=14,
14511c1c02c0SLois Curfman McInnes                MATOP_GETINFO=15,
14521c1c02c0SLois Curfman McInnes                MATOP_EQUAL=16,
14531c1c02c0SLois Curfman McInnes                MATOP_GET_DIAGONAL=17,
14541c1c02c0SLois Curfman McInnes                MATOP_DIAGONAL_SCALE=18,
14551c1c02c0SLois Curfman McInnes                MATOP_NORM=19,
14561c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_BEGIN=20,
14571c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_END=21,
1458d519adbfSMatthew Knepley                MATOP_SET_OPTION=22,
1459d519adbfSMatthew Knepley                MATOP_ZERO_ENTRIES=23,
1460d519adbfSMatthew Knepley                MATOP_ZERO_ROWS=24,
1461d519adbfSMatthew Knepley                MATOP_LUFACTOR_SYMBOLIC=25,
1462d519adbfSMatthew Knepley                MATOP_LUFACTOR_NUMERIC=26,
1463d519adbfSMatthew Knepley                MATOP_CHOLESKY_FACTOR_SYMBOLIC=27,
1464d519adbfSMatthew Knepley                MATOP_CHOLESKY_FACTOR_NUMERIC=28,
1465d519adbfSMatthew Knepley                MATOP_SETUP_PREALLOCATION=29,
1466d519adbfSMatthew Knepley                MATOP_ILUFACTOR_SYMBOLIC=30,
1467d519adbfSMatthew Knepley                MATOP_ICCFACTOR_SYMBOLIC=31,
1468d519adbfSMatthew Knepley                MATOP_GET_ARRAY=32,
1469d519adbfSMatthew Knepley                MATOP_RESTORE_ARRAY=33,
1470d519adbfSMatthew Knepley                MATOP_DUPLICATE=34,
1471d519adbfSMatthew Knepley                MATOP_FORWARD_SOLVE=35,
1472d519adbfSMatthew Knepley                MATOP_BACKWARD_SOLVE=36,
1473d519adbfSMatthew Knepley                MATOP_ILUFACTOR=37,
1474d519adbfSMatthew Knepley                MATOP_ICCFACTOR=38,
1475d519adbfSMatthew Knepley                MATOP_AXPY=39,
1476d519adbfSMatthew Knepley                MATOP_GET_SUBMATRICES=40,
1477d519adbfSMatthew Knepley                MATOP_INCREASE_OVERLAP=41,
1478d519adbfSMatthew Knepley                MATOP_GET_VALUES=42,
1479d519adbfSMatthew Knepley                MATOP_COPY=43,
1480d519adbfSMatthew Knepley                MATOP_GET_ROW_MAX=44,
1481d519adbfSMatthew Knepley                MATOP_SCALE=45,
1482d519adbfSMatthew Knepley                MATOP_SHIFT=46,
148335153367SBarry Smith                MATOP_DIAGONAL_SET=47,
1484d519adbfSMatthew Knepley                MATOP_ILUDT_FACTOR=48,
1485d519adbfSMatthew Knepley                MATOP_SET_BLOCK_SIZE=49,
1486d519adbfSMatthew Knepley                MATOP_GET_ROW_IJ=50,
1487d519adbfSMatthew Knepley                MATOP_RESTORE_ROW_IJ=51,
1488d519adbfSMatthew Knepley                MATOP_GET_COLUMN_IJ=52,
1489d519adbfSMatthew Knepley                MATOP_RESTORE_COLUMN_IJ=53,
1490d519adbfSMatthew Knepley                MATOP_FDCOLORING_CREATE=54,
1491d519adbfSMatthew Knepley                MATOP_COLORING_PATCH=55,
1492d519adbfSMatthew Knepley                MATOP_SET_UNFACTORED=56,
1493d519adbfSMatthew Knepley                MATOP_PERMUTE=57,
1494d519adbfSMatthew Knepley                MATOP_SET_VALUES_BLOCKED=58,
1495d519adbfSMatthew Knepley                MATOP_GET_SUBMATRIX=59,
1496d519adbfSMatthew Knepley                MATOP_DESTROY=60,
1497d519adbfSMatthew Knepley                MATOP_VIEW=61,
1498d519adbfSMatthew Knepley                MATOP_CONVERT_FROM=62,
1499d519adbfSMatthew Knepley                MATOP_USE_SCALED_FORM=63,
1500d519adbfSMatthew Knepley                MATOP_SCALE_SYSTEM=64,
1501d519adbfSMatthew Knepley                MATOP_UNSCALE_SYSTEM=65,
1502d519adbfSMatthew Knepley                MATOP_SET_LOCAL_TO_GLOBAL_MAP=66,
1503d519adbfSMatthew Knepley                MATOP_SET_VALUES_LOCAL=67,
1504d519adbfSMatthew Knepley                MATOP_ZERO_ROWS_LOCAL=68,
1505d519adbfSMatthew Knepley                MATOP_GET_ROW_MAX_ABS=69,
1506d519adbfSMatthew Knepley                MATOP_GET_ROW_MIN_ABS=70,
1507d519adbfSMatthew Knepley                MATOP_CONVERT=71,
1508d519adbfSMatthew Knepley                MATOP_SET_COLORING=72,
1509d519adbfSMatthew Knepley                MATOP_SET_VALUES_ADIC=73,
1510d519adbfSMatthew Knepley                MATOP_SET_VALUES_ADIFOR=74,
1511d519adbfSMatthew Knepley                MATOP_FD_COLORING_APPLY=75,
1512d519adbfSMatthew Knepley                MATOP_SET_FROM_OPTIONS=76,
1513d519adbfSMatthew Knepley                MATOP_MULT_CON=77,
1514d519adbfSMatthew Knepley                MATOP_MULT_TRANSPOSE_CON=78,
1515d519adbfSMatthew Knepley                MATOP_PERMUTE_SPARSIFY=79,
1516d519adbfSMatthew Knepley                MATOP_MULT_MULTIPLE=80,
1517d519adbfSMatthew Knepley                MATOP_SOLVE_MULTIPLE=81,
1518d519adbfSMatthew Knepley                MATOP_GET_INERTIA=82,
1519d519adbfSMatthew Knepley                MATOP_LOAD=83,
1520d519adbfSMatthew Knepley                MATOP_IS_SYMMETRIC=84,
1521d519adbfSMatthew Knepley                MATOP_IS_HERMITIAN=85,
1522d519adbfSMatthew Knepley                MATOP_IS_STRUCTURALLY_SYMMETRIC=86,
152341f059aeSBarry Smith                MATOP_DUMMY=87,
1524d519adbfSMatthew Knepley                MATOP_GET_VECS=88,
1525d519adbfSMatthew Knepley                MATOP_MAT_MULT=89,
1526d519adbfSMatthew Knepley                MATOP_MAT_MULT_SYMBOLIC=90,
1527d519adbfSMatthew Knepley                MATOP_MAT_MULT_NUMERIC=91,
1528d519adbfSMatthew Knepley                MATOP_PTAP=92,
1529d519adbfSMatthew Knepley                MATOP_PTAP_SYMBOLIC=93,
1530d519adbfSMatthew Knepley                MATOP_PTAP_NUMERIC=94,
1531d519adbfSMatthew Knepley                MATOP_MAT_MULTTRANSPOSE=95,
15321763ac4eSSatish Balay                MATOP_MAT_MULTTRANSPOSE_SYM=96,
15331763ac4eSSatish Balay                MATOP_MAT_MULTTRANSPOSE_NUM=97,
1534d519adbfSMatthew Knepley                MATOP_PTAP_SYMBOLIC_SEQAIJ=98,
1535d519adbfSMatthew Knepley                MATOP_PTAP_NUMERIC_SEQAIJ=99,
1536d519adbfSMatthew Knepley                MATOP_PTAP_SYMBOLIC_MPIAIJ=100,
1537d519adbfSMatthew Knepley                MATOP_PTAP_NUMERIC_MPIAIJ=101,
1538d519adbfSMatthew Knepley                MATOP_CONJUGATE=102,
1539d519adbfSMatthew Knepley                MATOP_SET_SIZES=103,
1540d519adbfSMatthew Knepley                MATOP_SET_VALUES_ROW=104,
1541d519adbfSMatthew Knepley                MATOP_REAL_PART=105,
1542d519adbfSMatthew Knepley                MATOP_IMAG_PART=106,
1543d519adbfSMatthew Knepley                MATOP_GET_ROW_UTRIANGULAR=107,
1544d519adbfSMatthew Knepley                MATOP_RESTORE_ROW_UTRIANGULAR=108,
1545d519adbfSMatthew Knepley                MATOP_MATSOLVE=109,
1546d519adbfSMatthew Knepley                MATOP_GET_REDUNDANTMATRIX=110,
1547d519adbfSMatthew Knepley                MATOP_GET_ROW_MIN=111,
1548d519adbfSMatthew Knepley                MATOP_GET_COLUMN_VEC=112,
1549d519adbfSMatthew Knepley                MATOP_MISSING_DIAGONAL=113,
1550d519adbfSMatthew Knepley                MATOP_MATGETSEQNONZEROSTRUCTURE=114,
155189c6957cSBarry Smith                MATOP_CREATE=115,
155289c6957cSBarry Smith                MATOP_GET_GHOSTS=116,
1553ea38565cSJed Brown                MATOP_GET_LOCALSUBMATRIX=117,
1554ea38565cSJed Brown                MATOP_RESTORE_LOCALSUBMATRIX=118,
1555eeffb40dSHong Zhang                MATOP_MULT_DIAGONAL_BLOCK=119,
1556547795f9SHong Zhang                MATOP_HERMITIANTRANSPOSE=120,
1557547795f9SHong Zhang                MATOP_MULTHERMITIANTRANSPOSE=121,
1558fbdbba38SShri Abhyankar                MATOP_MULTHERMITIANTRANSPOSEADD=122,
1559b9614d88SDmitry Karpeev                MATOP_GETMULTIPROCBLOCK=123,
1560b9614d88SDmitry Karpeev 	       MATOP_GET_SUBMATRICES_PARALLEL=128
1561fae171e0SBarry Smith              } MatOperation;
15627087cfbeSBarry Smith extern PetscErrorCode  MatHasOperation(Mat,MatOperation,PetscBool *);
15637087cfbeSBarry Smith extern PetscErrorCode  MatShellSetOperation(Mat,MatOperation,void(*)(void));
15647087cfbeSBarry Smith extern PetscErrorCode  MatShellGetOperation(Mat,MatOperation,void(**)(void));
15657087cfbeSBarry Smith extern PetscErrorCode  MatShellSetContext(Mat,void*);
1566112a2221SBarry Smith 
156790ace30eSBarry Smith /*
156890ace30eSBarry Smith    Codes for matrices stored on disk. By default they are
156990ace30eSBarry Smith    stored in a universal format. By changing the format with
15707973ad04SBarry Smith    PetscViewerSetFormat(viewer,PETSC_VIEWER_NATIVE); the matrices will
157190ace30eSBarry Smith    be stored in a way natural for the matrix, for example dense matrices
157290ace30eSBarry Smith    would be stored as dense. Matrices stored this way may only be
1573f4403165SShri Abhyankar    read into matrices of the same type.
157490ace30eSBarry Smith */
157590ace30eSBarry Smith #define MATRIX_BINARY_FORMAT_DENSE -1
157690ace30eSBarry Smith 
15777087cfbeSBarry Smith extern PetscErrorCode  MatMPIBAIJSetHashTableFactor(Mat,PetscReal);
15787087cfbeSBarry Smith extern PetscErrorCode  MatISGetLocalMat(Mat,Mat*);
15791f4e1ec7SBarry Smith 
1580d9274352SBarry Smith /*S
1581d9274352SBarry Smith      MatNullSpace - Object that removes a null space from a vector, i.e.
1582d9274352SBarry Smith          orthogonalizes the vector to a subsapce
1583d9274352SBarry Smith 
1584f7a9e4ceSBarry Smith    Level: advanced
1585d9274352SBarry Smith 
1586d9274352SBarry Smith   Concepts: matrix; linear operator, null space
1587d9274352SBarry Smith 
15886e1639daSBarry Smith   Users manual sections:
15896e1639daSBarry Smith .   sec_singular
15906e1639daSBarry Smith 
1591d9274352SBarry Smith .seealso:  MatNullSpaceCreate()
1592d9274352SBarry Smith S*/
159374637425SBarry Smith typedef struct _p_MatNullSpace* MatNullSpace;
1594d9274352SBarry Smith 
15957087cfbeSBarry Smith extern PetscErrorCode  MatNullSpaceCreate(MPI_Comm,PetscBool ,PetscInt,const Vec[],MatNullSpace*);
15967087cfbeSBarry Smith extern PetscErrorCode  MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(MatNullSpace,Vec,void*),void*);
15977087cfbeSBarry Smith extern PetscErrorCode  MatNullSpaceDestroy(MatNullSpace);
15987087cfbeSBarry Smith extern PetscErrorCode  MatNullSpaceRemove(MatNullSpace,Vec,Vec*);
15997087cfbeSBarry Smith extern PetscErrorCode  MatNullSpaceAttach(Mat,MatNullSpace);
16007087cfbeSBarry Smith extern PetscErrorCode  MatNullSpaceTest(MatNullSpace,Mat,PetscBool  *);
160174637425SBarry Smith 
16027087cfbeSBarry Smith extern PetscErrorCode  MatReorderingSeqSBAIJ(Mat,IS);
16037087cfbeSBarry Smith extern PetscErrorCode  MatMPISBAIJSetHashTableFactor(Mat,PetscReal);
16047087cfbeSBarry Smith extern PetscErrorCode  MatSeqSBAIJSetColumnIndices(Mat,PetscInt *);
16057087cfbeSBarry Smith extern PetscErrorCode  MatSeqBAIJInvertBlockDiagonal(Mat);
16063f1d51d7SBarry Smith 
16077087cfbeSBarry Smith extern PetscErrorCode  MatCreateMAIJ(Mat,PetscInt,Mat*);
16087087cfbeSBarry Smith extern PetscErrorCode  MatMAIJRedimension(Mat,PetscInt,Mat*);
16097087cfbeSBarry Smith extern PetscErrorCode  MatMAIJGetAIJ(Mat,Mat*);
1610c4f061fbSSatish Balay 
16117087cfbeSBarry Smith extern PetscErrorCode  MatComputeExplicitOperator(Mat,Mat*);
1612b0a32e0cSBarry Smith 
16137087cfbeSBarry Smith extern PetscErrorCode  MatDiagonalScaleLocal(Mat,Vec);
161404f1ad80SBarry Smith 
16157087cfbeSBarry Smith extern PetscErrorCode  MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*);
16167087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetBase(Mat,Vec,Vec);
16177087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*);
16187087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*));
16197087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec));
16207087cfbeSBarry Smith extern PetscErrorCode  MatMFFDAddNullSpace(Mat,MatNullSpace);
16217087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt);
16227087cfbeSBarry Smith extern PetscErrorCode  MatMFFDResetHHistory(Mat);
16237087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetFunctionError(Mat,PetscReal);
16247087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetPeriod(Mat,PetscInt);
16257087cfbeSBarry Smith extern PetscErrorCode  MatMFFDGetH(Mat,PetscScalar *);
16267087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetOptionsPrefix(Mat,const char[]);
16277087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetFromOptions(Mat);
16287087cfbeSBarry Smith extern PetscErrorCode  MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*);
16297087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*);
1630e884886fSBarry Smith 
16316370053bSBarry Smith /*S
16326370053bSBarry Smith     MatMFFD - A data structured used to manage the computation of the h differencing parameter for matrix-free
16336370053bSBarry Smith               Jacobian vector products
1634e884886fSBarry Smith 
16356370053bSBarry Smith     Notes: MATMFFD is a specific MatType which uses the MatMFFD data structure
16366370053bSBarry Smith 
16376370053bSBarry Smith            MatMFFD*() methods actually take the Mat as their first argument. Not a MatMFFD data structure
16386370053bSBarry Smith 
16396370053bSBarry Smith     Level: developer
16406370053bSBarry Smith 
16416370053bSBarry Smith .seealso: MATMFFD, MatCreateMFFD(), MatMFFDSetFuction(), MatMFFDSetType(), MatMFFDRegister()
16426370053bSBarry Smith S*/
1643e884886fSBarry Smith typedef struct _p_MatMFFD* MatMFFD;
1644e884886fSBarry Smith 
1645e884886fSBarry Smith /*E
1646e884886fSBarry Smith     MatMFFDType - algorithm used to compute the h used in computing matrix-vector products via differencing of the function
1647e884886fSBarry Smith 
1648e884886fSBarry Smith    Level: beginner
1649e884886fSBarry Smith 
1650e884886fSBarry Smith .seealso: MatMFFDSetType(), MatMFFDRegister()
1651e884886fSBarry Smith E*/
1652a313700dSBarry Smith #define MatMFFDType char*
1653e884886fSBarry Smith #define MATMFFD_DS  "ds"
1654e884886fSBarry Smith #define MATMFFD_WP  "wp"
1655e884886fSBarry Smith 
16567087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetType(Mat,const MatMFFDType);
16577087cfbeSBarry Smith extern PetscErrorCode  MatMFFDRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatMFFD));
1658e884886fSBarry Smith 
1659e884886fSBarry Smith /*MC
1660e884886fSBarry Smith    MatMFFDRegisterDynamic - Adds a method to the MatMFFD registry.
1661e884886fSBarry Smith 
1662e884886fSBarry Smith    Synopsis:
16631890ba74SBarry Smith    PetscErrorCode MatMFFDRegisterDynamic(const char *name_solver,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatMFFD))
1664e884886fSBarry Smith 
1665e884886fSBarry Smith    Not Collective
1666e884886fSBarry Smith 
1667e884886fSBarry Smith    Input Parameters:
1668e884886fSBarry Smith +  name_solver - name of a new user-defined compute-h module
1669e884886fSBarry Smith .  path - path (either absolute or relative) the library containing this solver
1670e884886fSBarry Smith .  name_create - name of routine to create method context
1671e884886fSBarry Smith -  routine_create - routine to create method context
1672e884886fSBarry Smith 
1673e884886fSBarry Smith    Level: developer
1674e884886fSBarry Smith 
1675e884886fSBarry Smith    Notes:
1676e884886fSBarry Smith    MatMFFDRegisterDynamic() may be called multiple times to add several user-defined solvers.
1677e884886fSBarry Smith 
1678e884886fSBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
1679e884886fSBarry Smith    is ignored.
1680e884886fSBarry Smith 
1681e884886fSBarry Smith    Sample usage:
1682e884886fSBarry Smith .vb
1683e884886fSBarry Smith    MatMFFDRegisterDynamic("my_h",/home/username/my_lib/lib/libO/solaris/mylib.a,
1684e884886fSBarry Smith                "MyHCreate",MyHCreate);
1685e884886fSBarry Smith .ve
1686e884886fSBarry Smith 
1687e884886fSBarry Smith    Then, your solver can be chosen with the procedural interface via
1688e884886fSBarry Smith $     MatMFFDSetType(mfctx,"my_h")
1689e884886fSBarry Smith    or at runtime via the option
1690e884886fSBarry Smith $     -snes_mf_type my_h
1691e884886fSBarry Smith 
1692e884886fSBarry Smith .keywords: MatMFFD, register
1693e884886fSBarry Smith 
1694e884886fSBarry Smith .seealso: MatMFFDRegisterAll(), MatMFFDRegisterDestroy()
1695e884886fSBarry Smith M*/
1696e884886fSBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1697e884886fSBarry Smith #define MatMFFDRegisterDynamic(a,b,c,d) MatMFFDRegister(a,b,c,0)
1698e884886fSBarry Smith #else
1699e884886fSBarry Smith #define MatMFFDRegisterDynamic(a,b,c,d) MatMFFDRegister(a,b,c,d)
1700e884886fSBarry Smith #endif
1701e884886fSBarry Smith 
17027087cfbeSBarry Smith extern PetscErrorCode  MatMFFDRegisterAll(const char[]);
17037087cfbeSBarry Smith extern PetscErrorCode  MatMFFDRegisterDestroy(void);
17047087cfbeSBarry Smith extern PetscErrorCode  MatMFFDDSSetUmin(Mat,PetscReal);
17057087cfbeSBarry Smith extern PetscErrorCode  MatMFFDWPSetComputeNormU(Mat,PetscBool );
1706e884886fSBarry Smith 
1707e884886fSBarry Smith 
17087087cfbeSBarry Smith extern PetscErrorCode  PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *);
17097087cfbeSBarry Smith extern PetscErrorCode  PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *);
17107dbadf16SMatthew Knepley 
171197969023SHong Zhang /*
171297969023SHong Zhang    PETSc interface to MUMPS
171397969023SHong Zhang */
171497969023SHong Zhang #ifdef PETSC_HAVE_MUMPS
1715f250808bSBarry Smith extern PetscErrorCode  MatMumpsSetIcntl(Mat,PetscInt,PetscInt);
171697969023SHong Zhang #endif
171723a5497aSJed Brown 
1718d954db57SHong Zhang /*
1719d954db57SHong Zhang    PETSc interface to SUPERLU
1720d954db57SHong Zhang */
1721d954db57SHong Zhang #ifdef PETSC_HAVE_SUPERLU
1722f250808bSBarry Smith extern PetscErrorCode  MatSuperluSetILUDropTol(Mat,PetscReal);
1723d954db57SHong Zhang #endif
1724d954db57SHong Zhang 
1725c8883902SJed Brown extern PetscErrorCode MatCreateNest(MPI_Comm,PetscInt,const IS[],PetscInt,const IS[],const Mat[],Mat*);
1726c8883902SJed Brown extern PetscErrorCode MatNestGetSize(Mat,PetscInt*,PetscInt*);
1727c8883902SJed Brown extern PetscErrorCode MatNestGetSubMats(Mat,PetscInt*,PetscInt*,Mat***);
1728c8883902SJed Brown extern PetscErrorCode MatNestGetSubMat(Mat,PetscInt,PetscInt,Mat*);
17297087cfbeSBarry Smith extern PetscErrorCode MatNestSetVecType(Mat,const VecType);
1730c8883902SJed Brown extern PetscErrorCode MatNestSetSubMats(Mat,PetscInt,const IS[],PetscInt,const IS[],const Mat[]);
1731d8588912SDave May 
173223a5497aSJed Brown PETSC_EXTERN_CXX_END
173323a5497aSJed Brown #endif
1734