xref: /petsc/include/petscmat.h (revision cd6b891e1504900ecbee365dbcdfc13910f27c49)
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"
415a11e1b2SBarry Smith #define MATAIJCUDA             "aijcuda"
425a11e1b2SBarry Smith #define MATSEQAIJCUDA            "seqaijcuda"
435a11e1b2SBarry Smith #define MATMPIAIJCUDA            "mpiaijcuda"
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"
655a7f1df3SHong Zhang #define MATSEQFFTW         "seqfftw"
66e133240eSMatthew G Knepley #define MATSEQCUFFT        "seqcufft"
67557cca28SSatish Balay #define MATTRANSPOSEMAT    "transpose"
6872ca8751SBarry Smith #define MATSCHURCOMPLEMENT "schurcomplement"
691d6018f0SLisandro Dalcin #define MATPYTHON          "python"
70f91d8e95SBarry Smith #define MATHYPRESTRUCT     "hyprestruct"
71a9e6138eSGlenn Hammond #define MATHYPRESSTRUCT    "hypresstruct"
72ee1cef2cSJed Brown #define MATSUBMATRIX       "submatrix"
73dbc6227fSDmitry Karpeev #define MATDD              "matdd"
74ba2f8784SDmitry Karpeev #define MATIM              "matim"
752c0dbf93SJed Brown #define MATLOCALREF        "localref"
76d8588912SDave May #define MATNEST            "nest"
77773941d6SBarry Smith 
789989ab13SBarry Smith /*E
79c7393fdbSBarry Smith     MatSolverPackage - String with the name of a PETSc matrix solver type.
809989ab13SBarry Smith 
819989ab13SBarry Smith     For example: "petsc" indicates what PETSc provides, "superlu" indicates either
829989ab13SBarry Smith        SuperLU or SuperLU_Dist etc.
839989ab13SBarry Smith 
849989ab13SBarry Smith 
859989ab13SBarry Smith    Level: beginner
869989ab13SBarry Smith 
875c9eb25fSBarry Smith .seealso: MatGetFactor(), Mat, MatSetType(), MatType
889989ab13SBarry Smith E*/
89c7393fdbSBarry Smith #define MatSolverPackage char*
902692d6eeSBarry Smith #define MATSOLVERSPOOLES      "spooles"
912692d6eeSBarry Smith #define MATSOLVERSUPERLU      "superlu"
922692d6eeSBarry Smith #define MATSOLVERSUPERLU_DIST "superlu_dist"
932692d6eeSBarry Smith #define MATSOLVERUMFPACK      "umfpack"
9420db9a53SJed Brown #define MATSOLVERCHOLMOD      "cholmod"
952692d6eeSBarry Smith #define MATSOLVERESSL         "essl"
962692d6eeSBarry Smith #define MATSOLVERLUSOL        "lusol"
972692d6eeSBarry Smith #define MATSOLVERMUMPS        "mumps"
982692d6eeSBarry Smith #define MATSOLVERPASTIX       "pastix"
992692d6eeSBarry Smith #define MATSOLVERDSCPACK      "dscpack"
1002692d6eeSBarry Smith #define MATSOLVERMATLAB       "matlab"
1012692d6eeSBarry Smith #define MATSOLVERPETSC        "petsc"
1022692d6eeSBarry Smith #define MATSOLVERPLAPACK      "plapack"
1032692d6eeSBarry Smith #define MATSOLVERBAS          "bas"
104773941d6SBarry Smith 
105b24902e0SBarry Smith /*E
106b24902e0SBarry Smith     MatFactorType - indicates what type of factorization is requested
107b24902e0SBarry Smith 
108b24902e0SBarry Smith     Level: beginner
109b24902e0SBarry Smith 
110b24902e0SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
111b24902e0SBarry Smith 
112c7393fdbSBarry Smith .seealso: MatSolverPackage, MatGetFactor()
113b24902e0SBarry Smith E*/
114599ef60dSHong Zhang typedef enum {MAT_FACTOR_NONE, MAT_FACTOR_LU, MAT_FACTOR_CHOLESKY, MAT_FACTOR_ILU, MAT_FACTOR_ICC,MAT_FACTOR_ILUDT} MatFactorType;
11534918c53SJed Brown extern const char *const MatFactorTypes[];
116e92e720dSBarry Smith 
1177087cfbeSBarry Smith extern PetscErrorCode  MatGetFactor(Mat,const MatSolverPackage,MatFactorType,Mat*);
1187087cfbeSBarry Smith extern PetscErrorCode  MatGetFactorAvailable(Mat,const MatSolverPackage,MatFactorType,PetscBool *);
1197087cfbeSBarry Smith extern PetscErrorCode  MatFactorGetSolverPackage(Mat,const MatSolverPackage*);
1207087cfbeSBarry Smith extern PetscErrorCode  MatGetFactorType(Mat,MatFactorType*);
1219989ab13SBarry Smith 
122c06d978dSMatthew Knepley /* Logging support */
1230700a824SBarry Smith #define    MAT_FILE_CLASSID 1211216    /* used to indicate matrices in binary files */
1247087cfbeSBarry Smith extern PetscClassId  MAT_CLASSID;
1257087cfbeSBarry Smith extern PetscClassId  MAT_FDCOLORING_CLASSID;
1267087cfbeSBarry Smith extern PetscClassId  MAT_PARTITIONING_CLASSID;
1277087cfbeSBarry Smith extern PetscClassId  MAT_NULLSPACE_CLASSID;
1287087cfbeSBarry Smith extern PetscClassId  MATMFFD_CLASSID;
129c06d978dSMatthew Knepley 
130ceb03754SKris Buschelman /*E
131ceb03754SKris Buschelman     MatReuse - Indicates if matrices obtained from a previous call to MatGetSubMatrices()
132d6eff37eSBarry Smith      or MatGetSubMatrix() are to be reused to store the new matrix values. For MatConvert() is used to indicate
133d6eff37eSBarry Smith      that the input matrix is to be replaced with the converted matrix.
134ceb03754SKris Buschelman 
135ceb03754SKris Buschelman     Level: beginner
136ceb03754SKris Buschelman 
137ceb03754SKris Buschelman    Any additions/changes here MUST also be made in include/finclude/petscmat.h
138ceb03754SKris Buschelman 
1390c451bc4SBarry Smith .seealso: MatGetSubMatrices(), MatGetSubMatrix(), MatDestroyMatrices(), MatConvert()
140ceb03754SKris Buschelman E*/
141dfe085dbSJed Brown typedef enum {MAT_INITIAL_MATRIX,MAT_REUSE_MATRIX,MAT_IGNORE_MATRIX} MatReuse;
142ceb03754SKris Buschelman 
1435494a064SHong Zhang /*E
1445494a064SHong Zhang     MatGetSubMatrixOption - Indicates if matrices obtained from a call to MatGetSubMatrices()
145829201f2SHong Zhang      include the matrix values. Currently it is only used by MatGetSeqNonzerostructure().
1465494a064SHong Zhang 
1475494a064SHong Zhang     Level: beginner
1485494a064SHong Zhang 
149829201f2SHong Zhang .seealso: MatGetSeqNonzerostructure()
1505494a064SHong Zhang E*/
1515494a064SHong Zhang typedef enum {MAT_DO_NOT_GET_VALUES,MAT_GET_VALUES} MatGetSubMatrixOption;
1525494a064SHong Zhang 
1537087cfbeSBarry Smith extern PetscErrorCode  MatInitializePackage(const char[]);
154c06d978dSMatthew Knepley 
1557087cfbeSBarry Smith extern PetscErrorCode  MatCreate(MPI_Comm,Mat*);
156a69a960dSMatthew Knepley PetscPolymorphicFunction(MatCreate,(MPI_Comm comm),(comm,&A),Mat,A)
157a69a960dSMatthew Knepley PetscPolymorphicFunction(MatCreate,(),(PETSC_COMM_WORLD,&A),Mat,A)
1587087cfbeSBarry Smith extern PetscErrorCode  MatSetSizes(Mat,PetscInt,PetscInt,PetscInt,PetscInt);
1597087cfbeSBarry Smith extern PetscErrorCode  MatSetType(Mat,const MatType);
1607087cfbeSBarry Smith extern PetscErrorCode  MatSetFromOptions(Mat);
1617087cfbeSBarry Smith extern PetscErrorCode  MatSetUpPreallocation(Mat);
1627087cfbeSBarry Smith extern PetscErrorCode  MatRegisterAll(const char[]);
1637087cfbeSBarry Smith extern PetscErrorCode  MatRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat));
16430de9b25SBarry Smith 
1657087cfbeSBarry Smith extern PetscErrorCode  MatSetOptionsPrefix(Mat,const char[]);
1667087cfbeSBarry Smith extern PetscErrorCode  MatAppendOptionsPrefix(Mat,const char[]);
1677087cfbeSBarry Smith extern PetscErrorCode  MatGetOptionsPrefix(Mat,const char*[]);
168f69a0ea3SMatthew Knepley 
16930de9b25SBarry Smith /*MC
17030de9b25SBarry Smith    MatRegisterDynamic - Adds a new matrix type
17130de9b25SBarry Smith 
17230de9b25SBarry Smith    Synopsis:
1731890ba74SBarry Smith    PetscErrorCode MatRegisterDynamic(const char *name,const char *path,const char *name_create,PetscErrorCode (*routine_create)(Mat))
17430de9b25SBarry Smith 
17530de9b25SBarry Smith    Not Collective
17630de9b25SBarry Smith 
17730de9b25SBarry Smith    Input Parameters:
17830de9b25SBarry Smith +  name - name of a new user-defined matrix type
17930de9b25SBarry Smith .  path - path (either absolute or relative) the library containing this solver
18030de9b25SBarry Smith .  name_create - name of routine to create method context
18130de9b25SBarry Smith -  routine_create - routine to create method context
18230de9b25SBarry Smith 
18330de9b25SBarry Smith    Notes:
18430de9b25SBarry Smith    MatRegisterDynamic() may be called multiple times to add several user-defined solvers.
18530de9b25SBarry Smith 
18630de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
18730de9b25SBarry Smith    is ignored.
18830de9b25SBarry Smith 
18930de9b25SBarry Smith    Sample usage:
19030de9b25SBarry Smith .vb
19130de9b25SBarry Smith    MatRegisterDynamic("my_mat",/home/username/my_lib/lib/libO/solaris/mylib.a,
19230de9b25SBarry Smith                "MyMatCreate",MyMatCreate);
19330de9b25SBarry Smith .ve
19430de9b25SBarry Smith 
19530de9b25SBarry Smith    Then, your solver can be chosen with the procedural interface via
19630de9b25SBarry Smith $     MatSetType(Mat,"my_mat")
19730de9b25SBarry Smith    or at runtime via the option
19830de9b25SBarry Smith $     -mat_type my_mat
19930de9b25SBarry Smith 
20030de9b25SBarry Smith    Level: advanced
20130de9b25SBarry Smith 
202ab901514SBarry Smith    Notes: ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values.
20330de9b25SBarry Smith          If your function is not being put into a shared library then use VecRegister() instead
20430de9b25SBarry Smith 
20530de9b25SBarry Smith .keywords: Mat, register
20630de9b25SBarry Smith 
20730de9b25SBarry Smith .seealso: MatRegisterAll(), MatRegisterDestroy()
20830de9b25SBarry Smith 
20930de9b25SBarry Smith M*/
210273d9f13SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
211273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,0)
212273d9f13SBarry Smith #else
213273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,d)
21430de9b25SBarry Smith #endif
21530de9b25SBarry Smith 
216ace3abfcSBarry Smith extern PetscBool  MatRegisterAllCalled;
217b0a32e0cSBarry Smith extern PetscFList MatList;
218b022a5c1SBarry Smith extern PetscFList MatColoringList;
219b022a5c1SBarry Smith extern PetscFList MatPartitioningList;
22028988994SBarry Smith 
2213b224e63SBarry Smith /*E
2223b224e63SBarry Smith     MatStructure - Indicates if the matrix has the same nonzero structure
2233b224e63SBarry Smith 
2243b224e63SBarry Smith     Level: beginner
2253b224e63SBarry Smith 
2263b224e63SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
2273b224e63SBarry Smith 
2283b224e63SBarry Smith .seealso: MatCopy(), KSPSetOperators(), PCSetOperators()
2293b224e63SBarry Smith E*/
2303b224e63SBarry Smith typedef enum {SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER,SUBSET_NONZERO_PATTERN} MatStructure;
2313b224e63SBarry Smith 
2327087cfbeSBarry Smith extern PetscErrorCode  MatCreateSeqDense(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*);
2337087cfbeSBarry Smith extern PetscErrorCode  MatCreateMPIDense(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*);
2347087cfbeSBarry Smith extern PetscErrorCode  MatCreateSeqAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
235ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,nz,nnz,&A),Mat,A)
236ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,&A),Mat,A)
237ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,0,nnz,&A),Mat,A)
238ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,&A),Mat,A)
239ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,A))
240ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,m,n,0,nnz,A))
241ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,A))
2427087cfbeSBarry Smith extern PetscErrorCode  MatCreateMPIAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
243ba966639SSatish 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)
244ba966639SSatish 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)
245ba966639SSatish 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)
246ba966639SSatish 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)
247ba966639SSatish 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))
248ba966639SSatish 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))
249ba966639SSatish 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))
250ba966639SSatish 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)
251ba966639SSatish 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)
252ba966639SSatish 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)
253ba966639SSatish 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)
254ba966639SSatish 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))
255ba966639SSatish 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))
256ba966639SSatish 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))
2577087cfbeSBarry Smith extern PetscErrorCode  MatCreateMPIAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *);
2587087cfbeSBarry Smith extern PetscErrorCode  MatCreateMPIAIJWithSplitArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],PetscInt[],PetscInt[],PetscScalar[],Mat*);
2598d7a6e47SBarry Smith 
2607087cfbeSBarry Smith extern PetscErrorCode  MatCreateSeqBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
261ba966639SSatish 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)
262ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A)
263ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A)
264ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A)
265ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A))
266ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A))
267ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A))
2687087cfbeSBarry Smith extern PetscErrorCode  MatCreateMPIBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
269ba966639SSatish 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)
270ba966639SSatish 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)
271ba966639SSatish 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)
272ba966639SSatish 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)
273ba966639SSatish 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))
274ba966639SSatish 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))
275ba966639SSatish 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))
276ba966639SSatish 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)
277ba966639SSatish 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)
278ba966639SSatish 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)
279ba966639SSatish 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)
280ba966639SSatish 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))
281ba966639SSatish 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))
282ba966639SSatish 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))
2837087cfbeSBarry Smith extern PetscErrorCode  MatCreateMPIBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat*);
284d21a29f3SJed Brown 
2857087cfbeSBarry Smith extern PetscErrorCode  MatCreateMPIAdj(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscInt[],Mat*);
2867087cfbeSBarry Smith extern PetscErrorCode  MatCreateSeqSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
287ba966639SSatish 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)
288ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A)
289ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A)
290ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A)
291ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A))
292ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A))
293ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A))
294d21a29f3SJed Brown 
2957087cfbeSBarry Smith extern PetscErrorCode  MatCreateMPISBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
296ba966639SSatish 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)
297ba966639SSatish 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)
298ba966639SSatish 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)
299ba966639SSatish 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)
300ba966639SSatish 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))
301ba966639SSatish 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))
302ba966639SSatish 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))
303ba966639SSatish 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)
304ba966639SSatish 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)
305ba966639SSatish 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)
306ba966639SSatish 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)
307ba966639SSatish 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))
308ba966639SSatish 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))
309ba966639SSatish 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))
3107087cfbeSBarry Smith extern PetscErrorCode  MatCreateMPISBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *);
3117087cfbeSBarry Smith extern PetscErrorCode  MatMPISBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
312dfb205c3SBarry Smith 
3137087cfbeSBarry Smith extern PetscErrorCode  MatCreateShell(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,void *,Mat*);
314ba966639SSatish 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)
315ba966639SSatish Balay PetscPolymorphicFunction(MatCreateShell,(PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx),(PETSC_COMM_WORLD,m,n,M,N,ctx,&A),Mat,A)
3167087cfbeSBarry Smith extern PetscErrorCode  MatCreateAdic(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,void (*)(void),Mat*);
3177087cfbeSBarry Smith extern PetscErrorCode  MatCreateNormal(Mat,Mat*);
318ba966639SSatish Balay PetscPolymorphicFunction(MatCreateNormal,(Mat mat),(mat,&A),Mat,A)
3197087cfbeSBarry Smith extern PetscErrorCode  MatCreateLRC(Mat,Mat,Mat,Mat*);
3207087cfbeSBarry Smith extern PetscErrorCode  MatCreateIS(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,ISLocalToGlobalMapping,Mat*);
3217087cfbeSBarry Smith extern PetscErrorCode  MatCreateSeqAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
3227087cfbeSBarry Smith extern PetscErrorCode  MatCreateMPIAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
3237087cfbeSBarry Smith extern PetscErrorCode  MatCreateScatter(MPI_Comm,VecScatter,Mat*);
3247087cfbeSBarry Smith extern PetscErrorCode  MatScatterSetVecScatter(Mat,VecScatter);
3257087cfbeSBarry Smith extern PetscErrorCode  MatScatterGetVecScatter(Mat,VecScatter*);
3267087cfbeSBarry Smith extern PetscErrorCode  MatCreateBlockMat(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,Mat*);
3277087cfbeSBarry Smith extern PetscErrorCode  MatCompositeAddMat(Mat,Mat);
3287087cfbeSBarry Smith extern PetscErrorCode  MatCompositeMerge(Mat);
3297087cfbeSBarry Smith extern PetscErrorCode  MatCreateComposite(MPI_Comm,PetscInt,const Mat*,Mat*);
3306d7c1e57SBarry Smith typedef enum {MAT_COMPOSITE_ADDITIVE,MAT_COMPOSITE_MULTIPLICATIVE} MatCompositeType;
3317087cfbeSBarry Smith extern PetscErrorCode  MatCompositeSetType(Mat,MatCompositeType);
3326d7c1e57SBarry Smith 
33345b63f25SDmitry Karpeev #if defined PETSC_HAVE_MATDD
3347087cfbeSBarry Smith extern PetscErrorCode  MatDDCreate(Mat A);
335d214e108SDmitry Karpeev typedef enum {MATDD_DOMAINS_COLUMN, MATDD_DOMAINS_ROW} MatDDDomainType;
3367087cfbeSBarry Smith extern PetscErrorCode  MatDDSetDomainsLocal(Mat A, MatDDDomainType type, PetscInt domain_count, const PetscInt *supported_domains, const PetscInt *domain_limits, PetscBool covering);
3377087cfbeSBarry Smith extern PetscErrorCode  MatDDSetDomainsLocalIS(Mat A, MatDDDomainType type, IS supported_domains, IS domain_limits, PetscBool covering);
3387087cfbeSBarry Smith extern PetscErrorCode  MatDDSetScatter(Mat A, Mat S);
3397087cfbeSBarry Smith extern PetscErrorCode  MatDDSetGather(Mat A,  Mat G);
340d214e108SDmitry Karpeev /**/
341dbc6227fSDmitry Karpeev typedef enum {MATDD_BLOCK_COMM_DEFAULT = 0, MATDD_BLOCK_COMM_SELF = -1, MATDD_BLOCK_COMM_DETERMINE = -2} MatDDBlockCommType;
3427087cfbeSBarry Smith extern PetscErrorCode  MatDDGetDefaltBlockType(Mat A, const MatType *type);
3437087cfbeSBarry Smith extern PetscErrorCode  MatDDSetDefaltBlockType(Mat A, const MatType type);
3447087cfbeSBarry Smith extern PetscErrorCode  MatDDAddBlockLocal(Mat A, PetscInt rowblock, PetscInt colblock, const MatType blockmattype,  MatDDBlockCommType blockcommtype, Mat *block);
3457087cfbeSBarry Smith extern PetscErrorCode  MatDDSetBlockLocal(Mat A, PetscInt rowblock, PetscInt colblock, Mat block);
3467087cfbeSBarry Smith extern PetscErrorCode  MatDDGetBlockLocal(Mat A, PetscInt rowblock, PetscInt colblock, Mat *block);
347d214e108SDmitry Karpeev /**/
3487087cfbeSBarry Smith extern PetscErrorCode  MatDDAIJSetPreallocation(Mat A,PetscInt nz,PetscInt *nnz);
34945b63f25SDmitry Karpeev #endif
35062ee7b53SDmitry Karpeev 
351ba2f8784SDmitry Karpeev #if defined PETSC_HAVE_MATIM
3527087cfbeSBarry Smith extern PetscErrorCode  MatIMSetIS(Mat A, IS in, IS out);
3537087cfbeSBarry Smith extern PetscErrorCode  MatIMGetIS(Mat A, IS *in, IS *out);
354ba2f8784SDmitry Karpeev #endif
355ba2f8784SDmitry Karpeev 
356ba2f8784SDmitry Karpeev 
357a6a5cd3fSDmitry Karpeev 
3587087cfbeSBarry Smith extern PetscErrorCode  MatCreateSeqFFTW(MPI_Comm,PetscInt,const PetscInt[],Mat*);
3597087cfbeSBarry Smith extern PetscErrorCode  MatCreateSeqCUFFT(MPI_Comm,PetscInt,const PetscInt[],Mat*);
3607087cfbeSBarry Smith extern PetscErrorCode  MatCreateTranspose(Mat,Mat*);
3617087cfbeSBarry Smith extern PetscErrorCode  MatCreateSubMatrix(Mat,IS,IS,Mat*);
3627087cfbeSBarry Smith extern PetscErrorCode  MatSubMatrixUpdate(Mat,Mat,IS,IS);
3637087cfbeSBarry Smith extern PetscErrorCode  MatCreateLocalRef(Mat,IS,IS,Mat*);
3641472f72bSBarry Smith 
3657087cfbeSBarry Smith extern PetscErrorCode  MatCreatePython(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const char[],Mat*);
3667087cfbeSBarry Smith extern PetscErrorCode  MatPythonSetType(Mat,const char[]);
3671d6018f0SLisandro Dalcin 
3681d6018f0SLisandro Dalcin 
3697087cfbeSBarry Smith extern PetscErrorCode  MatSetUp(Mat);
3707087cfbeSBarry Smith extern PetscErrorCode  MatDestroy(Mat);
37121c89e3eSBarry Smith 
3727087cfbeSBarry Smith extern PetscErrorCode  MatConjugate(Mat);
3737087cfbeSBarry Smith extern PetscErrorCode  MatRealPart(Mat);
3747087cfbeSBarry Smith extern PetscErrorCode  MatImaginaryPart(Mat);
3757087cfbeSBarry Smith extern PetscErrorCode  MatGetDiagonalBlock(Mat,PetscBool *,MatReuse,Mat*);
3767087cfbeSBarry Smith extern PetscErrorCode  MatGetTrace(Mat,PetscScalar*);
37799cafbc1SBarry Smith 
3788ed539a5SBarry Smith /* ------------------------------------------------------------*/
3797087cfbeSBarry Smith extern PetscErrorCode  MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
3807087cfbeSBarry Smith extern PetscErrorCode  MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
3817087cfbeSBarry Smith extern PetscErrorCode  MatSetValuesRow(Mat,PetscInt,const PetscScalar[]);
3827087cfbeSBarry Smith extern PetscErrorCode  MatSetValuesRowLocal(Mat,PetscInt,const PetscScalar[]);
38384cb2905SBarry Smith 
3842ef4de8bSBarry Smith /*S
3852ef4de8bSBarry Smith      MatStencil - Data structure (C struct) for storing information about a single row or
3862ef4de8bSBarry Smith         column of a matrix as index on an associated grid.
3872ef4de8bSBarry Smith 
3882ef4de8bSBarry Smith    Level: beginner
3892ef4de8bSBarry Smith 
3902ef4de8bSBarry Smith   Concepts: matrix; linear operator
3912ef4de8bSBarry Smith 
392d7bd93baSBarry Smith .seealso:  MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockStencil()
3932ef4de8bSBarry Smith S*/
394435da068SBarry Smith typedef struct {
395c1ac3661SBarry Smith   PetscInt k,j,i,c;
396435da068SBarry Smith } MatStencil;
3972ef4de8bSBarry Smith 
3987087cfbeSBarry Smith extern PetscErrorCode  MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
3997087cfbeSBarry Smith extern PetscErrorCode  MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
4007087cfbeSBarry Smith extern PetscErrorCode  MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt);
401435da068SBarry Smith 
4027087cfbeSBarry Smith extern PetscErrorCode  MatSetColoring(Mat,ISColoring);
4037087cfbeSBarry Smith extern PetscErrorCode  MatSetValuesAdic(Mat,void*);
4047087cfbeSBarry Smith extern PetscErrorCode  MatSetValuesAdifor(Mat,PetscInt,void*);
4053a7fca6bSBarry Smith 
406d91e6319SBarry Smith /*E
407d91e6319SBarry Smith     MatAssemblyType - Indicates if the matrix is now to be used, or if you plan
408d91e6319SBarry Smith      to continue to add values to it
409d91e6319SBarry Smith 
410d91e6319SBarry Smith     Level: beginner
411d91e6319SBarry Smith 
412d91e6319SBarry Smith .seealso: MatAssemblyBegin(), MatAssemblyEnd()
413d91e6319SBarry Smith E*/
4146d4a8577SBarry Smith typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType;
4157087cfbeSBarry Smith extern PetscErrorCode  MatAssemblyBegin(Mat,MatAssemblyType);
4167087cfbeSBarry Smith extern PetscErrorCode  MatAssemblyEnd(Mat,MatAssemblyType);
4177087cfbeSBarry Smith extern PetscErrorCode  MatAssembled(Mat,PetscBool *);
4184f9c727eSBarry Smith 
4191d73ed98SBarry Smith 
42030de9b25SBarry Smith 
421d91e6319SBarry Smith /*E
422d91e6319SBarry Smith     MatOption - Options that may be set for a matrix and its behavior or storage
423d91e6319SBarry Smith 
424d91e6319SBarry Smith     Level: beginner
425d91e6319SBarry Smith 
4260a835dfdSSatish Balay    Any additions/changes here MUST also be made in include/finclude/petscmat.h
427d91e6319SBarry Smith 
428d91e6319SBarry Smith .seealso: MatSetOption()
429d91e6319SBarry Smith E*/
430512a5fc5SBarry Smith typedef enum {MAT_ROW_ORIENTED,MAT_NEW_NONZERO_LOCATIONS,
4314e0d8c25SBarry Smith               MAT_SYMMETRIC,
4324e0d8c25SBarry Smith               MAT_STRUCTURALLY_SYMMETRIC,
4338047e77bSBarry Smith               MAT_NEW_DIAGONALS,MAT_IGNORE_OFF_PROC_ENTRIES,
4344e0d8c25SBarry Smith               MAT_NEW_NONZERO_LOCATION_ERR,
4354e0d8c25SBarry Smith               MAT_NEW_NONZERO_ALLOCATION_ERR,MAT_USE_HASH_TABLE,
436a9817697SBarry Smith               MAT_KEEP_NONZERO_PATTERN,MAT_IGNORE_ZERO_ENTRIES,
4374e0d8c25SBarry Smith               MAT_USE_INODES,
4384e0d8c25SBarry Smith               MAT_HERMITIAN,
4394e0d8c25SBarry Smith               MAT_SYMMETRY_ETERNAL,
440*cd6b891eSBarry Smith               MAT_CHECK_COMPRESSED_ROW,
4414e0d8c25SBarry Smith               MAT_IGNORE_LOWER_TRIANGULAR,MAT_ERROR_LOWER_TRIANGULAR,
44228b2fa4aSMatthew Knepley               MAT_GETROW_UPPERTRIANGULAR,MAT_UNUSED_NONZERO_LOCATION_ERR,
4434cb17eb5SBarry Smith               MAT_SPD,MAT_NO_OFF_PROC_ENTRIES,MAT_NO_OFF_PROC_ZERO_ROWS,
44428b2fa4aSMatthew Knepley               NUM_MAT_OPTIONS} MatOption;
445290bbb0aSBarry Smith extern const char *MatOptions[];
4467087cfbeSBarry Smith extern PetscErrorCode  MatSetOption(Mat,MatOption,PetscBool );
4477087cfbeSBarry Smith extern PetscErrorCode  MatGetType(Mat,const MatType*);
448a313700dSBarry Smith PetscPolymorphicFunction(MatGetType,(Mat mat),(mat,&t),const MatType,t)
44984cb2905SBarry Smith 
4507087cfbeSBarry Smith extern PetscErrorCode  MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]);
4517087cfbeSBarry Smith extern PetscErrorCode  MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
4527087cfbeSBarry Smith extern PetscErrorCode  MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
4537087cfbeSBarry Smith extern PetscErrorCode  MatGetRowUpperTriangular(Mat);
4547087cfbeSBarry Smith extern PetscErrorCode  MatRestoreRowUpperTriangular(Mat);
4557087cfbeSBarry Smith extern PetscErrorCode  MatGetColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
4567087cfbeSBarry Smith extern PetscErrorCode  MatRestoreColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
4577087cfbeSBarry Smith extern PetscErrorCode  MatGetColumnVector(Mat,Vec,PetscInt);
4587087cfbeSBarry Smith extern PetscErrorCode  MatGetArray(Mat,PetscScalar *[]);
459ba966639SSatish Balay PetscPolymorphicFunction(MatGetArray,(Mat mat),(mat,&a),PetscScalar*,a)
4607087cfbeSBarry Smith extern PetscErrorCode  MatRestoreArray(Mat,PetscScalar *[]);
4617087cfbeSBarry Smith extern PetscErrorCode  MatGetBlockSize(Mat,PetscInt *);
462ba966639SSatish Balay PetscPolymorphicFunction(MatGetBlockSize,(Mat mat),(mat,&a),PetscInt,a)
4637087cfbeSBarry Smith extern PetscErrorCode  MatSetBlockSize(Mat,PetscInt);
4647b80b807SBarry Smith 
4651620fd73SBarry Smith 
4667087cfbeSBarry Smith extern PetscErrorCode  MatMult(Mat,Vec,Vec);
4677087cfbeSBarry Smith extern PetscErrorCode  MatMultDiagonalBlock(Mat,Vec,Vec);
4687087cfbeSBarry Smith extern PetscErrorCode  MatMultAdd(Mat,Vec,Vec,Vec);
469ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
4707087cfbeSBarry Smith extern PetscErrorCode  MatMultTranspose(Mat,Vec,Vec);
4717087cfbeSBarry Smith extern PetscErrorCode  MatMultHermitianTranspose(Mat,Vec,Vec);
4727087cfbeSBarry Smith extern PetscErrorCode  MatIsTranspose(Mat,Mat,PetscReal,PetscBool *);
473ace3abfcSBarry Smith PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B,PetscReal tol),(A,B,tol,&t),PetscBool ,t)
474ace3abfcSBarry Smith PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B),(A,B,0,&t),PetscBool ,t)
4757087cfbeSBarry Smith extern PetscErrorCode  MatIsHermitianTranspose(Mat,Mat,PetscReal,PetscBool *);
4767087cfbeSBarry Smith extern PetscErrorCode  MatMultTransposeAdd(Mat,Vec,Vec,Vec);
4777087cfbeSBarry Smith extern PetscErrorCode  MatMultHermitianTransposeAdd(Mat,Vec,Vec,Vec);
478ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultTransposeAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
4797087cfbeSBarry Smith extern PetscErrorCode  MatMultConstrained(Mat,Vec,Vec);
4807087cfbeSBarry Smith extern PetscErrorCode  MatMultTransposeConstrained(Mat,Vec,Vec);
4817087cfbeSBarry Smith extern PetscErrorCode  MatMatSolve(Mat,Mat,Mat);
482f5edf698SHong Zhang 
483d91e6319SBarry Smith /*E
484d91e6319SBarry Smith     MatDuplicateOption - Indicates if a duplicated sparse matrix should have
485d91e6319SBarry Smith   its numerical values copied over or just its nonzero structure.
486d91e6319SBarry Smith 
487d91e6319SBarry Smith     Level: beginner
488d91e6319SBarry Smith 
489d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
490d91e6319SBarry Smith 
49170dcbbb9SBarry Smith $   MAT_SHARE_NONZERO_PATTERN - the i and j arrays in the new matrix will be shared with the original matrix
49270dcbbb9SBarry Smith $                               this also triggers the MAT_DO_NOT_COPY_VALUES option. This is used when you
49370dcbbb9SBarry Smith $                               have several matrices with the same nonzero pattern.
49470dcbbb9SBarry Smith 
495d91e6319SBarry Smith .seealso: MatDuplicate()
496d91e6319SBarry Smith E*/
49770dcbbb9SBarry Smith typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES,MAT_SHARE_NONZERO_PATTERN} MatDuplicateOption;
4982e8a6d31SBarry Smith 
4997087cfbeSBarry Smith extern PetscErrorCode  MatConvert(Mat,const MatType,MatReuse,Mat*);
500a313700dSBarry Smith PetscPolymorphicFunction(MatConvert,(Mat A,const MatType t),(A,t,MAT_INITIAL_MATRIX,&a),Mat,a)
5017087cfbeSBarry Smith extern PetscErrorCode  MatDuplicate(Mat,MatDuplicateOption,Mat*);
502ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A,MatDuplicateOption o),(A,o,&a),Mat,a)
503ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A),(A,MAT_COPY_VALUES,&a),Mat,a)
50494a9d846SBarry Smith 
505cb5b572fSBarry Smith 
5067087cfbeSBarry Smith extern PetscErrorCode  MatCopy(Mat,Mat,MatStructure);
5077087cfbeSBarry Smith extern PetscErrorCode  MatView(Mat,PetscViewer);
5087087cfbeSBarry Smith extern PetscErrorCode  MatIsSymmetric(Mat,PetscReal,PetscBool *);
509ace3abfcSBarry Smith PetscPolymorphicFunction(MatIsSymmetric,(Mat A,PetscReal tol),(A,tol,&t),PetscBool ,t)
510ace3abfcSBarry Smith PetscPolymorphicFunction(MatIsSymmetric,(Mat A),(A,0,&t),PetscBool ,t)
5117087cfbeSBarry Smith extern PetscErrorCode  MatIsStructurallySymmetric(Mat,PetscBool *);
512ace3abfcSBarry Smith PetscPolymorphicFunction(MatIsStructurallySymmetric,(Mat A),(A,&t),PetscBool ,t)
5137087cfbeSBarry Smith extern PetscErrorCode  MatIsHermitian(Mat,PetscReal,PetscBool *);
514ace3abfcSBarry Smith PetscPolymorphicFunction(MatIsHermitian,(Mat A),(A,0,&t),PetscBool ,t)
5157087cfbeSBarry Smith extern PetscErrorCode  MatIsSymmetricKnown(Mat,PetscBool *,PetscBool *);
5167087cfbeSBarry Smith extern PetscErrorCode  MatIsHermitianKnown(Mat,PetscBool *,PetscBool *);
5177087cfbeSBarry Smith extern PetscErrorCode  MatMissingDiagonal(Mat,PetscBool  *,PetscInt *);
5187087cfbeSBarry Smith extern PetscErrorCode  MatLoad(Mat, PetscViewer);
5197b80b807SBarry Smith 
5207087cfbeSBarry Smith extern PetscErrorCode  MatGetRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,PetscInt *[],PetscInt *[],PetscBool  *);
5217087cfbeSBarry Smith extern PetscErrorCode  MatRestoreRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,PetscInt *[],PetscInt *[],PetscBool  *);
5227087cfbeSBarry Smith extern PetscErrorCode  MatGetColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,PetscInt *[],PetscInt *[],PetscBool  *);
5237087cfbeSBarry Smith extern PetscErrorCode  MatRestoreColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,PetscInt *[],PetscInt *[],PetscBool  *);
524d4fbbf0eSBarry Smith 
525d91e6319SBarry Smith /*S
526d91e6319SBarry Smith      MatInfo - Context of matrix information, used with MatGetInfo()
527d91e6319SBarry Smith 
528d91e6319SBarry Smith    In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE
529d91e6319SBarry Smith 
530d91e6319SBarry Smith    Level: intermediate
531d91e6319SBarry Smith 
532d91e6319SBarry Smith   Concepts: matrix^nonzero information
533d91e6319SBarry Smith 
534d9274352SBarry Smith .seealso:  MatGetInfo(), MatInfoType
535d91e6319SBarry Smith S*/
5364e220ebcSLois Curfman McInnes typedef struct {
537b0a32e0cSBarry Smith   PetscLogDouble block_size;                         /* block size */
538b0a32e0cSBarry Smith   PetscLogDouble nz_allocated,nz_used,nz_unneeded;   /* number of nonzeros */
539b0a32e0cSBarry Smith   PetscLogDouble memory;                             /* memory allocated */
540b0a32e0cSBarry Smith   PetscLogDouble assemblies;                         /* number of matrix assemblies called */
541b0a32e0cSBarry Smith   PetscLogDouble mallocs;                            /* number of mallocs during MatSetValues() */
542b0a32e0cSBarry Smith   PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */
543b0a32e0cSBarry Smith   PetscLogDouble factor_mallocs;                     /* number of mallocs during factorization */
5444e220ebcSLois Curfman McInnes } MatInfo;
5454e220ebcSLois Curfman McInnes 
546d9274352SBarry Smith /*E
547d9274352SBarry Smith     MatInfoType - Indicates if you want information about the local part of the matrix,
548d9274352SBarry Smith      the entire parallel matrix or the maximum over all the local parts.
549d9274352SBarry Smith 
550d9274352SBarry Smith     Level: beginner
551d9274352SBarry Smith 
552d9274352SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
553d9274352SBarry Smith 
554d9274352SBarry Smith .seealso: MatGetInfo(), MatInfo
555d9274352SBarry Smith E*/
5567b80b807SBarry Smith typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType;
5577087cfbeSBarry Smith extern PetscErrorCode  MatGetInfo(Mat,MatInfoType,MatInfo*);
5587087cfbeSBarry Smith extern PetscErrorCode  MatGetDiagonal(Mat,Vec);
5597087cfbeSBarry Smith extern PetscErrorCode  MatGetRowMax(Mat,Vec,PetscInt[]);
5607087cfbeSBarry Smith extern PetscErrorCode  MatGetRowMin(Mat,Vec,PetscInt[]);
5617087cfbeSBarry Smith extern PetscErrorCode  MatGetRowMaxAbs(Mat,Vec,PetscInt[]);
5627087cfbeSBarry Smith extern PetscErrorCode  MatGetRowMinAbs(Mat,Vec,PetscInt[]);
5637087cfbeSBarry Smith extern PetscErrorCode  MatGetRowSum(Mat,Vec);
5647087cfbeSBarry Smith extern PetscErrorCode  MatTranspose(Mat,MatReuse,Mat*);
565fc4dec0aSBarry Smith PetscPolymorphicFunction(MatTranspose,(Mat A),(A,MAT_INITIAL_MATRIX,&t),Mat,t)
5667087cfbeSBarry Smith extern PetscErrorCode  MatHermitianTranspose(Mat,MatReuse,Mat*);
5677087cfbeSBarry Smith extern PetscErrorCode  MatPermute(Mat,IS,IS,Mat *);
568ba966639SSatish Balay PetscPolymorphicFunction(MatPermute,(Mat A,IS is1,IS is2),(A,is1,is2,&t),Mat,t)
5697087cfbeSBarry Smith extern PetscErrorCode  MatDiagonalScale(Mat,Vec,Vec);
5707087cfbeSBarry Smith extern PetscErrorCode  MatDiagonalSet(Mat,Vec,InsertMode);
5717087cfbeSBarry Smith extern PetscErrorCode  MatEqual(Mat,Mat,PetscBool *);
572ace3abfcSBarry Smith PetscPolymorphicFunction(MatEqual,(Mat A,Mat B),(A,B,&t),PetscBool ,t)
5737087cfbeSBarry Smith extern PetscErrorCode  MatMultEqual(Mat,Mat,PetscInt,PetscBool *);
5747087cfbeSBarry Smith extern PetscErrorCode  MatMultAddEqual(Mat,Mat,PetscInt,PetscBool *);
5757087cfbeSBarry Smith extern PetscErrorCode  MatMultTransposeEqual(Mat,Mat,PetscInt,PetscBool *);
5767087cfbeSBarry Smith extern PetscErrorCode  MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscBool *);
5777b80b807SBarry Smith 
5787087cfbeSBarry Smith extern PetscErrorCode  MatNorm(Mat,NormType,PetscReal *);
579ba966639SSatish Balay PetscPolymorphicFunction(MatNorm,(Mat A,NormType t),(A,t,&n),PetscReal,n)
58009573ac7SBarry Smith extern PetscErrorCode MatGetColumnNorms(Mat,NormType,PetscReal *);
5817087cfbeSBarry Smith extern PetscErrorCode  MatZeroEntries(Mat);
5827087cfbeSBarry Smith extern PetscErrorCode  MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
5837087cfbeSBarry Smith extern PetscErrorCode  MatZeroRowsIS(Mat,IS,PetscScalar,Vec,Vec);
5847087cfbeSBarry Smith extern PetscErrorCode  MatZeroRowsStencil(Mat,PetscInt,const MatStencil [],PetscScalar,Vec,Vec);
5857087cfbeSBarry Smith extern PetscErrorCode  MatZeroRowsColumns(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
5867087cfbeSBarry Smith extern PetscErrorCode  MatZeroRowsColumnsIS(Mat,IS,PetscScalar,Vec,Vec);
5877b80b807SBarry Smith 
5887087cfbeSBarry Smith extern PetscErrorCode  MatUseScaledForm(Mat,PetscBool );
5897087cfbeSBarry Smith extern PetscErrorCode  MatScaleSystem(Mat,Vec,Vec);
5907087cfbeSBarry Smith extern PetscErrorCode  MatUnScaleSystem(Mat,Vec,Vec);
5915ef9f2a5SBarry Smith 
5927087cfbeSBarry Smith extern PetscErrorCode  MatGetSize(Mat,PetscInt*,PetscInt*);
5937087cfbeSBarry Smith extern PetscErrorCode  MatGetLocalSize(Mat,PetscInt*,PetscInt*);
5947087cfbeSBarry Smith extern PetscErrorCode  MatGetOwnershipRange(Mat,PetscInt*,PetscInt*);
5957087cfbeSBarry Smith extern PetscErrorCode  MatGetOwnershipRanges(Mat,const PetscInt**);
5967087cfbeSBarry Smith extern PetscErrorCode  MatGetOwnershipRangeColumn(Mat,PetscInt*,PetscInt*);
5977087cfbeSBarry Smith extern PetscErrorCode  MatGetOwnershipRangesColumn(Mat,const PetscInt**);
5987b80b807SBarry Smith 
5997087cfbeSBarry Smith extern PetscErrorCode  MatGetSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
6007087cfbeSBarry Smith extern PetscErrorCode  MatGetSubMatricesParallel(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
6017087cfbeSBarry Smith extern PetscErrorCode  MatDestroyMatrices(PetscInt,Mat *[]);
6027087cfbeSBarry Smith extern PetscErrorCode  MatGetSubMatrix(Mat,IS,IS,MatReuse,Mat *);
6037087cfbeSBarry Smith extern PetscErrorCode  MatGetLocalSubMatrix(Mat,IS,IS,Mat*);
6047087cfbeSBarry Smith extern PetscErrorCode  MatRestoreLocalSubMatrix(Mat,IS,IS,Mat*);
6057087cfbeSBarry Smith extern PetscErrorCode  MatGetSeqNonzeroStructure(Mat,Mat*);
6067087cfbeSBarry Smith extern PetscErrorCode  MatDestroySeqNonzeroStructure(Mat*);
6075494a064SHong Zhang 
6087087cfbeSBarry Smith extern PetscErrorCode  MatMerge(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
6097087cfbeSBarry Smith extern PetscErrorCode  MatMerge_SeqsToMPI(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*);
6107087cfbeSBarry Smith extern PetscErrorCode  MatMerge_SeqsToMPISymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*);
6117087cfbeSBarry Smith extern PetscErrorCode  MatMerge_SeqsToMPINumeric(Mat,Mat);
6127087cfbeSBarry Smith extern PetscErrorCode  MatDestroy_MPIAIJ_SeqsToMPI(Mat);
6137087cfbeSBarry Smith extern PetscErrorCode  MatGetLocalMat(Mat,MatReuse,Mat*);
6147087cfbeSBarry Smith extern PetscErrorCode  MatGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*);
6157087cfbeSBarry Smith extern PetscErrorCode  MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,PetscInt*,Mat*);
6167087cfbeSBarry Smith extern PetscErrorCode  MatGetBrowsOfAoCols(Mat,Mat,MatReuse,PetscInt**,PetscInt**,MatScalar**,Mat*);
61743eb5e2fSMatthew Knepley #if defined (PETSC_USE_CTABLE)
618cd116e26SSatish Balay #include "petscctable.h"
6197087cfbeSBarry Smith extern PetscErrorCode  MatGetCommunicationStructs(Mat, Vec *, PetscTable *, VecScatter *);
62043eb5e2fSMatthew Knepley #else
6217087cfbeSBarry Smith extern PetscErrorCode  MatGetCommunicationStructs(Mat, Vec *, PetscInt *[], VecScatter *);
62243eb5e2fSMatthew Knepley #endif
6237087cfbeSBarry Smith extern PetscErrorCode  MatGetGhosts(Mat, PetscInt *,const PetscInt *[]);
6248efafbd8SBarry Smith 
6257087cfbeSBarry Smith extern PetscErrorCode  MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt);
6267b80b807SBarry Smith 
6277087cfbeSBarry Smith extern PetscErrorCode  MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*);
6287087cfbeSBarry Smith extern PetscErrorCode  MatMatMultSymbolic(Mat,Mat,PetscReal,Mat*);
6297087cfbeSBarry Smith extern PetscErrorCode  MatMatMultNumeric(Mat,Mat,Mat);
63022440eb1SKris Buschelman 
6317087cfbeSBarry Smith extern PetscErrorCode  MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*);
6327087cfbeSBarry Smith extern PetscErrorCode  MatPtAPSymbolic(Mat,Mat,PetscReal,Mat*);
6337087cfbeSBarry Smith extern PetscErrorCode  MatPtAPNumeric(Mat,Mat,Mat);
63422440eb1SKris Buschelman 
6357087cfbeSBarry Smith extern PetscErrorCode  MatMatMultTranspose(Mat,Mat,MatReuse,PetscReal,Mat*);
6367087cfbeSBarry Smith extern PetscErrorCode  MatMatMultTransposeSymbolic(Mat,Mat,PetscReal,Mat*);
6377087cfbeSBarry Smith extern PetscErrorCode  MatMatMultTransposeNumeric(Mat,Mat,Mat);
638bc011b1eSHong Zhang 
6397087cfbeSBarry Smith extern PetscErrorCode  MatAXPY(Mat,PetscScalar,Mat,MatStructure);
6407087cfbeSBarry Smith extern PetscErrorCode  MatAYPX(Mat,PetscScalar,Mat,MatStructure);
6411c741599SHong Zhang 
6427087cfbeSBarry Smith extern PetscErrorCode  MatScale(Mat,PetscScalar);
6437087cfbeSBarry Smith extern PetscErrorCode  MatShift(Mat,PetscScalar);
6447b80b807SBarry Smith 
6457087cfbeSBarry Smith extern PetscErrorCode  MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping);
6467087cfbeSBarry Smith extern PetscErrorCode  MatSetLocalToGlobalMappingBlock(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping);
6477087cfbeSBarry Smith extern PetscErrorCode  MatGetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*);
6487087cfbeSBarry Smith extern PetscErrorCode  MatGetLocalToGlobalMappingBlock(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*);
6497087cfbeSBarry Smith extern PetscErrorCode  MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
6507087cfbeSBarry Smith extern PetscErrorCode  MatZeroRowsLocalIS(Mat,IS,PetscScalar,Vec,Vec);
6517087cfbeSBarry Smith extern PetscErrorCode  MatZeroRowsColumnsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
6527087cfbeSBarry Smith extern PetscErrorCode  MatZeroRowsColumnsLocalIS(Mat,IS,PetscScalar,Vec,Vec);
6537087cfbeSBarry Smith extern PetscErrorCode  MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
6547087cfbeSBarry Smith extern PetscErrorCode  MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
655052efed2SBarry Smith 
6567087cfbeSBarry Smith extern PetscErrorCode  MatStashSetInitialSize(Mat,PetscInt,PetscInt);
6577087cfbeSBarry Smith extern PetscErrorCode  MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
65890f02eecSBarry Smith 
6597087cfbeSBarry Smith extern PetscErrorCode  MatInterpolate(Mat,Vec,Vec);
6607087cfbeSBarry Smith extern PetscErrorCode  MatInterpolateAdd(Mat,Vec,Vec,Vec);
661ba966639SSatish Balay PetscPolymorphicSubroutine(MatInterpolateAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
6627087cfbeSBarry Smith extern PetscErrorCode  MatRestrict(Mat,Vec,Vec);
6637087cfbeSBarry Smith extern PetscErrorCode  MatGetVecs(Mat,Vec*,Vec*);
6647087cfbeSBarry Smith extern PetscErrorCode  MatGetRedundantMatrix(Mat,PetscInt,MPI_Comm,PetscInt,MatReuse,Mat*);
6657087cfbeSBarry Smith extern PetscErrorCode  MatGetMultiProcBlock(Mat,MPI_Comm,Mat*);
666bd481603SBarry Smith 
667bd481603SBarry Smith /*MC
6681620fd73SBarry Smith    MatSetValue - Set a single entry into a matrix.
6691620fd73SBarry Smith 
6701620fd73SBarry Smith    Not collective
6711620fd73SBarry Smith 
6721620fd73SBarry Smith    Input Parameters:
6731620fd73SBarry Smith +  m - the matrix
6741620fd73SBarry Smith .  row - the row location of the entry
6751620fd73SBarry Smith .  col - the column location of the entry
6761620fd73SBarry Smith .  value - the value to insert
6771620fd73SBarry Smith -  mode - either INSERT_VALUES or ADD_VALUES
6781620fd73SBarry Smith 
6791620fd73SBarry Smith    Notes:
6801620fd73SBarry Smith    For efficiency one should use MatSetValues() and set several or many
6811620fd73SBarry Smith    values simultaneously if possible.
6821620fd73SBarry Smith 
6831620fd73SBarry Smith    Level: beginner
6841620fd73SBarry Smith 
6851620fd73SBarry Smith .seealso: MatSetValues(), MatSetValueLocal()
6861620fd73SBarry Smith M*/
6871620fd73SBarry 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);}
6881620fd73SBarry Smith 
6891620fd73SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatGetValue(Mat v,PetscInt i,PetscInt j,PetscScalar *va) {return MatGetValues(v,1,&i,1,&j,va);}
6901620fd73SBarry Smith 
6911620fd73SBarry 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);}
6921620fd73SBarry Smith 
6931620fd73SBarry Smith /*MC
6940d2a7ff2SSatish Balay    MatPreallocateInitialize - Begins the block of code that will count the number of nonzeros per
695bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
696bd481603SBarry Smith 
697bd481603SBarry Smith    Synopsis:
698c1ac3661SBarry Smith    PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)
699bd481603SBarry Smith 
700bd481603SBarry Smith    Collective on MPI_Comm
701bd481603SBarry Smith 
702bd481603SBarry Smith    Input Parameters:
703bd481603SBarry Smith +  comm - the communicator that will share the eventually allocated matrix
704859fcb39SBarry Smith .  nrows - the number of LOCAL rows in the matrix
705859fcb39SBarry Smith -  ncols - the number of LOCAL columns in the matrix
706bd481603SBarry Smith 
707bd481603SBarry Smith    Output Parameters:
708bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
709bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
710bd481603SBarry Smith 
711bd481603SBarry Smith 
712bd481603SBarry Smith    Level: intermediate
713bd481603SBarry Smith 
714bd481603SBarry Smith    Notes:
7150598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
716bd481603SBarry Smith 
7171d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
718bd481603SBarry Smith 
719bd481603SBarry Smith    Use MatPreallocateInitializeSymmetric() for symmetric matrices (MPISBAIJ matrices)
720bd481603SBarry Smith 
7211620fd73SBarry Smith    This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize().
7221620fd73SBarry Smith 
723bd481603SBarry Smith   Concepts: preallocation^Matrix
724bd481603SBarry Smith 
725bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
726bd481603SBarry Smith           MatPreallocateInitializeSymmetric(), MatPreallocateSymmetricSetLocal()
727bd481603SBarry Smith M*/
728c4f061fbSSatish Balay #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \
7297c922b88SBarry Smith { \
730a89bc211SBarry Smith   PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__start,__end; \
731a89bc211SBarry Smith   _4_ierr = PetscMalloc2(__nrows,PetscInt,&dnz,__nrows,PetscInt,&onz);CHKERRQ(_4_ierr); \
732a89bc211SBarry Smith   _4_ierr = PetscMemzero(dnz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
733a89bc211SBarry Smith   _4_ierr = PetscMemzero(onz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
734c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ctmp;\
735a89bc211SBarry Smith   _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows;
7367c922b88SBarry Smith 
737bd481603SBarry Smith /*MC
7380d2a7ff2SSatish Balay    MatPreallocateSymmetricInitialize - Begins the block of code that will count the number of nonzeros per
739bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
740bd481603SBarry Smith 
741bd481603SBarry Smith    Synopsis:
742c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)
743bd481603SBarry Smith 
744bd481603SBarry Smith    Collective on MPI_Comm
745bd481603SBarry Smith 
746bd481603SBarry Smith    Input Parameters:
747bd481603SBarry Smith +  comm - the communicator that will share the eventually allocated matrix
748859fcb39SBarry Smith .  nrows - the number of LOCAL rows in the matrix
749859fcb39SBarry Smith -  ncols - the number of LOCAL columns in the matrix
750bd481603SBarry Smith 
751bd481603SBarry Smith    Output Parameters:
752bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
753bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
754bd481603SBarry Smith 
755bd481603SBarry Smith 
756bd481603SBarry Smith    Level: intermediate
757bd481603SBarry Smith 
758bd481603SBarry Smith    Notes:
7590598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
760bd481603SBarry Smith 
7611d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
762bd481603SBarry Smith 
7631620fd73SBarry Smith    This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize().
7641620fd73SBarry Smith 
765bd481603SBarry Smith   Concepts: preallocation^Matrix
766bd481603SBarry Smith 
767bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
768bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal()
769bd481603SBarry Smith M*/
770222b16d4SBarry Smith #define MatPreallocateSymmetricInitialize(comm,nrows,ncols,dnz,onz) 0; \
771222b16d4SBarry Smith { \
772a89bc211SBarry Smith   PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__end; \
773a89bc211SBarry Smith   _4_ierr = PetscMalloc2(__nrows,PetscInt,&dnz,__nrows,PetscInt,&onz);CHKERRQ(_4_ierr); \
774a89bc211SBarry Smith   _4_ierr = PetscMemzero(dnz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
775a89bc211SBarry Smith   _4_ierr = PetscMemzero(onz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
776c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr);\
777a89bc211SBarry Smith   _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows;
778222b16d4SBarry Smith 
779bd481603SBarry Smith /*MC
780bd481603SBarry Smith    MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
781bd481603SBarry Smith        inserted using a local number of the rows and columns
782bd481603SBarry Smith 
783bd481603SBarry Smith    Synopsis:
784c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
785bd481603SBarry Smith 
786bd481603SBarry Smith    Not Collective
787bd481603SBarry Smith 
788bd481603SBarry Smith    Input Parameters:
789784ac674SJed Brown +  map - the row mapping from local numbering to global numbering
790bd481603SBarry Smith .  nrows - the number of rows indicated
7911d73ed98SBarry Smith .  rows - the indices of the rows
792784ac674SJed Brown .  cmap - the column mapping from local to global numbering
793bd481603SBarry Smith .  ncols - the number of columns in the matrix
794bd481603SBarry Smith .  cols - the columns indicated
795bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
796bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
797bd481603SBarry Smith 
798bd481603SBarry Smith 
799bd481603SBarry Smith    Level: intermediate
800bd481603SBarry Smith 
801bd481603SBarry Smith    Notes:
8020598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
803bd481603SBarry Smith 
8041d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
805bd481603SBarry Smith 
806bd481603SBarry Smith   Concepts: preallocation^Matrix
807bd481603SBarry Smith 
808bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
809bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal()
810bd481603SBarry Smith M*/
811784ac674SJed Brown #define MatPreallocateSetLocal(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \
812c4f061fbSSatish Balay {\
813c1ac3661SBarry Smith   PetscInt __l;\
814784ac674SJed Brown   _4_ierr = ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\
815784ac674SJed Brown   _4_ierr = ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\
816c4f061fbSSatish Balay   for (__l=0;__l<nrows;__l++) {\
817ef66eb69SBarry Smith     _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
818c4f061fbSSatish Balay   }\
819c4f061fbSSatish Balay }
820c4f061fbSSatish Balay 
821bd481603SBarry Smith /*MC
822bd481603SBarry Smith    MatPreallocateSymmetricSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
823bd481603SBarry Smith        inserted using a local number of the rows and columns
824bd481603SBarry Smith 
825bd481603SBarry Smith    Synopsis:
826c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
827bd481603SBarry Smith 
828bd481603SBarry Smith    Not Collective
829bd481603SBarry Smith 
830bd481603SBarry Smith    Input Parameters:
831bd481603SBarry Smith +  map - the mapping between local numbering and global numbering
832bd481603SBarry Smith .  nrows - the number of rows indicated
8331d73ed98SBarry Smith .  rows - the indices of the rows
834bd481603SBarry Smith .  ncols - the number of columns in the matrix
835bd481603SBarry Smith .  cols - the columns indicated
836bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
837bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
838bd481603SBarry Smith 
839bd481603SBarry Smith 
840bd481603SBarry Smith    Level: intermediate
841bd481603SBarry Smith 
842bd481603SBarry Smith    Notes:
8430598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
844bd481603SBarry Smith 
845bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
846bd481603SBarry Smith 
847bd481603SBarry Smith   Concepts: preallocation^Matrix
848bd481603SBarry Smith 
849bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
850bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
851bd481603SBarry Smith M*/
852d3d32019SBarry Smith #define MatPreallocateSymmetricSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\
853d3d32019SBarry Smith {\
854c1ac3661SBarry Smith   PetscInt __l;\
855d3d32019SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\
856d3d32019SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\
857d3d32019SBarry Smith   for (__l=0;__l<nrows;__l++) {\
858d3d32019SBarry Smith     _4_ierr = MatPreallocateSymmetricSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
859d3d32019SBarry Smith   }\
860d3d32019SBarry Smith }
861d3d32019SBarry Smith 
862bd481603SBarry Smith /*MC
863bd481603SBarry Smith    MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
864bd481603SBarry Smith        inserted using a local number of the rows and columns
865bd481603SBarry Smith 
866bd481603SBarry Smith    Synopsis:
867c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
868bd481603SBarry Smith 
869bd481603SBarry Smith    Not Collective
870bd481603SBarry Smith 
871bd481603SBarry Smith    Input Parameters:
87264075487SBarry Smith +  row - the row
873bd481603SBarry Smith .  ncols - the number of columns in the matrix
874a50f8bf6SHong Zhang -  cols - the columns indicated
875a50f8bf6SHong Zhang 
876a50f8bf6SHong Zhang    Output Parameters:
877a50f8bf6SHong Zhang +  dnz - the array that will be passed to the matrix preallocation routines
878bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
879bd481603SBarry Smith 
880bd481603SBarry Smith 
881bd481603SBarry Smith    Level: intermediate
882bd481603SBarry Smith 
883bd481603SBarry Smith    Notes:
8840598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
885bd481603SBarry Smith 
886bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
887bd481603SBarry Smith 
8881620fd73SBarry Smith    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
8891620fd73SBarry Smith 
890bd481603SBarry Smith   Concepts: preallocation^Matrix
891bd481603SBarry Smith 
892bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
893bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
894bd481603SBarry Smith M*/
895c4f061fbSSatish Balay #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\
896c1ac3661SBarry Smith { PetscInt __i; \
897e32f2f54SBarry 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);\
898e32f2f54SBarry 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);\
8997c922b88SBarry Smith   for (__i=0; __i<nc; __i++) {\
90064075487SBarry Smith     if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \
9017cc688d7SBarry Smith     else dnz[row - __rstart]++;\
9027c922b88SBarry Smith   }\
9037c922b88SBarry Smith }
9047c922b88SBarry Smith 
905bd481603SBarry Smith /*MC
906bd481603SBarry Smith    MatPreallocateSymmetricSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
907bd481603SBarry Smith        inserted using a local number of the rows and columns
908bd481603SBarry Smith 
909bd481603SBarry Smith    Synopsis:
910c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
911bd481603SBarry Smith 
912bd481603SBarry Smith    Not Collective
913bd481603SBarry Smith 
914bd481603SBarry Smith    Input Parameters:
915bd481603SBarry Smith +  nrows - the number of rows indicated
9161d73ed98SBarry Smith .  rows - the indices of the rows
917bd481603SBarry Smith .  ncols - the number of columns in the matrix
918bd481603SBarry Smith .  cols - the columns indicated
919bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
920bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
921bd481603SBarry Smith 
922bd481603SBarry Smith 
923bd481603SBarry Smith    Level: intermediate
924bd481603SBarry Smith 
925bd481603SBarry Smith    Notes:
9260598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
927bd481603SBarry Smith 
928bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
929bd481603SBarry Smith 
9301620fd73SBarry Smith    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
9311620fd73SBarry Smith 
932bd481603SBarry Smith   Concepts: preallocation^Matrix
933bd481603SBarry Smith 
934bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
935bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
936bd481603SBarry Smith M*/
937d3d32019SBarry Smith #define MatPreallocateSymmetricSet(row,nc,cols,dnz,onz) 0;\
938c1ac3661SBarry Smith { PetscInt __i; \
939d3d32019SBarry Smith   for (__i=0; __i<nc; __i++) {\
940d3d32019SBarry Smith     if (cols[__i] >= __end) onz[row - __rstart]++; \
941d3d32019SBarry Smith     else if (cols[__i] >= row) dnz[row - __rstart]++;\
942d3d32019SBarry Smith   }\
943d3d32019SBarry Smith }
944d3d32019SBarry Smith 
945bd481603SBarry Smith /*MC
94616371a99SBarry Smith    MatPreallocateLocation -  An alternative to MatPreallocationSet() that puts the nonzero locations into the matrix if it exists
94716371a99SBarry Smith 
94816371a99SBarry Smith    Synopsis:
94916371a99SBarry Smith    PetscErrorCode MatPreallocateLocations(Mat A,PetscInt row,PetscInt ncols,PetscInt *cols,PetscInt *dnz,PetscInt *onz)
95016371a99SBarry Smith 
95116371a99SBarry Smith    Not Collective
95216371a99SBarry Smith 
95316371a99SBarry Smith    Input Parameters:
95416371a99SBarry Smith .  A - matrix
95516371a99SBarry Smith .  row - row where values exist (must be local to this process)
95616371a99SBarry Smith .  ncols - number of columns
95716371a99SBarry Smith .  cols - columns with nonzeros
95816371a99SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
95916371a99SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
96016371a99SBarry Smith 
96116371a99SBarry Smith 
96216371a99SBarry Smith    Level: intermediate
96316371a99SBarry Smith 
96416371a99SBarry Smith    Notes:
9650598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
96616371a99SBarry Smith 
96716371a99SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
96816371a99SBarry Smith 
96916371a99SBarry Smith    This is a MACRO not a function because it uses a bunch of variables private to the MatPreallocation.... routines.
97016371a99SBarry Smith 
97116371a99SBarry Smith   Concepts: preallocation^Matrix
97216371a99SBarry Smith 
97316371a99SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
97416371a99SBarry Smith           MatPreallocateSymmetricInitialize(), MatPreallocateSymmetricSetLocal()
97516371a99SBarry Smith M*/
97616371a99SBarry 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);}
97716371a99SBarry Smith 
97816371a99SBarry Smith 
97916371a99SBarry Smith /*MC
9800d2a7ff2SSatish Balay    MatPreallocateFinalize - Ends the block of code that will count the number of nonzeros per
981bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
982bd481603SBarry Smith 
983bd481603SBarry Smith    Synopsis:
984c1ac3661SBarry Smith    PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz)
985bd481603SBarry Smith 
986bd481603SBarry Smith    Collective on MPI_Comm
987bd481603SBarry Smith 
988bd481603SBarry Smith    Input Parameters:
98916371a99SBarry Smith +  dnz - the array that was be passed to the matrix preallocation routines
990bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
991bd481603SBarry Smith 
992bd481603SBarry Smith 
993bd481603SBarry Smith    Level: intermediate
994bd481603SBarry Smith 
995bd481603SBarry Smith    Notes:
9960598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
997bd481603SBarry Smith 
998bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
999bd481603SBarry Smith 
10001620fd73SBarry Smith    This is a MACRO not a function because it closes the { started in MatPreallocateInitialize().
10011620fd73SBarry Smith 
1002bd481603SBarry Smith   Concepts: preallocation^Matrix
1003bd481603SBarry Smith 
1004bd481603SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
1005bd481603SBarry Smith           MatPreallocateSymmetricInitialize(), MatPreallocateSymmetricSetLocal()
1006bd481603SBarry Smith M*/
1007a89bc211SBarry Smith #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree2(dnz,onz);CHKERRQ(_4_ierr);}
10087c922b88SBarry Smith 
1009bd481603SBarry Smith 
1010bd481603SBarry Smith 
10117b80b807SBarry Smith /* Routines unique to particular data structures */
10127087cfbeSBarry Smith extern PetscErrorCode  MatShellGetContext(Mat,void **);
1013ba966639SSatish Balay PetscPolymorphicFunction(MatShellGetContext,(Mat A),(A,&t),void*,t)
1014698d4c6aSKris Buschelman 
10157087cfbeSBarry Smith extern PetscErrorCode  MatInodeAdjustForInodes(Mat,IS*,IS*);
10167087cfbeSBarry Smith extern PetscErrorCode  MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *);
1017698d4c6aSKris Buschelman 
10187087cfbeSBarry Smith extern PetscErrorCode  MatSeqAIJSetColumnIndices(Mat,PetscInt[]);
10197087cfbeSBarry Smith extern PetscErrorCode  MatSeqBAIJSetColumnIndices(Mat,PetscInt[]);
10207087cfbeSBarry Smith extern PetscErrorCode  MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
10217087cfbeSBarry Smith extern PetscErrorCode  MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
10227087cfbeSBarry Smith extern PetscErrorCode  MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
10237b80b807SBarry Smith 
1024a96a251dSBarry Smith #define MAT_SKIP_ALLOCATION -4
1025a96a251dSBarry Smith 
10267087cfbeSBarry Smith extern PetscErrorCode  MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
1027ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz))
10287087cfbeSBarry Smith extern PetscErrorCode  MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
1029ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqSBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz))
10307087cfbeSBarry Smith extern PetscErrorCode  MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]);
1031ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqAIJSetPreallocation,(Mat A,const PetscInt nnz[]),(A,0,nnz))
1032273d9f13SBarry Smith 
10337087cfbeSBarry Smith extern PetscErrorCode  MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
1034ba966639SSatish Balay PetscPolymorphicSubroutine(MatMPIBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[],const PetscInt onz[]),(A,bs,0,nnz,0,onz))
10357087cfbeSBarry Smith extern PetscErrorCode  MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
10367087cfbeSBarry Smith extern PetscErrorCode  MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
10377087cfbeSBarry Smith extern PetscErrorCode  MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []);
10387087cfbeSBarry Smith extern PetscErrorCode  MatSeqBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
10397087cfbeSBarry Smith extern PetscErrorCode  MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]);
10407087cfbeSBarry Smith extern PetscErrorCode  MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
10417087cfbeSBarry Smith extern PetscErrorCode  MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]);
10427087cfbeSBarry Smith extern PetscErrorCode  MatMPIDenseSetPreallocation(Mat,PetscScalar[]);
10437087cfbeSBarry Smith extern PetscErrorCode  MatSeqDenseSetPreallocation(Mat,PetscScalar[]);
10447087cfbeSBarry Smith extern PetscErrorCode  MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,PetscInt*[]);
10457087cfbeSBarry Smith extern PetscErrorCode  MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,PetscInt*[]);
10467087cfbeSBarry Smith extern PetscErrorCode  MatAdicSetLocalFunction(Mat,void (*)(void));
1047273d9f13SBarry Smith 
10487087cfbeSBarry Smith extern PetscErrorCode  MatSeqDenseSetLDA(Mat,PetscInt);
10497087cfbeSBarry Smith extern PetscErrorCode  MatDenseGetLocalMatrix(Mat,Mat*);
10501b807ce4Svictorle 
10517087cfbeSBarry Smith extern PetscErrorCode  MatStoreValues(Mat);
10527087cfbeSBarry Smith extern PetscErrorCode  MatRetrieveValues(Mat);
10532e8a6d31SBarry Smith 
10547087cfbeSBarry Smith extern PetscErrorCode  MatDAADSetCtx(Mat,void*);
10553a7fca6bSBarry Smith 
10567b80b807SBarry Smith /*
10577b80b807SBarry Smith   These routines are not usually accessed directly, rather solving is
105894b7f48cSBarry Smith   done through the KSP and PC interfaces.
10597b80b807SBarry Smith */
10607b80b807SBarry Smith 
1061d9274352SBarry Smith /*E
1062d9274352SBarry Smith     MatOrderingType - String with the name of a PETSc matrix ordering or the creation function
1063d9274352SBarry Smith        with an optional dynamic library name, for example
1064d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:orderingcreate()
1065d9274352SBarry Smith 
1066d9274352SBarry Smith    Level: beginner
1067d9274352SBarry Smith 
1068e5a9bf91SBarry Smith    Cannot use const because the PC objects manipulate the string
1069e5a9bf91SBarry Smith 
1070d9274352SBarry Smith .seealso: MatGetOrdering()
1071d9274352SBarry Smith E*/
10723eaad2caSSatish Balay #define MatOrderingType char*
10732692d6eeSBarry Smith #define MATORDERINGNATURAL     "natural"
10742692d6eeSBarry Smith #define MATORDERINGND          "nd"
10752692d6eeSBarry Smith #define MATORDERING1WD         "1wd"
10762692d6eeSBarry Smith #define MATORDERINGRCM         "rcm"
10772692d6eeSBarry Smith #define MATORDERINGQMD         "qmd"
10782692d6eeSBarry Smith #define MATORDERINGROWLENGTH   "rowlength"
10792692d6eeSBarry Smith #define MATORDERINGDSC_ND      "dsc_nd"         /* these three are only for DSCPACK, see its documentation for details */
10802692d6eeSBarry Smith #define MATORDERINGDSC_MMD     "dsc_mmd"
10812692d6eeSBarry Smith #define MATORDERINGDSC_MDF     "dsc_mdf"
1082312e372aSBarry Smith #define MATORDERINGAMD         "amd"            /* only works if UMFPACK is installed with PETSc */
1083b12f92e5SBarry Smith 
10847087cfbeSBarry Smith extern PetscErrorCode  MatGetOrdering(Mat,const MatOrderingType,IS*,IS*);
10857087cfbeSBarry Smith extern PetscErrorCode  MatGetOrderingList(PetscFList *list);
10867087cfbeSBarry Smith extern PetscErrorCode  MatOrderingRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,const MatOrderingType,IS*,IS*));
108730de9b25SBarry Smith 
108830de9b25SBarry Smith /*MC
10891890ba74SBarry Smith    MatOrderingRegisterDynamic - Adds a new sparse matrix ordering to the matrix package.
109030de9b25SBarry Smith 
109130de9b25SBarry Smith    Synopsis:
10921890ba74SBarry Smith    PetscErrorCode MatOrderingRegisterDynamic(const char *name_ordering,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatOrdering))
109330de9b25SBarry Smith 
109430de9b25SBarry Smith    Not Collective
109530de9b25SBarry Smith 
109630de9b25SBarry Smith    Input Parameters:
10972692d6eeSBarry Smith +  sname - name of ordering (for example MATORDERINGND)
109830de9b25SBarry Smith .  path - location of library where creation routine is
109930de9b25SBarry Smith .  name - name of function that creates the ordering type,a string
110030de9b25SBarry Smith -  function - function pointer that creates the ordering
110130de9b25SBarry Smith 
110230de9b25SBarry Smith    Level: developer
110330de9b25SBarry Smith 
110430de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
110530de9b25SBarry Smith    is ignored.
110630de9b25SBarry Smith 
110730de9b25SBarry Smith    Sample usage:
110830de9b25SBarry Smith .vb
110930de9b25SBarry Smith    MatOrderingRegisterDynamic("my_order",/home/username/my_lib/lib/libO/solaris/mylib.a,
111030de9b25SBarry Smith                "MyOrder",MyOrder);
111130de9b25SBarry Smith .ve
111230de9b25SBarry Smith 
111330de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
111430de9b25SBarry Smith $     MatOrderingSetType(part,"my_order)
111530de9b25SBarry Smith    or at runtime via the option
11162401956bSBarry Smith $     -pc_factor_mat_ordering_type my_order
111730de9b25SBarry Smith 
1118ab901514SBarry Smith    ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values.
111930de9b25SBarry Smith 
112030de9b25SBarry Smith .keywords: matrix, ordering, register
112130de9b25SBarry Smith 
112230de9b25SBarry Smith .seealso: MatOrderingRegisterDestroy(), MatOrderingRegisterAll()
112330de9b25SBarry Smith M*/
1124aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1125f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,0)
1126b12f92e5SBarry Smith #else
1127f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,d)
1128b12f92e5SBarry Smith #endif
112930de9b25SBarry Smith 
11307087cfbeSBarry Smith extern PetscErrorCode  MatOrderingRegisterDestroy(void);
11317087cfbeSBarry Smith extern PetscErrorCode  MatOrderingRegisterAll(const char[]);
1132ace3abfcSBarry Smith extern PetscBool  MatOrderingRegisterAllCalled;
1133b0a32e0cSBarry Smith extern PetscFList MatOrderingList;
1134d4fbbf0eSBarry Smith 
11357087cfbeSBarry Smith extern PetscErrorCode  MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS);
1136a2ce50c7SBarry Smith 
1137d91e6319SBarry Smith /*S
1138d90ac83dSHong Zhang     MatFactorShiftType - Numeric Shift.
1139d90ac83dSHong Zhang 
1140d90ac83dSHong Zhang    Level: beginner
1141d90ac83dSHong Zhang 
1142d90ac83dSHong Zhang S*/
1143d90ac83dSHong Zhang typedef enum {MAT_SHIFT_NONE,MAT_SHIFT_NONZERO,MAT_SHIFT_POSITIVE_DEFINITE,MAT_SHIFT_INBLOCKS} MatFactorShiftType;
1144d90ac83dSHong Zhang extern const char *MatFactorShiftTypes[];
1145d90ac83dSHong Zhang 
1146d90ac83dSHong Zhang /*S
11472401956bSBarry Smith    MatFactorInfo - Data passed into the matrix factorization routines
1148b00f7748SHong Zhang 
114961cad860SBarry Smith    In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE, that is use
115061cad860SBarry Smith $     MatFactorInfo  info(MAT_FACTORINFO_SIZE)
1151b00f7748SHong Zhang 
115215e8a5b3SHong Zhang    Notes: These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC.
1153b00f7748SHong Zhang 
115461cad860SBarry Smith       You can use MatFactorInfoInitialize() to set default values.
115561cad860SBarry Smith 
1156b00f7748SHong Zhang    Level: developer
1157b00f7748SHong Zhang 
1158d7d82daaSBarry Smith .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(),
1159d7d82daaSBarry Smith           MatFactorInfoInitialize()
1160b00f7748SHong Zhang 
1161b00f7748SHong Zhang S*/
1162b00f7748SHong Zhang typedef struct {
116315e8a5b3SHong Zhang   PetscReal     diagonal_fill;  /* force diagonal to fill in if initially not filled */
116485317021SBarry Smith   PetscReal     usedt;
116515e8a5b3SHong Zhang   PetscReal     dt;             /* drop tolerance */
1166b00f7748SHong Zhang   PetscReal     dtcol;          /* tolerance for pivoting */
116715e8a5b3SHong Zhang   PetscReal     dtcount;        /* maximum nonzeros to be allowed per row */
116867e28bfeSBarry Smith   PetscReal     fill;           /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */
1169348344bbSBarry Smith   PetscReal     levels;         /* ICC/ILU(levels) */
1170bcd9e38bSBarry Smith   PetscReal     pivotinblocks;  /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0
1171bcd9e38bSBarry Smith                                    factorization may be faster if do not pivot */
117215e8a5b3SHong Zhang   PetscReal     zeropivot;      /* pivot is called zero if less than this */
1173f4db908eSBarry Smith   PetscReal     shifttype;      /* type of shift added to matrix factor to prevent zero pivots */
1174f4db908eSBarry Smith   PetscReal     shiftamount;     /* how large the shift is */
117515e8a5b3SHong Zhang } MatFactorInfo;
1176ffa6d0a5SLois Curfman McInnes 
11777087cfbeSBarry Smith extern PetscErrorCode  MatFactorInfoInitialize(MatFactorInfo*);
11787087cfbeSBarry Smith extern PetscErrorCode  MatCholeskyFactor(Mat,IS,const MatFactorInfo*);
11797087cfbeSBarry Smith extern PetscErrorCode  MatCholeskyFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
11807087cfbeSBarry Smith extern PetscErrorCode  MatCholeskyFactorNumeric(Mat,Mat,const MatFactorInfo*);
11817087cfbeSBarry Smith extern PetscErrorCode  MatLUFactor(Mat,IS,IS,const MatFactorInfo*);
11827087cfbeSBarry Smith extern PetscErrorCode  MatILUFactor(Mat,IS,IS,const MatFactorInfo*);
11837087cfbeSBarry Smith extern PetscErrorCode  MatLUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
11847087cfbeSBarry Smith extern PetscErrorCode  MatILUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
11857087cfbeSBarry Smith extern PetscErrorCode  MatICCFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
11867087cfbeSBarry Smith extern PetscErrorCode  MatICCFactor(Mat,IS,const MatFactorInfo*);
11877087cfbeSBarry Smith extern PetscErrorCode  MatLUFactorNumeric(Mat,Mat,const MatFactorInfo*);
11887087cfbeSBarry Smith extern PetscErrorCode  MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*);
11897087cfbeSBarry Smith extern PetscErrorCode  MatSolve(Mat,Vec,Vec);
11907087cfbeSBarry Smith extern PetscErrorCode  MatForwardSolve(Mat,Vec,Vec);
11917087cfbeSBarry Smith extern PetscErrorCode  MatBackwardSolve(Mat,Vec,Vec);
11927087cfbeSBarry Smith extern PetscErrorCode  MatSolveAdd(Mat,Vec,Vec,Vec);
11937087cfbeSBarry Smith extern PetscErrorCode  MatSolveTranspose(Mat,Vec,Vec);
11947087cfbeSBarry Smith extern PetscErrorCode  MatSolveTransposeAdd(Mat,Vec,Vec,Vec);
11957087cfbeSBarry Smith extern PetscErrorCode  MatSolves(Mat,Vecs,Vecs);
11968ed539a5SBarry Smith 
11977087cfbeSBarry Smith extern PetscErrorCode  MatSetUnfactored(Mat);
1198bb5a7306SBarry Smith 
1199d91e6319SBarry Smith /*E
1200d91e6319SBarry Smith     MatSORType - What type of (S)SOR to perform
1201bb1eb677SSatish Balay 
1202d91e6319SBarry Smith     Level: beginner
1203d91e6319SBarry Smith 
1204d9274352SBarry Smith    May be bitwise ORd together
1205d9274352SBarry Smith 
1206d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
1207d91e6319SBarry Smith 
12084e7234bfSBarry Smith    MatSORType may be bitwise ORd together, so do not change the numbers
12094e7234bfSBarry Smith 
121041f059aeSBarry Smith .seealso: MatSOR()
1211d91e6319SBarry Smith E*/
1212ee50ffe9SBarry Smith typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3,
1213ee50ffe9SBarry Smith               SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8,
1214ee50ffe9SBarry Smith               SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16,
121584cb2905SBarry Smith               SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType;
12167087cfbeSBarry Smith extern PetscErrorCode  MatSOR(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
12178ed539a5SBarry Smith 
1218d4fbbf0eSBarry Smith /*
1219639f9d9dSBarry Smith     These routines are for efficiently computing Jacobians via finite differences.
1220639f9d9dSBarry Smith */
1221b12f92e5SBarry Smith 
1222d9274352SBarry Smith /*E
1223d9274352SBarry Smith     MatColoringType - String with the name of a PETSc matrix coloring or the creation function
1224d9274352SBarry Smith        with an optional dynamic library name, for example
1225d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:coloringcreate()
1226d9274352SBarry Smith 
1227d9274352SBarry Smith    Level: beginner
1228d9274352SBarry Smith 
1229d9274352SBarry Smith .seealso: MatGetColoring()
1230d9274352SBarry Smith E*/
1231a313700dSBarry Smith #define MatColoringType char*
12322692d6eeSBarry Smith #define MATCOLORINGNATURAL "natural"
12332692d6eeSBarry Smith #define MATCOLORINGSL      "sl"
12342692d6eeSBarry Smith #define MATCOLORINGLF      "lf"
12352692d6eeSBarry Smith #define MATCOLORINGID      "id"
1236b12f92e5SBarry Smith 
12377087cfbeSBarry Smith extern PetscErrorCode  MatGetColoring(Mat,const MatColoringType,ISColoring*);
12387087cfbeSBarry Smith extern PetscErrorCode  MatColoringRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,MatColoringType,ISColoring *));
123930de9b25SBarry Smith 
124030de9b25SBarry Smith /*MC
124130de9b25SBarry Smith    MatColoringRegisterDynamic - Adds a new sparse matrix coloring to the
124230de9b25SBarry Smith                                matrix package.
124330de9b25SBarry Smith 
124430de9b25SBarry Smith    Synopsis:
12451890ba74SBarry Smith    PetscErrorCode MatColoringRegisterDynamic(const char *name_coloring,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatColoring))
124630de9b25SBarry Smith 
124730de9b25SBarry Smith    Not Collective
124830de9b25SBarry Smith 
124930de9b25SBarry Smith    Input Parameters:
12502692d6eeSBarry Smith +  sname - name of Coloring (for example MATCOLORINGSL)
125130de9b25SBarry Smith .  path - location of library where creation routine is
125230de9b25SBarry Smith .  name - name of function that creates the Coloring type, a string
125330de9b25SBarry Smith -  function - function pointer that creates the coloring
125430de9b25SBarry Smith 
125530de9b25SBarry Smith    Level: developer
125630de9b25SBarry Smith 
125730de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
125830de9b25SBarry Smith    is ignored.
125930de9b25SBarry Smith 
126030de9b25SBarry Smith    Sample usage:
126130de9b25SBarry Smith .vb
126230de9b25SBarry Smith    MatColoringRegisterDynamic("my_color",/home/username/my_lib/lib/libO/solaris/mylib.a,
126330de9b25SBarry Smith                "MyColor",MyColor);
126430de9b25SBarry Smith .ve
126530de9b25SBarry Smith 
126630de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
126730de9b25SBarry Smith $     MatColoringSetType(part,"my_color")
126830de9b25SBarry Smith    or at runtime via the option
126930de9b25SBarry Smith $     -mat_coloring_type my_color
127030de9b25SBarry Smith 
1271ab901514SBarry Smith    $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
127230de9b25SBarry Smith 
127330de9b25SBarry Smith .keywords: matrix, Coloring, register
127430de9b25SBarry Smith 
127530de9b25SBarry Smith .seealso: MatColoringRegisterDestroy(), MatColoringRegisterAll()
127630de9b25SBarry Smith M*/
1277aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1278f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,0)
1279b12f92e5SBarry Smith #else
1280f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,d)
1281b12f92e5SBarry Smith #endif
128230de9b25SBarry Smith 
1283ace3abfcSBarry Smith extern PetscBool  MatColoringRegisterAllCalled;
1284f1144a30SSatish Balay 
12857087cfbeSBarry Smith extern PetscErrorCode  MatColoringRegisterAll(const char[]);
12867087cfbeSBarry Smith extern PetscErrorCode  MatColoringRegisterDestroy(void);
12877087cfbeSBarry Smith extern PetscErrorCode  MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*);
1288639f9d9dSBarry Smith 
1289d9274352SBarry Smith /*S
1290d9274352SBarry Smith      MatFDColoring - Object for computing a sparse Jacobian via finite differences
1291d9274352SBarry Smith         and coloring
1292639f9d9dSBarry Smith 
1293d9274352SBarry Smith    Level: beginner
1294d9274352SBarry Smith 
1295d9274352SBarry Smith   Concepts: coloring, sparse Jacobian, finite differences
1296d9274352SBarry Smith 
1297d9274352SBarry Smith .seealso:  MatFDColoringCreate()
1298d9274352SBarry Smith S*/
1299e2a1c21fSSatish Balay typedef struct _p_MatFDColoring* MatFDColoring;
1300639f9d9dSBarry Smith 
13017087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringCreate(Mat,ISColoring,MatFDColoring *);
13027087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringDestroy(MatFDColoring);
13037087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringView(MatFDColoring,PetscViewer);
13047087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*);
13057087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**);
13067087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal);
13077087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringSetFromOptions(MatFDColoring);
13087087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringApply(Mat,MatFDColoring,Vec,MatStructure*,void *);
13097087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringApplyTS(Mat,MatFDColoring,PetscReal,Vec,MatStructure*,void *);
13107087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringSetF(MatFDColoring,Vec);
13117087cfbeSBarry Smith extern PetscErrorCode  MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,PetscInt*[]);
1312639f9d9dSBarry Smith /*
13130752156aSBarry Smith     These routines are for partitioning matrices: currently used only
13143eda8832SBarry Smith   for adjacency matrix, MatCreateMPIAdj().
13150752156aSBarry Smith */
1316ca161407SBarry Smith 
1317d9274352SBarry Smith /*S
1318d9274352SBarry Smith      MatPartitioning - Object for managing the partitioning of a matrix or graph
1319d9274352SBarry Smith 
1320d9274352SBarry Smith    Level: beginner
1321d9274352SBarry Smith 
1322d9274352SBarry Smith   Concepts: partitioning
1323d9274352SBarry Smith 
1324743a1026Svictorle .seealso:  MatPartitioningCreate(), MatPartitioningType
1325d9274352SBarry Smith S*/
132691e9ee9fSBarry Smith typedef struct _p_MatPartitioning* MatPartitioning;
1327d9274352SBarry Smith 
1328d9274352SBarry Smith /*E
13295bc6e7ccSvictorle     MatPartitioningType - String with the name of a PETSc matrix partitioning or the creation function
1330d9274352SBarry Smith        with an optional dynamic library name, for example
1331d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:partitioningcreate()
1332d9274352SBarry Smith 
1333d9274352SBarry Smith    Level: beginner
1334d9274352SBarry Smith 
1335b547a13bSvictorle .seealso: MatPartitioningCreate(), MatPartitioning
1336d9274352SBarry Smith E*/
1337a313700dSBarry Smith #define MatPartitioningType char*
13382692d6eeSBarry Smith #define MATPARTITIONINGCURRENT  "current"
13392692d6eeSBarry Smith #define MATPARTITIONINGSQUARE   "square"
13402692d6eeSBarry Smith #define MATPARTITIONINGPARMETIS "parmetis"
13412692d6eeSBarry Smith #define MATPARTITIONINGCHACO    "chaco"
13422692d6eeSBarry Smith #define MATPARTITIONINGJOSTLE   "jostle"
13432692d6eeSBarry Smith #define MATPARTITIONINGPARTY    "party"
13442692d6eeSBarry Smith #define MATPARTITIONINGSCOTCH   "scotch"
134571306c60Sjroman 
1346ca161407SBarry Smith 
13477087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningCreate(MPI_Comm,MatPartitioning*);
13487087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningSetType(MatPartitioning,const MatPartitioningType);
13497087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningSetNParts(MatPartitioning,PetscInt);
13507087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningSetAdjacency(MatPartitioning,Mat);
13517087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]);
13527087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []);
13537087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningApply(MatPartitioning,IS*);
13547087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningDestroy(MatPartitioning);
13552aabb6bbSBarry Smith 
13567087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatPartitioning));
135730de9b25SBarry Smith 
135830de9b25SBarry Smith /*MC
135930de9b25SBarry Smith    MatPartitioningRegisterDynamic - Adds a new sparse matrix partitioning to the
136030de9b25SBarry Smith    matrix package.
136130de9b25SBarry Smith 
136230de9b25SBarry Smith    Synopsis:
13631890ba74SBarry Smith    PetscErrorCode MatPartitioningRegisterDynamic(const char *name_partitioning,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatPartitioning))
136430de9b25SBarry Smith 
136530de9b25SBarry Smith    Not Collective
136630de9b25SBarry Smith 
136730de9b25SBarry Smith    Input Parameters:
13682692d6eeSBarry Smith +  sname - name of partitioning (for example MATPARTITIONINGCURRENT) or parmetis
136930de9b25SBarry Smith .  path - location of library where creation routine is
137030de9b25SBarry Smith .  name - name of function that creates the partitioning type, a string
137130de9b25SBarry Smith -  function - function pointer that creates the partitioning type
137230de9b25SBarry Smith 
137330de9b25SBarry Smith    Level: developer
137430de9b25SBarry Smith 
137530de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
137630de9b25SBarry Smith    is ignored.
137730de9b25SBarry Smith 
137830de9b25SBarry Smith    Sample usage:
137930de9b25SBarry Smith .vb
138030de9b25SBarry Smith    MatPartitioningRegisterDynamic("my_part",/home/username/my_lib/lib/libO/solaris/mylib.a,
138130de9b25SBarry Smith                "MyPartCreate",MyPartCreate);
138230de9b25SBarry Smith .ve
138330de9b25SBarry Smith 
138430de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
138530de9b25SBarry Smith $     MatPartitioningSetType(part,"my_part")
138630de9b25SBarry Smith    or at runtime via the option
138730de9b25SBarry Smith $     -mat_partitioning_type my_part
138830de9b25SBarry Smith 
1389ab901514SBarry Smith    $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
139030de9b25SBarry Smith 
139130de9b25SBarry Smith .keywords: matrix, partitioning, register
139230de9b25SBarry Smith 
139330de9b25SBarry Smith .seealso: MatPartitioningRegisterDestroy(), MatPartitioningRegisterAll()
139430de9b25SBarry Smith M*/
1395aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1396f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,0)
13972aabb6bbSBarry Smith #else
1398f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,d)
13992aabb6bbSBarry Smith #endif
14002aabb6bbSBarry Smith 
1401ace3abfcSBarry Smith extern PetscBool  MatPartitioningRegisterAllCalled;
1402f1144a30SSatish Balay 
14037087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningRegisterAll(const char[]);
14047087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningRegisterDestroy(void);
14052bad1931SBarry Smith 
14067087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningView(MatPartitioning,PetscViewer);
14077087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningSetFromOptions(MatPartitioning);
14087087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningGetType(MatPartitioning,const MatPartitioningType*);
1409ca161407SBarry Smith 
14107087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningParmetisSetCoarseSequential(MatPartitioning);
14117087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningParmetisGetEdgeCut(MatPartitioning, PetscInt *);
14120752156aSBarry Smith 
14137087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningJostleSetCoarseLevel(MatPartitioning,PetscReal);
14147087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningJostleSetCoarseSequential(MatPartitioning);
141571306c60Sjroman 
1416954b3ec6SBarry Smith typedef enum {MP_CHACO_MULTILEVEL_KL,MP_CHACO_SPECTRAL,MP_CHACO_LINEAR,MP_CHACO_RANDOM, MP_CHACO_SCATTERED} MPChacoGlobalType;
14177087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningChacoSetGlobal(MatPartitioning, MPChacoGlobalType);
141871306c60Sjroman typedef enum { MP_CHACO_KERNIGHAN_LIN, MP_CHACO_NONE } MPChacoLocalType;
14197087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningChacoSetLocal(MatPartitioning, MPChacoLocalType);
14207087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal);
142171306c60Sjroman typedef enum { MP_CHACO_LANCZOS, MP_CHACO_RQI_SYMMLQ } MPChacoEigenType;
14227087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType);
14237087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningChacoSetEigenTol(MatPartitioning, PetscReal);
14247087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningChacoSetEigenNumber(MatPartitioning, PetscInt);
142571306c60Sjroman 
142671306c60Sjroman #define MP_PARTY_OPT "opt"
142771306c60Sjroman #define MP_PARTY_LIN "lin"
142871306c60Sjroman #define MP_PARTY_SCA "sca"
142971306c60Sjroman #define MP_PARTY_RAN "ran"
143071306c60Sjroman #define MP_PARTY_GBF "gbf"
143171306c60Sjroman #define MP_PARTY_GCF "gcf"
143271306c60Sjroman #define MP_PARTY_BUB "bub"
143371306c60Sjroman #define MP_PARTY_DEF "def"
14347087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningPartySetGlobal(MatPartitioning, const char*);
143571306c60Sjroman #define MP_PARTY_HELPFUL_SETS "hs"
143671306c60Sjroman #define MP_PARTY_KERNIGHAN_LIN "kl"
143771306c60Sjroman #define MP_PARTY_NONE "no"
14387087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningPartySetLocal(MatPartitioning, const char*);
14397087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal);
14407087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningPartySetBipart(MatPartitioning,PetscBool );
14417087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscBool );
144271306c60Sjroman 
144371306c60Sjroman typedef enum { MP_SCOTCH_GREEDY, MP_SCOTCH_GPS, MP_SCOTCH_GR_GPS } MPScotchGlobalType;
14447087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningScotchSetArch(MatPartitioning,const char*);
14457087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningScotchSetMultilevel(MatPartitioning);
14467087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningScotchSetGlobal(MatPartitioning,MPScotchGlobalType);
14477087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningScotchSetCoarseLevel(MatPartitioning,PetscReal);
14487087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningScotchSetHostList(MatPartitioning,const char*);
144971306c60Sjroman typedef enum { MP_SCOTCH_KERNIGHAN_LIN, MP_SCOTCH_NONE } MPScotchLocalType;
14507087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningScotchSetLocal(MatPartitioning,MPScotchLocalType);
14517087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningScotchSetMapping(MatPartitioning);
14527087cfbeSBarry Smith extern PetscErrorCode  MatPartitioningScotchSetStrategy(MatPartitioning,char*);
145371306c60Sjroman 
145409573ac7SBarry Smith extern PetscErrorCode MatMeshToVertexGraph(Mat,PetscInt,Mat*);
145509573ac7SBarry Smith extern PetscErrorCode MatMeshToCellGraph(Mat,PetscInt,Mat*);
1456591294e4SBarry Smith 
14570752156aSBarry Smith /*
14580a835dfdSSatish Balay     If you add entries here you must also add them to finclude/petscmat.h
1459d4fbbf0eSBarry Smith */
14601c1c02c0SLois Curfman McInnes typedef enum { MATOP_SET_VALUES=0,
14611c1c02c0SLois Curfman McInnes                MATOP_GET_ROW=1,
14621c1c02c0SLois Curfman McInnes                MATOP_RESTORE_ROW=2,
14631c1c02c0SLois Curfman McInnes                MATOP_MULT=3,
14641c1c02c0SLois Curfman McInnes                MATOP_MULT_ADD=4,
14657c922b88SBarry Smith                MATOP_MULT_TRANSPOSE=5,
14667c922b88SBarry Smith                MATOP_MULT_TRANSPOSE_ADD=6,
14671c1c02c0SLois Curfman McInnes                MATOP_SOLVE=7,
14681c1c02c0SLois Curfman McInnes                MATOP_SOLVE_ADD=8,
14697c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE=9,
14707c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE_ADD=10,
14711c1c02c0SLois Curfman McInnes                MATOP_LUFACTOR=11,
14721c1c02c0SLois Curfman McInnes                MATOP_CHOLESKYFACTOR=12,
1473a714c06dSBarry Smith                MATOP_SOR=13,
14741c1c02c0SLois Curfman McInnes                MATOP_TRANSPOSE=14,
14751c1c02c0SLois Curfman McInnes                MATOP_GETINFO=15,
14761c1c02c0SLois Curfman McInnes                MATOP_EQUAL=16,
14771c1c02c0SLois Curfman McInnes                MATOP_GET_DIAGONAL=17,
14781c1c02c0SLois Curfman McInnes                MATOP_DIAGONAL_SCALE=18,
14791c1c02c0SLois Curfman McInnes                MATOP_NORM=19,
14801c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_BEGIN=20,
14811c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_END=21,
1482d519adbfSMatthew Knepley                MATOP_SET_OPTION=22,
1483d519adbfSMatthew Knepley                MATOP_ZERO_ENTRIES=23,
1484d519adbfSMatthew Knepley                MATOP_ZERO_ROWS=24,
1485d519adbfSMatthew Knepley                MATOP_LUFACTOR_SYMBOLIC=25,
1486d519adbfSMatthew Knepley                MATOP_LUFACTOR_NUMERIC=26,
1487d519adbfSMatthew Knepley                MATOP_CHOLESKY_FACTOR_SYMBOLIC=27,
1488d519adbfSMatthew Knepley                MATOP_CHOLESKY_FACTOR_NUMERIC=28,
1489d519adbfSMatthew Knepley                MATOP_SETUP_PREALLOCATION=29,
1490d519adbfSMatthew Knepley                MATOP_ILUFACTOR_SYMBOLIC=30,
1491d519adbfSMatthew Knepley                MATOP_ICCFACTOR_SYMBOLIC=31,
1492d519adbfSMatthew Knepley                MATOP_GET_ARRAY=32,
1493d519adbfSMatthew Knepley                MATOP_RESTORE_ARRAY=33,
1494d519adbfSMatthew Knepley                MATOP_DUPLICATE=34,
1495d519adbfSMatthew Knepley                MATOP_FORWARD_SOLVE=35,
1496d519adbfSMatthew Knepley                MATOP_BACKWARD_SOLVE=36,
1497d519adbfSMatthew Knepley                MATOP_ILUFACTOR=37,
1498d519adbfSMatthew Knepley                MATOP_ICCFACTOR=38,
1499d519adbfSMatthew Knepley                MATOP_AXPY=39,
1500d519adbfSMatthew Knepley                MATOP_GET_SUBMATRICES=40,
1501d519adbfSMatthew Knepley                MATOP_INCREASE_OVERLAP=41,
1502d519adbfSMatthew Knepley                MATOP_GET_VALUES=42,
1503d519adbfSMatthew Knepley                MATOP_COPY=43,
1504d519adbfSMatthew Knepley                MATOP_GET_ROW_MAX=44,
1505d519adbfSMatthew Knepley                MATOP_SCALE=45,
1506d519adbfSMatthew Knepley                MATOP_SHIFT=46,
150735153367SBarry Smith                MATOP_DIAGONAL_SET=47,
1508d519adbfSMatthew Knepley                MATOP_ILUDT_FACTOR=48,
1509d519adbfSMatthew Knepley                MATOP_SET_BLOCK_SIZE=49,
1510d519adbfSMatthew Knepley                MATOP_GET_ROW_IJ=50,
1511d519adbfSMatthew Knepley                MATOP_RESTORE_ROW_IJ=51,
1512d519adbfSMatthew Knepley                MATOP_GET_COLUMN_IJ=52,
1513d519adbfSMatthew Knepley                MATOP_RESTORE_COLUMN_IJ=53,
1514d519adbfSMatthew Knepley                MATOP_FDCOLORING_CREATE=54,
1515d519adbfSMatthew Knepley                MATOP_COLORING_PATCH=55,
1516d519adbfSMatthew Knepley                MATOP_SET_UNFACTORED=56,
1517d519adbfSMatthew Knepley                MATOP_PERMUTE=57,
1518d519adbfSMatthew Knepley                MATOP_SET_VALUES_BLOCKED=58,
1519d519adbfSMatthew Knepley                MATOP_GET_SUBMATRIX=59,
1520d519adbfSMatthew Knepley                MATOP_DESTROY=60,
1521d519adbfSMatthew Knepley                MATOP_VIEW=61,
1522d519adbfSMatthew Knepley                MATOP_CONVERT_FROM=62,
1523d519adbfSMatthew Knepley                MATOP_USE_SCALED_FORM=63,
1524d519adbfSMatthew Knepley                MATOP_SCALE_SYSTEM=64,
1525d519adbfSMatthew Knepley                MATOP_UNSCALE_SYSTEM=65,
1526d519adbfSMatthew Knepley                MATOP_SET_LOCAL_TO_GLOBAL_MAP=66,
1527d519adbfSMatthew Knepley                MATOP_SET_VALUES_LOCAL=67,
1528d519adbfSMatthew Knepley                MATOP_ZERO_ROWS_LOCAL=68,
1529d519adbfSMatthew Knepley                MATOP_GET_ROW_MAX_ABS=69,
1530d519adbfSMatthew Knepley                MATOP_GET_ROW_MIN_ABS=70,
1531d519adbfSMatthew Knepley                MATOP_CONVERT=71,
1532d519adbfSMatthew Knepley                MATOP_SET_COLORING=72,
1533d519adbfSMatthew Knepley                MATOP_SET_VALUES_ADIC=73,
1534d519adbfSMatthew Knepley                MATOP_SET_VALUES_ADIFOR=74,
1535d519adbfSMatthew Knepley                MATOP_FD_COLORING_APPLY=75,
1536d519adbfSMatthew Knepley                MATOP_SET_FROM_OPTIONS=76,
1537d519adbfSMatthew Knepley                MATOP_MULT_CON=77,
1538d519adbfSMatthew Knepley                MATOP_MULT_TRANSPOSE_CON=78,
1539d519adbfSMatthew Knepley                MATOP_PERMUTE_SPARSIFY=79,
1540d519adbfSMatthew Knepley                MATOP_MULT_MULTIPLE=80,
1541d519adbfSMatthew Knepley                MATOP_SOLVE_MULTIPLE=81,
1542d519adbfSMatthew Knepley                MATOP_GET_INERTIA=82,
1543d519adbfSMatthew Knepley                MATOP_LOAD=83,
1544d519adbfSMatthew Knepley                MATOP_IS_SYMMETRIC=84,
1545d519adbfSMatthew Knepley                MATOP_IS_HERMITIAN=85,
1546d519adbfSMatthew Knepley                MATOP_IS_STRUCTURALLY_SYMMETRIC=86,
154741f059aeSBarry Smith                MATOP_DUMMY=87,
1548d519adbfSMatthew Knepley                MATOP_GET_VECS=88,
1549d519adbfSMatthew Knepley                MATOP_MAT_MULT=89,
1550d519adbfSMatthew Knepley                MATOP_MAT_MULT_SYMBOLIC=90,
1551d519adbfSMatthew Knepley                MATOP_MAT_MULT_NUMERIC=91,
1552d519adbfSMatthew Knepley                MATOP_PTAP=92,
1553d519adbfSMatthew Knepley                MATOP_PTAP_SYMBOLIC=93,
1554d519adbfSMatthew Knepley                MATOP_PTAP_NUMERIC=94,
1555d519adbfSMatthew Knepley                MATOP_MAT_MULTTRANSPOSE=95,
15561763ac4eSSatish Balay                MATOP_MAT_MULTTRANSPOSE_SYM=96,
15571763ac4eSSatish Balay                MATOP_MAT_MULTTRANSPOSE_NUM=97,
1558d519adbfSMatthew Knepley                MATOP_PTAP_SYMBOLIC_SEQAIJ=98,
1559d519adbfSMatthew Knepley                MATOP_PTAP_NUMERIC_SEQAIJ=99,
1560d519adbfSMatthew Knepley                MATOP_PTAP_SYMBOLIC_MPIAIJ=100,
1561d519adbfSMatthew Knepley                MATOP_PTAP_NUMERIC_MPIAIJ=101,
1562d519adbfSMatthew Knepley                MATOP_CONJUGATE=102,
1563d519adbfSMatthew Knepley                MATOP_SET_SIZES=103,
1564d519adbfSMatthew Knepley                MATOP_SET_VALUES_ROW=104,
1565d519adbfSMatthew Knepley                MATOP_REAL_PART=105,
1566d519adbfSMatthew Knepley                MATOP_IMAG_PART=106,
1567d519adbfSMatthew Knepley                MATOP_GET_ROW_UTRIANGULAR=107,
1568d519adbfSMatthew Knepley                MATOP_RESTORE_ROW_UTRIANGULAR=108,
1569d519adbfSMatthew Knepley                MATOP_MATSOLVE=109,
1570d519adbfSMatthew Knepley                MATOP_GET_REDUNDANTMATRIX=110,
1571d519adbfSMatthew Knepley                MATOP_GET_ROW_MIN=111,
1572d519adbfSMatthew Knepley                MATOP_GET_COLUMN_VEC=112,
1573d519adbfSMatthew Knepley                MATOP_MISSING_DIAGONAL=113,
1574d519adbfSMatthew Knepley                MATOP_MATGETSEQNONZEROSTRUCTURE=114,
157589c6957cSBarry Smith                MATOP_CREATE=115,
157689c6957cSBarry Smith                MATOP_GET_GHOSTS=116,
1577ea38565cSJed Brown                MATOP_GET_LOCALSUBMATRIX=117,
1578ea38565cSJed Brown                MATOP_RESTORE_LOCALSUBMATRIX=118,
1579eeffb40dSHong Zhang                MATOP_MULT_DIAGONAL_BLOCK=119,
1580547795f9SHong Zhang                MATOP_HERMITIANTRANSPOSE=120,
1581547795f9SHong Zhang                MATOP_MULTHERMITIANTRANSPOSE=121,
1582fbdbba38SShri Abhyankar                MATOP_MULTHERMITIANTRANSPOSEADD=122,
1583b9614d88SDmitry Karpeev                MATOP_GETMULTIPROCBLOCK=123,
1584b9614d88SDmitry Karpeev 	       MATOP_GET_SUBMATRICES_PARALLEL=128
1585fae171e0SBarry Smith              } MatOperation;
15867087cfbeSBarry Smith extern PetscErrorCode  MatHasOperation(Mat,MatOperation,PetscBool *);
15877087cfbeSBarry Smith extern PetscErrorCode  MatShellSetOperation(Mat,MatOperation,void(*)(void));
15887087cfbeSBarry Smith extern PetscErrorCode  MatShellGetOperation(Mat,MatOperation,void(**)(void));
15897087cfbeSBarry Smith extern PetscErrorCode  MatShellSetContext(Mat,void*);
1590112a2221SBarry Smith 
159190ace30eSBarry Smith /*
159290ace30eSBarry Smith    Codes for matrices stored on disk. By default they are
159390ace30eSBarry Smith    stored in a universal format. By changing the format with
15947973ad04SBarry Smith    PetscViewerSetFormat(viewer,PETSC_VIEWER_NATIVE); the matrices will
159590ace30eSBarry Smith    be stored in a way natural for the matrix, for example dense matrices
159690ace30eSBarry Smith    would be stored as dense. Matrices stored this way may only be
1597f4403165SShri Abhyankar    read into matrices of the same type.
159890ace30eSBarry Smith */
159990ace30eSBarry Smith #define MATRIX_BINARY_FORMAT_DENSE -1
160090ace30eSBarry Smith 
16017087cfbeSBarry Smith extern PetscErrorCode  MatMPIBAIJSetHashTableFactor(Mat,PetscReal);
16027087cfbeSBarry Smith extern PetscErrorCode  MatISGetLocalMat(Mat,Mat*);
16031f4e1ec7SBarry Smith 
1604d9274352SBarry Smith /*S
1605d9274352SBarry Smith      MatNullSpace - Object that removes a null space from a vector, i.e.
1606d9274352SBarry Smith          orthogonalizes the vector to a subsapce
1607d9274352SBarry Smith 
1608f7a9e4ceSBarry Smith    Level: advanced
1609d9274352SBarry Smith 
1610d9274352SBarry Smith   Concepts: matrix; linear operator, null space
1611d9274352SBarry Smith 
16126e1639daSBarry Smith   Users manual sections:
16136e1639daSBarry Smith .   sec_singular
16146e1639daSBarry Smith 
1615d9274352SBarry Smith .seealso:  MatNullSpaceCreate()
1616d9274352SBarry Smith S*/
161774637425SBarry Smith typedef struct _p_MatNullSpace* MatNullSpace;
1618d9274352SBarry Smith 
16197087cfbeSBarry Smith extern PetscErrorCode  MatNullSpaceCreate(MPI_Comm,PetscBool ,PetscInt,const Vec[],MatNullSpace*);
16207087cfbeSBarry Smith extern PetscErrorCode  MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(MatNullSpace,Vec,void*),void*);
16217087cfbeSBarry Smith extern PetscErrorCode  MatNullSpaceDestroy(MatNullSpace);
16227087cfbeSBarry Smith extern PetscErrorCode  MatNullSpaceRemove(MatNullSpace,Vec,Vec*);
16237087cfbeSBarry Smith extern PetscErrorCode  MatNullSpaceAttach(Mat,MatNullSpace);
16247087cfbeSBarry Smith extern PetscErrorCode  MatNullSpaceTest(MatNullSpace,Mat,PetscBool  *);
162574637425SBarry Smith 
16267087cfbeSBarry Smith extern PetscErrorCode  MatReorderingSeqSBAIJ(Mat,IS);
16277087cfbeSBarry Smith extern PetscErrorCode  MatMPISBAIJSetHashTableFactor(Mat,PetscReal);
16287087cfbeSBarry Smith extern PetscErrorCode  MatSeqSBAIJSetColumnIndices(Mat,PetscInt *);
16297087cfbeSBarry Smith extern PetscErrorCode  MatSeqBAIJInvertBlockDiagonal(Mat);
16303f1d51d7SBarry Smith 
16317087cfbeSBarry Smith extern PetscErrorCode  MatCreateMAIJ(Mat,PetscInt,Mat*);
16327087cfbeSBarry Smith extern PetscErrorCode  MatMAIJRedimension(Mat,PetscInt,Mat*);
16337087cfbeSBarry Smith extern PetscErrorCode  MatMAIJGetAIJ(Mat,Mat*);
1634c4f061fbSSatish Balay 
16357087cfbeSBarry Smith extern PetscErrorCode  MatComputeExplicitOperator(Mat,Mat*);
1636b0a32e0cSBarry Smith 
16377087cfbeSBarry Smith extern PetscErrorCode  MatDiagonalScaleLocal(Mat,Vec);
163804f1ad80SBarry Smith 
16397087cfbeSBarry Smith extern PetscErrorCode  MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*);
16407087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetBase(Mat,Vec,Vec);
16417087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*);
16427087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*));
16437087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec));
16447087cfbeSBarry Smith extern PetscErrorCode  MatMFFDAddNullSpace(Mat,MatNullSpace);
16457087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt);
16467087cfbeSBarry Smith extern PetscErrorCode  MatMFFDResetHHistory(Mat);
16477087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetFunctionError(Mat,PetscReal);
16487087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetPeriod(Mat,PetscInt);
16497087cfbeSBarry Smith extern PetscErrorCode  MatMFFDGetH(Mat,PetscScalar *);
16507087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetOptionsPrefix(Mat,const char[]);
16517087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetFromOptions(Mat);
16527087cfbeSBarry Smith extern PetscErrorCode  MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*);
16537087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*);
1654e884886fSBarry Smith 
16556370053bSBarry Smith /*S
16566370053bSBarry Smith     MatMFFD - A data structured used to manage the computation of the h differencing parameter for matrix-free
16576370053bSBarry Smith               Jacobian vector products
1658e884886fSBarry Smith 
16596370053bSBarry Smith     Notes: MATMFFD is a specific MatType which uses the MatMFFD data structure
16606370053bSBarry Smith 
16616370053bSBarry Smith            MatMFFD*() methods actually take the Mat as their first argument. Not a MatMFFD data structure
16626370053bSBarry Smith 
16636370053bSBarry Smith     Level: developer
16646370053bSBarry Smith 
16656370053bSBarry Smith .seealso: MATMFFD, MatCreateMFFD(), MatMFFDSetFuction(), MatMFFDSetType(), MatMFFDRegister()
16666370053bSBarry Smith S*/
1667e884886fSBarry Smith typedef struct _p_MatMFFD* MatMFFD;
1668e884886fSBarry Smith 
1669e884886fSBarry Smith /*E
1670e884886fSBarry Smith     MatMFFDType - algorithm used to compute the h used in computing matrix-vector products via differencing of the function
1671e884886fSBarry Smith 
1672e884886fSBarry Smith    Level: beginner
1673e884886fSBarry Smith 
1674e884886fSBarry Smith .seealso: MatMFFDSetType(), MatMFFDRegister()
1675e884886fSBarry Smith E*/
1676a313700dSBarry Smith #define MatMFFDType char*
1677e884886fSBarry Smith #define MATMFFD_DS  "ds"
1678e884886fSBarry Smith #define MATMFFD_WP  "wp"
1679e884886fSBarry Smith 
16807087cfbeSBarry Smith extern PetscErrorCode  MatMFFDSetType(Mat,const MatMFFDType);
16817087cfbeSBarry Smith extern PetscErrorCode  MatMFFDRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatMFFD));
1682e884886fSBarry Smith 
1683e884886fSBarry Smith /*MC
1684e884886fSBarry Smith    MatMFFDRegisterDynamic - Adds a method to the MatMFFD registry.
1685e884886fSBarry Smith 
1686e884886fSBarry Smith    Synopsis:
16871890ba74SBarry Smith    PetscErrorCode MatMFFDRegisterDynamic(const char *name_solver,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatMFFD))
1688e884886fSBarry Smith 
1689e884886fSBarry Smith    Not Collective
1690e884886fSBarry Smith 
1691e884886fSBarry Smith    Input Parameters:
1692e884886fSBarry Smith +  name_solver - name of a new user-defined compute-h module
1693e884886fSBarry Smith .  path - path (either absolute or relative) the library containing this solver
1694e884886fSBarry Smith .  name_create - name of routine to create method context
1695e884886fSBarry Smith -  routine_create - routine to create method context
1696e884886fSBarry Smith 
1697e884886fSBarry Smith    Level: developer
1698e884886fSBarry Smith 
1699e884886fSBarry Smith    Notes:
1700e884886fSBarry Smith    MatMFFDRegisterDynamic() may be called multiple times to add several user-defined solvers.
1701e884886fSBarry Smith 
1702e884886fSBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
1703e884886fSBarry Smith    is ignored.
1704e884886fSBarry Smith 
1705e884886fSBarry Smith    Sample usage:
1706e884886fSBarry Smith .vb
1707e884886fSBarry Smith    MatMFFDRegisterDynamic("my_h",/home/username/my_lib/lib/libO/solaris/mylib.a,
1708e884886fSBarry Smith                "MyHCreate",MyHCreate);
1709e884886fSBarry Smith .ve
1710e884886fSBarry Smith 
1711e884886fSBarry Smith    Then, your solver can be chosen with the procedural interface via
1712e884886fSBarry Smith $     MatMFFDSetType(mfctx,"my_h")
1713e884886fSBarry Smith    or at runtime via the option
1714e884886fSBarry Smith $     -snes_mf_type my_h
1715e884886fSBarry Smith 
1716e884886fSBarry Smith .keywords: MatMFFD, register
1717e884886fSBarry Smith 
1718e884886fSBarry Smith .seealso: MatMFFDRegisterAll(), MatMFFDRegisterDestroy()
1719e884886fSBarry Smith M*/
1720e884886fSBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1721e884886fSBarry Smith #define MatMFFDRegisterDynamic(a,b,c,d) MatMFFDRegister(a,b,c,0)
1722e884886fSBarry Smith #else
1723e884886fSBarry Smith #define MatMFFDRegisterDynamic(a,b,c,d) MatMFFDRegister(a,b,c,d)
1724e884886fSBarry Smith #endif
1725e884886fSBarry Smith 
17267087cfbeSBarry Smith extern PetscErrorCode  MatMFFDRegisterAll(const char[]);
17277087cfbeSBarry Smith extern PetscErrorCode  MatMFFDRegisterDestroy(void);
17287087cfbeSBarry Smith extern PetscErrorCode  MatMFFDDSSetUmin(Mat,PetscReal);
17297087cfbeSBarry Smith extern PetscErrorCode  MatMFFDWPSetComputeNormU(Mat,PetscBool );
1730e884886fSBarry Smith 
1731e884886fSBarry Smith 
17327087cfbeSBarry Smith extern PetscErrorCode  PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *);
17337087cfbeSBarry Smith extern PetscErrorCode  PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *);
17347dbadf16SMatthew Knepley 
173597969023SHong Zhang /*
173697969023SHong Zhang    PETSc interface to MUMPS
173797969023SHong Zhang */
173897969023SHong Zhang #ifdef PETSC_HAVE_MUMPS
17397087cfbeSBarry Smith extern PetscErrorCode  MatSetMumpsIcntl(Mat,PetscInt,PetscInt);
174097969023SHong Zhang #endif
174123a5497aSJed Brown 
17427087cfbeSBarry Smith extern PetscErrorCode  MatCreateNest(MPI_Comm comm,PetscInt,const IS[],PetscInt,const IS[],const Mat[],Mat*);
17437087cfbeSBarry Smith extern PetscErrorCode  MatNestGetSize(Mat A,PetscInt *M,PetscInt *N);
17447087cfbeSBarry Smith extern PetscErrorCode  MatNestGetSubMats(Mat A,PetscInt *M,PetscInt *N,Mat ***mat);
17457087cfbeSBarry Smith extern PetscErrorCode  MatNestGetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat *sub);
17467087cfbeSBarry Smith extern PetscErrorCode  MatNestSetVecType(Mat,const VecType);
1747d8588912SDave May 
174823a5497aSJed Brown PETSC_EXTERN_CXX_END
174923a5497aSJed Brown #endif
1750