xref: /petsc/include/petscmat.h (revision 599ef60d4b3dadfc876bc133d4ddc9cfd927b79c)
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 MATMPIROWBS        "mpirowbs"
36273d9f13SBarry Smith #define MATSEQAIJ          "seqaij"
37273d9f13SBarry Smith #define MATMPIAIJ          "mpiaij"
38209238afSKris Buschelman #define MATAIJ             "aij"
39273d9f13SBarry Smith #define MATSHELL           "shell"
40209238afSKris Buschelman #define MATSEQDENSE        "seqdense"
41273d9f13SBarry Smith #define MATMPIDENSE        "mpidense"
42209238afSKris Buschelman #define MATDENSE           "dense"
43273d9f13SBarry Smith #define MATSEQBAIJ         "seqbaij"
44273d9f13SBarry Smith #define MATMPIBAIJ         "mpibaij"
45209238afSKris Buschelman #define MATBAIJ            "baij"
46273d9f13SBarry Smith #define MATMPIADJ          "mpiadj"
47273d9f13SBarry Smith #define MATSEQSBAIJ        "seqsbaij"
48273d9f13SBarry Smith #define MATMPISBAIJ        "mpisbaij"
49209238afSKris Buschelman #define MATSBAIJ           "sbaij"
50cebc7f6cSBarry Smith #define MATDAAD            "daad"
51cebc7f6cSBarry Smith #define MATMFFD            "mffd"
52c8a8475eSBarry Smith #define MATNORMAL          "normal"
53ab92ecdeSBarry Smith #define MATLRC             "lrc"
544b8d99adSRichard Tran Mills #define MATSEQCSRPERM      "seqcsrperm"
5518def556SRichard Tran Mills #define MATMPICSRPERM      "mpicsrperm"
56814655b8SBarry Smith #define MATCSRPERM         "csrperm"
5781824310SBarry Smith #define MATSEQCRL          "seqcrl"
58c02b24c5SBarry Smith #define MATMPICRL          "mpicrl"
59c02b24c5SBarry Smith #define MATCRL             "crl"
602a6744ebSBarry Smith #define MATSCATTER         "scatter"
61421e10b8SBarry Smith #define MATBLOCKMAT        "blockmat"
62793850ffSBarry Smith #define MATCOMPOSITE       "composite"
635a7f1df3SHong Zhang #define MATSEQFFTW         "seqfftw"
64557cca28SSatish Balay #define MATTRANSPOSEMAT    "transpose"
6572ca8751SBarry Smith #define MATSCHURCOMPLEMENT "schurcomplement"
661d6018f0SLisandro Dalcin #define MATPYTHON          "python"
67d91e6319SBarry Smith 
689989ab13SBarry Smith /*E
69c7393fdbSBarry Smith     MatSolverPackage - String with the name of a PETSc matrix solver type.
709989ab13SBarry Smith 
719989ab13SBarry Smith     For example: "petsc" indicates what PETSc provides, "superlu" indicates either
729989ab13SBarry Smith        SuperLU or SuperLU_Dist etc.
739989ab13SBarry Smith 
749989ab13SBarry Smith 
759989ab13SBarry Smith    Level: beginner
769989ab13SBarry Smith 
775c9eb25fSBarry Smith .seealso: MatGetFactor(), Mat, MatSetType(), MatType
789989ab13SBarry Smith E*/
79c7393fdbSBarry Smith #define MatSolverPackage char*
80b24902e0SBarry Smith #define MAT_SOLVER_SPOOLES      "spooles"
81b24902e0SBarry Smith #define MAT_SOLVER_SUPERLU      "superlu"
82b24902e0SBarry Smith #define MAT_SOLVER_SUPERLU_DIST "superlu_dist"
83b24902e0SBarry Smith #define MAT_SOLVER_UMFPACK      "umfpack"
84b24902e0SBarry Smith #define MAT_SOLVER_ESSL         "essl"
85b24902e0SBarry Smith #define MAT_SOLVER_LUSOL        "lusol"
86b24902e0SBarry Smith #define MAT_SOLVER_MUMPS        "mumps"
873bf14a46SMatthew Knepley #define MAT_SOLVER_PASTIX       "pastix"
88b24902e0SBarry Smith #define MAT_SOLVER_DSCPACK      "dscpack"
89b24902e0SBarry Smith #define MAT_SOLVER_MATLAB       "matlab"
9043244d56SBarry Smith #define MAT_SOLVER_PETSC        "petsc"
91b6806ab0SHong Zhang #define MAT_SOLVER_PLAPACK      "plapack"
92b24902e0SBarry Smith 
93b24902e0SBarry Smith /*E
94b24902e0SBarry Smith     MatFactorType - indicates what type of factorization is requested
95b24902e0SBarry Smith 
96b24902e0SBarry Smith     Level: beginner
97b24902e0SBarry Smith 
98b24902e0SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
99b24902e0SBarry Smith 
100c7393fdbSBarry Smith .seealso: MatSolverPackage, MatGetFactor()
101b24902e0SBarry Smith E*/
102*599ef60dSHong Zhang typedef enum {MAT_FACTOR_NONE, MAT_FACTOR_LU, MAT_FACTOR_CHOLESKY, MAT_FACTOR_ILU, MAT_FACTOR_ICC,MAT_FACTOR_ILUDT} MatFactorType;
103c7393fdbSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetFactor(Mat,const MatSolverPackage,MatFactorType,Mat*);
104db4efbfdSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetFactorAvailable(Mat,const MatSolverPackage,MatFactorType,PetscTruth*);
10535bd34faSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFactorGetSolverPackage(Mat,const MatSolverPackage*);
106b24902e0SBarry Smith 
1079989ab13SBarry Smith 
108c06d978dSMatthew Knepley /* Logging support */
109552e946dSBarry Smith #define    MAT_FILE_COOKIE 1211216    /* used to indicate matrices in binary files */
110be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_COOKIE;
111be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_FDCOLORING_COOKIE;
112be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_PARTITIONING_COOKIE;
113be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_NULLSPACE_COOKIE;
114e884886fSBarry Smith extern PetscCookie PETSCMAT_DLLEXPORT MATMFFD_COOKIE;
115c06d978dSMatthew Knepley 
116ceb03754SKris Buschelman /*E
117ceb03754SKris Buschelman     MatReuse - Indicates if matrices obtained from a previous call to MatGetSubMatrices()
118ceb03754SKris Buschelman      or MatGetSubMatrix() are to be reused to store the new matrix values.
119ceb03754SKris Buschelman 
120ceb03754SKris Buschelman     Level: beginner
121ceb03754SKris Buschelman 
122ceb03754SKris Buschelman    Any additions/changes here MUST also be made in include/finclude/petscmat.h
123ceb03754SKris Buschelman 
1240c451bc4SBarry Smith .seealso: MatGetSubMatrices(), MatGetSubMatrix(), MatDestroyMatrices(), MatConvert()
125ceb03754SKris Buschelman E*/
126ceb03754SKris Buschelman typedef enum {MAT_INITIAL_MATRIX,MAT_REUSE_MATRIX} MatReuse;
127ceb03754SKris Buschelman 
1285494a064SHong Zhang /*E
1295494a064SHong Zhang     MatGetSubMatrixOption - Indicates if matrices obtained from a call to MatGetSubMatrices()
130829201f2SHong Zhang      include the matrix values. Currently it is only used by MatGetSeqNonzerostructure().
1315494a064SHong Zhang 
1325494a064SHong Zhang     Level: beginner
1335494a064SHong Zhang 
134829201f2SHong Zhang .seealso: MatGetSeqNonzerostructure()
1355494a064SHong Zhang E*/
1365494a064SHong Zhang typedef enum {MAT_DO_NOT_GET_VALUES,MAT_GET_VALUES} MatGetSubMatrixOption;
1375494a064SHong Zhang 
138e5bd5246SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInitializePackage(const char[]);
139c06d978dSMatthew Knepley 
140f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreate(MPI_Comm,Mat*);
141a69a960dSMatthew Knepley PetscPolymorphicFunction(MatCreate,(MPI_Comm comm),(comm,&A),Mat,A)
142a69a960dSMatthew Knepley PetscPolymorphicFunction(MatCreate,(),(PETSC_COMM_WORLD,&A),Mat,A)
143f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetSizes(Mat,PetscInt,PetscInt,PetscInt,PetscInt);
144a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetType(Mat,const MatType);
145be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetFromOptions(Mat);
146be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUpPreallocation(Mat);
147be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRegisterAll(const char[]);
148be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat));
14930de9b25SBarry Smith 
150f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetOptionsPrefix(Mat,const char[]);
151f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAppendOptionsPrefix(Mat,const char[]);
152f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOptionsPrefix(Mat,const char*[]);
153f69a0ea3SMatthew Knepley 
15430de9b25SBarry Smith /*MC
15530de9b25SBarry Smith    MatRegisterDynamic - Adds a new matrix type
15630de9b25SBarry Smith 
15730de9b25SBarry Smith    Synopsis:
158c1ac3661SBarry Smith    PetscErrorCode MatRegisterDynamic(char *name,char *path,char *name_create,PetscErrorCode (*routine_create)(Mat))
15930de9b25SBarry Smith 
16030de9b25SBarry Smith    Not Collective
16130de9b25SBarry Smith 
16230de9b25SBarry Smith    Input Parameters:
16330de9b25SBarry Smith +  name - name of a new user-defined matrix type
16430de9b25SBarry Smith .  path - path (either absolute or relative) the library containing this solver
16530de9b25SBarry Smith .  name_create - name of routine to create method context
16630de9b25SBarry Smith -  routine_create - routine to create method context
16730de9b25SBarry Smith 
16830de9b25SBarry Smith    Notes:
16930de9b25SBarry Smith    MatRegisterDynamic() may be called multiple times to add several user-defined solvers.
17030de9b25SBarry Smith 
17130de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
17230de9b25SBarry Smith    is ignored.
17330de9b25SBarry Smith 
17430de9b25SBarry Smith    Sample usage:
17530de9b25SBarry Smith .vb
17630de9b25SBarry Smith    MatRegisterDynamic("my_mat",/home/username/my_lib/lib/libO/solaris/mylib.a,
17730de9b25SBarry Smith                "MyMatCreate",MyMatCreate);
17830de9b25SBarry Smith .ve
17930de9b25SBarry Smith 
18030de9b25SBarry Smith    Then, your solver can be chosen with the procedural interface via
18130de9b25SBarry Smith $     MatSetType(Mat,"my_mat")
18230de9b25SBarry Smith    or at runtime via the option
18330de9b25SBarry Smith $     -mat_type my_mat
18430de9b25SBarry Smith 
18530de9b25SBarry Smith    Level: advanced
18630de9b25SBarry Smith 
187ab901514SBarry Smith    Notes: ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values.
18830de9b25SBarry Smith          If your function is not being put into a shared library then use VecRegister() instead
18930de9b25SBarry Smith 
19030de9b25SBarry Smith .keywords: Mat, register
19130de9b25SBarry Smith 
19230de9b25SBarry Smith .seealso: MatRegisterAll(), MatRegisterDestroy()
19330de9b25SBarry Smith 
19430de9b25SBarry Smith M*/
195273d9f13SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
196273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,0)
197273d9f13SBarry Smith #else
198273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,d)
19930de9b25SBarry Smith #endif
20030de9b25SBarry Smith 
201273d9f13SBarry Smith extern PetscTruth MatRegisterAllCalled;
202b0a32e0cSBarry Smith extern PetscFList MatList;
20328988994SBarry Smith 
2043b224e63SBarry Smith /*E
2053b224e63SBarry Smith     MatStructure - Indicates if the matrix has the same nonzero structure
2063b224e63SBarry Smith 
2073b224e63SBarry Smith     Level: beginner
2083b224e63SBarry Smith 
2093b224e63SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
2103b224e63SBarry Smith 
2113b224e63SBarry Smith .seealso: MatCopy(), KSPSetOperators(), PCSetOperators()
2123b224e63SBarry Smith E*/
2133b224e63SBarry Smith typedef enum {SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER,SUBSET_NONZERO_PATTERN} MatStructure;
2143b224e63SBarry Smith 
215be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqDense(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*);
216be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIDense(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*);
217be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
218ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,nz,nnz,&A),Mat,A)
219ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,&A),Mat,A)
220ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,0,nnz,&A),Mat,A)
221ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,&A),Mat,A)
222ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,A))
223ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,m,n,0,nnz,A))
224ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,A))
225be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
226ba966639SSatish 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)
227ba966639SSatish 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)
228ba966639SSatish 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)
229ba966639SSatish 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)
230ba966639SSatish 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))
231ba966639SSatish 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))
232ba966639SSatish 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))
233ba966639SSatish 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)
234ba966639SSatish 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)
235ba966639SSatish 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)
236ba966639SSatish 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)
237ba966639SSatish 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))
238ba966639SSatish 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))
239ba966639SSatish 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))
24063ab93bfSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *);
2418d7a6e47SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJWithSplitArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],PetscInt[],PetscInt[],PetscScalar[],Mat*);
2428d7a6e47SBarry Smith 
243be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIRowbs(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
244be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
245ba966639SSatish 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)
246ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A)
247ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A)
248ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A)
249ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A))
250ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A))
251ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A))
252be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
253ba966639SSatish 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)
254ba966639SSatish 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)
255ba966639SSatish 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)
256ba966639SSatish 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)
257ba966639SSatish 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))
258ba966639SSatish 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))
259ba966639SSatish 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))
260ba966639SSatish 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)
261ba966639SSatish 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)
262ba966639SSatish 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)
263ba966639SSatish 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)
264ba966639SSatish 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))
265ba966639SSatish 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))
266ba966639SSatish 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))
267be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAdj(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscInt[],Mat*);
268be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
269ba966639SSatish 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)
270ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A)
271ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A)
272ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A)
273ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A))
274ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A))
275ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A))
276be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPISBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
277ba966639SSatish 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)
278ba966639SSatish 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)
279ba966639SSatish 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)
280ba966639SSatish 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)
281ba966639SSatish 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))
282ba966639SSatish 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))
283ba966639SSatish 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))
284ba966639SSatish 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)
285ba966639SSatish 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)
286ba966639SSatish 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)
287ba966639SSatish 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)
288ba966639SSatish 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))
289ba966639SSatish 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))
290ba966639SSatish 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))
291be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateShell(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,void *,Mat*);
292ba966639SSatish 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)
293ba966639SSatish Balay PetscPolymorphicFunction(MatCreateShell,(PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx),(PETSC_COMM_WORLD,m,n,M,N,ctx,&A),Mat,A)
294be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateAdic(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,void (*)(void),Mat*);
295be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateNormal(Mat,Mat*);
296ba966639SSatish Balay PetscPolymorphicFunction(MatCreateNormal,(Mat mat),(mat,&A),Mat,A)
297ab92ecdeSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateLRC(Mat,Mat,Mat,Mat*);
298284134d9SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateIS(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,ISLocalToGlobalMapping,Mat*);
2991472f72bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
3001472f72bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPICRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
3012a6744ebSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateScatter(MPI_Comm,VecScatter,Mat*);
3022a6744ebSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScatterSetVecScatter(Mat,VecScatter);
3032a6744ebSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScatterGetVecScatter(Mat,VecScatter*);
3048cdf2d9bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateBlockMat(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,Mat*);
305793850ffSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompositeAddMat(Mat,Mat);
306b52f573bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompositeMerge(Mat);
307793850ffSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateComposite(MPI_Comm,PetscInt,const Mat*,Mat*);
3086d7c1e57SBarry Smith typedef enum {MAT_COMPOSITE_ADDITIVE,MAT_COMPOSITE_MULTIPLICATIVE} MatCompositeType;
3096d7c1e57SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompositeSetType(Mat,MatCompositeType);
3106d7c1e57SBarry Smith 
3115a7f1df3SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqFFTW(MPI_Comm,PetscInt,const PetscInt[],Mat*);
312f20085c4SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateTranspose(Mat,Mat*);
3131472f72bSBarry Smith 
3141d6018f0SLisandro Dalcin EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreatePython(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const char[],Mat*);
3151d6018f0SLisandro Dalcin EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPythonSetType(Mat,const char[]);
3161d6018f0SLisandro Dalcin 
3171d6018f0SLisandro Dalcin 
318f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUp(Mat);
319be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroy(Mat);
32021c89e3eSBarry Smith 
3210c451bc4SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConjugate(Mat);
32299cafbc1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRealPart(Mat);
32399cafbc1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatImaginaryPart(Mat);
324bfc7d00aSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetDiagonalBlock(Mat,PetscTruth*,MatReuse,Mat*);
32599cafbc1SBarry Smith 
3268ed539a5SBarry Smith /* ------------------------------------------------------------*/
327be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
328be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
32987d4246cSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesRow(Mat,PetscInt,const PetscScalar[]);
330f259bd47SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesRowLocal(Mat,PetscInt,const PetscScalar[]);
33184cb2905SBarry Smith 
3322ef4de8bSBarry Smith /*S
3332ef4de8bSBarry Smith      MatStencil - Data structure (C struct) for storing information about a single row or
3342ef4de8bSBarry Smith         column of a matrix as index on an associated grid.
3352ef4de8bSBarry Smith 
3362ef4de8bSBarry Smith    Level: beginner
3372ef4de8bSBarry Smith 
3382ef4de8bSBarry Smith   Concepts: matrix; linear operator
3392ef4de8bSBarry Smith 
340d7bd93baSBarry Smith .seealso:  MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockStencil()
3412ef4de8bSBarry Smith S*/
342435da068SBarry Smith typedef struct {
343c1ac3661SBarry Smith   PetscInt k,j,i,c;
344435da068SBarry Smith } MatStencil;
3452ef4de8bSBarry Smith 
346be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
347be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
348be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt);
349435da068SBarry Smith 
350be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetColoring(Mat,ISColoring);
351be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesAdic(Mat,void*);
352be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesAdifor(Mat,PetscInt,void*);
3533a7fca6bSBarry Smith 
354d91e6319SBarry Smith /*E
355d91e6319SBarry Smith     MatAssemblyType - Indicates if the matrix is now to be used, or if you plan
356d91e6319SBarry Smith      to continue to add values to it
357d91e6319SBarry Smith 
358d91e6319SBarry Smith     Level: beginner
359d91e6319SBarry Smith 
360d91e6319SBarry Smith .seealso: MatAssemblyBegin(), MatAssemblyEnd()
361d91e6319SBarry Smith E*/
3626d4a8577SBarry Smith typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType;
363be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssemblyBegin(Mat,MatAssemblyType);
364be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssemblyEnd(Mat,MatAssemblyType);
365be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssembled(Mat,PetscTruth*);
3664f9c727eSBarry Smith 
3671d73ed98SBarry Smith 
36830de9b25SBarry Smith 
369d91e6319SBarry Smith /*E
370d91e6319SBarry Smith     MatOption - Options that may be set for a matrix and its behavior or storage
371d91e6319SBarry Smith 
372d91e6319SBarry Smith     Level: beginner
373d91e6319SBarry Smith 
3740a835dfdSSatish Balay    Any additions/changes here MUST also be made in include/finclude/petscmat.h
375d91e6319SBarry Smith 
376d91e6319SBarry Smith .seealso: MatSetOption()
377d91e6319SBarry Smith E*/
378512a5fc5SBarry Smith typedef enum {MAT_ROW_ORIENTED,MAT_NEW_NONZERO_LOCATIONS,
3794e0d8c25SBarry Smith               MAT_SYMMETRIC,
3804e0d8c25SBarry Smith               MAT_STRUCTURALLY_SYMMETRIC,
3818047e77bSBarry Smith               MAT_NEW_DIAGONALS,MAT_IGNORE_OFF_PROC_ENTRIES,
3824e0d8c25SBarry Smith               MAT_NEW_NONZERO_LOCATION_ERR,
3834e0d8c25SBarry Smith               MAT_NEW_NONZERO_ALLOCATION_ERR,MAT_USE_HASH_TABLE,
3844e0d8c25SBarry Smith               MAT_KEEP_ZEROED_ROWS,MAT_IGNORE_ZERO_ENTRIES,
3854e0d8c25SBarry Smith               MAT_USE_INODES,
3864e0d8c25SBarry Smith               MAT_HERMITIAN,
3874e0d8c25SBarry Smith               MAT_SYMMETRY_ETERNAL,
3884e0d8c25SBarry Smith               MAT_USE_COMPRESSEDROW,
3894e0d8c25SBarry Smith               MAT_IGNORE_LOWER_TRIANGULAR,MAT_ERROR_LOWER_TRIANGULAR,
39028b2fa4aSMatthew Knepley               MAT_GETROW_UPPERTRIANGULAR,MAT_UNUSED_NONZERO_LOCATION_ERR,
39128b2fa4aSMatthew Knepley               NUM_MAT_OPTIONS} MatOption;
392290bbb0aSBarry Smith extern const char *MatOptions[];
3934e0d8c25SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetOption(Mat,MatOption,PetscTruth);
394a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetType(Mat,const MatType*);
395a313700dSBarry Smith PetscPolymorphicFunction(MatGetType,(Mat mat),(mat,&t),const MatType,t)
39684cb2905SBarry Smith 
397be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]);
398be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
399be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
400f5edf698SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowUpperTriangular(Mat);
401f5edf698SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRowUpperTriangular(Mat);
402be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
403be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
404be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumnVector(Mat,Vec,PetscInt);
405be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetArray(Mat,PetscScalar *[]);
406ba966639SSatish Balay PetscPolymorphicFunction(MatGetArray,(Mat mat),(mat,&a),PetscScalar*,a)
407be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreArray(Mat,PetscScalar *[]);
408be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBlockSize(Mat,PetscInt *);
409ba966639SSatish Balay PetscPolymorphicFunction(MatGetBlockSize,(Mat mat),(mat,&a),PetscInt,a)
410be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetBlockSize(Mat,PetscInt);
4117b80b807SBarry Smith 
4121620fd73SBarry Smith 
413be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMult(Mat,Vec,Vec);
414be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultAdd(Mat,Vec,Vec,Vec);
415ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
416be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTranspose(Mat,Vec,Vec);
417be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsTranspose(Mat,Mat,PetscReal,PetscTruth*);
418ba966639SSatish Balay PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B,PetscReal tol),(A,B,tol,&t),PetscTruth,t)
419e37f64cdSBarry Smith PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B),(A,B,0,&t),PetscTruth,t)
4201cbb95d3SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitianTranspose(Mat,Mat,PetscReal,PetscTruth*);
421be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeAdd(Mat,Vec,Vec,Vec);
422ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultTransposeAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
423be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultConstrained(Mat,Vec,Vec);
424be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeConstrained(Mat,Vec,Vec);
4252bebee5dSHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatSolve(Mat,Mat,Mat);
426f5edf698SHong Zhang 
427d91e6319SBarry Smith /*E
428d91e6319SBarry Smith     MatDuplicateOption - Indicates if a duplicated sparse matrix should have
429d91e6319SBarry Smith   its numerical values copied over or just its nonzero structure.
430d91e6319SBarry Smith 
431d91e6319SBarry Smith     Level: beginner
432d91e6319SBarry Smith 
433d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
434d91e6319SBarry Smith 
435d91e6319SBarry Smith .seealso: MatDuplicate()
436d91e6319SBarry Smith E*/
4372e8a6d31SBarry Smith typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES} MatDuplicateOption;
4382e8a6d31SBarry Smith 
439a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert(Mat,const MatType,MatReuse,Mat*);
440a313700dSBarry Smith PetscPolymorphicFunction(MatConvert,(Mat A,const MatType t),(A,t,MAT_INITIAL_MATRIX,&a),Mat,a)
441be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDuplicate(Mat,MatDuplicateOption,Mat*);
442ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A,MatDuplicateOption o),(A,o,&a),Mat,a)
443ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A),(A,MAT_COPY_VALUES,&a),Mat,a)
44494a9d846SBarry Smith 
445cb5b572fSBarry Smith 
446be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCopy(Mat,Mat,MatStructure);
447be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatView(Mat,PetscViewer);
448be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsSymmetric(Mat,PetscReal,PetscTruth*);
449ba966639SSatish Balay PetscPolymorphicFunction(MatIsSymmetric,(Mat A,PetscReal tol),(A,tol,&t),PetscTruth,t)
450e37f64cdSBarry Smith PetscPolymorphicFunction(MatIsSymmetric,(Mat A),(A,0,&t),PetscTruth,t)
451be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsStructurallySymmetric(Mat,PetscTruth*);
452ba966639SSatish Balay PetscPolymorphicFunction(MatIsStructurallySymmetric,(Mat A),(A,&t),PetscTruth,t)
4531cbb95d3SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitian(Mat,PetscReal,PetscTruth*);
4541cbb95d3SBarry Smith PetscPolymorphicFunction(MatIsHermitian,(Mat A),(A,0,&t),PetscTruth,t)
455be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsSymmetricKnown(Mat,PetscTruth*,PetscTruth*);
456be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitianKnown(Mat,PetscTruth*,PetscTruth*);
45764c62002SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMissingDiagonal(Mat,PetscTruth *,PetscInt *);
458a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLoad(PetscViewer,const MatType,Mat*);
459a313700dSBarry Smith PetscPolymorphicFunction(MatLoad,(PetscViewer v,const MatType t),(v,t,&a),Mat,a)
4607b80b807SBarry Smith 
4618f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *);
4628f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRowIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *);
4638f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumnIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *);
4648f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreColumnIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *);
465d4fbbf0eSBarry Smith 
466d91e6319SBarry Smith /*S
467d91e6319SBarry Smith      MatInfo - Context of matrix information, used with MatGetInfo()
468d91e6319SBarry Smith 
469d91e6319SBarry Smith    In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE
470d91e6319SBarry Smith 
471d91e6319SBarry Smith    Level: intermediate
472d91e6319SBarry Smith 
473d91e6319SBarry Smith   Concepts: matrix^nonzero information
474d91e6319SBarry Smith 
475d9274352SBarry Smith .seealso:  MatGetInfo(), MatInfoType
476d91e6319SBarry Smith S*/
4774e220ebcSLois Curfman McInnes typedef struct {
478b0a32e0cSBarry Smith   PetscLogDouble block_size;                         /* block size */
479b0a32e0cSBarry Smith   PetscLogDouble nz_allocated,nz_used,nz_unneeded;   /* number of nonzeros */
480b0a32e0cSBarry Smith   PetscLogDouble memory;                             /* memory allocated */
481b0a32e0cSBarry Smith   PetscLogDouble assemblies;                         /* number of matrix assemblies called */
482b0a32e0cSBarry Smith   PetscLogDouble mallocs;                            /* number of mallocs during MatSetValues() */
483b0a32e0cSBarry Smith   PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */
484b0a32e0cSBarry Smith   PetscLogDouble factor_mallocs;                     /* number of mallocs during factorization */
4854e220ebcSLois Curfman McInnes } MatInfo;
4864e220ebcSLois Curfman McInnes 
487d9274352SBarry Smith /*E
488d9274352SBarry Smith     MatInfoType - Indicates if you want information about the local part of the matrix,
489d9274352SBarry Smith      the entire parallel matrix or the maximum over all the local parts.
490d9274352SBarry Smith 
491d9274352SBarry Smith     Level: beginner
492d9274352SBarry Smith 
493d9274352SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
494d9274352SBarry Smith 
495d9274352SBarry Smith .seealso: MatGetInfo(), MatInfo
496d9274352SBarry Smith E*/
4977b80b807SBarry Smith typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType;
498be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetInfo(Mat,MatInfoType,MatInfo*);
499be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatValid(Mat,PetscTruth*);
500be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetDiagonal(Mat,Vec);
501985db425SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMax(Mat,Vec,PetscInt[]);
502985db425SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMin(Mat,Vec,PetscInt[]);
503985db425SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMaxAbs(Mat,Vec,PetscInt[]);
504c87e5d42SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMinAbs(Mat,Vec,PetscInt[]);
50586697f06SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowSum(Mat,Vec);
506fc4dec0aSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatTranspose(Mat,MatReuse,Mat*);
507fc4dec0aSBarry Smith PetscPolymorphicFunction(MatTranspose,(Mat A),(A,MAT_INITIAL_MATRIX,&t),Mat,t)
508be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPermute(Mat,IS,IS,Mat *);
509ba966639SSatish Balay PetscPolymorphicFunction(MatPermute,(Mat A,IS is1,IS is2),(A,is1,is2,&t),Mat,t)
510be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPermuteSparsify(Mat,PetscInt,PetscReal,PetscReal,IS,IS,Mat *);
511be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalScale(Mat,Vec,Vec);
512be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalSet(Mat,Vec,InsertMode);
513be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatEqual(Mat,Mat,PetscTruth*);
514ba966639SSatish Balay PetscPolymorphicFunction(MatEqual,(Mat A,Mat B),(A,B,&t),PetscTruth,t)
515be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultEqual(Mat,Mat,PetscInt,PetscTruth*);
516be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultAddEqual(Mat,Mat,PetscInt,PetscTruth*);
517be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeEqual(Mat,Mat,PetscInt,PetscTruth*);
518be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscTruth*);
5197b80b807SBarry Smith 
520be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNorm(Mat,NormType,PetscReal *);
521ba966639SSatish Balay PetscPolymorphicFunction(MatNorm,(Mat A,NormType t),(A,t,&n),PetscReal,n)
522be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroEntries(Mat);
523f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar);
524f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsIS(Mat,IS,PetscScalar);
525f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroColumns(Mat,PetscInt,const PetscInt [],const PetscScalar*);
526f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroColumnsIS(Mat,IS,const PetscScalar*);
5277b80b807SBarry Smith 
528be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatUseScaledForm(Mat,PetscTruth);
529be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScaleSystem(Mat,Vec,Vec);
530be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatUnScaleSystem(Mat,Vec,Vec);
5315ef9f2a5SBarry Smith 
532be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSize(Mat,PetscInt*,PetscInt*);
533be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalSize(Mat,PetscInt*,PetscInt*);
534be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRange(Mat,PetscInt*,PetscInt*);
5358ee2e534SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRanges(Mat,const PetscInt**);
536f4b8409dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRangeColumn(Mat,PetscInt*,PetscInt*);
53772ca8751SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRangesColumn(Mat,const PetscInt**);
5387b80b807SBarry Smith 
539be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
540be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroyMatrices(PetscInt,Mat *[]);
541be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrix(Mat,IS,IS,PetscInt,MatReuse,Mat *);
542abf0e977SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrixRaw(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscInt,MatReuse,Mat *);
543829201f2SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSeqNonzeroStructure(Mat,Mat *[]);
544829201f2SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroySeqNonzeroStructure(Mat *[]);
5455494a064SHong Zhang 
546be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
547be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPI(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*);
548be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPISymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*);
549be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPINumeric(Mat,Mat);
550be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroy_MPIAIJ_SeqsToMPI(Mat);
551be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalMat(Mat,MatReuse,Mat*);
552be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*);
553be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,PetscInt*,Mat*);
554dd6ea824SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBrowsOfAoCols(Mat,Mat,MatReuse,PetscInt**,MatScalar**,Mat*);
55543eb5e2fSMatthew Knepley #if defined (PETSC_USE_CTABLE)
556cd116e26SSatish Balay #include "petscctable.h"
55743eb5e2fSMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetCommunicationStructs(Mat, Vec *, PetscTable *, VecScatter *);
55843eb5e2fSMatthew Knepley #else
55943eb5e2fSMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetCommunicationStructs(Mat, Vec *, PetscInt *[], VecScatter *);
56043eb5e2fSMatthew Knepley #endif
5618c7482ecSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetGhosts(Mat, PetscInt *,const PetscInt *[]);
5628efafbd8SBarry Smith 
563be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt);
5647b80b807SBarry Smith 
565be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*);
566be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultSymbolic(Mat,Mat,PetscReal,Mat*);
567be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultNumeric(Mat,Mat,Mat);
56822440eb1SKris Buschelman 
569be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*);
570be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAPSymbolic(Mat,Mat,PetscReal,Mat*);
571be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAPNumeric(Mat,Mat,Mat);
57222440eb1SKris Buschelman 
573be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTranspose(Mat,Mat,MatReuse,PetscReal,Mat*);
574be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTransposeSymbolic(Mat,Mat,PetscReal,Mat*);
575be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTransposeNumeric(Mat,Mat,Mat);
576bc011b1eSHong Zhang 
577f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAXPY(Mat,PetscScalar,Mat,MatStructure);
57804aac2b0SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAYPX(Mat,PetscScalar,Mat,MatStructure);
579be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompress(Mat);
5801c741599SHong Zhang 
581f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScale(Mat,PetscScalar);
582f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShift(Mat,PetscScalar);
5837b80b807SBarry Smith 
584be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping);
585be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetLocalToGlobalMappingBlock(Mat,ISLocalToGlobalMapping);
586f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar);
587f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsLocalIS(Mat,IS,PetscScalar);
588be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
589be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
590052efed2SBarry Smith 
591be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStashSetInitialSize(Mat,PetscInt,PetscInt);
592be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
59390f02eecSBarry Smith 
594be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInterpolate(Mat,Vec,Vec);
5950c451bc4SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInterpolateAdd(Mat,Vec,Vec,Vec);
596ba966639SSatish Balay PetscPolymorphicSubroutine(MatInterpolateAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
597be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestrict(Mat,Vec,Vec);
598be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetVecs(Mat,Vec*,Vec*);
59969db28dcSHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRedundantMatrix(Mat,PetscInt,MPI_Comm,PetscInt,MatReuse,Mat*);
600bd481603SBarry Smith 
601bd481603SBarry Smith /*MC
6021620fd73SBarry Smith    MatSetValue - Set a single entry into a matrix.
6031620fd73SBarry Smith 
6041620fd73SBarry Smith    Not collective
6051620fd73SBarry Smith 
6061620fd73SBarry Smith    Input Parameters:
6071620fd73SBarry Smith +  m - the matrix
6081620fd73SBarry Smith .  row - the row location of the entry
6091620fd73SBarry Smith .  col - the column location of the entry
6101620fd73SBarry Smith .  value - the value to insert
6111620fd73SBarry Smith -  mode - either INSERT_VALUES or ADD_VALUES
6121620fd73SBarry Smith 
6131620fd73SBarry Smith    Notes:
6141620fd73SBarry Smith    For efficiency one should use MatSetValues() and set several or many
6151620fd73SBarry Smith    values simultaneously if possible.
6161620fd73SBarry Smith 
6171620fd73SBarry Smith    Level: beginner
6181620fd73SBarry Smith 
6191620fd73SBarry Smith .seealso: MatSetValues(), MatSetValueLocal()
6201620fd73SBarry Smith M*/
6211620fd73SBarry 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);}
6221620fd73SBarry Smith 
6231620fd73SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatGetValue(Mat v,PetscInt i,PetscInt j,PetscScalar *va) {return MatGetValues(v,1,&i,1,&j,va);}
6241620fd73SBarry Smith 
6251620fd73SBarry 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);}
6261620fd73SBarry Smith 
6271620fd73SBarry Smith /*MC
6280d2a7ff2SSatish Balay    MatPreallocateInitialize - Begins the block of code that will count the number of nonzeros per
629bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
630bd481603SBarry Smith 
631bd481603SBarry Smith    Synopsis:
632c1ac3661SBarry Smith    PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)
633bd481603SBarry Smith 
634bd481603SBarry Smith    Collective on MPI_Comm
635bd481603SBarry Smith 
636bd481603SBarry Smith    Input Parameters:
637bd481603SBarry Smith +  comm - the communicator that will share the eventually allocated matrix
638859fcb39SBarry Smith .  nrows - the number of LOCAL rows in the matrix
639859fcb39SBarry Smith -  ncols - the number of LOCAL columns in the matrix
640bd481603SBarry Smith 
641bd481603SBarry Smith    Output Parameters:
642bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
643bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
644bd481603SBarry Smith 
645bd481603SBarry Smith 
646bd481603SBarry Smith    Level: intermediate
647bd481603SBarry Smith 
648bd481603SBarry Smith    Notes:
649bd481603SBarry Smith    See the chapter in the users manual on performance for more details
650bd481603SBarry Smith 
6511d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
652bd481603SBarry Smith 
653bd481603SBarry Smith    Use MatPreallocateInitializeSymmetric() for symmetric matrices (MPISBAIJ matrices)
654bd481603SBarry Smith 
6551620fd73SBarry Smith    This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize().
6561620fd73SBarry Smith 
657bd481603SBarry Smith   Concepts: preallocation^Matrix
658bd481603SBarry Smith 
659bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
660bd481603SBarry Smith           MatPreallocateInitializeSymmetric(), MatPreallocateSymmetricSetLocal()
661bd481603SBarry Smith M*/
662c4f061fbSSatish Balay #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \
6637c922b88SBarry Smith { \
664a89bc211SBarry Smith   PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__start,__end; \
665a89bc211SBarry Smith   _4_ierr = PetscMalloc2(__nrows,PetscInt,&dnz,__nrows,PetscInt,&onz);CHKERRQ(_4_ierr); \
666a89bc211SBarry Smith   _4_ierr = PetscMemzero(dnz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
667a89bc211SBarry Smith   _4_ierr = PetscMemzero(onz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
668c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ctmp;\
669a89bc211SBarry Smith   _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows;
6707c922b88SBarry Smith 
671bd481603SBarry Smith /*MC
6720d2a7ff2SSatish Balay    MatPreallocateSymmetricInitialize - Begins the block of code that will count the number of nonzeros per
673bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
674bd481603SBarry Smith 
675bd481603SBarry Smith    Synopsis:
676c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)
677bd481603SBarry Smith 
678bd481603SBarry Smith    Collective on MPI_Comm
679bd481603SBarry Smith 
680bd481603SBarry Smith    Input Parameters:
681bd481603SBarry Smith +  comm - the communicator that will share the eventually allocated matrix
682859fcb39SBarry Smith .  nrows - the number of LOCAL rows in the matrix
683859fcb39SBarry Smith -  ncols - the number of LOCAL columns in the matrix
684bd481603SBarry Smith 
685bd481603SBarry Smith    Output Parameters:
686bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
687bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
688bd481603SBarry Smith 
689bd481603SBarry Smith 
690bd481603SBarry Smith    Level: intermediate
691bd481603SBarry Smith 
692bd481603SBarry Smith    Notes:
693bd481603SBarry Smith    See the chapter in the users manual on performance for more details
694bd481603SBarry Smith 
6951d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
696bd481603SBarry Smith 
6971620fd73SBarry Smith    This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize().
6981620fd73SBarry Smith 
699bd481603SBarry Smith   Concepts: preallocation^Matrix
700bd481603SBarry Smith 
701bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
702bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal()
703bd481603SBarry Smith M*/
704222b16d4SBarry Smith #define MatPreallocateSymmetricInitialize(comm,nrows,ncols,dnz,onz) 0; \
705222b16d4SBarry Smith { \
706a89bc211SBarry Smith   PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__end; \
707a89bc211SBarry Smith   _4_ierr = PetscMalloc2(__nrows,PetscInt,&dnz,__nrows,PetscInt,&onz);CHKERRQ(_4_ierr); \
708a89bc211SBarry Smith   _4_ierr = PetscMemzero(dnz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
709a89bc211SBarry Smith   _4_ierr = PetscMemzero(onz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
710c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr);\
711a89bc211SBarry Smith   _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows;
712222b16d4SBarry Smith 
713bd481603SBarry Smith /*MC
714bd481603SBarry Smith    MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
715bd481603SBarry Smith        inserted using a local number of the rows and columns
716bd481603SBarry Smith 
717bd481603SBarry Smith    Synopsis:
718c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
719bd481603SBarry Smith 
720bd481603SBarry Smith    Not Collective
721bd481603SBarry Smith 
722bd481603SBarry Smith    Input Parameters:
723bd481603SBarry Smith +  map - the mapping between local numbering and global numbering
724bd481603SBarry Smith .  nrows - the number of rows indicated
7251d73ed98SBarry Smith .  rows - the indices of the rows
726bd481603SBarry Smith .  ncols - the number of columns in the matrix
727bd481603SBarry Smith .  cols - the columns indicated
728bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
729bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
730bd481603SBarry Smith 
731bd481603SBarry Smith 
732bd481603SBarry Smith    Level: intermediate
733bd481603SBarry Smith 
734bd481603SBarry Smith    Notes:
735bd481603SBarry Smith    See the chapter in the users manual on performance for more details
736bd481603SBarry Smith 
7371d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
738bd481603SBarry Smith 
739bd481603SBarry Smith   Concepts: preallocation^Matrix
740bd481603SBarry Smith 
741bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
742bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal()
743bd481603SBarry Smith M*/
744c4f061fbSSatish Balay #define MatPreallocateSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\
745c4f061fbSSatish Balay {\
746c1ac3661SBarry Smith   PetscInt __l;\
747ef66eb69SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\
748ef66eb69SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\
749c4f061fbSSatish Balay   for (__l=0;__l<nrows;__l++) {\
750ef66eb69SBarry Smith     _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
751c4f061fbSSatish Balay   }\
752c4f061fbSSatish Balay }
753c4f061fbSSatish Balay 
754bd481603SBarry Smith /*MC
755bd481603SBarry Smith    MatPreallocateSymmetricSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
756bd481603SBarry Smith        inserted using a local number of the rows and columns
757bd481603SBarry Smith 
758bd481603SBarry Smith    Synopsis:
759c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
760bd481603SBarry Smith 
761bd481603SBarry Smith    Not Collective
762bd481603SBarry Smith 
763bd481603SBarry Smith    Input Parameters:
764bd481603SBarry Smith +  map - the mapping between local numbering and global numbering
765bd481603SBarry Smith .  nrows - the number of rows indicated
7661d73ed98SBarry Smith .  rows - the indices of the rows
767bd481603SBarry Smith .  ncols - the number of columns in the matrix
768bd481603SBarry Smith .  cols - the columns indicated
769bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
770bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
771bd481603SBarry Smith 
772bd481603SBarry Smith 
773bd481603SBarry Smith    Level: intermediate
774bd481603SBarry Smith 
775bd481603SBarry Smith    Notes:
776bd481603SBarry Smith    See the chapter in the users manual on performance for more details
777bd481603SBarry Smith 
778bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
779bd481603SBarry Smith 
780bd481603SBarry Smith   Concepts: preallocation^Matrix
781bd481603SBarry Smith 
782bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
783bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
784bd481603SBarry Smith M*/
785d3d32019SBarry Smith #define MatPreallocateSymmetricSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\
786d3d32019SBarry Smith {\
787c1ac3661SBarry Smith   PetscInt __l;\
788d3d32019SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\
789d3d32019SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\
790d3d32019SBarry Smith   for (__l=0;__l<nrows;__l++) {\
791d3d32019SBarry Smith     _4_ierr = MatPreallocateSymmetricSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
792d3d32019SBarry Smith   }\
793d3d32019SBarry Smith }
794d3d32019SBarry Smith 
795bd481603SBarry Smith /*MC
796bd481603SBarry Smith    MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
797bd481603SBarry Smith        inserted using a local number of the rows and columns
798bd481603SBarry Smith 
799bd481603SBarry Smith    Synopsis:
800c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
801bd481603SBarry Smith 
802bd481603SBarry Smith    Not Collective
803bd481603SBarry Smith 
804bd481603SBarry Smith    Input Parameters:
80564075487SBarry Smith +  row - the row
806bd481603SBarry Smith .  ncols - the number of columns in the matrix
807a50f8bf6SHong Zhang -  cols - the columns indicated
808a50f8bf6SHong Zhang 
809a50f8bf6SHong Zhang    Output Parameters:
810a50f8bf6SHong Zhang +  dnz - the array that will be passed to the matrix preallocation routines
811bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
812bd481603SBarry Smith 
813bd481603SBarry Smith 
814bd481603SBarry Smith    Level: intermediate
815bd481603SBarry Smith 
816bd481603SBarry Smith    Notes:
817bd481603SBarry Smith    See the chapter in the users manual on performance for more details
818bd481603SBarry Smith 
819bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
820bd481603SBarry Smith 
8211620fd73SBarry Smith    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
8221620fd73SBarry Smith 
823bd481603SBarry Smith   Concepts: preallocation^Matrix
824bd481603SBarry Smith 
825bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
826bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
827bd481603SBarry Smith M*/
828c4f061fbSSatish Balay #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\
829c1ac3661SBarry Smith { PetscInt __i; \
8308ee2e534SBarry Smith   if (row < __rstart) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Trying to set preallocation for row %D less than first local row %D",row,__rstart);\
831a89bc211SBarry 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);\
8327c922b88SBarry Smith   for (__i=0; __i<nc; __i++) {\
83364075487SBarry Smith     if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \
8347cc688d7SBarry Smith     else dnz[row - __rstart]++;\
8357c922b88SBarry Smith   }\
8367c922b88SBarry Smith }
8377c922b88SBarry Smith 
838bd481603SBarry Smith /*MC
839bd481603SBarry Smith    MatPreallocateSymmetricSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
840bd481603SBarry Smith        inserted using a local number of the rows and columns
841bd481603SBarry Smith 
842bd481603SBarry Smith    Synopsis:
843c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
844bd481603SBarry Smith 
845bd481603SBarry Smith    Not Collective
846bd481603SBarry Smith 
847bd481603SBarry Smith    Input Parameters:
848bd481603SBarry Smith +  nrows - the number of rows indicated
8491d73ed98SBarry Smith .  rows - the indices of the rows
850bd481603SBarry Smith .  ncols - the number of columns in the matrix
851bd481603SBarry Smith .  cols - the columns indicated
852bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
853bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
854bd481603SBarry Smith 
855bd481603SBarry Smith 
856bd481603SBarry Smith    Level: intermediate
857bd481603SBarry Smith 
858bd481603SBarry Smith    Notes:
859bd481603SBarry Smith    See the chapter in the users manual on performance for more details
860bd481603SBarry Smith 
861bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
862bd481603SBarry Smith 
8631620fd73SBarry Smith    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
8641620fd73SBarry Smith 
865bd481603SBarry Smith   Concepts: preallocation^Matrix
866bd481603SBarry Smith 
867bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
868bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
869bd481603SBarry Smith M*/
870d3d32019SBarry Smith #define MatPreallocateSymmetricSet(row,nc,cols,dnz,onz) 0;\
871c1ac3661SBarry Smith { PetscInt __i; \
872d3d32019SBarry Smith   for (__i=0; __i<nc; __i++) {\
873d3d32019SBarry Smith     if (cols[__i] >= __end) onz[row - __rstart]++; \
874d3d32019SBarry Smith     else if (cols[__i] >= row) dnz[row - __rstart]++;\
875d3d32019SBarry Smith   }\
876d3d32019SBarry Smith }
877d3d32019SBarry Smith 
878bd481603SBarry Smith /*MC
87916371a99SBarry Smith    MatPreallocateLocation -  An alternative to MatPreallocationSet() that puts the nonzero locations into the matrix if it exists
88016371a99SBarry Smith 
88116371a99SBarry Smith    Synopsis:
88216371a99SBarry Smith    PetscErrorCode MatPreallocateLocations(Mat A,PetscInt row,PetscInt ncols,PetscInt *cols,PetscInt *dnz,PetscInt *onz)
88316371a99SBarry Smith 
88416371a99SBarry Smith    Not Collective
88516371a99SBarry Smith 
88616371a99SBarry Smith    Input Parameters:
88716371a99SBarry Smith .  A - matrix
88816371a99SBarry Smith .  row - row where values exist (must be local to this process)
88916371a99SBarry Smith .  ncols - number of columns
89016371a99SBarry Smith .  cols - columns with nonzeros
89116371a99SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
89216371a99SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
89316371a99SBarry Smith 
89416371a99SBarry Smith 
89516371a99SBarry Smith    Level: intermediate
89616371a99SBarry Smith 
89716371a99SBarry Smith    Notes:
89816371a99SBarry Smith    See the chapter in the users manual on performance for more details
89916371a99SBarry Smith 
90016371a99SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
90116371a99SBarry Smith 
90216371a99SBarry Smith    This is a MACRO not a function because it uses a bunch of variables private to the MatPreallocation.... routines.
90316371a99SBarry Smith 
90416371a99SBarry Smith   Concepts: preallocation^Matrix
90516371a99SBarry Smith 
90616371a99SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
90716371a99SBarry Smith           MatPreallocateSymmetricInitialize(), MatPreallocateSymmetricSetLocal()
90816371a99SBarry Smith M*/
90916371a99SBarry 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);}
91016371a99SBarry Smith 
91116371a99SBarry Smith 
91216371a99SBarry Smith /*MC
9130d2a7ff2SSatish Balay    MatPreallocateFinalize - Ends the block of code that will count the number of nonzeros per
914bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
915bd481603SBarry Smith 
916bd481603SBarry Smith    Synopsis:
917c1ac3661SBarry Smith    PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz)
918bd481603SBarry Smith 
919bd481603SBarry Smith    Collective on MPI_Comm
920bd481603SBarry Smith 
921bd481603SBarry Smith    Input Parameters:
92216371a99SBarry Smith +  dnz - the array that was be passed to the matrix preallocation routines
923bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
924bd481603SBarry Smith 
925bd481603SBarry Smith 
926bd481603SBarry Smith    Level: intermediate
927bd481603SBarry Smith 
928bd481603SBarry Smith    Notes:
929bd481603SBarry Smith    See the chapter in the users manual on performance for more details
930bd481603SBarry Smith 
931bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
932bd481603SBarry Smith 
9331620fd73SBarry Smith    This is a MACRO not a function because it closes the { started in MatPreallocateInitialize().
9341620fd73SBarry Smith 
935bd481603SBarry Smith   Concepts: preallocation^Matrix
936bd481603SBarry Smith 
937bd481603SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
938bd481603SBarry Smith           MatPreallocateSymmetricInitialize(), MatPreallocateSymmetricSetLocal()
939bd481603SBarry Smith M*/
940a89bc211SBarry Smith #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree2(dnz,onz);CHKERRQ(_4_ierr);}
9417c922b88SBarry Smith 
942bd481603SBarry Smith 
943bd481603SBarry Smith 
9447b80b807SBarry Smith /* Routines unique to particular data structures */
945be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellGetContext(Mat,void **);
946ba966639SSatish Balay PetscPolymorphicFunction(MatShellGetContext,(Mat A),(A,&t),void*,t)
947698d4c6aSKris Buschelman 
948be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInodeAdjustForInodes(Mat,IS*,IS*);
949be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *);
950698d4c6aSKris Buschelman 
951be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetColumnIndices(Mat,PetscInt[]);
952be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetColumnIndices(Mat,PetscInt[]);
953be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
954c75a6043SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
955c75a6043SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
9567b80b807SBarry Smith 
957a96a251dSBarry Smith #define MAT_SKIP_ALLOCATION -4
958a96a251dSBarry Smith 
959be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
960ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz))
961be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
962ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqSBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz))
963be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]);
964ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqAIJSetPreallocation,(Mat A,const PetscInt nnz[]),(A,0,nnz))
965273d9f13SBarry Smith 
966be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
967ba966639SSatish Balay PetscPolymorphicSubroutine(MatMPIBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[],const PetscInt onz[]),(A,bs,0,nnz,0,onz))
968be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
969be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
97053803f26SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []);
971725b52f3SLisandro Dalcin EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
972be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]);
973be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
974be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]);
975be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIDenseSetPreallocation(Mat,PetscScalar[]);
976284134d9SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation(Mat,PetscScalar[]);
977be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIRowbsSetPreallocation(Mat,PetscInt,const PetscInt[]);
978be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,PetscInt*[]);
979be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,PetscInt*[]);
980be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAdicSetLocalFunction(Mat,void (*)(void));
981273d9f13SBarry Smith 
982be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetLDA(Mat,PetscInt);
983ab92ecdeSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDenseGetLocalMatrix(Mat,Mat*);
9841b807ce4Svictorle 
985be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStoreValues(Mat);
986be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRetrieveValues(Mat);
9872e8a6d31SBarry Smith 
988be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDAADSetCtx(Mat,void*);
9893a7fca6bSBarry Smith 
9907b80b807SBarry Smith /*
9917b80b807SBarry Smith   These routines are not usually accessed directly, rather solving is
99294b7f48cSBarry Smith   done through the KSP and PC interfaces.
9937b80b807SBarry Smith */
9947b80b807SBarry Smith 
995d9274352SBarry Smith /*E
996d9274352SBarry Smith     MatOrderingType - String with the name of a PETSc matrix ordering or the creation function
997d9274352SBarry Smith        with an optional dynamic library name, for example
998d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:orderingcreate()
999d9274352SBarry Smith 
1000d9274352SBarry Smith    Level: beginner
1001d9274352SBarry Smith 
1002e5a9bf91SBarry Smith    Cannot use const because the PC objects manipulate the string
1003e5a9bf91SBarry Smith 
1004d9274352SBarry Smith .seealso: MatGetOrdering()
1005d9274352SBarry Smith E*/
10063eaad2caSSatish Balay #define MatOrderingType char*
1007b12f92e5SBarry Smith #define MATORDERING_NATURAL   "natural"
1008b12f92e5SBarry Smith #define MATORDERING_ND        "nd"
1009b12f92e5SBarry Smith #define MATORDERING_1WD       "1wd"
1010b12f92e5SBarry Smith #define MATORDERING_RCM       "rcm"
1011b12f92e5SBarry Smith #define MATORDERING_QMD       "qmd"
1012b12f92e5SBarry Smith #define MATORDERING_ROWLENGTH "rowlength"
101362152c8bSBarry Smith #define MATORDERING_DSC_ND    "dsc_nd"
101462152c8bSBarry Smith #define MATORDERING_DSC_MMD   "dsc_mmd"
101562152c8bSBarry Smith #define MATORDERING_DSC_MDF   "dsc_mdf"
1016c06d978dSMatthew Knepley #define MATORDERING_CONSTRAINED "constrained"
1017c06d978dSMatthew Knepley #define MATORDERING_IDENTITY  "identity"
1018c06d978dSMatthew Knepley #define MATORDERING_REVERSE   "reverse"
1019b12f92e5SBarry Smith 
1020f2416520SSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOrdering(Mat,const MatOrderingType,IS*,IS*);
1021be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOrderingList(PetscFList *list);
1022f2416520SSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,const MatOrderingType,IS*,IS*));
102330de9b25SBarry Smith 
102430de9b25SBarry Smith /*MC
102530de9b25SBarry Smith    MatOrderingRegisterDynamic - Adds a new sparse matrix ordering to the
102630de9b25SBarry Smith                                matrix package.
102730de9b25SBarry Smith 
102830de9b25SBarry Smith    Synopsis:
1029c1ac3661SBarry Smith    PetscErrorCode MatOrderingRegisterDynamic(char *name_ordering,char *path,char *name_create,PetscErrorCode (*routine_create)(MatOrdering))
103030de9b25SBarry Smith 
103130de9b25SBarry Smith    Not Collective
103230de9b25SBarry Smith 
103330de9b25SBarry Smith    Input Parameters:
103430de9b25SBarry Smith +  sname - name of ordering (for example MATORDERING_ND)
103530de9b25SBarry Smith .  path - location of library where creation routine is
103630de9b25SBarry Smith .  name - name of function that creates the ordering type,a string
103730de9b25SBarry Smith -  function - function pointer that creates the ordering
103830de9b25SBarry Smith 
103930de9b25SBarry Smith    Level: developer
104030de9b25SBarry Smith 
104130de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
104230de9b25SBarry Smith    is ignored.
104330de9b25SBarry Smith 
104430de9b25SBarry Smith    Sample usage:
104530de9b25SBarry Smith .vb
104630de9b25SBarry Smith    MatOrderingRegisterDynamic("my_order",/home/username/my_lib/lib/libO/solaris/mylib.a,
104730de9b25SBarry Smith                "MyOrder",MyOrder);
104830de9b25SBarry Smith .ve
104930de9b25SBarry Smith 
105030de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
105130de9b25SBarry Smith $     MatOrderingSetType(part,"my_order)
105230de9b25SBarry Smith    or at runtime via the option
10532401956bSBarry Smith $     -pc_factor_mat_ordering_type my_order
105430de9b25SBarry Smith 
1055ab901514SBarry Smith    ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values.
105630de9b25SBarry Smith 
105730de9b25SBarry Smith .keywords: matrix, ordering, register
105830de9b25SBarry Smith 
105930de9b25SBarry Smith .seealso: MatOrderingRegisterDestroy(), MatOrderingRegisterAll()
106030de9b25SBarry Smith M*/
1061aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1062f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,0)
1063b12f92e5SBarry Smith #else
1064f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,d)
1065b12f92e5SBarry Smith #endif
106630de9b25SBarry Smith 
1067be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegisterDestroy(void);
1068be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegisterAll(const char[]);
10692bad1931SBarry Smith extern PetscTruth MatOrderingRegisterAllCalled;
1070b0a32e0cSBarry Smith extern PetscFList MatOrderingList;
1071d4fbbf0eSBarry Smith 
1072be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS);
1073a2ce50c7SBarry Smith 
1074d91e6319SBarry Smith /*S
10752401956bSBarry Smith    MatFactorInfo - Data passed into the matrix factorization routines
1076b00f7748SHong Zhang 
107761cad860SBarry Smith    In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE, that is use
107861cad860SBarry Smith $     MatFactorInfo  info(MAT_FACTORINFO_SIZE)
1079b00f7748SHong Zhang 
108015e8a5b3SHong Zhang    Notes: These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC.
1081b00f7748SHong Zhang 
108261cad860SBarry Smith       You can use MatFactorInfoInitialize() to set default values.
108361cad860SBarry Smith 
1084b00f7748SHong Zhang    Level: developer
1085b00f7748SHong Zhang 
1086d7d82daaSBarry Smith .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(),
1087d7d82daaSBarry Smith           MatFactorInfoInitialize()
1088b00f7748SHong Zhang 
1089b00f7748SHong Zhang S*/
1090b00f7748SHong Zhang typedef struct {
10910a29876aSHong Zhang   PetscReal     shiftnz;        /* scaling of identity added to matrix to prevent zero pivots */
1092fbf22428SSatish Balay   PetscReal     shiftpd;        /* if true, shift until positive pivots */
109315e8a5b3SHong Zhang   PetscReal     diagonal_fill;  /* force diagonal to fill in if initially not filled */
109485317021SBarry Smith   PetscReal     usedt;
109515e8a5b3SHong Zhang   PetscReal     dt;             /* drop tolerance */
1096b00f7748SHong Zhang   PetscReal     dtcol;          /* tolerance for pivoting */
109715e8a5b3SHong Zhang   PetscReal     dtcount;        /* maximum nonzeros to be allowed per row */
109867e28bfeSBarry Smith   PetscReal     fill;           /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */
1099348344bbSBarry Smith   PetscReal     levels;         /* ICC/ILU(levels) */
1100bcd9e38bSBarry Smith   PetscReal     pivotinblocks;  /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0
1101bcd9e38bSBarry Smith                                    factorization may be faster if do not pivot */
110262bba022SBarry Smith   PetscReal     shiftinblocks;  /* if block in block factorization has zero pivot then shift diagonal until non-singular */
110315e8a5b3SHong Zhang   PetscReal     zeropivot;      /* pivot is called zero if less than this */
110415e8a5b3SHong Zhang } MatFactorInfo;
1105ffa6d0a5SLois Curfman McInnes 
1106be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFactorInfoInitialize(MatFactorInfo*);
11070481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactor(Mat,IS,const MatFactorInfo*);
11080481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
11090481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactorNumeric(Mat,Mat,const MatFactorInfo*);
11100481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactor(Mat,IS,IS,const MatFactorInfo*);
11110481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUFactor(Mat,IS,IS,const MatFactorInfo*);
11120481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
11130481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
11140481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatICCFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
11150481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatICCFactor(Mat,IS,const MatFactorInfo*);
11160481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactorNumeric(Mat,Mat,const MatFactorInfo*);
11170481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUDTFactor(Mat,IS,IS,const MatFactorInfo*,Mat *);
1118*599ef60dSHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUDTFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo *);
1119*599ef60dSHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUDTFactorNumeric(Mat,Mat,const MatFactorInfo*);
1120be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*);
1121be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolve(Mat,Vec,Vec);
1122be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatForwardSolve(Mat,Vec,Vec);
1123be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatBackwardSolve(Mat,Vec,Vec);
1124be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveAdd(Mat,Vec,Vec,Vec);
1125be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveTranspose(Mat,Vec,Vec);
1126be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveTransposeAdd(Mat,Vec,Vec,Vec);
1127be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolves(Mat,Vecs,Vecs);
11288ed539a5SBarry Smith 
1129be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUnfactored(Mat);
1130bb5a7306SBarry Smith 
1131d91e6319SBarry Smith /*E
1132d91e6319SBarry Smith     MatSORType - What type of (S)SOR to perform
1133bb1eb677SSatish Balay 
1134d91e6319SBarry Smith     Level: beginner
1135d91e6319SBarry Smith 
1136d9274352SBarry Smith    May be bitwise ORd together
1137d9274352SBarry Smith 
1138d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
1139d91e6319SBarry Smith 
11404e7234bfSBarry Smith    MatSORType may be bitwise ORd together, so do not change the numbers
11414e7234bfSBarry Smith 
1142a1d92eedSBarry Smith .seealso: MatRelax(), MatPBRelax()
1143d91e6319SBarry Smith E*/
1144ee50ffe9SBarry Smith typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3,
1145ee50ffe9SBarry Smith               SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8,
1146ee50ffe9SBarry Smith               SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16,
114784cb2905SBarry Smith               SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType;
1148be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRelax(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
1149be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPBRelax(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
11508ed539a5SBarry Smith 
1151d4fbbf0eSBarry Smith /*
1152639f9d9dSBarry Smith     These routines are for efficiently computing Jacobians via finite differences.
1153639f9d9dSBarry Smith */
1154b12f92e5SBarry Smith 
1155d9274352SBarry Smith /*E
1156d9274352SBarry Smith     MatColoringType - String with the name of a PETSc matrix coloring or the creation function
1157d9274352SBarry Smith        with an optional dynamic library name, for example
1158d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:coloringcreate()
1159d9274352SBarry Smith 
1160d9274352SBarry Smith    Level: beginner
1161d9274352SBarry Smith 
1162d9274352SBarry Smith .seealso: MatGetColoring()
1163d9274352SBarry Smith E*/
1164a313700dSBarry Smith #define MatColoringType char*
1165b12f92e5SBarry Smith #define MATCOLORING_NATURAL "natural"
1166b12f92e5SBarry Smith #define MATCOLORING_SL      "sl"
1167b12f92e5SBarry Smith #define MATCOLORING_LF      "lf"
1168b12f92e5SBarry Smith #define MATCOLORING_ID      "id"
1169b12f92e5SBarry Smith 
1170ab970d9bSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColoring(Mat,const MatColoringType,ISColoring*);
11712e90c967SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,MatColoringType,ISColoring *));
117230de9b25SBarry Smith 
117330de9b25SBarry Smith /*MC
117430de9b25SBarry Smith    MatColoringRegisterDynamic - Adds a new sparse matrix coloring to the
117530de9b25SBarry Smith                                matrix package.
117630de9b25SBarry Smith 
117730de9b25SBarry Smith    Synopsis:
1178c1ac3661SBarry Smith    PetscErrorCode MatColoringRegisterDynamic(char *name_coloring,char *path,char *name_create,PetscErrorCode (*routine_create)(MatColoring))
117930de9b25SBarry Smith 
118030de9b25SBarry Smith    Not Collective
118130de9b25SBarry Smith 
118230de9b25SBarry Smith    Input Parameters:
118330de9b25SBarry Smith +  sname - name of Coloring (for example MATCOLORING_SL)
118430de9b25SBarry Smith .  path - location of library where creation routine is
118530de9b25SBarry Smith .  name - name of function that creates the Coloring type, a string
118630de9b25SBarry Smith -  function - function pointer that creates the coloring
118730de9b25SBarry Smith 
118830de9b25SBarry Smith    Level: developer
118930de9b25SBarry Smith 
119030de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
119130de9b25SBarry Smith    is ignored.
119230de9b25SBarry Smith 
119330de9b25SBarry Smith    Sample usage:
119430de9b25SBarry Smith .vb
119530de9b25SBarry Smith    MatColoringRegisterDynamic("my_color",/home/username/my_lib/lib/libO/solaris/mylib.a,
119630de9b25SBarry Smith                "MyColor",MyColor);
119730de9b25SBarry Smith .ve
119830de9b25SBarry Smith 
119930de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
120030de9b25SBarry Smith $     MatColoringSetType(part,"my_color")
120130de9b25SBarry Smith    or at runtime via the option
120230de9b25SBarry Smith $     -mat_coloring_type my_color
120330de9b25SBarry Smith 
1204ab901514SBarry Smith    $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
120530de9b25SBarry Smith 
120630de9b25SBarry Smith .keywords: matrix, Coloring, register
120730de9b25SBarry Smith 
120830de9b25SBarry Smith .seealso: MatColoringRegisterDestroy(), MatColoringRegisterAll()
120930de9b25SBarry Smith M*/
1210aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1211f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,0)
1212b12f92e5SBarry Smith #else
1213f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,d)
1214b12f92e5SBarry Smith #endif
121530de9b25SBarry Smith 
12162bad1931SBarry Smith extern PetscTruth MatColoringRegisterAllCalled;
1217f1144a30SSatish Balay 
1218be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegisterAll(const char[]);
1219be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegisterDestroy(void);
1220be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*);
1221639f9d9dSBarry Smith 
1222d9274352SBarry Smith /*S
1223d9274352SBarry Smith      MatFDColoring - Object for computing a sparse Jacobian via finite differences
1224d9274352SBarry Smith         and coloring
1225639f9d9dSBarry Smith 
1226d9274352SBarry Smith    Level: beginner
1227d9274352SBarry Smith 
1228d9274352SBarry Smith   Concepts: coloring, sparse Jacobian, finite differences
1229d9274352SBarry Smith 
1230d9274352SBarry Smith .seealso:  MatFDColoringCreate()
1231d9274352SBarry Smith S*/
1232e2a1c21fSSatish Balay typedef struct _p_MatFDColoring* MatFDColoring;
1233639f9d9dSBarry Smith 
1234be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringCreate(Mat,ISColoring,MatFDColoring *);
1235be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringDestroy(MatFDColoring);
1236be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringView(MatFDColoring,PetscViewer);
1237be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*);
12384a9d489dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**);
1239be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal);
1240be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFromOptions(MatFDColoring);
1241be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApply(Mat,MatFDColoring,Vec,MatStructure*,void *);
1242be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApplyTS(Mat,MatFDColoring,PetscReal,Vec,MatStructure*,void *);
1243be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetF(MatFDColoring,Vec);
1244be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,PetscInt*[]);
1245639f9d9dSBarry Smith /*
12460752156aSBarry Smith     These routines are for partitioning matrices: currently used only
12473eda8832SBarry Smith   for adjacency matrix, MatCreateMPIAdj().
12480752156aSBarry Smith */
1249ca161407SBarry Smith 
1250d9274352SBarry Smith /*S
1251d9274352SBarry Smith      MatPartitioning - Object for managing the partitioning of a matrix or graph
1252d9274352SBarry Smith 
1253d9274352SBarry Smith    Level: beginner
1254d9274352SBarry Smith 
1255d9274352SBarry Smith   Concepts: partitioning
1256d9274352SBarry Smith 
1257743a1026Svictorle .seealso:  MatPartitioningCreate(), MatPartitioningType
1258d9274352SBarry Smith S*/
125991e9ee9fSBarry Smith typedef struct _p_MatPartitioning* MatPartitioning;
1260d9274352SBarry Smith 
1261d9274352SBarry Smith /*E
12625bc6e7ccSvictorle     MatPartitioningType - String with the name of a PETSc matrix partitioning or the creation function
1263d9274352SBarry Smith        with an optional dynamic library name, for example
1264d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:partitioningcreate()
1265d9274352SBarry Smith 
1266d9274352SBarry Smith    Level: beginner
1267d9274352SBarry Smith 
1268b547a13bSvictorle .seealso: MatPartitioningCreate(), MatPartitioning
1269d9274352SBarry Smith E*/
1270a313700dSBarry Smith #define MatPartitioningType char*
12718ba1e511SMatthew Knepley #define MAT_PARTITIONING_CURRENT  "current"
1272c2ec2a73SHong Zhang #define MAT_PARTITIONING_SQUARE   "square"
12738ba1e511SMatthew Knepley #define MAT_PARTITIONING_PARMETIS "parmetis"
127471306c60Sjroman #define MAT_PARTITIONING_CHACO    "chaco"
127571306c60Sjroman #define MAT_PARTITIONING_JOSTLE   "jostle"
127671306c60Sjroman #define MAT_PARTITIONING_PARTY    "party"
127771306c60Sjroman #define MAT_PARTITIONING_SCOTCH   "scotch"
127871306c60Sjroman 
1279ca161407SBarry Smith 
1280be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningCreate(MPI_Comm,MatPartitioning*);
1281a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetType(MatPartitioning,const MatPartitioningType);
1282be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetNParts(MatPartitioning,PetscInt);
1283be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetAdjacency(MatPartitioning,Mat);
1284be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]);
1285be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []);
1286be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningApply(MatPartitioning,IS*);
1287be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningDestroy(MatPartitioning);
12882aabb6bbSBarry Smith 
1289be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatPartitioning));
129030de9b25SBarry Smith 
129130de9b25SBarry Smith /*MC
129230de9b25SBarry Smith    MatPartitioningRegisterDynamic - Adds a new sparse matrix partitioning to the
129330de9b25SBarry Smith    matrix package.
129430de9b25SBarry Smith 
129530de9b25SBarry Smith    Synopsis:
1296c1ac3661SBarry Smith    PetscErrorCode MatPartitioningRegisterDynamic(char *name_partitioning,char *path,char *name_create,PetscErrorCode (*routine_create)(MatPartitioning))
129730de9b25SBarry Smith 
129830de9b25SBarry Smith    Not Collective
129930de9b25SBarry Smith 
130030de9b25SBarry Smith    Input Parameters:
130130de9b25SBarry Smith +  sname - name of partitioning (for example MAT_PARTITIONING_CURRENT) or parmetis
130230de9b25SBarry Smith .  path - location of library where creation routine is
130330de9b25SBarry Smith .  name - name of function that creates the partitioning type, a string
130430de9b25SBarry Smith -  function - function pointer that creates the partitioning type
130530de9b25SBarry Smith 
130630de9b25SBarry Smith    Level: developer
130730de9b25SBarry Smith 
130830de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
130930de9b25SBarry Smith    is ignored.
131030de9b25SBarry Smith 
131130de9b25SBarry Smith    Sample usage:
131230de9b25SBarry Smith .vb
131330de9b25SBarry Smith    MatPartitioningRegisterDynamic("my_part",/home/username/my_lib/lib/libO/solaris/mylib.a,
131430de9b25SBarry Smith                "MyPartCreate",MyPartCreate);
131530de9b25SBarry Smith .ve
131630de9b25SBarry Smith 
131730de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
131830de9b25SBarry Smith $     MatPartitioningSetType(part,"my_part")
131930de9b25SBarry Smith    or at runtime via the option
132030de9b25SBarry Smith $     -mat_partitioning_type my_part
132130de9b25SBarry Smith 
1322ab901514SBarry Smith    $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
132330de9b25SBarry Smith 
132430de9b25SBarry Smith .keywords: matrix, partitioning, register
132530de9b25SBarry Smith 
132630de9b25SBarry Smith .seealso: MatPartitioningRegisterDestroy(), MatPartitioningRegisterAll()
132730de9b25SBarry Smith M*/
1328aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1329f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,0)
13302aabb6bbSBarry Smith #else
1331f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,d)
13322aabb6bbSBarry Smith #endif
13332aabb6bbSBarry Smith 
13342bad1931SBarry Smith extern PetscTruth MatPartitioningRegisterAllCalled;
1335f1144a30SSatish Balay 
1336be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegisterAll(const char[]);
1337be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegisterDestroy(void);
13382bad1931SBarry Smith 
1339be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningView(MatPartitioning,PetscViewer);
1340be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetFromOptions(MatPartitioning);
1341a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningGetType(MatPartitioning,const MatPartitioningType*);
1342ca161407SBarry Smith 
1343be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningParmetisSetCoarseSequential(MatPartitioning);
13440e8af3caSMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningParmetisGetEdgeCut(MatPartitioning, PetscInt *);
13450752156aSBarry Smith 
1346be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningJostleSetCoarseLevel(MatPartitioning,PetscReal);
1347be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningJostleSetCoarseSequential(MatPartitioning);
134871306c60Sjroman 
1349954b3ec6SBarry Smith typedef enum {MP_CHACO_MULTILEVEL_KL,MP_CHACO_SPECTRAL,MP_CHACO_LINEAR,MP_CHACO_RANDOM, MP_CHACO_SCATTERED} MPChacoGlobalType;
1350be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetGlobal(MatPartitioning, MPChacoGlobalType);
135171306c60Sjroman typedef enum { MP_CHACO_KERNIGHAN_LIN, MP_CHACO_NONE } MPChacoLocalType;
1352be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetLocal(MatPartitioning, MPChacoLocalType);
1353be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal);
135471306c60Sjroman typedef enum { MP_CHACO_LANCZOS, MP_CHACO_RQI_SYMMLQ } MPChacoEigenType;
1355be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType);
1356be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenTol(MatPartitioning, PetscReal);
1357be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenNumber(MatPartitioning, PetscInt);
135871306c60Sjroman 
135971306c60Sjroman #define MP_PARTY_OPT "opt"
136071306c60Sjroman #define MP_PARTY_LIN "lin"
136171306c60Sjroman #define MP_PARTY_SCA "sca"
136271306c60Sjroman #define MP_PARTY_RAN "ran"
136371306c60Sjroman #define MP_PARTY_GBF "gbf"
136471306c60Sjroman #define MP_PARTY_GCF "gcf"
136571306c60Sjroman #define MP_PARTY_BUB "bub"
136671306c60Sjroman #define MP_PARTY_DEF "def"
1367be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetGlobal(MatPartitioning, const char*);
136871306c60Sjroman #define MP_PARTY_HELPFUL_SETS "hs"
136971306c60Sjroman #define MP_PARTY_KERNIGHAN_LIN "kl"
137071306c60Sjroman #define MP_PARTY_NONE "no"
1371be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetLocal(MatPartitioning, const char*);
1372be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal);
1373be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetBipart(MatPartitioning,PetscTruth);
1374be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscTruth);
137571306c60Sjroman 
137671306c60Sjroman typedef enum { MP_SCOTCH_GREEDY, MP_SCOTCH_GPS, MP_SCOTCH_GR_GPS } MPScotchGlobalType;
1377be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetArch(MatPartitioning,const char*);
1378be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetMultilevel(MatPartitioning);
1379be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetGlobal(MatPartitioning,MPScotchGlobalType);
1380be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetCoarseLevel(MatPartitioning,PetscReal);
1381be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetHostList(MatPartitioning,const char*);
138271306c60Sjroman typedef enum { MP_SCOTCH_KERNIGHAN_LIN, MP_SCOTCH_NONE } MPScotchLocalType;
1383be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetLocal(MatPartitioning,MPScotchLocalType);
1384be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetMapping(MatPartitioning);
1385be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetStrategy(MatPartitioning,char*);
138671306c60Sjroman 
13870752156aSBarry Smith /*
13880a835dfdSSatish Balay     If you add entries here you must also add them to finclude/petscmat.h
1389d4fbbf0eSBarry Smith */
13901c1c02c0SLois Curfman McInnes typedef enum { MATOP_SET_VALUES=0,
13911c1c02c0SLois Curfman McInnes                MATOP_GET_ROW=1,
13921c1c02c0SLois Curfman McInnes                MATOP_RESTORE_ROW=2,
13931c1c02c0SLois Curfman McInnes                MATOP_MULT=3,
13941c1c02c0SLois Curfman McInnes                MATOP_MULT_ADD=4,
13957c922b88SBarry Smith                MATOP_MULT_TRANSPOSE=5,
13967c922b88SBarry Smith                MATOP_MULT_TRANSPOSE_ADD=6,
13971c1c02c0SLois Curfman McInnes                MATOP_SOLVE=7,
13981c1c02c0SLois Curfman McInnes                MATOP_SOLVE_ADD=8,
13997c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE=9,
14007c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE_ADD=10,
14011c1c02c0SLois Curfman McInnes                MATOP_LUFACTOR=11,
14021c1c02c0SLois Curfman McInnes                MATOP_CHOLESKYFACTOR=12,
14031c1c02c0SLois Curfman McInnes                MATOP_RELAX=13,
14041c1c02c0SLois Curfman McInnes                MATOP_TRANSPOSE=14,
14051c1c02c0SLois Curfman McInnes                MATOP_GETINFO=15,
14061c1c02c0SLois Curfman McInnes                MATOP_EQUAL=16,
14071c1c02c0SLois Curfman McInnes                MATOP_GET_DIAGONAL=17,
14081c1c02c0SLois Curfman McInnes                MATOP_DIAGONAL_SCALE=18,
14091c1c02c0SLois Curfman McInnes                MATOP_NORM=19,
14101c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_BEGIN=20,
14111c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_END=21,
14121c1c02c0SLois Curfman McInnes                MATOP_COMPRESS=22,
14131c1c02c0SLois Curfman McInnes                MATOP_SET_OPTION=23,
14141c1c02c0SLois Curfman McInnes                MATOP_ZERO_ENTRIES=24,
14151c1c02c0SLois Curfman McInnes                MATOP_ZERO_ROWS=25,
14161c1c02c0SLois Curfman McInnes                MATOP_LUFACTOR_SYMBOLIC=26,
14171c1c02c0SLois Curfman McInnes                MATOP_LUFACTOR_NUMERIC=27,
14181c1c02c0SLois Curfman McInnes                MATOP_CHOLESKY_FACTOR_SYMBOLIC=28,
14191c1c02c0SLois Curfman McInnes                MATOP_CHOLESKY_FACTOR_NUMERIC=29,
1420d643ce63SMatthew Knepley                MATOP_SETUP_PREALLOCATION=30,
1421d643ce63SMatthew Knepley                MATOP_ILUFACTOR_SYMBOLIC=31,
1422d643ce63SMatthew Knepley                MATOP_ICCFACTOR_SYMBOLIC=32,
1423d643ce63SMatthew Knepley                MATOP_GET_ARRAY=33,
1424d643ce63SMatthew Knepley                MATOP_RESTORE_ARRAY=34,
1425e26f1c3cSBarry Smith                MATOP_DUPLICATE=35,
1426d643ce63SMatthew Knepley                MATOP_FORWARD_SOLVE=36,
1427d643ce63SMatthew Knepley                MATOP_BACKWARD_SOLVE=37,
1428d643ce63SMatthew Knepley                MATOP_ILUFACTOR=38,
1429d643ce63SMatthew Knepley                MATOP_ICCFACTOR=39,
1430d643ce63SMatthew Knepley                MATOP_AXPY=40,
1431d643ce63SMatthew Knepley                MATOP_GET_SUBMATRICES=41,
1432d643ce63SMatthew Knepley                MATOP_INCREASE_OVERLAP=42,
1433d643ce63SMatthew Knepley                MATOP_GET_VALUES=43,
1434d643ce63SMatthew Knepley                MATOP_COPY=44,
1435ff1cced0SMatthew Knepley                MATOP_GET_ROW_MAX=45,
1436d643ce63SMatthew Knepley                MATOP_SCALE=46,
1437d643ce63SMatthew Knepley                MATOP_SHIFT=47,
1438d643ce63SMatthew Knepley                MATOP_DIAGONAL_SHIFT=48,
1439d643ce63SMatthew Knepley                MATOP_ILUDT_FACTOR=49,
1440ff1cced0SMatthew Knepley                MATOP_SET_BLOCK_SIZE=50,
1441d643ce63SMatthew Knepley                MATOP_GET_ROW_IJ=51,
1442d643ce63SMatthew Knepley                MATOP_RESTORE_ROW_IJ=52,
1443d643ce63SMatthew Knepley                MATOP_GET_COLUMN_IJ=53,
1444d643ce63SMatthew Knepley                MATOP_RESTORE_COLUMN_IJ=54,
1445d643ce63SMatthew Knepley                MATOP_FDCOLORING_CREATE=55,
1446d643ce63SMatthew Knepley                MATOP_COLORING_PATCH=56,
1447d643ce63SMatthew Knepley                MATOP_SET_UNFACTORED=57,
1448d643ce63SMatthew Knepley                MATOP_PERMUTE=58,
1449d643ce63SMatthew Knepley                MATOP_SET_VALUES_BLOCKED=59,
1450d643ce63SMatthew Knepley                MATOP_GET_SUBMATRIX=60,
1451d643ce63SMatthew Knepley                MATOP_DESTROY=61,
1452d643ce63SMatthew Knepley                MATOP_VIEW=62,
1453ff1cced0SMatthew Knepley                MATOP_CONVERT_FROM=63,
1454d643ce63SMatthew Knepley                MATOP_USE_SCALED_FORM=64,
1455d643ce63SMatthew Knepley                MATOP_SCALE_SYSTEM=65,
1456d643ce63SMatthew Knepley                MATOP_UNSCALE_SYSTEM=66,
1457ba16a8cbSKris Buschelman                MATOP_SET_LOCAL_TO_GLOBAL_MAP=67,
1458d643ce63SMatthew Knepley                MATOP_SET_VALUES_LOCAL=68,
1459d643ce63SMatthew Knepley                MATOP_ZERO_ROWS_LOCAL=69,
1460ff1cced0SMatthew Knepley                MATOP_GET_ROW_MAX_ABS=70,
1461c87e5d42SMatthew Knepley                MATOP_GET_ROW_MIN_ABS=71,
1462c87e5d42SMatthew Knepley                MATOP_CONVERT=72,
1463c87e5d42SMatthew Knepley                MATOP_SET_COLORING=73,
1464c87e5d42SMatthew Knepley                MATOP_SET_VALUES_ADIC=74,
1465c87e5d42SMatthew Knepley                MATOP_SET_VALUES_ADIFOR=75,
1466c87e5d42SMatthew Knepley                MATOP_FD_COLORING_APPLY=76,
1467c87e5d42SMatthew Knepley                MATOP_SET_FROM_OPTIONS=77,
1468c87e5d42SMatthew Knepley                MATOP_MULT_CON=78,
1469c87e5d42SMatthew Knepley                MATOP_MULT_TRANSPOSE_CON=79,
1470d643ce63SMatthew Knepley                MATOP_PERMUTE_SPARSIFY=80,
1471d643ce63SMatthew Knepley                MATOP_MULT_MULTIPLE=81,
147241acf15aSKris Buschelman                MATOP_SOLVE_MULTIPLE=82,
147341acf15aSKris Buschelman                MATOP_GET_INERTIA=83,
147441acf15aSKris Buschelman                MATOP_LOAD=84,
147541acf15aSKris Buschelman                MATOP_IS_SYMMETRIC=85,
147641acf15aSKris Buschelman                MATOP_IS_HERMITIAN=86,
147741acf15aSKris Buschelman                MATOP_IS_STRUCTURALLY_SYMMETRIC=87,
147841acf15aSKris Buschelman                MATOP_PB_RELAX=88,
147941acf15aSKris Buschelman                MATOP_GET_VECS=89,
148041acf15aSKris Buschelman                MATOP_MAT_MULT=90,
148141acf15aSKris Buschelman                MATOP_MAT_MULT_SYMBOLIC=91,
148241acf15aSKris Buschelman                MATOP_MAT_MULT_NUMERIC=92,
148341acf15aSKris Buschelman                MATOP_PTAP=93,
148441acf15aSKris Buschelman                MATOP_PTAP_SYMBOLIC=94,
1485bc011b1eSHong Zhang                MATOP_PTAP_NUMERIC=95,
1486bc011b1eSHong Zhang                MATOP_MAT_MULTTRANSPOSE=96,
1487bc011b1eSHong Zhang                MATOP_MAT_MULTTRANSPOSE_SYMBOLIC=97,
14887ba1a0bfSKris Buschelman                MATOP_MAT_MULTTRANSPOSE_NUMERIC=98,
14897ba1a0bfSKris Buschelman                MATOP_PTAP_SYMBOLIC_SEQAIJ=99,
14907ba1a0bfSKris Buschelman                MATOP_PTAP_NUMERIC_SEQAIJ=100,
14917ba1a0bfSKris Buschelman                MATOP_PTAP_SYMBOLIC_MPIAIJ=101,
149287d4246cSBarry Smith                MATOP_PTAP_NUMERIC_MPIAIJ=102,
1493ff1cced0SMatthew Knepley                MATOP_CONJUGATE=103,
1494ff1cced0SMatthew Knepley                MATOP_SET_SIZES=104,
1495f5edf698SHong Zhang                MATOP_SET_VALUES_ROW = 105,
1496ff1cced0SMatthew Knepley                MATOP_REAL_PART=106,
1497ff1cced0SMatthew Knepley                MATOP_IMAG_PART=107,
1498d5c63f83SSatish Balay                MATOP_GET_ROW_UTRIANGULAR=108,
14992bebee5dSHong Zhang                MATOP_RESTORE_ROW_UTRIANGULAR=109,
150069db28dcSHong Zhang                MATOP_MATSOLVE=110,
1501829201f2SHong Zhang                MATOP_GET_REDUNDANTMATRIX=111,
1502ff1cced0SMatthew Knepley                MATOP_GET_ROW_MIN=112,
1503ff1cced0SMatthew Knepley                MATOP_GET_COLUMN_VEC=113,
1504ff1cced0SMatthew Knepley                MATOP_MISSING_DIAGONAL=114,
1505ff1cced0SMatthew Knepley                MATOP_MATGETSEQNONZEROSTRUCTURE=115,
1506ff1cced0SMatthew Knepley                MATOP_DESTROY_SOLVER=116
1507fae171e0SBarry Smith              } MatOperation;
1508be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatHasOperation(Mat,MatOperation,PetscTruth*);
1509be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellSetOperation(Mat,MatOperation,void(*)(void));
1510be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellGetOperation(Mat,MatOperation,void(**)(void));
1511be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellSetContext(Mat,void*);
1512112a2221SBarry Smith 
151390ace30eSBarry Smith /*
151490ace30eSBarry Smith    Codes for matrices stored on disk. By default they are
151590ace30eSBarry Smith  stored in a universal format. By changing the format with
15167973ad04SBarry Smith  PetscViewerSetFormat(viewer,PETSC_VIEWER_NATIVE); the matrices will
151790ace30eSBarry Smith  be stored in a way natural for the matrix, for example dense matrices
151890ace30eSBarry Smith  would be stored as dense. Matrices stored this way may only be
151990ace30eSBarry Smith  read into matrices of the same time.
152090ace30eSBarry Smith */
152190ace30eSBarry Smith #define MATRIX_BINARY_FORMAT_DENSE -1
152290ace30eSBarry Smith 
1523be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetHashTableFactor(Mat,PetscReal);
1524be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatISGetLocalMat(Mat,Mat*);
15251f4e1ec7SBarry Smith 
1526d9274352SBarry Smith /*S
1527d9274352SBarry Smith      MatNullSpace - Object that removes a null space from a vector, i.e.
1528d9274352SBarry Smith          orthogonalizes the vector to a subsapce
1529d9274352SBarry Smith 
1530f7a9e4ceSBarry Smith    Level: advanced
1531d9274352SBarry Smith 
1532d9274352SBarry Smith   Concepts: matrix; linear operator, null space
1533d9274352SBarry Smith 
15346e1639daSBarry Smith   Users manual sections:
15356e1639daSBarry Smith .   sec_singular
15366e1639daSBarry Smith 
1537d9274352SBarry Smith .seealso:  MatNullSpaceCreate()
1538d9274352SBarry Smith S*/
153974637425SBarry Smith typedef struct _p_MatNullSpace* MatNullSpace;
1540d9274352SBarry Smith 
1541be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceCreate(MPI_Comm,PetscTruth,PetscInt,const Vec[],MatNullSpace*);
1542281d1b2eSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(Vec,void*),void*);
1543be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceDestroy(MatNullSpace);
1544be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceRemove(MatNullSpace,Vec,Vec*);
1545be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceAttach(Mat,MatNullSpace);
154695902228SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceTest(MatNullSpace,Mat,PetscTruth *);
154774637425SBarry Smith 
1548be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatReorderingSeqSBAIJ(Mat,IS);
1549be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPISBAIJSetHashTableFactor(Mat,PetscReal);
1550be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetColumnIndices(Mat,PetscInt *);
155146129b97SKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJInvertBlockDiagonal(Mat);
15523f1d51d7SBarry Smith 
1553be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMAIJ(Mat,PetscInt,Mat*);
1554be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMAIJRedimension(Mat,PetscInt,Mat*);
1555be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMAIJGetAIJ(Mat,Mat*);
1556c4f061fbSSatish Balay 
1557be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatComputeExplicitOperator(Mat,Mat*);
1558b0a32e0cSBarry Smith 
1559be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalScaleLocal(Mat,Vec);
156004f1ad80SBarry Smith 
1561fef1beadSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*);
15623ec795f1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetBase(Mat,Vec,Vec);
15633ec795f1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*);
1564e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*));
1565e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec));
1566e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDAddNullSpace(Mat,MatNullSpace);
1567e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt);
1568e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDResetHHistory(Mat);
1569e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctionError(Mat,PetscReal);
1570e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetPeriod(Mat,PetscInt);
1571e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDGetH(Mat,PetscScalar *);
1572e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFromOptions(Mat);
1573e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*);
1574e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*);
1575e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDDSSetUmin(Mat,PetscReal);
1576e884886fSBarry Smith 
1577e884886fSBarry Smith 
1578e884886fSBarry Smith typedef struct _p_MatMFFD* MatMFFD;
1579e884886fSBarry Smith 
1580e884886fSBarry Smith /*E
1581e884886fSBarry Smith     MatMFFDType - algorithm used to compute the h used in computing matrix-vector products via differencing of the function
1582e884886fSBarry Smith 
1583e884886fSBarry Smith    Level: beginner
1584e884886fSBarry Smith 
1585e884886fSBarry Smith .seealso: MatMFFDSetType(), MatMFFDRegister()
1586e884886fSBarry Smith E*/
1587a313700dSBarry Smith #define MatMFFDType char*
1588e884886fSBarry Smith #define MATMFFD_DS  "ds"
1589e884886fSBarry Smith #define MATMFFD_WP  "wp"
1590e884886fSBarry Smith 
1591a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetType(Mat,const MatMFFDType);
1592e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatMFFD));
1593e884886fSBarry Smith 
1594e884886fSBarry Smith /*MC
1595e884886fSBarry Smith    MatMFFDRegisterDynamic - Adds a method to the MatMFFD registry.
1596e884886fSBarry Smith 
1597e884886fSBarry Smith    Synopsis:
1598e884886fSBarry Smith    PetscErrorCode MatMFFDRegisterDynamic(char *name_solver,char *path,char *name_create,PetscErrorCode (*routine_create)(MatMFFD))
1599e884886fSBarry Smith 
1600e884886fSBarry Smith    Not Collective
1601e884886fSBarry Smith 
1602e884886fSBarry Smith    Input Parameters:
1603e884886fSBarry Smith +  name_solver - name of a new user-defined compute-h module
1604e884886fSBarry Smith .  path - path (either absolute or relative) the library containing this solver
1605e884886fSBarry Smith .  name_create - name of routine to create method context
1606e884886fSBarry Smith -  routine_create - routine to create method context
1607e884886fSBarry Smith 
1608e884886fSBarry Smith    Level: developer
1609e884886fSBarry Smith 
1610e884886fSBarry Smith    Notes:
1611e884886fSBarry Smith    MatMFFDRegisterDynamic() may be called multiple times to add several user-defined solvers.
1612e884886fSBarry Smith 
1613e884886fSBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
1614e884886fSBarry Smith    is ignored.
1615e884886fSBarry Smith 
1616e884886fSBarry Smith    Sample usage:
1617e884886fSBarry Smith .vb
1618e884886fSBarry Smith    MatMFFDRegisterDynamic("my_h",/home/username/my_lib/lib/libO/solaris/mylib.a,
1619e884886fSBarry Smith                "MyHCreate",MyHCreate);
1620e884886fSBarry Smith .ve
1621e884886fSBarry Smith 
1622e884886fSBarry Smith    Then, your solver can be chosen with the procedural interface via
1623e884886fSBarry Smith $     MatMFFDSetType(mfctx,"my_h")
1624e884886fSBarry Smith    or at runtime via the option
1625e884886fSBarry Smith $     -snes_mf_type my_h
1626e884886fSBarry Smith 
1627e884886fSBarry Smith .keywords: MatMFFD, register
1628e884886fSBarry Smith 
1629e884886fSBarry Smith .seealso: MatMFFDRegisterAll(), MatMFFDRegisterDestroy()
1630e884886fSBarry Smith M*/
1631e884886fSBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1632e884886fSBarry Smith #define MatMFFDRegisterDynamic(a,b,c,d) MatMFFDRegister(a,b,c,0)
1633e884886fSBarry Smith #else
1634e884886fSBarry Smith #define MatMFFDRegisterDynamic(a,b,c,d) MatMFFDRegister(a,b,c,d)
1635e884886fSBarry Smith #endif
1636e884886fSBarry Smith 
1637e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDRegisterAll(const char[]);
1638e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDRegisterDestroy(void);
1639e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDDefaultSetUmin(Mat,PetscReal);
1640e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDWPSetComputeNormU(Mat,PetscTruth);
1641e884886fSBarry Smith 
1642e884886fSBarry Smith 
1643be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *);
1644be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *);
16457dbadf16SMatthew Knepley 
1646e9fa29b7SSatish Balay PETSC_EXTERN_CXX_END
16472eac72dbSBarry Smith #endif
1648