xref: /petsc/include/petscmat.h (revision a714c06d2e1c14eac20fc9e0765d1d6303f28782)
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"
31273d9f13SBarry Smith #define MATSEQMAIJ         "seqmaij"
32273d9f13SBarry Smith #define MATMPIMAIJ         "mpimaij"
33209238afSKris Buschelman #define MATMAIJ            "maij"
34273d9f13SBarry Smith #define MATIS              "is"
35273d9f13SBarry Smith #define MATSEQAIJ          "seqaij"
36273d9f13SBarry Smith #define MATMPIAIJ          "mpiaij"
37209238afSKris Buschelman #define MATAIJ             "aij"
38273d9f13SBarry Smith #define MATSHELL           "shell"
39209238afSKris Buschelman #define MATSEQDENSE        "seqdense"
40273d9f13SBarry Smith #define MATMPIDENSE        "mpidense"
41209238afSKris Buschelman #define MATDENSE           "dense"
42273d9f13SBarry Smith #define MATSEQBAIJ         "seqbaij"
43273d9f13SBarry Smith #define MATMPIBAIJ         "mpibaij"
44209238afSKris Buschelman #define MATBAIJ            "baij"
45273d9f13SBarry Smith #define MATMPIADJ          "mpiadj"
46273d9f13SBarry Smith #define MATSEQSBAIJ        "seqsbaij"
47273d9f13SBarry Smith #define MATMPISBAIJ        "mpisbaij"
48209238afSKris Buschelman #define MATSBAIJ           "sbaij"
49cebc7f6cSBarry Smith #define MATDAAD            "daad"
50cebc7f6cSBarry Smith #define MATMFFD            "mffd"
51c8a8475eSBarry Smith #define MATNORMAL          "normal"
52ab92ecdeSBarry Smith #define MATLRC             "lrc"
534b8d99adSRichard Tran Mills #define MATSEQCSRPERM      "seqcsrperm"
5418def556SRichard Tran Mills #define MATMPICSRPERM      "mpicsrperm"
55814655b8SBarry Smith #define MATCSRPERM         "csrperm"
5681824310SBarry Smith #define MATSEQCRL          "seqcrl"
57c02b24c5SBarry Smith #define MATMPICRL          "mpicrl"
58c02b24c5SBarry Smith #define MATCRL             "crl"
592a6744ebSBarry Smith #define MATSCATTER         "scatter"
60421e10b8SBarry Smith #define MATBLOCKMAT        "blockmat"
61793850ffSBarry Smith #define MATCOMPOSITE       "composite"
625a7f1df3SHong Zhang #define MATSEQFFTW         "seqfftw"
63557cca28SSatish Balay #define MATTRANSPOSEMAT    "transpose"
6472ca8751SBarry Smith #define MATSCHURCOMPLEMENT "schurcomplement"
651d6018f0SLisandro Dalcin #define MATPYTHON          "python"
66f91d8e95SBarry Smith #define MATHYPRESTRUCT     "hyprestruct"
67ee1cef2cSJed Brown #define MATSUBMATRIX       "submatrix"
68d91e6319SBarry Smith 
69a6a5cd3fSDmitry Karpeev #ifdef PETSC_USE_MATFRAMEWORK
70a6a5cd3fSDmitry Karpeev #define MATFRAMEWORK           "framework"
71a6a5cd3fSDmitry Karpeev #endif
72773941d6SBarry Smith 
739989ab13SBarry Smith /*E
74c7393fdbSBarry Smith     MatSolverPackage - String with the name of a PETSc matrix solver type.
759989ab13SBarry Smith 
769989ab13SBarry Smith     For example: "petsc" indicates what PETSc provides, "superlu" indicates either
779989ab13SBarry Smith        SuperLU or SuperLU_Dist etc.
789989ab13SBarry Smith 
799989ab13SBarry Smith 
809989ab13SBarry Smith    Level: beginner
819989ab13SBarry Smith 
825c9eb25fSBarry Smith .seealso: MatGetFactor(), Mat, MatSetType(), MatType
839989ab13SBarry Smith E*/
84c7393fdbSBarry Smith #define MatSolverPackage char*
85b24902e0SBarry Smith #define MAT_SOLVER_SPOOLES      "spooles"
86b24902e0SBarry Smith #define MAT_SOLVER_SUPERLU      "superlu"
87b24902e0SBarry Smith #define MAT_SOLVER_SUPERLU_DIST "superlu_dist"
88b24902e0SBarry Smith #define MAT_SOLVER_UMFPACK      "umfpack"
89b24902e0SBarry Smith #define MAT_SOLVER_ESSL         "essl"
90b24902e0SBarry Smith #define MAT_SOLVER_LUSOL        "lusol"
91b24902e0SBarry Smith #define MAT_SOLVER_MUMPS        "mumps"
923bf14a46SMatthew Knepley #define MAT_SOLVER_PASTIX       "pastix"
93b24902e0SBarry Smith #define MAT_SOLVER_DSCPACK      "dscpack"
94b24902e0SBarry Smith #define MAT_SOLVER_MATLAB       "matlab"
9543244d56SBarry Smith #define MAT_SOLVER_PETSC        "petsc"
96b6806ab0SHong Zhang #define MAT_SOLVER_PLAPACK      "plapack"
97b24902e0SBarry Smith 
98773941d6SBarry Smith 
99b24902e0SBarry Smith /*E
100b24902e0SBarry Smith     MatFactorType - indicates what type of factorization is requested
101b24902e0SBarry Smith 
102b24902e0SBarry Smith     Level: beginner
103b24902e0SBarry Smith 
104b24902e0SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
105b24902e0SBarry Smith 
106c7393fdbSBarry Smith .seealso: MatSolverPackage, MatGetFactor()
107b24902e0SBarry Smith E*/
108599ef60dSHong Zhang typedef enum {MAT_FACTOR_NONE, MAT_FACTOR_LU, MAT_FACTOR_CHOLESKY, MAT_FACTOR_ILU, MAT_FACTOR_ICC,MAT_FACTOR_ILUDT} MatFactorType;
109c7393fdbSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetFactor(Mat,const MatSolverPackage,MatFactorType,Mat*);
110db4efbfdSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetFactorAvailable(Mat,const MatSolverPackage,MatFactorType,PetscTruth*);
11135bd34faSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFactorGetSolverPackage(Mat,const MatSolverPackage*);
112b24902e0SBarry Smith 
1139989ab13SBarry Smith 
114c06d978dSMatthew Knepley /* Logging support */
115552e946dSBarry Smith #define    MAT_FILE_COOKIE 1211216    /* used to indicate matrices in binary files */
116be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_COOKIE;
117be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_FDCOLORING_COOKIE;
118be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_PARTITIONING_COOKIE;
119be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_NULLSPACE_COOKIE;
120e884886fSBarry Smith extern PetscCookie PETSCMAT_DLLEXPORT MATMFFD_COOKIE;
121c06d978dSMatthew Knepley 
122ceb03754SKris Buschelman /*E
123ceb03754SKris Buschelman     MatReuse - Indicates if matrices obtained from a previous call to MatGetSubMatrices()
124ceb03754SKris Buschelman      or MatGetSubMatrix() are to be reused to store the new matrix values.
125ceb03754SKris Buschelman 
126ceb03754SKris Buschelman     Level: beginner
127ceb03754SKris Buschelman 
128ceb03754SKris Buschelman    Any additions/changes here MUST also be made in include/finclude/petscmat.h
129ceb03754SKris Buschelman 
1300c451bc4SBarry Smith .seealso: MatGetSubMatrices(), MatGetSubMatrix(), MatDestroyMatrices(), MatConvert()
131ceb03754SKris Buschelman E*/
132dfe085dbSJed Brown typedef enum {MAT_INITIAL_MATRIX,MAT_REUSE_MATRIX,MAT_IGNORE_MATRIX} MatReuse;
133ceb03754SKris Buschelman 
1345494a064SHong Zhang /*E
1355494a064SHong Zhang     MatGetSubMatrixOption - Indicates if matrices obtained from a call to MatGetSubMatrices()
136829201f2SHong Zhang      include the matrix values. Currently it is only used by MatGetSeqNonzerostructure().
1375494a064SHong Zhang 
1385494a064SHong Zhang     Level: beginner
1395494a064SHong Zhang 
140829201f2SHong Zhang .seealso: MatGetSeqNonzerostructure()
1415494a064SHong Zhang E*/
1425494a064SHong Zhang typedef enum {MAT_DO_NOT_GET_VALUES,MAT_GET_VALUES} MatGetSubMatrixOption;
1435494a064SHong Zhang 
144e5bd5246SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInitializePackage(const char[]);
145c06d978dSMatthew Knepley 
146f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreate(MPI_Comm,Mat*);
147a69a960dSMatthew Knepley PetscPolymorphicFunction(MatCreate,(MPI_Comm comm),(comm,&A),Mat,A)
148a69a960dSMatthew Knepley PetscPolymorphicFunction(MatCreate,(),(PETSC_COMM_WORLD,&A),Mat,A)
149f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetSizes(Mat,PetscInt,PetscInt,PetscInt,PetscInt);
150a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetType(Mat,const MatType);
151be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetFromOptions(Mat);
152be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUpPreallocation(Mat);
153be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRegisterAll(const char[]);
154be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat));
15530de9b25SBarry Smith 
156f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetOptionsPrefix(Mat,const char[]);
157f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAppendOptionsPrefix(Mat,const char[]);
158f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOptionsPrefix(Mat,const char*[]);
159f69a0ea3SMatthew Knepley 
16030de9b25SBarry Smith /*MC
16130de9b25SBarry Smith    MatRegisterDynamic - Adds a new matrix type
16230de9b25SBarry Smith 
16330de9b25SBarry Smith    Synopsis:
164c1ac3661SBarry Smith    PetscErrorCode MatRegisterDynamic(char *name,char *path,char *name_create,PetscErrorCode (*routine_create)(Mat))
16530de9b25SBarry Smith 
16630de9b25SBarry Smith    Not Collective
16730de9b25SBarry Smith 
16830de9b25SBarry Smith    Input Parameters:
16930de9b25SBarry Smith +  name - name of a new user-defined matrix type
17030de9b25SBarry Smith .  path - path (either absolute or relative) the library containing this solver
17130de9b25SBarry Smith .  name_create - name of routine to create method context
17230de9b25SBarry Smith -  routine_create - routine to create method context
17330de9b25SBarry Smith 
17430de9b25SBarry Smith    Notes:
17530de9b25SBarry Smith    MatRegisterDynamic() may be called multiple times to add several user-defined solvers.
17630de9b25SBarry Smith 
17730de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
17830de9b25SBarry Smith    is ignored.
17930de9b25SBarry Smith 
18030de9b25SBarry Smith    Sample usage:
18130de9b25SBarry Smith .vb
18230de9b25SBarry Smith    MatRegisterDynamic("my_mat",/home/username/my_lib/lib/libO/solaris/mylib.a,
18330de9b25SBarry Smith                "MyMatCreate",MyMatCreate);
18430de9b25SBarry Smith .ve
18530de9b25SBarry Smith 
18630de9b25SBarry Smith    Then, your solver can be chosen with the procedural interface via
18730de9b25SBarry Smith $     MatSetType(Mat,"my_mat")
18830de9b25SBarry Smith    or at runtime via the option
18930de9b25SBarry Smith $     -mat_type my_mat
19030de9b25SBarry Smith 
19130de9b25SBarry Smith    Level: advanced
19230de9b25SBarry Smith 
193ab901514SBarry Smith    Notes: ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values.
19430de9b25SBarry Smith          If your function is not being put into a shared library then use VecRegister() instead
19530de9b25SBarry Smith 
19630de9b25SBarry Smith .keywords: Mat, register
19730de9b25SBarry Smith 
19830de9b25SBarry Smith .seealso: MatRegisterAll(), MatRegisterDestroy()
19930de9b25SBarry Smith 
20030de9b25SBarry Smith M*/
201273d9f13SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
202273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,0)
203273d9f13SBarry Smith #else
204273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,d)
20530de9b25SBarry Smith #endif
20630de9b25SBarry Smith 
207273d9f13SBarry Smith extern PetscTruth MatRegisterAllCalled;
208b0a32e0cSBarry Smith extern PetscFList MatList;
209b022a5c1SBarry Smith extern PetscFList MatColoringList;
210b022a5c1SBarry Smith extern PetscFList MatPartitioningList;
21128988994SBarry Smith 
2123b224e63SBarry Smith /*E
2133b224e63SBarry Smith     MatStructure - Indicates if the matrix has the same nonzero structure
2143b224e63SBarry Smith 
2153b224e63SBarry Smith     Level: beginner
2163b224e63SBarry Smith 
2173b224e63SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
2183b224e63SBarry Smith 
2193b224e63SBarry Smith .seealso: MatCopy(), KSPSetOperators(), PCSetOperators()
2203b224e63SBarry Smith E*/
2213b224e63SBarry Smith typedef enum {SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER,SUBSET_NONZERO_PATTERN} MatStructure;
2223b224e63SBarry Smith 
223be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqDense(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*);
224be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIDense(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*);
225be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
226ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,nz,nnz,&A),Mat,A)
227ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,&A),Mat,A)
228ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,0,nnz,&A),Mat,A)
229ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,&A),Mat,A)
230ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,A))
231ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,m,n,0,nnz,A))
232ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,A))
233be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
234ba966639SSatish 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)
235ba966639SSatish 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)
236ba966639SSatish 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)
237ba966639SSatish 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)
238ba966639SSatish 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))
239ba966639SSatish 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))
240ba966639SSatish 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))
241ba966639SSatish 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)
242ba966639SSatish 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)
243ba966639SSatish 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)
244ba966639SSatish 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)
245ba966639SSatish 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))
246ba966639SSatish 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))
247ba966639SSatish 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))
24863ab93bfSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *);
2498d7a6e47SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJWithSplitArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],PetscInt[],PetscInt[],PetscScalar[],Mat*);
2508d7a6e47SBarry Smith 
251be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
252ba966639SSatish 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)
253ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A)
254ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A)
255ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A)
256ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A))
257ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A))
258ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A))
259be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
260ba966639SSatish 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)
261ba966639SSatish 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)
262ba966639SSatish 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)
263ba966639SSatish 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)
264ba966639SSatish 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))
265ba966639SSatish 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))
266ba966639SSatish 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))
267ba966639SSatish 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)
268ba966639SSatish 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)
269ba966639SSatish 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)
270ba966639SSatish 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)
271ba966639SSatish 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))
272ba966639SSatish 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))
273ba966639SSatish 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))
274be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAdj(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscInt[],Mat*);
275be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
276ba966639SSatish 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)
277ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A)
278ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A)
279ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A)
280ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A))
281ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A))
282ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A))
283be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPISBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
284ba966639SSatish 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)
285ba966639SSatish 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)
286ba966639SSatish 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)
287ba966639SSatish 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)
288ba966639SSatish 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))
289ba966639SSatish 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))
290ba966639SSatish 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))
291ba966639SSatish 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)
292ba966639SSatish 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)
293ba966639SSatish 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)
294ba966639SSatish 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)
295ba966639SSatish 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))
296ba966639SSatish 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))
297ba966639SSatish 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))
298be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateShell(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,void *,Mat*);
299ba966639SSatish 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)
300ba966639SSatish Balay PetscPolymorphicFunction(MatCreateShell,(PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx),(PETSC_COMM_WORLD,m,n,M,N,ctx,&A),Mat,A)
301be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateAdic(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,void (*)(void),Mat*);
302be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateNormal(Mat,Mat*);
303ba966639SSatish Balay PetscPolymorphicFunction(MatCreateNormal,(Mat mat),(mat,&A),Mat,A)
304ab92ecdeSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateLRC(Mat,Mat,Mat,Mat*);
305284134d9SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateIS(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,ISLocalToGlobalMapping,Mat*);
3061472f72bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
3071472f72bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPICRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
3082a6744ebSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateScatter(MPI_Comm,VecScatter,Mat*);
3092a6744ebSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScatterSetVecScatter(Mat,VecScatter);
3102a6744ebSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScatterGetVecScatter(Mat,VecScatter*);
3118cdf2d9bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateBlockMat(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,Mat*);
312793850ffSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompositeAddMat(Mat,Mat);
313b52f573bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompositeMerge(Mat);
314793850ffSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateComposite(MPI_Comm,PetscInt,const Mat*,Mat*);
3156d7c1e57SBarry Smith typedef enum {MAT_COMPOSITE_ADDITIVE,MAT_COMPOSITE_MULTIPLICATIVE} MatCompositeType;
3166d7c1e57SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompositeSetType(Mat,MatCompositeType);
3176d7c1e57SBarry Smith 
318a6a5cd3fSDmitry Karpeev #ifdef PETSC_USE_MATFRAMEWORK
319a6a5cd3fSDmitry Karpeev typedef enum {MAT_FRAMEWORK_ADDITIVE, MAT_FRAMEWORK_MULTIPLICATIVE} MatFrameworkType;
320a6a5cd3fSDmitry Karpeev typedef enum {MAT_FRAMEWORK_ASSEMBLED, MAT_FRAMEWORK_DISASSEMBLED} MatFrameworkMode;
321a6a5cd3fSDmitry Karpeev 
322a6a5cd3fSDmitry Karpeev EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateFramework(MPI_Comm,Mat *);
323a6a5cd3fSDmitry Karpeev EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFrameworkSetType(Mat,MatFrameworkType);
324a6a5cd3fSDmitry Karpeev EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFrameworkGetType(Mat,MatFrameworkType*);
325a6a5cd3fSDmitry Karpeev EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFrameworkSetMode(Mat,MatFrameworkMode);
326a6a5cd3fSDmitry Karpeev EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFrameworkGetMode(Mat,MatFrameworkMode*);
327a6a5cd3fSDmitry Karpeev EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFrameworkGetBlockCount(Mat, PetscInt*);
328a6a5cd3fSDmitry Karpeev EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFrameworkAddBlock(Mat, PetscInt*);
329a6a5cd3fSDmitry Karpeev EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFrameworkSetBlockEmbedding(Mat, PetscInt, IS,IS);
330a6a5cd3fSDmitry Karpeev EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFrameworkGetBlockEmbedding(Mat, PetscInt, IS, IS);
331a6a5cd3fSDmitry Karpeev EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFrameworkSetBlock(Mat, PetscInt, Mat);
332a6a5cd3fSDmitry Karpeev EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFrameworkGetBlock(Mat, PetscInt, Mat*);
333a6a5cd3fSDmitry Karpeev #endif
334a6a5cd3fSDmitry Karpeev 
3355a7f1df3SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqFFTW(MPI_Comm,PetscInt,const PetscInt[],Mat*);
336f20085c4SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateTranspose(Mat,Mat*);
337ee1cef2cSJed Brown EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSubMatrix(Mat,IS,IS,Mat*);
338ee1cef2cSJed Brown EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSubMatrixUpdate(Mat,Mat,IS,IS);
3391472f72bSBarry Smith 
3401d6018f0SLisandro Dalcin EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreatePython(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const char[],Mat*);
3411d6018f0SLisandro Dalcin EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPythonSetType(Mat,const char[]);
3421d6018f0SLisandro Dalcin 
3431d6018f0SLisandro Dalcin 
344f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUp(Mat);
345be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroy(Mat);
34621c89e3eSBarry Smith 
3470c451bc4SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConjugate(Mat);
34899cafbc1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRealPart(Mat);
34999cafbc1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatImaginaryPart(Mat);
350bfc7d00aSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetDiagonalBlock(Mat,PetscTruth*,MatReuse,Mat*);
35199cafbc1SBarry Smith 
3528ed539a5SBarry Smith /* ------------------------------------------------------------*/
353be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
354be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
35587d4246cSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesRow(Mat,PetscInt,const PetscScalar[]);
356f259bd47SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesRowLocal(Mat,PetscInt,const PetscScalar[]);
35784cb2905SBarry Smith 
3582ef4de8bSBarry Smith /*S
3592ef4de8bSBarry Smith      MatStencil - Data structure (C struct) for storing information about a single row or
3602ef4de8bSBarry Smith         column of a matrix as index on an associated grid.
3612ef4de8bSBarry Smith 
3622ef4de8bSBarry Smith    Level: beginner
3632ef4de8bSBarry Smith 
3642ef4de8bSBarry Smith   Concepts: matrix; linear operator
3652ef4de8bSBarry Smith 
366d7bd93baSBarry Smith .seealso:  MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockStencil()
3672ef4de8bSBarry Smith S*/
368435da068SBarry Smith typedef struct {
369c1ac3661SBarry Smith   PetscInt k,j,i,c;
370435da068SBarry Smith } MatStencil;
3712ef4de8bSBarry Smith 
372be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
373be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
374be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt);
375435da068SBarry Smith 
376be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetColoring(Mat,ISColoring);
377be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesAdic(Mat,void*);
378be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesAdifor(Mat,PetscInt,void*);
3793a7fca6bSBarry Smith 
380d91e6319SBarry Smith /*E
381d91e6319SBarry Smith     MatAssemblyType - Indicates if the matrix is now to be used, or if you plan
382d91e6319SBarry Smith      to continue to add values to it
383d91e6319SBarry Smith 
384d91e6319SBarry Smith     Level: beginner
385d91e6319SBarry Smith 
386d91e6319SBarry Smith .seealso: MatAssemblyBegin(), MatAssemblyEnd()
387d91e6319SBarry Smith E*/
3886d4a8577SBarry Smith typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType;
389be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssemblyBegin(Mat,MatAssemblyType);
390be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssemblyEnd(Mat,MatAssemblyType);
391be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssembled(Mat,PetscTruth*);
3924f9c727eSBarry Smith 
3931d73ed98SBarry Smith 
39430de9b25SBarry Smith 
395d91e6319SBarry Smith /*E
396d91e6319SBarry Smith     MatOption - Options that may be set for a matrix and its behavior or storage
397d91e6319SBarry Smith 
398d91e6319SBarry Smith     Level: beginner
399d91e6319SBarry Smith 
4000a835dfdSSatish Balay    Any additions/changes here MUST also be made in include/finclude/petscmat.h
401d91e6319SBarry Smith 
402d91e6319SBarry Smith .seealso: MatSetOption()
403d91e6319SBarry Smith E*/
404512a5fc5SBarry Smith typedef enum {MAT_ROW_ORIENTED,MAT_NEW_NONZERO_LOCATIONS,
4054e0d8c25SBarry Smith               MAT_SYMMETRIC,
4064e0d8c25SBarry Smith               MAT_STRUCTURALLY_SYMMETRIC,
4078047e77bSBarry Smith               MAT_NEW_DIAGONALS,MAT_IGNORE_OFF_PROC_ENTRIES,
4084e0d8c25SBarry Smith               MAT_NEW_NONZERO_LOCATION_ERR,
4094e0d8c25SBarry Smith               MAT_NEW_NONZERO_ALLOCATION_ERR,MAT_USE_HASH_TABLE,
410a9817697SBarry Smith               MAT_KEEP_NONZERO_PATTERN,MAT_IGNORE_ZERO_ENTRIES,
4114e0d8c25SBarry Smith               MAT_USE_INODES,
4124e0d8c25SBarry Smith               MAT_HERMITIAN,
4134e0d8c25SBarry Smith               MAT_SYMMETRY_ETERNAL,
4144e0d8c25SBarry Smith               MAT_USE_COMPRESSEDROW,
4154e0d8c25SBarry Smith               MAT_IGNORE_LOWER_TRIANGULAR,MAT_ERROR_LOWER_TRIANGULAR,
41628b2fa4aSMatthew Knepley               MAT_GETROW_UPPERTRIANGULAR,MAT_UNUSED_NONZERO_LOCATION_ERR,
41728b2fa4aSMatthew Knepley               NUM_MAT_OPTIONS} MatOption;
418290bbb0aSBarry Smith extern const char *MatOptions[];
4194e0d8c25SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetOption(Mat,MatOption,PetscTruth);
420a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetType(Mat,const MatType*);
421a313700dSBarry Smith PetscPolymorphicFunction(MatGetType,(Mat mat),(mat,&t),const MatType,t)
42284cb2905SBarry Smith 
423be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]);
424be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
425be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
426f5edf698SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowUpperTriangular(Mat);
427f5edf698SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRowUpperTriangular(Mat);
428be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
429be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
430be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumnVector(Mat,Vec,PetscInt);
431be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetArray(Mat,PetscScalar *[]);
432ba966639SSatish Balay PetscPolymorphicFunction(MatGetArray,(Mat mat),(mat,&a),PetscScalar*,a)
433be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreArray(Mat,PetscScalar *[]);
434be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBlockSize(Mat,PetscInt *);
435ba966639SSatish Balay PetscPolymorphicFunction(MatGetBlockSize,(Mat mat),(mat,&a),PetscInt,a)
436be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetBlockSize(Mat,PetscInt);
4377b80b807SBarry Smith 
4381620fd73SBarry Smith 
439be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMult(Mat,Vec,Vec);
44089c6957cSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultDiagonalBlock(Mat,Vec,Vec);
441be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultAdd(Mat,Vec,Vec,Vec);
442ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
443be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTranspose(Mat,Vec,Vec);
444547795f9SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultHermitianTranspose(Mat,Vec,Vec);
445be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsTranspose(Mat,Mat,PetscReal,PetscTruth*);
446ba966639SSatish Balay PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B,PetscReal tol),(A,B,tol,&t),PetscTruth,t)
447e37f64cdSBarry Smith PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B),(A,B,0,&t),PetscTruth,t)
4481cbb95d3SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitianTranspose(Mat,Mat,PetscReal,PetscTruth*);
449be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeAdd(Mat,Vec,Vec,Vec);
4505cb05578SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultHermitianTransposeAdd(Mat,Vec,Vec,Vec);
451ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultTransposeAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
452be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultConstrained(Mat,Vec,Vec);
453be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeConstrained(Mat,Vec,Vec);
4542bebee5dSHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatSolve(Mat,Mat,Mat);
455f5edf698SHong Zhang 
456d91e6319SBarry Smith /*E
457d91e6319SBarry Smith     MatDuplicateOption - Indicates if a duplicated sparse matrix should have
458d91e6319SBarry Smith   its numerical values copied over or just its nonzero structure.
459d91e6319SBarry Smith 
460d91e6319SBarry Smith     Level: beginner
461d91e6319SBarry Smith 
462d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
463d91e6319SBarry Smith 
46470dcbbb9SBarry Smith $   MAT_SHARE_NONZERO_PATTERN - the i and j arrays in the new matrix will be shared with the original matrix
46570dcbbb9SBarry Smith $                               this also triggers the MAT_DO_NOT_COPY_VALUES option. This is used when you
46670dcbbb9SBarry Smith $                               have several matrices with the same nonzero pattern.
46770dcbbb9SBarry Smith 
468d91e6319SBarry Smith .seealso: MatDuplicate()
469d91e6319SBarry Smith E*/
47070dcbbb9SBarry Smith typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES,MAT_SHARE_NONZERO_PATTERN} MatDuplicateOption;
4712e8a6d31SBarry Smith 
472a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert(Mat,const MatType,MatReuse,Mat*);
473a313700dSBarry Smith PetscPolymorphicFunction(MatConvert,(Mat A,const MatType t),(A,t,MAT_INITIAL_MATRIX,&a),Mat,a)
474be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDuplicate(Mat,MatDuplicateOption,Mat*);
475ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A,MatDuplicateOption o),(A,o,&a),Mat,a)
476ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A),(A,MAT_COPY_VALUES,&a),Mat,a)
47794a9d846SBarry Smith 
478cb5b572fSBarry Smith 
479be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCopy(Mat,Mat,MatStructure);
480be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatView(Mat,PetscViewer);
481be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsSymmetric(Mat,PetscReal,PetscTruth*);
482ba966639SSatish Balay PetscPolymorphicFunction(MatIsSymmetric,(Mat A,PetscReal tol),(A,tol,&t),PetscTruth,t)
483e37f64cdSBarry Smith PetscPolymorphicFunction(MatIsSymmetric,(Mat A),(A,0,&t),PetscTruth,t)
484be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsStructurallySymmetric(Mat,PetscTruth*);
485ba966639SSatish Balay PetscPolymorphicFunction(MatIsStructurallySymmetric,(Mat A),(A,&t),PetscTruth,t)
4861cbb95d3SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitian(Mat,PetscReal,PetscTruth*);
4871cbb95d3SBarry Smith PetscPolymorphicFunction(MatIsHermitian,(Mat A),(A,0,&t),PetscTruth,t)
488be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsSymmetricKnown(Mat,PetscTruth*,PetscTruth*);
489be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitianKnown(Mat,PetscTruth*,PetscTruth*);
49064c62002SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMissingDiagonal(Mat,PetscTruth *,PetscInt *);
491a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLoad(PetscViewer,const MatType,Mat*);
492a313700dSBarry Smith PetscPolymorphicFunction(MatLoad,(PetscViewer v,const MatType t),(v,t,&a),Mat,a)
4937b80b807SBarry Smith 
4948f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *);
4958f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRowIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *);
4968f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumnIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *);
4978f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreColumnIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *);
498d4fbbf0eSBarry Smith 
499d91e6319SBarry Smith /*S
500d91e6319SBarry Smith      MatInfo - Context of matrix information, used with MatGetInfo()
501d91e6319SBarry Smith 
502d91e6319SBarry Smith    In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE
503d91e6319SBarry Smith 
504d91e6319SBarry Smith    Level: intermediate
505d91e6319SBarry Smith 
506d91e6319SBarry Smith   Concepts: matrix^nonzero information
507d91e6319SBarry Smith 
508d9274352SBarry Smith .seealso:  MatGetInfo(), MatInfoType
509d91e6319SBarry Smith S*/
5104e220ebcSLois Curfman McInnes typedef struct {
511b0a32e0cSBarry Smith   PetscLogDouble block_size;                         /* block size */
512b0a32e0cSBarry Smith   PetscLogDouble nz_allocated,nz_used,nz_unneeded;   /* number of nonzeros */
513b0a32e0cSBarry Smith   PetscLogDouble memory;                             /* memory allocated */
514b0a32e0cSBarry Smith   PetscLogDouble assemblies;                         /* number of matrix assemblies called */
515b0a32e0cSBarry Smith   PetscLogDouble mallocs;                            /* number of mallocs during MatSetValues() */
516b0a32e0cSBarry Smith   PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */
517b0a32e0cSBarry Smith   PetscLogDouble factor_mallocs;                     /* number of mallocs during factorization */
5184e220ebcSLois Curfman McInnes } MatInfo;
5194e220ebcSLois Curfman McInnes 
520d9274352SBarry Smith /*E
521d9274352SBarry Smith     MatInfoType - Indicates if you want information about the local part of the matrix,
522d9274352SBarry Smith      the entire parallel matrix or the maximum over all the local parts.
523d9274352SBarry Smith 
524d9274352SBarry Smith     Level: beginner
525d9274352SBarry Smith 
526d9274352SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
527d9274352SBarry Smith 
528d9274352SBarry Smith .seealso: MatGetInfo(), MatInfo
529d9274352SBarry Smith E*/
5307b80b807SBarry Smith typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType;
531be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetInfo(Mat,MatInfoType,MatInfo*);
532be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatValid(Mat,PetscTruth*);
533be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetDiagonal(Mat,Vec);
534985db425SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMax(Mat,Vec,PetscInt[]);
535985db425SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMin(Mat,Vec,PetscInt[]);
536985db425SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMaxAbs(Mat,Vec,PetscInt[]);
537c87e5d42SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMinAbs(Mat,Vec,PetscInt[]);
53886697f06SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowSum(Mat,Vec);
539fc4dec0aSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatTranspose(Mat,MatReuse,Mat*);
540fc4dec0aSBarry Smith PetscPolymorphicFunction(MatTranspose,(Mat A),(A,MAT_INITIAL_MATRIX,&t),Mat,t)
541eeffb40dSHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatHermitianTranspose(Mat,MatReuse,Mat*);
542be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPermute(Mat,IS,IS,Mat *);
543ba966639SSatish Balay PetscPolymorphicFunction(MatPermute,(Mat A,IS is1,IS is2),(A,is1,is2,&t),Mat,t)
544be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPermuteSparsify(Mat,PetscInt,PetscReal,PetscReal,IS,IS,Mat *);
545be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalScale(Mat,Vec,Vec);
546be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalSet(Mat,Vec,InsertMode);
547be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatEqual(Mat,Mat,PetscTruth*);
548ba966639SSatish Balay PetscPolymorphicFunction(MatEqual,(Mat A,Mat B),(A,B,&t),PetscTruth,t)
549be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultEqual(Mat,Mat,PetscInt,PetscTruth*);
550be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultAddEqual(Mat,Mat,PetscInt,PetscTruth*);
551be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeEqual(Mat,Mat,PetscInt,PetscTruth*);
552be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscTruth*);
5537b80b807SBarry Smith 
554be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNorm(Mat,NormType,PetscReal *);
555ba966639SSatish Balay PetscPolymorphicFunction(MatNorm,(Mat A,NormType t),(A,t,&n),PetscReal,n)
556242f1d38SBarry Smith EXTERN PetscErrorCode MatGetColumnNorms(Mat,NormType,PetscReal *);
557be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroEntries(Mat);
558f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar);
559f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsIS(Mat,IS,PetscScalar);
560f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroColumns(Mat,PetscInt,const PetscInt [],const PetscScalar*);
561f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroColumnsIS(Mat,IS,const PetscScalar*);
5627b80b807SBarry Smith 
563be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatUseScaledForm(Mat,PetscTruth);
564be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScaleSystem(Mat,Vec,Vec);
565be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatUnScaleSystem(Mat,Vec,Vec);
5665ef9f2a5SBarry Smith 
567be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSize(Mat,PetscInt*,PetscInt*);
568be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalSize(Mat,PetscInt*,PetscInt*);
569be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRange(Mat,PetscInt*,PetscInt*);
5708ee2e534SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRanges(Mat,const PetscInt**);
571f4b8409dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRangeColumn(Mat,PetscInt*,PetscInt*);
57272ca8751SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRangesColumn(Mat,const PetscInt**);
5737b80b807SBarry Smith 
574be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
575be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroyMatrices(PetscInt,Mat *[]);
5764aa3045dSJed Brown EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrix(Mat,IS,IS,MatReuse,Mat *);
577f6d58c54SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSeqNonzeroStructure(Mat,Mat*);
578f6d58c54SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroySeqNonzeroStructure(Mat*);
5795494a064SHong Zhang 
580be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
581be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPI(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*);
582be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPISymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*);
583be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPINumeric(Mat,Mat);
584be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroy_MPIAIJ_SeqsToMPI(Mat);
585be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalMat(Mat,MatReuse,Mat*);
586be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*);
587be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,PetscInt*,Mat*);
5881d79065fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBrowsOfAoCols(Mat,Mat,MatReuse,PetscInt**,PetscInt**,MatScalar**,Mat*);
58943eb5e2fSMatthew Knepley #if defined (PETSC_USE_CTABLE)
590cd116e26SSatish Balay #include "petscctable.h"
59143eb5e2fSMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetCommunicationStructs(Mat, Vec *, PetscTable *, VecScatter *);
59243eb5e2fSMatthew Knepley #else
59343eb5e2fSMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetCommunicationStructs(Mat, Vec *, PetscInt *[], VecScatter *);
59443eb5e2fSMatthew Knepley #endif
5958c7482ecSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetGhosts(Mat, PetscInt *,const PetscInt *[]);
5968efafbd8SBarry Smith 
597be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt);
5987b80b807SBarry Smith 
599be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*);
600be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultSymbolic(Mat,Mat,PetscReal,Mat*);
601be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultNumeric(Mat,Mat,Mat);
60222440eb1SKris Buschelman 
603be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*);
604be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAPSymbolic(Mat,Mat,PetscReal,Mat*);
605be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAPNumeric(Mat,Mat,Mat);
60622440eb1SKris Buschelman 
607be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTranspose(Mat,Mat,MatReuse,PetscReal,Mat*);
608be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTransposeSymbolic(Mat,Mat,PetscReal,Mat*);
609be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTransposeNumeric(Mat,Mat,Mat);
610bc011b1eSHong Zhang 
611f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAXPY(Mat,PetscScalar,Mat,MatStructure);
61204aac2b0SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAYPX(Mat,PetscScalar,Mat,MatStructure);
6131c741599SHong Zhang 
614f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScale(Mat,PetscScalar);
615f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShift(Mat,PetscScalar);
6167b80b807SBarry Smith 
617be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping);
618be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetLocalToGlobalMappingBlock(Mat,ISLocalToGlobalMapping);
619f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar);
620f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsLocalIS(Mat,IS,PetscScalar);
621be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
622be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
623052efed2SBarry Smith 
624be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStashSetInitialSize(Mat,PetscInt,PetscInt);
625be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
62690f02eecSBarry Smith 
627be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInterpolate(Mat,Vec,Vec);
6280c451bc4SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInterpolateAdd(Mat,Vec,Vec,Vec);
629ba966639SSatish Balay PetscPolymorphicSubroutine(MatInterpolateAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
630be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestrict(Mat,Vec,Vec);
631be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetVecs(Mat,Vec*,Vec*);
63269db28dcSHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRedundantMatrix(Mat,PetscInt,MPI_Comm,PetscInt,MatReuse,Mat*);
633bd481603SBarry Smith 
634bd481603SBarry Smith /*MC
6351620fd73SBarry Smith    MatSetValue - Set a single entry into a matrix.
6361620fd73SBarry Smith 
6371620fd73SBarry Smith    Not collective
6381620fd73SBarry Smith 
6391620fd73SBarry Smith    Input Parameters:
6401620fd73SBarry Smith +  m - the matrix
6411620fd73SBarry Smith .  row - the row location of the entry
6421620fd73SBarry Smith .  col - the column location of the entry
6431620fd73SBarry Smith .  value - the value to insert
6441620fd73SBarry Smith -  mode - either INSERT_VALUES or ADD_VALUES
6451620fd73SBarry Smith 
6461620fd73SBarry Smith    Notes:
6471620fd73SBarry Smith    For efficiency one should use MatSetValues() and set several or many
6481620fd73SBarry Smith    values simultaneously if possible.
6491620fd73SBarry Smith 
6501620fd73SBarry Smith    Level: beginner
6511620fd73SBarry Smith 
6521620fd73SBarry Smith .seealso: MatSetValues(), MatSetValueLocal()
6531620fd73SBarry Smith M*/
6541620fd73SBarry 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);}
6551620fd73SBarry Smith 
6561620fd73SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatGetValue(Mat v,PetscInt i,PetscInt j,PetscScalar *va) {return MatGetValues(v,1,&i,1,&j,va);}
6571620fd73SBarry Smith 
6581620fd73SBarry 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);}
6591620fd73SBarry Smith 
6601620fd73SBarry Smith /*MC
6610d2a7ff2SSatish Balay    MatPreallocateInitialize - Begins the block of code that will count the number of nonzeros per
662bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
663bd481603SBarry Smith 
664bd481603SBarry Smith    Synopsis:
665c1ac3661SBarry Smith    PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)
666bd481603SBarry Smith 
667bd481603SBarry Smith    Collective on MPI_Comm
668bd481603SBarry Smith 
669bd481603SBarry Smith    Input Parameters:
670bd481603SBarry Smith +  comm - the communicator that will share the eventually allocated matrix
671859fcb39SBarry Smith .  nrows - the number of LOCAL rows in the matrix
672859fcb39SBarry Smith -  ncols - the number of LOCAL columns in the matrix
673bd481603SBarry Smith 
674bd481603SBarry Smith    Output Parameters:
675bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
676bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
677bd481603SBarry Smith 
678bd481603SBarry Smith 
679bd481603SBarry Smith    Level: intermediate
680bd481603SBarry Smith 
681bd481603SBarry Smith    Notes:
682bd481603SBarry Smith    See the chapter in the users manual on performance for more details
683bd481603SBarry Smith 
6841d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
685bd481603SBarry Smith 
686bd481603SBarry Smith    Use MatPreallocateInitializeSymmetric() for symmetric matrices (MPISBAIJ matrices)
687bd481603SBarry Smith 
6881620fd73SBarry Smith    This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize().
6891620fd73SBarry Smith 
690bd481603SBarry Smith   Concepts: preallocation^Matrix
691bd481603SBarry Smith 
692bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
693bd481603SBarry Smith           MatPreallocateInitializeSymmetric(), MatPreallocateSymmetricSetLocal()
694bd481603SBarry Smith M*/
695c4f061fbSSatish Balay #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \
6967c922b88SBarry Smith { \
697a89bc211SBarry Smith   PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__start,__end; \
698a89bc211SBarry Smith   _4_ierr = PetscMalloc2(__nrows,PetscInt,&dnz,__nrows,PetscInt,&onz);CHKERRQ(_4_ierr); \
699a89bc211SBarry Smith   _4_ierr = PetscMemzero(dnz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
700a89bc211SBarry Smith   _4_ierr = PetscMemzero(onz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
701c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ctmp;\
702a89bc211SBarry Smith   _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows;
7037c922b88SBarry Smith 
704bd481603SBarry Smith /*MC
7050d2a7ff2SSatish Balay    MatPreallocateSymmetricInitialize - Begins the block of code that will count the number of nonzeros per
706bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
707bd481603SBarry Smith 
708bd481603SBarry Smith    Synopsis:
709c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)
710bd481603SBarry Smith 
711bd481603SBarry Smith    Collective on MPI_Comm
712bd481603SBarry Smith 
713bd481603SBarry Smith    Input Parameters:
714bd481603SBarry Smith +  comm - the communicator that will share the eventually allocated matrix
715859fcb39SBarry Smith .  nrows - the number of LOCAL rows in the matrix
716859fcb39SBarry Smith -  ncols - the number of LOCAL columns in the matrix
717bd481603SBarry Smith 
718bd481603SBarry Smith    Output Parameters:
719bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
720bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
721bd481603SBarry Smith 
722bd481603SBarry Smith 
723bd481603SBarry Smith    Level: intermediate
724bd481603SBarry Smith 
725bd481603SBarry Smith    Notes:
726bd481603SBarry Smith    See the chapter in the users manual on performance for more details
727bd481603SBarry Smith 
7281d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
729bd481603SBarry Smith 
7301620fd73SBarry Smith    This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize().
7311620fd73SBarry Smith 
732bd481603SBarry Smith   Concepts: preallocation^Matrix
733bd481603SBarry Smith 
734bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
735bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal()
736bd481603SBarry Smith M*/
737222b16d4SBarry Smith #define MatPreallocateSymmetricInitialize(comm,nrows,ncols,dnz,onz) 0; \
738222b16d4SBarry Smith { \
739a89bc211SBarry Smith   PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__end; \
740a89bc211SBarry Smith   _4_ierr = PetscMalloc2(__nrows,PetscInt,&dnz,__nrows,PetscInt,&onz);CHKERRQ(_4_ierr); \
741a89bc211SBarry Smith   _4_ierr = PetscMemzero(dnz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
742a89bc211SBarry Smith   _4_ierr = PetscMemzero(onz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
743c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr);\
744a89bc211SBarry Smith   _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows;
745222b16d4SBarry Smith 
746bd481603SBarry Smith /*MC
747bd481603SBarry Smith    MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
748bd481603SBarry Smith        inserted using a local number of the rows and columns
749bd481603SBarry Smith 
750bd481603SBarry Smith    Synopsis:
751c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
752bd481603SBarry Smith 
753bd481603SBarry Smith    Not Collective
754bd481603SBarry Smith 
755bd481603SBarry Smith    Input Parameters:
756bd481603SBarry Smith +  map - the mapping between local numbering and global numbering
757bd481603SBarry Smith .  nrows - the number of rows indicated
7581d73ed98SBarry Smith .  rows - the indices of the rows
759bd481603SBarry Smith .  ncols - the number of columns in the matrix
760bd481603SBarry Smith .  cols - the columns indicated
761bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
762bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
763bd481603SBarry Smith 
764bd481603SBarry Smith 
765bd481603SBarry Smith    Level: intermediate
766bd481603SBarry Smith 
767bd481603SBarry Smith    Notes:
768bd481603SBarry Smith    See the chapter in the users manual on performance for more details
769bd481603SBarry Smith 
7701d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
771bd481603SBarry Smith 
772bd481603SBarry Smith   Concepts: preallocation^Matrix
773bd481603SBarry Smith 
774bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
775bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal()
776bd481603SBarry Smith M*/
777c4f061fbSSatish Balay #define MatPreallocateSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\
778c4f061fbSSatish Balay {\
779c1ac3661SBarry Smith   PetscInt __l;\
780ef66eb69SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\
781ef66eb69SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\
782c4f061fbSSatish Balay   for (__l=0;__l<nrows;__l++) {\
783ef66eb69SBarry Smith     _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
784c4f061fbSSatish Balay   }\
785c4f061fbSSatish Balay }
786c4f061fbSSatish Balay 
787bd481603SBarry Smith /*MC
788bd481603SBarry Smith    MatPreallocateSymmetricSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
789bd481603SBarry Smith        inserted using a local number of the rows and columns
790bd481603SBarry Smith 
791bd481603SBarry Smith    Synopsis:
792c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
793bd481603SBarry Smith 
794bd481603SBarry Smith    Not Collective
795bd481603SBarry Smith 
796bd481603SBarry Smith    Input Parameters:
797bd481603SBarry Smith +  map - the mapping between local numbering and global numbering
798bd481603SBarry Smith .  nrows - the number of rows indicated
7991d73ed98SBarry Smith .  rows - the indices of the rows
800bd481603SBarry Smith .  ncols - the number of columns in the matrix
801bd481603SBarry Smith .  cols - the columns indicated
802bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
803bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
804bd481603SBarry Smith 
805bd481603SBarry Smith 
806bd481603SBarry Smith    Level: intermediate
807bd481603SBarry Smith 
808bd481603SBarry Smith    Notes:
809bd481603SBarry Smith    See the chapter in the users manual on performance for more details
810bd481603SBarry Smith 
811bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
812bd481603SBarry Smith 
813bd481603SBarry Smith   Concepts: preallocation^Matrix
814bd481603SBarry Smith 
815bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
816bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
817bd481603SBarry Smith M*/
818d3d32019SBarry Smith #define MatPreallocateSymmetricSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\
819d3d32019SBarry Smith {\
820c1ac3661SBarry Smith   PetscInt __l;\
821d3d32019SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\
822d3d32019SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\
823d3d32019SBarry Smith   for (__l=0;__l<nrows;__l++) {\
824d3d32019SBarry Smith     _4_ierr = MatPreallocateSymmetricSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
825d3d32019SBarry Smith   }\
826d3d32019SBarry Smith }
827d3d32019SBarry Smith 
828bd481603SBarry Smith /*MC
829bd481603SBarry Smith    MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
830bd481603SBarry Smith        inserted using a local number of the rows and columns
831bd481603SBarry Smith 
832bd481603SBarry Smith    Synopsis:
833c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
834bd481603SBarry Smith 
835bd481603SBarry Smith    Not Collective
836bd481603SBarry Smith 
837bd481603SBarry Smith    Input Parameters:
83864075487SBarry Smith +  row - the row
839bd481603SBarry Smith .  ncols - the number of columns in the matrix
840a50f8bf6SHong Zhang -  cols - the columns indicated
841a50f8bf6SHong Zhang 
842a50f8bf6SHong Zhang    Output Parameters:
843a50f8bf6SHong Zhang +  dnz - the array that will be passed to the matrix preallocation routines
844bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
845bd481603SBarry Smith 
846bd481603SBarry Smith 
847bd481603SBarry Smith    Level: intermediate
848bd481603SBarry Smith 
849bd481603SBarry Smith    Notes:
850bd481603SBarry Smith    See the chapter in the users manual on performance for more details
851bd481603SBarry Smith 
852bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
853bd481603SBarry Smith 
8541620fd73SBarry Smith    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
8551620fd73SBarry Smith 
856bd481603SBarry Smith   Concepts: preallocation^Matrix
857bd481603SBarry Smith 
858bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
859bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
860bd481603SBarry Smith M*/
861c4f061fbSSatish Balay #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\
862c1ac3661SBarry Smith { PetscInt __i; \
8638ee2e534SBarry Smith   if (row < __rstart) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Trying to set preallocation for row %D less than first local row %D",row,__rstart);\
864a89bc211SBarry Smith   if (row >= __rstart+__nrows) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Trying to set preallocation for row %D greater than last local row %D",row,__rstart+__nrows-1);\
8657c922b88SBarry Smith   for (__i=0; __i<nc; __i++) {\
86664075487SBarry Smith     if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \
8677cc688d7SBarry Smith     else dnz[row - __rstart]++;\
8687c922b88SBarry Smith   }\
8697c922b88SBarry Smith }
8707c922b88SBarry Smith 
871bd481603SBarry Smith /*MC
872bd481603SBarry Smith    MatPreallocateSymmetricSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
873bd481603SBarry Smith        inserted using a local number of the rows and columns
874bd481603SBarry Smith 
875bd481603SBarry Smith    Synopsis:
876c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
877bd481603SBarry Smith 
878bd481603SBarry Smith    Not Collective
879bd481603SBarry Smith 
880bd481603SBarry Smith    Input Parameters:
881bd481603SBarry Smith +  nrows - the number of rows indicated
8821d73ed98SBarry Smith .  rows - the indices of the rows
883bd481603SBarry Smith .  ncols - the number of columns in the matrix
884bd481603SBarry Smith .  cols - the columns indicated
885bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
886bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
887bd481603SBarry Smith 
888bd481603SBarry Smith 
889bd481603SBarry Smith    Level: intermediate
890bd481603SBarry Smith 
891bd481603SBarry Smith    Notes:
892bd481603SBarry Smith    See the chapter in the users manual on performance for more details
893bd481603SBarry Smith 
894bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
895bd481603SBarry Smith 
8961620fd73SBarry Smith    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
8971620fd73SBarry Smith 
898bd481603SBarry Smith   Concepts: preallocation^Matrix
899bd481603SBarry Smith 
900bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
901bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
902bd481603SBarry Smith M*/
903d3d32019SBarry Smith #define MatPreallocateSymmetricSet(row,nc,cols,dnz,onz) 0;\
904c1ac3661SBarry Smith { PetscInt __i; \
905d3d32019SBarry Smith   for (__i=0; __i<nc; __i++) {\
906d3d32019SBarry Smith     if (cols[__i] >= __end) onz[row - __rstart]++; \
907d3d32019SBarry Smith     else if (cols[__i] >= row) dnz[row - __rstart]++;\
908d3d32019SBarry Smith   }\
909d3d32019SBarry Smith }
910d3d32019SBarry Smith 
911bd481603SBarry Smith /*MC
91216371a99SBarry Smith    MatPreallocateLocation -  An alternative to MatPreallocationSet() that puts the nonzero locations into the matrix if it exists
91316371a99SBarry Smith 
91416371a99SBarry Smith    Synopsis:
91516371a99SBarry Smith    PetscErrorCode MatPreallocateLocations(Mat A,PetscInt row,PetscInt ncols,PetscInt *cols,PetscInt *dnz,PetscInt *onz)
91616371a99SBarry Smith 
91716371a99SBarry Smith    Not Collective
91816371a99SBarry Smith 
91916371a99SBarry Smith    Input Parameters:
92016371a99SBarry Smith .  A - matrix
92116371a99SBarry Smith .  row - row where values exist (must be local to this process)
92216371a99SBarry Smith .  ncols - number of columns
92316371a99SBarry Smith .  cols - columns with nonzeros
92416371a99SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
92516371a99SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
92616371a99SBarry Smith 
92716371a99SBarry Smith 
92816371a99SBarry Smith    Level: intermediate
92916371a99SBarry Smith 
93016371a99SBarry Smith    Notes:
93116371a99SBarry Smith    See the chapter in the users manual on performance for more details
93216371a99SBarry Smith 
93316371a99SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
93416371a99SBarry Smith 
93516371a99SBarry Smith    This is a MACRO not a function because it uses a bunch of variables private to the MatPreallocation.... routines.
93616371a99SBarry Smith 
93716371a99SBarry Smith   Concepts: preallocation^Matrix
93816371a99SBarry Smith 
93916371a99SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
94016371a99SBarry Smith           MatPreallocateSymmetricInitialize(), MatPreallocateSymmetricSetLocal()
94116371a99SBarry Smith M*/
94216371a99SBarry 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);}
94316371a99SBarry Smith 
94416371a99SBarry Smith 
94516371a99SBarry Smith /*MC
9460d2a7ff2SSatish Balay    MatPreallocateFinalize - Ends the block of code that will count the number of nonzeros per
947bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
948bd481603SBarry Smith 
949bd481603SBarry Smith    Synopsis:
950c1ac3661SBarry Smith    PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz)
951bd481603SBarry Smith 
952bd481603SBarry Smith    Collective on MPI_Comm
953bd481603SBarry Smith 
954bd481603SBarry Smith    Input Parameters:
95516371a99SBarry Smith +  dnz - the array that was be passed to the matrix preallocation routines
956bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
957bd481603SBarry Smith 
958bd481603SBarry Smith 
959bd481603SBarry Smith    Level: intermediate
960bd481603SBarry Smith 
961bd481603SBarry Smith    Notes:
962bd481603SBarry Smith    See the chapter in the users manual on performance for more details
963bd481603SBarry Smith 
964bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
965bd481603SBarry Smith 
9661620fd73SBarry Smith    This is a MACRO not a function because it closes the { started in MatPreallocateInitialize().
9671620fd73SBarry Smith 
968bd481603SBarry Smith   Concepts: preallocation^Matrix
969bd481603SBarry Smith 
970bd481603SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
971bd481603SBarry Smith           MatPreallocateSymmetricInitialize(), MatPreallocateSymmetricSetLocal()
972bd481603SBarry Smith M*/
973a89bc211SBarry Smith #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree2(dnz,onz);CHKERRQ(_4_ierr);}
9747c922b88SBarry Smith 
975bd481603SBarry Smith 
976bd481603SBarry Smith 
9777b80b807SBarry Smith /* Routines unique to particular data structures */
978be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellGetContext(Mat,void **);
979ba966639SSatish Balay PetscPolymorphicFunction(MatShellGetContext,(Mat A),(A,&t),void*,t)
980698d4c6aSKris Buschelman 
981be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInodeAdjustForInodes(Mat,IS*,IS*);
982be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *);
983698d4c6aSKris Buschelman 
984be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetColumnIndices(Mat,PetscInt[]);
985be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetColumnIndices(Mat,PetscInt[]);
986be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
987c75a6043SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
988c75a6043SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
9897b80b807SBarry Smith 
990a96a251dSBarry Smith #define MAT_SKIP_ALLOCATION -4
991a96a251dSBarry Smith 
992be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
993ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz))
994be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
995ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqSBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz))
996be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]);
997ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqAIJSetPreallocation,(Mat A,const PetscInt nnz[]),(A,0,nnz))
998273d9f13SBarry Smith 
999be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
1000ba966639SSatish Balay PetscPolymorphicSubroutine(MatMPIBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[],const PetscInt onz[]),(A,bs,0,nnz,0,onz))
1001be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
1002be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
100353803f26SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []);
1004725b52f3SLisandro Dalcin EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
1005be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]);
1006be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
1007be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]);
1008be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIDenseSetPreallocation(Mat,PetscScalar[]);
1009284134d9SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation(Mat,PetscScalar[]);
1010be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,PetscInt*[]);
1011be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,PetscInt*[]);
1012be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAdicSetLocalFunction(Mat,void (*)(void));
1013273d9f13SBarry Smith 
1014be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetLDA(Mat,PetscInt);
1015ab92ecdeSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDenseGetLocalMatrix(Mat,Mat*);
10161b807ce4Svictorle 
1017be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStoreValues(Mat);
1018be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRetrieveValues(Mat);
10192e8a6d31SBarry Smith 
1020be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDAADSetCtx(Mat,void*);
10213a7fca6bSBarry Smith 
10227b80b807SBarry Smith /*
10237b80b807SBarry Smith   These routines are not usually accessed directly, rather solving is
102494b7f48cSBarry Smith   done through the KSP and PC interfaces.
10257b80b807SBarry Smith */
10267b80b807SBarry Smith 
1027d9274352SBarry Smith /*E
1028d9274352SBarry Smith     MatOrderingType - String with the name of a PETSc matrix ordering or the creation function
1029d9274352SBarry Smith        with an optional dynamic library name, for example
1030d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:orderingcreate()
1031d9274352SBarry Smith 
1032d9274352SBarry Smith    Level: beginner
1033d9274352SBarry Smith 
1034e5a9bf91SBarry Smith    Cannot use const because the PC objects manipulate the string
1035e5a9bf91SBarry Smith 
1036d9274352SBarry Smith .seealso: MatGetOrdering()
1037d9274352SBarry Smith E*/
10383eaad2caSSatish Balay #define MatOrderingType char*
1039b12f92e5SBarry Smith #define MATORDERING_NATURAL     "natural"
1040b12f92e5SBarry Smith #define MATORDERING_ND          "nd"
1041b12f92e5SBarry Smith #define MATORDERING_1WD         "1wd"
1042b12f92e5SBarry Smith #define MATORDERING_RCM         "rcm"
1043b12f92e5SBarry Smith #define MATORDERING_QMD         "qmd"
1044b12f92e5SBarry Smith #define MATORDERING_ROWLENGTH   "rowlength"
1045e106eecfSBarry Smith #define MATORDERING_DSC_ND      "dsc_nd"         /* these three are only for DSCPACK, see its documentation for details */
104662152c8bSBarry Smith #define MATORDERING_DSC_MMD     "dsc_mmd"
104762152c8bSBarry Smith #define MATORDERING_DSC_MDF     "dsc_mdf"
1048c06d978dSMatthew Knepley #define MATORDERING_CONSTRAINED "constrained"
1049c06d978dSMatthew Knepley #define MATORDERING_IDENTITY    "identity"
1050c06d978dSMatthew Knepley #define MATORDERING_REVERSE     "reverse"
10518fa24674SBarry Smith #define MATORDERING_FLOW        "flow"
1052b12f92e5SBarry Smith 
1053f2416520SSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOrdering(Mat,const MatOrderingType,IS*,IS*);
1054be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOrderingList(PetscFList *list);
1055f2416520SSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,const MatOrderingType,IS*,IS*));
105630de9b25SBarry Smith 
105730de9b25SBarry Smith /*MC
105830de9b25SBarry Smith    MatOrderingRegisterDynamic - Adds a new sparse matrix ordering to the
105930de9b25SBarry Smith                                matrix package.
106030de9b25SBarry Smith 
106130de9b25SBarry Smith    Synopsis:
1062c1ac3661SBarry Smith    PetscErrorCode MatOrderingRegisterDynamic(char *name_ordering,char *path,char *name_create,PetscErrorCode (*routine_create)(MatOrdering))
106330de9b25SBarry Smith 
106430de9b25SBarry Smith    Not Collective
106530de9b25SBarry Smith 
106630de9b25SBarry Smith    Input Parameters:
106730de9b25SBarry Smith +  sname - name of ordering (for example MATORDERING_ND)
106830de9b25SBarry Smith .  path - location of library where creation routine is
106930de9b25SBarry Smith .  name - name of function that creates the ordering type,a string
107030de9b25SBarry Smith -  function - function pointer that creates the ordering
107130de9b25SBarry Smith 
107230de9b25SBarry Smith    Level: developer
107330de9b25SBarry Smith 
107430de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
107530de9b25SBarry Smith    is ignored.
107630de9b25SBarry Smith 
107730de9b25SBarry Smith    Sample usage:
107830de9b25SBarry Smith .vb
107930de9b25SBarry Smith    MatOrderingRegisterDynamic("my_order",/home/username/my_lib/lib/libO/solaris/mylib.a,
108030de9b25SBarry Smith                "MyOrder",MyOrder);
108130de9b25SBarry Smith .ve
108230de9b25SBarry Smith 
108330de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
108430de9b25SBarry Smith $     MatOrderingSetType(part,"my_order)
108530de9b25SBarry Smith    or at runtime via the option
10862401956bSBarry Smith $     -pc_factor_mat_ordering_type my_order
108730de9b25SBarry Smith 
1088ab901514SBarry Smith    ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values.
108930de9b25SBarry Smith 
109030de9b25SBarry Smith .keywords: matrix, ordering, register
109130de9b25SBarry Smith 
109230de9b25SBarry Smith .seealso: MatOrderingRegisterDestroy(), MatOrderingRegisterAll()
109330de9b25SBarry Smith M*/
1094aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1095f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,0)
1096b12f92e5SBarry Smith #else
1097f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,d)
1098b12f92e5SBarry Smith #endif
109930de9b25SBarry Smith 
1100be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegisterDestroy(void);
1101be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegisterAll(const char[]);
11022bad1931SBarry Smith extern PetscTruth MatOrderingRegisterAllCalled;
1103b0a32e0cSBarry Smith extern PetscFList MatOrderingList;
1104d4fbbf0eSBarry Smith 
1105be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS);
1106a2ce50c7SBarry Smith 
1107d91e6319SBarry Smith /*S
11082401956bSBarry Smith    MatFactorInfo - Data passed into the matrix factorization routines
1109b00f7748SHong Zhang 
111061cad860SBarry Smith    In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE, that is use
111161cad860SBarry Smith $     MatFactorInfo  info(MAT_FACTORINFO_SIZE)
1112b00f7748SHong Zhang 
111315e8a5b3SHong Zhang    Notes: These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC.
1114b00f7748SHong Zhang 
111561cad860SBarry Smith       You can use MatFactorInfoInitialize() to set default values.
111661cad860SBarry Smith 
1117b00f7748SHong Zhang    Level: developer
1118b00f7748SHong Zhang 
1119d7d82daaSBarry Smith .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(),
1120d7d82daaSBarry Smith           MatFactorInfoInitialize()
1121b00f7748SHong Zhang 
1122b00f7748SHong Zhang S*/
1123b00f7748SHong Zhang typedef struct {
11240a29876aSHong Zhang   PetscReal     shiftnz;        /* scaling of identity added to matrix to prevent zero pivots */
1125fbf22428SSatish Balay   PetscReal     shiftpd;        /* if true, shift until positive pivots */
112615e8a5b3SHong Zhang   PetscReal     diagonal_fill;  /* force diagonal to fill in if initially not filled */
112785317021SBarry Smith   PetscReal     usedt;
112815e8a5b3SHong Zhang   PetscReal     dt;             /* drop tolerance */
1129b00f7748SHong Zhang   PetscReal     dtcol;          /* tolerance for pivoting */
113015e8a5b3SHong Zhang   PetscReal     dtcount;        /* maximum nonzeros to be allowed per row */
113167e28bfeSBarry Smith   PetscReal     fill;           /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */
1132348344bbSBarry Smith   PetscReal     levels;         /* ICC/ILU(levels) */
1133bcd9e38bSBarry Smith   PetscReal     pivotinblocks;  /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0
1134bcd9e38bSBarry Smith                                    factorization may be faster if do not pivot */
113562bba022SBarry Smith   PetscReal     shiftinblocks;  /* if block in block factorization has zero pivot then shift diagonal until non-singular */
113615e8a5b3SHong Zhang   PetscReal     zeropivot;      /* pivot is called zero if less than this */
113715e8a5b3SHong Zhang } MatFactorInfo;
1138ffa6d0a5SLois Curfman McInnes 
1139be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFactorInfoInitialize(MatFactorInfo*);
11400481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactor(Mat,IS,const MatFactorInfo*);
11410481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
11420481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactorNumeric(Mat,Mat,const MatFactorInfo*);
11430481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactor(Mat,IS,IS,const MatFactorInfo*);
11440481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUFactor(Mat,IS,IS,const MatFactorInfo*);
11450481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
11460481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
11470481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatICCFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
11480481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatICCFactor(Mat,IS,const MatFactorInfo*);
11490481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactorNumeric(Mat,Mat,const MatFactorInfo*);
11500481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUDTFactor(Mat,IS,IS,const MatFactorInfo*,Mat *);
11513c2a7987SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUDTFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
11523c2a7987SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUDTFactorNumeric(Mat,Mat,const MatFactorInfo*);
1153be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*);
1154be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolve(Mat,Vec,Vec);
1155be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatForwardSolve(Mat,Vec,Vec);
1156be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatBackwardSolve(Mat,Vec,Vec);
1157be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveAdd(Mat,Vec,Vec,Vec);
1158be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveTranspose(Mat,Vec,Vec);
1159be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveTransposeAdd(Mat,Vec,Vec,Vec);
1160be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolves(Mat,Vecs,Vecs);
11618ed539a5SBarry Smith 
1162be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUnfactored(Mat);
1163bb5a7306SBarry Smith 
1164d91e6319SBarry Smith /*E
1165d91e6319SBarry Smith     MatSORType - What type of (S)SOR to perform
1166bb1eb677SSatish Balay 
1167d91e6319SBarry Smith     Level: beginner
1168d91e6319SBarry Smith 
1169d9274352SBarry Smith    May be bitwise ORd together
1170d9274352SBarry Smith 
1171d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
1172d91e6319SBarry Smith 
11734e7234bfSBarry Smith    MatSORType may be bitwise ORd together, so do not change the numbers
11744e7234bfSBarry Smith 
117541f059aeSBarry Smith .seealso: MatSOR()
1176d91e6319SBarry Smith E*/
1177ee50ffe9SBarry Smith typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3,
1178ee50ffe9SBarry Smith               SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8,
1179ee50ffe9SBarry Smith               SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16,
118084cb2905SBarry Smith               SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType;
118141f059aeSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSOR(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
11828ed539a5SBarry Smith 
1183d4fbbf0eSBarry Smith /*
1184639f9d9dSBarry Smith     These routines are for efficiently computing Jacobians via finite differences.
1185639f9d9dSBarry Smith */
1186b12f92e5SBarry Smith 
1187d9274352SBarry Smith /*E
1188d9274352SBarry Smith     MatColoringType - String with the name of a PETSc matrix coloring or the creation function
1189d9274352SBarry Smith        with an optional dynamic library name, for example
1190d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:coloringcreate()
1191d9274352SBarry Smith 
1192d9274352SBarry Smith    Level: beginner
1193d9274352SBarry Smith 
1194d9274352SBarry Smith .seealso: MatGetColoring()
1195d9274352SBarry Smith E*/
1196a313700dSBarry Smith #define MatColoringType char*
1197b12f92e5SBarry Smith #define MATCOLORING_NATURAL "natural"
1198b12f92e5SBarry Smith #define MATCOLORING_SL      "sl"
1199b12f92e5SBarry Smith #define MATCOLORING_LF      "lf"
1200b12f92e5SBarry Smith #define MATCOLORING_ID      "id"
1201b12f92e5SBarry Smith 
1202ab970d9bSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColoring(Mat,const MatColoringType,ISColoring*);
12032e90c967SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,MatColoringType,ISColoring *));
120430de9b25SBarry Smith 
120530de9b25SBarry Smith /*MC
120630de9b25SBarry Smith    MatColoringRegisterDynamic - Adds a new sparse matrix coloring to the
120730de9b25SBarry Smith                                matrix package.
120830de9b25SBarry Smith 
120930de9b25SBarry Smith    Synopsis:
1210c1ac3661SBarry Smith    PetscErrorCode MatColoringRegisterDynamic(char *name_coloring,char *path,char *name_create,PetscErrorCode (*routine_create)(MatColoring))
121130de9b25SBarry Smith 
121230de9b25SBarry Smith    Not Collective
121330de9b25SBarry Smith 
121430de9b25SBarry Smith    Input Parameters:
121530de9b25SBarry Smith +  sname - name of Coloring (for example MATCOLORING_SL)
121630de9b25SBarry Smith .  path - location of library where creation routine is
121730de9b25SBarry Smith .  name - name of function that creates the Coloring type, a string
121830de9b25SBarry Smith -  function - function pointer that creates the coloring
121930de9b25SBarry Smith 
122030de9b25SBarry Smith    Level: developer
122130de9b25SBarry Smith 
122230de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
122330de9b25SBarry Smith    is ignored.
122430de9b25SBarry Smith 
122530de9b25SBarry Smith    Sample usage:
122630de9b25SBarry Smith .vb
122730de9b25SBarry Smith    MatColoringRegisterDynamic("my_color",/home/username/my_lib/lib/libO/solaris/mylib.a,
122830de9b25SBarry Smith                "MyColor",MyColor);
122930de9b25SBarry Smith .ve
123030de9b25SBarry Smith 
123130de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
123230de9b25SBarry Smith $     MatColoringSetType(part,"my_color")
123330de9b25SBarry Smith    or at runtime via the option
123430de9b25SBarry Smith $     -mat_coloring_type my_color
123530de9b25SBarry Smith 
1236ab901514SBarry Smith    $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
123730de9b25SBarry Smith 
123830de9b25SBarry Smith .keywords: matrix, Coloring, register
123930de9b25SBarry Smith 
124030de9b25SBarry Smith .seealso: MatColoringRegisterDestroy(), MatColoringRegisterAll()
124130de9b25SBarry Smith M*/
1242aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1243f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,0)
1244b12f92e5SBarry Smith #else
1245f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,d)
1246b12f92e5SBarry Smith #endif
124730de9b25SBarry Smith 
12482bad1931SBarry Smith extern PetscTruth MatColoringRegisterAllCalled;
1249f1144a30SSatish Balay 
1250be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegisterAll(const char[]);
1251be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegisterDestroy(void);
1252be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*);
1253639f9d9dSBarry Smith 
1254d9274352SBarry Smith /*S
1255d9274352SBarry Smith      MatFDColoring - Object for computing a sparse Jacobian via finite differences
1256d9274352SBarry Smith         and coloring
1257639f9d9dSBarry Smith 
1258d9274352SBarry Smith    Level: beginner
1259d9274352SBarry Smith 
1260d9274352SBarry Smith   Concepts: coloring, sparse Jacobian, finite differences
1261d9274352SBarry Smith 
1262d9274352SBarry Smith .seealso:  MatFDColoringCreate()
1263d9274352SBarry Smith S*/
1264e2a1c21fSSatish Balay typedef struct _p_MatFDColoring* MatFDColoring;
1265639f9d9dSBarry Smith 
1266be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringCreate(Mat,ISColoring,MatFDColoring *);
1267be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringDestroy(MatFDColoring);
1268be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringView(MatFDColoring,PetscViewer);
1269be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*);
12704a9d489dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**);
1271be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal);
1272be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFromOptions(MatFDColoring);
1273be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApply(Mat,MatFDColoring,Vec,MatStructure*,void *);
1274be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApplyTS(Mat,MatFDColoring,PetscReal,Vec,MatStructure*,void *);
1275be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetF(MatFDColoring,Vec);
1276be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,PetscInt*[]);
1277639f9d9dSBarry Smith /*
12780752156aSBarry Smith     These routines are for partitioning matrices: currently used only
12793eda8832SBarry Smith   for adjacency matrix, MatCreateMPIAdj().
12800752156aSBarry Smith */
1281ca161407SBarry Smith 
1282d9274352SBarry Smith /*S
1283d9274352SBarry Smith      MatPartitioning - Object for managing the partitioning of a matrix or graph
1284d9274352SBarry Smith 
1285d9274352SBarry Smith    Level: beginner
1286d9274352SBarry Smith 
1287d9274352SBarry Smith   Concepts: partitioning
1288d9274352SBarry Smith 
1289743a1026Svictorle .seealso:  MatPartitioningCreate(), MatPartitioningType
1290d9274352SBarry Smith S*/
129191e9ee9fSBarry Smith typedef struct _p_MatPartitioning* MatPartitioning;
1292d9274352SBarry Smith 
1293d9274352SBarry Smith /*E
12945bc6e7ccSvictorle     MatPartitioningType - String with the name of a PETSc matrix partitioning or the creation function
1295d9274352SBarry Smith        with an optional dynamic library name, for example
1296d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:partitioningcreate()
1297d9274352SBarry Smith 
1298d9274352SBarry Smith    Level: beginner
1299d9274352SBarry Smith 
1300b547a13bSvictorle .seealso: MatPartitioningCreate(), MatPartitioning
1301d9274352SBarry Smith E*/
1302a313700dSBarry Smith #define MatPartitioningType char*
13038ba1e511SMatthew Knepley #define MAT_PARTITIONING_CURRENT  "current"
1304c2ec2a73SHong Zhang #define MAT_PARTITIONING_SQUARE   "square"
13058ba1e511SMatthew Knepley #define MAT_PARTITIONING_PARMETIS "parmetis"
130671306c60Sjroman #define MAT_PARTITIONING_CHACO    "chaco"
130771306c60Sjroman #define MAT_PARTITIONING_JOSTLE   "jostle"
130871306c60Sjroman #define MAT_PARTITIONING_PARTY    "party"
130971306c60Sjroman #define MAT_PARTITIONING_SCOTCH   "scotch"
131071306c60Sjroman 
1311ca161407SBarry Smith 
1312be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningCreate(MPI_Comm,MatPartitioning*);
1313a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetType(MatPartitioning,const MatPartitioningType);
1314be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetNParts(MatPartitioning,PetscInt);
1315be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetAdjacency(MatPartitioning,Mat);
1316be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]);
1317be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []);
1318be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningApply(MatPartitioning,IS*);
1319be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningDestroy(MatPartitioning);
13202aabb6bbSBarry Smith 
1321be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatPartitioning));
132230de9b25SBarry Smith 
132330de9b25SBarry Smith /*MC
132430de9b25SBarry Smith    MatPartitioningRegisterDynamic - Adds a new sparse matrix partitioning to the
132530de9b25SBarry Smith    matrix package.
132630de9b25SBarry Smith 
132730de9b25SBarry Smith    Synopsis:
1328c1ac3661SBarry Smith    PetscErrorCode MatPartitioningRegisterDynamic(char *name_partitioning,char *path,char *name_create,PetscErrorCode (*routine_create)(MatPartitioning))
132930de9b25SBarry Smith 
133030de9b25SBarry Smith    Not Collective
133130de9b25SBarry Smith 
133230de9b25SBarry Smith    Input Parameters:
133330de9b25SBarry Smith +  sname - name of partitioning (for example MAT_PARTITIONING_CURRENT) or parmetis
133430de9b25SBarry Smith .  path - location of library where creation routine is
133530de9b25SBarry Smith .  name - name of function that creates the partitioning type, a string
133630de9b25SBarry Smith -  function - function pointer that creates the partitioning type
133730de9b25SBarry Smith 
133830de9b25SBarry Smith    Level: developer
133930de9b25SBarry Smith 
134030de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
134130de9b25SBarry Smith    is ignored.
134230de9b25SBarry Smith 
134330de9b25SBarry Smith    Sample usage:
134430de9b25SBarry Smith .vb
134530de9b25SBarry Smith    MatPartitioningRegisterDynamic("my_part",/home/username/my_lib/lib/libO/solaris/mylib.a,
134630de9b25SBarry Smith                "MyPartCreate",MyPartCreate);
134730de9b25SBarry Smith .ve
134830de9b25SBarry Smith 
134930de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
135030de9b25SBarry Smith $     MatPartitioningSetType(part,"my_part")
135130de9b25SBarry Smith    or at runtime via the option
135230de9b25SBarry Smith $     -mat_partitioning_type my_part
135330de9b25SBarry Smith 
1354ab901514SBarry Smith    $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
135530de9b25SBarry Smith 
135630de9b25SBarry Smith .keywords: matrix, partitioning, register
135730de9b25SBarry Smith 
135830de9b25SBarry Smith .seealso: MatPartitioningRegisterDestroy(), MatPartitioningRegisterAll()
135930de9b25SBarry Smith M*/
1360aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1361f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,0)
13622aabb6bbSBarry Smith #else
1363f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,d)
13642aabb6bbSBarry Smith #endif
13652aabb6bbSBarry Smith 
13662bad1931SBarry Smith extern PetscTruth MatPartitioningRegisterAllCalled;
1367f1144a30SSatish Balay 
1368be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegisterAll(const char[]);
1369be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegisterDestroy(void);
13702bad1931SBarry Smith 
1371be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningView(MatPartitioning,PetscViewer);
1372be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetFromOptions(MatPartitioning);
1373a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningGetType(MatPartitioning,const MatPartitioningType*);
1374ca161407SBarry Smith 
1375be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningParmetisSetCoarseSequential(MatPartitioning);
13760e8af3caSMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningParmetisGetEdgeCut(MatPartitioning, PetscInt *);
13770752156aSBarry Smith 
1378be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningJostleSetCoarseLevel(MatPartitioning,PetscReal);
1379be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningJostleSetCoarseSequential(MatPartitioning);
138071306c60Sjroman 
1381954b3ec6SBarry Smith typedef enum {MP_CHACO_MULTILEVEL_KL,MP_CHACO_SPECTRAL,MP_CHACO_LINEAR,MP_CHACO_RANDOM, MP_CHACO_SCATTERED} MPChacoGlobalType;
1382be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetGlobal(MatPartitioning, MPChacoGlobalType);
138371306c60Sjroman typedef enum { MP_CHACO_KERNIGHAN_LIN, MP_CHACO_NONE } MPChacoLocalType;
1384be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetLocal(MatPartitioning, MPChacoLocalType);
1385be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal);
138671306c60Sjroman typedef enum { MP_CHACO_LANCZOS, MP_CHACO_RQI_SYMMLQ } MPChacoEigenType;
1387be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType);
1388be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenTol(MatPartitioning, PetscReal);
1389be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenNumber(MatPartitioning, PetscInt);
139071306c60Sjroman 
139171306c60Sjroman #define MP_PARTY_OPT "opt"
139271306c60Sjroman #define MP_PARTY_LIN "lin"
139371306c60Sjroman #define MP_PARTY_SCA "sca"
139471306c60Sjroman #define MP_PARTY_RAN "ran"
139571306c60Sjroman #define MP_PARTY_GBF "gbf"
139671306c60Sjroman #define MP_PARTY_GCF "gcf"
139771306c60Sjroman #define MP_PARTY_BUB "bub"
139871306c60Sjroman #define MP_PARTY_DEF "def"
1399be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetGlobal(MatPartitioning, const char*);
140071306c60Sjroman #define MP_PARTY_HELPFUL_SETS "hs"
140171306c60Sjroman #define MP_PARTY_KERNIGHAN_LIN "kl"
140271306c60Sjroman #define MP_PARTY_NONE "no"
1403be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetLocal(MatPartitioning, const char*);
1404be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal);
1405be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetBipart(MatPartitioning,PetscTruth);
1406be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscTruth);
140771306c60Sjroman 
140871306c60Sjroman typedef enum { MP_SCOTCH_GREEDY, MP_SCOTCH_GPS, MP_SCOTCH_GR_GPS } MPScotchGlobalType;
1409be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetArch(MatPartitioning,const char*);
1410be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetMultilevel(MatPartitioning);
1411be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetGlobal(MatPartitioning,MPScotchGlobalType);
1412be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetCoarseLevel(MatPartitioning,PetscReal);
1413be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetHostList(MatPartitioning,const char*);
141471306c60Sjroman typedef enum { MP_SCOTCH_KERNIGHAN_LIN, MP_SCOTCH_NONE } MPScotchLocalType;
1415be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetLocal(MatPartitioning,MPScotchLocalType);
1416be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetMapping(MatPartitioning);
1417be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetStrategy(MatPartitioning,char*);
141871306c60Sjroman 
1419591294e4SBarry Smith EXTERN PetscErrorCode MatMeshToVertexGraph(Mat,PetscInt,Mat*);
142014196391SBarry Smith EXTERN PetscErrorCode MatMeshToCellGraph(Mat,PetscInt,Mat*);
1421591294e4SBarry Smith 
14220752156aSBarry Smith /*
14230a835dfdSSatish Balay     If you add entries here you must also add them to finclude/petscmat.h
1424d4fbbf0eSBarry Smith */
14251c1c02c0SLois Curfman McInnes typedef enum { MATOP_SET_VALUES=0,
14261c1c02c0SLois Curfman McInnes                MATOP_GET_ROW=1,
14271c1c02c0SLois Curfman McInnes                MATOP_RESTORE_ROW=2,
14281c1c02c0SLois Curfman McInnes                MATOP_MULT=3,
14291c1c02c0SLois Curfman McInnes                MATOP_MULT_ADD=4,
14307c922b88SBarry Smith                MATOP_MULT_TRANSPOSE=5,
14317c922b88SBarry Smith                MATOP_MULT_TRANSPOSE_ADD=6,
14321c1c02c0SLois Curfman McInnes                MATOP_SOLVE=7,
14331c1c02c0SLois Curfman McInnes                MATOP_SOLVE_ADD=8,
14347c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE=9,
14357c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE_ADD=10,
14361c1c02c0SLois Curfman McInnes                MATOP_LUFACTOR=11,
14371c1c02c0SLois Curfman McInnes                MATOP_CHOLESKYFACTOR=12,
1438*a714c06dSBarry Smith                MATOP_SOR=13,
14391c1c02c0SLois Curfman McInnes                MATOP_TRANSPOSE=14,
14401c1c02c0SLois Curfman McInnes                MATOP_GETINFO=15,
14411c1c02c0SLois Curfman McInnes                MATOP_EQUAL=16,
14421c1c02c0SLois Curfman McInnes                MATOP_GET_DIAGONAL=17,
14431c1c02c0SLois Curfman McInnes                MATOP_DIAGONAL_SCALE=18,
14441c1c02c0SLois Curfman McInnes                MATOP_NORM=19,
14451c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_BEGIN=20,
14461c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_END=21,
1447d519adbfSMatthew Knepley                MATOP_SET_OPTION=22,
1448d519adbfSMatthew Knepley                MATOP_ZERO_ENTRIES=23,
1449d519adbfSMatthew Knepley                MATOP_ZERO_ROWS=24,
1450d519adbfSMatthew Knepley                MATOP_LUFACTOR_SYMBOLIC=25,
1451d519adbfSMatthew Knepley                MATOP_LUFACTOR_NUMERIC=26,
1452d519adbfSMatthew Knepley                MATOP_CHOLESKY_FACTOR_SYMBOLIC=27,
1453d519adbfSMatthew Knepley                MATOP_CHOLESKY_FACTOR_NUMERIC=28,
1454d519adbfSMatthew Knepley                MATOP_SETUP_PREALLOCATION=29,
1455d519adbfSMatthew Knepley                MATOP_ILUFACTOR_SYMBOLIC=30,
1456d519adbfSMatthew Knepley                MATOP_ICCFACTOR_SYMBOLIC=31,
1457d519adbfSMatthew Knepley                MATOP_GET_ARRAY=32,
1458d519adbfSMatthew Knepley                MATOP_RESTORE_ARRAY=33,
1459d519adbfSMatthew Knepley                MATOP_DUPLICATE=34,
1460d519adbfSMatthew Knepley                MATOP_FORWARD_SOLVE=35,
1461d519adbfSMatthew Knepley                MATOP_BACKWARD_SOLVE=36,
1462d519adbfSMatthew Knepley                MATOP_ILUFACTOR=37,
1463d519adbfSMatthew Knepley                MATOP_ICCFACTOR=38,
1464d519adbfSMatthew Knepley                MATOP_AXPY=39,
1465d519adbfSMatthew Knepley                MATOP_GET_SUBMATRICES=40,
1466d519adbfSMatthew Knepley                MATOP_INCREASE_OVERLAP=41,
1467d519adbfSMatthew Knepley                MATOP_GET_VALUES=42,
1468d519adbfSMatthew Knepley                MATOP_COPY=43,
1469d519adbfSMatthew Knepley                MATOP_GET_ROW_MAX=44,
1470d519adbfSMatthew Knepley                MATOP_SCALE=45,
1471d519adbfSMatthew Knepley                MATOP_SHIFT=46,
147235153367SBarry Smith                MATOP_DIAGONAL_SET=47,
1473d519adbfSMatthew Knepley                MATOP_ILUDT_FACTOR=48,
1474d519adbfSMatthew Knepley                MATOP_SET_BLOCK_SIZE=49,
1475d519adbfSMatthew Knepley                MATOP_GET_ROW_IJ=50,
1476d519adbfSMatthew Knepley                MATOP_RESTORE_ROW_IJ=51,
1477d519adbfSMatthew Knepley                MATOP_GET_COLUMN_IJ=52,
1478d519adbfSMatthew Knepley                MATOP_RESTORE_COLUMN_IJ=53,
1479d519adbfSMatthew Knepley                MATOP_FDCOLORING_CREATE=54,
1480d519adbfSMatthew Knepley                MATOP_COLORING_PATCH=55,
1481d519adbfSMatthew Knepley                MATOP_SET_UNFACTORED=56,
1482d519adbfSMatthew Knepley                MATOP_PERMUTE=57,
1483d519adbfSMatthew Knepley                MATOP_SET_VALUES_BLOCKED=58,
1484d519adbfSMatthew Knepley                MATOP_GET_SUBMATRIX=59,
1485d519adbfSMatthew Knepley                MATOP_DESTROY=60,
1486d519adbfSMatthew Knepley                MATOP_VIEW=61,
1487d519adbfSMatthew Knepley                MATOP_CONVERT_FROM=62,
1488d519adbfSMatthew Knepley                MATOP_USE_SCALED_FORM=63,
1489d519adbfSMatthew Knepley                MATOP_SCALE_SYSTEM=64,
1490d519adbfSMatthew Knepley                MATOP_UNSCALE_SYSTEM=65,
1491d519adbfSMatthew Knepley                MATOP_SET_LOCAL_TO_GLOBAL_MAP=66,
1492d519adbfSMatthew Knepley                MATOP_SET_VALUES_LOCAL=67,
1493d519adbfSMatthew Knepley                MATOP_ZERO_ROWS_LOCAL=68,
1494d519adbfSMatthew Knepley                MATOP_GET_ROW_MAX_ABS=69,
1495d519adbfSMatthew Knepley                MATOP_GET_ROW_MIN_ABS=70,
1496d519adbfSMatthew Knepley                MATOP_CONVERT=71,
1497d519adbfSMatthew Knepley                MATOP_SET_COLORING=72,
1498d519adbfSMatthew Knepley                MATOP_SET_VALUES_ADIC=73,
1499d519adbfSMatthew Knepley                MATOP_SET_VALUES_ADIFOR=74,
1500d519adbfSMatthew Knepley                MATOP_FD_COLORING_APPLY=75,
1501d519adbfSMatthew Knepley                MATOP_SET_FROM_OPTIONS=76,
1502d519adbfSMatthew Knepley                MATOP_MULT_CON=77,
1503d519adbfSMatthew Knepley                MATOP_MULT_TRANSPOSE_CON=78,
1504d519adbfSMatthew Knepley                MATOP_PERMUTE_SPARSIFY=79,
1505d519adbfSMatthew Knepley                MATOP_MULT_MULTIPLE=80,
1506d519adbfSMatthew Knepley                MATOP_SOLVE_MULTIPLE=81,
1507d519adbfSMatthew Knepley                MATOP_GET_INERTIA=82,
1508d519adbfSMatthew Knepley                MATOP_LOAD=83,
1509d519adbfSMatthew Knepley                MATOP_IS_SYMMETRIC=84,
1510d519adbfSMatthew Knepley                MATOP_IS_HERMITIAN=85,
1511d519adbfSMatthew Knepley                MATOP_IS_STRUCTURALLY_SYMMETRIC=86,
151241f059aeSBarry Smith                MATOP_DUMMY=87,
1513d519adbfSMatthew Knepley                MATOP_GET_VECS=88,
1514d519adbfSMatthew Knepley                MATOP_MAT_MULT=89,
1515d519adbfSMatthew Knepley                MATOP_MAT_MULT_SYMBOLIC=90,
1516d519adbfSMatthew Knepley                MATOP_MAT_MULT_NUMERIC=91,
1517d519adbfSMatthew Knepley                MATOP_PTAP=92,
1518d519adbfSMatthew Knepley                MATOP_PTAP_SYMBOLIC=93,
1519d519adbfSMatthew Knepley                MATOP_PTAP_NUMERIC=94,
1520d519adbfSMatthew Knepley                MATOP_MAT_MULTTRANSPOSE=95,
15211763ac4eSSatish Balay                MATOP_MAT_MULTTRANSPOSE_SYM=96,
15221763ac4eSSatish Balay                MATOP_MAT_MULTTRANSPOSE_NUM=97,
1523d519adbfSMatthew Knepley                MATOP_PTAP_SYMBOLIC_SEQAIJ=98,
1524d519adbfSMatthew Knepley                MATOP_PTAP_NUMERIC_SEQAIJ=99,
1525d519adbfSMatthew Knepley                MATOP_PTAP_SYMBOLIC_MPIAIJ=100,
1526d519adbfSMatthew Knepley                MATOP_PTAP_NUMERIC_MPIAIJ=101,
1527d519adbfSMatthew Knepley                MATOP_CONJUGATE=102,
1528d519adbfSMatthew Knepley                MATOP_SET_SIZES=103,
1529d519adbfSMatthew Knepley                MATOP_SET_VALUES_ROW=104,
1530d519adbfSMatthew Knepley                MATOP_REAL_PART=105,
1531d519adbfSMatthew Knepley                MATOP_IMAG_PART=106,
1532d519adbfSMatthew Knepley                MATOP_GET_ROW_UTRIANGULAR=107,
1533d519adbfSMatthew Knepley                MATOP_RESTORE_ROW_UTRIANGULAR=108,
1534d519adbfSMatthew Knepley                MATOP_MATSOLVE=109,
1535d519adbfSMatthew Knepley                MATOP_GET_REDUNDANTMATRIX=110,
1536d519adbfSMatthew Knepley                MATOP_GET_ROW_MIN=111,
1537d519adbfSMatthew Knepley                MATOP_GET_COLUMN_VEC=112,
1538d519adbfSMatthew Knepley                MATOP_MISSING_DIAGONAL=113,
1539d519adbfSMatthew Knepley                MATOP_MATGETSEQNONZEROSTRUCTURE=114,
154089c6957cSBarry Smith                MATOP_CREATE=115,
154189c6957cSBarry Smith                MATOP_GET_GHOSTS=116,
154289c6957cSBarry Smith                MATOP_ILUDTFACTORSYMBOLIC=117,
154389c6957cSBarry Smith                MATOP_ILUDTFACTORNUMERIC=118,
1544eeffb40dSHong Zhang                MATOP_MULT_DIAGONAL_BLOCK=119,
1545547795f9SHong Zhang                MATOP_HERMITIANTRANSPOSE=120,
1546547795f9SHong Zhang                MATOP_MULTHERMITIANTRANSPOSE=121,
1547547795f9SHong Zhang                MATOP_MULTHERMITIANTRANSPOSEADD=122
1548fae171e0SBarry Smith              } MatOperation;
1549be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatHasOperation(Mat,MatOperation,PetscTruth*);
1550be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellSetOperation(Mat,MatOperation,void(*)(void));
1551be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellGetOperation(Mat,MatOperation,void(**)(void));
1552be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellSetContext(Mat,void*);
1553112a2221SBarry Smith 
155490ace30eSBarry Smith /*
155590ace30eSBarry Smith    Codes for matrices stored on disk. By default they are
155690ace30eSBarry Smith  stored in a universal format. By changing the format with
15577973ad04SBarry Smith  PetscViewerSetFormat(viewer,PETSC_VIEWER_NATIVE); the matrices will
155890ace30eSBarry Smith  be stored in a way natural for the matrix, for example dense matrices
155990ace30eSBarry Smith  would be stored as dense. Matrices stored this way may only be
156090ace30eSBarry Smith  read into matrices of the same time.
156190ace30eSBarry Smith */
156290ace30eSBarry Smith #define MATRIX_BINARY_FORMAT_DENSE -1
156390ace30eSBarry Smith 
1564be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetHashTableFactor(Mat,PetscReal);
1565be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatISGetLocalMat(Mat,Mat*);
15661f4e1ec7SBarry Smith 
1567d9274352SBarry Smith /*S
1568d9274352SBarry Smith      MatNullSpace - Object that removes a null space from a vector, i.e.
1569d9274352SBarry Smith          orthogonalizes the vector to a subsapce
1570d9274352SBarry Smith 
1571f7a9e4ceSBarry Smith    Level: advanced
1572d9274352SBarry Smith 
1573d9274352SBarry Smith   Concepts: matrix; linear operator, null space
1574d9274352SBarry Smith 
15756e1639daSBarry Smith   Users manual sections:
15766e1639daSBarry Smith .   sec_singular
15776e1639daSBarry Smith 
1578d9274352SBarry Smith .seealso:  MatNullSpaceCreate()
1579d9274352SBarry Smith S*/
158074637425SBarry Smith typedef struct _p_MatNullSpace* MatNullSpace;
1581d9274352SBarry Smith 
1582be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceCreate(MPI_Comm,PetscTruth,PetscInt,const Vec[],MatNullSpace*);
1583b22b330cSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(MatNullSpace,Vec,void*),void*);
1584be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceDestroy(MatNullSpace);
1585be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceRemove(MatNullSpace,Vec,Vec*);
1586be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceAttach(Mat,MatNullSpace);
158795902228SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceTest(MatNullSpace,Mat,PetscTruth *);
158874637425SBarry Smith 
1589be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatReorderingSeqSBAIJ(Mat,IS);
1590be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPISBAIJSetHashTableFactor(Mat,PetscReal);
1591be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetColumnIndices(Mat,PetscInt *);
159246129b97SKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJInvertBlockDiagonal(Mat);
15933f1d51d7SBarry Smith 
1594be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMAIJ(Mat,PetscInt,Mat*);
1595be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMAIJRedimension(Mat,PetscInt,Mat*);
1596be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMAIJGetAIJ(Mat,Mat*);
1597c4f061fbSSatish Balay 
1598be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatComputeExplicitOperator(Mat,Mat*);
1599b0a32e0cSBarry Smith 
1600be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalScaleLocal(Mat,Vec);
160104f1ad80SBarry Smith 
1602fef1beadSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*);
16033ec795f1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetBase(Mat,Vec,Vec);
16043ec795f1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*);
1605e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*));
1606e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec));
1607e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDAddNullSpace(Mat,MatNullSpace);
1608e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt);
1609e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDResetHHistory(Mat);
1610e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctionError(Mat,PetscReal);
1611e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetPeriod(Mat,PetscInt);
1612e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDGetH(Mat,PetscScalar *);
16136aa9148fSLisandro Dalcin EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetOptionsPrefix(Mat,const char[]);
1614e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFromOptions(Mat);
1615e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*);
1616e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*);
1617e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDDSSetUmin(Mat,PetscReal);
1618e884886fSBarry Smith 
1619e884886fSBarry Smith 
1620e884886fSBarry Smith typedef struct _p_MatMFFD* MatMFFD;
1621e884886fSBarry Smith 
1622e884886fSBarry Smith /*E
1623e884886fSBarry Smith     MatMFFDType - algorithm used to compute the h used in computing matrix-vector products via differencing of the function
1624e884886fSBarry Smith 
1625e884886fSBarry Smith    Level: beginner
1626e884886fSBarry Smith 
1627e884886fSBarry Smith .seealso: MatMFFDSetType(), MatMFFDRegister()
1628e884886fSBarry Smith E*/
1629a313700dSBarry Smith #define MatMFFDType char*
1630e884886fSBarry Smith #define MATMFFD_DS  "ds"
1631e884886fSBarry Smith #define MATMFFD_WP  "wp"
1632e884886fSBarry Smith 
1633a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetType(Mat,const MatMFFDType);
1634e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatMFFD));
1635e884886fSBarry Smith 
1636e884886fSBarry Smith /*MC
1637e884886fSBarry Smith    MatMFFDRegisterDynamic - Adds a method to the MatMFFD registry.
1638e884886fSBarry Smith 
1639e884886fSBarry Smith    Synopsis:
1640e884886fSBarry Smith    PetscErrorCode MatMFFDRegisterDynamic(char *name_solver,char *path,char *name_create,PetscErrorCode (*routine_create)(MatMFFD))
1641e884886fSBarry Smith 
1642e884886fSBarry Smith    Not Collective
1643e884886fSBarry Smith 
1644e884886fSBarry Smith    Input Parameters:
1645e884886fSBarry Smith +  name_solver - name of a new user-defined compute-h module
1646e884886fSBarry Smith .  path - path (either absolute or relative) the library containing this solver
1647e884886fSBarry Smith .  name_create - name of routine to create method context
1648e884886fSBarry Smith -  routine_create - routine to create method context
1649e884886fSBarry Smith 
1650e884886fSBarry Smith    Level: developer
1651e884886fSBarry Smith 
1652e884886fSBarry Smith    Notes:
1653e884886fSBarry Smith    MatMFFDRegisterDynamic() may be called multiple times to add several user-defined solvers.
1654e884886fSBarry Smith 
1655e884886fSBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
1656e884886fSBarry Smith    is ignored.
1657e884886fSBarry Smith 
1658e884886fSBarry Smith    Sample usage:
1659e884886fSBarry Smith .vb
1660e884886fSBarry Smith    MatMFFDRegisterDynamic("my_h",/home/username/my_lib/lib/libO/solaris/mylib.a,
1661e884886fSBarry Smith                "MyHCreate",MyHCreate);
1662e884886fSBarry Smith .ve
1663e884886fSBarry Smith 
1664e884886fSBarry Smith    Then, your solver can be chosen with the procedural interface via
1665e884886fSBarry Smith $     MatMFFDSetType(mfctx,"my_h")
1666e884886fSBarry Smith    or at runtime via the option
1667e884886fSBarry Smith $     -snes_mf_type my_h
1668e884886fSBarry Smith 
1669e884886fSBarry Smith .keywords: MatMFFD, register
1670e884886fSBarry Smith 
1671e884886fSBarry Smith .seealso: MatMFFDRegisterAll(), MatMFFDRegisterDestroy()
1672e884886fSBarry Smith M*/
1673e884886fSBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1674e884886fSBarry Smith #define MatMFFDRegisterDynamic(a,b,c,d) MatMFFDRegister(a,b,c,0)
1675e884886fSBarry Smith #else
1676e884886fSBarry Smith #define MatMFFDRegisterDynamic(a,b,c,d) MatMFFDRegister(a,b,c,d)
1677e884886fSBarry Smith #endif
1678e884886fSBarry Smith 
1679e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDRegisterAll(const char[]);
1680e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDRegisterDestroy(void);
1681e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDDefaultSetUmin(Mat,PetscReal);
1682e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDWPSetComputeNormU(Mat,PetscTruth);
1683e884886fSBarry Smith 
1684e884886fSBarry Smith 
1685be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *);
1686be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *);
16877dbadf16SMatthew Knepley 
168897969023SHong Zhang /*
168997969023SHong Zhang    PETSc interface to MUMPS
169097969023SHong Zhang */
169197969023SHong Zhang #ifdef PETSC_HAVE_MUMPS
169297969023SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMumpsSetIcntl(Mat,PetscInt,PetscInt);
169397969023SHong Zhang #endif
169423a5497aSJed Brown 
169523a5497aSJed Brown PETSC_EXTERN_CXX_END
169623a5497aSJed Brown #endif
1697