xref: /petsc/include/petscmat.h (revision d34fcf5f3c1110c5dc2577f606d52f97e4a96ba3)
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 
332dedccee8SHong Zhang extern PetscErrorCode  MatCreateFFT(MPI_Comm,PetscInt,const PetscInt[],const MatType,Mat*);
3337087cfbeSBarry Smith extern PetscErrorCode  MatCreateSeqCUFFT(MPI_Comm,PetscInt,const PetscInt[],Mat*);
334dedccee8SHong 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);
344dce485f0SBarry Smith extern PetscErrorCode  MatDestroy_(Mat);
345dce485f0SBarry Smith #define MatDestroy(a)  (MatDestroy_(a) || (((a) = 0),0))
34621c89e3eSBarry Smith 
3477087cfbeSBarry Smith extern PetscErrorCode  MatConjugate(Mat);
3487087cfbeSBarry Smith extern PetscErrorCode  MatRealPart(Mat);
3497087cfbeSBarry Smith extern PetscErrorCode  MatImaginaryPart(Mat);
3507087cfbeSBarry Smith extern PetscErrorCode  MatGetDiagonalBlock(Mat,PetscBool *,MatReuse,Mat*);
3517087cfbeSBarry Smith extern PetscErrorCode  MatGetTrace(Mat,PetscScalar*);
35299cafbc1SBarry Smith 
3538ed539a5SBarry Smith /* ------------------------------------------------------------*/
3547087cfbeSBarry Smith extern PetscErrorCode  MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
3557087cfbeSBarry Smith extern PetscErrorCode  MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
3567087cfbeSBarry Smith extern PetscErrorCode  MatSetValuesRow(Mat,PetscInt,const PetscScalar[]);
3577087cfbeSBarry Smith extern PetscErrorCode  MatSetValuesRowLocal(Mat,PetscInt,const PetscScalar[]);
35884cb2905SBarry Smith 
3592ef4de8bSBarry Smith /*S
3602ef4de8bSBarry Smith      MatStencil - Data structure (C struct) for storing information about a single row or
3612ef4de8bSBarry Smith         column of a matrix as index on an associated grid.
3622ef4de8bSBarry Smith 
3632ef4de8bSBarry Smith    Level: beginner
3642ef4de8bSBarry Smith 
3652ef4de8bSBarry Smith   Concepts: matrix; linear operator
3662ef4de8bSBarry Smith 
367d7bd93baSBarry Smith .seealso:  MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockStencil()
3682ef4de8bSBarry Smith S*/
369435da068SBarry Smith typedef struct {
370c1ac3661SBarry Smith   PetscInt k,j,i,c;
371435da068SBarry Smith } MatStencil;
3722ef4de8bSBarry Smith 
3737087cfbeSBarry Smith extern PetscErrorCode  MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
3747087cfbeSBarry Smith extern PetscErrorCode  MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
3757087cfbeSBarry Smith extern PetscErrorCode  MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt);
376435da068SBarry Smith 
3777087cfbeSBarry Smith extern PetscErrorCode  MatSetColoring(Mat,ISColoring);
3787087cfbeSBarry Smith extern PetscErrorCode  MatSetValuesAdic(Mat,void*);
3797087cfbeSBarry Smith extern PetscErrorCode  MatSetValuesAdifor(Mat,PetscInt,void*);
3803a7fca6bSBarry Smith 
381d91e6319SBarry Smith /*E
382d91e6319SBarry Smith     MatAssemblyType - Indicates if the matrix is now to be used, or if you plan
383d91e6319SBarry Smith      to continue to add values to it
384d91e6319SBarry Smith 
385d91e6319SBarry Smith     Level: beginner
386d91e6319SBarry Smith 
387d91e6319SBarry Smith .seealso: MatAssemblyBegin(), MatAssemblyEnd()
388d91e6319SBarry Smith E*/
3896d4a8577SBarry Smith typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType;
3907087cfbeSBarry Smith extern PetscErrorCode  MatAssemblyBegin(Mat,MatAssemblyType);
3917087cfbeSBarry Smith extern PetscErrorCode  MatAssemblyEnd(Mat,MatAssemblyType);
3927087cfbeSBarry Smith extern PetscErrorCode  MatAssembled(Mat,PetscBool *);
3934f9c727eSBarry Smith 
3941d73ed98SBarry Smith 
39530de9b25SBarry Smith 
396d91e6319SBarry Smith /*E
397d91e6319SBarry Smith     MatOption - Options that may be set for a matrix and its behavior or storage
398d91e6319SBarry Smith 
399d91e6319SBarry Smith     Level: beginner
400d91e6319SBarry Smith 
4010a835dfdSSatish Balay    Any additions/changes here MUST also be made in include/finclude/petscmat.h
402d91e6319SBarry Smith 
403d91e6319SBarry Smith .seealso: MatSetOption()
404d91e6319SBarry Smith E*/
405512a5fc5SBarry Smith typedef enum {MAT_ROW_ORIENTED,MAT_NEW_NONZERO_LOCATIONS,
4064e0d8c25SBarry Smith               MAT_SYMMETRIC,
4074e0d8c25SBarry Smith               MAT_STRUCTURALLY_SYMMETRIC,
4088047e77bSBarry Smith               MAT_NEW_DIAGONALS,MAT_IGNORE_OFF_PROC_ENTRIES,
4094e0d8c25SBarry Smith               MAT_NEW_NONZERO_LOCATION_ERR,
4104e0d8c25SBarry Smith               MAT_NEW_NONZERO_ALLOCATION_ERR,MAT_USE_HASH_TABLE,
411a9817697SBarry Smith               MAT_KEEP_NONZERO_PATTERN,MAT_IGNORE_ZERO_ENTRIES,
4124e0d8c25SBarry Smith               MAT_USE_INODES,
4134e0d8c25SBarry Smith               MAT_HERMITIAN,
4144e0d8c25SBarry Smith               MAT_SYMMETRY_ETERNAL,
415cd6b891eSBarry Smith               MAT_CHECK_COMPRESSED_ROW,
4164e0d8c25SBarry Smith               MAT_IGNORE_LOWER_TRIANGULAR,MAT_ERROR_LOWER_TRIANGULAR,
41728b2fa4aSMatthew Knepley               MAT_GETROW_UPPERTRIANGULAR,MAT_UNUSED_NONZERO_LOCATION_ERR,
4184cb17eb5SBarry Smith               MAT_SPD,MAT_NO_OFF_PROC_ENTRIES,MAT_NO_OFF_PROC_ZERO_ROWS,
41928b2fa4aSMatthew Knepley               NUM_MAT_OPTIONS} MatOption;
420290bbb0aSBarry Smith extern const char *MatOptions[];
4217087cfbeSBarry Smith extern PetscErrorCode  MatSetOption(Mat,MatOption,PetscBool );
4227087cfbeSBarry Smith extern PetscErrorCode  MatGetType(Mat,const MatType*);
423a313700dSBarry Smith PetscPolymorphicFunction(MatGetType,(Mat mat),(mat,&t),const MatType,t)
42484cb2905SBarry Smith 
4257087cfbeSBarry Smith extern PetscErrorCode  MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]);
4267087cfbeSBarry Smith extern PetscErrorCode  MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
4277087cfbeSBarry Smith extern PetscErrorCode  MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
4287087cfbeSBarry Smith extern PetscErrorCode  MatGetRowUpperTriangular(Mat);
4297087cfbeSBarry Smith extern PetscErrorCode  MatRestoreRowUpperTriangular(Mat);
4307087cfbeSBarry Smith extern PetscErrorCode  MatGetColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
4317087cfbeSBarry Smith extern PetscErrorCode  MatRestoreColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
4327087cfbeSBarry Smith extern PetscErrorCode  MatGetColumnVector(Mat,Vec,PetscInt);
4337087cfbeSBarry Smith extern PetscErrorCode  MatGetArray(Mat,PetscScalar *[]);
434ba966639SSatish Balay PetscPolymorphicFunction(MatGetArray,(Mat mat),(mat,&a),PetscScalar*,a)
4357087cfbeSBarry Smith extern PetscErrorCode  MatRestoreArray(Mat,PetscScalar *[]);
4367087cfbeSBarry Smith extern PetscErrorCode  MatGetBlockSize(Mat,PetscInt *);
437ba966639SSatish Balay PetscPolymorphicFunction(MatGetBlockSize,(Mat mat),(mat,&a),PetscInt,a)
4387087cfbeSBarry Smith extern PetscErrorCode  MatSetBlockSize(Mat,PetscInt);
4397b80b807SBarry Smith 
4401620fd73SBarry Smith 
4417087cfbeSBarry Smith extern PetscErrorCode  MatMult(Mat,Vec,Vec);
4427087cfbeSBarry Smith extern PetscErrorCode  MatMultDiagonalBlock(Mat,Vec,Vec);
4437087cfbeSBarry Smith extern PetscErrorCode  MatMultAdd(Mat,Vec,Vec,Vec);
444ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
4457087cfbeSBarry Smith extern PetscErrorCode  MatMultTranspose(Mat,Vec,Vec);
4467087cfbeSBarry Smith extern PetscErrorCode  MatMultHermitianTranspose(Mat,Vec,Vec);
4477087cfbeSBarry Smith extern PetscErrorCode  MatIsTranspose(Mat,Mat,PetscReal,PetscBool *);
448ace3abfcSBarry Smith PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B,PetscReal tol),(A,B,tol,&t),PetscBool ,t)
449ace3abfcSBarry Smith PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B),(A,B,0,&t),PetscBool ,t)
4507087cfbeSBarry Smith extern PetscErrorCode  MatIsHermitianTranspose(Mat,Mat,PetscReal,PetscBool *);
4517087cfbeSBarry Smith extern PetscErrorCode  MatMultTransposeAdd(Mat,Vec,Vec,Vec);
4527087cfbeSBarry Smith extern PetscErrorCode  MatMultHermitianTransposeAdd(Mat,Vec,Vec,Vec);
453ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultTransposeAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
4547087cfbeSBarry Smith extern PetscErrorCode  MatMultConstrained(Mat,Vec,Vec);
4557087cfbeSBarry Smith extern PetscErrorCode  MatMultTransposeConstrained(Mat,Vec,Vec);
4567087cfbeSBarry Smith extern PetscErrorCode  MatMatSolve(Mat,Mat,Mat);
457f5edf698SHong Zhang 
458d91e6319SBarry Smith /*E
459d91e6319SBarry Smith     MatDuplicateOption - Indicates if a duplicated sparse matrix should have
460d91e6319SBarry Smith   its numerical values copied over or just its nonzero structure.
461d91e6319SBarry Smith 
462d91e6319SBarry Smith     Level: beginner
463d91e6319SBarry Smith 
464d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
465d91e6319SBarry Smith 
46670dcbbb9SBarry Smith $   MAT_SHARE_NONZERO_PATTERN - the i and j arrays in the new matrix will be shared with the original matrix
46770dcbbb9SBarry Smith $                               this also triggers the MAT_DO_NOT_COPY_VALUES option. This is used when you
46870dcbbb9SBarry Smith $                               have several matrices with the same nonzero pattern.
46970dcbbb9SBarry Smith 
470d91e6319SBarry Smith .seealso: MatDuplicate()
471d91e6319SBarry Smith E*/
47270dcbbb9SBarry Smith typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES,MAT_SHARE_NONZERO_PATTERN} MatDuplicateOption;
4732e8a6d31SBarry Smith 
4747087cfbeSBarry Smith extern PetscErrorCode  MatConvert(Mat,const MatType,MatReuse,Mat*);
475a313700dSBarry Smith PetscPolymorphicFunction(MatConvert,(Mat A,const MatType t),(A,t,MAT_INITIAL_MATRIX,&a),Mat,a)
4767087cfbeSBarry Smith extern PetscErrorCode  MatDuplicate(Mat,MatDuplicateOption,Mat*);
477ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A,MatDuplicateOption o),(A,o,&a),Mat,a)
478ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A),(A,MAT_COPY_VALUES,&a),Mat,a)
47994a9d846SBarry Smith 
480cb5b572fSBarry Smith 
4817087cfbeSBarry Smith extern PetscErrorCode  MatCopy(Mat,Mat,MatStructure);
4827087cfbeSBarry Smith extern PetscErrorCode  MatView(Mat,PetscViewer);
4837087cfbeSBarry Smith extern PetscErrorCode  MatIsSymmetric(Mat,PetscReal,PetscBool *);
484ace3abfcSBarry Smith PetscPolymorphicFunction(MatIsSymmetric,(Mat A,PetscReal tol),(A,tol,&t),PetscBool ,t)
485ace3abfcSBarry Smith PetscPolymorphicFunction(MatIsSymmetric,(Mat A),(A,0,&t),PetscBool ,t)
4867087cfbeSBarry Smith extern PetscErrorCode  MatIsStructurallySymmetric(Mat,PetscBool *);
487ace3abfcSBarry Smith PetscPolymorphicFunction(MatIsStructurallySymmetric,(Mat A),(A,&t),PetscBool ,t)
4887087cfbeSBarry Smith extern PetscErrorCode  MatIsHermitian(Mat,PetscReal,PetscBool *);
489ace3abfcSBarry Smith PetscPolymorphicFunction(MatIsHermitian,(Mat A),(A,0,&t),PetscBool ,t)
4907087cfbeSBarry Smith extern PetscErrorCode  MatIsSymmetricKnown(Mat,PetscBool *,PetscBool *);
4917087cfbeSBarry Smith extern PetscErrorCode  MatIsHermitianKnown(Mat,PetscBool *,PetscBool *);
4927087cfbeSBarry Smith extern PetscErrorCode  MatMissingDiagonal(Mat,PetscBool  *,PetscInt *);
4937087cfbeSBarry Smith extern PetscErrorCode  MatLoad(Mat, PetscViewer);
4947b80b807SBarry Smith 
4957087cfbeSBarry Smith extern PetscErrorCode  MatGetRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,PetscInt *[],PetscInt *[],PetscBool  *);
4967087cfbeSBarry Smith extern PetscErrorCode  MatRestoreRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,PetscInt *[],PetscInt *[],PetscBool  *);
4977087cfbeSBarry Smith extern PetscErrorCode  MatGetColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,PetscInt *[],PetscInt *[],PetscBool  *);
4987087cfbeSBarry Smith extern PetscErrorCode  MatRestoreColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,PetscInt *[],PetscInt *[],PetscBool  *);
499d4fbbf0eSBarry Smith 
500d91e6319SBarry Smith /*S
501d91e6319SBarry Smith      MatInfo - Context of matrix information, used with MatGetInfo()
502d91e6319SBarry Smith 
503d91e6319SBarry Smith    In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE
504d91e6319SBarry Smith 
505d91e6319SBarry Smith    Level: intermediate
506d91e6319SBarry Smith 
507d91e6319SBarry Smith   Concepts: matrix^nonzero information
508d91e6319SBarry Smith 
509d9274352SBarry Smith .seealso:  MatGetInfo(), MatInfoType
510d91e6319SBarry Smith S*/
5114e220ebcSLois Curfman McInnes typedef struct {
512b0a32e0cSBarry Smith   PetscLogDouble block_size;                         /* block size */
513b0a32e0cSBarry Smith   PetscLogDouble nz_allocated,nz_used,nz_unneeded;   /* number of nonzeros */
514b0a32e0cSBarry Smith   PetscLogDouble memory;                             /* memory allocated */
515b0a32e0cSBarry Smith   PetscLogDouble assemblies;                         /* number of matrix assemblies called */
516b0a32e0cSBarry Smith   PetscLogDouble mallocs;                            /* number of mallocs during MatSetValues() */
517b0a32e0cSBarry Smith   PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */
518b0a32e0cSBarry Smith   PetscLogDouble factor_mallocs;                     /* number of mallocs during factorization */
5194e220ebcSLois Curfman McInnes } MatInfo;
5204e220ebcSLois Curfman McInnes 
521d9274352SBarry Smith /*E
522d9274352SBarry Smith     MatInfoType - Indicates if you want information about the local part of the matrix,
523d9274352SBarry Smith      the entire parallel matrix or the maximum over all the local parts.
524d9274352SBarry Smith 
525d9274352SBarry Smith     Level: beginner
526d9274352SBarry Smith 
527d9274352SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
528d9274352SBarry Smith 
529d9274352SBarry Smith .seealso: MatGetInfo(), MatInfo
530d9274352SBarry Smith E*/
5317b80b807SBarry Smith typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType;
5327087cfbeSBarry Smith extern PetscErrorCode  MatGetInfo(Mat,MatInfoType,MatInfo*);
5337087cfbeSBarry Smith extern PetscErrorCode  MatGetDiagonal(Mat,Vec);
5347087cfbeSBarry Smith extern PetscErrorCode  MatGetRowMax(Mat,Vec,PetscInt[]);
5357087cfbeSBarry Smith extern PetscErrorCode  MatGetRowMin(Mat,Vec,PetscInt[]);
5367087cfbeSBarry Smith extern PetscErrorCode  MatGetRowMaxAbs(Mat,Vec,PetscInt[]);
5377087cfbeSBarry Smith extern PetscErrorCode  MatGetRowMinAbs(Mat,Vec,PetscInt[]);
5387087cfbeSBarry Smith extern PetscErrorCode  MatGetRowSum(Mat,Vec);
5397087cfbeSBarry Smith extern PetscErrorCode  MatTranspose(Mat,MatReuse,Mat*);
540fc4dec0aSBarry Smith PetscPolymorphicFunction(MatTranspose,(Mat A),(A,MAT_INITIAL_MATRIX,&t),Mat,t)
5417087cfbeSBarry Smith extern PetscErrorCode  MatHermitianTranspose(Mat,MatReuse,Mat*);
5427087cfbeSBarry Smith extern PetscErrorCode  MatPermute(Mat,IS,IS,Mat *);
543ba966639SSatish Balay PetscPolymorphicFunction(MatPermute,(Mat A,IS is1,IS is2),(A,is1,is2,&t),Mat,t)
5447087cfbeSBarry Smith extern PetscErrorCode  MatDiagonalScale(Mat,Vec,Vec);
5457087cfbeSBarry Smith extern PetscErrorCode  MatDiagonalSet(Mat,Vec,InsertMode);
5467087cfbeSBarry Smith extern PetscErrorCode  MatEqual(Mat,Mat,PetscBool *);
547ace3abfcSBarry Smith PetscPolymorphicFunction(MatEqual,(Mat A,Mat B),(A,B,&t),PetscBool ,t)
5487087cfbeSBarry Smith extern PetscErrorCode  MatMultEqual(Mat,Mat,PetscInt,PetscBool *);
5497087cfbeSBarry Smith extern PetscErrorCode  MatMultAddEqual(Mat,Mat,PetscInt,PetscBool *);
5507087cfbeSBarry Smith extern PetscErrorCode  MatMultTransposeEqual(Mat,Mat,PetscInt,PetscBool *);
5517087cfbeSBarry Smith extern PetscErrorCode  MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscBool *);
5527b80b807SBarry Smith 
5537087cfbeSBarry Smith extern PetscErrorCode  MatNorm(Mat,NormType,PetscReal *);
554ba966639SSatish Balay PetscPolymorphicFunction(MatNorm,(Mat A,NormType t),(A,t,&n),PetscReal,n)
55509573ac7SBarry Smith extern PetscErrorCode MatGetColumnNorms(Mat,NormType,PetscReal *);
5567087cfbeSBarry Smith extern PetscErrorCode  MatZeroEntries(Mat);
5577087cfbeSBarry Smith extern PetscErrorCode  MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
5587087cfbeSBarry Smith extern PetscErrorCode  MatZeroRowsIS(Mat,IS,PetscScalar,Vec,Vec);
5597087cfbeSBarry Smith extern PetscErrorCode  MatZeroRowsStencil(Mat,PetscInt,const MatStencil [],PetscScalar,Vec,Vec);
5607087cfbeSBarry Smith extern PetscErrorCode  MatZeroRowsColumns(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
5617087cfbeSBarry Smith extern PetscErrorCode  MatZeroRowsColumnsIS(Mat,IS,PetscScalar,Vec,Vec);
5627b80b807SBarry Smith 
5637087cfbeSBarry Smith extern PetscErrorCode  MatUseScaledForm(Mat,PetscBool );
5647087cfbeSBarry Smith extern PetscErrorCode  MatScaleSystem(Mat,Vec,Vec);
5657087cfbeSBarry Smith extern PetscErrorCode  MatUnScaleSystem(Mat,Vec,Vec);
5665ef9f2a5SBarry Smith 
5677087cfbeSBarry Smith extern PetscErrorCode  MatGetSize(Mat,PetscInt*,PetscInt*);
5687087cfbeSBarry Smith extern PetscErrorCode  MatGetLocalSize(Mat,PetscInt*,PetscInt*);
5697087cfbeSBarry Smith extern PetscErrorCode  MatGetOwnershipRange(Mat,PetscInt*,PetscInt*);
5707087cfbeSBarry Smith extern PetscErrorCode  MatGetOwnershipRanges(Mat,const PetscInt**);
5717087cfbeSBarry Smith extern PetscErrorCode  MatGetOwnershipRangeColumn(Mat,PetscInt*,PetscInt*);
5727087cfbeSBarry Smith extern PetscErrorCode  MatGetOwnershipRangesColumn(Mat,const PetscInt**);
5737b80b807SBarry Smith 
5747087cfbeSBarry Smith extern PetscErrorCode  MatGetSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
5757087cfbeSBarry Smith extern PetscErrorCode  MatGetSubMatricesParallel(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
5767087cfbeSBarry Smith extern PetscErrorCode  MatDestroyMatrices(PetscInt,Mat *[]);
5777087cfbeSBarry Smith extern PetscErrorCode  MatGetSubMatrix(Mat,IS,IS,MatReuse,Mat *);
5787087cfbeSBarry Smith extern PetscErrorCode  MatGetLocalSubMatrix(Mat,IS,IS,Mat*);
5797087cfbeSBarry Smith extern PetscErrorCode  MatRestoreLocalSubMatrix(Mat,IS,IS,Mat*);
5807087cfbeSBarry Smith extern PetscErrorCode  MatGetSeqNonzeroStructure(Mat,Mat*);
5817087cfbeSBarry Smith extern PetscErrorCode  MatDestroySeqNonzeroStructure(Mat*);
5825494a064SHong Zhang 
5837087cfbeSBarry Smith extern PetscErrorCode  MatMerge(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
5847087cfbeSBarry Smith extern PetscErrorCode  MatMerge_SeqsToMPI(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*);
5857087cfbeSBarry Smith extern PetscErrorCode  MatMerge_SeqsToMPISymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*);
5867087cfbeSBarry Smith extern PetscErrorCode  MatMerge_SeqsToMPINumeric(Mat,Mat);
5877087cfbeSBarry Smith extern PetscErrorCode  MatDestroy_MPIAIJ_SeqsToMPI(Mat);
5887087cfbeSBarry Smith extern PetscErrorCode  MatGetLocalMat(Mat,MatReuse,Mat*);
5897087cfbeSBarry Smith extern PetscErrorCode  MatGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*);
5907087cfbeSBarry Smith extern PetscErrorCode  MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,PetscInt*,Mat*);
5917087cfbeSBarry Smith extern PetscErrorCode  MatGetBrowsOfAoCols(Mat,Mat,MatReuse,PetscInt**,PetscInt**,MatScalar**,Mat*);
59243eb5e2fSMatthew Knepley #if defined (PETSC_USE_CTABLE)
593cd116e26SSatish Balay #include "petscctable.h"
5947087cfbeSBarry Smith extern PetscErrorCode  MatGetCommunicationStructs(Mat, Vec *, PetscTable *, VecScatter *);
59543eb5e2fSMatthew Knepley #else
5967087cfbeSBarry Smith extern PetscErrorCode  MatGetCommunicationStructs(Mat, Vec *, PetscInt *[], VecScatter *);
59743eb5e2fSMatthew Knepley #endif
5987087cfbeSBarry Smith extern PetscErrorCode  MatGetGhosts(Mat, PetscInt *,const PetscInt *[]);
5998efafbd8SBarry Smith 
6007087cfbeSBarry Smith extern PetscErrorCode  MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt);
6017b80b807SBarry Smith 
6027087cfbeSBarry Smith extern PetscErrorCode  MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*);
6037087cfbeSBarry Smith extern PetscErrorCode  MatMatMultSymbolic(Mat,Mat,PetscReal,Mat*);
6047087cfbeSBarry Smith extern PetscErrorCode  MatMatMultNumeric(Mat,Mat,Mat);
60522440eb1SKris Buschelman 
6067087cfbeSBarry Smith extern PetscErrorCode  MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*);
6077087cfbeSBarry Smith extern PetscErrorCode  MatPtAPSymbolic(Mat,Mat,PetscReal,Mat*);
6087087cfbeSBarry Smith extern PetscErrorCode  MatPtAPNumeric(Mat,Mat,Mat);
60922440eb1SKris Buschelman 
6107087cfbeSBarry Smith extern PetscErrorCode  MatMatMultTranspose(Mat,Mat,MatReuse,PetscReal,Mat*);
6117087cfbeSBarry Smith extern PetscErrorCode  MatMatMultTransposeSymbolic(Mat,Mat,PetscReal,Mat*);
6127087cfbeSBarry Smith extern PetscErrorCode  MatMatMultTransposeNumeric(Mat,Mat,Mat);
613bc011b1eSHong Zhang 
6147087cfbeSBarry Smith extern PetscErrorCode  MatAXPY(Mat,PetscScalar,Mat,MatStructure);
6157087cfbeSBarry Smith extern PetscErrorCode  MatAYPX(Mat,PetscScalar,Mat,MatStructure);
6161c741599SHong Zhang 
6177087cfbeSBarry Smith extern PetscErrorCode  MatScale(Mat,PetscScalar);
6187087cfbeSBarry Smith extern PetscErrorCode  MatShift(Mat,PetscScalar);
6197b80b807SBarry Smith 
6207087cfbeSBarry Smith extern PetscErrorCode  MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping);
6217087cfbeSBarry Smith extern PetscErrorCode  MatSetLocalToGlobalMappingBlock(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping);
6227087cfbeSBarry Smith extern PetscErrorCode  MatGetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*);
6237087cfbeSBarry Smith extern PetscErrorCode  MatGetLocalToGlobalMappingBlock(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*);
6247087cfbeSBarry Smith extern PetscErrorCode  MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
6257087cfbeSBarry Smith extern PetscErrorCode  MatZeroRowsLocalIS(Mat,IS,PetscScalar,Vec,Vec);
6267087cfbeSBarry Smith extern PetscErrorCode  MatZeroRowsColumnsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
6277087cfbeSBarry Smith extern PetscErrorCode  MatZeroRowsColumnsLocalIS(Mat,IS,PetscScalar,Vec,Vec);
6287087cfbeSBarry Smith extern PetscErrorCode  MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
6297087cfbeSBarry Smith extern PetscErrorCode  MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
630052efed2SBarry Smith 
6317087cfbeSBarry Smith extern PetscErrorCode  MatStashSetInitialSize(Mat,PetscInt,PetscInt);
6327087cfbeSBarry Smith extern PetscErrorCode  MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
63390f02eecSBarry Smith 
6347087cfbeSBarry Smith extern PetscErrorCode  MatInterpolate(Mat,Vec,Vec);
6357087cfbeSBarry Smith extern PetscErrorCode  MatInterpolateAdd(Mat,Vec,Vec,Vec);
636ba966639SSatish Balay PetscPolymorphicSubroutine(MatInterpolateAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
6377087cfbeSBarry Smith extern PetscErrorCode  MatRestrict(Mat,Vec,Vec);
6387087cfbeSBarry Smith extern PetscErrorCode  MatGetVecs(Mat,Vec*,Vec*);
6397087cfbeSBarry Smith extern PetscErrorCode  MatGetRedundantMatrix(Mat,PetscInt,MPI_Comm,PetscInt,MatReuse,Mat*);
6407087cfbeSBarry Smith extern PetscErrorCode  MatGetMultiProcBlock(Mat,MPI_Comm,Mat*);
641fc9bc008SSatish Balay extern PetscErrorCode  MatFindZeroDiagonals(Mat,IS*);
642bd481603SBarry Smith 
643bd481603SBarry Smith /*MC
6441620fd73SBarry Smith    MatSetValue - Set a single entry into a matrix.
6451620fd73SBarry Smith 
6461620fd73SBarry Smith    Not collective
6471620fd73SBarry Smith 
6481620fd73SBarry Smith    Input Parameters:
6491620fd73SBarry Smith +  m - the matrix
6501620fd73SBarry Smith .  row - the row location of the entry
6511620fd73SBarry Smith .  col - the column location of the entry
6521620fd73SBarry Smith .  value - the value to insert
6531620fd73SBarry Smith -  mode - either INSERT_VALUES or ADD_VALUES
6541620fd73SBarry Smith 
6551620fd73SBarry Smith    Notes:
6561620fd73SBarry Smith    For efficiency one should use MatSetValues() and set several or many
6571620fd73SBarry Smith    values simultaneously if possible.
6581620fd73SBarry Smith 
6591620fd73SBarry Smith    Level: beginner
6601620fd73SBarry Smith 
6611620fd73SBarry Smith .seealso: MatSetValues(), MatSetValueLocal()
6621620fd73SBarry Smith M*/
6631620fd73SBarry 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);}
6641620fd73SBarry Smith 
6651620fd73SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatGetValue(Mat v,PetscInt i,PetscInt j,PetscScalar *va) {return MatGetValues(v,1,&i,1,&j,va);}
6661620fd73SBarry Smith 
6671620fd73SBarry 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);}
6681620fd73SBarry Smith 
6691620fd73SBarry Smith /*MC
6700d2a7ff2SSatish Balay    MatPreallocateInitialize - Begins the block of code that will count the number of nonzeros per
671bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
672bd481603SBarry Smith 
673bd481603SBarry Smith    Synopsis:
674c1ac3661SBarry Smith    PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)
675bd481603SBarry Smith 
676bd481603SBarry Smith    Collective on MPI_Comm
677bd481603SBarry Smith 
678bd481603SBarry Smith    Input Parameters:
679bd481603SBarry Smith +  comm - the communicator that will share the eventually allocated matrix
680859fcb39SBarry Smith .  nrows - the number of LOCAL rows in the matrix
681859fcb39SBarry Smith -  ncols - the number of LOCAL columns in the matrix
682bd481603SBarry Smith 
683bd481603SBarry Smith    Output Parameters:
684bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
685bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
686bd481603SBarry Smith 
687bd481603SBarry Smith 
688bd481603SBarry Smith    Level: intermediate
689bd481603SBarry Smith 
690bd481603SBarry Smith    Notes:
6910598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
692bd481603SBarry Smith 
6931d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
694bd481603SBarry Smith 
695bd481603SBarry Smith    Use MatPreallocateInitializeSymmetric() for symmetric matrices (MPISBAIJ matrices)
696bd481603SBarry Smith 
6971620fd73SBarry Smith    This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize().
6981620fd73SBarry Smith 
699bd481603SBarry Smith   Concepts: preallocation^Matrix
700bd481603SBarry Smith 
701bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
702bd481603SBarry Smith           MatPreallocateInitializeSymmetric(), MatPreallocateSymmetricSetLocal()
703bd481603SBarry Smith M*/
704c4f061fbSSatish Balay #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \
7057c922b88SBarry Smith { \
706a89bc211SBarry Smith   PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__start,__end; \
707a89bc211SBarry Smith   _4_ierr = PetscMalloc2(__nrows,PetscInt,&dnz,__nrows,PetscInt,&onz);CHKERRQ(_4_ierr); \
708a89bc211SBarry Smith   _4_ierr = PetscMemzero(dnz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
709a89bc211SBarry Smith   _4_ierr = PetscMemzero(onz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
710c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ctmp;\
711a89bc211SBarry Smith   _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows;
7127c922b88SBarry Smith 
713bd481603SBarry Smith /*MC
7140d2a7ff2SSatish Balay    MatPreallocateSymmetricInitialize - Begins the block of code that will count the number of nonzeros per
715bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
716bd481603SBarry Smith 
717bd481603SBarry Smith    Synopsis:
718c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)
719bd481603SBarry Smith 
720bd481603SBarry Smith    Collective on MPI_Comm
721bd481603SBarry Smith 
722bd481603SBarry Smith    Input Parameters:
723bd481603SBarry Smith +  comm - the communicator that will share the eventually allocated matrix
724859fcb39SBarry Smith .  nrows - the number of LOCAL rows in the matrix
725859fcb39SBarry Smith -  ncols - the number of LOCAL columns in the matrix
726bd481603SBarry Smith 
727bd481603SBarry Smith    Output Parameters:
728bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
729bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
730bd481603SBarry Smith 
731bd481603SBarry Smith 
732bd481603SBarry Smith    Level: intermediate
733bd481603SBarry Smith 
734bd481603SBarry Smith    Notes:
7350598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
736bd481603SBarry Smith 
7371d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
738bd481603SBarry Smith 
7391620fd73SBarry Smith    This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize().
7401620fd73SBarry Smith 
741bd481603SBarry Smith   Concepts: preallocation^Matrix
742bd481603SBarry Smith 
743bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
744bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal()
745bd481603SBarry Smith M*/
746222b16d4SBarry Smith #define MatPreallocateSymmetricInitialize(comm,nrows,ncols,dnz,onz) 0; \
747222b16d4SBarry Smith { \
748a89bc211SBarry Smith   PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__end; \
749a89bc211SBarry Smith   _4_ierr = PetscMalloc2(__nrows,PetscInt,&dnz,__nrows,PetscInt,&onz);CHKERRQ(_4_ierr); \
750a89bc211SBarry Smith   _4_ierr = PetscMemzero(dnz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
751a89bc211SBarry Smith   _4_ierr = PetscMemzero(onz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
752c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr);\
753a89bc211SBarry Smith   _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows;
754222b16d4SBarry Smith 
755bd481603SBarry Smith /*MC
756bd481603SBarry Smith    MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
757bd481603SBarry Smith        inserted using a local number of the rows and columns
758bd481603SBarry Smith 
759bd481603SBarry Smith    Synopsis:
760c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
761bd481603SBarry Smith 
762bd481603SBarry Smith    Not Collective
763bd481603SBarry Smith 
764bd481603SBarry Smith    Input Parameters:
765784ac674SJed Brown +  map - the row mapping from local numbering to global numbering
766bd481603SBarry Smith .  nrows - the number of rows indicated
7671d73ed98SBarry Smith .  rows - the indices of the rows
768784ac674SJed Brown .  cmap - the column mapping from local to global numbering
769bd481603SBarry Smith .  ncols - the number of columns in the matrix
770bd481603SBarry Smith .  cols - the columns indicated
771bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
772bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
773bd481603SBarry Smith 
774bd481603SBarry Smith 
775bd481603SBarry Smith    Level: intermediate
776bd481603SBarry Smith 
777bd481603SBarry Smith    Notes:
7780598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
779bd481603SBarry Smith 
7801d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
781bd481603SBarry Smith 
782bd481603SBarry Smith   Concepts: preallocation^Matrix
783bd481603SBarry Smith 
784bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
785bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal()
786bd481603SBarry Smith M*/
787784ac674SJed Brown #define MatPreallocateSetLocal(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \
788c4f061fbSSatish Balay {\
789c1ac3661SBarry Smith   PetscInt __l;\
790784ac674SJed Brown   _4_ierr = ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\
791784ac674SJed Brown   _4_ierr = ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\
792c4f061fbSSatish Balay   for (__l=0;__l<nrows;__l++) {\
793ef66eb69SBarry Smith     _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
794c4f061fbSSatish Balay   }\
795c4f061fbSSatish Balay }
796c4f061fbSSatish Balay 
797bd481603SBarry Smith /*MC
798bd481603SBarry Smith    MatPreallocateSymmetricSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
799bd481603SBarry Smith        inserted using a local number of the rows and columns
800bd481603SBarry Smith 
801bd481603SBarry Smith    Synopsis:
802c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
803bd481603SBarry Smith 
804bd481603SBarry Smith    Not Collective
805bd481603SBarry Smith 
806bd481603SBarry Smith    Input Parameters:
807bd481603SBarry Smith +  map - the mapping between local numbering and global numbering
808bd481603SBarry Smith .  nrows - the number of rows indicated
8091d73ed98SBarry Smith .  rows - the indices of the rows
810bd481603SBarry Smith .  ncols - the number of columns in the matrix
811bd481603SBarry Smith .  cols - the columns indicated
812bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
813bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
814bd481603SBarry Smith 
815bd481603SBarry Smith 
816bd481603SBarry Smith    Level: intermediate
817bd481603SBarry Smith 
818bd481603SBarry Smith    Notes:
8190598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
820bd481603SBarry Smith 
821bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
822bd481603SBarry Smith 
823bd481603SBarry Smith   Concepts: preallocation^Matrix
824bd481603SBarry Smith 
825bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
826bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
827bd481603SBarry Smith M*/
828d3d32019SBarry Smith #define MatPreallocateSymmetricSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\
829d3d32019SBarry Smith {\
830c1ac3661SBarry Smith   PetscInt __l;\
831d3d32019SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\
832d3d32019SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\
833d3d32019SBarry Smith   for (__l=0;__l<nrows;__l++) {\
834d3d32019SBarry Smith     _4_ierr = MatPreallocateSymmetricSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
835d3d32019SBarry Smith   }\
836d3d32019SBarry Smith }
837d3d32019SBarry Smith 
838bd481603SBarry Smith /*MC
839bd481603SBarry Smith    MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
840bd481603SBarry Smith        inserted using a local number of the rows and columns
841bd481603SBarry Smith 
842bd481603SBarry Smith    Synopsis:
843c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
844bd481603SBarry Smith 
845bd481603SBarry Smith    Not Collective
846bd481603SBarry Smith 
847bd481603SBarry Smith    Input Parameters:
84864075487SBarry Smith +  row - the row
849bd481603SBarry Smith .  ncols - the number of columns in the matrix
850a50f8bf6SHong Zhang -  cols - the columns indicated
851a50f8bf6SHong Zhang 
852a50f8bf6SHong Zhang    Output Parameters:
853a50f8bf6SHong Zhang +  dnz - the array that will be passed to the matrix preallocation routines
854bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
855bd481603SBarry Smith 
856bd481603SBarry Smith 
857bd481603SBarry Smith    Level: intermediate
858bd481603SBarry Smith 
859bd481603SBarry Smith    Notes:
8600598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
861bd481603SBarry Smith 
862bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
863bd481603SBarry Smith 
8641620fd73SBarry Smith    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
8651620fd73SBarry Smith 
866bd481603SBarry Smith   Concepts: preallocation^Matrix
867bd481603SBarry Smith 
868bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
869bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
870bd481603SBarry Smith M*/
871c4f061fbSSatish Balay #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\
872c1ac3661SBarry Smith { PetscInt __i; \
873e32f2f54SBarry 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);\
874e32f2f54SBarry 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);\
8757c922b88SBarry Smith   for (__i=0; __i<nc; __i++) {\
87664075487SBarry Smith     if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \
8777cc688d7SBarry Smith     else dnz[row - __rstart]++;\
8787c922b88SBarry Smith   }\
8797c922b88SBarry Smith }
8807c922b88SBarry Smith 
881bd481603SBarry Smith /*MC
882bd481603SBarry Smith    MatPreallocateSymmetricSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
883bd481603SBarry Smith        inserted using a local number of the rows and columns
884bd481603SBarry Smith 
885bd481603SBarry Smith    Synopsis:
886c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
887bd481603SBarry Smith 
888bd481603SBarry Smith    Not Collective
889bd481603SBarry Smith 
890bd481603SBarry Smith    Input Parameters:
891bd481603SBarry Smith +  nrows - the number of rows indicated
8921d73ed98SBarry Smith .  rows - the indices of the rows
893bd481603SBarry Smith .  ncols - the number of columns in the matrix
894bd481603SBarry Smith .  cols - the columns indicated
895bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
896bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
897bd481603SBarry Smith 
898bd481603SBarry Smith 
899bd481603SBarry Smith    Level: intermediate
900bd481603SBarry Smith 
901bd481603SBarry Smith    Notes:
9020598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
903bd481603SBarry Smith 
904bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
905bd481603SBarry Smith 
9061620fd73SBarry Smith    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
9071620fd73SBarry Smith 
908bd481603SBarry Smith   Concepts: preallocation^Matrix
909bd481603SBarry Smith 
910bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
911bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
912bd481603SBarry Smith M*/
913d3d32019SBarry Smith #define MatPreallocateSymmetricSet(row,nc,cols,dnz,onz) 0;\
914c1ac3661SBarry Smith { PetscInt __i; \
915d3d32019SBarry Smith   for (__i=0; __i<nc; __i++) {\
916d3d32019SBarry Smith     if (cols[__i] >= __end) onz[row - __rstart]++; \
917d3d32019SBarry Smith     else if (cols[__i] >= row) dnz[row - __rstart]++;\
918d3d32019SBarry Smith   }\
919d3d32019SBarry Smith }
920d3d32019SBarry Smith 
921bd481603SBarry Smith /*MC
92216371a99SBarry Smith    MatPreallocateLocation -  An alternative to MatPreallocationSet() that puts the nonzero locations into the matrix if it exists
92316371a99SBarry Smith 
92416371a99SBarry Smith    Synopsis:
92516371a99SBarry Smith    PetscErrorCode MatPreallocateLocations(Mat A,PetscInt row,PetscInt ncols,PetscInt *cols,PetscInt *dnz,PetscInt *onz)
92616371a99SBarry Smith 
92716371a99SBarry Smith    Not Collective
92816371a99SBarry Smith 
92916371a99SBarry Smith    Input Parameters:
93016371a99SBarry Smith .  A - matrix
93116371a99SBarry Smith .  row - row where values exist (must be local to this process)
93216371a99SBarry Smith .  ncols - number of columns
93316371a99SBarry Smith .  cols - columns with nonzeros
93416371a99SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
93516371a99SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
93616371a99SBarry Smith 
93716371a99SBarry Smith 
93816371a99SBarry Smith    Level: intermediate
93916371a99SBarry Smith 
94016371a99SBarry Smith    Notes:
9410598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
94216371a99SBarry Smith 
94316371a99SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
94416371a99SBarry Smith 
94516371a99SBarry Smith    This is a MACRO not a function because it uses a bunch of variables private to the MatPreallocation.... routines.
94616371a99SBarry Smith 
94716371a99SBarry Smith   Concepts: preallocation^Matrix
94816371a99SBarry Smith 
94916371a99SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
95016371a99SBarry Smith           MatPreallocateSymmetricInitialize(), MatPreallocateSymmetricSetLocal()
95116371a99SBarry Smith M*/
95216371a99SBarry 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);}
95316371a99SBarry Smith 
95416371a99SBarry Smith 
95516371a99SBarry Smith /*MC
9560d2a7ff2SSatish Balay    MatPreallocateFinalize - Ends the block of code that will count the number of nonzeros per
957bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
958bd481603SBarry Smith 
959bd481603SBarry Smith    Synopsis:
960c1ac3661SBarry Smith    PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz)
961bd481603SBarry Smith 
962bd481603SBarry Smith    Collective on MPI_Comm
963bd481603SBarry Smith 
964bd481603SBarry Smith    Input Parameters:
96516371a99SBarry Smith +  dnz - the array that was be passed to the matrix preallocation routines
966bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
967bd481603SBarry Smith 
968bd481603SBarry Smith 
969bd481603SBarry Smith    Level: intermediate
970bd481603SBarry Smith 
971bd481603SBarry Smith    Notes:
9720598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
973bd481603SBarry Smith 
974bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
975bd481603SBarry Smith 
9761620fd73SBarry Smith    This is a MACRO not a function because it closes the { started in MatPreallocateInitialize().
9771620fd73SBarry Smith 
978bd481603SBarry Smith   Concepts: preallocation^Matrix
979bd481603SBarry Smith 
980bd481603SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
981bd481603SBarry Smith           MatPreallocateSymmetricInitialize(), MatPreallocateSymmetricSetLocal()
982bd481603SBarry Smith M*/
983a89bc211SBarry Smith #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree2(dnz,onz);CHKERRQ(_4_ierr);}
9847c922b88SBarry Smith 
985bd481603SBarry Smith 
986bd481603SBarry Smith 
9877b80b807SBarry Smith /* Routines unique to particular data structures */
9887087cfbeSBarry Smith extern PetscErrorCode  MatShellGetContext(Mat,void **);
989ba966639SSatish Balay PetscPolymorphicFunction(MatShellGetContext,(Mat A),(A,&t),void*,t)
990698d4c6aSKris Buschelman 
9917087cfbeSBarry Smith extern PetscErrorCode  MatInodeAdjustForInodes(Mat,IS*,IS*);
9927087cfbeSBarry Smith extern PetscErrorCode  MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *);
993698d4c6aSKris Buschelman 
9947087cfbeSBarry Smith extern PetscErrorCode  MatSeqAIJSetColumnIndices(Mat,PetscInt[]);
9957087cfbeSBarry Smith extern PetscErrorCode  MatSeqBAIJSetColumnIndices(Mat,PetscInt[]);
9967087cfbeSBarry Smith extern PetscErrorCode  MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
9977087cfbeSBarry Smith extern PetscErrorCode  MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
9987087cfbeSBarry Smith extern PetscErrorCode  MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
9997b80b807SBarry Smith 
1000a96a251dSBarry Smith #define MAT_SKIP_ALLOCATION -4
1001a96a251dSBarry Smith 
10027087cfbeSBarry Smith extern PetscErrorCode  MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
1003ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz))
10047087cfbeSBarry Smith extern PetscErrorCode  MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
1005ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqSBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz))
10067087cfbeSBarry Smith extern PetscErrorCode  MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]);
1007ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqAIJSetPreallocation,(Mat A,const PetscInt nnz[]),(A,0,nnz))
1008273d9f13SBarry Smith 
10097087cfbeSBarry Smith extern PetscErrorCode  MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
1010ba966639SSatish Balay PetscPolymorphicSubroutine(MatMPIBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[],const PetscInt onz[]),(A,bs,0,nnz,0,onz))
10117087cfbeSBarry Smith extern PetscErrorCode  MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
10127087cfbeSBarry Smith extern PetscErrorCode  MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
10137087cfbeSBarry Smith extern PetscErrorCode  MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []);
10147087cfbeSBarry Smith extern PetscErrorCode  MatSeqBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
10157087cfbeSBarry Smith extern PetscErrorCode  MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]);
10167087cfbeSBarry Smith extern PetscErrorCode  MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
10177087cfbeSBarry Smith extern PetscErrorCode  MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]);
10187087cfbeSBarry Smith extern PetscErrorCode  MatMPIDenseSetPreallocation(Mat,PetscScalar[]);
10197087cfbeSBarry Smith extern PetscErrorCode  MatSeqDenseSetPreallocation(Mat,PetscScalar[]);
10207087cfbeSBarry Smith extern PetscErrorCode  MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,PetscInt*[]);
10217087cfbeSBarry Smith extern PetscErrorCode  MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,PetscInt*[]);
10227087cfbeSBarry Smith extern PetscErrorCode  MatAdicSetLocalFunction(Mat,void (*)(void));
1023273d9f13SBarry Smith 
10247087cfbeSBarry Smith extern PetscErrorCode  MatSeqDenseSetLDA(Mat,PetscInt);
10257087cfbeSBarry Smith extern PetscErrorCode  MatDenseGetLocalMatrix(Mat,Mat*);
10261b807ce4Svictorle 
10277087cfbeSBarry Smith extern PetscErrorCode  MatStoreValues(Mat);
10287087cfbeSBarry Smith extern PetscErrorCode  MatRetrieveValues(Mat);
10292e8a6d31SBarry Smith 
10307087cfbeSBarry Smith extern PetscErrorCode  MatDAADSetCtx(Mat,void*);
10313a7fca6bSBarry Smith 
1032b3a44c85SBarry Smith extern PetscErrorCode  MatFindNonzeroRows(Mat,IS*);
10337b80b807SBarry Smith /*
10347b80b807SBarry Smith   These routines are not usually accessed directly, rather solving is
103594b7f48cSBarry Smith   done through the KSP and PC interfaces.
10367b80b807SBarry Smith */
10377b80b807SBarry Smith 
1038d9274352SBarry Smith /*E
1039d9274352SBarry Smith     MatOrderingType - String with the name of a PETSc matrix ordering or the creation function
1040d9274352SBarry Smith        with an optional dynamic library name, for example
1041d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:orderingcreate()
1042d9274352SBarry Smith 
1043d9274352SBarry Smith    Level: beginner
1044d9274352SBarry Smith 
1045e5a9bf91SBarry Smith    Cannot use const because the PC objects manipulate the string
1046e5a9bf91SBarry Smith 
1047d9274352SBarry Smith .seealso: MatGetOrdering()
1048d9274352SBarry Smith E*/
10493eaad2caSSatish Balay #define MatOrderingType char*
10502692d6eeSBarry Smith #define MATORDERINGNATURAL     "natural"
10512692d6eeSBarry Smith #define MATORDERINGND          "nd"
10522692d6eeSBarry Smith #define MATORDERING1WD         "1wd"
10532692d6eeSBarry Smith #define MATORDERINGRCM         "rcm"
10542692d6eeSBarry Smith #define MATORDERINGQMD         "qmd"
10552692d6eeSBarry Smith #define MATORDERINGROWLENGTH   "rowlength"
10562692d6eeSBarry Smith #define MATORDERINGDSC_ND      "dsc_nd"         /* these three are only for DSCPACK, see its documentation for details */
10572692d6eeSBarry Smith #define MATORDERINGDSC_MMD     "dsc_mmd"
10582692d6eeSBarry Smith #define MATORDERINGDSC_MDF     "dsc_mdf"
1059312e372aSBarry Smith #define MATORDERINGAMD         "amd"            /* only works if UMFPACK is installed with PETSc */
1060b12f92e5SBarry Smith 
10617087cfbeSBarry Smith extern PetscErrorCode  MatGetOrdering(Mat,const MatOrderingType,IS*,IS*);
10627087cfbeSBarry Smith extern PetscErrorCode  MatGetOrderingList(PetscFList *list);
10637087cfbeSBarry Smith extern PetscErrorCode  MatOrderingRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,const MatOrderingType,IS*,IS*));
106430de9b25SBarry Smith 
106530de9b25SBarry Smith /*MC
10661890ba74SBarry Smith    MatOrderingRegisterDynamic - Adds a new sparse matrix ordering to the matrix package.
106730de9b25SBarry Smith 
106830de9b25SBarry Smith    Synopsis:
10691890ba74SBarry Smith    PetscErrorCode MatOrderingRegisterDynamic(const char *name_ordering,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatOrdering))
107030de9b25SBarry Smith 
107130de9b25SBarry Smith    Not Collective
107230de9b25SBarry Smith 
107330de9b25SBarry Smith    Input Parameters:
10742692d6eeSBarry Smith +  sname - name of ordering (for example MATORDERINGND)
107530de9b25SBarry Smith .  path - location of library where creation routine is
107630de9b25SBarry Smith .  name - name of function that creates the ordering type,a string
107730de9b25SBarry Smith -  function - function pointer that creates the ordering
107830de9b25SBarry Smith 
107930de9b25SBarry Smith    Level: developer
108030de9b25SBarry Smith 
108130de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
108230de9b25SBarry Smith    is ignored.
108330de9b25SBarry Smith 
108430de9b25SBarry Smith    Sample usage:
108530de9b25SBarry Smith .vb
108630de9b25SBarry Smith    MatOrderingRegisterDynamic("my_order",/home/username/my_lib/lib/libO/solaris/mylib.a,
108730de9b25SBarry Smith                "MyOrder",MyOrder);
108830de9b25SBarry Smith .ve
108930de9b25SBarry Smith 
109030de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
109130de9b25SBarry Smith $     MatOrderingSetType(part,"my_order)
109230de9b25SBarry Smith    or at runtime via the option
10932401956bSBarry Smith $     -pc_factor_mat_ordering_type my_order
109430de9b25SBarry Smith 
1095ab901514SBarry Smith    ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values.
109630de9b25SBarry Smith 
109730de9b25SBarry Smith .keywords: matrix, ordering, register
109830de9b25SBarry Smith 
109930de9b25SBarry Smith .seealso: MatOrderingRegisterDestroy(), MatOrderingRegisterAll()
110030de9b25SBarry Smith M*/
1101aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1102f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,0)
1103b12f92e5SBarry Smith #else
1104f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,d)
1105b12f92e5SBarry Smith #endif
110630de9b25SBarry Smith 
11077087cfbeSBarry Smith extern PetscErrorCode  MatOrderingRegisterDestroy(void);
11087087cfbeSBarry Smith extern PetscErrorCode  MatOrderingRegisterAll(const char[]);
1109ace3abfcSBarry Smith extern PetscBool  MatOrderingRegisterAllCalled;
1110b0a32e0cSBarry Smith extern PetscFList MatOrderingList;
1111d4fbbf0eSBarry Smith 
11127087cfbeSBarry Smith extern PetscErrorCode  MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS);
1113a2ce50c7SBarry Smith 
1114d91e6319SBarry Smith /*S
1115d90ac83dSHong Zhang     MatFactorShiftType - Numeric Shift.
1116d90ac83dSHong Zhang 
1117d90ac83dSHong Zhang    Level: beginner
1118d90ac83dSHong Zhang 
1119d90ac83dSHong Zhang S*/
1120d90ac83dSHong Zhang typedef enum {MAT_SHIFT_NONE,MAT_SHIFT_NONZERO,MAT_SHIFT_POSITIVE_DEFINITE,MAT_SHIFT_INBLOCKS} MatFactorShiftType;
1121d90ac83dSHong Zhang extern const char *MatFactorShiftTypes[];
1122d90ac83dSHong Zhang 
1123d90ac83dSHong Zhang /*S
11242401956bSBarry Smith    MatFactorInfo - Data passed into the matrix factorization routines
1125b00f7748SHong Zhang 
112661cad860SBarry Smith    In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE, that is use
112761cad860SBarry Smith $     MatFactorInfo  info(MAT_FACTORINFO_SIZE)
1128b00f7748SHong Zhang 
112915e8a5b3SHong Zhang    Notes: These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC.
1130b00f7748SHong Zhang 
113161cad860SBarry Smith       You can use MatFactorInfoInitialize() to set default values.
113261cad860SBarry Smith 
1133b00f7748SHong Zhang    Level: developer
1134b00f7748SHong Zhang 
1135d7d82daaSBarry Smith .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(),
1136d7d82daaSBarry Smith           MatFactorInfoInitialize()
1137b00f7748SHong Zhang 
1138b00f7748SHong Zhang S*/
1139b00f7748SHong Zhang typedef struct {
114015e8a5b3SHong Zhang   PetscReal     diagonal_fill;  /* force diagonal to fill in if initially not filled */
114185317021SBarry Smith   PetscReal     usedt;
114215e8a5b3SHong Zhang   PetscReal     dt;             /* drop tolerance */
1143b00f7748SHong Zhang   PetscReal     dtcol;          /* tolerance for pivoting */
114415e8a5b3SHong Zhang   PetscReal     dtcount;        /* maximum nonzeros to be allowed per row */
114567e28bfeSBarry Smith   PetscReal     fill;           /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */
1146348344bbSBarry Smith   PetscReal     levels;         /* ICC/ILU(levels) */
1147bcd9e38bSBarry Smith   PetscReal     pivotinblocks;  /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0
1148bcd9e38bSBarry Smith                                    factorization may be faster if do not pivot */
114915e8a5b3SHong Zhang   PetscReal     zeropivot;      /* pivot is called zero if less than this */
1150f4db908eSBarry Smith   PetscReal     shifttype;      /* type of shift added to matrix factor to prevent zero pivots */
1151f4db908eSBarry Smith   PetscReal     shiftamount;     /* how large the shift is */
115215e8a5b3SHong Zhang } MatFactorInfo;
1153ffa6d0a5SLois Curfman McInnes 
11547087cfbeSBarry Smith extern PetscErrorCode  MatFactorInfoInitialize(MatFactorInfo*);
11557087cfbeSBarry Smith extern PetscErrorCode  MatCholeskyFactor(Mat,IS,const MatFactorInfo*);
11567087cfbeSBarry Smith extern PetscErrorCode  MatCholeskyFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
11577087cfbeSBarry Smith extern PetscErrorCode  MatCholeskyFactorNumeric(Mat,Mat,const MatFactorInfo*);
11587087cfbeSBarry Smith extern PetscErrorCode  MatLUFactor(Mat,IS,IS,const MatFactorInfo*);
11597087cfbeSBarry Smith extern PetscErrorCode  MatILUFactor(Mat,IS,IS,const MatFactorInfo*);
11607087cfbeSBarry Smith extern PetscErrorCode  MatLUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
11617087cfbeSBarry Smith extern PetscErrorCode  MatILUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
11627087cfbeSBarry Smith extern PetscErrorCode  MatICCFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
11637087cfbeSBarry Smith extern PetscErrorCode  MatICCFactor(Mat,IS,const MatFactorInfo*);
11647087cfbeSBarry Smith extern PetscErrorCode  MatLUFactorNumeric(Mat,Mat,const MatFactorInfo*);
11657087cfbeSBarry Smith extern PetscErrorCode  MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*);
11667087cfbeSBarry Smith extern PetscErrorCode  MatSolve(Mat,Vec,Vec);
11677087cfbeSBarry Smith extern PetscErrorCode  MatForwardSolve(Mat,Vec,Vec);
11687087cfbeSBarry Smith extern PetscErrorCode  MatBackwardSolve(Mat,Vec,Vec);
11697087cfbeSBarry Smith extern PetscErrorCode  MatSolveAdd(Mat,Vec,Vec,Vec);
11707087cfbeSBarry Smith extern PetscErrorCode  MatSolveTranspose(Mat,Vec,Vec);
11717087cfbeSBarry Smith extern PetscErrorCode  MatSolveTransposeAdd(Mat,Vec,Vec,Vec);
11727087cfbeSBarry Smith extern PetscErrorCode  MatSolves(Mat,Vecs,Vecs);
11738ed539a5SBarry Smith 
11747087cfbeSBarry Smith extern PetscErrorCode  MatSetUnfactored(Mat);
1175bb5a7306SBarry Smith 
1176d91e6319SBarry Smith /*E
1177d91e6319SBarry Smith     MatSORType - What type of (S)SOR to perform
1178bb1eb677SSatish Balay 
1179d91e6319SBarry Smith     Level: beginner
1180d91e6319SBarry Smith 
1181d9274352SBarry Smith    May be bitwise ORd together
1182d9274352SBarry Smith 
1183d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
1184d91e6319SBarry Smith 
11854e7234bfSBarry Smith    MatSORType may be bitwise ORd together, so do not change the numbers
11864e7234bfSBarry Smith 
118741f059aeSBarry Smith .seealso: MatSOR()
1188d91e6319SBarry Smith E*/
1189ee50ffe9SBarry Smith typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3,
1190ee50ffe9SBarry Smith               SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8,
1191ee50ffe9SBarry Smith               SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16,
119284cb2905SBarry Smith               SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType;
11937087cfbeSBarry Smith extern PetscErrorCode  MatSOR(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
11948ed539a5SBarry Smith 
1195d4fbbf0eSBarry Smith /*
1196639f9d9dSBarry Smith     These routines are for efficiently computing Jacobians via finite differences.
1197639f9d9dSBarry Smith */
1198b12f92e5SBarry Smith 
1199d9274352SBarry Smith /*E
1200d9274352SBarry Smith     MatColoringType - String with the name of a PETSc matrix coloring or the creation function
1201d9274352SBarry Smith        with an optional dynamic library name, for example
1202d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:coloringcreate()
1203d9274352SBarry Smith 
1204d9274352SBarry Smith    Level: beginner
1205d9274352SBarry Smith 
1206d9274352SBarry Smith .seealso: MatGetColoring()
1207d9274352SBarry Smith E*/
1208a313700dSBarry Smith #define MatColoringType char*
12092692d6eeSBarry Smith #define MATCOLORINGNATURAL "natural"
12102692d6eeSBarry Smith #define MATCOLORINGSL      "sl"
12112692d6eeSBarry Smith #define MATCOLORINGLF      "lf"
12122692d6eeSBarry Smith #define MATCOLORINGID      "id"
1213b12f92e5SBarry Smith 
12147087cfbeSBarry Smith extern PetscErrorCode  MatGetColoring(Mat,const MatColoringType,ISColoring*);
12157087cfbeSBarry Smith extern PetscErrorCode  MatColoringRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,MatColoringType,ISColoring *));
121630de9b25SBarry Smith 
121730de9b25SBarry Smith /*MC
121830de9b25SBarry Smith    MatColoringRegisterDynamic - Adds a new sparse matrix coloring to the
121930de9b25SBarry Smith                                matrix package.
122030de9b25SBarry Smith 
122130de9b25SBarry Smith    Synopsis:
12221890ba74SBarry Smith    PetscErrorCode MatColoringRegisterDynamic(const char *name_coloring,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatColoring))
122330de9b25SBarry Smith 
122430de9b25SBarry Smith    Not Collective
122530de9b25SBarry Smith 
122630de9b25SBarry Smith    Input Parameters:
12272692d6eeSBarry Smith +  sname - name of Coloring (for example MATCOLORINGSL)
122830de9b25SBarry Smith .  path - location of library where creation routine is
122930de9b25SBarry Smith .  name - name of function that creates the Coloring type, a string
123030de9b25SBarry Smith -  function - function pointer that creates the coloring
123130de9b25SBarry Smith 
123230de9b25SBarry Smith    Level: developer
123330de9b25SBarry Smith 
123430de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
123530de9b25SBarry Smith    is ignored.
123630de9b25SBarry Smith 
123730de9b25SBarry Smith    Sample usage:
123830de9b25SBarry Smith .vb
123930de9b25SBarry Smith    MatColoringRegisterDynamic("my_color",/home/username/my_lib/lib/libO/solaris/mylib.a,
124030de9b25SBarry Smith                "MyColor",MyColor);
124130de9b25SBarry Smith .ve
124230de9b25SBarry Smith 
124330de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
124430de9b25SBarry Smith $     MatColoringSetType(part,"my_color")
124530de9b25SBarry Smith    or at runtime via the option
124630de9b25SBarry Smith $     -mat_coloring_type my_color
124730de9b25SBarry Smith 
1248ab901514SBarry Smith    $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
124930de9b25SBarry Smith 
125030de9b25SBarry Smith .keywords: matrix, Coloring, register
125130de9b25SBarry Smith 
125230de9b25SBarry Smith .seealso: MatColoringRegisterDestroy(), MatColoringRegisterAll()
125330de9b25SBarry Smith M*/
1254aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1255f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,0)
1256b12f92e5SBarry Smith #else
1257f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,d)
1258b12f92e5SBarry Smith #endif
125930de9b25SBarry Smith 
1260ace3abfcSBarry Smith extern PetscBool  MatColoringRegisterAllCalled;
1261f1144a30SSatish Balay 
12627087cfbeSBarry Smith extern PetscErrorCode  MatColoringRegisterAll(const char[]);
12637087cfbeSBarry Smith extern PetscErrorCode  MatColoringRegisterDestroy(void);
12647087cfbeSBarry Smith extern PetscErrorCode  MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*);
1265639f9d9dSBarry Smith 
1266d9274352SBarry Smith /*S
1267d9274352SBarry Smith      MatFDColoring - Object for computing a sparse Jacobian via finite differences
1268d9274352SBarry Smith         and coloring
1269639f9d9dSBarry Smith 
1270d9274352SBarry Smith    Level: beginner
1271d9274352SBarry Smith 
1272d9274352SBarry Smith   Concepts: coloring, sparse Jacobian, finite differences
1273d9274352SBarry Smith 
1274d9274352SBarry Smith .seealso:  MatFDColoringCreate()
1275d9274352SBarry Smith S*/
1276e2a1c21fSSatish Balay typedef struct _p_MatFDColoring* MatFDColoring;
1277639f9d9dSBarry Smith 
12787087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringCreate(Mat,ISColoring,MatFDColoring *);
12797087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringDestroy(MatFDColoring);
12807087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringView(MatFDColoring,PetscViewer);
12817087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*);
12827087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**);
12837087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal);
12847087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringSetFromOptions(MatFDColoring);
12857087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringApply(Mat,MatFDColoring,Vec,MatStructure*,void *);
12867087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringApplyTS(Mat,MatFDColoring,PetscReal,Vec,MatStructure*,void *);
12877087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringSetF(MatFDColoring,Vec);
12887087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,PetscInt*[]);
1289639f9d9dSBarry Smith /*
12900752156aSBarry Smith     These routines are for partitioning matrices: currently used only
12913eda8832SBarry Smith   for adjacency matrix, MatCreateMPIAdj().
12920752156aSBarry Smith */
1293ca161407SBarry Smith 
1294d9274352SBarry Smith /*S
1295d9274352SBarry Smith      MatPartitioning - Object for managing the partitioning of a matrix or graph
1296d9274352SBarry Smith 
1297d9274352SBarry Smith    Level: beginner
1298d9274352SBarry Smith 
1299d9274352SBarry Smith   Concepts: partitioning
1300d9274352SBarry Smith 
1301743a1026Svictorle .seealso:  MatPartitioningCreate(), MatPartitioningType
1302d9274352SBarry Smith S*/
130391e9ee9fSBarry Smith typedef struct _p_MatPartitioning* MatPartitioning;
1304d9274352SBarry Smith 
1305d9274352SBarry Smith /*E
13065bc6e7ccSvictorle     MatPartitioningType - String with the name of a PETSc matrix partitioning or the creation function
1307d9274352SBarry Smith        with an optional dynamic library name, for example
1308d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:partitioningcreate()
1309d9274352SBarry Smith 
1310d9274352SBarry Smith    Level: beginner
1311d9274352SBarry Smith 
1312b547a13bSvictorle .seealso: MatPartitioningCreate(), MatPartitioning
1313d9274352SBarry Smith E*/
1314a313700dSBarry Smith #define MatPartitioningType char*
13152692d6eeSBarry Smith #define MATPARTITIONINGCURRENT  "current"
13162692d6eeSBarry Smith #define MATPARTITIONINGSQUARE   "square"
13172692d6eeSBarry Smith #define MATPARTITIONINGPARMETIS "parmetis"
13182692d6eeSBarry Smith #define MATPARTITIONINGCHACO    "chaco"
13192692d6eeSBarry Smith #define MATPARTITIONINGJOSTLE   "jostle"
13202692d6eeSBarry Smith #define MATPARTITIONINGPARTY    "party"
13212692d6eeSBarry Smith #define MATPARTITIONINGSCOTCH   "scotch"
132271306c60Sjroman 
1323ca161407SBarry Smith 
13247087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningCreate(MPI_Comm,MatPartitioning*);
13257087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningSetType(MatPartitioning,const MatPartitioningType);
13267087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningSetNParts(MatPartitioning,PetscInt);
13277087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningSetAdjacency(MatPartitioning,Mat);
13287087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]);
13297087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []);
13307087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningApply(MatPartitioning,IS*);
13317087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningDestroy(MatPartitioning);
13322aabb6bbSBarry Smith 
13337087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatPartitioning));
133430de9b25SBarry Smith 
133530de9b25SBarry Smith /*MC
133630de9b25SBarry Smith    MatPartitioningRegisterDynamic - Adds a new sparse matrix partitioning to the
133730de9b25SBarry Smith    matrix package.
133830de9b25SBarry Smith 
133930de9b25SBarry Smith    Synopsis:
13401890ba74SBarry Smith    PetscErrorCode MatPartitioningRegisterDynamic(const char *name_partitioning,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatPartitioning))
134130de9b25SBarry Smith 
134230de9b25SBarry Smith    Not Collective
134330de9b25SBarry Smith 
134430de9b25SBarry Smith    Input Parameters:
13452692d6eeSBarry Smith +  sname - name of partitioning (for example MATPARTITIONINGCURRENT) or parmetis
134630de9b25SBarry Smith .  path - location of library where creation routine is
134730de9b25SBarry Smith .  name - name of function that creates the partitioning type, a string
134830de9b25SBarry Smith -  function - function pointer that creates the partitioning type
134930de9b25SBarry Smith 
135030de9b25SBarry Smith    Level: developer
135130de9b25SBarry Smith 
135230de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
135330de9b25SBarry Smith    is ignored.
135430de9b25SBarry Smith 
135530de9b25SBarry Smith    Sample usage:
135630de9b25SBarry Smith .vb
135730de9b25SBarry Smith    MatPartitioningRegisterDynamic("my_part",/home/username/my_lib/lib/libO/solaris/mylib.a,
135830de9b25SBarry Smith                "MyPartCreate",MyPartCreate);
135930de9b25SBarry Smith .ve
136030de9b25SBarry Smith 
136130de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
136230de9b25SBarry Smith $     MatPartitioningSetType(part,"my_part")
136330de9b25SBarry Smith    or at runtime via the option
136430de9b25SBarry Smith $     -mat_partitioning_type my_part
136530de9b25SBarry Smith 
1366ab901514SBarry Smith    $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
136730de9b25SBarry Smith 
136830de9b25SBarry Smith .keywords: matrix, partitioning, register
136930de9b25SBarry Smith 
137030de9b25SBarry Smith .seealso: MatPartitioningRegisterDestroy(), MatPartitioningRegisterAll()
137130de9b25SBarry Smith M*/
1372aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1373f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,0)
13742aabb6bbSBarry Smith #else
1375f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,d)
13762aabb6bbSBarry Smith #endif
13772aabb6bbSBarry Smith 
1378ace3abfcSBarry Smith extern PetscBool  MatPartitioningRegisterAllCalled;
1379f1144a30SSatish Balay 
13807087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningRegisterAll(const char[]);
13817087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningRegisterDestroy(void);
13822bad1931SBarry Smith 
13837087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningView(MatPartitioning,PetscViewer);
13847087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningSetFromOptions(MatPartitioning);
13857087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningGetType(MatPartitioning,const MatPartitioningType*);
1386ca161407SBarry Smith 
13877087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningParmetisSetCoarseSequential(MatPartitioning);
13887087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningParmetisGetEdgeCut(MatPartitioning, PetscInt *);
13890752156aSBarry Smith 
13907087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningJostleSetCoarseLevel(MatPartitioning,PetscReal);
13917087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningJostleSetCoarseSequential(MatPartitioning);
139271306c60Sjroman 
1393954b3ec6SBarry Smith typedef enum {MP_CHACO_MULTILEVEL_KL,MP_CHACO_SPECTRAL,MP_CHACO_LINEAR,MP_CHACO_RANDOM, MP_CHACO_SCATTERED} MPChacoGlobalType;
13947087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningChacoSetGlobal(MatPartitioning, MPChacoGlobalType);
139571306c60Sjroman typedef enum { MP_CHACO_KERNIGHAN_LIN, MP_CHACO_NONE } MPChacoLocalType;
13967087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningChacoSetLocal(MatPartitioning, MPChacoLocalType);
13977087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal);
139871306c60Sjroman typedef enum { MP_CHACO_LANCZOS, MP_CHACO_RQI_SYMMLQ } MPChacoEigenType;
13997087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType);
14007087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningChacoSetEigenTol(MatPartitioning, PetscReal);
14017087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningChacoSetEigenNumber(MatPartitioning, PetscInt);
140271306c60Sjroman 
140371306c60Sjroman #define MP_PARTY_OPT "opt"
140471306c60Sjroman #define MP_PARTY_LIN "lin"
140571306c60Sjroman #define MP_PARTY_SCA "sca"
140671306c60Sjroman #define MP_PARTY_RAN "ran"
140771306c60Sjroman #define MP_PARTY_GBF "gbf"
140871306c60Sjroman #define MP_PARTY_GCF "gcf"
140971306c60Sjroman #define MP_PARTY_BUB "bub"
141071306c60Sjroman #define MP_PARTY_DEF "def"
14117087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningPartySetGlobal(MatPartitioning, const char*);
141271306c60Sjroman #define MP_PARTY_HELPFUL_SETS "hs"
141371306c60Sjroman #define MP_PARTY_KERNIGHAN_LIN "kl"
141471306c60Sjroman #define MP_PARTY_NONE "no"
14157087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningPartySetLocal(MatPartitioning, const char*);
14167087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal);
14177087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningPartySetBipart(MatPartitioning,PetscBool );
14187087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscBool );
141971306c60Sjroman 
142071306c60Sjroman typedef enum { MP_SCOTCH_GREEDY, MP_SCOTCH_GPS, MP_SCOTCH_GR_GPS } MPScotchGlobalType;
14217087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningScotchSetArch(MatPartitioning,const char*);
14227087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningScotchSetMultilevel(MatPartitioning);
14237087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningScotchSetGlobal(MatPartitioning,MPScotchGlobalType);
14247087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningScotchSetCoarseLevel(MatPartitioning,PetscReal);
14257087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningScotchSetHostList(MatPartitioning,const char*);
142671306c60Sjroman typedef enum { MP_SCOTCH_KERNIGHAN_LIN, MP_SCOTCH_NONE } MPScotchLocalType;
14277087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningScotchSetLocal(MatPartitioning,MPScotchLocalType);
14287087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningScotchSetMapping(MatPartitioning);
14297087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningScotchSetStrategy(MatPartitioning,char*);
143071306c60Sjroman 
143109573ac7SBarry Smith extern PetscErrorCode MatMeshToVertexGraph(Mat,PetscInt,Mat*);
143209573ac7SBarry Smith extern PetscErrorCode MatMeshToCellGraph(Mat,PetscInt,Mat*);
1433591294e4SBarry Smith 
14340752156aSBarry Smith /*
14350a835dfdSSatish Balay     If you add entries here you must also add them to finclude/petscmat.h
1436d4fbbf0eSBarry Smith */
14371c1c02c0SLois Curfman McInnes typedef enum { MATOP_SET_VALUES=0,
14381c1c02c0SLois Curfman McInnes                MATOP_GET_ROW=1,
14391c1c02c0SLois Curfman McInnes                MATOP_RESTORE_ROW=2,
14401c1c02c0SLois Curfman McInnes                MATOP_MULT=3,
14411c1c02c0SLois Curfman McInnes                MATOP_MULT_ADD=4,
14427c922b88SBarry Smith                MATOP_MULT_TRANSPOSE=5,
14437c922b88SBarry Smith                MATOP_MULT_TRANSPOSE_ADD=6,
14441c1c02c0SLois Curfman McInnes                MATOP_SOLVE=7,
14451c1c02c0SLois Curfman McInnes                MATOP_SOLVE_ADD=8,
14467c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE=9,
14477c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE_ADD=10,
14481c1c02c0SLois Curfman McInnes                MATOP_LUFACTOR=11,
14491c1c02c0SLois Curfman McInnes                MATOP_CHOLESKYFACTOR=12,
1450a714c06dSBarry Smith                MATOP_SOR=13,
14511c1c02c0SLois Curfman McInnes                MATOP_TRANSPOSE=14,
14521c1c02c0SLois Curfman McInnes                MATOP_GETINFO=15,
14531c1c02c0SLois Curfman McInnes                MATOP_EQUAL=16,
14541c1c02c0SLois Curfman McInnes                MATOP_GET_DIAGONAL=17,
14551c1c02c0SLois Curfman McInnes                MATOP_DIAGONAL_SCALE=18,
14561c1c02c0SLois Curfman McInnes                MATOP_NORM=19,
14571c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_BEGIN=20,
14581c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_END=21,
1459d519adbfSMatthew Knepley                MATOP_SET_OPTION=22,
1460d519adbfSMatthew Knepley                MATOP_ZERO_ENTRIES=23,
1461d519adbfSMatthew Knepley                MATOP_ZERO_ROWS=24,
1462d519adbfSMatthew Knepley                MATOP_LUFACTOR_SYMBOLIC=25,
1463d519adbfSMatthew Knepley                MATOP_LUFACTOR_NUMERIC=26,
1464d519adbfSMatthew Knepley                MATOP_CHOLESKY_FACTOR_SYMBOLIC=27,
1465d519adbfSMatthew Knepley                MATOP_CHOLESKY_FACTOR_NUMERIC=28,
1466d519adbfSMatthew Knepley                MATOP_SETUP_PREALLOCATION=29,
1467d519adbfSMatthew Knepley                MATOP_ILUFACTOR_SYMBOLIC=30,
1468d519adbfSMatthew Knepley                MATOP_ICCFACTOR_SYMBOLIC=31,
1469d519adbfSMatthew Knepley                MATOP_GET_ARRAY=32,
1470d519adbfSMatthew Knepley                MATOP_RESTORE_ARRAY=33,
1471d519adbfSMatthew Knepley                MATOP_DUPLICATE=34,
1472d519adbfSMatthew Knepley                MATOP_FORWARD_SOLVE=35,
1473d519adbfSMatthew Knepley                MATOP_BACKWARD_SOLVE=36,
1474d519adbfSMatthew Knepley                MATOP_ILUFACTOR=37,
1475d519adbfSMatthew Knepley                MATOP_ICCFACTOR=38,
1476d519adbfSMatthew Knepley                MATOP_AXPY=39,
1477d519adbfSMatthew Knepley                MATOP_GET_SUBMATRICES=40,
1478d519adbfSMatthew Knepley                MATOP_INCREASE_OVERLAP=41,
1479d519adbfSMatthew Knepley                MATOP_GET_VALUES=42,
1480d519adbfSMatthew Knepley                MATOP_COPY=43,
1481d519adbfSMatthew Knepley                MATOP_GET_ROW_MAX=44,
1482d519adbfSMatthew Knepley                MATOP_SCALE=45,
1483d519adbfSMatthew Knepley                MATOP_SHIFT=46,
148435153367SBarry Smith                MATOP_DIAGONAL_SET=47,
1485d519adbfSMatthew Knepley                MATOP_ILUDT_FACTOR=48,
1486d519adbfSMatthew Knepley                MATOP_SET_BLOCK_SIZE=49,
1487d519adbfSMatthew Knepley                MATOP_GET_ROW_IJ=50,
1488d519adbfSMatthew Knepley                MATOP_RESTORE_ROW_IJ=51,
1489d519adbfSMatthew Knepley                MATOP_GET_COLUMN_IJ=52,
1490d519adbfSMatthew Knepley                MATOP_RESTORE_COLUMN_IJ=53,
1491d519adbfSMatthew Knepley                MATOP_FDCOLORING_CREATE=54,
1492d519adbfSMatthew Knepley                MATOP_COLORING_PATCH=55,
1493d519adbfSMatthew Knepley                MATOP_SET_UNFACTORED=56,
1494d519adbfSMatthew Knepley                MATOP_PERMUTE=57,
1495d519adbfSMatthew Knepley                MATOP_SET_VALUES_BLOCKED=58,
1496d519adbfSMatthew Knepley                MATOP_GET_SUBMATRIX=59,
1497d519adbfSMatthew Knepley                MATOP_DESTROY=60,
1498d519adbfSMatthew Knepley                MATOP_VIEW=61,
1499d519adbfSMatthew Knepley                MATOP_CONVERT_FROM=62,
1500d519adbfSMatthew Knepley                MATOP_USE_SCALED_FORM=63,
1501d519adbfSMatthew Knepley                MATOP_SCALE_SYSTEM=64,
1502d519adbfSMatthew Knepley                MATOP_UNSCALE_SYSTEM=65,
1503d519adbfSMatthew Knepley                MATOP_SET_LOCAL_TO_GLOBAL_MAP=66,
1504d519adbfSMatthew Knepley                MATOP_SET_VALUES_LOCAL=67,
1505d519adbfSMatthew Knepley                MATOP_ZERO_ROWS_LOCAL=68,
1506d519adbfSMatthew Knepley                MATOP_GET_ROW_MAX_ABS=69,
1507d519adbfSMatthew Knepley                MATOP_GET_ROW_MIN_ABS=70,
1508d519adbfSMatthew Knepley                MATOP_CONVERT=71,
1509d519adbfSMatthew Knepley                MATOP_SET_COLORING=72,
1510d519adbfSMatthew Knepley                MATOP_SET_VALUES_ADIC=73,
1511d519adbfSMatthew Knepley                MATOP_SET_VALUES_ADIFOR=74,
1512d519adbfSMatthew Knepley                MATOP_FD_COLORING_APPLY=75,
1513d519adbfSMatthew Knepley                MATOP_SET_FROM_OPTIONS=76,
1514d519adbfSMatthew Knepley                MATOP_MULT_CON=77,
1515d519adbfSMatthew Knepley                MATOP_MULT_TRANSPOSE_CON=78,
1516d519adbfSMatthew Knepley                MATOP_PERMUTE_SPARSIFY=79,
1517d519adbfSMatthew Knepley                MATOP_MULT_MULTIPLE=80,
1518d519adbfSMatthew Knepley                MATOP_SOLVE_MULTIPLE=81,
1519d519adbfSMatthew Knepley                MATOP_GET_INERTIA=82,
1520d519adbfSMatthew Knepley                MATOP_LOAD=83,
1521d519adbfSMatthew Knepley                MATOP_IS_SYMMETRIC=84,
1522d519adbfSMatthew Knepley                MATOP_IS_HERMITIAN=85,
1523d519adbfSMatthew Knepley                MATOP_IS_STRUCTURALLY_SYMMETRIC=86,
152441f059aeSBarry Smith                MATOP_DUMMY=87,
1525d519adbfSMatthew Knepley                MATOP_GET_VECS=88,
1526d519adbfSMatthew Knepley                MATOP_MAT_MULT=89,
1527d519adbfSMatthew Knepley                MATOP_MAT_MULT_SYMBOLIC=90,
1528d519adbfSMatthew Knepley                MATOP_MAT_MULT_NUMERIC=91,
1529d519adbfSMatthew Knepley                MATOP_PTAP=92,
1530d519adbfSMatthew Knepley                MATOP_PTAP_SYMBOLIC=93,
1531d519adbfSMatthew Knepley                MATOP_PTAP_NUMERIC=94,
1532d519adbfSMatthew Knepley                MATOP_MAT_MULTTRANSPOSE=95,
15331763ac4eSSatish Balay                MATOP_MAT_MULTTRANSPOSE_SYM=96,
15341763ac4eSSatish Balay                MATOP_MAT_MULTTRANSPOSE_NUM=97,
1535d519adbfSMatthew Knepley                MATOP_PTAP_SYMBOLIC_SEQAIJ=98,
1536d519adbfSMatthew Knepley                MATOP_PTAP_NUMERIC_SEQAIJ=99,
1537d519adbfSMatthew Knepley                MATOP_PTAP_SYMBOLIC_MPIAIJ=100,
1538d519adbfSMatthew Knepley                MATOP_PTAP_NUMERIC_MPIAIJ=101,
1539d519adbfSMatthew Knepley                MATOP_CONJUGATE=102,
1540d519adbfSMatthew Knepley                MATOP_SET_SIZES=103,
1541d519adbfSMatthew Knepley                MATOP_SET_VALUES_ROW=104,
1542d519adbfSMatthew Knepley                MATOP_REAL_PART=105,
1543d519adbfSMatthew Knepley                MATOP_IMAG_PART=106,
1544d519adbfSMatthew Knepley                MATOP_GET_ROW_UTRIANGULAR=107,
1545d519adbfSMatthew Knepley                MATOP_RESTORE_ROW_UTRIANGULAR=108,
1546d519adbfSMatthew Knepley                MATOP_MATSOLVE=109,
1547d519adbfSMatthew Knepley                MATOP_GET_REDUNDANTMATRIX=110,
1548d519adbfSMatthew Knepley                MATOP_GET_ROW_MIN=111,
1549d519adbfSMatthew Knepley                MATOP_GET_COLUMN_VEC=112,
1550d519adbfSMatthew Knepley                MATOP_MISSING_DIAGONAL=113,
1551d519adbfSMatthew Knepley                MATOP_MATGETSEQNONZEROSTRUCTURE=114,
155289c6957cSBarry Smith                MATOP_CREATE=115,
155389c6957cSBarry Smith                MATOP_GET_GHOSTS=116,
1554ea38565cSJed Brown                MATOP_GET_LOCALSUBMATRIX=117,
1555ea38565cSJed Brown                MATOP_RESTORE_LOCALSUBMATRIX=118,
1556eeffb40dSHong Zhang                MATOP_MULT_DIAGONAL_BLOCK=119,
1557547795f9SHong Zhang                MATOP_HERMITIANTRANSPOSE=120,
1558547795f9SHong Zhang                MATOP_MULTHERMITIANTRANSPOSE=121,
1559fbdbba38SShri Abhyankar                MATOP_MULTHERMITIANTRANSPOSEADD=122,
1560b9614d88SDmitry Karpeev                MATOP_GETMULTIPROCBLOCK=123,
1561b9614d88SDmitry Karpeev 	       MATOP_GET_SUBMATRICES_PARALLEL=128
1562fae171e0SBarry Smith              } MatOperation;
15637087cfbeSBarry Smith extern PetscErrorCode  MatHasOperation(Mat,MatOperation,PetscBool *);
15647087cfbeSBarry Smith extern PetscErrorCode  MatShellSetOperation(Mat,MatOperation,void(*)(void));
15657087cfbeSBarry Smith extern PetscErrorCode  MatShellGetOperation(Mat,MatOperation,void(**)(void));
15667087cfbeSBarry Smith extern PetscErrorCode  MatShellSetContext(Mat,void*);
1567112a2221SBarry Smith 
156890ace30eSBarry Smith /*
156990ace30eSBarry Smith    Codes for matrices stored on disk. By default they are
157090ace30eSBarry Smith    stored in a universal format. By changing the format with
15717973ad04SBarry Smith    PetscViewerSetFormat(viewer,PETSC_VIEWER_NATIVE); the matrices will
157290ace30eSBarry Smith    be stored in a way natural for the matrix, for example dense matrices
157390ace30eSBarry Smith    would be stored as dense. Matrices stored this way may only be
1574f4403165SShri Abhyankar    read into matrices of the same type.
157590ace30eSBarry Smith */
157690ace30eSBarry Smith #define MATRIX_BINARY_FORMAT_DENSE -1
157790ace30eSBarry Smith 
15787087cfbeSBarry Smith extern PetscErrorCode  MatMPIBAIJSetHashTableFactor(Mat,PetscReal);
15797087cfbeSBarry Smith extern PetscErrorCode  MatISGetLocalMat(Mat,Mat*);
15801f4e1ec7SBarry Smith 
1581d9274352SBarry Smith /*S
1582d9274352SBarry Smith      MatNullSpace - Object that removes a null space from a vector, i.e.
1583d9274352SBarry Smith          orthogonalizes the vector to a subsapce
1584d9274352SBarry Smith 
1585f7a9e4ceSBarry Smith    Level: advanced
1586d9274352SBarry Smith 
1587d9274352SBarry Smith   Concepts: matrix; linear operator, null space
1588d9274352SBarry Smith 
15896e1639daSBarry Smith   Users manual sections:
15906e1639daSBarry Smith .   sec_singular
15916e1639daSBarry Smith 
1592d9274352SBarry Smith .seealso:  MatNullSpaceCreate()
1593d9274352SBarry Smith S*/
159474637425SBarry Smith typedef struct _p_MatNullSpace* MatNullSpace;
1595d9274352SBarry Smith 
15967087cfbeSBarry Smith extern PetscErrorCode  MatNullSpaceCreate(MPI_Comm,PetscBool ,PetscInt,const Vec[],MatNullSpace*);
15977087cfbeSBarry Smith extern PetscErrorCode  MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(MatNullSpace,Vec,void*),void*);
1598*d34fcf5fSBarry Smith extern PetscErrorCode  MatNullSpaceDestroy(MatNullSpace*);
15997087cfbeSBarry Smith extern PetscErrorCode  MatNullSpaceRemove(MatNullSpace,Vec,Vec*);
16007087cfbeSBarry Smith extern PetscErrorCode  MatNullSpaceAttach(Mat,MatNullSpace);
16017087cfbeSBarry Smith extern PetscErrorCode  MatNullSpaceTest(MatNullSpace,Mat,PetscBool  *);
160274637425SBarry Smith 
16037087cfbeSBarry Smith extern PetscErrorCode  MatReorderingSeqSBAIJ(Mat,IS);
16047087cfbeSBarry Smith extern PetscErrorCode  MatMPISBAIJSetHashTableFactor(Mat,PetscReal);
16057087cfbeSBarry Smith extern PetscErrorCode  MatSeqSBAIJSetColumnIndices(Mat,PetscInt *);
16067087cfbeSBarry Smith extern PetscErrorCode  MatSeqBAIJInvertBlockDiagonal(Mat);
16073f1d51d7SBarry Smith 
16087087cfbeSBarry Smith extern PetscErrorCode  MatCreateMAIJ(Mat,PetscInt,Mat*);
16097087cfbeSBarry Smith extern PetscErrorCode  MatMAIJRedimension(Mat,PetscInt,Mat*);
16107087cfbeSBarry Smith extern PetscErrorCode  MatMAIJGetAIJ(Mat,Mat*);
1611c4f061fbSSatish Balay 
16127087cfbeSBarry Smith extern PetscErrorCode  MatComputeExplicitOperator(Mat,Mat*);
1613b0a32e0cSBarry Smith 
16147087cfbeSBarry Smith extern PetscErrorCode  MatDiagonalScaleLocal(Mat,Vec);
161504f1ad80SBarry Smith 
16167087cfbeSBarry Smith extern PetscErrorCode  MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*);
16177087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetBase(Mat,Vec,Vec);
16187087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*);
16197087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*));
16207087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec));
16217087cfbeSBarry Smith extern PetscErrorCode  MatMFFDAddNullSpace(Mat,MatNullSpace);
16227087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt);
16237087cfbeSBarry Smith extern PetscErrorCode  MatMFFDResetHHistory(Mat);
16247087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetFunctionError(Mat,PetscReal);
16257087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetPeriod(Mat,PetscInt);
16267087cfbeSBarry Smith extern PetscErrorCode  MatMFFDGetH(Mat,PetscScalar *);
16277087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetOptionsPrefix(Mat,const char[]);
16287087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetFromOptions(Mat);
16297087cfbeSBarry Smith extern PetscErrorCode  MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*);
16307087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*);
1631e884886fSBarry Smith 
16326370053bSBarry Smith /*S
16336370053bSBarry Smith     MatMFFD - A data structured used to manage the computation of the h differencing parameter for matrix-free
16346370053bSBarry Smith               Jacobian vector products
1635e884886fSBarry Smith 
16366370053bSBarry Smith     Notes: MATMFFD is a specific MatType which uses the MatMFFD data structure
16376370053bSBarry Smith 
16386370053bSBarry Smith            MatMFFD*() methods actually take the Mat as their first argument. Not a MatMFFD data structure
16396370053bSBarry Smith 
16406370053bSBarry Smith     Level: developer
16416370053bSBarry Smith 
16426370053bSBarry Smith .seealso: MATMFFD, MatCreateMFFD(), MatMFFDSetFuction(), MatMFFDSetType(), MatMFFDRegister()
16436370053bSBarry Smith S*/
1644e884886fSBarry Smith typedef struct _p_MatMFFD* MatMFFD;
1645e884886fSBarry Smith 
1646e884886fSBarry Smith /*E
1647e884886fSBarry Smith     MatMFFDType - algorithm used to compute the h used in computing matrix-vector products via differencing of the function
1648e884886fSBarry Smith 
1649e884886fSBarry Smith    Level: beginner
1650e884886fSBarry Smith 
1651e884886fSBarry Smith .seealso: MatMFFDSetType(), MatMFFDRegister()
1652e884886fSBarry Smith E*/
1653a313700dSBarry Smith #define MatMFFDType char*
1654e884886fSBarry Smith #define MATMFFD_DS  "ds"
1655e884886fSBarry Smith #define MATMFFD_WP  "wp"
1656e884886fSBarry Smith 
16577087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetType(Mat,const MatMFFDType);
16587087cfbeSBarry Smith extern PetscErrorCode  MatMFFDRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatMFFD));
1659e884886fSBarry Smith 
1660e884886fSBarry Smith /*MC
1661e884886fSBarry Smith    MatMFFDRegisterDynamic - Adds a method to the MatMFFD registry.
1662e884886fSBarry Smith 
1663e884886fSBarry Smith    Synopsis:
16641890ba74SBarry Smith    PetscErrorCode MatMFFDRegisterDynamic(const char *name_solver,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatMFFD))
1665e884886fSBarry Smith 
1666e884886fSBarry Smith    Not Collective
1667e884886fSBarry Smith 
1668e884886fSBarry Smith    Input Parameters:
1669e884886fSBarry Smith +  name_solver - name of a new user-defined compute-h module
1670e884886fSBarry Smith .  path - path (either absolute or relative) the library containing this solver
1671e884886fSBarry Smith .  name_create - name of routine to create method context
1672e884886fSBarry Smith -  routine_create - routine to create method context
1673e884886fSBarry Smith 
1674e884886fSBarry Smith    Level: developer
1675e884886fSBarry Smith 
1676e884886fSBarry Smith    Notes:
1677e884886fSBarry Smith    MatMFFDRegisterDynamic() may be called multiple times to add several user-defined solvers.
1678e884886fSBarry Smith 
1679e884886fSBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
1680e884886fSBarry Smith    is ignored.
1681e884886fSBarry Smith 
1682e884886fSBarry Smith    Sample usage:
1683e884886fSBarry Smith .vb
1684e884886fSBarry Smith    MatMFFDRegisterDynamic("my_h",/home/username/my_lib/lib/libO/solaris/mylib.a,
1685e884886fSBarry Smith                "MyHCreate",MyHCreate);
1686e884886fSBarry Smith .ve
1687e884886fSBarry Smith 
1688e884886fSBarry Smith    Then, your solver can be chosen with the procedural interface via
1689e884886fSBarry Smith $     MatMFFDSetType(mfctx,"my_h")
1690e884886fSBarry Smith    or at runtime via the option
1691e884886fSBarry Smith $     -snes_mf_type my_h
1692e884886fSBarry Smith 
1693e884886fSBarry Smith .keywords: MatMFFD, register
1694e884886fSBarry Smith 
1695e884886fSBarry Smith .seealso: MatMFFDRegisterAll(), MatMFFDRegisterDestroy()
1696e884886fSBarry Smith M*/
1697e884886fSBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1698e884886fSBarry Smith #define MatMFFDRegisterDynamic(a,b,c,d) MatMFFDRegister(a,b,c,0)
1699e884886fSBarry Smith #else
1700e884886fSBarry Smith #define MatMFFDRegisterDynamic(a,b,c,d) MatMFFDRegister(a,b,c,d)
1701e884886fSBarry Smith #endif
1702e884886fSBarry Smith 
17037087cfbeSBarry Smith extern PetscErrorCode  MatMFFDRegisterAll(const char[]);
17047087cfbeSBarry Smith extern PetscErrorCode  MatMFFDRegisterDestroy(void);
17057087cfbeSBarry Smith extern PetscErrorCode  MatMFFDDSSetUmin(Mat,PetscReal);
17067087cfbeSBarry Smith extern PetscErrorCode  MatMFFDWPSetComputeNormU(Mat,PetscBool );
1707e884886fSBarry Smith 
1708e884886fSBarry Smith 
17097087cfbeSBarry Smith extern PetscErrorCode  PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *);
17107087cfbeSBarry Smith extern PetscErrorCode  PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *);
17117dbadf16SMatthew Knepley 
171297969023SHong Zhang /*
171397969023SHong Zhang    PETSc interface to MUMPS
171497969023SHong Zhang */
171597969023SHong Zhang #ifdef PETSC_HAVE_MUMPS
1716f250808bSBarry Smith extern PetscErrorCode  MatMumpsSetIcntl(Mat,PetscInt,PetscInt);
171797969023SHong Zhang #endif
171823a5497aSJed Brown 
1719d954db57SHong Zhang /*
1720d954db57SHong Zhang    PETSc interface to SUPERLU
1721d954db57SHong Zhang */
1722d954db57SHong Zhang #ifdef PETSC_HAVE_SUPERLU
1723f250808bSBarry Smith extern PetscErrorCode  MatSuperluSetILUDropTol(Mat,PetscReal);
1724d954db57SHong Zhang #endif
1725d954db57SHong Zhang 
1726c8883902SJed Brown extern PetscErrorCode MatCreateNest(MPI_Comm,PetscInt,const IS[],PetscInt,const IS[],const Mat[],Mat*);
1727c8883902SJed Brown extern PetscErrorCode MatNestGetSize(Mat,PetscInt*,PetscInt*);
1728c8883902SJed Brown extern PetscErrorCode MatNestGetSubMats(Mat,PetscInt*,PetscInt*,Mat***);
1729c8883902SJed Brown extern PetscErrorCode MatNestGetSubMat(Mat,PetscInt,PetscInt,Mat*);
17307087cfbeSBarry Smith extern PetscErrorCode MatNestSetVecType(Mat,const VecType);
1731c8883902SJed Brown extern PetscErrorCode MatNestSetSubMats(Mat,PetscInt,const IS[],PetscInt,const IS[],const Mat[]);
1732d8588912SDave May 
173323a5497aSJed Brown PETSC_EXTERN_CXX_END
173423a5497aSJed Brown #endif
1735