xref: /petsc/include/petscmat.h (revision a981769725c6d790eb1ad7452c1aa2c8ba06b9bd)
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"
67f91d8e95SBarry Smith #define MATHYPRESTRUCT     "hyprestruct"
68d91e6319SBarry Smith 
699989ab13SBarry Smith /*E
70c7393fdbSBarry Smith     MatSolverPackage - String with the name of a PETSc matrix solver type.
719989ab13SBarry Smith 
729989ab13SBarry Smith     For example: "petsc" indicates what PETSc provides, "superlu" indicates either
739989ab13SBarry Smith        SuperLU or SuperLU_Dist etc.
749989ab13SBarry Smith 
759989ab13SBarry Smith 
769989ab13SBarry Smith    Level: beginner
779989ab13SBarry Smith 
785c9eb25fSBarry Smith .seealso: MatGetFactor(), Mat, MatSetType(), MatType
799989ab13SBarry Smith E*/
80c7393fdbSBarry Smith #define MatSolverPackage char*
81b24902e0SBarry Smith #define MAT_SOLVER_SPOOLES      "spooles"
82b24902e0SBarry Smith #define MAT_SOLVER_SUPERLU      "superlu"
83b24902e0SBarry Smith #define MAT_SOLVER_SUPERLU_DIST "superlu_dist"
84b24902e0SBarry Smith #define MAT_SOLVER_UMFPACK      "umfpack"
85b24902e0SBarry Smith #define MAT_SOLVER_ESSL         "essl"
86b24902e0SBarry Smith #define MAT_SOLVER_LUSOL        "lusol"
87b24902e0SBarry Smith #define MAT_SOLVER_MUMPS        "mumps"
883bf14a46SMatthew Knepley #define MAT_SOLVER_PASTIX       "pastix"
89b24902e0SBarry Smith #define MAT_SOLVER_DSCPACK      "dscpack"
90b24902e0SBarry Smith #define MAT_SOLVER_MATLAB       "matlab"
9143244d56SBarry Smith #define MAT_SOLVER_PETSC        "petsc"
92b6806ab0SHong Zhang #define MAT_SOLVER_PLAPACK      "plapack"
93b24902e0SBarry Smith 
94b24902e0SBarry Smith /*E
95b24902e0SBarry Smith     MatFactorType - indicates what type of factorization is requested
96b24902e0SBarry Smith 
97b24902e0SBarry Smith     Level: beginner
98b24902e0SBarry Smith 
99b24902e0SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
100b24902e0SBarry Smith 
101c7393fdbSBarry Smith .seealso: MatSolverPackage, MatGetFactor()
102b24902e0SBarry Smith E*/
103599ef60dSHong Zhang typedef enum {MAT_FACTOR_NONE, MAT_FACTOR_LU, MAT_FACTOR_CHOLESKY, MAT_FACTOR_ILU, MAT_FACTOR_ICC,MAT_FACTOR_ILUDT} MatFactorType;
104c7393fdbSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetFactor(Mat,const MatSolverPackage,MatFactorType,Mat*);
105db4efbfdSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetFactorAvailable(Mat,const MatSolverPackage,MatFactorType,PetscTruth*);
10635bd34faSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFactorGetSolverPackage(Mat,const MatSolverPackage*);
107b24902e0SBarry Smith 
1089989ab13SBarry Smith 
109c06d978dSMatthew Knepley /* Logging support */
110552e946dSBarry Smith #define    MAT_FILE_COOKIE 1211216    /* used to indicate matrices in binary files */
111be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_COOKIE;
112be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_FDCOLORING_COOKIE;
113be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_PARTITIONING_COOKIE;
114be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_NULLSPACE_COOKIE;
115e884886fSBarry Smith extern PetscCookie PETSCMAT_DLLEXPORT MATMFFD_COOKIE;
116c06d978dSMatthew Knepley 
117ceb03754SKris Buschelman /*E
118ceb03754SKris Buschelman     MatReuse - Indicates if matrices obtained from a previous call to MatGetSubMatrices()
119ceb03754SKris Buschelman      or MatGetSubMatrix() are to be reused to store the new matrix values.
120ceb03754SKris Buschelman 
121ceb03754SKris Buschelman     Level: beginner
122ceb03754SKris Buschelman 
123ceb03754SKris Buschelman    Any additions/changes here MUST also be made in include/finclude/petscmat.h
124ceb03754SKris Buschelman 
1250c451bc4SBarry Smith .seealso: MatGetSubMatrices(), MatGetSubMatrix(), MatDestroyMatrices(), MatConvert()
126ceb03754SKris Buschelman E*/
127ceb03754SKris Buschelman typedef enum {MAT_INITIAL_MATRIX,MAT_REUSE_MATRIX} MatReuse;
128ceb03754SKris Buschelman 
1295494a064SHong Zhang /*E
1305494a064SHong Zhang     MatGetSubMatrixOption - Indicates if matrices obtained from a call to MatGetSubMatrices()
131829201f2SHong Zhang      include the matrix values. Currently it is only used by MatGetSeqNonzerostructure().
1325494a064SHong Zhang 
1335494a064SHong Zhang     Level: beginner
1345494a064SHong Zhang 
135829201f2SHong Zhang .seealso: MatGetSeqNonzerostructure()
1365494a064SHong Zhang E*/
1375494a064SHong Zhang typedef enum {MAT_DO_NOT_GET_VALUES,MAT_GET_VALUES} MatGetSubMatrixOption;
1385494a064SHong Zhang 
139e5bd5246SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInitializePackage(const char[]);
140c06d978dSMatthew Knepley 
141f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreate(MPI_Comm,Mat*);
142a69a960dSMatthew Knepley PetscPolymorphicFunction(MatCreate,(MPI_Comm comm),(comm,&A),Mat,A)
143a69a960dSMatthew Knepley PetscPolymorphicFunction(MatCreate,(),(PETSC_COMM_WORLD,&A),Mat,A)
144f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetSizes(Mat,PetscInt,PetscInt,PetscInt,PetscInt);
145a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetType(Mat,const MatType);
146be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetFromOptions(Mat);
147be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUpPreallocation(Mat);
148be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRegisterAll(const char[]);
149be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat));
15030de9b25SBarry Smith 
151f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetOptionsPrefix(Mat,const char[]);
152f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAppendOptionsPrefix(Mat,const char[]);
153f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOptionsPrefix(Mat,const char*[]);
154f69a0ea3SMatthew Knepley 
15530de9b25SBarry Smith /*MC
15630de9b25SBarry Smith    MatRegisterDynamic - Adds a new matrix type
15730de9b25SBarry Smith 
15830de9b25SBarry Smith    Synopsis:
159c1ac3661SBarry Smith    PetscErrorCode MatRegisterDynamic(char *name,char *path,char *name_create,PetscErrorCode (*routine_create)(Mat))
16030de9b25SBarry Smith 
16130de9b25SBarry Smith    Not Collective
16230de9b25SBarry Smith 
16330de9b25SBarry Smith    Input Parameters:
16430de9b25SBarry Smith +  name - name of a new user-defined matrix type
16530de9b25SBarry Smith .  path - path (either absolute or relative) the library containing this solver
16630de9b25SBarry Smith .  name_create - name of routine to create method context
16730de9b25SBarry Smith -  routine_create - routine to create method context
16830de9b25SBarry Smith 
16930de9b25SBarry Smith    Notes:
17030de9b25SBarry Smith    MatRegisterDynamic() may be called multiple times to add several user-defined solvers.
17130de9b25SBarry Smith 
17230de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
17330de9b25SBarry Smith    is ignored.
17430de9b25SBarry Smith 
17530de9b25SBarry Smith    Sample usage:
17630de9b25SBarry Smith .vb
17730de9b25SBarry Smith    MatRegisterDynamic("my_mat",/home/username/my_lib/lib/libO/solaris/mylib.a,
17830de9b25SBarry Smith                "MyMatCreate",MyMatCreate);
17930de9b25SBarry Smith .ve
18030de9b25SBarry Smith 
18130de9b25SBarry Smith    Then, your solver can be chosen with the procedural interface via
18230de9b25SBarry Smith $     MatSetType(Mat,"my_mat")
18330de9b25SBarry Smith    or at runtime via the option
18430de9b25SBarry Smith $     -mat_type my_mat
18530de9b25SBarry Smith 
18630de9b25SBarry Smith    Level: advanced
18730de9b25SBarry Smith 
188ab901514SBarry Smith    Notes: ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values.
18930de9b25SBarry Smith          If your function is not being put into a shared library then use VecRegister() instead
19030de9b25SBarry Smith 
19130de9b25SBarry Smith .keywords: Mat, register
19230de9b25SBarry Smith 
19330de9b25SBarry Smith .seealso: MatRegisterAll(), MatRegisterDestroy()
19430de9b25SBarry Smith 
19530de9b25SBarry Smith M*/
196273d9f13SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
197273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,0)
198273d9f13SBarry Smith #else
199273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,d)
20030de9b25SBarry Smith #endif
20130de9b25SBarry Smith 
202273d9f13SBarry Smith extern PetscTruth MatRegisterAllCalled;
203b0a32e0cSBarry Smith extern PetscFList MatList;
20428988994SBarry Smith 
2053b224e63SBarry Smith /*E
2063b224e63SBarry Smith     MatStructure - Indicates if the matrix has the same nonzero structure
2073b224e63SBarry Smith 
2083b224e63SBarry Smith     Level: beginner
2093b224e63SBarry Smith 
2103b224e63SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
2113b224e63SBarry Smith 
2123b224e63SBarry Smith .seealso: MatCopy(), KSPSetOperators(), PCSetOperators()
2133b224e63SBarry Smith E*/
2143b224e63SBarry Smith typedef enum {SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER,SUBSET_NONZERO_PATTERN} MatStructure;
2153b224e63SBarry Smith 
216be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqDense(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*);
217be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIDense(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*);
218be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
219ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,nz,nnz,&A),Mat,A)
220ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,&A),Mat,A)
221ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,0,nnz,&A),Mat,A)
222ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,&A),Mat,A)
223ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,A))
224ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,m,n,0,nnz,A))
225ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,A))
226be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
227ba966639SSatish 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)
228ba966639SSatish 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)
229ba966639SSatish 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)
230ba966639SSatish 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)
231ba966639SSatish 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))
232ba966639SSatish 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))
233ba966639SSatish 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))
234ba966639SSatish 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)
235ba966639SSatish 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)
236ba966639SSatish 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)
237ba966639SSatish 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)
238ba966639SSatish 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))
239ba966639SSatish 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))
240ba966639SSatish 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))
24163ab93bfSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *);
2428d7a6e47SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJWithSplitArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],PetscInt[],PetscInt[],PetscScalar[],Mat*);
2438d7a6e47SBarry Smith 
244be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIRowbs(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
245be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
246ba966639SSatish 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)
247ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A)
248ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A)
249ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A)
250ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A))
251ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A))
252ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A))
253be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
254ba966639SSatish 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)
255ba966639SSatish 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)
256ba966639SSatish 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)
257ba966639SSatish 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)
258ba966639SSatish 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))
259ba966639SSatish 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))
260ba966639SSatish 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))
261ba966639SSatish 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)
262ba966639SSatish 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)
263ba966639SSatish 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)
264ba966639SSatish 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)
265ba966639SSatish 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))
266ba966639SSatish 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))
267ba966639SSatish 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))
268be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAdj(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscInt[],Mat*);
269be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
270ba966639SSatish 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)
271ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A)
272ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A)
273ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A)
274ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A))
275ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A))
276ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A))
277be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPISBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
278ba966639SSatish 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)
279ba966639SSatish 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)
280ba966639SSatish 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)
281ba966639SSatish 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)
282ba966639SSatish 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))
283ba966639SSatish 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))
284ba966639SSatish 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))
285ba966639SSatish 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)
286ba966639SSatish 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)
287ba966639SSatish 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)
288ba966639SSatish 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)
289ba966639SSatish 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))
290ba966639SSatish 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))
291ba966639SSatish 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))
292be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateShell(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,void *,Mat*);
293ba966639SSatish 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)
294ba966639SSatish Balay PetscPolymorphicFunction(MatCreateShell,(PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx),(PETSC_COMM_WORLD,m,n,M,N,ctx,&A),Mat,A)
295be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateAdic(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,void (*)(void),Mat*);
296be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateNormal(Mat,Mat*);
297ba966639SSatish Balay PetscPolymorphicFunction(MatCreateNormal,(Mat mat),(mat,&A),Mat,A)
298ab92ecdeSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateLRC(Mat,Mat,Mat,Mat*);
299284134d9SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateIS(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,ISLocalToGlobalMapping,Mat*);
3001472f72bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
3011472f72bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPICRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
3022a6744ebSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateScatter(MPI_Comm,VecScatter,Mat*);
3032a6744ebSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScatterSetVecScatter(Mat,VecScatter);
3042a6744ebSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScatterGetVecScatter(Mat,VecScatter*);
3058cdf2d9bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateBlockMat(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,Mat*);
306793850ffSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompositeAddMat(Mat,Mat);
307b52f573bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompositeMerge(Mat);
308793850ffSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateComposite(MPI_Comm,PetscInt,const Mat*,Mat*);
3096d7c1e57SBarry Smith typedef enum {MAT_COMPOSITE_ADDITIVE,MAT_COMPOSITE_MULTIPLICATIVE} MatCompositeType;
3106d7c1e57SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompositeSetType(Mat,MatCompositeType);
3116d7c1e57SBarry Smith 
3125a7f1df3SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqFFTW(MPI_Comm,PetscInt,const PetscInt[],Mat*);
313f20085c4SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateTranspose(Mat,Mat*);
3141472f72bSBarry Smith 
3151d6018f0SLisandro Dalcin EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreatePython(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const char[],Mat*);
3161d6018f0SLisandro Dalcin EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPythonSetType(Mat,const char[]);
3171d6018f0SLisandro Dalcin 
3181d6018f0SLisandro Dalcin 
319f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUp(Mat);
320be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroy(Mat);
32121c89e3eSBarry Smith 
3220c451bc4SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConjugate(Mat);
32399cafbc1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRealPart(Mat);
32499cafbc1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatImaginaryPart(Mat);
325bfc7d00aSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetDiagonalBlock(Mat,PetscTruth*,MatReuse,Mat*);
32699cafbc1SBarry Smith 
3278ed539a5SBarry Smith /* ------------------------------------------------------------*/
328be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
329be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
33087d4246cSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesRow(Mat,PetscInt,const PetscScalar[]);
331f259bd47SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesRowLocal(Mat,PetscInt,const PetscScalar[]);
33284cb2905SBarry Smith 
3332ef4de8bSBarry Smith /*S
3342ef4de8bSBarry Smith      MatStencil - Data structure (C struct) for storing information about a single row or
3352ef4de8bSBarry Smith         column of a matrix as index on an associated grid.
3362ef4de8bSBarry Smith 
3372ef4de8bSBarry Smith    Level: beginner
3382ef4de8bSBarry Smith 
3392ef4de8bSBarry Smith   Concepts: matrix; linear operator
3402ef4de8bSBarry Smith 
341d7bd93baSBarry Smith .seealso:  MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockStencil()
3422ef4de8bSBarry Smith S*/
343435da068SBarry Smith typedef struct {
344c1ac3661SBarry Smith   PetscInt k,j,i,c;
345435da068SBarry Smith } MatStencil;
3462ef4de8bSBarry Smith 
347be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
348be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
349be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt);
350435da068SBarry Smith 
351be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetColoring(Mat,ISColoring);
352be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesAdic(Mat,void*);
353be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesAdifor(Mat,PetscInt,void*);
3543a7fca6bSBarry Smith 
355d91e6319SBarry Smith /*E
356d91e6319SBarry Smith     MatAssemblyType - Indicates if the matrix is now to be used, or if you plan
357d91e6319SBarry Smith      to continue to add values to it
358d91e6319SBarry Smith 
359d91e6319SBarry Smith     Level: beginner
360d91e6319SBarry Smith 
361d91e6319SBarry Smith .seealso: MatAssemblyBegin(), MatAssemblyEnd()
362d91e6319SBarry Smith E*/
3636d4a8577SBarry Smith typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType;
364be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssemblyBegin(Mat,MatAssemblyType);
365be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssemblyEnd(Mat,MatAssemblyType);
366be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssembled(Mat,PetscTruth*);
3674f9c727eSBarry Smith 
3681d73ed98SBarry Smith 
36930de9b25SBarry Smith 
370d91e6319SBarry Smith /*E
371d91e6319SBarry Smith     MatOption - Options that may be set for a matrix and its behavior or storage
372d91e6319SBarry Smith 
373d91e6319SBarry Smith     Level: beginner
374d91e6319SBarry Smith 
3750a835dfdSSatish Balay    Any additions/changes here MUST also be made in include/finclude/petscmat.h
376d91e6319SBarry Smith 
377d91e6319SBarry Smith .seealso: MatSetOption()
378d91e6319SBarry Smith E*/
379512a5fc5SBarry Smith typedef enum {MAT_ROW_ORIENTED,MAT_NEW_NONZERO_LOCATIONS,
3804e0d8c25SBarry Smith               MAT_SYMMETRIC,
3814e0d8c25SBarry Smith               MAT_STRUCTURALLY_SYMMETRIC,
3828047e77bSBarry Smith               MAT_NEW_DIAGONALS,MAT_IGNORE_OFF_PROC_ENTRIES,
3834e0d8c25SBarry Smith               MAT_NEW_NONZERO_LOCATION_ERR,
3844e0d8c25SBarry Smith               MAT_NEW_NONZERO_ALLOCATION_ERR,MAT_USE_HASH_TABLE,
385*a9817697SBarry Smith               MAT_KEEP_NONZERO_PATTERN,MAT_IGNORE_ZERO_ENTRIES,
3864e0d8c25SBarry Smith               MAT_USE_INODES,
3874e0d8c25SBarry Smith               MAT_HERMITIAN,
3884e0d8c25SBarry Smith               MAT_SYMMETRY_ETERNAL,
3894e0d8c25SBarry Smith               MAT_USE_COMPRESSEDROW,
3904e0d8c25SBarry Smith               MAT_IGNORE_LOWER_TRIANGULAR,MAT_ERROR_LOWER_TRIANGULAR,
39128b2fa4aSMatthew Knepley               MAT_GETROW_UPPERTRIANGULAR,MAT_UNUSED_NONZERO_LOCATION_ERR,
39228b2fa4aSMatthew Knepley               NUM_MAT_OPTIONS} MatOption;
393290bbb0aSBarry Smith extern const char *MatOptions[];
3944e0d8c25SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetOption(Mat,MatOption,PetscTruth);
395a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetType(Mat,const MatType*);
396a313700dSBarry Smith PetscPolymorphicFunction(MatGetType,(Mat mat),(mat,&t),const MatType,t)
39784cb2905SBarry Smith 
398be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]);
399be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
400be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
401f5edf698SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowUpperTriangular(Mat);
402f5edf698SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRowUpperTriangular(Mat);
403be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
404be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
405be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumnVector(Mat,Vec,PetscInt);
406be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetArray(Mat,PetscScalar *[]);
407ba966639SSatish Balay PetscPolymorphicFunction(MatGetArray,(Mat mat),(mat,&a),PetscScalar*,a)
408be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreArray(Mat,PetscScalar *[]);
409be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBlockSize(Mat,PetscInt *);
410ba966639SSatish Balay PetscPolymorphicFunction(MatGetBlockSize,(Mat mat),(mat,&a),PetscInt,a)
411be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetBlockSize(Mat,PetscInt);
4127b80b807SBarry Smith 
4131620fd73SBarry Smith 
414be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMult(Mat,Vec,Vec);
415be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultAdd(Mat,Vec,Vec,Vec);
416ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
417be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTranspose(Mat,Vec,Vec);
418be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsTranspose(Mat,Mat,PetscReal,PetscTruth*);
419ba966639SSatish Balay PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B,PetscReal tol),(A,B,tol,&t),PetscTruth,t)
420e37f64cdSBarry Smith PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B),(A,B,0,&t),PetscTruth,t)
4211cbb95d3SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitianTranspose(Mat,Mat,PetscReal,PetscTruth*);
422be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeAdd(Mat,Vec,Vec,Vec);
423ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultTransposeAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
424be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultConstrained(Mat,Vec,Vec);
425be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeConstrained(Mat,Vec,Vec);
4262bebee5dSHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatSolve(Mat,Mat,Mat);
427f5edf698SHong Zhang 
428d91e6319SBarry Smith /*E
429d91e6319SBarry Smith     MatDuplicateOption - Indicates if a duplicated sparse matrix should have
430d91e6319SBarry Smith   its numerical values copied over or just its nonzero structure.
431d91e6319SBarry Smith 
432d91e6319SBarry Smith     Level: beginner
433d91e6319SBarry Smith 
434d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
435d91e6319SBarry Smith 
436d91e6319SBarry Smith .seealso: MatDuplicate()
437d91e6319SBarry Smith E*/
4382e8a6d31SBarry Smith typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES} MatDuplicateOption;
4392e8a6d31SBarry Smith 
440a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert(Mat,const MatType,MatReuse,Mat*);
441a313700dSBarry Smith PetscPolymorphicFunction(MatConvert,(Mat A,const MatType t),(A,t,MAT_INITIAL_MATRIX,&a),Mat,a)
442be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDuplicate(Mat,MatDuplicateOption,Mat*);
443ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A,MatDuplicateOption o),(A,o,&a),Mat,a)
444ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A),(A,MAT_COPY_VALUES,&a),Mat,a)
44594a9d846SBarry Smith 
446cb5b572fSBarry Smith 
447be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCopy(Mat,Mat,MatStructure);
448be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatView(Mat,PetscViewer);
449be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsSymmetric(Mat,PetscReal,PetscTruth*);
450ba966639SSatish Balay PetscPolymorphicFunction(MatIsSymmetric,(Mat A,PetscReal tol),(A,tol,&t),PetscTruth,t)
451e37f64cdSBarry Smith PetscPolymorphicFunction(MatIsSymmetric,(Mat A),(A,0,&t),PetscTruth,t)
452be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsStructurallySymmetric(Mat,PetscTruth*);
453ba966639SSatish Balay PetscPolymorphicFunction(MatIsStructurallySymmetric,(Mat A),(A,&t),PetscTruth,t)
4541cbb95d3SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitian(Mat,PetscReal,PetscTruth*);
4551cbb95d3SBarry Smith PetscPolymorphicFunction(MatIsHermitian,(Mat A),(A,0,&t),PetscTruth,t)
456be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsSymmetricKnown(Mat,PetscTruth*,PetscTruth*);
457be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitianKnown(Mat,PetscTruth*,PetscTruth*);
45864c62002SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMissingDiagonal(Mat,PetscTruth *,PetscInt *);
459a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLoad(PetscViewer,const MatType,Mat*);
460a313700dSBarry Smith PetscPolymorphicFunction(MatLoad,(PetscViewer v,const MatType t),(v,t,&a),Mat,a)
4617b80b807SBarry Smith 
4628f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *);
4638f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRowIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *);
4648f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumnIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *);
4658f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreColumnIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *);
466d4fbbf0eSBarry Smith 
467d91e6319SBarry Smith /*S
468d91e6319SBarry Smith      MatInfo - Context of matrix information, used with MatGetInfo()
469d91e6319SBarry Smith 
470d91e6319SBarry Smith    In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE
471d91e6319SBarry Smith 
472d91e6319SBarry Smith    Level: intermediate
473d91e6319SBarry Smith 
474d91e6319SBarry Smith   Concepts: matrix^nonzero information
475d91e6319SBarry Smith 
476d9274352SBarry Smith .seealso:  MatGetInfo(), MatInfoType
477d91e6319SBarry Smith S*/
4784e220ebcSLois Curfman McInnes typedef struct {
479b0a32e0cSBarry Smith   PetscLogDouble block_size;                         /* block size */
480b0a32e0cSBarry Smith   PetscLogDouble nz_allocated,nz_used,nz_unneeded;   /* number of nonzeros */
481b0a32e0cSBarry Smith   PetscLogDouble memory;                             /* memory allocated */
482b0a32e0cSBarry Smith   PetscLogDouble assemblies;                         /* number of matrix assemblies called */
483b0a32e0cSBarry Smith   PetscLogDouble mallocs;                            /* number of mallocs during MatSetValues() */
484b0a32e0cSBarry Smith   PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */
485b0a32e0cSBarry Smith   PetscLogDouble factor_mallocs;                     /* number of mallocs during factorization */
4864e220ebcSLois Curfman McInnes } MatInfo;
4874e220ebcSLois Curfman McInnes 
488d9274352SBarry Smith /*E
489d9274352SBarry Smith     MatInfoType - Indicates if you want information about the local part of the matrix,
490d9274352SBarry Smith      the entire parallel matrix or the maximum over all the local parts.
491d9274352SBarry Smith 
492d9274352SBarry Smith     Level: beginner
493d9274352SBarry Smith 
494d9274352SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
495d9274352SBarry Smith 
496d9274352SBarry Smith .seealso: MatGetInfo(), MatInfo
497d9274352SBarry Smith E*/
4987b80b807SBarry Smith typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType;
499be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetInfo(Mat,MatInfoType,MatInfo*);
500be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatValid(Mat,PetscTruth*);
501be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetDiagonal(Mat,Vec);
502985db425SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMax(Mat,Vec,PetscInt[]);
503985db425SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMin(Mat,Vec,PetscInt[]);
504985db425SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMaxAbs(Mat,Vec,PetscInt[]);
505c87e5d42SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMinAbs(Mat,Vec,PetscInt[]);
50686697f06SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowSum(Mat,Vec);
507fc4dec0aSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatTranspose(Mat,MatReuse,Mat*);
508fc4dec0aSBarry Smith PetscPolymorphicFunction(MatTranspose,(Mat A),(A,MAT_INITIAL_MATRIX,&t),Mat,t)
509be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPermute(Mat,IS,IS,Mat *);
510ba966639SSatish Balay PetscPolymorphicFunction(MatPermute,(Mat A,IS is1,IS is2),(A,is1,is2,&t),Mat,t)
511be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPermuteSparsify(Mat,PetscInt,PetscReal,PetscReal,IS,IS,Mat *);
512be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalScale(Mat,Vec,Vec);
513be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalSet(Mat,Vec,InsertMode);
514be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatEqual(Mat,Mat,PetscTruth*);
515ba966639SSatish Balay PetscPolymorphicFunction(MatEqual,(Mat A,Mat B),(A,B,&t),PetscTruth,t)
516be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultEqual(Mat,Mat,PetscInt,PetscTruth*);
517be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultAddEqual(Mat,Mat,PetscInt,PetscTruth*);
518be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeEqual(Mat,Mat,PetscInt,PetscTruth*);
519be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscTruth*);
5207b80b807SBarry Smith 
521be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNorm(Mat,NormType,PetscReal *);
522ba966639SSatish Balay PetscPolymorphicFunction(MatNorm,(Mat A,NormType t),(A,t,&n),PetscReal,n)
523be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroEntries(Mat);
524f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar);
525f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsIS(Mat,IS,PetscScalar);
526f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroColumns(Mat,PetscInt,const PetscInt [],const PetscScalar*);
527f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroColumnsIS(Mat,IS,const PetscScalar*);
5287b80b807SBarry Smith 
529be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatUseScaledForm(Mat,PetscTruth);
530be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScaleSystem(Mat,Vec,Vec);
531be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatUnScaleSystem(Mat,Vec,Vec);
5325ef9f2a5SBarry Smith 
533be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSize(Mat,PetscInt*,PetscInt*);
534be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalSize(Mat,PetscInt*,PetscInt*);
535be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRange(Mat,PetscInt*,PetscInt*);
5368ee2e534SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRanges(Mat,const PetscInt**);
537f4b8409dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRangeColumn(Mat,PetscInt*,PetscInt*);
53872ca8751SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRangesColumn(Mat,const PetscInt**);
5397b80b807SBarry Smith 
540be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
541be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroyMatrices(PetscInt,Mat *[]);
5424aa3045dSJed Brown EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrix(Mat,IS,IS,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);
5791c741599SHong Zhang 
580f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScale(Mat,PetscScalar);
581f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShift(Mat,PetscScalar);
5827b80b807SBarry Smith 
583be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping);
584be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetLocalToGlobalMappingBlock(Mat,ISLocalToGlobalMapping);
585f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar);
586f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsLocalIS(Mat,IS,PetscScalar);
587be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
588be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
589052efed2SBarry Smith 
590be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStashSetInitialSize(Mat,PetscInt,PetscInt);
591be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
59290f02eecSBarry Smith 
593be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInterpolate(Mat,Vec,Vec);
5940c451bc4SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInterpolateAdd(Mat,Vec,Vec,Vec);
595ba966639SSatish Balay PetscPolymorphicSubroutine(MatInterpolateAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
596be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestrict(Mat,Vec,Vec);
597be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetVecs(Mat,Vec*,Vec*);
59869db28dcSHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRedundantMatrix(Mat,PetscInt,MPI_Comm,PetscInt,MatReuse,Mat*);
599bd481603SBarry Smith 
600bd481603SBarry Smith /*MC
6011620fd73SBarry Smith    MatSetValue - Set a single entry into a matrix.
6021620fd73SBarry Smith 
6031620fd73SBarry Smith    Not collective
6041620fd73SBarry Smith 
6051620fd73SBarry Smith    Input Parameters:
6061620fd73SBarry Smith +  m - the matrix
6071620fd73SBarry Smith .  row - the row location of the entry
6081620fd73SBarry Smith .  col - the column location of the entry
6091620fd73SBarry Smith .  value - the value to insert
6101620fd73SBarry Smith -  mode - either INSERT_VALUES or ADD_VALUES
6111620fd73SBarry Smith 
6121620fd73SBarry Smith    Notes:
6131620fd73SBarry Smith    For efficiency one should use MatSetValues() and set several or many
6141620fd73SBarry Smith    values simultaneously if possible.
6151620fd73SBarry Smith 
6161620fd73SBarry Smith    Level: beginner
6171620fd73SBarry Smith 
6181620fd73SBarry Smith .seealso: MatSetValues(), MatSetValueLocal()
6191620fd73SBarry Smith M*/
6201620fd73SBarry 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);}
6211620fd73SBarry Smith 
6221620fd73SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatGetValue(Mat v,PetscInt i,PetscInt j,PetscScalar *va) {return MatGetValues(v,1,&i,1,&j,va);}
6231620fd73SBarry Smith 
6241620fd73SBarry 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);}
6251620fd73SBarry Smith 
6261620fd73SBarry Smith /*MC
6270d2a7ff2SSatish Balay    MatPreallocateInitialize - Begins the block of code that will count the number of nonzeros per
628bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
629bd481603SBarry Smith 
630bd481603SBarry Smith    Synopsis:
631c1ac3661SBarry Smith    PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)
632bd481603SBarry Smith 
633bd481603SBarry Smith    Collective on MPI_Comm
634bd481603SBarry Smith 
635bd481603SBarry Smith    Input Parameters:
636bd481603SBarry Smith +  comm - the communicator that will share the eventually allocated matrix
637859fcb39SBarry Smith .  nrows - the number of LOCAL rows in the matrix
638859fcb39SBarry Smith -  ncols - the number of LOCAL columns in the matrix
639bd481603SBarry Smith 
640bd481603SBarry Smith    Output Parameters:
641bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
642bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
643bd481603SBarry Smith 
644bd481603SBarry Smith 
645bd481603SBarry Smith    Level: intermediate
646bd481603SBarry Smith 
647bd481603SBarry Smith    Notes:
648bd481603SBarry Smith    See the chapter in the users manual on performance for more details
649bd481603SBarry Smith 
6501d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
651bd481603SBarry Smith 
652bd481603SBarry Smith    Use MatPreallocateInitializeSymmetric() for symmetric matrices (MPISBAIJ matrices)
653bd481603SBarry Smith 
6541620fd73SBarry Smith    This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize().
6551620fd73SBarry Smith 
656bd481603SBarry Smith   Concepts: preallocation^Matrix
657bd481603SBarry Smith 
658bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
659bd481603SBarry Smith           MatPreallocateInitializeSymmetric(), MatPreallocateSymmetricSetLocal()
660bd481603SBarry Smith M*/
661c4f061fbSSatish Balay #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \
6627c922b88SBarry Smith { \
663a89bc211SBarry Smith   PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__start,__end; \
664a89bc211SBarry Smith   _4_ierr = PetscMalloc2(__nrows,PetscInt,&dnz,__nrows,PetscInt,&onz);CHKERRQ(_4_ierr); \
665a89bc211SBarry Smith   _4_ierr = PetscMemzero(dnz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
666a89bc211SBarry Smith   _4_ierr = PetscMemzero(onz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
667c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ctmp;\
668a89bc211SBarry Smith   _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows;
6697c922b88SBarry Smith 
670bd481603SBarry Smith /*MC
6710d2a7ff2SSatish Balay    MatPreallocateSymmetricInitialize - Begins the block of code that will count the number of nonzeros per
672bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
673bd481603SBarry Smith 
674bd481603SBarry Smith    Synopsis:
675c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)
676bd481603SBarry Smith 
677bd481603SBarry Smith    Collective on MPI_Comm
678bd481603SBarry Smith 
679bd481603SBarry Smith    Input Parameters:
680bd481603SBarry Smith +  comm - the communicator that will share the eventually allocated matrix
681859fcb39SBarry Smith .  nrows - the number of LOCAL rows in the matrix
682859fcb39SBarry Smith -  ncols - the number of LOCAL columns in the matrix
683bd481603SBarry Smith 
684bd481603SBarry Smith    Output Parameters:
685bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
686bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
687bd481603SBarry Smith 
688bd481603SBarry Smith 
689bd481603SBarry Smith    Level: intermediate
690bd481603SBarry Smith 
691bd481603SBarry Smith    Notes:
692bd481603SBarry Smith    See the chapter in the users manual on performance for more details
693bd481603SBarry Smith 
6941d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
695bd481603SBarry Smith 
6961620fd73SBarry Smith    This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize().
6971620fd73SBarry Smith 
698bd481603SBarry Smith   Concepts: preallocation^Matrix
699bd481603SBarry Smith 
700bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
701bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal()
702bd481603SBarry Smith M*/
703222b16d4SBarry Smith #define MatPreallocateSymmetricInitialize(comm,nrows,ncols,dnz,onz) 0; \
704222b16d4SBarry Smith { \
705a89bc211SBarry Smith   PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__end; \
706a89bc211SBarry Smith   _4_ierr = PetscMalloc2(__nrows,PetscInt,&dnz,__nrows,PetscInt,&onz);CHKERRQ(_4_ierr); \
707a89bc211SBarry Smith   _4_ierr = PetscMemzero(dnz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
708a89bc211SBarry Smith   _4_ierr = PetscMemzero(onz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
709c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr);\
710a89bc211SBarry Smith   _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows;
711222b16d4SBarry Smith 
712bd481603SBarry Smith /*MC
713bd481603SBarry Smith    MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
714bd481603SBarry Smith        inserted using a local number of the rows and columns
715bd481603SBarry Smith 
716bd481603SBarry Smith    Synopsis:
717c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
718bd481603SBarry Smith 
719bd481603SBarry Smith    Not Collective
720bd481603SBarry Smith 
721bd481603SBarry Smith    Input Parameters:
722bd481603SBarry Smith +  map - the mapping between local numbering and global numbering
723bd481603SBarry Smith .  nrows - the number of rows indicated
7241d73ed98SBarry Smith .  rows - the indices of the rows
725bd481603SBarry Smith .  ncols - the number of columns in the matrix
726bd481603SBarry Smith .  cols - the columns indicated
727bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
728bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
729bd481603SBarry Smith 
730bd481603SBarry Smith 
731bd481603SBarry Smith    Level: intermediate
732bd481603SBarry Smith 
733bd481603SBarry Smith    Notes:
734bd481603SBarry Smith    See the chapter in the users manual on performance for more details
735bd481603SBarry Smith 
7361d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
737bd481603SBarry Smith 
738bd481603SBarry Smith   Concepts: preallocation^Matrix
739bd481603SBarry Smith 
740bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
741bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal()
742bd481603SBarry Smith M*/
743c4f061fbSSatish Balay #define MatPreallocateSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\
744c4f061fbSSatish Balay {\
745c1ac3661SBarry Smith   PetscInt __l;\
746ef66eb69SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\
747ef66eb69SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\
748c4f061fbSSatish Balay   for (__l=0;__l<nrows;__l++) {\
749ef66eb69SBarry Smith     _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
750c4f061fbSSatish Balay   }\
751c4f061fbSSatish Balay }
752c4f061fbSSatish Balay 
753bd481603SBarry Smith /*MC
754bd481603SBarry Smith    MatPreallocateSymmetricSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
755bd481603SBarry Smith        inserted using a local number of the rows and columns
756bd481603SBarry Smith 
757bd481603SBarry Smith    Synopsis:
758c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
759bd481603SBarry Smith 
760bd481603SBarry Smith    Not Collective
761bd481603SBarry Smith 
762bd481603SBarry Smith    Input Parameters:
763bd481603SBarry Smith +  map - the mapping between local numbering and global numbering
764bd481603SBarry Smith .  nrows - the number of rows indicated
7651d73ed98SBarry Smith .  rows - the indices of the rows
766bd481603SBarry Smith .  ncols - the number of columns in the matrix
767bd481603SBarry Smith .  cols - the columns indicated
768bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
769bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
770bd481603SBarry Smith 
771bd481603SBarry Smith 
772bd481603SBarry Smith    Level: intermediate
773bd481603SBarry Smith 
774bd481603SBarry Smith    Notes:
775bd481603SBarry Smith    See the chapter in the users manual on performance for more details
776bd481603SBarry Smith 
777bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
778bd481603SBarry Smith 
779bd481603SBarry Smith   Concepts: preallocation^Matrix
780bd481603SBarry Smith 
781bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
782bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
783bd481603SBarry Smith M*/
784d3d32019SBarry Smith #define MatPreallocateSymmetricSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\
785d3d32019SBarry Smith {\
786c1ac3661SBarry Smith   PetscInt __l;\
787d3d32019SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\
788d3d32019SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\
789d3d32019SBarry Smith   for (__l=0;__l<nrows;__l++) {\
790d3d32019SBarry Smith     _4_ierr = MatPreallocateSymmetricSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
791d3d32019SBarry Smith   }\
792d3d32019SBarry Smith }
793d3d32019SBarry Smith 
794bd481603SBarry Smith /*MC
795bd481603SBarry Smith    MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
796bd481603SBarry Smith        inserted using a local number of the rows and columns
797bd481603SBarry Smith 
798bd481603SBarry Smith    Synopsis:
799c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
800bd481603SBarry Smith 
801bd481603SBarry Smith    Not Collective
802bd481603SBarry Smith 
803bd481603SBarry Smith    Input Parameters:
80464075487SBarry Smith +  row - the row
805bd481603SBarry Smith .  ncols - the number of columns in the matrix
806a50f8bf6SHong Zhang -  cols - the columns indicated
807a50f8bf6SHong Zhang 
808a50f8bf6SHong Zhang    Output Parameters:
809a50f8bf6SHong Zhang +  dnz - the array that will be passed to the matrix preallocation routines
810bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
811bd481603SBarry Smith 
812bd481603SBarry Smith 
813bd481603SBarry Smith    Level: intermediate
814bd481603SBarry Smith 
815bd481603SBarry Smith    Notes:
816bd481603SBarry Smith    See the chapter in the users manual on performance for more details
817bd481603SBarry Smith 
818bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
819bd481603SBarry Smith 
8201620fd73SBarry Smith    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
8211620fd73SBarry Smith 
822bd481603SBarry Smith   Concepts: preallocation^Matrix
823bd481603SBarry Smith 
824bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
825bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
826bd481603SBarry Smith M*/
827c4f061fbSSatish Balay #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\
828c1ac3661SBarry Smith { PetscInt __i; \
8298ee2e534SBarry Smith   if (row < __rstart) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Trying to set preallocation for row %D less than first local row %D",row,__rstart);\
830a89bc211SBarry 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);\
8317c922b88SBarry Smith   for (__i=0; __i<nc; __i++) {\
83264075487SBarry Smith     if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \
8337cc688d7SBarry Smith     else dnz[row - __rstart]++;\
8347c922b88SBarry Smith   }\
8357c922b88SBarry Smith }
8367c922b88SBarry Smith 
837bd481603SBarry Smith /*MC
838bd481603SBarry Smith    MatPreallocateSymmetricSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
839bd481603SBarry Smith        inserted using a local number of the rows and columns
840bd481603SBarry Smith 
841bd481603SBarry Smith    Synopsis:
842c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
843bd481603SBarry Smith 
844bd481603SBarry Smith    Not Collective
845bd481603SBarry Smith 
846bd481603SBarry Smith    Input Parameters:
847bd481603SBarry Smith +  nrows - the number of rows indicated
8481d73ed98SBarry Smith .  rows - the indices of the rows
849bd481603SBarry Smith .  ncols - the number of columns in the matrix
850bd481603SBarry Smith .  cols - the columns indicated
851bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
852bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
853bd481603SBarry Smith 
854bd481603SBarry Smith 
855bd481603SBarry Smith    Level: intermediate
856bd481603SBarry Smith 
857bd481603SBarry Smith    Notes:
858bd481603SBarry Smith    See the chapter in the users manual on performance for more details
859bd481603SBarry Smith 
860bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
861bd481603SBarry Smith 
8621620fd73SBarry Smith    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
8631620fd73SBarry Smith 
864bd481603SBarry Smith   Concepts: preallocation^Matrix
865bd481603SBarry Smith 
866bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
867bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
868bd481603SBarry Smith M*/
869d3d32019SBarry Smith #define MatPreallocateSymmetricSet(row,nc,cols,dnz,onz) 0;\
870c1ac3661SBarry Smith { PetscInt __i; \
871d3d32019SBarry Smith   for (__i=0; __i<nc; __i++) {\
872d3d32019SBarry Smith     if (cols[__i] >= __end) onz[row - __rstart]++; \
873d3d32019SBarry Smith     else if (cols[__i] >= row) dnz[row - __rstart]++;\
874d3d32019SBarry Smith   }\
875d3d32019SBarry Smith }
876d3d32019SBarry Smith 
877bd481603SBarry Smith /*MC
87816371a99SBarry Smith    MatPreallocateLocation -  An alternative to MatPreallocationSet() that puts the nonzero locations into the matrix if it exists
87916371a99SBarry Smith 
88016371a99SBarry Smith    Synopsis:
88116371a99SBarry Smith    PetscErrorCode MatPreallocateLocations(Mat A,PetscInt row,PetscInt ncols,PetscInt *cols,PetscInt *dnz,PetscInt *onz)
88216371a99SBarry Smith 
88316371a99SBarry Smith    Not Collective
88416371a99SBarry Smith 
88516371a99SBarry Smith    Input Parameters:
88616371a99SBarry Smith .  A - matrix
88716371a99SBarry Smith .  row - row where values exist (must be local to this process)
88816371a99SBarry Smith .  ncols - number of columns
88916371a99SBarry Smith .  cols - columns with nonzeros
89016371a99SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
89116371a99SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
89216371a99SBarry Smith 
89316371a99SBarry Smith 
89416371a99SBarry Smith    Level: intermediate
89516371a99SBarry Smith 
89616371a99SBarry Smith    Notes:
89716371a99SBarry Smith    See the chapter in the users manual on performance for more details
89816371a99SBarry Smith 
89916371a99SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
90016371a99SBarry Smith 
90116371a99SBarry Smith    This is a MACRO not a function because it uses a bunch of variables private to the MatPreallocation.... routines.
90216371a99SBarry Smith 
90316371a99SBarry Smith   Concepts: preallocation^Matrix
90416371a99SBarry Smith 
90516371a99SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
90616371a99SBarry Smith           MatPreallocateSymmetricInitialize(), MatPreallocateSymmetricSetLocal()
90716371a99SBarry Smith M*/
90816371a99SBarry 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);}
90916371a99SBarry Smith 
91016371a99SBarry Smith 
91116371a99SBarry Smith /*MC
9120d2a7ff2SSatish Balay    MatPreallocateFinalize - Ends the block of code that will count the number of nonzeros per
913bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
914bd481603SBarry Smith 
915bd481603SBarry Smith    Synopsis:
916c1ac3661SBarry Smith    PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz)
917bd481603SBarry Smith 
918bd481603SBarry Smith    Collective on MPI_Comm
919bd481603SBarry Smith 
920bd481603SBarry Smith    Input Parameters:
92116371a99SBarry Smith +  dnz - the array that was be passed to the matrix preallocation routines
922bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
923bd481603SBarry Smith 
924bd481603SBarry Smith 
925bd481603SBarry Smith    Level: intermediate
926bd481603SBarry Smith 
927bd481603SBarry Smith    Notes:
928bd481603SBarry Smith    See the chapter in the users manual on performance for more details
929bd481603SBarry Smith 
930bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
931bd481603SBarry Smith 
9321620fd73SBarry Smith    This is a MACRO not a function because it closes the { started in MatPreallocateInitialize().
9331620fd73SBarry Smith 
934bd481603SBarry Smith   Concepts: preallocation^Matrix
935bd481603SBarry Smith 
936bd481603SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
937bd481603SBarry Smith           MatPreallocateSymmetricInitialize(), MatPreallocateSymmetricSetLocal()
938bd481603SBarry Smith M*/
939a89bc211SBarry Smith #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree2(dnz,onz);CHKERRQ(_4_ierr);}
9407c922b88SBarry Smith 
941bd481603SBarry Smith 
942bd481603SBarry Smith 
9437b80b807SBarry Smith /* Routines unique to particular data structures */
944be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellGetContext(Mat,void **);
945ba966639SSatish Balay PetscPolymorphicFunction(MatShellGetContext,(Mat A),(A,&t),void*,t)
946698d4c6aSKris Buschelman 
947be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInodeAdjustForInodes(Mat,IS*,IS*);
948be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *);
949698d4c6aSKris Buschelman 
950be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetColumnIndices(Mat,PetscInt[]);
951be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetColumnIndices(Mat,PetscInt[]);
952be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
953c75a6043SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
954c75a6043SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
9557b80b807SBarry Smith 
956a96a251dSBarry Smith #define MAT_SKIP_ALLOCATION -4
957a96a251dSBarry Smith 
958be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
959ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz))
960be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
961ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqSBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz))
962be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]);
963ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqAIJSetPreallocation,(Mat A,const PetscInt nnz[]),(A,0,nnz))
964273d9f13SBarry Smith 
965be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
966ba966639SSatish Balay PetscPolymorphicSubroutine(MatMPIBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[],const PetscInt onz[]),(A,bs,0,nnz,0,onz))
967be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
968be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
96953803f26SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []);
970725b52f3SLisandro Dalcin EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
971be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]);
972be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
973be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]);
974be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIDenseSetPreallocation(Mat,PetscScalar[]);
975284134d9SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation(Mat,PetscScalar[]);
976be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIRowbsSetPreallocation(Mat,PetscInt,const PetscInt[]);
977be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,PetscInt*[]);
978be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,PetscInt*[]);
979be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAdicSetLocalFunction(Mat,void (*)(void));
980273d9f13SBarry Smith 
981be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetLDA(Mat,PetscInt);
982ab92ecdeSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDenseGetLocalMatrix(Mat,Mat*);
9831b807ce4Svictorle 
984be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStoreValues(Mat);
985be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRetrieveValues(Mat);
9862e8a6d31SBarry Smith 
987be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDAADSetCtx(Mat,void*);
9883a7fca6bSBarry Smith 
9897b80b807SBarry Smith /*
9907b80b807SBarry Smith   These routines are not usually accessed directly, rather solving is
99194b7f48cSBarry Smith   done through the KSP and PC interfaces.
9927b80b807SBarry Smith */
9937b80b807SBarry Smith 
994d9274352SBarry Smith /*E
995d9274352SBarry Smith     MatOrderingType - String with the name of a PETSc matrix ordering or the creation function
996d9274352SBarry Smith        with an optional dynamic library name, for example
997d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:orderingcreate()
998d9274352SBarry Smith 
999d9274352SBarry Smith    Level: beginner
1000d9274352SBarry Smith 
1001e5a9bf91SBarry Smith    Cannot use const because the PC objects manipulate the string
1002e5a9bf91SBarry Smith 
1003d9274352SBarry Smith .seealso: MatGetOrdering()
1004d9274352SBarry Smith E*/
10053eaad2caSSatish Balay #define MatOrderingType char*
1006b12f92e5SBarry Smith #define MATORDERING_NATURAL     "natural"
1007b12f92e5SBarry Smith #define MATORDERING_ND          "nd"
1008b12f92e5SBarry Smith #define MATORDERING_1WD         "1wd"
1009b12f92e5SBarry Smith #define MATORDERING_RCM         "rcm"
1010b12f92e5SBarry Smith #define MATORDERING_QMD         "qmd"
1011b12f92e5SBarry Smith #define MATORDERING_ROWLENGTH   "rowlength"
1012e106eecfSBarry Smith #define MATORDERING_DSC_ND      "dsc_nd"         /* these three are only for DSCPACK, see its documentation for details */
101362152c8bSBarry Smith #define MATORDERING_DSC_MMD     "dsc_mmd"
101462152c8bSBarry Smith #define MATORDERING_DSC_MDF     "dsc_mdf"
1015c06d978dSMatthew Knepley #define MATORDERING_CONSTRAINED "constrained"
1016c06d978dSMatthew Knepley #define MATORDERING_IDENTITY    "identity"
1017c06d978dSMatthew Knepley #define MATORDERING_REVERSE     "reverse"
1018b12f92e5SBarry Smith 
1019f2416520SSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOrdering(Mat,const MatOrderingType,IS*,IS*);
1020be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOrderingList(PetscFList *list);
1021f2416520SSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,const MatOrderingType,IS*,IS*));
102230de9b25SBarry Smith 
102330de9b25SBarry Smith /*MC
102430de9b25SBarry Smith    MatOrderingRegisterDynamic - Adds a new sparse matrix ordering to the
102530de9b25SBarry Smith                                matrix package.
102630de9b25SBarry Smith 
102730de9b25SBarry Smith    Synopsis:
1028c1ac3661SBarry Smith    PetscErrorCode MatOrderingRegisterDynamic(char *name_ordering,char *path,char *name_create,PetscErrorCode (*routine_create)(MatOrdering))
102930de9b25SBarry Smith 
103030de9b25SBarry Smith    Not Collective
103130de9b25SBarry Smith 
103230de9b25SBarry Smith    Input Parameters:
103330de9b25SBarry Smith +  sname - name of ordering (for example MATORDERING_ND)
103430de9b25SBarry Smith .  path - location of library where creation routine is
103530de9b25SBarry Smith .  name - name of function that creates the ordering type,a string
103630de9b25SBarry Smith -  function - function pointer that creates the ordering
103730de9b25SBarry Smith 
103830de9b25SBarry Smith    Level: developer
103930de9b25SBarry Smith 
104030de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
104130de9b25SBarry Smith    is ignored.
104230de9b25SBarry Smith 
104330de9b25SBarry Smith    Sample usage:
104430de9b25SBarry Smith .vb
104530de9b25SBarry Smith    MatOrderingRegisterDynamic("my_order",/home/username/my_lib/lib/libO/solaris/mylib.a,
104630de9b25SBarry Smith                "MyOrder",MyOrder);
104730de9b25SBarry Smith .ve
104830de9b25SBarry Smith 
104930de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
105030de9b25SBarry Smith $     MatOrderingSetType(part,"my_order)
105130de9b25SBarry Smith    or at runtime via the option
10522401956bSBarry Smith $     -pc_factor_mat_ordering_type my_order
105330de9b25SBarry Smith 
1054ab901514SBarry Smith    ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values.
105530de9b25SBarry Smith 
105630de9b25SBarry Smith .keywords: matrix, ordering, register
105730de9b25SBarry Smith 
105830de9b25SBarry Smith .seealso: MatOrderingRegisterDestroy(), MatOrderingRegisterAll()
105930de9b25SBarry Smith M*/
1060aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1061f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,0)
1062b12f92e5SBarry Smith #else
1063f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,d)
1064b12f92e5SBarry Smith #endif
106530de9b25SBarry Smith 
1066be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegisterDestroy(void);
1067be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegisterAll(const char[]);
10682bad1931SBarry Smith extern PetscTruth MatOrderingRegisterAllCalled;
1069b0a32e0cSBarry Smith extern PetscFList MatOrderingList;
1070d4fbbf0eSBarry Smith 
1071be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS);
1072a2ce50c7SBarry Smith 
1073d91e6319SBarry Smith /*S
10742401956bSBarry Smith    MatFactorInfo - Data passed into the matrix factorization routines
1075b00f7748SHong Zhang 
107661cad860SBarry Smith    In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE, that is use
107761cad860SBarry Smith $     MatFactorInfo  info(MAT_FACTORINFO_SIZE)
1078b00f7748SHong Zhang 
107915e8a5b3SHong Zhang    Notes: These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC.
1080b00f7748SHong Zhang 
108161cad860SBarry Smith       You can use MatFactorInfoInitialize() to set default values.
108261cad860SBarry Smith 
1083b00f7748SHong Zhang    Level: developer
1084b00f7748SHong Zhang 
1085d7d82daaSBarry Smith .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(),
1086d7d82daaSBarry Smith           MatFactorInfoInitialize()
1087b00f7748SHong Zhang 
1088b00f7748SHong Zhang S*/
1089b00f7748SHong Zhang typedef struct {
10900a29876aSHong Zhang   PetscReal     shiftnz;        /* scaling of identity added to matrix to prevent zero pivots */
1091fbf22428SSatish Balay   PetscReal     shiftpd;        /* if true, shift until positive pivots */
109215e8a5b3SHong Zhang   PetscReal     diagonal_fill;  /* force diagonal to fill in if initially not filled */
109385317021SBarry Smith   PetscReal     usedt;
109415e8a5b3SHong Zhang   PetscReal     dt;             /* drop tolerance */
1095b00f7748SHong Zhang   PetscReal     dtcol;          /* tolerance for pivoting */
109615e8a5b3SHong Zhang   PetscReal     dtcount;        /* maximum nonzeros to be allowed per row */
109767e28bfeSBarry Smith   PetscReal     fill;           /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */
1098348344bbSBarry Smith   PetscReal     levels;         /* ICC/ILU(levels) */
1099bcd9e38bSBarry Smith   PetscReal     pivotinblocks;  /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0
1100bcd9e38bSBarry Smith                                    factorization may be faster if do not pivot */
110162bba022SBarry Smith   PetscReal     shiftinblocks;  /* if block in block factorization has zero pivot then shift diagonal until non-singular */
110215e8a5b3SHong Zhang   PetscReal     zeropivot;      /* pivot is called zero if less than this */
110315e8a5b3SHong Zhang } MatFactorInfo;
1104ffa6d0a5SLois Curfman McInnes 
1105be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFactorInfoInitialize(MatFactorInfo*);
11060481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactor(Mat,IS,const MatFactorInfo*);
11070481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
11080481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactorNumeric(Mat,Mat,const MatFactorInfo*);
11090481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactor(Mat,IS,IS,const MatFactorInfo*);
11100481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUFactor(Mat,IS,IS,const MatFactorInfo*);
11110481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
11120481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
11130481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatICCFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
11140481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatICCFactor(Mat,IS,const MatFactorInfo*);
11150481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactorNumeric(Mat,Mat,const MatFactorInfo*);
11160481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUDTFactor(Mat,IS,IS,const MatFactorInfo*,Mat *);
11173c2a7987SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUDTFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
11183c2a7987SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUDTFactorNumeric(Mat,Mat,const MatFactorInfo*);
1119be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*);
1120be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolve(Mat,Vec,Vec);
1121be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatForwardSolve(Mat,Vec,Vec);
1122be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatBackwardSolve(Mat,Vec,Vec);
1123be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveAdd(Mat,Vec,Vec,Vec);
1124be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveTranspose(Mat,Vec,Vec);
1125be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveTransposeAdd(Mat,Vec,Vec,Vec);
1126be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolves(Mat,Vecs,Vecs);
11278ed539a5SBarry Smith 
1128be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUnfactored(Mat);
1129bb5a7306SBarry Smith 
1130d91e6319SBarry Smith /*E
1131d91e6319SBarry Smith     MatSORType - What type of (S)SOR to perform
1132bb1eb677SSatish Balay 
1133d91e6319SBarry Smith     Level: beginner
1134d91e6319SBarry Smith 
1135d9274352SBarry Smith    May be bitwise ORd together
1136d9274352SBarry Smith 
1137d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
1138d91e6319SBarry Smith 
11394e7234bfSBarry Smith    MatSORType may be bitwise ORd together, so do not change the numbers
11404e7234bfSBarry Smith 
1141a1d92eedSBarry Smith .seealso: MatRelax(), MatPBRelax()
1142d91e6319SBarry Smith E*/
1143ee50ffe9SBarry Smith typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3,
1144ee50ffe9SBarry Smith               SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8,
1145ee50ffe9SBarry Smith               SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16,
114684cb2905SBarry Smith               SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType;
1147be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRelax(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
1148be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPBRelax(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
11498ed539a5SBarry Smith 
1150d4fbbf0eSBarry Smith /*
1151639f9d9dSBarry Smith     These routines are for efficiently computing Jacobians via finite differences.
1152639f9d9dSBarry Smith */
1153b12f92e5SBarry Smith 
1154d9274352SBarry Smith /*E
1155d9274352SBarry Smith     MatColoringType - String with the name of a PETSc matrix coloring or the creation function
1156d9274352SBarry Smith        with an optional dynamic library name, for example
1157d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:coloringcreate()
1158d9274352SBarry Smith 
1159d9274352SBarry Smith    Level: beginner
1160d9274352SBarry Smith 
1161d9274352SBarry Smith .seealso: MatGetColoring()
1162d9274352SBarry Smith E*/
1163a313700dSBarry Smith #define MatColoringType char*
1164b12f92e5SBarry Smith #define MATCOLORING_NATURAL "natural"
1165b12f92e5SBarry Smith #define MATCOLORING_SL      "sl"
1166b12f92e5SBarry Smith #define MATCOLORING_LF      "lf"
1167b12f92e5SBarry Smith #define MATCOLORING_ID      "id"
1168b12f92e5SBarry Smith 
1169ab970d9bSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColoring(Mat,const MatColoringType,ISColoring*);
11702e90c967SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,MatColoringType,ISColoring *));
117130de9b25SBarry Smith 
117230de9b25SBarry Smith /*MC
117330de9b25SBarry Smith    MatColoringRegisterDynamic - Adds a new sparse matrix coloring to the
117430de9b25SBarry Smith                                matrix package.
117530de9b25SBarry Smith 
117630de9b25SBarry Smith    Synopsis:
1177c1ac3661SBarry Smith    PetscErrorCode MatColoringRegisterDynamic(char *name_coloring,char *path,char *name_create,PetscErrorCode (*routine_create)(MatColoring))
117830de9b25SBarry Smith 
117930de9b25SBarry Smith    Not Collective
118030de9b25SBarry Smith 
118130de9b25SBarry Smith    Input Parameters:
118230de9b25SBarry Smith +  sname - name of Coloring (for example MATCOLORING_SL)
118330de9b25SBarry Smith .  path - location of library where creation routine is
118430de9b25SBarry Smith .  name - name of function that creates the Coloring type, a string
118530de9b25SBarry Smith -  function - function pointer that creates the coloring
118630de9b25SBarry Smith 
118730de9b25SBarry Smith    Level: developer
118830de9b25SBarry Smith 
118930de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
119030de9b25SBarry Smith    is ignored.
119130de9b25SBarry Smith 
119230de9b25SBarry Smith    Sample usage:
119330de9b25SBarry Smith .vb
119430de9b25SBarry Smith    MatColoringRegisterDynamic("my_color",/home/username/my_lib/lib/libO/solaris/mylib.a,
119530de9b25SBarry Smith                "MyColor",MyColor);
119630de9b25SBarry Smith .ve
119730de9b25SBarry Smith 
119830de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
119930de9b25SBarry Smith $     MatColoringSetType(part,"my_color")
120030de9b25SBarry Smith    or at runtime via the option
120130de9b25SBarry Smith $     -mat_coloring_type my_color
120230de9b25SBarry Smith 
1203ab901514SBarry Smith    $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
120430de9b25SBarry Smith 
120530de9b25SBarry Smith .keywords: matrix, Coloring, register
120630de9b25SBarry Smith 
120730de9b25SBarry Smith .seealso: MatColoringRegisterDestroy(), MatColoringRegisterAll()
120830de9b25SBarry Smith M*/
1209aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1210f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,0)
1211b12f92e5SBarry Smith #else
1212f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,d)
1213b12f92e5SBarry Smith #endif
121430de9b25SBarry Smith 
12152bad1931SBarry Smith extern PetscTruth MatColoringRegisterAllCalled;
1216f1144a30SSatish Balay 
1217be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegisterAll(const char[]);
1218be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegisterDestroy(void);
1219be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*);
1220639f9d9dSBarry Smith 
1221d9274352SBarry Smith /*S
1222d9274352SBarry Smith      MatFDColoring - Object for computing a sparse Jacobian via finite differences
1223d9274352SBarry Smith         and coloring
1224639f9d9dSBarry Smith 
1225d9274352SBarry Smith    Level: beginner
1226d9274352SBarry Smith 
1227d9274352SBarry Smith   Concepts: coloring, sparse Jacobian, finite differences
1228d9274352SBarry Smith 
1229d9274352SBarry Smith .seealso:  MatFDColoringCreate()
1230d9274352SBarry Smith S*/
1231e2a1c21fSSatish Balay typedef struct _p_MatFDColoring* MatFDColoring;
1232639f9d9dSBarry Smith 
1233be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringCreate(Mat,ISColoring,MatFDColoring *);
1234be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringDestroy(MatFDColoring);
1235be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringView(MatFDColoring,PetscViewer);
1236be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*);
12374a9d489dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**);
1238be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal);
1239be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFromOptions(MatFDColoring);
1240be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApply(Mat,MatFDColoring,Vec,MatStructure*,void *);
1241be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApplyTS(Mat,MatFDColoring,PetscReal,Vec,MatStructure*,void *);
1242be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetF(MatFDColoring,Vec);
1243be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,PetscInt*[]);
1244639f9d9dSBarry Smith /*
12450752156aSBarry Smith     These routines are for partitioning matrices: currently used only
12463eda8832SBarry Smith   for adjacency matrix, MatCreateMPIAdj().
12470752156aSBarry Smith */
1248ca161407SBarry Smith 
1249d9274352SBarry Smith /*S
1250d9274352SBarry Smith      MatPartitioning - Object for managing the partitioning of a matrix or graph
1251d9274352SBarry Smith 
1252d9274352SBarry Smith    Level: beginner
1253d9274352SBarry Smith 
1254d9274352SBarry Smith   Concepts: partitioning
1255d9274352SBarry Smith 
1256743a1026Svictorle .seealso:  MatPartitioningCreate(), MatPartitioningType
1257d9274352SBarry Smith S*/
125891e9ee9fSBarry Smith typedef struct _p_MatPartitioning* MatPartitioning;
1259d9274352SBarry Smith 
1260d9274352SBarry Smith /*E
12615bc6e7ccSvictorle     MatPartitioningType - String with the name of a PETSc matrix partitioning or the creation function
1262d9274352SBarry Smith        with an optional dynamic library name, for example
1263d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:partitioningcreate()
1264d9274352SBarry Smith 
1265d9274352SBarry Smith    Level: beginner
1266d9274352SBarry Smith 
1267b547a13bSvictorle .seealso: MatPartitioningCreate(), MatPartitioning
1268d9274352SBarry Smith E*/
1269a313700dSBarry Smith #define MatPartitioningType char*
12708ba1e511SMatthew Knepley #define MAT_PARTITIONING_CURRENT  "current"
1271c2ec2a73SHong Zhang #define MAT_PARTITIONING_SQUARE   "square"
12728ba1e511SMatthew Knepley #define MAT_PARTITIONING_PARMETIS "parmetis"
127371306c60Sjroman #define MAT_PARTITIONING_CHACO    "chaco"
127471306c60Sjroman #define MAT_PARTITIONING_JOSTLE   "jostle"
127571306c60Sjroman #define MAT_PARTITIONING_PARTY    "party"
127671306c60Sjroman #define MAT_PARTITIONING_SCOTCH   "scotch"
127771306c60Sjroman 
1278ca161407SBarry Smith 
1279be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningCreate(MPI_Comm,MatPartitioning*);
1280a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetType(MatPartitioning,const MatPartitioningType);
1281be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetNParts(MatPartitioning,PetscInt);
1282be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetAdjacency(MatPartitioning,Mat);
1283be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]);
1284be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []);
1285be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningApply(MatPartitioning,IS*);
1286be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningDestroy(MatPartitioning);
12872aabb6bbSBarry Smith 
1288be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatPartitioning));
128930de9b25SBarry Smith 
129030de9b25SBarry Smith /*MC
129130de9b25SBarry Smith    MatPartitioningRegisterDynamic - Adds a new sparse matrix partitioning to the
129230de9b25SBarry Smith    matrix package.
129330de9b25SBarry Smith 
129430de9b25SBarry Smith    Synopsis:
1295c1ac3661SBarry Smith    PetscErrorCode MatPartitioningRegisterDynamic(char *name_partitioning,char *path,char *name_create,PetscErrorCode (*routine_create)(MatPartitioning))
129630de9b25SBarry Smith 
129730de9b25SBarry Smith    Not Collective
129830de9b25SBarry Smith 
129930de9b25SBarry Smith    Input Parameters:
130030de9b25SBarry Smith +  sname - name of partitioning (for example MAT_PARTITIONING_CURRENT) or parmetis
130130de9b25SBarry Smith .  path - location of library where creation routine is
130230de9b25SBarry Smith .  name - name of function that creates the partitioning type, a string
130330de9b25SBarry Smith -  function - function pointer that creates the partitioning type
130430de9b25SBarry Smith 
130530de9b25SBarry Smith    Level: developer
130630de9b25SBarry Smith 
130730de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
130830de9b25SBarry Smith    is ignored.
130930de9b25SBarry Smith 
131030de9b25SBarry Smith    Sample usage:
131130de9b25SBarry Smith .vb
131230de9b25SBarry Smith    MatPartitioningRegisterDynamic("my_part",/home/username/my_lib/lib/libO/solaris/mylib.a,
131330de9b25SBarry Smith                "MyPartCreate",MyPartCreate);
131430de9b25SBarry Smith .ve
131530de9b25SBarry Smith 
131630de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
131730de9b25SBarry Smith $     MatPartitioningSetType(part,"my_part")
131830de9b25SBarry Smith    or at runtime via the option
131930de9b25SBarry Smith $     -mat_partitioning_type my_part
132030de9b25SBarry Smith 
1321ab901514SBarry Smith    $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
132230de9b25SBarry Smith 
132330de9b25SBarry Smith .keywords: matrix, partitioning, register
132430de9b25SBarry Smith 
132530de9b25SBarry Smith .seealso: MatPartitioningRegisterDestroy(), MatPartitioningRegisterAll()
132630de9b25SBarry Smith M*/
1327aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1328f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,0)
13292aabb6bbSBarry Smith #else
1330f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,d)
13312aabb6bbSBarry Smith #endif
13322aabb6bbSBarry Smith 
13332bad1931SBarry Smith extern PetscTruth MatPartitioningRegisterAllCalled;
1334f1144a30SSatish Balay 
1335be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegisterAll(const char[]);
1336be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegisterDestroy(void);
13372bad1931SBarry Smith 
1338be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningView(MatPartitioning,PetscViewer);
1339be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetFromOptions(MatPartitioning);
1340a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningGetType(MatPartitioning,const MatPartitioningType*);
1341ca161407SBarry Smith 
1342be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningParmetisSetCoarseSequential(MatPartitioning);
13430e8af3caSMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningParmetisGetEdgeCut(MatPartitioning, PetscInt *);
13440752156aSBarry Smith 
1345be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningJostleSetCoarseLevel(MatPartitioning,PetscReal);
1346be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningJostleSetCoarseSequential(MatPartitioning);
134771306c60Sjroman 
1348954b3ec6SBarry Smith typedef enum {MP_CHACO_MULTILEVEL_KL,MP_CHACO_SPECTRAL,MP_CHACO_LINEAR,MP_CHACO_RANDOM, MP_CHACO_SCATTERED} MPChacoGlobalType;
1349be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetGlobal(MatPartitioning, MPChacoGlobalType);
135071306c60Sjroman typedef enum { MP_CHACO_KERNIGHAN_LIN, MP_CHACO_NONE } MPChacoLocalType;
1351be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetLocal(MatPartitioning, MPChacoLocalType);
1352be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal);
135371306c60Sjroman typedef enum { MP_CHACO_LANCZOS, MP_CHACO_RQI_SYMMLQ } MPChacoEigenType;
1354be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType);
1355be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenTol(MatPartitioning, PetscReal);
1356be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenNumber(MatPartitioning, PetscInt);
135771306c60Sjroman 
135871306c60Sjroman #define MP_PARTY_OPT "opt"
135971306c60Sjroman #define MP_PARTY_LIN "lin"
136071306c60Sjroman #define MP_PARTY_SCA "sca"
136171306c60Sjroman #define MP_PARTY_RAN "ran"
136271306c60Sjroman #define MP_PARTY_GBF "gbf"
136371306c60Sjroman #define MP_PARTY_GCF "gcf"
136471306c60Sjroman #define MP_PARTY_BUB "bub"
136571306c60Sjroman #define MP_PARTY_DEF "def"
1366be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetGlobal(MatPartitioning, const char*);
136771306c60Sjroman #define MP_PARTY_HELPFUL_SETS "hs"
136871306c60Sjroman #define MP_PARTY_KERNIGHAN_LIN "kl"
136971306c60Sjroman #define MP_PARTY_NONE "no"
1370be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetLocal(MatPartitioning, const char*);
1371be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal);
1372be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetBipart(MatPartitioning,PetscTruth);
1373be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscTruth);
137471306c60Sjroman 
137571306c60Sjroman typedef enum { MP_SCOTCH_GREEDY, MP_SCOTCH_GPS, MP_SCOTCH_GR_GPS } MPScotchGlobalType;
1376be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetArch(MatPartitioning,const char*);
1377be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetMultilevel(MatPartitioning);
1378be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetGlobal(MatPartitioning,MPScotchGlobalType);
1379be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetCoarseLevel(MatPartitioning,PetscReal);
1380be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetHostList(MatPartitioning,const char*);
138171306c60Sjroman typedef enum { MP_SCOTCH_KERNIGHAN_LIN, MP_SCOTCH_NONE } MPScotchLocalType;
1382be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetLocal(MatPartitioning,MPScotchLocalType);
1383be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetMapping(MatPartitioning);
1384be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetStrategy(MatPartitioning,char*);
138571306c60Sjroman 
13860752156aSBarry Smith /*
13870a835dfdSSatish Balay     If you add entries here you must also add them to finclude/petscmat.h
1388d4fbbf0eSBarry Smith */
13891c1c02c0SLois Curfman McInnes typedef enum { MATOP_SET_VALUES=0,
13901c1c02c0SLois Curfman McInnes                MATOP_GET_ROW=1,
13911c1c02c0SLois Curfman McInnes                MATOP_RESTORE_ROW=2,
13921c1c02c0SLois Curfman McInnes                MATOP_MULT=3,
13931c1c02c0SLois Curfman McInnes                MATOP_MULT_ADD=4,
13947c922b88SBarry Smith                MATOP_MULT_TRANSPOSE=5,
13957c922b88SBarry Smith                MATOP_MULT_TRANSPOSE_ADD=6,
13961c1c02c0SLois Curfman McInnes                MATOP_SOLVE=7,
13971c1c02c0SLois Curfman McInnes                MATOP_SOLVE_ADD=8,
13987c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE=9,
13997c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE_ADD=10,
14001c1c02c0SLois Curfman McInnes                MATOP_LUFACTOR=11,
14011c1c02c0SLois Curfman McInnes                MATOP_CHOLESKYFACTOR=12,
14021c1c02c0SLois Curfman McInnes                MATOP_RELAX=13,
14031c1c02c0SLois Curfman McInnes                MATOP_TRANSPOSE=14,
14041c1c02c0SLois Curfman McInnes                MATOP_GETINFO=15,
14051c1c02c0SLois Curfman McInnes                MATOP_EQUAL=16,
14061c1c02c0SLois Curfman McInnes                MATOP_GET_DIAGONAL=17,
14071c1c02c0SLois Curfman McInnes                MATOP_DIAGONAL_SCALE=18,
14081c1c02c0SLois Curfman McInnes                MATOP_NORM=19,
14091c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_BEGIN=20,
14101c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_END=21,
1411d519adbfSMatthew Knepley                MATOP_SET_OPTION=22,
1412d519adbfSMatthew Knepley                MATOP_ZERO_ENTRIES=23,
1413d519adbfSMatthew Knepley                MATOP_ZERO_ROWS=24,
1414d519adbfSMatthew Knepley                MATOP_LUFACTOR_SYMBOLIC=25,
1415d519adbfSMatthew Knepley                MATOP_LUFACTOR_NUMERIC=26,
1416d519adbfSMatthew Knepley                MATOP_CHOLESKY_FACTOR_SYMBOLIC=27,
1417d519adbfSMatthew Knepley                MATOP_CHOLESKY_FACTOR_NUMERIC=28,
1418d519adbfSMatthew Knepley                MATOP_SETUP_PREALLOCATION=29,
1419d519adbfSMatthew Knepley                MATOP_ILUFACTOR_SYMBOLIC=30,
1420d519adbfSMatthew Knepley                MATOP_ICCFACTOR_SYMBOLIC=31,
1421d519adbfSMatthew Knepley                MATOP_GET_ARRAY=32,
1422d519adbfSMatthew Knepley                MATOP_RESTORE_ARRAY=33,
1423d519adbfSMatthew Knepley                MATOP_DUPLICATE=34,
1424d519adbfSMatthew Knepley                MATOP_FORWARD_SOLVE=35,
1425d519adbfSMatthew Knepley                MATOP_BACKWARD_SOLVE=36,
1426d519adbfSMatthew Knepley                MATOP_ILUFACTOR=37,
1427d519adbfSMatthew Knepley                MATOP_ICCFACTOR=38,
1428d519adbfSMatthew Knepley                MATOP_AXPY=39,
1429d519adbfSMatthew Knepley                MATOP_GET_SUBMATRICES=40,
1430d519adbfSMatthew Knepley                MATOP_INCREASE_OVERLAP=41,
1431d519adbfSMatthew Knepley                MATOP_GET_VALUES=42,
1432d519adbfSMatthew Knepley                MATOP_COPY=43,
1433d519adbfSMatthew Knepley                MATOP_GET_ROW_MAX=44,
1434d519adbfSMatthew Knepley                MATOP_SCALE=45,
1435d519adbfSMatthew Knepley                MATOP_SHIFT=46,
1436d519adbfSMatthew Knepley                MATOP_DIAGONAL_SHIFT=47,
1437d519adbfSMatthew Knepley                MATOP_ILUDT_FACTOR=48,
1438d519adbfSMatthew Knepley                MATOP_SET_BLOCK_SIZE=49,
1439d519adbfSMatthew Knepley                MATOP_GET_ROW_IJ=50,
1440d519adbfSMatthew Knepley                MATOP_RESTORE_ROW_IJ=51,
1441d519adbfSMatthew Knepley                MATOP_GET_COLUMN_IJ=52,
1442d519adbfSMatthew Knepley                MATOP_RESTORE_COLUMN_IJ=53,
1443d519adbfSMatthew Knepley                MATOP_FDCOLORING_CREATE=54,
1444d519adbfSMatthew Knepley                MATOP_COLORING_PATCH=55,
1445d519adbfSMatthew Knepley                MATOP_SET_UNFACTORED=56,
1446d519adbfSMatthew Knepley                MATOP_PERMUTE=57,
1447d519adbfSMatthew Knepley                MATOP_SET_VALUES_BLOCKED=58,
1448d519adbfSMatthew Knepley                MATOP_GET_SUBMATRIX=59,
1449d519adbfSMatthew Knepley                MATOP_DESTROY=60,
1450d519adbfSMatthew Knepley                MATOP_VIEW=61,
1451d519adbfSMatthew Knepley                MATOP_CONVERT_FROM=62,
1452d519adbfSMatthew Knepley                MATOP_USE_SCALED_FORM=63,
1453d519adbfSMatthew Knepley                MATOP_SCALE_SYSTEM=64,
1454d519adbfSMatthew Knepley                MATOP_UNSCALE_SYSTEM=65,
1455d519adbfSMatthew Knepley                MATOP_SET_LOCAL_TO_GLOBAL_MAP=66,
1456d519adbfSMatthew Knepley                MATOP_SET_VALUES_LOCAL=67,
1457d519adbfSMatthew Knepley                MATOP_ZERO_ROWS_LOCAL=68,
1458d519adbfSMatthew Knepley                MATOP_GET_ROW_MAX_ABS=69,
1459d519adbfSMatthew Knepley                MATOP_GET_ROW_MIN_ABS=70,
1460d519adbfSMatthew Knepley                MATOP_CONVERT=71,
1461d519adbfSMatthew Knepley                MATOP_SET_COLORING=72,
1462d519adbfSMatthew Knepley                MATOP_SET_VALUES_ADIC=73,
1463d519adbfSMatthew Knepley                MATOP_SET_VALUES_ADIFOR=74,
1464d519adbfSMatthew Knepley                MATOP_FD_COLORING_APPLY=75,
1465d519adbfSMatthew Knepley                MATOP_SET_FROM_OPTIONS=76,
1466d519adbfSMatthew Knepley                MATOP_MULT_CON=77,
1467d519adbfSMatthew Knepley                MATOP_MULT_TRANSPOSE_CON=78,
1468d519adbfSMatthew Knepley                MATOP_PERMUTE_SPARSIFY=79,
1469d519adbfSMatthew Knepley                MATOP_MULT_MULTIPLE=80,
1470d519adbfSMatthew Knepley                MATOP_SOLVE_MULTIPLE=81,
1471d519adbfSMatthew Knepley                MATOP_GET_INERTIA=82,
1472d519adbfSMatthew Knepley                MATOP_LOAD=83,
1473d519adbfSMatthew Knepley                MATOP_IS_SYMMETRIC=84,
1474d519adbfSMatthew Knepley                MATOP_IS_HERMITIAN=85,
1475d519adbfSMatthew Knepley                MATOP_IS_STRUCTURALLY_SYMMETRIC=86,
1476d519adbfSMatthew Knepley                MATOP_PB_RELAX=87,
1477d519adbfSMatthew Knepley                MATOP_GET_VECS=88,
1478d519adbfSMatthew Knepley                MATOP_MAT_MULT=89,
1479d519adbfSMatthew Knepley                MATOP_MAT_MULT_SYMBOLIC=90,
1480d519adbfSMatthew Knepley                MATOP_MAT_MULT_NUMERIC=91,
1481d519adbfSMatthew Knepley                MATOP_PTAP=92,
1482d519adbfSMatthew Knepley                MATOP_PTAP_SYMBOLIC=93,
1483d519adbfSMatthew Knepley                MATOP_PTAP_NUMERIC=94,
1484d519adbfSMatthew Knepley                MATOP_MAT_MULTTRANSPOSE=95,
1485d519adbfSMatthew Knepley                MATOP_MAT_MULTTRANSPOSE_SYMBOLIC=96,
1486d519adbfSMatthew Knepley                MATOP_MAT_MULTTRANSPOSE_NUMERIC=97,
1487d519adbfSMatthew Knepley                MATOP_PTAP_SYMBOLIC_SEQAIJ=98,
1488d519adbfSMatthew Knepley                MATOP_PTAP_NUMERIC_SEQAIJ=99,
1489d519adbfSMatthew Knepley                MATOP_PTAP_SYMBOLIC_MPIAIJ=100,
1490d519adbfSMatthew Knepley                MATOP_PTAP_NUMERIC_MPIAIJ=101,
1491d519adbfSMatthew Knepley                MATOP_CONJUGATE=102,
1492d519adbfSMatthew Knepley                MATOP_SET_SIZES=103,
1493d519adbfSMatthew Knepley                MATOP_SET_VALUES_ROW = 104,
1494d519adbfSMatthew Knepley                MATOP_REAL_PART=105,
1495d519adbfSMatthew Knepley                MATOP_IMAG_PART=106,
1496d519adbfSMatthew Knepley                MATOP_GET_ROW_UTRIANGULAR=107,
1497d519adbfSMatthew Knepley                MATOP_RESTORE_ROW_UTRIANGULAR=108,
1498d519adbfSMatthew Knepley                MATOP_MATSOLVE=109,
1499d519adbfSMatthew Knepley                MATOP_GET_REDUNDANTMATRIX=110,
1500d519adbfSMatthew Knepley                MATOP_GET_ROW_MIN=111,
1501d519adbfSMatthew Knepley                MATOP_GET_COLUMN_VEC=112,
1502d519adbfSMatthew Knepley                MATOP_MISSING_DIAGONAL=113,
1503d519adbfSMatthew Knepley                MATOP_MATGETSEQNONZEROSTRUCTURE=114,
1504d519adbfSMatthew Knepley                MATOP_DESTROY_SOLVER=115
1505fae171e0SBarry Smith              } MatOperation;
1506be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatHasOperation(Mat,MatOperation,PetscTruth*);
1507be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellSetOperation(Mat,MatOperation,void(*)(void));
1508be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellGetOperation(Mat,MatOperation,void(**)(void));
1509be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellSetContext(Mat,void*);
1510112a2221SBarry Smith 
151190ace30eSBarry Smith /*
151290ace30eSBarry Smith    Codes for matrices stored on disk. By default they are
151390ace30eSBarry Smith  stored in a universal format. By changing the format with
15147973ad04SBarry Smith  PetscViewerSetFormat(viewer,PETSC_VIEWER_NATIVE); the matrices will
151590ace30eSBarry Smith  be stored in a way natural for the matrix, for example dense matrices
151690ace30eSBarry Smith  would be stored as dense. Matrices stored this way may only be
151790ace30eSBarry Smith  read into matrices of the same time.
151890ace30eSBarry Smith */
151990ace30eSBarry Smith #define MATRIX_BINARY_FORMAT_DENSE -1
152090ace30eSBarry Smith 
1521be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetHashTableFactor(Mat,PetscReal);
1522be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatISGetLocalMat(Mat,Mat*);
15231f4e1ec7SBarry Smith 
1524d9274352SBarry Smith /*S
1525d9274352SBarry Smith      MatNullSpace - Object that removes a null space from a vector, i.e.
1526d9274352SBarry Smith          orthogonalizes the vector to a subsapce
1527d9274352SBarry Smith 
1528f7a9e4ceSBarry Smith    Level: advanced
1529d9274352SBarry Smith 
1530d9274352SBarry Smith   Concepts: matrix; linear operator, null space
1531d9274352SBarry Smith 
15326e1639daSBarry Smith   Users manual sections:
15336e1639daSBarry Smith .   sec_singular
15346e1639daSBarry Smith 
1535d9274352SBarry Smith .seealso:  MatNullSpaceCreate()
1536d9274352SBarry Smith S*/
153774637425SBarry Smith typedef struct _p_MatNullSpace* MatNullSpace;
1538d9274352SBarry Smith 
1539be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceCreate(MPI_Comm,PetscTruth,PetscInt,const Vec[],MatNullSpace*);
1540281d1b2eSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(Vec,void*),void*);
1541be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceDestroy(MatNullSpace);
1542be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceRemove(MatNullSpace,Vec,Vec*);
1543be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceAttach(Mat,MatNullSpace);
154495902228SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceTest(MatNullSpace,Mat,PetscTruth *);
154574637425SBarry Smith 
1546be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatReorderingSeqSBAIJ(Mat,IS);
1547be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPISBAIJSetHashTableFactor(Mat,PetscReal);
1548be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetColumnIndices(Mat,PetscInt *);
154946129b97SKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJInvertBlockDiagonal(Mat);
15503f1d51d7SBarry Smith 
1551be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMAIJ(Mat,PetscInt,Mat*);
1552be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMAIJRedimension(Mat,PetscInt,Mat*);
1553be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMAIJGetAIJ(Mat,Mat*);
1554c4f061fbSSatish Balay 
1555be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatComputeExplicitOperator(Mat,Mat*);
1556b0a32e0cSBarry Smith 
1557be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalScaleLocal(Mat,Vec);
155804f1ad80SBarry Smith 
1559fef1beadSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*);
15603ec795f1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetBase(Mat,Vec,Vec);
15613ec795f1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*);
1562e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*));
1563e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec));
1564e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDAddNullSpace(Mat,MatNullSpace);
1565e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt);
1566e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDResetHHistory(Mat);
1567e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctionError(Mat,PetscReal);
1568e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetPeriod(Mat,PetscInt);
1569e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDGetH(Mat,PetscScalar *);
1570e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFromOptions(Mat);
1571e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*);
1572e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*);
1573e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDDSSetUmin(Mat,PetscReal);
1574e884886fSBarry Smith 
1575e884886fSBarry Smith 
1576e884886fSBarry Smith typedef struct _p_MatMFFD* MatMFFD;
1577e884886fSBarry Smith 
1578e884886fSBarry Smith /*E
1579e884886fSBarry Smith     MatMFFDType - algorithm used to compute the h used in computing matrix-vector products via differencing of the function
1580e884886fSBarry Smith 
1581e884886fSBarry Smith    Level: beginner
1582e884886fSBarry Smith 
1583e884886fSBarry Smith .seealso: MatMFFDSetType(), MatMFFDRegister()
1584e884886fSBarry Smith E*/
1585a313700dSBarry Smith #define MatMFFDType char*
1586e884886fSBarry Smith #define MATMFFD_DS  "ds"
1587e884886fSBarry Smith #define MATMFFD_WP  "wp"
1588e884886fSBarry Smith 
1589a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetType(Mat,const MatMFFDType);
1590e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatMFFD));
1591e884886fSBarry Smith 
1592e884886fSBarry Smith /*MC
1593e884886fSBarry Smith    MatMFFDRegisterDynamic - Adds a method to the MatMFFD registry.
1594e884886fSBarry Smith 
1595e884886fSBarry Smith    Synopsis:
1596e884886fSBarry Smith    PetscErrorCode MatMFFDRegisterDynamic(char *name_solver,char *path,char *name_create,PetscErrorCode (*routine_create)(MatMFFD))
1597e884886fSBarry Smith 
1598e884886fSBarry Smith    Not Collective
1599e884886fSBarry Smith 
1600e884886fSBarry Smith    Input Parameters:
1601e884886fSBarry Smith +  name_solver - name of a new user-defined compute-h module
1602e884886fSBarry Smith .  path - path (either absolute or relative) the library containing this solver
1603e884886fSBarry Smith .  name_create - name of routine to create method context
1604e884886fSBarry Smith -  routine_create - routine to create method context
1605e884886fSBarry Smith 
1606e884886fSBarry Smith    Level: developer
1607e884886fSBarry Smith 
1608e884886fSBarry Smith    Notes:
1609e884886fSBarry Smith    MatMFFDRegisterDynamic() may be called multiple times to add several user-defined solvers.
1610e884886fSBarry Smith 
1611e884886fSBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
1612e884886fSBarry Smith    is ignored.
1613e884886fSBarry Smith 
1614e884886fSBarry Smith    Sample usage:
1615e884886fSBarry Smith .vb
1616e884886fSBarry Smith    MatMFFDRegisterDynamic("my_h",/home/username/my_lib/lib/libO/solaris/mylib.a,
1617e884886fSBarry Smith                "MyHCreate",MyHCreate);
1618e884886fSBarry Smith .ve
1619e884886fSBarry Smith 
1620e884886fSBarry Smith    Then, your solver can be chosen with the procedural interface via
1621e884886fSBarry Smith $     MatMFFDSetType(mfctx,"my_h")
1622e884886fSBarry Smith    or at runtime via the option
1623e884886fSBarry Smith $     -snes_mf_type my_h
1624e884886fSBarry Smith 
1625e884886fSBarry Smith .keywords: MatMFFD, register
1626e884886fSBarry Smith 
1627e884886fSBarry Smith .seealso: MatMFFDRegisterAll(), MatMFFDRegisterDestroy()
1628e884886fSBarry Smith M*/
1629e884886fSBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1630e884886fSBarry Smith #define MatMFFDRegisterDynamic(a,b,c,d) MatMFFDRegister(a,b,c,0)
1631e884886fSBarry Smith #else
1632e884886fSBarry Smith #define MatMFFDRegisterDynamic(a,b,c,d) MatMFFDRegister(a,b,c,d)
1633e884886fSBarry Smith #endif
1634e884886fSBarry Smith 
1635e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDRegisterAll(const char[]);
1636e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDRegisterDestroy(void);
1637e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDDefaultSetUmin(Mat,PetscReal);
1638e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDWPSetComputeNormU(Mat,PetscTruth);
1639e884886fSBarry Smith 
1640e884886fSBarry Smith 
1641be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *);
1642be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *);
16437dbadf16SMatthew Knepley 
1644e9fa29b7SSatish Balay PETSC_EXTERN_CXX_END
16452eac72dbSBarry Smith #endif
1646