xref: /petsc/include/petscmat.h (revision 72ca8751900d6af5fbc172fe35212f0f69bfacc6)
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"
6485e3dda7SBarry Smith #define MATTRANSPOSE       "transpose"
65*72ca8751SBarry Smith #define MATSCHURCOMPLEMENT "schurcomplement"
66d91e6319SBarry Smith 
679989ab13SBarry Smith /*E
68c7393fdbSBarry Smith     MatSolverPackage - String with the name of a PETSc matrix solver type.
699989ab13SBarry Smith 
709989ab13SBarry Smith     For example: "petsc" indicates what PETSc provides, "superlu" indicates either
719989ab13SBarry Smith        SuperLU or SuperLU_Dist etc.
729989ab13SBarry Smith 
739989ab13SBarry Smith 
749989ab13SBarry Smith    Level: beginner
759989ab13SBarry Smith 
765c9eb25fSBarry Smith .seealso: MatGetFactor(), Mat, MatSetType(), MatType
779989ab13SBarry Smith E*/
78c7393fdbSBarry Smith #define MatSolverPackage char*
79b24902e0SBarry Smith #define MAT_SOLVER_SPOOLES      "spooles"
80b24902e0SBarry Smith #define MAT_SOLVER_SUPERLU      "superlu"
81b24902e0SBarry Smith #define MAT_SOLVER_SUPERLU_DIST "superlu_dist"
82b24902e0SBarry Smith #define MAT_SOLVER_UMFPACK      "umfpack"
83b24902e0SBarry Smith #define MAT_SOLVER_ESSL         "essl"
84b24902e0SBarry Smith #define MAT_SOLVER_LUSOL        "lusol"
85b24902e0SBarry Smith #define MAT_SOLVER_MUMPS        "mumps"
86b24902e0SBarry Smith #define MAT_SOLVER_DSCPACK      "dscpack"
87b24902e0SBarry Smith #define MAT_SOLVER_MATLAB       "matlab"
8843244d56SBarry Smith #define MAT_SOLVER_PETSC        "petsc"
89b24902e0SBarry Smith 
90b24902e0SBarry Smith /*E
91b24902e0SBarry Smith     MatFactorType - indicates what type of factorization is requested
92b24902e0SBarry Smith 
93b24902e0SBarry Smith     Level: beginner
94b24902e0SBarry Smith 
95b24902e0SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
96b24902e0SBarry Smith 
97c7393fdbSBarry Smith .seealso: MatSolverPackage, MatGetFactor()
98b24902e0SBarry Smith E*/
995c9eb25fSBarry Smith typedef enum {MAT_FACTOR_NONE, MAT_FACTOR_LU, MAT_FACTOR_CHOLESKY, MAT_FACTOR_ILU, MAT_FACTOR_ICC} MatFactorType;
100c7393fdbSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetFactor(Mat,const MatSolverPackage,MatFactorType,Mat*);
101b24902e0SBarry Smith 
1029989ab13SBarry Smith 
103c06d978dSMatthew Knepley /* Logging support */
104552e946dSBarry Smith #define    MAT_FILE_COOKIE 1211216    /* used to indicate matrices in binary files */
105be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_COOKIE;
106be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_FDCOLORING_COOKIE;
107be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_PARTITIONING_COOKIE;
108be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_NULLSPACE_COOKIE;
109e884886fSBarry Smith extern PetscCookie PETSCMAT_DLLEXPORT MATMFFD_COOKIE;
110c06d978dSMatthew Knepley 
111ceb03754SKris Buschelman /*E
112ceb03754SKris Buschelman     MatReuse - Indicates if matrices obtained from a previous call to MatGetSubMatrices()
113ceb03754SKris Buschelman      or MatGetSubMatrix() are to be reused to store the new matrix values.
114ceb03754SKris Buschelman 
115ceb03754SKris Buschelman     Level: beginner
116ceb03754SKris Buschelman 
117ceb03754SKris Buschelman    Any additions/changes here MUST also be made in include/finclude/petscmat.h
118ceb03754SKris Buschelman 
1190c451bc4SBarry Smith .seealso: MatGetSubMatrices(), MatGetSubMatrix(), MatDestroyMatrices(), MatConvert()
120ceb03754SKris Buschelman E*/
121ceb03754SKris Buschelman typedef enum {MAT_INITIAL_MATRIX,MAT_REUSE_MATRIX} MatReuse;
122ceb03754SKris Buschelman 
1235494a064SHong Zhang /*E
1245494a064SHong Zhang     MatGetSubMatrixOption - Indicates if matrices obtained from a call to MatGetSubMatrices()
125829201f2SHong Zhang      include the matrix values. Currently it is only used by MatGetSeqNonzerostructure().
1265494a064SHong Zhang 
1275494a064SHong Zhang     Level: beginner
1285494a064SHong Zhang 
129829201f2SHong Zhang .seealso: MatGetSeqNonzerostructure()
1305494a064SHong Zhang E*/
1315494a064SHong Zhang typedef enum {MAT_DO_NOT_GET_VALUES,MAT_GET_VALUES} MatGetSubMatrixOption;
1325494a064SHong Zhang 
133e5bd5246SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInitializePackage(const char[]);
134c06d978dSMatthew Knepley 
135f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreate(MPI_Comm,Mat*);
136a69a960dSMatthew Knepley PetscPolymorphicFunction(MatCreate,(MPI_Comm comm),(comm,&A),Mat,A)
137a69a960dSMatthew Knepley PetscPolymorphicFunction(MatCreate,(),(PETSC_COMM_WORLD,&A),Mat,A)
138f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetSizes(Mat,PetscInt,PetscInt,PetscInt,PetscInt);
139a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetType(Mat,const MatType);
140be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetFromOptions(Mat);
141be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUpPreallocation(Mat);
142be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRegisterAll(const char[]);
143be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat));
14430de9b25SBarry Smith 
145f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetOptionsPrefix(Mat,const char[]);
146f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAppendOptionsPrefix(Mat,const char[]);
147f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOptionsPrefix(Mat,const char*[]);
148f69a0ea3SMatthew Knepley 
14930de9b25SBarry Smith /*MC
15030de9b25SBarry Smith    MatRegisterDynamic - Adds a new matrix type
15130de9b25SBarry Smith 
15230de9b25SBarry Smith    Synopsis:
153c1ac3661SBarry Smith    PetscErrorCode MatRegisterDynamic(char *name,char *path,char *name_create,PetscErrorCode (*routine_create)(Mat))
15430de9b25SBarry Smith 
15530de9b25SBarry Smith    Not Collective
15630de9b25SBarry Smith 
15730de9b25SBarry Smith    Input Parameters:
15830de9b25SBarry Smith +  name - name of a new user-defined matrix type
15930de9b25SBarry Smith .  path - path (either absolute or relative) the library containing this solver
16030de9b25SBarry Smith .  name_create - name of routine to create method context
16130de9b25SBarry Smith -  routine_create - routine to create method context
16230de9b25SBarry Smith 
16330de9b25SBarry Smith    Notes:
16430de9b25SBarry Smith    MatRegisterDynamic() may be called multiple times to add several user-defined solvers.
16530de9b25SBarry Smith 
16630de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
16730de9b25SBarry Smith    is ignored.
16830de9b25SBarry Smith 
16930de9b25SBarry Smith    Sample usage:
17030de9b25SBarry Smith .vb
17130de9b25SBarry Smith    MatRegisterDynamic("my_mat",/home/username/my_lib/lib/libO/solaris/mylib.a,
17230de9b25SBarry Smith                "MyMatCreate",MyMatCreate);
17330de9b25SBarry Smith .ve
17430de9b25SBarry Smith 
17530de9b25SBarry Smith    Then, your solver can be chosen with the procedural interface via
17630de9b25SBarry Smith $     MatSetType(Mat,"my_mat")
17730de9b25SBarry Smith    or at runtime via the option
17830de9b25SBarry Smith $     -mat_type my_mat
17930de9b25SBarry Smith 
18030de9b25SBarry Smith    Level: advanced
18130de9b25SBarry Smith 
182ab901514SBarry Smith    Notes: ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values.
18330de9b25SBarry Smith          If your function is not being put into a shared library then use VecRegister() instead
18430de9b25SBarry Smith 
18530de9b25SBarry Smith .keywords: Mat, register
18630de9b25SBarry Smith 
18730de9b25SBarry Smith .seealso: MatRegisterAll(), MatRegisterDestroy()
18830de9b25SBarry Smith 
18930de9b25SBarry Smith M*/
190273d9f13SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
191273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,0)
192273d9f13SBarry Smith #else
193273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,d)
19430de9b25SBarry Smith #endif
19530de9b25SBarry Smith 
1969c666560SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRegisterAll(const char[]);
197273d9f13SBarry Smith extern PetscTruth MatRegisterAllCalled;
198b0a32e0cSBarry Smith extern PetscFList MatList;
19928988994SBarry Smith 
200be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqDense(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*);
201be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIDense(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*);
202be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
203ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,nz,nnz,&A),Mat,A)
204ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,&A),Mat,A)
205ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,0,nnz,&A),Mat,A)
206ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,&A),Mat,A)
207ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,A))
208ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,m,n,0,nnz,A))
209ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,A))
210be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
211ba966639SSatish 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)
212ba966639SSatish 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)
213ba966639SSatish 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)
214ba966639SSatish 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)
215ba966639SSatish 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))
216ba966639SSatish 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))
217ba966639SSatish 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))
218ba966639SSatish 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)
219ba966639SSatish 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)
220ba966639SSatish 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)
221ba966639SSatish 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)
222ba966639SSatish 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))
223ba966639SSatish 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))
224ba966639SSatish 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))
22563ab93bfSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *);
2268d7a6e47SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJWithSplitArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],PetscInt[],PetscInt[],PetscScalar[],Mat*);
2278d7a6e47SBarry Smith 
228be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIRowbs(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
229be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
230ba966639SSatish 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)
231ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A)
232ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A)
233ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A)
234ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A))
235ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A))
236ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A))
237be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
238ba966639SSatish 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)
239ba966639SSatish 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)
240ba966639SSatish 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)
241ba966639SSatish 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)
242ba966639SSatish 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))
243ba966639SSatish 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))
244ba966639SSatish 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))
245ba966639SSatish 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)
246ba966639SSatish 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)
247ba966639SSatish 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)
248ba966639SSatish 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)
249ba966639SSatish 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))
250ba966639SSatish 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))
251ba966639SSatish 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))
252be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAdj(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscInt[],Mat*);
253be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
254ba966639SSatish 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)
255ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A)
256ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A)
257ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A)
258ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A))
259ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A))
260ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A))
261be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPISBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
262ba966639SSatish 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)
263ba966639SSatish 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)
264ba966639SSatish 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)
265ba966639SSatish 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)
266ba966639SSatish 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))
267ba966639SSatish 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))
268ba966639SSatish 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))
269ba966639SSatish 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)
270ba966639SSatish 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)
271ba966639SSatish 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)
272ba966639SSatish 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)
273ba966639SSatish 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))
274ba966639SSatish 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))
275ba966639SSatish 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))
276be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateShell(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,void *,Mat*);
277ba966639SSatish 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)
278ba966639SSatish Balay PetscPolymorphicFunction(MatCreateShell,(PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx),(PETSC_COMM_WORLD,m,n,M,N,ctx,&A),Mat,A)
279be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateAdic(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,void (*)(void),Mat*);
280be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateNormal(Mat,Mat*);
281ba966639SSatish Balay PetscPolymorphicFunction(MatCreateNormal,(Mat mat),(mat,&A),Mat,A)
282ab92ecdeSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateLRC(Mat,Mat,Mat,Mat*);
283284134d9SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateIS(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,ISLocalToGlobalMapping,Mat*);
2841472f72bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
2851472f72bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPICRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
2862a6744ebSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateScatter(MPI_Comm,VecScatter,Mat*);
2872a6744ebSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScatterSetVecScatter(Mat,VecScatter);
2882a6744ebSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScatterGetVecScatter(Mat,VecScatter*);
2898cdf2d9bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateBlockMat(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,Mat*);
290793850ffSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompositeAddMat(Mat,Mat);
291b52f573bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompositeMerge(Mat);
292793850ffSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateComposite(MPI_Comm,PetscInt,const Mat*,Mat*);
2935a7f1df3SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqFFTW(MPI_Comm,PetscInt,const PetscInt[],Mat*);
2941472f72bSBarry Smith 
295f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUp(Mat);
296be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroy(Mat);
29721c89e3eSBarry Smith 
2980c451bc4SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConjugate(Mat);
29999cafbc1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRealPart(Mat);
30099cafbc1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatImaginaryPart(Mat);
30199cafbc1SBarry Smith 
3028ed539a5SBarry Smith /* ------------------------------------------------------------*/
303be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
304be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
30587d4246cSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesRow(Mat,PetscInt,const PetscScalar[]);
306f259bd47SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesRowLocal(Mat,PetscInt,const PetscScalar[]);
30784cb2905SBarry Smith 
3082ef4de8bSBarry Smith /*S
3092ef4de8bSBarry Smith      MatStencil - Data structure (C struct) for storing information about a single row or
3102ef4de8bSBarry Smith         column of a matrix as index on an associated grid.
3112ef4de8bSBarry Smith 
3122ef4de8bSBarry Smith    Level: beginner
3132ef4de8bSBarry Smith 
3142ef4de8bSBarry Smith   Concepts: matrix; linear operator
3152ef4de8bSBarry Smith 
316d7bd93baSBarry Smith .seealso:  MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockStencil()
3172ef4de8bSBarry Smith S*/
318435da068SBarry Smith typedef struct {
319c1ac3661SBarry Smith   PetscInt k,j,i,c;
320435da068SBarry Smith } MatStencil;
3212ef4de8bSBarry Smith 
322be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
323be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
324be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt);
325435da068SBarry Smith 
326be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetColoring(Mat,ISColoring);
327be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesAdic(Mat,void*);
328be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesAdifor(Mat,PetscInt,void*);
3293a7fca6bSBarry Smith 
330d91e6319SBarry Smith /*E
331d91e6319SBarry Smith     MatAssemblyType - Indicates if the matrix is now to be used, or if you plan
332d91e6319SBarry Smith      to continue to add values to it
333d91e6319SBarry Smith 
334d91e6319SBarry Smith     Level: beginner
335d91e6319SBarry Smith 
336d91e6319SBarry Smith .seealso: MatAssemblyBegin(), MatAssemblyEnd()
337d91e6319SBarry Smith E*/
3386d4a8577SBarry Smith typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType;
339be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssemblyBegin(Mat,MatAssemblyType);
340be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssemblyEnd(Mat,MatAssemblyType);
341be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssembled(Mat,PetscTruth*);
3424f9c727eSBarry Smith 
3431d73ed98SBarry Smith 
34430de9b25SBarry Smith 
345d91e6319SBarry Smith /*E
346d91e6319SBarry Smith     MatOption - Options that may be set for a matrix and its behavior or storage
347d91e6319SBarry Smith 
348d91e6319SBarry Smith     Level: beginner
349d91e6319SBarry Smith 
3500a835dfdSSatish Balay    Any additions/changes here MUST also be made in include/finclude/petscmat.h
351d91e6319SBarry Smith 
352d91e6319SBarry Smith .seealso: MatSetOption()
353d91e6319SBarry Smith E*/
354512a5fc5SBarry Smith typedef enum {MAT_ROW_ORIENTED,MAT_NEW_NONZERO_LOCATIONS,
3554e0d8c25SBarry Smith               MAT_SYMMETRIC,
3564e0d8c25SBarry Smith               MAT_STRUCTURALLY_SYMMETRIC,
3578047e77bSBarry Smith               MAT_NEW_DIAGONALS,MAT_IGNORE_OFF_PROC_ENTRIES,
3584e0d8c25SBarry Smith               MAT_NEW_NONZERO_LOCATION_ERR,
3594e0d8c25SBarry Smith               MAT_NEW_NONZERO_ALLOCATION_ERR,MAT_USE_HASH_TABLE,
3604e0d8c25SBarry Smith               MAT_KEEP_ZEROED_ROWS,MAT_IGNORE_ZERO_ENTRIES,
3614e0d8c25SBarry Smith               MAT_USE_INODES,
3624e0d8c25SBarry Smith               MAT_HERMITIAN,
3634e0d8c25SBarry Smith               MAT_SYMMETRY_ETERNAL,
3644e0d8c25SBarry Smith               MAT_USE_COMPRESSEDROW,
3654e0d8c25SBarry Smith               MAT_IGNORE_LOWER_TRIANGULAR,MAT_ERROR_LOWER_TRIANGULAR,
3664e0d8c25SBarry Smith               MAT_GETROW_UPPERTRIANGULAR} MatOption;
367290bbb0aSBarry Smith extern const char *MatOptions[];
3684e0d8c25SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetOption(Mat,MatOption,PetscTruth);
369a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetType(Mat,const MatType*);
370a313700dSBarry Smith PetscPolymorphicFunction(MatGetType,(Mat mat),(mat,&t),const MatType,t)
37184cb2905SBarry Smith 
372be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]);
373be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
374be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
375f5edf698SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowUpperTriangular(Mat);
376f5edf698SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRowUpperTriangular(Mat);
377be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
378be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
379be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumnVector(Mat,Vec,PetscInt);
380be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetArray(Mat,PetscScalar *[]);
381ba966639SSatish Balay PetscPolymorphicFunction(MatGetArray,(Mat mat),(mat,&a),PetscScalar*,a)
382be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreArray(Mat,PetscScalar *[]);
383be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBlockSize(Mat,PetscInt *);
384ba966639SSatish Balay PetscPolymorphicFunction(MatGetBlockSize,(Mat mat),(mat,&a),PetscInt,a)
385be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetBlockSize(Mat,PetscInt);
3867b80b807SBarry Smith 
3871620fd73SBarry Smith 
388be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMult(Mat,Vec,Vec);
389be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultAdd(Mat,Vec,Vec,Vec);
390ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
391be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTranspose(Mat,Vec,Vec);
392be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsTranspose(Mat,Mat,PetscReal,PetscTruth*);
393ba966639SSatish Balay PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B,PetscReal tol),(A,B,tol,&t),PetscTruth,t)
394e37f64cdSBarry Smith PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B),(A,B,0,&t),PetscTruth,t)
3951cbb95d3SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitianTranspose(Mat,Mat,PetscReal,PetscTruth*);
396be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeAdd(Mat,Vec,Vec,Vec);
397ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultTransposeAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
398be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultConstrained(Mat,Vec,Vec);
399be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeConstrained(Mat,Vec,Vec);
4002bebee5dSHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatSolve(Mat,Mat,Mat);
401f5edf698SHong Zhang 
402d91e6319SBarry Smith /*E
403d91e6319SBarry Smith     MatDuplicateOption - Indicates if a duplicated sparse matrix should have
404d91e6319SBarry Smith   its numerical values copied over or just its nonzero structure.
405d91e6319SBarry Smith 
406d91e6319SBarry Smith     Level: beginner
407d91e6319SBarry Smith 
408d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
409d91e6319SBarry Smith 
410d91e6319SBarry Smith .seealso: MatDuplicate()
411d91e6319SBarry Smith E*/
4122e8a6d31SBarry Smith typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES} MatDuplicateOption;
4132e8a6d31SBarry Smith 
414a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert(Mat,const MatType,MatReuse,Mat*);
415a313700dSBarry Smith PetscPolymorphicFunction(MatConvert,(Mat A,const MatType t),(A,t,MAT_INITIAL_MATRIX,&a),Mat,a)
416be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDuplicate(Mat,MatDuplicateOption,Mat*);
417ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A,MatDuplicateOption o),(A,o,&a),Mat,a)
418ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A),(A,MAT_COPY_VALUES,&a),Mat,a)
41994a9d846SBarry Smith 
420d91e6319SBarry Smith /*E
421d91e6319SBarry Smith     MatStructure - Indicates if the matrix has the same nonzero structure
422d91e6319SBarry Smith 
423d91e6319SBarry Smith     Level: beginner
424d91e6319SBarry Smith 
425d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
426d91e6319SBarry Smith 
42794b7f48cSBarry Smith .seealso: MatCopy(), KSPSetOperators(), PCSetOperators()
428d91e6319SBarry Smith E*/
429c537a176SHong Zhang typedef enum {SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER,SUBSET_NONZERO_PATTERN} MatStructure;
430cb5b572fSBarry Smith 
431be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCopy(Mat,Mat,MatStructure);
432be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatView(Mat,PetscViewer);
433be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsSymmetric(Mat,PetscReal,PetscTruth*);
434ba966639SSatish Balay PetscPolymorphicFunction(MatIsSymmetric,(Mat A,PetscReal tol),(A,tol,&t),PetscTruth,t)
435e37f64cdSBarry Smith PetscPolymorphicFunction(MatIsSymmetric,(Mat A),(A,0,&t),PetscTruth,t)
436be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsStructurallySymmetric(Mat,PetscTruth*);
437ba966639SSatish Balay PetscPolymorphicFunction(MatIsStructurallySymmetric,(Mat A),(A,&t),PetscTruth,t)
4381cbb95d3SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitian(Mat,PetscReal,PetscTruth*);
4391cbb95d3SBarry Smith PetscPolymorphicFunction(MatIsHermitian,(Mat A),(A,0,&t),PetscTruth,t)
440be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsSymmetricKnown(Mat,PetscTruth*,PetscTruth*);
441be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitianKnown(Mat,PetscTruth*,PetscTruth*);
44264c62002SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMissingDiagonal(Mat,PetscTruth *,PetscInt *);
443a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLoad(PetscViewer,const MatType,Mat*);
444a313700dSBarry Smith PetscPolymorphicFunction(MatLoad,(PetscViewer v,const MatType t),(v,t,&a),Mat,a)
4457b80b807SBarry Smith 
4468f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *);
4478f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRowIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *);
4488f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumnIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *);
4498f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreColumnIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *);
450d4fbbf0eSBarry Smith 
451d91e6319SBarry Smith /*S
452d91e6319SBarry Smith      MatInfo - Context of matrix information, used with MatGetInfo()
453d91e6319SBarry Smith 
454d91e6319SBarry Smith    In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE
455d91e6319SBarry Smith 
456d91e6319SBarry Smith    Level: intermediate
457d91e6319SBarry Smith 
458d91e6319SBarry Smith   Concepts: matrix^nonzero information
459d91e6319SBarry Smith 
460d9274352SBarry Smith .seealso:  MatGetInfo(), MatInfoType
461d91e6319SBarry Smith S*/
4624e220ebcSLois Curfman McInnes typedef struct {
463b0a32e0cSBarry Smith   PetscLogDouble rows_global,columns_global;         /* number of global rows and columns */
464b0a32e0cSBarry Smith   PetscLogDouble rows_local,columns_local;           /* number of local rows and columns */
465b0a32e0cSBarry Smith   PetscLogDouble block_size;                         /* block size */
466b0a32e0cSBarry Smith   PetscLogDouble nz_allocated,nz_used,nz_unneeded;   /* number of nonzeros */
467b0a32e0cSBarry Smith   PetscLogDouble memory;                             /* memory allocated */
468b0a32e0cSBarry Smith   PetscLogDouble assemblies;                         /* number of matrix assemblies called */
469b0a32e0cSBarry Smith   PetscLogDouble mallocs;                            /* number of mallocs during MatSetValues() */
470b0a32e0cSBarry Smith   PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */
471b0a32e0cSBarry Smith   PetscLogDouble factor_mallocs;                     /* number of mallocs during factorization */
4724e220ebcSLois Curfman McInnes } MatInfo;
4734e220ebcSLois Curfman McInnes 
474d9274352SBarry Smith /*E
475d9274352SBarry Smith     MatInfoType - Indicates if you want information about the local part of the matrix,
476d9274352SBarry Smith      the entire parallel matrix or the maximum over all the local parts.
477d9274352SBarry Smith 
478d9274352SBarry Smith     Level: beginner
479d9274352SBarry Smith 
480d9274352SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
481d9274352SBarry Smith 
482d9274352SBarry Smith .seealso: MatGetInfo(), MatInfo
483d9274352SBarry Smith E*/
4847b80b807SBarry Smith typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType;
485be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetInfo(Mat,MatInfoType,MatInfo*);
486be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatValid(Mat,PetscTruth*);
487be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetDiagonal(Mat,Vec);
488985db425SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMax(Mat,Vec,PetscInt[]);
489985db425SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMin(Mat,Vec,PetscInt[]);
490985db425SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMaxAbs(Mat,Vec,PetscInt[]);
49186697f06SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowSum(Mat,Vec);
492fc4dec0aSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatTranspose(Mat,MatReuse,Mat*);
493fc4dec0aSBarry Smith PetscPolymorphicFunction(MatTranspose,(Mat A),(A,MAT_INITIAL_MATRIX,&t),Mat,t)
494be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPermute(Mat,IS,IS,Mat *);
495ba966639SSatish Balay PetscPolymorphicFunction(MatPermute,(Mat A,IS is1,IS is2),(A,is1,is2,&t),Mat,t)
496be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPermuteSparsify(Mat,PetscInt,PetscReal,PetscReal,IS,IS,Mat *);
497be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalScale(Mat,Vec,Vec);
498be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalSet(Mat,Vec,InsertMode);
499be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatEqual(Mat,Mat,PetscTruth*);
500ba966639SSatish Balay PetscPolymorphicFunction(MatEqual,(Mat A,Mat B),(A,B,&t),PetscTruth,t)
501be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultEqual(Mat,Mat,PetscInt,PetscTruth*);
502be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultAddEqual(Mat,Mat,PetscInt,PetscTruth*);
503be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeEqual(Mat,Mat,PetscInt,PetscTruth*);
504be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscTruth*);
5057b80b807SBarry Smith 
506be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNorm(Mat,NormType,PetscReal *);
507ba966639SSatish Balay PetscPolymorphicFunction(MatNorm,(Mat A,NormType t),(A,t,&n),PetscReal,n)
508be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroEntries(Mat);
509f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar);
510f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsIS(Mat,IS,PetscScalar);
511f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroColumns(Mat,PetscInt,const PetscInt [],const PetscScalar*);
512f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroColumnsIS(Mat,IS,const PetscScalar*);
5137b80b807SBarry Smith 
514be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatUseScaledForm(Mat,PetscTruth);
515be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScaleSystem(Mat,Vec,Vec);
516be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatUnScaleSystem(Mat,Vec,Vec);
5175ef9f2a5SBarry Smith 
518be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSize(Mat,PetscInt*,PetscInt*);
519be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalSize(Mat,PetscInt*,PetscInt*);
520be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRange(Mat,PetscInt*,PetscInt*);
5218ee2e534SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRanges(Mat,const PetscInt**);
522f4b8409dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRangeColumn(Mat,PetscInt*,PetscInt*);
523*72ca8751SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRangesColumn(Mat,const PetscInt**);
5247b80b807SBarry Smith 
525be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
526be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroyMatrices(PetscInt,Mat *[]);
527be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrix(Mat,IS,IS,PetscInt,MatReuse,Mat *);
528abf0e977SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrixRaw(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscInt,MatReuse,Mat *);
529829201f2SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSeqNonzeroStructure(Mat,Mat *[]);
530829201f2SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroySeqNonzeroStructure(Mat *[]);
5315494a064SHong Zhang 
532be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
533be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPI(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*);
534be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPISymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*);
535be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPINumeric(Mat,Mat);
536be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroy_MPIAIJ_SeqsToMPI(Mat);
537be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalMat(Mat,MatReuse,Mat*);
538be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*);
539be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,PetscInt*,Mat*);
540dd6ea824SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBrowsOfAoCols(Mat,Mat,MatReuse,PetscInt**,MatScalar**,Mat*);
54143eb5e2fSMatthew Knepley #if defined (PETSC_USE_CTABLE)
542cd116e26SSatish Balay #include "petscctable.h"
54343eb5e2fSMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetCommunicationStructs(Mat, Vec *, PetscTable *, VecScatter *);
54443eb5e2fSMatthew Knepley #else
54543eb5e2fSMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetCommunicationStructs(Mat, Vec *, PetscInt *[], VecScatter *);
54643eb5e2fSMatthew Knepley #endif
5478efafbd8SBarry Smith 
548be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt);
5497b80b807SBarry Smith 
550be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*);
551be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultSymbolic(Mat,Mat,PetscReal,Mat*);
552be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultNumeric(Mat,Mat,Mat);
55322440eb1SKris Buschelman 
554be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*);
555be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAPSymbolic(Mat,Mat,PetscReal,Mat*);
556be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAPNumeric(Mat,Mat,Mat);
55722440eb1SKris Buschelman 
558be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTranspose(Mat,Mat,MatReuse,PetscReal,Mat*);
559be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTransposeSymbolic(Mat,Mat,PetscReal,Mat*);
560be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTransposeNumeric(Mat,Mat,Mat);
561bc011b1eSHong Zhang 
562f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAXPY(Mat,PetscScalar,Mat,MatStructure);
56304aac2b0SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAYPX(Mat,PetscScalar,Mat,MatStructure);
564be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompress(Mat);
5651c741599SHong Zhang 
566f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScale(Mat,PetscScalar);
567f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShift(Mat,PetscScalar);
5687b80b807SBarry Smith 
569be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping);
570be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetLocalToGlobalMappingBlock(Mat,ISLocalToGlobalMapping);
571f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar);
572f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsLocalIS(Mat,IS,PetscScalar);
573be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
574be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
575052efed2SBarry Smith 
576be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStashSetInitialSize(Mat,PetscInt,PetscInt);
577be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
57890f02eecSBarry Smith 
579be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInterpolate(Mat,Vec,Vec);
5800c451bc4SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInterpolateAdd(Mat,Vec,Vec,Vec);
581ba966639SSatish Balay PetscPolymorphicSubroutine(MatInterpolateAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
582be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestrict(Mat,Vec,Vec);
583be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetVecs(Mat,Vec*,Vec*);
58469db28dcSHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRedundantMatrix(Mat,PetscInt,MPI_Comm,PetscInt,MatReuse,Mat*);
585bd481603SBarry Smith 
586bd481603SBarry Smith /*MC
5871620fd73SBarry Smith    MatSetValue - Set a single entry into a matrix.
5881620fd73SBarry Smith 
5891620fd73SBarry Smith    Not collective
5901620fd73SBarry Smith 
5911620fd73SBarry Smith    Input Parameters:
5921620fd73SBarry Smith +  m - the matrix
5931620fd73SBarry Smith .  row - the row location of the entry
5941620fd73SBarry Smith .  col - the column location of the entry
5951620fd73SBarry Smith .  value - the value to insert
5961620fd73SBarry Smith -  mode - either INSERT_VALUES or ADD_VALUES
5971620fd73SBarry Smith 
5981620fd73SBarry Smith    Notes:
5991620fd73SBarry Smith    For efficiency one should use MatSetValues() and set several or many
6001620fd73SBarry Smith    values simultaneously if possible.
6011620fd73SBarry Smith 
6021620fd73SBarry Smith    Level: beginner
6031620fd73SBarry Smith 
6041620fd73SBarry Smith .seealso: MatSetValues(), MatSetValueLocal()
6051620fd73SBarry Smith M*/
6061620fd73SBarry 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);}
6071620fd73SBarry Smith 
6081620fd73SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatGetValue(Mat v,PetscInt i,PetscInt j,PetscScalar *va) {return MatGetValues(v,1,&i,1,&j,va);}
6091620fd73SBarry Smith 
6101620fd73SBarry 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);}
6111620fd73SBarry Smith 
6121620fd73SBarry Smith /*MC
6130d2a7ff2SSatish Balay    MatPreallocateInitialize - Begins the block of code that will count the number of nonzeros per
614bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
615bd481603SBarry Smith 
616bd481603SBarry Smith    Synopsis:
617c1ac3661SBarry Smith    PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)
618bd481603SBarry Smith 
619bd481603SBarry Smith    Collective on MPI_Comm
620bd481603SBarry Smith 
621bd481603SBarry Smith    Input Parameters:
622bd481603SBarry Smith +  comm - the communicator that will share the eventually allocated matrix
623859fcb39SBarry Smith .  nrows - the number of LOCAL rows in the matrix
624859fcb39SBarry Smith -  ncols - the number of LOCAL columns in the matrix
625bd481603SBarry Smith 
626bd481603SBarry Smith    Output Parameters:
627bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
628bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
629bd481603SBarry Smith 
630bd481603SBarry Smith 
631bd481603SBarry Smith    Level: intermediate
632bd481603SBarry Smith 
633bd481603SBarry Smith    Notes:
634bd481603SBarry Smith    See the chapter in the users manual on performance for more details
635bd481603SBarry Smith 
6361d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
637bd481603SBarry Smith 
638bd481603SBarry Smith    Use MatPreallocateInitializeSymmetric() for symmetric matrices (MPISBAIJ matrices)
639bd481603SBarry Smith 
6401620fd73SBarry Smith    This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize().
6411620fd73SBarry Smith 
642bd481603SBarry Smith   Concepts: preallocation^Matrix
643bd481603SBarry Smith 
644bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
645bd481603SBarry Smith           MatPreallocateInitializeSymmetric(), MatPreallocateSymmetricSetLocal()
646bd481603SBarry Smith M*/
647c4f061fbSSatish Balay #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \
6487c922b88SBarry Smith { \
649a89bc211SBarry Smith   PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__start,__end; \
650a89bc211SBarry Smith   _4_ierr = PetscMalloc2(__nrows,PetscInt,&dnz,__nrows,PetscInt,&onz);CHKERRQ(_4_ierr); \
651a89bc211SBarry Smith   _4_ierr = PetscMemzero(dnz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
652a89bc211SBarry Smith   _4_ierr = PetscMemzero(onz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
653c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ctmp;\
654a89bc211SBarry Smith   _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows;
6557c922b88SBarry Smith 
656bd481603SBarry Smith /*MC
6570d2a7ff2SSatish Balay    MatPreallocateSymmetricInitialize - Begins the block of code that will count the number of nonzeros per
658bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
659bd481603SBarry Smith 
660bd481603SBarry Smith    Synopsis:
661c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)
662bd481603SBarry Smith 
663bd481603SBarry Smith    Collective on MPI_Comm
664bd481603SBarry Smith 
665bd481603SBarry Smith    Input Parameters:
666bd481603SBarry Smith +  comm - the communicator that will share the eventually allocated matrix
667859fcb39SBarry Smith .  nrows - the number of LOCAL rows in the matrix
668859fcb39SBarry Smith -  ncols - the number of LOCAL columns in the matrix
669bd481603SBarry Smith 
670bd481603SBarry Smith    Output Parameters:
671bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
672bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
673bd481603SBarry Smith 
674bd481603SBarry Smith 
675bd481603SBarry Smith    Level: intermediate
676bd481603SBarry Smith 
677bd481603SBarry Smith    Notes:
678bd481603SBarry Smith    See the chapter in the users manual on performance for more details
679bd481603SBarry Smith 
6801d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
681bd481603SBarry Smith 
6821620fd73SBarry Smith    This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize().
6831620fd73SBarry Smith 
684bd481603SBarry Smith   Concepts: preallocation^Matrix
685bd481603SBarry Smith 
686bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
687bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal()
688bd481603SBarry Smith M*/
689222b16d4SBarry Smith #define MatPreallocateSymmetricInitialize(comm,nrows,ncols,dnz,onz) 0; \
690222b16d4SBarry Smith { \
691a89bc211SBarry Smith   PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__end; \
692a89bc211SBarry Smith   _4_ierr = PetscMalloc2(__nrows,PetscInt,&dnz,__nrows,PetscInt,&onz);CHKERRQ(_4_ierr); \
693a89bc211SBarry Smith   _4_ierr = PetscMemzero(dnz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
694a89bc211SBarry Smith   _4_ierr = PetscMemzero(onz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
695c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr);\
696a89bc211SBarry Smith   _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows;
697222b16d4SBarry Smith 
698bd481603SBarry Smith /*MC
699bd481603SBarry Smith    MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
700bd481603SBarry Smith        inserted using a local number of the rows and columns
701bd481603SBarry Smith 
702bd481603SBarry Smith    Synopsis:
703c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
704bd481603SBarry Smith 
705bd481603SBarry Smith    Not Collective
706bd481603SBarry Smith 
707bd481603SBarry Smith    Input Parameters:
708bd481603SBarry Smith +  map - the mapping between local numbering and global numbering
709bd481603SBarry Smith .  nrows - the number of rows indicated
7101d73ed98SBarry Smith .  rows - the indices of the rows
711bd481603SBarry Smith .  ncols - the number of columns in the matrix
712bd481603SBarry Smith .  cols - the columns indicated
713bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
714bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
715bd481603SBarry Smith 
716bd481603SBarry Smith 
717bd481603SBarry Smith    Level: intermediate
718bd481603SBarry Smith 
719bd481603SBarry Smith    Notes:
720bd481603SBarry Smith    See the chapter in the users manual on performance for more details
721bd481603SBarry Smith 
7221d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
723bd481603SBarry Smith 
724bd481603SBarry Smith   Concepts: preallocation^Matrix
725bd481603SBarry Smith 
726bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
727bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal()
728bd481603SBarry Smith M*/
729c4f061fbSSatish Balay #define MatPreallocateSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\
730c4f061fbSSatish Balay {\
731c1ac3661SBarry Smith   PetscInt __l;\
732ef66eb69SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\
733ef66eb69SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\
734c4f061fbSSatish Balay   for (__l=0;__l<nrows;__l++) {\
735ef66eb69SBarry Smith     _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
736c4f061fbSSatish Balay   }\
737c4f061fbSSatish Balay }
738c4f061fbSSatish Balay 
739bd481603SBarry Smith /*MC
740bd481603SBarry Smith    MatPreallocateSymmetricSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
741bd481603SBarry Smith        inserted using a local number of the rows and columns
742bd481603SBarry Smith 
743bd481603SBarry Smith    Synopsis:
744c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
745bd481603SBarry Smith 
746bd481603SBarry Smith    Not Collective
747bd481603SBarry Smith 
748bd481603SBarry Smith    Input Parameters:
749bd481603SBarry Smith +  map - the mapping between local numbering and global numbering
750bd481603SBarry Smith .  nrows - the number of rows indicated
7511d73ed98SBarry Smith .  rows - the indices of the rows
752bd481603SBarry Smith .  ncols - the number of columns in the matrix
753bd481603SBarry Smith .  cols - the columns indicated
754bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
755bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
756bd481603SBarry Smith 
757bd481603SBarry Smith 
758bd481603SBarry Smith    Level: intermediate
759bd481603SBarry Smith 
760bd481603SBarry Smith    Notes:
761bd481603SBarry Smith    See the chapter in the users manual on performance for more details
762bd481603SBarry Smith 
763bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
764bd481603SBarry Smith 
765bd481603SBarry Smith   Concepts: preallocation^Matrix
766bd481603SBarry Smith 
767bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
768bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
769bd481603SBarry Smith M*/
770d3d32019SBarry Smith #define MatPreallocateSymmetricSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\
771d3d32019SBarry Smith {\
772c1ac3661SBarry Smith   PetscInt __l;\
773d3d32019SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\
774d3d32019SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\
775d3d32019SBarry Smith   for (__l=0;__l<nrows;__l++) {\
776d3d32019SBarry Smith     _4_ierr = MatPreallocateSymmetricSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
777d3d32019SBarry Smith   }\
778d3d32019SBarry Smith }
779d3d32019SBarry Smith 
780bd481603SBarry Smith /*MC
781bd481603SBarry Smith    MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
782bd481603SBarry Smith        inserted using a local number of the rows and columns
783bd481603SBarry Smith 
784bd481603SBarry Smith    Synopsis:
785c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
786bd481603SBarry Smith 
787bd481603SBarry Smith    Not Collective
788bd481603SBarry Smith 
789bd481603SBarry Smith    Input Parameters:
79064075487SBarry Smith +  row - the row
791bd481603SBarry Smith .  ncols - the number of columns in the matrix
792a50f8bf6SHong Zhang -  cols - the columns indicated
793a50f8bf6SHong Zhang 
794a50f8bf6SHong Zhang    Output Parameters:
795a50f8bf6SHong Zhang +  dnz - the array that will be passed to the matrix preallocation routines
796bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
797bd481603SBarry Smith 
798bd481603SBarry Smith 
799bd481603SBarry Smith    Level: intermediate
800bd481603SBarry Smith 
801bd481603SBarry Smith    Notes:
802bd481603SBarry Smith    See the chapter in the users manual on performance for more details
803bd481603SBarry Smith 
804bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
805bd481603SBarry Smith 
8061620fd73SBarry Smith    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
8071620fd73SBarry Smith 
808bd481603SBarry Smith   Concepts: preallocation^Matrix
809bd481603SBarry Smith 
810bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
811bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
812bd481603SBarry Smith M*/
813c4f061fbSSatish Balay #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\
814c1ac3661SBarry Smith { PetscInt __i; \
8158ee2e534SBarry Smith   if (row < __rstart) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Trying to set preallocation for row %D less than first local row %D",row,__rstart);\
816a89bc211SBarry 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);\
8177c922b88SBarry Smith   for (__i=0; __i<nc; __i++) {\
81864075487SBarry Smith     if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \
8197cc688d7SBarry Smith     else dnz[row - __rstart]++;\
8207c922b88SBarry Smith   }\
8217c922b88SBarry Smith }
8227c922b88SBarry Smith 
823bd481603SBarry Smith /*MC
824bd481603SBarry Smith    MatPreallocateSymmetricSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
825bd481603SBarry Smith        inserted using a local number of the rows and columns
826bd481603SBarry Smith 
827bd481603SBarry Smith    Synopsis:
828c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
829bd481603SBarry Smith 
830bd481603SBarry Smith    Not Collective
831bd481603SBarry Smith 
832bd481603SBarry Smith    Input Parameters:
833bd481603SBarry Smith +  nrows - the number of rows indicated
8341d73ed98SBarry Smith .  rows - the indices of the rows
835bd481603SBarry Smith .  ncols - the number of columns in the matrix
836bd481603SBarry Smith .  cols - the columns indicated
837bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
838bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
839bd481603SBarry Smith 
840bd481603SBarry Smith 
841bd481603SBarry Smith    Level: intermediate
842bd481603SBarry Smith 
843bd481603SBarry Smith    Notes:
844bd481603SBarry Smith    See the chapter in the users manual on performance for more details
845bd481603SBarry Smith 
846bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
847bd481603SBarry Smith 
8481620fd73SBarry Smith    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
8491620fd73SBarry Smith 
850bd481603SBarry Smith   Concepts: preallocation^Matrix
851bd481603SBarry Smith 
852bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
853bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
854bd481603SBarry Smith M*/
855d3d32019SBarry Smith #define MatPreallocateSymmetricSet(row,nc,cols,dnz,onz) 0;\
856c1ac3661SBarry Smith { PetscInt __i; \
857d3d32019SBarry Smith   for (__i=0; __i<nc; __i++) {\
858d3d32019SBarry Smith     if (cols[__i] >= __end) onz[row - __rstart]++; \
859d3d32019SBarry Smith     else if (cols[__i] >= row) dnz[row - __rstart]++;\
860d3d32019SBarry Smith   }\
861d3d32019SBarry Smith }
862d3d32019SBarry Smith 
863bd481603SBarry Smith /*MC
86416371a99SBarry Smith    MatPreallocateLocation -  An alternative to MatPreallocationSet() that puts the nonzero locations into the matrix if it exists
86516371a99SBarry Smith 
86616371a99SBarry Smith    Synopsis:
86716371a99SBarry Smith    PetscErrorCode MatPreallocateLocations(Mat A,PetscInt row,PetscInt ncols,PetscInt *cols,PetscInt *dnz,PetscInt *onz)
86816371a99SBarry Smith 
86916371a99SBarry Smith    Not Collective
87016371a99SBarry Smith 
87116371a99SBarry Smith    Input Parameters:
87216371a99SBarry Smith .  A - matrix
87316371a99SBarry Smith .  row - row where values exist (must be local to this process)
87416371a99SBarry Smith .  ncols - number of columns
87516371a99SBarry Smith .  cols - columns with nonzeros
87616371a99SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
87716371a99SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
87816371a99SBarry Smith 
87916371a99SBarry Smith 
88016371a99SBarry Smith    Level: intermediate
88116371a99SBarry Smith 
88216371a99SBarry Smith    Notes:
88316371a99SBarry Smith    See the chapter in the users manual on performance for more details
88416371a99SBarry Smith 
88516371a99SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
88616371a99SBarry Smith 
88716371a99SBarry Smith    This is a MACRO not a function because it uses a bunch of variables private to the MatPreallocation.... routines.
88816371a99SBarry Smith 
88916371a99SBarry Smith   Concepts: preallocation^Matrix
89016371a99SBarry Smith 
89116371a99SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
89216371a99SBarry Smith           MatPreallocateSymmetricInitialize(), MatPreallocateSymmetricSetLocal()
89316371a99SBarry Smith M*/
89416371a99SBarry 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);}
89516371a99SBarry Smith 
89616371a99SBarry Smith 
89716371a99SBarry Smith /*MC
8980d2a7ff2SSatish Balay    MatPreallocateFinalize - Ends the block of code that will count the number of nonzeros per
899bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
900bd481603SBarry Smith 
901bd481603SBarry Smith    Synopsis:
902c1ac3661SBarry Smith    PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz)
903bd481603SBarry Smith 
904bd481603SBarry Smith    Collective on MPI_Comm
905bd481603SBarry Smith 
906bd481603SBarry Smith    Input Parameters:
90716371a99SBarry Smith +  dnz - the array that was be passed to the matrix preallocation routines
908bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
909bd481603SBarry Smith 
910bd481603SBarry Smith 
911bd481603SBarry Smith    Level: intermediate
912bd481603SBarry Smith 
913bd481603SBarry Smith    Notes:
914bd481603SBarry Smith    See the chapter in the users manual on performance for more details
915bd481603SBarry Smith 
916bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
917bd481603SBarry Smith 
9181620fd73SBarry Smith    This is a MACRO not a function because it closes the { started in MatPreallocateInitialize().
9191620fd73SBarry Smith 
920bd481603SBarry Smith   Concepts: preallocation^Matrix
921bd481603SBarry Smith 
922bd481603SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
923bd481603SBarry Smith           MatPreallocateSymmetricInitialize(), MatPreallocateSymmetricSetLocal()
924bd481603SBarry Smith M*/
925a89bc211SBarry Smith #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree2(dnz,onz);CHKERRQ(_4_ierr);}
9267c922b88SBarry Smith 
927bd481603SBarry Smith 
928bd481603SBarry Smith 
9297b80b807SBarry Smith /* Routines unique to particular data structures */
930be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellGetContext(Mat,void **);
931ba966639SSatish Balay PetscPolymorphicFunction(MatShellGetContext,(Mat A),(A,&t),void*,t)
932698d4c6aSKris Buschelman 
933be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInodeAdjustForInodes(Mat,IS*,IS*);
934be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *);
935698d4c6aSKris Buschelman 
936be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetColumnIndices(Mat,PetscInt[]);
937be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetColumnIndices(Mat,PetscInt[]);
938be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
939c75a6043SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
940c75a6043SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
9417b80b807SBarry Smith 
942a96a251dSBarry Smith #define MAT_SKIP_ALLOCATION -4
943a96a251dSBarry Smith 
944be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
945ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz))
946be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
947ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqSBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz))
948be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]);
949ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqAIJSetPreallocation,(Mat A,const PetscInt nnz[]),(A,0,nnz))
950be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDensePreallocation(Mat,PetscScalar[]);
951273d9f13SBarry Smith 
952be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
953ba966639SSatish Balay PetscPolymorphicSubroutine(MatMPIBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[],const PetscInt onz[]),(A,bs,0,nnz,0,onz))
954be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
955be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
95653803f26SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []);
957725b52f3SLisandro Dalcin EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
958be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]);
959be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
960be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]);
961be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIDenseSetPreallocation(Mat,PetscScalar[]);
962284134d9SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation(Mat,PetscScalar[]);
963be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIRowbsSetPreallocation(Mat,PetscInt,const PetscInt[]);
964be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,PetscInt*[]);
965be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,PetscInt*[]);
966be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAdicSetLocalFunction(Mat,void (*)(void));
967273d9f13SBarry Smith 
968be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetLDA(Mat,PetscInt);
969ab92ecdeSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDenseGetLocalMatrix(Mat,Mat*);
9701b807ce4Svictorle 
971be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStoreValues(Mat);
972be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRetrieveValues(Mat);
9732e8a6d31SBarry Smith 
974be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDAADSetCtx(Mat,void*);
9753a7fca6bSBarry Smith 
9767b80b807SBarry Smith /*
9777b80b807SBarry Smith   These routines are not usually accessed directly, rather solving is
97894b7f48cSBarry Smith   done through the KSP and PC interfaces.
9797b80b807SBarry Smith */
9807b80b807SBarry Smith 
981d9274352SBarry Smith /*E
982d9274352SBarry Smith     MatOrderingType - String with the name of a PETSc matrix ordering or the creation function
983d9274352SBarry Smith        with an optional dynamic library name, for example
984d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:orderingcreate()
985d9274352SBarry Smith 
986d9274352SBarry Smith    Level: beginner
987d9274352SBarry Smith 
988e5a9bf91SBarry Smith    Cannot use const because the PC objects manipulate the string
989e5a9bf91SBarry Smith 
990d9274352SBarry Smith .seealso: MatGetOrdering()
991d9274352SBarry Smith E*/
9923eaad2caSSatish Balay #define MatOrderingType char*
993b12f92e5SBarry Smith #define MATORDERING_NATURAL   "natural"
994b12f92e5SBarry Smith #define MATORDERING_ND        "nd"
995b12f92e5SBarry Smith #define MATORDERING_1WD       "1wd"
996b12f92e5SBarry Smith #define MATORDERING_RCM       "rcm"
997b12f92e5SBarry Smith #define MATORDERING_QMD       "qmd"
998b12f92e5SBarry Smith #define MATORDERING_ROWLENGTH "rowlength"
99962152c8bSBarry Smith #define MATORDERING_DSC_ND    "dsc_nd"
100062152c8bSBarry Smith #define MATORDERING_DSC_MMD   "dsc_mmd"
100162152c8bSBarry Smith #define MATORDERING_DSC_MDF   "dsc_mdf"
1002c06d978dSMatthew Knepley #define MATORDERING_CONSTRAINED "constrained"
1003c06d978dSMatthew Knepley #define MATORDERING_IDENTITY  "identity"
1004c06d978dSMatthew Knepley #define MATORDERING_REVERSE   "reverse"
1005b12f92e5SBarry Smith 
1006f2416520SSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOrdering(Mat,const MatOrderingType,IS*,IS*);
1007be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOrderingList(PetscFList *list);
1008f2416520SSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,const MatOrderingType,IS*,IS*));
100930de9b25SBarry Smith 
101030de9b25SBarry Smith /*MC
101130de9b25SBarry Smith    MatOrderingRegisterDynamic - Adds a new sparse matrix ordering to the
101230de9b25SBarry Smith                                matrix package.
101330de9b25SBarry Smith 
101430de9b25SBarry Smith    Synopsis:
1015c1ac3661SBarry Smith    PetscErrorCode MatOrderingRegisterDynamic(char *name_ordering,char *path,char *name_create,PetscErrorCode (*routine_create)(MatOrdering))
101630de9b25SBarry Smith 
101730de9b25SBarry Smith    Not Collective
101830de9b25SBarry Smith 
101930de9b25SBarry Smith    Input Parameters:
102030de9b25SBarry Smith +  sname - name of ordering (for example MATORDERING_ND)
102130de9b25SBarry Smith .  path - location of library where creation routine is
102230de9b25SBarry Smith .  name - name of function that creates the ordering type,a string
102330de9b25SBarry Smith -  function - function pointer that creates the ordering
102430de9b25SBarry Smith 
102530de9b25SBarry Smith    Level: developer
102630de9b25SBarry Smith 
102730de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
102830de9b25SBarry Smith    is ignored.
102930de9b25SBarry Smith 
103030de9b25SBarry Smith    Sample usage:
103130de9b25SBarry Smith .vb
103230de9b25SBarry Smith    MatOrderingRegisterDynamic("my_order",/home/username/my_lib/lib/libO/solaris/mylib.a,
103330de9b25SBarry Smith                "MyOrder",MyOrder);
103430de9b25SBarry Smith .ve
103530de9b25SBarry Smith 
103630de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
103730de9b25SBarry Smith $     MatOrderingSetType(part,"my_order)
103830de9b25SBarry Smith    or at runtime via the option
10392401956bSBarry Smith $     -pc_factor_mat_ordering_type my_order
104030de9b25SBarry Smith 
1041ab901514SBarry Smith    ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values.
104230de9b25SBarry Smith 
104330de9b25SBarry Smith .keywords: matrix, ordering, register
104430de9b25SBarry Smith 
104530de9b25SBarry Smith .seealso: MatOrderingRegisterDestroy(), MatOrderingRegisterAll()
104630de9b25SBarry Smith M*/
1047aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1048f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,0)
1049b12f92e5SBarry Smith #else
1050f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,d)
1051b12f92e5SBarry Smith #endif
105230de9b25SBarry Smith 
1053be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegisterDestroy(void);
1054be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegisterAll(const char[]);
10552bad1931SBarry Smith extern PetscTruth MatOrderingRegisterAllCalled;
1056b0a32e0cSBarry Smith extern PetscFList MatOrderingList;
1057d4fbbf0eSBarry Smith 
1058be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS);
1059a2ce50c7SBarry Smith 
1060d91e6319SBarry Smith /*S
10612401956bSBarry Smith    MatFactorInfo - Data passed into the matrix factorization routines
1062b00f7748SHong Zhang 
106361cad860SBarry Smith    In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE, that is use
106461cad860SBarry Smith $     MatFactorInfo  info(MAT_FACTORINFO_SIZE)
1065b00f7748SHong Zhang 
106615e8a5b3SHong Zhang    Notes: These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC.
1067b00f7748SHong Zhang 
106861cad860SBarry Smith       You can use MatFactorInfoInitialize() to set default values.
106961cad860SBarry Smith 
1070b00f7748SHong Zhang    Level: developer
1071b00f7748SHong Zhang 
1072d7d82daaSBarry Smith .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(),
1073d7d82daaSBarry Smith           MatFactorInfoInitialize()
1074b00f7748SHong Zhang 
1075b00f7748SHong Zhang S*/
1076b00f7748SHong Zhang typedef struct {
10770a29876aSHong Zhang   PetscReal     shiftnz;        /* scaling of identity added to matrix to prevent zero pivots */
1078fbf22428SSatish Balay   PetscReal     shiftpd;        /* if true, shift until positive pivots */
10792cea7109SBarry Smith   PetscReal     shift_fraction; /* record shift fraction taken */
108015e8a5b3SHong Zhang   PetscReal     diagonal_fill;  /* force diagonal to fill in if initially not filled */
108115e8a5b3SHong Zhang   PetscReal     dt;             /* drop tolerance */
1082b00f7748SHong Zhang   PetscReal     dtcol;          /* tolerance for pivoting */
108315e8a5b3SHong Zhang   PetscReal     dtcount;        /* maximum nonzeros to be allowed per row */
108467e28bfeSBarry Smith   PetscReal     fill;           /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */
1085348344bbSBarry Smith   PetscReal     levels;         /* ICC/ILU(levels) */
1086bcd9e38bSBarry Smith   PetscReal     pivotinblocks;  /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0
1087bcd9e38bSBarry Smith                                    factorization may be faster if do not pivot */
108862bba022SBarry Smith   PetscReal     shiftinblocks;  /* if block in block factorization has zero pivot then shift diagonal until non-singular */
108915e8a5b3SHong Zhang   PetscReal     zeropivot;      /* pivot is called zero if less than this */
109015e8a5b3SHong Zhang } MatFactorInfo;
1091ffa6d0a5SLois Curfman McInnes 
1092be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFactorInfoInitialize(MatFactorInfo*);
1093be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactor(Mat,IS,MatFactorInfo*);
1094be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactorSymbolic(Mat,IS,MatFactorInfo*,Mat*);
1095be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactorNumeric(Mat,MatFactorInfo*,Mat*);
1096be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactor(Mat,IS,IS,MatFactorInfo*);
1097be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUFactor(Mat,IS,IS,MatFactorInfo*);
1098be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactorSymbolic(Mat,IS,IS,MatFactorInfo*,Mat*);
1099be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUFactorSymbolic(Mat,IS,IS,MatFactorInfo*,Mat*);
1100be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatICCFactorSymbolic(Mat,IS,MatFactorInfo*,Mat*);
1101be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatICCFactor(Mat,IS,MatFactorInfo*);
1102be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactorNumeric(Mat,MatFactorInfo*,Mat*);
1103be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUDTFactor(Mat,IS,IS,MatFactorInfo*,Mat *);
1104be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*);
1105be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolve(Mat,Vec,Vec);
1106be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatForwardSolve(Mat,Vec,Vec);
1107be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatBackwardSolve(Mat,Vec,Vec);
1108be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveAdd(Mat,Vec,Vec,Vec);
1109be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveTranspose(Mat,Vec,Vec);
1110be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveTransposeAdd(Mat,Vec,Vec,Vec);
1111be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolves(Mat,Vecs,Vecs);
11128ed539a5SBarry Smith 
1113be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUnfactored(Mat);
1114bb5a7306SBarry Smith 
1115d91e6319SBarry Smith /*E
1116d91e6319SBarry Smith     MatSORType - What type of (S)SOR to perform
1117bb1eb677SSatish Balay 
1118d91e6319SBarry Smith     Level: beginner
1119d91e6319SBarry Smith 
1120d9274352SBarry Smith    May be bitwise ORd together
1121d9274352SBarry Smith 
1122d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
1123d91e6319SBarry Smith 
11244e7234bfSBarry Smith    MatSORType may be bitwise ORd together, so do not change the numbers
11254e7234bfSBarry Smith 
1126a1d92eedSBarry Smith .seealso: MatRelax(), MatPBRelax()
1127d91e6319SBarry Smith E*/
1128ee50ffe9SBarry Smith typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3,
1129ee50ffe9SBarry Smith               SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8,
1130ee50ffe9SBarry Smith               SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16,
113184cb2905SBarry Smith               SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType;
1132be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRelax(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
1133be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPBRelax(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
11348ed539a5SBarry Smith 
1135d4fbbf0eSBarry Smith /*
1136639f9d9dSBarry Smith     These routines are for efficiently computing Jacobians via finite differences.
1137639f9d9dSBarry Smith */
1138b12f92e5SBarry Smith 
1139d9274352SBarry Smith /*E
1140d9274352SBarry Smith     MatColoringType - String with the name of a PETSc matrix coloring or the creation function
1141d9274352SBarry Smith        with an optional dynamic library name, for example
1142d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:coloringcreate()
1143d9274352SBarry Smith 
1144d9274352SBarry Smith    Level: beginner
1145d9274352SBarry Smith 
1146d9274352SBarry Smith .seealso: MatGetColoring()
1147d9274352SBarry Smith E*/
1148a313700dSBarry Smith #define MatColoringType char*
1149b12f92e5SBarry Smith #define MATCOLORING_NATURAL "natural"
1150b12f92e5SBarry Smith #define MATCOLORING_SL      "sl"
1151b12f92e5SBarry Smith #define MATCOLORING_LF      "lf"
1152b12f92e5SBarry Smith #define MATCOLORING_ID      "id"
1153b12f92e5SBarry Smith 
1154ab970d9bSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColoring(Mat,const MatColoringType,ISColoring*);
11552e90c967SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,MatColoringType,ISColoring *));
115630de9b25SBarry Smith 
115730de9b25SBarry Smith /*MC
115830de9b25SBarry Smith    MatColoringRegisterDynamic - Adds a new sparse matrix coloring to the
115930de9b25SBarry Smith                                matrix package.
116030de9b25SBarry Smith 
116130de9b25SBarry Smith    Synopsis:
1162c1ac3661SBarry Smith    PetscErrorCode MatColoringRegisterDynamic(char *name_coloring,char *path,char *name_create,PetscErrorCode (*routine_create)(MatColoring))
116330de9b25SBarry Smith 
116430de9b25SBarry Smith    Not Collective
116530de9b25SBarry Smith 
116630de9b25SBarry Smith    Input Parameters:
116730de9b25SBarry Smith +  sname - name of Coloring (for example MATCOLORING_SL)
116830de9b25SBarry Smith .  path - location of library where creation routine is
116930de9b25SBarry Smith .  name - name of function that creates the Coloring type, a string
117030de9b25SBarry Smith -  function - function pointer that creates the coloring
117130de9b25SBarry Smith 
117230de9b25SBarry Smith    Level: developer
117330de9b25SBarry Smith 
117430de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
117530de9b25SBarry Smith    is ignored.
117630de9b25SBarry Smith 
117730de9b25SBarry Smith    Sample usage:
117830de9b25SBarry Smith .vb
117930de9b25SBarry Smith    MatColoringRegisterDynamic("my_color",/home/username/my_lib/lib/libO/solaris/mylib.a,
118030de9b25SBarry Smith                "MyColor",MyColor);
118130de9b25SBarry Smith .ve
118230de9b25SBarry Smith 
118330de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
118430de9b25SBarry Smith $     MatColoringSetType(part,"my_color")
118530de9b25SBarry Smith    or at runtime via the option
118630de9b25SBarry Smith $     -mat_coloring_type my_color
118730de9b25SBarry Smith 
1188ab901514SBarry Smith    $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
118930de9b25SBarry Smith 
119030de9b25SBarry Smith .keywords: matrix, Coloring, register
119130de9b25SBarry Smith 
119230de9b25SBarry Smith .seealso: MatColoringRegisterDestroy(), MatColoringRegisterAll()
119330de9b25SBarry Smith M*/
1194aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1195f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,0)
1196b12f92e5SBarry Smith #else
1197f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,d)
1198b12f92e5SBarry Smith #endif
119930de9b25SBarry Smith 
12002bad1931SBarry Smith extern PetscTruth MatColoringRegisterAllCalled;
1201f1144a30SSatish Balay 
1202be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegisterAll(const char[]);
1203be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegisterDestroy(void);
1204be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*);
1205639f9d9dSBarry Smith 
1206d9274352SBarry Smith /*S
1207d9274352SBarry Smith      MatFDColoring - Object for computing a sparse Jacobian via finite differences
1208d9274352SBarry Smith         and coloring
1209639f9d9dSBarry Smith 
1210d9274352SBarry Smith    Level: beginner
1211d9274352SBarry Smith 
1212d9274352SBarry Smith   Concepts: coloring, sparse Jacobian, finite differences
1213d9274352SBarry Smith 
1214d9274352SBarry Smith .seealso:  MatFDColoringCreate()
1215d9274352SBarry Smith S*/
1216e2a1c21fSSatish Balay typedef struct _p_MatFDColoring* MatFDColoring;
1217639f9d9dSBarry Smith 
1218be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringCreate(Mat,ISColoring,MatFDColoring *);
1219be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringDestroy(MatFDColoring);
1220be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringView(MatFDColoring,PetscViewer);
1221be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*);
12224a9d489dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**);
1223be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal);
1224be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFrequency(MatFDColoring,PetscInt);
1225be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetFrequency(MatFDColoring,PetscInt*);
1226be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFromOptions(MatFDColoring);
1227be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApply(Mat,MatFDColoring,Vec,MatStructure*,void *);
1228be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApplyTS(Mat,MatFDColoring,PetscReal,Vec,MatStructure*,void *);
1229be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetRecompute(MatFDColoring);
1230be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetF(MatFDColoring,Vec);
1231be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,PetscInt*[]);
1232639f9d9dSBarry Smith /*
12330752156aSBarry Smith     These routines are for partitioning matrices: currently used only
12343eda8832SBarry Smith   for adjacency matrix, MatCreateMPIAdj().
12350752156aSBarry Smith */
1236ca161407SBarry Smith 
1237d9274352SBarry Smith /*S
1238d9274352SBarry Smith      MatPartitioning - Object for managing the partitioning of a matrix or graph
1239d9274352SBarry Smith 
1240d9274352SBarry Smith    Level: beginner
1241d9274352SBarry Smith 
1242d9274352SBarry Smith   Concepts: partitioning
1243d9274352SBarry Smith 
1244743a1026Svictorle .seealso:  MatPartitioningCreate(), MatPartitioningType
1245d9274352SBarry Smith S*/
124691e9ee9fSBarry Smith typedef struct _p_MatPartitioning* MatPartitioning;
1247d9274352SBarry Smith 
1248d9274352SBarry Smith /*E
12495bc6e7ccSvictorle     MatPartitioningType - String with the name of a PETSc matrix partitioning or the creation function
1250d9274352SBarry Smith        with an optional dynamic library name, for example
1251d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:partitioningcreate()
1252d9274352SBarry Smith 
1253d9274352SBarry Smith    Level: beginner
1254d9274352SBarry Smith 
1255b547a13bSvictorle .seealso: MatPartitioningCreate(), MatPartitioning
1256d9274352SBarry Smith E*/
1257a313700dSBarry Smith #define MatPartitioningType char*
12588ba1e511SMatthew Knepley #define MAT_PARTITIONING_CURRENT  "current"
1259c2ec2a73SHong Zhang #define MAT_PARTITIONING_SQUARE   "square"
12608ba1e511SMatthew Knepley #define MAT_PARTITIONING_PARMETIS "parmetis"
126171306c60Sjroman #define MAT_PARTITIONING_CHACO    "chaco"
126271306c60Sjroman #define MAT_PARTITIONING_JOSTLE   "jostle"
126371306c60Sjroman #define MAT_PARTITIONING_PARTY    "party"
126471306c60Sjroman #define MAT_PARTITIONING_SCOTCH   "scotch"
126571306c60Sjroman 
1266ca161407SBarry Smith 
1267be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningCreate(MPI_Comm,MatPartitioning*);
1268a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetType(MatPartitioning,const MatPartitioningType);
1269be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetNParts(MatPartitioning,PetscInt);
1270be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetAdjacency(MatPartitioning,Mat);
1271be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]);
1272be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []);
1273be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningApply(MatPartitioning,IS*);
1274be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningDestroy(MatPartitioning);
12752aabb6bbSBarry Smith 
1276be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatPartitioning));
127730de9b25SBarry Smith 
127830de9b25SBarry Smith /*MC
127930de9b25SBarry Smith    MatPartitioningRegisterDynamic - Adds a new sparse matrix partitioning to the
128030de9b25SBarry Smith    matrix package.
128130de9b25SBarry Smith 
128230de9b25SBarry Smith    Synopsis:
1283c1ac3661SBarry Smith    PetscErrorCode MatPartitioningRegisterDynamic(char *name_partitioning,char *path,char *name_create,PetscErrorCode (*routine_create)(MatPartitioning))
128430de9b25SBarry Smith 
128530de9b25SBarry Smith    Not Collective
128630de9b25SBarry Smith 
128730de9b25SBarry Smith    Input Parameters:
128830de9b25SBarry Smith +  sname - name of partitioning (for example MAT_PARTITIONING_CURRENT) or parmetis
128930de9b25SBarry Smith .  path - location of library where creation routine is
129030de9b25SBarry Smith .  name - name of function that creates the partitioning type, a string
129130de9b25SBarry Smith -  function - function pointer that creates the partitioning type
129230de9b25SBarry Smith 
129330de9b25SBarry Smith    Level: developer
129430de9b25SBarry Smith 
129530de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
129630de9b25SBarry Smith    is ignored.
129730de9b25SBarry Smith 
129830de9b25SBarry Smith    Sample usage:
129930de9b25SBarry Smith .vb
130030de9b25SBarry Smith    MatPartitioningRegisterDynamic("my_part",/home/username/my_lib/lib/libO/solaris/mylib.a,
130130de9b25SBarry Smith                "MyPartCreate",MyPartCreate);
130230de9b25SBarry Smith .ve
130330de9b25SBarry Smith 
130430de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
130530de9b25SBarry Smith $     MatPartitioningSetType(part,"my_part")
130630de9b25SBarry Smith    or at runtime via the option
130730de9b25SBarry Smith $     -mat_partitioning_type my_part
130830de9b25SBarry Smith 
1309ab901514SBarry Smith    $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
131030de9b25SBarry Smith 
131130de9b25SBarry Smith .keywords: matrix, partitioning, register
131230de9b25SBarry Smith 
131330de9b25SBarry Smith .seealso: MatPartitioningRegisterDestroy(), MatPartitioningRegisterAll()
131430de9b25SBarry Smith M*/
1315aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1316f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,0)
13172aabb6bbSBarry Smith #else
1318f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,d)
13192aabb6bbSBarry Smith #endif
13202aabb6bbSBarry Smith 
13212bad1931SBarry Smith extern PetscTruth MatPartitioningRegisterAllCalled;
1322f1144a30SSatish Balay 
1323be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegisterAll(const char[]);
1324be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegisterDestroy(void);
13252bad1931SBarry Smith 
1326be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningView(MatPartitioning,PetscViewer);
1327be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetFromOptions(MatPartitioning);
1328a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningGetType(MatPartitioning,const MatPartitioningType*);
1329ca161407SBarry Smith 
1330be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningParmetisSetCoarseSequential(MatPartitioning);
13310752156aSBarry Smith 
1332be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningJostleSetCoarseLevel(MatPartitioning,PetscReal);
1333be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningJostleSetCoarseSequential(MatPartitioning);
133471306c60Sjroman 
1335954b3ec6SBarry Smith typedef enum {MP_CHACO_MULTILEVEL_KL,MP_CHACO_SPECTRAL,MP_CHACO_LINEAR,MP_CHACO_RANDOM, MP_CHACO_SCATTERED} MPChacoGlobalType;
1336be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetGlobal(MatPartitioning, MPChacoGlobalType);
133771306c60Sjroman typedef enum { MP_CHACO_KERNIGHAN_LIN, MP_CHACO_NONE } MPChacoLocalType;
1338be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetLocal(MatPartitioning, MPChacoLocalType);
1339be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal);
134071306c60Sjroman typedef enum { MP_CHACO_LANCZOS, MP_CHACO_RQI_SYMMLQ } MPChacoEigenType;
1341be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType);
1342be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenTol(MatPartitioning, PetscReal);
1343be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenNumber(MatPartitioning, PetscInt);
134471306c60Sjroman 
134571306c60Sjroman #define MP_PARTY_OPT "opt"
134671306c60Sjroman #define MP_PARTY_LIN "lin"
134771306c60Sjroman #define MP_PARTY_SCA "sca"
134871306c60Sjroman #define MP_PARTY_RAN "ran"
134971306c60Sjroman #define MP_PARTY_GBF "gbf"
135071306c60Sjroman #define MP_PARTY_GCF "gcf"
135171306c60Sjroman #define MP_PARTY_BUB "bub"
135271306c60Sjroman #define MP_PARTY_DEF "def"
1353be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetGlobal(MatPartitioning, const char*);
135471306c60Sjroman #define MP_PARTY_HELPFUL_SETS "hs"
135571306c60Sjroman #define MP_PARTY_KERNIGHAN_LIN "kl"
135671306c60Sjroman #define MP_PARTY_NONE "no"
1357be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetLocal(MatPartitioning, const char*);
1358be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal);
1359be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetBipart(MatPartitioning,PetscTruth);
1360be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscTruth);
136171306c60Sjroman 
136271306c60Sjroman typedef enum { MP_SCOTCH_GREEDY, MP_SCOTCH_GPS, MP_SCOTCH_GR_GPS } MPScotchGlobalType;
1363be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetArch(MatPartitioning,const char*);
1364be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetMultilevel(MatPartitioning);
1365be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetGlobal(MatPartitioning,MPScotchGlobalType);
1366be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetCoarseLevel(MatPartitioning,PetscReal);
1367be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetHostList(MatPartitioning,const char*);
136871306c60Sjroman typedef enum { MP_SCOTCH_KERNIGHAN_LIN, MP_SCOTCH_NONE } MPScotchLocalType;
1369be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetLocal(MatPartitioning,MPScotchLocalType);
1370be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetMapping(MatPartitioning);
1371be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetStrategy(MatPartitioning,char*);
137271306c60Sjroman 
13730752156aSBarry Smith /*
13740a835dfdSSatish Balay     If you add entries here you must also add them to finclude/petscmat.h
1375d4fbbf0eSBarry Smith */
13761c1c02c0SLois Curfman McInnes typedef enum { MATOP_SET_VALUES=0,
13771c1c02c0SLois Curfman McInnes                MATOP_GET_ROW=1,
13781c1c02c0SLois Curfman McInnes                MATOP_RESTORE_ROW=2,
13791c1c02c0SLois Curfman McInnes                MATOP_MULT=3,
13801c1c02c0SLois Curfman McInnes                MATOP_MULT_ADD=4,
13817c922b88SBarry Smith                MATOP_MULT_TRANSPOSE=5,
13827c922b88SBarry Smith                MATOP_MULT_TRANSPOSE_ADD=6,
13831c1c02c0SLois Curfman McInnes                MATOP_SOLVE=7,
13841c1c02c0SLois Curfman McInnes                MATOP_SOLVE_ADD=8,
13857c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE=9,
13867c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE_ADD=10,
13871c1c02c0SLois Curfman McInnes                MATOP_LUFACTOR=11,
13881c1c02c0SLois Curfman McInnes                MATOP_CHOLESKYFACTOR=12,
13891c1c02c0SLois Curfman McInnes                MATOP_RELAX=13,
13901c1c02c0SLois Curfman McInnes                MATOP_TRANSPOSE=14,
13911c1c02c0SLois Curfman McInnes                MATOP_GETINFO=15,
13921c1c02c0SLois Curfman McInnes                MATOP_EQUAL=16,
13931c1c02c0SLois Curfman McInnes                MATOP_GET_DIAGONAL=17,
13941c1c02c0SLois Curfman McInnes                MATOP_DIAGONAL_SCALE=18,
13951c1c02c0SLois Curfman McInnes                MATOP_NORM=19,
13961c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_BEGIN=20,
13971c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_END=21,
13981c1c02c0SLois Curfman McInnes                MATOP_COMPRESS=22,
13991c1c02c0SLois Curfman McInnes                MATOP_SET_OPTION=23,
14001c1c02c0SLois Curfman McInnes                MATOP_ZERO_ENTRIES=24,
14011c1c02c0SLois Curfman McInnes                MATOP_ZERO_ROWS=25,
14021c1c02c0SLois Curfman McInnes                MATOP_LUFACTOR_SYMBOLIC=26,
14031c1c02c0SLois Curfman McInnes                MATOP_LUFACTOR_NUMERIC=27,
14041c1c02c0SLois Curfman McInnes                MATOP_CHOLESKY_FACTOR_SYMBOLIC=28,
14051c1c02c0SLois Curfman McInnes                MATOP_CHOLESKY_FACTOR_NUMERIC=29,
1406d643ce63SMatthew Knepley                MATOP_SETUP_PREALLOCATION=30,
1407d643ce63SMatthew Knepley                MATOP_ILUFACTOR_SYMBOLIC=31,
1408d643ce63SMatthew Knepley                MATOP_ICCFACTOR_SYMBOLIC=32,
1409d643ce63SMatthew Knepley                MATOP_GET_ARRAY=33,
1410d643ce63SMatthew Knepley                MATOP_RESTORE_ARRAY=34,
1411e26f1c3cSBarry Smith                MATOP_DUPLICATE=35,
1412d643ce63SMatthew Knepley                MATOP_FORWARD_SOLVE=36,
1413d643ce63SMatthew Knepley                MATOP_BACKWARD_SOLVE=37,
1414d643ce63SMatthew Knepley                MATOP_ILUFACTOR=38,
1415d643ce63SMatthew Knepley                MATOP_ICCFACTOR=39,
1416d643ce63SMatthew Knepley                MATOP_AXPY=40,
1417d643ce63SMatthew Knepley                MATOP_GET_SUBMATRICES=41,
1418d643ce63SMatthew Knepley                MATOP_INCREASE_OVERLAP=42,
1419d643ce63SMatthew Knepley                MATOP_GET_VALUES=43,
1420d643ce63SMatthew Knepley                MATOP_COPY=44,
1421ff1cced0SMatthew Knepley                MATOP_GET_ROW_MAX=45,
1422d643ce63SMatthew Knepley                MATOP_SCALE=46,
1423d643ce63SMatthew Knepley                MATOP_SHIFT=47,
1424d643ce63SMatthew Knepley                MATOP_DIAGONAL_SHIFT=48,
1425d643ce63SMatthew Knepley                MATOP_ILUDT_FACTOR=49,
1426ff1cced0SMatthew Knepley                MATOP_SET_BLOCK_SIZE=50,
1427d643ce63SMatthew Knepley                MATOP_GET_ROW_IJ=51,
1428d643ce63SMatthew Knepley                MATOP_RESTORE_ROW_IJ=52,
1429d643ce63SMatthew Knepley                MATOP_GET_COLUMN_IJ=53,
1430d643ce63SMatthew Knepley                MATOP_RESTORE_COLUMN_IJ=54,
1431d643ce63SMatthew Knepley                MATOP_FDCOLORING_CREATE=55,
1432d643ce63SMatthew Knepley                MATOP_COLORING_PATCH=56,
1433d643ce63SMatthew Knepley                MATOP_SET_UNFACTORED=57,
1434d643ce63SMatthew Knepley                MATOP_PERMUTE=58,
1435d643ce63SMatthew Knepley                MATOP_SET_VALUES_BLOCKED=59,
1436d643ce63SMatthew Knepley                MATOP_GET_SUBMATRIX=60,
1437d643ce63SMatthew Knepley                MATOP_DESTROY=61,
1438d643ce63SMatthew Knepley                MATOP_VIEW=62,
1439ff1cced0SMatthew Knepley                MATOP_CONVERT_FROM=63,
1440d643ce63SMatthew Knepley                MATOP_USE_SCALED_FORM=64,
1441d643ce63SMatthew Knepley                MATOP_SCALE_SYSTEM=65,
1442d643ce63SMatthew Knepley                MATOP_UNSCALE_SYSTEM=66,
1443ba16a8cbSKris Buschelman                MATOP_SET_LOCAL_TO_GLOBAL_MAP=67,
1444d643ce63SMatthew Knepley                MATOP_SET_VALUES_LOCAL=68,
1445d643ce63SMatthew Knepley                MATOP_ZERO_ROWS_LOCAL=69,
1446ff1cced0SMatthew Knepley                MATOP_GET_ROW_MAX_ABS=70,
1447d643ce63SMatthew Knepley                MATOP_CONVERT=71,
1448d643ce63SMatthew Knepley                MATOP_SET_COLORING=72,
1449d643ce63SMatthew Knepley                MATOP_SET_VALUES_ADIC=73,
1450d643ce63SMatthew Knepley                MATOP_SET_VALUES_ADIFOR=74,
1451d643ce63SMatthew Knepley                MATOP_FD_COLORING_APPLY=75,
1452d643ce63SMatthew Knepley                MATOP_SET_FROM_OPTIONS=76,
1453ba16a8cbSKris Buschelman                MATOP_MULT_CON=77,
1454ba16a8cbSKris Buschelman                MATOP_MULT_TRANSPOSE_CON=78,
1455ba16a8cbSKris Buschelman                MATOP_ILU_FACTOR_SYMBOLIC_CON=79,
1456d643ce63SMatthew Knepley                MATOP_PERMUTE_SPARSIFY=80,
1457d643ce63SMatthew Knepley                MATOP_MULT_MULTIPLE=81,
145841acf15aSKris Buschelman                MATOP_SOLVE_MULTIPLE=82,
145941acf15aSKris Buschelman                MATOP_GET_INERTIA=83,
146041acf15aSKris Buschelman                MATOP_LOAD=84,
146141acf15aSKris Buschelman                MATOP_IS_SYMMETRIC=85,
146241acf15aSKris Buschelman                MATOP_IS_HERMITIAN=86,
146341acf15aSKris Buschelman                MATOP_IS_STRUCTURALLY_SYMMETRIC=87,
146441acf15aSKris Buschelman                MATOP_PB_RELAX=88,
146541acf15aSKris Buschelman                MATOP_GET_VECS=89,
146641acf15aSKris Buschelman                MATOP_MAT_MULT=90,
146741acf15aSKris Buschelman                MATOP_MAT_MULT_SYMBOLIC=91,
146841acf15aSKris Buschelman                MATOP_MAT_MULT_NUMERIC=92,
146941acf15aSKris Buschelman                MATOP_PTAP=93,
147041acf15aSKris Buschelman                MATOP_PTAP_SYMBOLIC=94,
1471bc011b1eSHong Zhang                MATOP_PTAP_NUMERIC=95,
1472bc011b1eSHong Zhang                MATOP_MAT_MULTTRANSPOSE=96,
1473bc011b1eSHong Zhang                MATOP_MAT_MULTTRANSPOSE_SYMBOLIC=97,
14747ba1a0bfSKris Buschelman                MATOP_MAT_MULTTRANSPOSE_NUMERIC=98,
14757ba1a0bfSKris Buschelman                MATOP_PTAP_SYMBOLIC_SEQAIJ=99,
14767ba1a0bfSKris Buschelman                MATOP_PTAP_NUMERIC_SEQAIJ=100,
14777ba1a0bfSKris Buschelman                MATOP_PTAP_SYMBOLIC_MPIAIJ=101,
147887d4246cSBarry Smith                MATOP_PTAP_NUMERIC_MPIAIJ=102,
1479ff1cced0SMatthew Knepley                MATOP_CONJUGATE=103,
1480ff1cced0SMatthew Knepley                MATOP_SET_SIZES=104,
1481f5edf698SHong Zhang                MATOP_SET_VALUES_ROW = 105,
1482ff1cced0SMatthew Knepley                MATOP_REAL_PART=106,
1483ff1cced0SMatthew Knepley                MATOP_IMAG_PART=107,
1484d5c63f83SSatish Balay                MATOP_GET_ROW_UTRIANGULAR=108,
14852bebee5dSHong Zhang                MATOP_RESTORE_ROW_UTRIANGULAR=109,
148669db28dcSHong Zhang                MATOP_MATSOLVE=110,
1487829201f2SHong Zhang                MATOP_GET_REDUNDANTMATRIX=111,
1488ff1cced0SMatthew Knepley                MATOP_GET_ROW_MIN=112,
1489ff1cced0SMatthew Knepley                MATOP_GET_COLUMN_VEC=113,
1490ff1cced0SMatthew Knepley                MATOP_MISSING_DIAGONAL=114,
1491ff1cced0SMatthew Knepley                MATOP_MATGETSEQNONZEROSTRUCTURE=115,
1492ff1cced0SMatthew Knepley                MATOP_DESTROY_SOLVER=116
1493fae171e0SBarry Smith              } MatOperation;
1494be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatHasOperation(Mat,MatOperation,PetscTruth*);
1495be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellSetOperation(Mat,MatOperation,void(*)(void));
1496be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellGetOperation(Mat,MatOperation,void(**)(void));
1497be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellSetContext(Mat,void*);
1498112a2221SBarry Smith 
149990ace30eSBarry Smith /*
150090ace30eSBarry Smith    Codes for matrices stored on disk. By default they are
150190ace30eSBarry Smith  stored in a universal format. By changing the format with
1502fb9695e5SSatish Balay  PetscViewerSetFormat(viewer,PETSC_VIEWER_BINARY_NATIVE); the matrices will
150390ace30eSBarry Smith  be stored in a way natural for the matrix, for example dense matrices
150490ace30eSBarry Smith  would be stored as dense. Matrices stored this way may only be
150590ace30eSBarry Smith  read into matrices of the same time.
150690ace30eSBarry Smith */
150790ace30eSBarry Smith #define MATRIX_BINARY_FORMAT_DENSE -1
150890ace30eSBarry Smith 
1509be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetHashTableFactor(Mat,PetscReal);
1510be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatISGetLocalMat(Mat,Mat*);
15111f4e1ec7SBarry Smith 
1512d9274352SBarry Smith /*S
1513d9274352SBarry Smith      MatNullSpace - Object that removes a null space from a vector, i.e.
1514d9274352SBarry Smith          orthogonalizes the vector to a subsapce
1515d9274352SBarry Smith 
1516f7a9e4ceSBarry Smith    Level: advanced
1517d9274352SBarry Smith 
1518d9274352SBarry Smith   Concepts: matrix; linear operator, null space
1519d9274352SBarry Smith 
15206e1639daSBarry Smith   Users manual sections:
15216e1639daSBarry Smith .   sec_singular
15226e1639daSBarry Smith 
1523d9274352SBarry Smith .seealso:  MatNullSpaceCreate()
1524d9274352SBarry Smith S*/
152574637425SBarry Smith typedef struct _p_MatNullSpace* MatNullSpace;
1526d9274352SBarry Smith 
1527be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceCreate(MPI_Comm,PetscTruth,PetscInt,const Vec[],MatNullSpace*);
1528281d1b2eSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(Vec,void*),void*);
1529be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceDestroy(MatNullSpace);
1530be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceRemove(MatNullSpace,Vec,Vec*);
1531be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceAttach(Mat,MatNullSpace);
1532be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceTest(MatNullSpace,Mat);
153374637425SBarry Smith 
1534be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatReorderingSeqSBAIJ(Mat,IS);
1535be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPISBAIJSetHashTableFactor(Mat,PetscReal);
1536be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetColumnIndices(Mat,PetscInt *);
153746129b97SKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJInvertBlockDiagonal(Mat);
15383f1d51d7SBarry Smith 
1539be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMAIJ(Mat,PetscInt,Mat*);
1540be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMAIJRedimension(Mat,PetscInt,Mat*);
1541be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMAIJGetAIJ(Mat,Mat*);
1542c4f061fbSSatish Balay 
1543be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatComputeExplicitOperator(Mat,Mat*);
1544b0a32e0cSBarry Smith 
1545be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalScaleLocal(Mat,Vec);
154604f1ad80SBarry Smith 
1547fef1beadSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*);
15483ec795f1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetBase(Mat,Vec,Vec);
15493ec795f1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*);
1550e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*));
1551e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec));
1552e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDAddNullSpace(Mat,MatNullSpace);
1553e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt);
1554e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDResetHHistory(Mat);
1555e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctionError(Mat,PetscReal);
1556e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetPeriod(Mat,PetscInt);
1557e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDGetH(Mat,PetscScalar *);
1558e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFromOptions(Mat);
1559e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*);
1560e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*);
1561e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDDSSetUmin(Mat,PetscReal);
1562e884886fSBarry Smith 
1563e884886fSBarry Smith 
1564e884886fSBarry Smith typedef struct _p_MatMFFD* MatMFFD;
1565e884886fSBarry Smith 
1566e884886fSBarry Smith /*E
1567e884886fSBarry Smith     MatMFFDType - algorithm used to compute the h used in computing matrix-vector products via differencing of the function
1568e884886fSBarry Smith 
1569e884886fSBarry Smith    Level: beginner
1570e884886fSBarry Smith 
1571e884886fSBarry Smith .seealso: MatMFFDSetType(), MatMFFDRegister()
1572e884886fSBarry Smith E*/
1573a313700dSBarry Smith #define MatMFFDType char*
1574e884886fSBarry Smith #define MATMFFD_DS  "ds"
1575e884886fSBarry Smith #define MATMFFD_WP  "wp"
1576e884886fSBarry Smith 
1577a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetType(Mat,const MatMFFDType);
1578e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatMFFD));
1579e884886fSBarry Smith 
1580e884886fSBarry Smith /*MC
1581e884886fSBarry Smith    MatMFFDRegisterDynamic - Adds a method to the MatMFFD registry.
1582e884886fSBarry Smith 
1583e884886fSBarry Smith    Synopsis:
1584e884886fSBarry Smith    PetscErrorCode MatMFFDRegisterDynamic(char *name_solver,char *path,char *name_create,PetscErrorCode (*routine_create)(MatMFFD))
1585e884886fSBarry Smith 
1586e884886fSBarry Smith    Not Collective
1587e884886fSBarry Smith 
1588e884886fSBarry Smith    Input Parameters:
1589e884886fSBarry Smith +  name_solver - name of a new user-defined compute-h module
1590e884886fSBarry Smith .  path - path (either absolute or relative) the library containing this solver
1591e884886fSBarry Smith .  name_create - name of routine to create method context
1592e884886fSBarry Smith -  routine_create - routine to create method context
1593e884886fSBarry Smith 
1594e884886fSBarry Smith    Level: developer
1595e884886fSBarry Smith 
1596e884886fSBarry Smith    Notes:
1597e884886fSBarry Smith    MatMFFDRegisterDynamic() may be called multiple times to add several user-defined solvers.
1598e884886fSBarry Smith 
1599e884886fSBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
1600e884886fSBarry Smith    is ignored.
1601e884886fSBarry Smith 
1602e884886fSBarry Smith    Sample usage:
1603e884886fSBarry Smith .vb
1604e884886fSBarry Smith    MatMFFDRegisterDynamic("my_h",/home/username/my_lib/lib/libO/solaris/mylib.a,
1605e884886fSBarry Smith                "MyHCreate",MyHCreate);
1606e884886fSBarry Smith .ve
1607e884886fSBarry Smith 
1608e884886fSBarry Smith    Then, your solver can be chosen with the procedural interface via
1609e884886fSBarry Smith $     MatMFFDSetType(mfctx,"my_h")
1610e884886fSBarry Smith    or at runtime via the option
1611e884886fSBarry Smith $     -snes_mf_type my_h
1612e884886fSBarry Smith 
1613e884886fSBarry Smith .keywords: MatMFFD, register
1614e884886fSBarry Smith 
1615e884886fSBarry Smith .seealso: MatMFFDRegisterAll(), MatMFFDRegisterDestroy()
1616e884886fSBarry Smith M*/
1617e884886fSBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1618e884886fSBarry Smith #define MatMFFDRegisterDynamic(a,b,c,d) MatMFFDRegister(a,b,c,0)
1619e884886fSBarry Smith #else
1620e884886fSBarry Smith #define MatMFFDRegisterDynamic(a,b,c,d) MatMFFDRegister(a,b,c,d)
1621e884886fSBarry Smith #endif
1622e884886fSBarry Smith 
1623e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDRegisterAll(const char[]);
1624e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDRegisterDestroy(void);
1625e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDDefaultSetUmin(Mat,PetscReal);
1626e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDWPSetComputeNormU(Mat,PetscTruth);
1627e884886fSBarry Smith 
1628e884886fSBarry Smith 
1629be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *);
1630be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *);
16317dbadf16SMatthew Knepley 
1632e9fa29b7SSatish Balay PETSC_EXTERN_CXX_END
16332eac72dbSBarry Smith #endif
1634