xref: /petsc/include/petscmat.h (revision ab970d9bea782d29202c9bc5dbfd4088487bf7e5)
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 
27d9274352SBarry Smith .seealso: MatSetType(), Mat
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"
54b3a1e11cSKris Buschelman #define MATSEQAIJSPOOLES   "seqaijspooles"
55d10c748bSKris Buschelman #define MATMPIAIJSPOOLES   "mpiaijspooles"
569abb65ffSKris Buschelman #define MATSEQSBAIJSPOOLES "seqsbaijspooles"
5722191285SKris Buschelman #define MATMPISBAIJSPOOLES "mpisbaijspooles"
58bb4d4166SHong Zhang #define MATAIJSPOOLES      "aijspooles"
59bb4d4166SHong Zhang #define MATSBAIJSPOOLES    "sbaijspooles"
6014b396bbSKris Buschelman #define MATSUPERLU         "superlu"
611d96aa28SKris Buschelman #define MATSUPERLU_DIST    "superlu_dist"
621677d0b8SKris Buschelman #define MATUMFPACK         "umfpack"
63e8aa55a4SKris Buschelman #define MATESSL            "essl"
644eb8e494SKris Buschelman #define MATLUSOL           "lusol"
65397b6df1SKris Buschelman #define MATAIJMUMPS        "aijmumps"
66397b6df1SKris Buschelman #define MATSBAIJMUMPS      "sbaijmumps"
678da957c5SKris Buschelman #define MATDSCPACK         "dscpack"
687065e2aaSKris Buschelman #define MATMATLAB          "matlab"
694b8d99adSRichard Tran Mills #define MATSEQCSRPERM      "seqcsrperm"
7018def556SRichard Tran Mills #define MATMPICSRPERM      "mpicsrperm"
71814655b8SBarry Smith #define MATCSRPERM         "csrperm"
7281824310SBarry Smith #define MATSEQCRL          "seqcrl"
73c02b24c5SBarry Smith #define MATMPICRL          "mpicrl"
74c02b24c5SBarry Smith #define MATCRL             "crl"
752a6744ebSBarry Smith #define MATSCATTER         "scatter"
76421e10b8SBarry Smith #define MATBLOCKMAT        "blockmat"
77793850ffSBarry Smith #define MATCOMPOSITE       "composite"
785a7f1df3SHong Zhang #define MATSEQFFTW         "seqfftw"
79d91e6319SBarry Smith 
809989ab13SBarry Smith /*E
819989ab13SBarry Smith     MatSolverType - String with the name of a PETSc matrix solver type.
829989ab13SBarry Smith 
839989ab13SBarry Smith     For example: "petsc" indicates what PETSc provides, "superlu" indicates either
849989ab13SBarry Smith        SuperLU or SuperLU_Dist etc.
859989ab13SBarry Smith 
869989ab13SBarry Smith 
879989ab13SBarry Smith    Level: beginner
889989ab13SBarry Smith 
899989ab13SBarry Smith .seealso: MatSetSolverType(), Mat, MatSetType(), MatType
909989ab13SBarry Smith E*/
919989ab13SBarry Smith #define MatSolverType char*
929989ab13SBarry Smith 
93c06d978dSMatthew Knepley /* Logging support */
94552e946dSBarry Smith #define    MAT_FILE_COOKIE 1211216    /* used to indicate matrices in binary files */
95be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_COOKIE;
96be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_FDCOLORING_COOKIE;
97be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_PARTITIONING_COOKIE;
98be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_NULLSPACE_COOKIE;
99e884886fSBarry Smith extern PetscCookie PETSCMAT_DLLEXPORT MATMFFD_COOKIE;
100c06d978dSMatthew Knepley 
101ceb03754SKris Buschelman /*E
102ceb03754SKris Buschelman     MatReuse - Indicates if matrices obtained from a previous call to MatGetSubMatrices()
103ceb03754SKris Buschelman      or MatGetSubMatrix() are to be reused to store the new matrix values.
104ceb03754SKris Buschelman 
105ceb03754SKris Buschelman     Level: beginner
106ceb03754SKris Buschelman 
107ceb03754SKris Buschelman    Any additions/changes here MUST also be made in include/finclude/petscmat.h
108ceb03754SKris Buschelman 
1090c451bc4SBarry Smith .seealso: MatGetSubMatrices(), MatGetSubMatrix(), MatDestroyMatrices(), MatConvert()
110ceb03754SKris Buschelman E*/
111ceb03754SKris Buschelman typedef enum {MAT_INITIAL_MATRIX,MAT_REUSE_MATRIX} MatReuse;
112ceb03754SKris Buschelman 
1135494a064SHong Zhang /*E
1145494a064SHong Zhang     MatGetSubMatrixOption - Indicates if matrices obtained from a call to MatGetSubMatrices()
115829201f2SHong Zhang      include the matrix values. Currently it is only used by MatGetSeqNonzerostructure().
1165494a064SHong Zhang 
1175494a064SHong Zhang     Level: beginner
1185494a064SHong Zhang 
119829201f2SHong Zhang .seealso: MatGetSeqNonzerostructure()
1205494a064SHong Zhang E*/
1215494a064SHong Zhang typedef enum {MAT_DO_NOT_GET_VALUES,MAT_GET_VALUES} MatGetSubMatrixOption;
1225494a064SHong Zhang 
123e5bd5246SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInitializePackage(const char[]);
124c06d978dSMatthew Knepley 
125f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreate(MPI_Comm,Mat*);
126a69a960dSMatthew Knepley PetscPolymorphicFunction(MatCreate,(MPI_Comm comm),(comm,&A),Mat,A)
127a69a960dSMatthew Knepley PetscPolymorphicFunction(MatCreate,(),(PETSC_COMM_WORLD,&A),Mat,A)
128f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetSizes(Mat,PetscInt,PetscInt,PetscInt,PetscInt);
129a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetType(Mat,const MatType);
1309989ab13SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetSolverType(Mat,MatSolverType);
1319989ab13SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSolverType(Mat,MatSolverType*);
132be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetFromOptions(Mat);
133be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUpPreallocation(Mat);
134be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRegisterAll(const char[]);
135be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat));
13630de9b25SBarry Smith 
137f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetOptionsPrefix(Mat,const char[]);
138f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAppendOptionsPrefix(Mat,const char[]);
139f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOptionsPrefix(Mat,const char*[]);
140f69a0ea3SMatthew Knepley 
14130de9b25SBarry Smith /*MC
14230de9b25SBarry Smith    MatRegisterDynamic - Adds a new matrix type
14330de9b25SBarry Smith 
14430de9b25SBarry Smith    Synopsis:
145c1ac3661SBarry Smith    PetscErrorCode MatRegisterDynamic(char *name,char *path,char *name_create,PetscErrorCode (*routine_create)(Mat))
14630de9b25SBarry Smith 
14730de9b25SBarry Smith    Not Collective
14830de9b25SBarry Smith 
14930de9b25SBarry Smith    Input Parameters:
15030de9b25SBarry Smith +  name - name of a new user-defined matrix type
15130de9b25SBarry Smith .  path - path (either absolute or relative) the library containing this solver
15230de9b25SBarry Smith .  name_create - name of routine to create method context
15330de9b25SBarry Smith -  routine_create - routine to create method context
15430de9b25SBarry Smith 
15530de9b25SBarry Smith    Notes:
15630de9b25SBarry Smith    MatRegisterDynamic() may be called multiple times to add several user-defined solvers.
15730de9b25SBarry Smith 
15830de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
15930de9b25SBarry Smith    is ignored.
16030de9b25SBarry Smith 
16130de9b25SBarry Smith    Sample usage:
16230de9b25SBarry Smith .vb
16330de9b25SBarry Smith    MatRegisterDynamic("my_mat",/home/username/my_lib/lib/libO/solaris/mylib.a,
16430de9b25SBarry Smith                "MyMatCreate",MyMatCreate);
16530de9b25SBarry Smith .ve
16630de9b25SBarry Smith 
16730de9b25SBarry Smith    Then, your solver can be chosen with the procedural interface via
16830de9b25SBarry Smith $     MatSetType(Mat,"my_mat")
16930de9b25SBarry Smith    or at runtime via the option
17030de9b25SBarry Smith $     -mat_type my_mat
17130de9b25SBarry Smith 
17230de9b25SBarry Smith    Level: advanced
17330de9b25SBarry Smith 
174ab901514SBarry Smith    Notes: ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values.
17530de9b25SBarry Smith          If your function is not being put into a shared library then use VecRegister() instead
17630de9b25SBarry Smith 
17730de9b25SBarry Smith .keywords: Mat, register
17830de9b25SBarry Smith 
17930de9b25SBarry Smith .seealso: MatRegisterAll(), MatRegisterDestroy()
18030de9b25SBarry Smith 
18130de9b25SBarry Smith M*/
182273d9f13SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
183273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,0)
184273d9f13SBarry Smith #else
185273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,d)
18630de9b25SBarry Smith #endif
18730de9b25SBarry Smith 
1889c666560SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRegisterAll(const char[]);
189273d9f13SBarry Smith extern PetscTruth MatRegisterAllCalled;
190b0a32e0cSBarry Smith extern PetscFList MatList;
19128988994SBarry Smith 
192be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqDense(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*);
193be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIDense(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*);
194be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
195ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,nz,nnz,&A),Mat,A)
196ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,&A),Mat,A)
197ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,0,nnz,&A),Mat,A)
198ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,&A),Mat,A)
199ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,A))
200ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,m,n,0,nnz,A))
201ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,A))
202be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
203ba966639SSatish 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)
204ba966639SSatish 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)
205ba966639SSatish 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)
206ba966639SSatish 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)
207ba966639SSatish 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))
208ba966639SSatish 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))
209ba966639SSatish 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))
210ba966639SSatish 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)
211ba966639SSatish 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)
212ba966639SSatish 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)
213ba966639SSatish 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)
214ba966639SSatish 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))
215ba966639SSatish 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))
216ba966639SSatish 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))
21763ab93bfSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *);
2188d7a6e47SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJWithSplitArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],PetscInt[],PetscInt[],PetscScalar[],Mat*);
2198d7a6e47SBarry Smith 
220be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIRowbs(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
221be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
222ba966639SSatish 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)
223ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A)
224ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A)
225ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A)
226ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A))
227ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A))
228ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A))
229be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
230ba966639SSatish 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)
231ba966639SSatish 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)
232ba966639SSatish 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)
233ba966639SSatish 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)
234ba966639SSatish 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))
235ba966639SSatish 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))
236ba966639SSatish 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))
237ba966639SSatish 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)
238ba966639SSatish 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)
239ba966639SSatish 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)
240ba966639SSatish 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)
241ba966639SSatish 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))
242ba966639SSatish 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))
243ba966639SSatish 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))
244be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAdj(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscInt[],Mat*);
245be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
246ba966639SSatish 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)
247ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A)
248ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A)
249ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A)
250ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A))
251ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A))
252ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A))
253be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPISBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
254ba966639SSatish 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)
255ba966639SSatish 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)
256ba966639SSatish 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)
257ba966639SSatish 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)
258ba966639SSatish 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))
259ba966639SSatish 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))
260ba966639SSatish 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))
261ba966639SSatish 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)
262ba966639SSatish 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)
263ba966639SSatish 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)
264ba966639SSatish 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)
265ba966639SSatish 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))
266ba966639SSatish 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))
267ba966639SSatish 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))
268be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateShell(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,void *,Mat*);
269ba966639SSatish 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)
270ba966639SSatish Balay PetscPolymorphicFunction(MatCreateShell,(PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx),(PETSC_COMM_WORLD,m,n,M,N,ctx,&A),Mat,A)
271be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateAdic(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,void (*)(void),Mat*);
272be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateNormal(Mat,Mat*);
273ba966639SSatish Balay PetscPolymorphicFunction(MatCreateNormal,(Mat mat),(mat,&A),Mat,A)
274ab92ecdeSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateLRC(Mat,Mat,Mat,Mat*);
275284134d9SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateIS(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,ISLocalToGlobalMapping,Mat*);
2761472f72bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
2771472f72bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPICRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
2782a6744ebSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateScatter(MPI_Comm,VecScatter,Mat*);
2792a6744ebSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScatterSetVecScatter(Mat,VecScatter);
2802a6744ebSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScatterGetVecScatter(Mat,VecScatter*);
2818cdf2d9bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateBlockMat(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,Mat*);
282793850ffSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompositeAddMat(Mat,Mat);
283b52f573bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompositeMerge(Mat);
284793850ffSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateComposite(MPI_Comm,PetscInt,const Mat*,Mat*);
2855a7f1df3SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqFFTW(MPI_Comm,PetscInt,const PetscInt[],Mat*);
2861472f72bSBarry Smith 
287f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUp(Mat);
288be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroy(Mat);
28921c89e3eSBarry Smith 
2900c451bc4SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConjugate(Mat);
29199cafbc1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRealPart(Mat);
29299cafbc1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatImaginaryPart(Mat);
29399cafbc1SBarry Smith 
2948ed539a5SBarry Smith /* ------------------------------------------------------------*/
295be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
296be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
29787d4246cSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesRow(Mat,PetscInt,const PetscScalar[]);
298f259bd47SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesRowLocal(Mat,PetscInt,const PetscScalar[]);
29984cb2905SBarry Smith 
3002ef4de8bSBarry Smith /*S
3012ef4de8bSBarry Smith      MatStencil - Data structure (C struct) for storing information about a single row or
3022ef4de8bSBarry Smith         column of a matrix as index on an associated grid.
3032ef4de8bSBarry Smith 
3042ef4de8bSBarry Smith    Level: beginner
3052ef4de8bSBarry Smith 
3062ef4de8bSBarry Smith   Concepts: matrix; linear operator
3072ef4de8bSBarry Smith 
308d7bd93baSBarry Smith .seealso:  MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockStencil()
3092ef4de8bSBarry Smith S*/
310435da068SBarry Smith typedef struct {
311c1ac3661SBarry Smith   PetscInt k,j,i,c;
312435da068SBarry Smith } MatStencil;
3132ef4de8bSBarry Smith 
314be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
315be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
316be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt);
317435da068SBarry Smith 
318be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetColoring(Mat,ISColoring);
319be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesAdic(Mat,void*);
320be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesAdifor(Mat,PetscInt,void*);
3213a7fca6bSBarry Smith 
322d91e6319SBarry Smith /*E
323d91e6319SBarry Smith     MatAssemblyType - Indicates if the matrix is now to be used, or if you plan
324d91e6319SBarry Smith      to continue to add values to it
325d91e6319SBarry Smith 
326d91e6319SBarry Smith     Level: beginner
327d91e6319SBarry Smith 
328d91e6319SBarry Smith .seealso: MatAssemblyBegin(), MatAssemblyEnd()
329d91e6319SBarry Smith E*/
3306d4a8577SBarry Smith typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType;
331be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssemblyBegin(Mat,MatAssemblyType);
332be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssemblyEnd(Mat,MatAssemblyType);
333be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssembled(Mat,PetscTruth*);
3344f9c727eSBarry Smith 
3351d73ed98SBarry Smith 
33630de9b25SBarry Smith 
337d91e6319SBarry Smith /*E
338d91e6319SBarry Smith     MatOption - Options that may be set for a matrix and its behavior or storage
339d91e6319SBarry Smith 
340d91e6319SBarry Smith     Level: beginner
341d91e6319SBarry Smith 
3420a835dfdSSatish Balay    Any additions/changes here MUST also be made in include/finclude/petscmat.h
343d91e6319SBarry Smith 
344d91e6319SBarry Smith .seealso: MatSetOption()
345d91e6319SBarry Smith E*/
346512a5fc5SBarry Smith typedef enum {MAT_ROW_ORIENTED,MAT_NEW_NONZERO_LOCATIONS,
3474e0d8c25SBarry Smith               MAT_SYMMETRIC,
3484e0d8c25SBarry Smith               MAT_STRUCTURALLY_SYMMETRIC,
3498047e77bSBarry Smith               MAT_NEW_DIAGONALS,MAT_IGNORE_OFF_PROC_ENTRIES,
3504e0d8c25SBarry Smith               MAT_NEW_NONZERO_LOCATION_ERR,
3514e0d8c25SBarry Smith               MAT_NEW_NONZERO_ALLOCATION_ERR,MAT_USE_HASH_TABLE,
3524e0d8c25SBarry Smith               MAT_KEEP_ZEROED_ROWS,MAT_IGNORE_ZERO_ENTRIES,
3534e0d8c25SBarry Smith               MAT_USE_INODES,
3544e0d8c25SBarry Smith               MAT_HERMITIAN,
3554e0d8c25SBarry Smith               MAT_SYMMETRY_ETERNAL,
3564e0d8c25SBarry Smith               MAT_USE_COMPRESSEDROW,
3574e0d8c25SBarry Smith               MAT_IGNORE_LOWER_TRIANGULAR,MAT_ERROR_LOWER_TRIANGULAR,
3584e0d8c25SBarry Smith               MAT_GETROW_UPPERTRIANGULAR} MatOption;
359290bbb0aSBarry Smith extern const char *MatOptions[];
3604e0d8c25SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetOption(Mat,MatOption,PetscTruth);
361a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetType(Mat,const MatType*);
362a313700dSBarry Smith PetscPolymorphicFunction(MatGetType,(Mat mat),(mat,&t),const MatType,t)
36384cb2905SBarry Smith 
364be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]);
365be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
366be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
367f5edf698SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowUpperTriangular(Mat);
368f5edf698SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRowUpperTriangular(Mat);
369be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
370be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
371be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumnVector(Mat,Vec,PetscInt);
372be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetArray(Mat,PetscScalar *[]);
373ba966639SSatish Balay PetscPolymorphicFunction(MatGetArray,(Mat mat),(mat,&a),PetscScalar*,a)
374be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreArray(Mat,PetscScalar *[]);
375be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBlockSize(Mat,PetscInt *);
376ba966639SSatish Balay PetscPolymorphicFunction(MatGetBlockSize,(Mat mat),(mat,&a),PetscInt,a)
377be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetBlockSize(Mat,PetscInt);
3787b80b807SBarry Smith 
3791620fd73SBarry Smith 
380be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMult(Mat,Vec,Vec);
381be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultAdd(Mat,Vec,Vec,Vec);
382ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
383be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTranspose(Mat,Vec,Vec);
384be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsTranspose(Mat,Mat,PetscReal,PetscTruth*);
385ba966639SSatish Balay PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B,PetscReal tol),(A,B,tol,&t),PetscTruth,t)
386e37f64cdSBarry Smith PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B),(A,B,0,&t),PetscTruth,t)
3871cbb95d3SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitianTranspose(Mat,Mat,PetscReal,PetscTruth*);
388be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeAdd(Mat,Vec,Vec,Vec);
389ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultTransposeAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
390be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultConstrained(Mat,Vec,Vec);
391be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeConstrained(Mat,Vec,Vec);
3922bebee5dSHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatSolve(Mat,Mat,Mat);
393f5edf698SHong Zhang 
394d91e6319SBarry Smith /*E
395d91e6319SBarry Smith     MatDuplicateOption - Indicates if a duplicated sparse matrix should have
396d91e6319SBarry Smith   its numerical values copied over or just its nonzero structure.
397d91e6319SBarry Smith 
398d91e6319SBarry Smith     Level: beginner
399d91e6319SBarry Smith 
400d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
401d91e6319SBarry Smith 
402d91e6319SBarry Smith .seealso: MatDuplicate()
403d91e6319SBarry Smith E*/
4042e8a6d31SBarry Smith typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES} MatDuplicateOption;
4052e8a6d31SBarry Smith 
406a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert(Mat,const MatType,MatReuse,Mat*);
407a313700dSBarry Smith PetscPolymorphicFunction(MatConvert,(Mat A,const MatType t),(A,t,MAT_INITIAL_MATRIX,&a),Mat,a)
408be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDuplicate(Mat,MatDuplicateOption,Mat*);
409ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A,MatDuplicateOption o),(A,o,&a),Mat,a)
410ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A),(A,MAT_COPY_VALUES,&a),Mat,a)
41194a9d846SBarry Smith 
412d91e6319SBarry Smith /*E
413d91e6319SBarry Smith     MatStructure - Indicates if the matrix has the same nonzero structure
414d91e6319SBarry Smith 
415d91e6319SBarry Smith     Level: beginner
416d91e6319SBarry Smith 
417d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
418d91e6319SBarry Smith 
41994b7f48cSBarry Smith .seealso: MatCopy(), KSPSetOperators(), PCSetOperators()
420d91e6319SBarry Smith E*/
421c537a176SHong Zhang typedef enum {SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER,SUBSET_NONZERO_PATTERN} MatStructure;
422cb5b572fSBarry Smith 
423be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCopy(Mat,Mat,MatStructure);
424be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatView(Mat,PetscViewer);
425be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsSymmetric(Mat,PetscReal,PetscTruth*);
426ba966639SSatish Balay PetscPolymorphicFunction(MatIsSymmetric,(Mat A,PetscReal tol),(A,tol,&t),PetscTruth,t)
427e37f64cdSBarry Smith PetscPolymorphicFunction(MatIsSymmetric,(Mat A),(A,0,&t),PetscTruth,t)
428be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsStructurallySymmetric(Mat,PetscTruth*);
429ba966639SSatish Balay PetscPolymorphicFunction(MatIsStructurallySymmetric,(Mat A),(A,&t),PetscTruth,t)
4301cbb95d3SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitian(Mat,PetscReal,PetscTruth*);
4311cbb95d3SBarry Smith PetscPolymorphicFunction(MatIsHermitian,(Mat A),(A,0,&t),PetscTruth,t)
432be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsSymmetricKnown(Mat,PetscTruth*,PetscTruth*);
433be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitianKnown(Mat,PetscTruth*,PetscTruth*);
43464c62002SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMissingDiagonal(Mat,PetscTruth *,PetscInt *);
435a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLoad(PetscViewer,const MatType,Mat*);
436a313700dSBarry Smith PetscPolymorphicFunction(MatLoad,(PetscViewer v,const MatType t),(v,t,&a),Mat,a)
4377b80b807SBarry Smith 
4388f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *);
4398f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRowIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *);
4408f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumnIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *);
4418f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreColumnIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *);
442d4fbbf0eSBarry Smith 
443d91e6319SBarry Smith /*S
444d91e6319SBarry Smith      MatInfo - Context of matrix information, used with MatGetInfo()
445d91e6319SBarry Smith 
446d91e6319SBarry Smith    In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE
447d91e6319SBarry Smith 
448d91e6319SBarry Smith    Level: intermediate
449d91e6319SBarry Smith 
450d91e6319SBarry Smith   Concepts: matrix^nonzero information
451d91e6319SBarry Smith 
452d9274352SBarry Smith .seealso:  MatGetInfo(), MatInfoType
453d91e6319SBarry Smith S*/
4544e220ebcSLois Curfman McInnes typedef struct {
455b0a32e0cSBarry Smith   PetscLogDouble rows_global,columns_global;         /* number of global rows and columns */
456b0a32e0cSBarry Smith   PetscLogDouble rows_local,columns_local;           /* number of local rows and columns */
457b0a32e0cSBarry Smith   PetscLogDouble block_size;                         /* block size */
458b0a32e0cSBarry Smith   PetscLogDouble nz_allocated,nz_used,nz_unneeded;   /* number of nonzeros */
459b0a32e0cSBarry Smith   PetscLogDouble memory;                             /* memory allocated */
460b0a32e0cSBarry Smith   PetscLogDouble assemblies;                         /* number of matrix assemblies called */
461b0a32e0cSBarry Smith   PetscLogDouble mallocs;                            /* number of mallocs during MatSetValues() */
462b0a32e0cSBarry Smith   PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */
463b0a32e0cSBarry Smith   PetscLogDouble factor_mallocs;                     /* number of mallocs during factorization */
4644e220ebcSLois Curfman McInnes } MatInfo;
4654e220ebcSLois Curfman McInnes 
466d9274352SBarry Smith /*E
467d9274352SBarry Smith     MatInfoType - Indicates if you want information about the local part of the matrix,
468d9274352SBarry Smith      the entire parallel matrix or the maximum over all the local parts.
469d9274352SBarry Smith 
470d9274352SBarry Smith     Level: beginner
471d9274352SBarry Smith 
472d9274352SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
473d9274352SBarry Smith 
474d9274352SBarry Smith .seealso: MatGetInfo(), MatInfo
475d9274352SBarry Smith E*/
4767b80b807SBarry Smith typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType;
477be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetInfo(Mat,MatInfoType,MatInfo*);
478be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatValid(Mat,PetscTruth*);
479be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetDiagonal(Mat,Vec);
480985db425SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMax(Mat,Vec,PetscInt[]);
481985db425SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMin(Mat,Vec,PetscInt[]);
482985db425SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMaxAbs(Mat,Vec,PetscInt[]);
48386697f06SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowSum(Mat,Vec);
484fc4dec0aSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatTranspose(Mat,MatReuse,Mat*);
485fc4dec0aSBarry Smith PetscPolymorphicFunction(MatTranspose,(Mat A),(A,MAT_INITIAL_MATRIX,&t),Mat,t)
486be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPermute(Mat,IS,IS,Mat *);
487ba966639SSatish Balay PetscPolymorphicFunction(MatPermute,(Mat A,IS is1,IS is2),(A,is1,is2,&t),Mat,t)
488be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPermuteSparsify(Mat,PetscInt,PetscReal,PetscReal,IS,IS,Mat *);
489be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalScale(Mat,Vec,Vec);
490be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalSet(Mat,Vec,InsertMode);
491be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatEqual(Mat,Mat,PetscTruth*);
492ba966639SSatish Balay PetscPolymorphicFunction(MatEqual,(Mat A,Mat B),(A,B,&t),PetscTruth,t)
493be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultEqual(Mat,Mat,PetscInt,PetscTruth*);
494be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultAddEqual(Mat,Mat,PetscInt,PetscTruth*);
495be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeEqual(Mat,Mat,PetscInt,PetscTruth*);
496be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscTruth*);
4977b80b807SBarry Smith 
498be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNorm(Mat,NormType,PetscReal *);
499ba966639SSatish Balay PetscPolymorphicFunction(MatNorm,(Mat A,NormType t),(A,t,&n),PetscReal,n)
500be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroEntries(Mat);
501f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar);
502f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsIS(Mat,IS,PetscScalar);
503f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroColumns(Mat,PetscInt,const PetscInt [],const PetscScalar*);
504f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroColumnsIS(Mat,IS,const PetscScalar*);
5057b80b807SBarry Smith 
506be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatUseScaledForm(Mat,PetscTruth);
507be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScaleSystem(Mat,Vec,Vec);
508be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatUnScaleSystem(Mat,Vec,Vec);
5095ef9f2a5SBarry Smith 
510be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSize(Mat,PetscInt*,PetscInt*);
511be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalSize(Mat,PetscInt*,PetscInt*);
512be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRange(Mat,PetscInt*,PetscInt*);
5138ee2e534SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRanges(Mat,const PetscInt**);
514f4b8409dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRangeColumn(Mat,PetscInt*,PetscInt*);
5157b80b807SBarry Smith 
516be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
517be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroyMatrices(PetscInt,Mat *[]);
518be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrix(Mat,IS,IS,PetscInt,MatReuse,Mat *);
519abf0e977SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrixRaw(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscInt,MatReuse,Mat *);
520829201f2SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSeqNonzeroStructure(Mat,Mat *[]);
521829201f2SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroySeqNonzeroStructure(Mat *[]);
5225494a064SHong Zhang 
523be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
524be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPI(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*);
525be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPISymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*);
526be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPINumeric(Mat,Mat);
527be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroy_MPIAIJ_SeqsToMPI(Mat);
528be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalMat(Mat,MatReuse,Mat*);
529be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*);
530be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,PetscInt*,Mat*);
531dd6ea824SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBrowsOfAoCols(Mat,Mat,MatReuse,PetscInt**,MatScalar**,Mat*);
53243eb5e2fSMatthew Knepley #if defined (PETSC_USE_CTABLE)
533cd116e26SSatish Balay #include "petscctable.h"
53443eb5e2fSMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetCommunicationStructs(Mat, Vec *, PetscTable *, VecScatter *);
53543eb5e2fSMatthew Knepley #else
53643eb5e2fSMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetCommunicationStructs(Mat, Vec *, PetscInt *[], VecScatter *);
53743eb5e2fSMatthew Knepley #endif
5388efafbd8SBarry Smith 
539be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt);
5407b80b807SBarry Smith 
541be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*);
542be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultSymbolic(Mat,Mat,PetscReal,Mat*);
543be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultNumeric(Mat,Mat,Mat);
54422440eb1SKris Buschelman 
545be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*);
546be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAPSymbolic(Mat,Mat,PetscReal,Mat*);
547be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAPNumeric(Mat,Mat,Mat);
54822440eb1SKris Buschelman 
549be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTranspose(Mat,Mat,MatReuse,PetscReal,Mat*);
550be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTransposeSymbolic(Mat,Mat,PetscReal,Mat*);
551be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTransposeNumeric(Mat,Mat,Mat);
552bc011b1eSHong Zhang 
553f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAXPY(Mat,PetscScalar,Mat,MatStructure);
55404aac2b0SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAYPX(Mat,PetscScalar,Mat,MatStructure);
555be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompress(Mat);
5561c741599SHong Zhang 
557f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScale(Mat,PetscScalar);
558f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShift(Mat,PetscScalar);
5597b80b807SBarry Smith 
560be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping);
561be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetLocalToGlobalMappingBlock(Mat,ISLocalToGlobalMapping);
562f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar);
563f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsLocalIS(Mat,IS,PetscScalar);
564be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
565be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
566052efed2SBarry Smith 
567be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStashSetInitialSize(Mat,PetscInt,PetscInt);
568be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
56990f02eecSBarry Smith 
570be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInterpolate(Mat,Vec,Vec);
5710c451bc4SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInterpolateAdd(Mat,Vec,Vec,Vec);
572ba966639SSatish Balay PetscPolymorphicSubroutine(MatInterpolateAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
573be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestrict(Mat,Vec,Vec);
574be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetVecs(Mat,Vec*,Vec*);
57569db28dcSHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRedundantMatrix(Mat,PetscInt,MPI_Comm,PetscInt,MatReuse,Mat*);
576bd481603SBarry Smith 
577bd481603SBarry Smith /*MC
5781620fd73SBarry Smith    MatSetValue - Set a single entry into a matrix.
5791620fd73SBarry Smith 
5801620fd73SBarry Smith    Not collective
5811620fd73SBarry Smith 
5821620fd73SBarry Smith    Input Parameters:
5831620fd73SBarry Smith +  m - the matrix
5841620fd73SBarry Smith .  row - the row location of the entry
5851620fd73SBarry Smith .  col - the column location of the entry
5861620fd73SBarry Smith .  value - the value to insert
5871620fd73SBarry Smith -  mode - either INSERT_VALUES or ADD_VALUES
5881620fd73SBarry Smith 
5891620fd73SBarry Smith    Notes:
5901620fd73SBarry Smith    For efficiency one should use MatSetValues() and set several or many
5911620fd73SBarry Smith    values simultaneously if possible.
5921620fd73SBarry Smith 
5931620fd73SBarry Smith    Level: beginner
5941620fd73SBarry Smith 
5951620fd73SBarry Smith .seealso: MatSetValues(), MatSetValueLocal()
5961620fd73SBarry Smith M*/
5971620fd73SBarry 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);}
5981620fd73SBarry Smith 
5991620fd73SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatGetValue(Mat v,PetscInt i,PetscInt j,PetscScalar *va) {return MatGetValues(v,1,&i,1,&j,va);}
6001620fd73SBarry Smith 
6011620fd73SBarry 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);}
6021620fd73SBarry Smith 
6031620fd73SBarry Smith /*MC
6040d2a7ff2SSatish Balay    MatPreallocateInitialize - Begins the block of code that will count the number of nonzeros per
605bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
606bd481603SBarry Smith 
607bd481603SBarry Smith    Synopsis:
608c1ac3661SBarry Smith    PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)
609bd481603SBarry Smith 
610bd481603SBarry Smith    Collective on MPI_Comm
611bd481603SBarry Smith 
612bd481603SBarry Smith    Input Parameters:
613bd481603SBarry Smith +  comm - the communicator that will share the eventually allocated matrix
614859fcb39SBarry Smith .  nrows - the number of LOCAL rows in the matrix
615859fcb39SBarry Smith -  ncols - the number of LOCAL columns in the matrix
616bd481603SBarry Smith 
617bd481603SBarry Smith    Output Parameters:
618bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
619bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
620bd481603SBarry Smith 
621bd481603SBarry Smith 
622bd481603SBarry Smith    Level: intermediate
623bd481603SBarry Smith 
624bd481603SBarry Smith    Notes:
625bd481603SBarry Smith    See the chapter in the users manual on performance for more details
626bd481603SBarry Smith 
6271d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
628bd481603SBarry Smith 
629bd481603SBarry Smith    Use MatPreallocateInitializeSymmetric() for symmetric matrices (MPISBAIJ matrices)
630bd481603SBarry Smith 
6311620fd73SBarry Smith    This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize().
6321620fd73SBarry Smith 
633bd481603SBarry Smith   Concepts: preallocation^Matrix
634bd481603SBarry Smith 
635bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
636bd481603SBarry Smith           MatPreallocateInitializeSymmetric(), MatPreallocateSymmetricSetLocal()
637bd481603SBarry Smith M*/
638c4f061fbSSatish Balay #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \
6397c922b88SBarry Smith { \
640a89bc211SBarry Smith   PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__start,__end; \
641a89bc211SBarry Smith   _4_ierr = PetscMalloc2(__nrows,PetscInt,&dnz,__nrows,PetscInt,&onz);CHKERRQ(_4_ierr); \
642a89bc211SBarry Smith   _4_ierr = PetscMemzero(dnz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
643a89bc211SBarry Smith   _4_ierr = PetscMemzero(onz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
644c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ctmp;\
645a89bc211SBarry Smith   _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows;
6467c922b88SBarry Smith 
647bd481603SBarry Smith /*MC
6480d2a7ff2SSatish Balay    MatPreallocateSymmetricInitialize - Begins the block of code that will count the number of nonzeros per
649bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
650bd481603SBarry Smith 
651bd481603SBarry Smith    Synopsis:
652c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)
653bd481603SBarry Smith 
654bd481603SBarry Smith    Collective on MPI_Comm
655bd481603SBarry Smith 
656bd481603SBarry Smith    Input Parameters:
657bd481603SBarry Smith +  comm - the communicator that will share the eventually allocated matrix
658859fcb39SBarry Smith .  nrows - the number of LOCAL rows in the matrix
659859fcb39SBarry Smith -  ncols - the number of LOCAL columns in the matrix
660bd481603SBarry Smith 
661bd481603SBarry Smith    Output Parameters:
662bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
663bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
664bd481603SBarry Smith 
665bd481603SBarry Smith 
666bd481603SBarry Smith    Level: intermediate
667bd481603SBarry Smith 
668bd481603SBarry Smith    Notes:
669bd481603SBarry Smith    See the chapter in the users manual on performance for more details
670bd481603SBarry Smith 
6711d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
672bd481603SBarry Smith 
6731620fd73SBarry Smith    This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize().
6741620fd73SBarry Smith 
675bd481603SBarry Smith   Concepts: preallocation^Matrix
676bd481603SBarry Smith 
677bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
678bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal()
679bd481603SBarry Smith M*/
680222b16d4SBarry Smith #define MatPreallocateSymmetricInitialize(comm,nrows,ncols,dnz,onz) 0; \
681222b16d4SBarry Smith { \
682a89bc211SBarry Smith   PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__end; \
683a89bc211SBarry Smith   _4_ierr = PetscMalloc2(__nrows,PetscInt,&dnz,__nrows,PetscInt,&onz);CHKERRQ(_4_ierr); \
684a89bc211SBarry Smith   _4_ierr = PetscMemzero(dnz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
685a89bc211SBarry Smith   _4_ierr = PetscMemzero(onz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
686c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr);\
687a89bc211SBarry Smith   _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows;
688222b16d4SBarry Smith 
689bd481603SBarry Smith /*MC
690bd481603SBarry Smith    MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
691bd481603SBarry Smith        inserted using a local number of the rows and columns
692bd481603SBarry Smith 
693bd481603SBarry Smith    Synopsis:
694c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
695bd481603SBarry Smith 
696bd481603SBarry Smith    Not Collective
697bd481603SBarry Smith 
698bd481603SBarry Smith    Input Parameters:
699bd481603SBarry Smith +  map - the mapping between local numbering and global numbering
700bd481603SBarry Smith .  nrows - the number of rows indicated
7011d73ed98SBarry Smith .  rows - the indices of the rows
702bd481603SBarry Smith .  ncols - the number of columns in the matrix
703bd481603SBarry Smith .  cols - the columns indicated
704bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
705bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
706bd481603SBarry Smith 
707bd481603SBarry Smith 
708bd481603SBarry Smith    Level: intermediate
709bd481603SBarry Smith 
710bd481603SBarry Smith    Notes:
711bd481603SBarry Smith    See the chapter in the users manual on performance for more details
712bd481603SBarry Smith 
7131d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
714bd481603SBarry Smith 
715bd481603SBarry Smith   Concepts: preallocation^Matrix
716bd481603SBarry Smith 
717bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
718bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal()
719bd481603SBarry Smith M*/
720c4f061fbSSatish Balay #define MatPreallocateSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\
721c4f061fbSSatish Balay {\
722c1ac3661SBarry Smith   PetscInt __l;\
723ef66eb69SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\
724ef66eb69SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\
725c4f061fbSSatish Balay   for (__l=0;__l<nrows;__l++) {\
726ef66eb69SBarry Smith     _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
727c4f061fbSSatish Balay   }\
728c4f061fbSSatish Balay }
729c4f061fbSSatish Balay 
730bd481603SBarry Smith /*MC
731bd481603SBarry Smith    MatPreallocateSymmetricSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
732bd481603SBarry Smith        inserted using a local number of the rows and columns
733bd481603SBarry Smith 
734bd481603SBarry Smith    Synopsis:
735c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
736bd481603SBarry Smith 
737bd481603SBarry Smith    Not Collective
738bd481603SBarry Smith 
739bd481603SBarry Smith    Input Parameters:
740bd481603SBarry Smith +  map - the mapping between local numbering and global numbering
741bd481603SBarry Smith .  nrows - the number of rows indicated
7421d73ed98SBarry Smith .  rows - the indices of the rows
743bd481603SBarry Smith .  ncols - the number of columns in the matrix
744bd481603SBarry Smith .  cols - the columns indicated
745bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
746bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
747bd481603SBarry Smith 
748bd481603SBarry Smith 
749bd481603SBarry Smith    Level: intermediate
750bd481603SBarry Smith 
751bd481603SBarry Smith    Notes:
752bd481603SBarry Smith    See the chapter in the users manual on performance for more details
753bd481603SBarry Smith 
754bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
755bd481603SBarry Smith 
756bd481603SBarry Smith   Concepts: preallocation^Matrix
757bd481603SBarry Smith 
758bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
759bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
760bd481603SBarry Smith M*/
761d3d32019SBarry Smith #define MatPreallocateSymmetricSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\
762d3d32019SBarry Smith {\
763c1ac3661SBarry Smith   PetscInt __l;\
764d3d32019SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\
765d3d32019SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\
766d3d32019SBarry Smith   for (__l=0;__l<nrows;__l++) {\
767d3d32019SBarry Smith     _4_ierr = MatPreallocateSymmetricSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
768d3d32019SBarry Smith   }\
769d3d32019SBarry Smith }
770d3d32019SBarry Smith 
771bd481603SBarry Smith /*MC
772bd481603SBarry Smith    MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
773bd481603SBarry Smith        inserted using a local number of the rows and columns
774bd481603SBarry Smith 
775bd481603SBarry Smith    Synopsis:
776c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
777bd481603SBarry Smith 
778bd481603SBarry Smith    Not Collective
779bd481603SBarry Smith 
780bd481603SBarry Smith    Input Parameters:
78164075487SBarry Smith +  row - the row
782bd481603SBarry Smith .  ncols - the number of columns in the matrix
783a50f8bf6SHong Zhang -  cols - the columns indicated
784a50f8bf6SHong Zhang 
785a50f8bf6SHong Zhang    Output Parameters:
786a50f8bf6SHong Zhang +  dnz - the array that will be passed to the matrix preallocation routines
787bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
788bd481603SBarry Smith 
789bd481603SBarry Smith 
790bd481603SBarry Smith    Level: intermediate
791bd481603SBarry Smith 
792bd481603SBarry Smith    Notes:
793bd481603SBarry Smith    See the chapter in the users manual on performance for more details
794bd481603SBarry Smith 
795bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
796bd481603SBarry Smith 
7971620fd73SBarry Smith    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
7981620fd73SBarry Smith 
799bd481603SBarry Smith   Concepts: preallocation^Matrix
800bd481603SBarry Smith 
801bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
802bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
803bd481603SBarry Smith M*/
804c4f061fbSSatish Balay #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\
805c1ac3661SBarry Smith { PetscInt __i; \
8068ee2e534SBarry Smith   if (row < __rstart) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Trying to set preallocation for row %D less than first local row %D",row,__rstart);\
807a89bc211SBarry 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);\
8087c922b88SBarry Smith   for (__i=0; __i<nc; __i++) {\
80964075487SBarry Smith     if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \
8107cc688d7SBarry Smith     else dnz[row - __rstart]++;\
8117c922b88SBarry Smith   }\
8127c922b88SBarry Smith }
8137c922b88SBarry Smith 
814bd481603SBarry Smith /*MC
815bd481603SBarry Smith    MatPreallocateSymmetricSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
816bd481603SBarry Smith        inserted using a local number of the rows and columns
817bd481603SBarry Smith 
818bd481603SBarry Smith    Synopsis:
819c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
820bd481603SBarry Smith 
821bd481603SBarry Smith    Not Collective
822bd481603SBarry Smith 
823bd481603SBarry Smith    Input Parameters:
824bd481603SBarry Smith +  nrows - the number of rows indicated
8251d73ed98SBarry Smith .  rows - the indices of the rows
826bd481603SBarry Smith .  ncols - the number of columns in the matrix
827bd481603SBarry Smith .  cols - the columns indicated
828bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
829bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
830bd481603SBarry Smith 
831bd481603SBarry Smith 
832bd481603SBarry Smith    Level: intermediate
833bd481603SBarry Smith 
834bd481603SBarry Smith    Notes:
835bd481603SBarry Smith    See the chapter in the users manual on performance for more details
836bd481603SBarry Smith 
837bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
838bd481603SBarry Smith 
8391620fd73SBarry Smith    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
8401620fd73SBarry Smith 
841bd481603SBarry Smith   Concepts: preallocation^Matrix
842bd481603SBarry Smith 
843bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
844bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
845bd481603SBarry Smith M*/
846d3d32019SBarry Smith #define MatPreallocateSymmetricSet(row,nc,cols,dnz,onz) 0;\
847c1ac3661SBarry Smith { PetscInt __i; \
848d3d32019SBarry Smith   for (__i=0; __i<nc; __i++) {\
849d3d32019SBarry Smith     if (cols[__i] >= __end) onz[row - __rstart]++; \
850d3d32019SBarry Smith     else if (cols[__i] >= row) dnz[row - __rstart]++;\
851d3d32019SBarry Smith   }\
852d3d32019SBarry Smith }
853d3d32019SBarry Smith 
854bd481603SBarry Smith /*MC
85516371a99SBarry Smith    MatPreallocateLocation -  An alternative to MatPreallocationSet() that puts the nonzero locations into the matrix if it exists
85616371a99SBarry Smith 
85716371a99SBarry Smith    Synopsis:
85816371a99SBarry Smith    PetscErrorCode MatPreallocateLocations(Mat A,PetscInt row,PetscInt ncols,PetscInt *cols,PetscInt *dnz,PetscInt *onz)
85916371a99SBarry Smith 
86016371a99SBarry Smith    Not Collective
86116371a99SBarry Smith 
86216371a99SBarry Smith    Input Parameters:
86316371a99SBarry Smith .  A - matrix
86416371a99SBarry Smith .  row - row where values exist (must be local to this process)
86516371a99SBarry Smith .  ncols - number of columns
86616371a99SBarry Smith .  cols - columns with nonzeros
86716371a99SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
86816371a99SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
86916371a99SBarry Smith 
87016371a99SBarry Smith 
87116371a99SBarry Smith    Level: intermediate
87216371a99SBarry Smith 
87316371a99SBarry Smith    Notes:
87416371a99SBarry Smith    See the chapter in the users manual on performance for more details
87516371a99SBarry Smith 
87616371a99SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
87716371a99SBarry Smith 
87816371a99SBarry Smith    This is a MACRO not a function because it uses a bunch of variables private to the MatPreallocation.... routines.
87916371a99SBarry Smith 
88016371a99SBarry Smith   Concepts: preallocation^Matrix
88116371a99SBarry Smith 
88216371a99SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
88316371a99SBarry Smith           MatPreallocateSymmetricInitialize(), MatPreallocateSymmetricSetLocal()
88416371a99SBarry Smith M*/
88516371a99SBarry 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);}
88616371a99SBarry Smith 
88716371a99SBarry Smith 
88816371a99SBarry Smith /*MC
8890d2a7ff2SSatish Balay    MatPreallocateFinalize - Ends the block of code that will count the number of nonzeros per
890bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
891bd481603SBarry Smith 
892bd481603SBarry Smith    Synopsis:
893c1ac3661SBarry Smith    PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz)
894bd481603SBarry Smith 
895bd481603SBarry Smith    Collective on MPI_Comm
896bd481603SBarry Smith 
897bd481603SBarry Smith    Input Parameters:
89816371a99SBarry Smith +  dnz - the array that was be passed to the matrix preallocation routines
899bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
900bd481603SBarry Smith 
901bd481603SBarry Smith 
902bd481603SBarry Smith    Level: intermediate
903bd481603SBarry Smith 
904bd481603SBarry Smith    Notes:
905bd481603SBarry Smith    See the chapter in the users manual on performance for more details
906bd481603SBarry Smith 
907bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
908bd481603SBarry Smith 
9091620fd73SBarry Smith    This is a MACRO not a function because it closes the { started in MatPreallocateInitialize().
9101620fd73SBarry Smith 
911bd481603SBarry Smith   Concepts: preallocation^Matrix
912bd481603SBarry Smith 
913bd481603SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
914bd481603SBarry Smith           MatPreallocateSymmetricInitialize(), MatPreallocateSymmetricSetLocal()
915bd481603SBarry Smith M*/
916a89bc211SBarry Smith #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree2(dnz,onz);CHKERRQ(_4_ierr);}
9177c922b88SBarry Smith 
918bd481603SBarry Smith 
919bd481603SBarry Smith 
9207b80b807SBarry Smith /* Routines unique to particular data structures */
921be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellGetContext(Mat,void **);
922ba966639SSatish Balay PetscPolymorphicFunction(MatShellGetContext,(Mat A),(A,&t),void*,t)
923698d4c6aSKris Buschelman 
924be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInodeAdjustForInodes(Mat,IS*,IS*);
925be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *);
926698d4c6aSKris Buschelman 
927be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetColumnIndices(Mat,PetscInt[]);
928be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetColumnIndices(Mat,PetscInt[]);
929be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
930c75a6043SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
931c75a6043SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
9327b80b807SBarry Smith 
933a96a251dSBarry Smith #define MAT_SKIP_ALLOCATION -4
934a96a251dSBarry Smith 
935be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
936ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz))
937be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
938ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqSBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz))
939be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]);
940ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqAIJSetPreallocation,(Mat A,const PetscInt nnz[]),(A,0,nnz))
941be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDensePreallocation(Mat,PetscScalar[]);
942273d9f13SBarry Smith 
943be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
944ba966639SSatish Balay PetscPolymorphicSubroutine(MatMPIBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[],const PetscInt onz[]),(A,bs,0,nnz,0,onz))
945be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
946be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
94753803f26SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []);
948725b52f3SLisandro Dalcin EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
949be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]);
950be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
951be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]);
952be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIDenseSetPreallocation(Mat,PetscScalar[]);
953284134d9SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation(Mat,PetscScalar[]);
954be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIRowbsSetPreallocation(Mat,PetscInt,const PetscInt[]);
955be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,PetscInt*[]);
956be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,PetscInt*[]);
957be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAdicSetLocalFunction(Mat,void (*)(void));
958273d9f13SBarry Smith 
959be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetLDA(Mat,PetscInt);
960ab92ecdeSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDenseGetLocalMatrix(Mat,Mat*);
9611b807ce4Svictorle 
962be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStoreValues(Mat);
963be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRetrieveValues(Mat);
9642e8a6d31SBarry Smith 
965be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDAADSetCtx(Mat,void*);
9663a7fca6bSBarry Smith 
9677b80b807SBarry Smith /*
9687b80b807SBarry Smith   These routines are not usually accessed directly, rather solving is
96994b7f48cSBarry Smith   done through the KSP and PC interfaces.
9707b80b807SBarry Smith */
9717b80b807SBarry Smith 
972d9274352SBarry Smith /*E
973d9274352SBarry Smith     MatOrderingType - String with the name of a PETSc matrix ordering or the creation function
974d9274352SBarry Smith        with an optional dynamic library name, for example
975d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:orderingcreate()
976d9274352SBarry Smith 
977d9274352SBarry Smith    Level: beginner
978d9274352SBarry Smith 
979e5a9bf91SBarry Smith    Cannot use const because the PC objects manipulate the string
980e5a9bf91SBarry Smith 
981d9274352SBarry Smith .seealso: MatGetOrdering()
982d9274352SBarry Smith E*/
9833eaad2caSSatish Balay #define MatOrderingType char*
984b12f92e5SBarry Smith #define MATORDERING_NATURAL   "natural"
985b12f92e5SBarry Smith #define MATORDERING_ND        "nd"
986b12f92e5SBarry Smith #define MATORDERING_1WD       "1wd"
987b12f92e5SBarry Smith #define MATORDERING_RCM       "rcm"
988b12f92e5SBarry Smith #define MATORDERING_QMD       "qmd"
989b12f92e5SBarry Smith #define MATORDERING_ROWLENGTH "rowlength"
99062152c8bSBarry Smith #define MATORDERING_DSC_ND    "dsc_nd"
99162152c8bSBarry Smith #define MATORDERING_DSC_MMD   "dsc_mmd"
99262152c8bSBarry Smith #define MATORDERING_DSC_MDF   "dsc_mdf"
993c06d978dSMatthew Knepley #define MATORDERING_CONSTRAINED "constrained"
994c06d978dSMatthew Knepley #define MATORDERING_IDENTITY  "identity"
995c06d978dSMatthew Knepley #define MATORDERING_REVERSE   "reverse"
996b12f92e5SBarry Smith 
997f2416520SSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOrdering(Mat,const MatOrderingType,IS*,IS*);
998be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOrderingList(PetscFList *list);
999f2416520SSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,const MatOrderingType,IS*,IS*));
100030de9b25SBarry Smith 
100130de9b25SBarry Smith /*MC
100230de9b25SBarry Smith    MatOrderingRegisterDynamic - Adds a new sparse matrix ordering to the
100330de9b25SBarry Smith                                matrix package.
100430de9b25SBarry Smith 
100530de9b25SBarry Smith    Synopsis:
1006c1ac3661SBarry Smith    PetscErrorCode MatOrderingRegisterDynamic(char *name_ordering,char *path,char *name_create,PetscErrorCode (*routine_create)(MatOrdering))
100730de9b25SBarry Smith 
100830de9b25SBarry Smith    Not Collective
100930de9b25SBarry Smith 
101030de9b25SBarry Smith    Input Parameters:
101130de9b25SBarry Smith +  sname - name of ordering (for example MATORDERING_ND)
101230de9b25SBarry Smith .  path - location of library where creation routine is
101330de9b25SBarry Smith .  name - name of function that creates the ordering type,a string
101430de9b25SBarry Smith -  function - function pointer that creates the ordering
101530de9b25SBarry Smith 
101630de9b25SBarry Smith    Level: developer
101730de9b25SBarry Smith 
101830de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
101930de9b25SBarry Smith    is ignored.
102030de9b25SBarry Smith 
102130de9b25SBarry Smith    Sample usage:
102230de9b25SBarry Smith .vb
102330de9b25SBarry Smith    MatOrderingRegisterDynamic("my_order",/home/username/my_lib/lib/libO/solaris/mylib.a,
102430de9b25SBarry Smith                "MyOrder",MyOrder);
102530de9b25SBarry Smith .ve
102630de9b25SBarry Smith 
102730de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
102830de9b25SBarry Smith $     MatOrderingSetType(part,"my_order)
102930de9b25SBarry Smith    or at runtime via the option
10302401956bSBarry Smith $     -pc_factor_mat_ordering_type my_order
103130de9b25SBarry Smith 
1032ab901514SBarry Smith    ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values.
103330de9b25SBarry Smith 
103430de9b25SBarry Smith .keywords: matrix, ordering, register
103530de9b25SBarry Smith 
103630de9b25SBarry Smith .seealso: MatOrderingRegisterDestroy(), MatOrderingRegisterAll()
103730de9b25SBarry Smith M*/
1038aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1039f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,0)
1040b12f92e5SBarry Smith #else
1041f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,d)
1042b12f92e5SBarry Smith #endif
104330de9b25SBarry Smith 
1044be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegisterDestroy(void);
1045be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegisterAll(const char[]);
10462bad1931SBarry Smith extern PetscTruth MatOrderingRegisterAllCalled;
1047b0a32e0cSBarry Smith extern PetscFList MatOrderingList;
1048d4fbbf0eSBarry Smith 
1049be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS);
1050a2ce50c7SBarry Smith 
1051d91e6319SBarry Smith /*S
10522401956bSBarry Smith    MatFactorInfo - Data passed into the matrix factorization routines
1053b00f7748SHong Zhang 
105461cad860SBarry Smith    In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE, that is use
105561cad860SBarry Smith $     MatFactorInfo  info(MAT_FACTORINFO_SIZE)
1056b00f7748SHong Zhang 
105715e8a5b3SHong Zhang    Notes: These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC.
1058b00f7748SHong Zhang 
105961cad860SBarry Smith       You can use MatFactorInfoInitialize() to set default values.
106061cad860SBarry Smith 
1061b00f7748SHong Zhang    Level: developer
1062b00f7748SHong Zhang 
1063d7d82daaSBarry Smith .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(),
1064d7d82daaSBarry Smith           MatFactorInfoInitialize()
1065b00f7748SHong Zhang 
1066b00f7748SHong Zhang S*/
1067b00f7748SHong Zhang typedef struct {
10680a29876aSHong Zhang   PetscReal     shiftnz;        /* scaling of identity added to matrix to prevent zero pivots */
1069fbf22428SSatish Balay   PetscReal     shiftpd;        /* if true, shift until positive pivots */
10702cea7109SBarry Smith   PetscReal     shift_fraction; /* record shift fraction taken */
107115e8a5b3SHong Zhang   PetscReal     diagonal_fill;  /* force diagonal to fill in if initially not filled */
107215e8a5b3SHong Zhang   PetscReal     dt;             /* drop tolerance */
1073b00f7748SHong Zhang   PetscReal     dtcol;          /* tolerance for pivoting */
107415e8a5b3SHong Zhang   PetscReal     dtcount;        /* maximum nonzeros to be allowed per row */
107567e28bfeSBarry Smith   PetscReal     fill;           /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */
1076348344bbSBarry Smith   PetscReal     levels;         /* ICC/ILU(levels) */
1077bcd9e38bSBarry Smith   PetscReal     pivotinblocks;  /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0
1078bcd9e38bSBarry Smith                                    factorization may be faster if do not pivot */
107962bba022SBarry Smith   PetscReal     shiftinblocks;  /* if block in block factorization has zero pivot then shift diagonal until non-singular */
108015e8a5b3SHong Zhang   PetscReal     zeropivot;      /* pivot is called zero if less than this */
108115e8a5b3SHong Zhang } MatFactorInfo;
1082ffa6d0a5SLois Curfman McInnes 
1083be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFactorInfoInitialize(MatFactorInfo*);
1084be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactor(Mat,IS,MatFactorInfo*);
1085be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactorSymbolic(Mat,IS,MatFactorInfo*,Mat*);
1086be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactorNumeric(Mat,MatFactorInfo*,Mat*);
1087be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactor(Mat,IS,IS,MatFactorInfo*);
1088be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUFactor(Mat,IS,IS,MatFactorInfo*);
1089be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactorSymbolic(Mat,IS,IS,MatFactorInfo*,Mat*);
1090be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUFactorSymbolic(Mat,IS,IS,MatFactorInfo*,Mat*);
1091be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatICCFactorSymbolic(Mat,IS,MatFactorInfo*,Mat*);
1092be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatICCFactor(Mat,IS,MatFactorInfo*);
1093be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactorNumeric(Mat,MatFactorInfo*,Mat*);
1094be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUDTFactor(Mat,IS,IS,MatFactorInfo*,Mat *);
1095be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*);
1096be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolve(Mat,Vec,Vec);
1097be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatForwardSolve(Mat,Vec,Vec);
1098be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatBackwardSolve(Mat,Vec,Vec);
1099be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveAdd(Mat,Vec,Vec,Vec);
1100be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveTranspose(Mat,Vec,Vec);
1101be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveTransposeAdd(Mat,Vec,Vec,Vec);
1102be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolves(Mat,Vecs,Vecs);
11038ed539a5SBarry Smith 
1104be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUnfactored(Mat);
1105bb5a7306SBarry Smith 
1106d91e6319SBarry Smith /*E
1107d91e6319SBarry Smith     MatSORType - What type of (S)SOR to perform
1108bb1eb677SSatish Balay 
1109d91e6319SBarry Smith     Level: beginner
1110d91e6319SBarry Smith 
1111d9274352SBarry Smith    May be bitwise ORd together
1112d9274352SBarry Smith 
1113d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
1114d91e6319SBarry Smith 
11154e7234bfSBarry Smith    MatSORType may be bitwise ORd together, so do not change the numbers
11164e7234bfSBarry Smith 
1117a1d92eedSBarry Smith .seealso: MatRelax(), MatPBRelax()
1118d91e6319SBarry Smith E*/
1119ee50ffe9SBarry Smith typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3,
1120ee50ffe9SBarry Smith               SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8,
1121ee50ffe9SBarry Smith               SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16,
112284cb2905SBarry Smith               SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType;
1123be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRelax(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
1124be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPBRelax(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
11258ed539a5SBarry Smith 
1126d4fbbf0eSBarry Smith /*
1127639f9d9dSBarry Smith     These routines are for efficiently computing Jacobians via finite differences.
1128639f9d9dSBarry Smith */
1129b12f92e5SBarry Smith 
1130d9274352SBarry Smith /*E
1131d9274352SBarry Smith     MatColoringType - String with the name of a PETSc matrix coloring or the creation function
1132d9274352SBarry Smith        with an optional dynamic library name, for example
1133d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:coloringcreate()
1134d9274352SBarry Smith 
1135d9274352SBarry Smith    Level: beginner
1136d9274352SBarry Smith 
1137d9274352SBarry Smith .seealso: MatGetColoring()
1138d9274352SBarry Smith E*/
1139a313700dSBarry Smith #define MatColoringType char*
1140b12f92e5SBarry Smith #define MATCOLORING_NATURAL "natural"
1141b12f92e5SBarry Smith #define MATCOLORING_SL      "sl"
1142b12f92e5SBarry Smith #define MATCOLORING_LF      "lf"
1143b12f92e5SBarry Smith #define MATCOLORING_ID      "id"
1144b12f92e5SBarry Smith 
1145*ab970d9bSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColoring(Mat,const MatColoringType,ISColoring*);
11462e90c967SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,MatColoringType,ISColoring *));
114730de9b25SBarry Smith 
114830de9b25SBarry Smith /*MC
114930de9b25SBarry Smith    MatColoringRegisterDynamic - Adds a new sparse matrix coloring to the
115030de9b25SBarry Smith                                matrix package.
115130de9b25SBarry Smith 
115230de9b25SBarry Smith    Synopsis:
1153c1ac3661SBarry Smith    PetscErrorCode MatColoringRegisterDynamic(char *name_coloring,char *path,char *name_create,PetscErrorCode (*routine_create)(MatColoring))
115430de9b25SBarry Smith 
115530de9b25SBarry Smith    Not Collective
115630de9b25SBarry Smith 
115730de9b25SBarry Smith    Input Parameters:
115830de9b25SBarry Smith +  sname - name of Coloring (for example MATCOLORING_SL)
115930de9b25SBarry Smith .  path - location of library where creation routine is
116030de9b25SBarry Smith .  name - name of function that creates the Coloring type, a string
116130de9b25SBarry Smith -  function - function pointer that creates the coloring
116230de9b25SBarry Smith 
116330de9b25SBarry Smith    Level: developer
116430de9b25SBarry Smith 
116530de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
116630de9b25SBarry Smith    is ignored.
116730de9b25SBarry Smith 
116830de9b25SBarry Smith    Sample usage:
116930de9b25SBarry Smith .vb
117030de9b25SBarry Smith    MatColoringRegisterDynamic("my_color",/home/username/my_lib/lib/libO/solaris/mylib.a,
117130de9b25SBarry Smith                "MyColor",MyColor);
117230de9b25SBarry Smith .ve
117330de9b25SBarry Smith 
117430de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
117530de9b25SBarry Smith $     MatColoringSetType(part,"my_color")
117630de9b25SBarry Smith    or at runtime via the option
117730de9b25SBarry Smith $     -mat_coloring_type my_color
117830de9b25SBarry Smith 
1179ab901514SBarry Smith    $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
118030de9b25SBarry Smith 
118130de9b25SBarry Smith .keywords: matrix, Coloring, register
118230de9b25SBarry Smith 
118330de9b25SBarry Smith .seealso: MatColoringRegisterDestroy(), MatColoringRegisterAll()
118430de9b25SBarry Smith M*/
1185aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1186f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,0)
1187b12f92e5SBarry Smith #else
1188f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,d)
1189b12f92e5SBarry Smith #endif
119030de9b25SBarry Smith 
11912bad1931SBarry Smith extern PetscTruth MatColoringRegisterAllCalled;
1192f1144a30SSatish Balay 
1193be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegisterAll(const char[]);
1194be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegisterDestroy(void);
1195be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*);
1196639f9d9dSBarry Smith 
1197d9274352SBarry Smith /*S
1198d9274352SBarry Smith      MatFDColoring - Object for computing a sparse Jacobian via finite differences
1199d9274352SBarry Smith         and coloring
1200639f9d9dSBarry Smith 
1201d9274352SBarry Smith    Level: beginner
1202d9274352SBarry Smith 
1203d9274352SBarry Smith   Concepts: coloring, sparse Jacobian, finite differences
1204d9274352SBarry Smith 
1205d9274352SBarry Smith .seealso:  MatFDColoringCreate()
1206d9274352SBarry Smith S*/
1207e2a1c21fSSatish Balay typedef struct _p_MatFDColoring* MatFDColoring;
1208639f9d9dSBarry Smith 
1209be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringCreate(Mat,ISColoring,MatFDColoring *);
1210be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringDestroy(MatFDColoring);
1211be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringView(MatFDColoring,PetscViewer);
1212be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*);
12134a9d489dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**);
1214be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal);
1215be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFrequency(MatFDColoring,PetscInt);
1216be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetFrequency(MatFDColoring,PetscInt*);
1217be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFromOptions(MatFDColoring);
1218be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApply(Mat,MatFDColoring,Vec,MatStructure*,void *);
1219be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApplyTS(Mat,MatFDColoring,PetscReal,Vec,MatStructure*,void *);
1220be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetRecompute(MatFDColoring);
1221be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetF(MatFDColoring,Vec);
1222be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,PetscInt*[]);
1223639f9d9dSBarry Smith /*
12240752156aSBarry Smith     These routines are for partitioning matrices: currently used only
12253eda8832SBarry Smith   for adjacency matrix, MatCreateMPIAdj().
12260752156aSBarry Smith */
1227ca161407SBarry Smith 
1228d9274352SBarry Smith /*S
1229d9274352SBarry Smith      MatPartitioning - Object for managing the partitioning of a matrix or graph
1230d9274352SBarry Smith 
1231d9274352SBarry Smith    Level: beginner
1232d9274352SBarry Smith 
1233d9274352SBarry Smith   Concepts: partitioning
1234d9274352SBarry Smith 
1235743a1026Svictorle .seealso:  MatPartitioningCreate(), MatPartitioningType
1236d9274352SBarry Smith S*/
123791e9ee9fSBarry Smith typedef struct _p_MatPartitioning* MatPartitioning;
1238d9274352SBarry Smith 
1239d9274352SBarry Smith /*E
12405bc6e7ccSvictorle     MatPartitioningType - String with the name of a PETSc matrix partitioning or the creation function
1241d9274352SBarry Smith        with an optional dynamic library name, for example
1242d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:partitioningcreate()
1243d9274352SBarry Smith 
1244d9274352SBarry Smith    Level: beginner
1245d9274352SBarry Smith 
1246b547a13bSvictorle .seealso: MatPartitioningCreate(), MatPartitioning
1247d9274352SBarry Smith E*/
1248a313700dSBarry Smith #define MatPartitioningType char*
12498ba1e511SMatthew Knepley #define MAT_PARTITIONING_CURRENT  "current"
1250c2ec2a73SHong Zhang #define MAT_PARTITIONING_SQUARE   "square"
12518ba1e511SMatthew Knepley #define MAT_PARTITIONING_PARMETIS "parmetis"
125271306c60Sjroman #define MAT_PARTITIONING_CHACO    "chaco"
125371306c60Sjroman #define MAT_PARTITIONING_JOSTLE   "jostle"
125471306c60Sjroman #define MAT_PARTITIONING_PARTY    "party"
125571306c60Sjroman #define MAT_PARTITIONING_SCOTCH   "scotch"
125671306c60Sjroman 
1257ca161407SBarry Smith 
1258be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningCreate(MPI_Comm,MatPartitioning*);
1259a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetType(MatPartitioning,const MatPartitioningType);
1260be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetNParts(MatPartitioning,PetscInt);
1261be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetAdjacency(MatPartitioning,Mat);
1262be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]);
1263be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []);
1264be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningApply(MatPartitioning,IS*);
1265be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningDestroy(MatPartitioning);
12662aabb6bbSBarry Smith 
1267be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatPartitioning));
126830de9b25SBarry Smith 
126930de9b25SBarry Smith /*MC
127030de9b25SBarry Smith    MatPartitioningRegisterDynamic - Adds a new sparse matrix partitioning to the
127130de9b25SBarry Smith    matrix package.
127230de9b25SBarry Smith 
127330de9b25SBarry Smith    Synopsis:
1274c1ac3661SBarry Smith    PetscErrorCode MatPartitioningRegisterDynamic(char *name_partitioning,char *path,char *name_create,PetscErrorCode (*routine_create)(MatPartitioning))
127530de9b25SBarry Smith 
127630de9b25SBarry Smith    Not Collective
127730de9b25SBarry Smith 
127830de9b25SBarry Smith    Input Parameters:
127930de9b25SBarry Smith +  sname - name of partitioning (for example MAT_PARTITIONING_CURRENT) or parmetis
128030de9b25SBarry Smith .  path - location of library where creation routine is
128130de9b25SBarry Smith .  name - name of function that creates the partitioning type, a string
128230de9b25SBarry Smith -  function - function pointer that creates the partitioning type
128330de9b25SBarry Smith 
128430de9b25SBarry Smith    Level: developer
128530de9b25SBarry Smith 
128630de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
128730de9b25SBarry Smith    is ignored.
128830de9b25SBarry Smith 
128930de9b25SBarry Smith    Sample usage:
129030de9b25SBarry Smith .vb
129130de9b25SBarry Smith    MatPartitioningRegisterDynamic("my_part",/home/username/my_lib/lib/libO/solaris/mylib.a,
129230de9b25SBarry Smith                "MyPartCreate",MyPartCreate);
129330de9b25SBarry Smith .ve
129430de9b25SBarry Smith 
129530de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
129630de9b25SBarry Smith $     MatPartitioningSetType(part,"my_part")
129730de9b25SBarry Smith    or at runtime via the option
129830de9b25SBarry Smith $     -mat_partitioning_type my_part
129930de9b25SBarry Smith 
1300ab901514SBarry Smith    $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
130130de9b25SBarry Smith 
130230de9b25SBarry Smith .keywords: matrix, partitioning, register
130330de9b25SBarry Smith 
130430de9b25SBarry Smith .seealso: MatPartitioningRegisterDestroy(), MatPartitioningRegisterAll()
130530de9b25SBarry Smith M*/
1306aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1307f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,0)
13082aabb6bbSBarry Smith #else
1309f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,d)
13102aabb6bbSBarry Smith #endif
13112aabb6bbSBarry Smith 
13122bad1931SBarry Smith extern PetscTruth MatPartitioningRegisterAllCalled;
1313f1144a30SSatish Balay 
1314be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegisterAll(const char[]);
1315be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegisterDestroy(void);
13162bad1931SBarry Smith 
1317be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningView(MatPartitioning,PetscViewer);
1318be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetFromOptions(MatPartitioning);
1319a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningGetType(MatPartitioning,const MatPartitioningType*);
1320ca161407SBarry Smith 
1321be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningParmetisSetCoarseSequential(MatPartitioning);
13220752156aSBarry Smith 
1323be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningJostleSetCoarseLevel(MatPartitioning,PetscReal);
1324be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningJostleSetCoarseSequential(MatPartitioning);
132571306c60Sjroman 
1326954b3ec6SBarry Smith typedef enum {MP_CHACO_MULTILEVEL_KL,MP_CHACO_SPECTRAL,MP_CHACO_LINEAR,MP_CHACO_RANDOM, MP_CHACO_SCATTERED} MPChacoGlobalType;
1327be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetGlobal(MatPartitioning, MPChacoGlobalType);
132871306c60Sjroman typedef enum { MP_CHACO_KERNIGHAN_LIN, MP_CHACO_NONE } MPChacoLocalType;
1329be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetLocal(MatPartitioning, MPChacoLocalType);
1330be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal);
133171306c60Sjroman typedef enum { MP_CHACO_LANCZOS, MP_CHACO_RQI_SYMMLQ } MPChacoEigenType;
1332be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType);
1333be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenTol(MatPartitioning, PetscReal);
1334be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenNumber(MatPartitioning, PetscInt);
133571306c60Sjroman 
133671306c60Sjroman #define MP_PARTY_OPT "opt"
133771306c60Sjroman #define MP_PARTY_LIN "lin"
133871306c60Sjroman #define MP_PARTY_SCA "sca"
133971306c60Sjroman #define MP_PARTY_RAN "ran"
134071306c60Sjroman #define MP_PARTY_GBF "gbf"
134171306c60Sjroman #define MP_PARTY_GCF "gcf"
134271306c60Sjroman #define MP_PARTY_BUB "bub"
134371306c60Sjroman #define MP_PARTY_DEF "def"
1344be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetGlobal(MatPartitioning, const char*);
134571306c60Sjroman #define MP_PARTY_HELPFUL_SETS "hs"
134671306c60Sjroman #define MP_PARTY_KERNIGHAN_LIN "kl"
134771306c60Sjroman #define MP_PARTY_NONE "no"
1348be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetLocal(MatPartitioning, const char*);
1349be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal);
1350be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetBipart(MatPartitioning,PetscTruth);
1351be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscTruth);
135271306c60Sjroman 
135371306c60Sjroman typedef enum { MP_SCOTCH_GREEDY, MP_SCOTCH_GPS, MP_SCOTCH_GR_GPS } MPScotchGlobalType;
1354be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetArch(MatPartitioning,const char*);
1355be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetMultilevel(MatPartitioning);
1356be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetGlobal(MatPartitioning,MPScotchGlobalType);
1357be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetCoarseLevel(MatPartitioning,PetscReal);
1358be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetHostList(MatPartitioning,const char*);
135971306c60Sjroman typedef enum { MP_SCOTCH_KERNIGHAN_LIN, MP_SCOTCH_NONE } MPScotchLocalType;
1360be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetLocal(MatPartitioning,MPScotchLocalType);
1361be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetMapping(MatPartitioning);
1362be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetStrategy(MatPartitioning,char*);
136371306c60Sjroman 
13640752156aSBarry Smith /*
13650a835dfdSSatish Balay     If you add entries here you must also add them to finclude/petscmat.h
1366d4fbbf0eSBarry Smith */
13671c1c02c0SLois Curfman McInnes typedef enum { MATOP_SET_VALUES=0,
13681c1c02c0SLois Curfman McInnes                MATOP_GET_ROW=1,
13691c1c02c0SLois Curfman McInnes                MATOP_RESTORE_ROW=2,
13701c1c02c0SLois Curfman McInnes                MATOP_MULT=3,
13711c1c02c0SLois Curfman McInnes                MATOP_MULT_ADD=4,
13727c922b88SBarry Smith                MATOP_MULT_TRANSPOSE=5,
13737c922b88SBarry Smith                MATOP_MULT_TRANSPOSE_ADD=6,
13741c1c02c0SLois Curfman McInnes                MATOP_SOLVE=7,
13751c1c02c0SLois Curfman McInnes                MATOP_SOLVE_ADD=8,
13767c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE=9,
13777c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE_ADD=10,
13781c1c02c0SLois Curfman McInnes                MATOP_LUFACTOR=11,
13791c1c02c0SLois Curfman McInnes                MATOP_CHOLESKYFACTOR=12,
13801c1c02c0SLois Curfman McInnes                MATOP_RELAX=13,
13811c1c02c0SLois Curfman McInnes                MATOP_TRANSPOSE=14,
13821c1c02c0SLois Curfman McInnes                MATOP_GETINFO=15,
13831c1c02c0SLois Curfman McInnes                MATOP_EQUAL=16,
13841c1c02c0SLois Curfman McInnes                MATOP_GET_DIAGONAL=17,
13851c1c02c0SLois Curfman McInnes                MATOP_DIAGONAL_SCALE=18,
13861c1c02c0SLois Curfman McInnes                MATOP_NORM=19,
13871c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_BEGIN=20,
13881c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_END=21,
13891c1c02c0SLois Curfman McInnes                MATOP_COMPRESS=22,
13901c1c02c0SLois Curfman McInnes                MATOP_SET_OPTION=23,
13911c1c02c0SLois Curfman McInnes                MATOP_ZERO_ENTRIES=24,
13921c1c02c0SLois Curfman McInnes                MATOP_ZERO_ROWS=25,
13931c1c02c0SLois Curfman McInnes                MATOP_LUFACTOR_SYMBOLIC=26,
13941c1c02c0SLois Curfman McInnes                MATOP_LUFACTOR_NUMERIC=27,
13951c1c02c0SLois Curfman McInnes                MATOP_CHOLESKY_FACTOR_SYMBOLIC=28,
13961c1c02c0SLois Curfman McInnes                MATOP_CHOLESKY_FACTOR_NUMERIC=29,
1397d643ce63SMatthew Knepley                MATOP_SETUP_PREALLOCATION=30,
1398d643ce63SMatthew Knepley                MATOP_ILUFACTOR_SYMBOLIC=31,
1399d643ce63SMatthew Knepley                MATOP_ICCFACTOR_SYMBOLIC=32,
1400d643ce63SMatthew Knepley                MATOP_GET_ARRAY=33,
1401d643ce63SMatthew Knepley                MATOP_RESTORE_ARRAY=34,
1402e26f1c3cSBarry Smith                MATOP_DUPLICATE=35,
1403d643ce63SMatthew Knepley                MATOP_FORWARD_SOLVE=36,
1404d643ce63SMatthew Knepley                MATOP_BACKWARD_SOLVE=37,
1405d643ce63SMatthew Knepley                MATOP_ILUFACTOR=38,
1406d643ce63SMatthew Knepley                MATOP_ICCFACTOR=39,
1407d643ce63SMatthew Knepley                MATOP_AXPY=40,
1408d643ce63SMatthew Knepley                MATOP_GET_SUBMATRICES=41,
1409d643ce63SMatthew Knepley                MATOP_INCREASE_OVERLAP=42,
1410d643ce63SMatthew Knepley                MATOP_GET_VALUES=43,
1411d643ce63SMatthew Knepley                MATOP_COPY=44,
1412d643ce63SMatthew Knepley                MATOP_PRINT_HELP=45,
1413d643ce63SMatthew Knepley                MATOP_SCALE=46,
1414d643ce63SMatthew Knepley                MATOP_SHIFT=47,
1415d643ce63SMatthew Knepley                MATOP_DIAGONAL_SHIFT=48,
1416d643ce63SMatthew Knepley                MATOP_ILUDT_FACTOR=49,
1417d643ce63SMatthew Knepley                MATOP_GET_BLOCK_SIZE=50,
1418d643ce63SMatthew Knepley                MATOP_GET_ROW_IJ=51,
1419d643ce63SMatthew Knepley                MATOP_RESTORE_ROW_IJ=52,
1420d643ce63SMatthew Knepley                MATOP_GET_COLUMN_IJ=53,
1421d643ce63SMatthew Knepley                MATOP_RESTORE_COLUMN_IJ=54,
1422d643ce63SMatthew Knepley                MATOP_FDCOLORING_CREATE=55,
1423d643ce63SMatthew Knepley                MATOP_COLORING_PATCH=56,
1424d643ce63SMatthew Knepley                MATOP_SET_UNFACTORED=57,
1425d643ce63SMatthew Knepley                MATOP_PERMUTE=58,
1426d643ce63SMatthew Knepley                MATOP_SET_VALUES_BLOCKED=59,
1427d643ce63SMatthew Knepley                MATOP_GET_SUBMATRIX=60,
1428d643ce63SMatthew Knepley                MATOP_DESTROY=61,
1429d643ce63SMatthew Knepley                MATOP_VIEW=62,
1430d643ce63SMatthew Knepley                MATOP_GET_MAPS=63,
1431d643ce63SMatthew Knepley                MATOP_USE_SCALED_FORM=64,
1432d643ce63SMatthew Knepley                MATOP_SCALE_SYSTEM=65,
1433d643ce63SMatthew Knepley                MATOP_UNSCALE_SYSTEM=66,
1434ba16a8cbSKris Buschelman                MATOP_SET_LOCAL_TO_GLOBAL_MAP=67,
1435d643ce63SMatthew Knepley                MATOP_SET_VALUES_LOCAL=68,
1436d643ce63SMatthew Knepley                MATOP_ZERO_ROWS_LOCAL=69,
1437d643ce63SMatthew Knepley                MATOP_GET_ROW_MAX=70,
1438d643ce63SMatthew Knepley                MATOP_CONVERT=71,
1439d643ce63SMatthew Knepley                MATOP_SET_COLORING=72,
1440d643ce63SMatthew Knepley                MATOP_SET_VALUES_ADIC=73,
1441d643ce63SMatthew Knepley                MATOP_SET_VALUES_ADIFOR=74,
1442d643ce63SMatthew Knepley                MATOP_FD_COLORING_APPLY=75,
1443d643ce63SMatthew Knepley                MATOP_SET_FROM_OPTIONS=76,
1444ba16a8cbSKris Buschelman                MATOP_MULT_CON=77,
1445ba16a8cbSKris Buschelman                MATOP_MULT_TRANSPOSE_CON=78,
1446ba16a8cbSKris Buschelman                MATOP_ILU_FACTOR_SYMBOLIC_CON=79,
1447d643ce63SMatthew Knepley                MATOP_PERMUTE_SPARSIFY=80,
1448d643ce63SMatthew Knepley                MATOP_MULT_MULTIPLE=81,
144941acf15aSKris Buschelman                MATOP_SOLVE_MULTIPLE=82,
145041acf15aSKris Buschelman                MATOP_GET_INERTIA=83,
145141acf15aSKris Buschelman                MATOP_LOAD=84,
145241acf15aSKris Buschelman                MATOP_IS_SYMMETRIC=85,
145341acf15aSKris Buschelman                MATOP_IS_HERMITIAN=86,
145441acf15aSKris Buschelman                MATOP_IS_STRUCTURALLY_SYMMETRIC=87,
145541acf15aSKris Buschelman                MATOP_PB_RELAX=88,
145641acf15aSKris Buschelman                MATOP_GET_VECS=89,
145741acf15aSKris Buschelman                MATOP_MAT_MULT=90,
145841acf15aSKris Buschelman                MATOP_MAT_MULT_SYMBOLIC=91,
145941acf15aSKris Buschelman                MATOP_MAT_MULT_NUMERIC=92,
146041acf15aSKris Buschelman                MATOP_PTAP=93,
146141acf15aSKris Buschelman                MATOP_PTAP_SYMBOLIC=94,
1462bc011b1eSHong Zhang                MATOP_PTAP_NUMERIC=95,
1463bc011b1eSHong Zhang                MATOP_MAT_MULTTRANSPOSE=96,
1464bc011b1eSHong Zhang                MATOP_MAT_MULTTRANSPOSE_SYMBOLIC=97,
14657ba1a0bfSKris Buschelman                MATOP_MAT_MULTTRANSPOSE_NUMERIC=98,
14667ba1a0bfSKris Buschelman                MATOP_PTAP_SYMBOLIC_SEQAIJ=99,
14677ba1a0bfSKris Buschelman                MATOP_PTAP_NUMERIC_SEQAIJ=100,
14687ba1a0bfSKris Buschelman                MATOP_PTAP_SYMBOLIC_MPIAIJ=101,
146987d4246cSBarry Smith                MATOP_PTAP_NUMERIC_MPIAIJ=102,
1470f5edf698SHong Zhang                MATOP_SET_VALUES_ROW = 105,
1471d5c63f83SSatish Balay                MATOP_GET_ROW_UTRIANGULAR=108,
14722bebee5dSHong Zhang                MATOP_RESTORE_ROW_UTRIANGULAR=109,
147369db28dcSHong Zhang                MATOP_MATSOLVE=110,
1474829201f2SHong Zhang                MATOP_GET_REDUNDANTMATRIX=111,
1475829201f2SHong Zhang                MATOP_MATGETSEQNONZEROSTRUCTURE=115
1476fae171e0SBarry Smith              } MatOperation;
1477be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatHasOperation(Mat,MatOperation,PetscTruth*);
1478be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellSetOperation(Mat,MatOperation,void(*)(void));
1479be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellGetOperation(Mat,MatOperation,void(**)(void));
1480be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellSetContext(Mat,void*);
1481112a2221SBarry Smith 
148290ace30eSBarry Smith /*
148390ace30eSBarry Smith    Codes for matrices stored on disk. By default they are
148490ace30eSBarry Smith  stored in a universal format. By changing the format with
1485fb9695e5SSatish Balay  PetscViewerSetFormat(viewer,PETSC_VIEWER_BINARY_NATIVE); the matrices will
148690ace30eSBarry Smith  be stored in a way natural for the matrix, for example dense matrices
148790ace30eSBarry Smith  would be stored as dense. Matrices stored this way may only be
148890ace30eSBarry Smith  read into matrices of the same time.
148990ace30eSBarry Smith */
149090ace30eSBarry Smith #define MATRIX_BINARY_FORMAT_DENSE -1
149190ace30eSBarry Smith 
1492be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetHashTableFactor(Mat,PetscReal);
1493be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatISGetLocalMat(Mat,Mat*);
14941f4e1ec7SBarry Smith 
1495d9274352SBarry Smith /*S
1496d9274352SBarry Smith      MatNullSpace - Object that removes a null space from a vector, i.e.
1497d9274352SBarry Smith          orthogonalizes the vector to a subsapce
1498d9274352SBarry Smith 
1499f7a9e4ceSBarry Smith    Level: advanced
1500d9274352SBarry Smith 
1501d9274352SBarry Smith   Concepts: matrix; linear operator, null space
1502d9274352SBarry Smith 
15036e1639daSBarry Smith   Users manual sections:
15046e1639daSBarry Smith .   sec_singular
15056e1639daSBarry Smith 
1506d9274352SBarry Smith .seealso:  MatNullSpaceCreate()
1507d9274352SBarry Smith S*/
150874637425SBarry Smith typedef struct _p_MatNullSpace* MatNullSpace;
1509d9274352SBarry Smith 
1510be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceCreate(MPI_Comm,PetscTruth,PetscInt,const Vec[],MatNullSpace*);
1511281d1b2eSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(Vec,void*),void*);
1512be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceDestroy(MatNullSpace);
1513be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceRemove(MatNullSpace,Vec,Vec*);
1514be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceAttach(Mat,MatNullSpace);
1515be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceTest(MatNullSpace,Mat);
151674637425SBarry Smith 
1517be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatReorderingSeqSBAIJ(Mat,IS);
1518be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPISBAIJSetHashTableFactor(Mat,PetscReal);
1519be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetColumnIndices(Mat,PetscInt *);
152046129b97SKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJInvertBlockDiagonal(Mat);
15213f1d51d7SBarry Smith 
1522be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMAIJ(Mat,PetscInt,Mat*);
1523be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMAIJRedimension(Mat,PetscInt,Mat*);
1524be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMAIJGetAIJ(Mat,Mat*);
1525c4f061fbSSatish Balay 
1526be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatComputeExplicitOperator(Mat,Mat*);
1527b0a32e0cSBarry Smith 
1528be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalScaleLocal(Mat,Vec);
152904f1ad80SBarry Smith 
1530fef1beadSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*);
15313ec795f1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetBase(Mat,Vec,Vec);
15323ec795f1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*);
1533e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*));
1534e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec));
1535e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDAddNullSpace(Mat,MatNullSpace);
1536e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt);
1537e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDResetHHistory(Mat);
1538e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctionError(Mat,PetscReal);
1539e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetPeriod(Mat,PetscInt);
1540e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDGetH(Mat,PetscScalar *);
1541e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFromOptions(Mat);
1542e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*);
1543e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*);
1544e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDDSSetUmin(Mat,PetscReal);
1545e884886fSBarry Smith 
1546e884886fSBarry Smith 
1547e884886fSBarry Smith typedef struct _p_MatMFFD* MatMFFD;
1548e884886fSBarry Smith 
1549e884886fSBarry Smith /*E
1550e884886fSBarry Smith     MatMFFDType - algorithm used to compute the h used in computing matrix-vector products via differencing of the function
1551e884886fSBarry Smith 
1552e884886fSBarry Smith    Level: beginner
1553e884886fSBarry Smith 
1554e884886fSBarry Smith .seealso: MatMFFDSetType(), MatMFFDRegister()
1555e884886fSBarry Smith E*/
1556a313700dSBarry Smith #define MatMFFDType char*
1557e884886fSBarry Smith #define MATMFFD_DS  "ds"
1558e884886fSBarry Smith #define MATMFFD_WP  "wp"
1559e884886fSBarry Smith 
1560a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetType(Mat,const MatMFFDType);
1561e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatMFFD));
1562e884886fSBarry Smith 
1563e884886fSBarry Smith /*MC
1564e884886fSBarry Smith    MatMFFDRegisterDynamic - Adds a method to the MatMFFD registry.
1565e884886fSBarry Smith 
1566e884886fSBarry Smith    Synopsis:
1567e884886fSBarry Smith    PetscErrorCode MatMFFDRegisterDynamic(char *name_solver,char *path,char *name_create,PetscErrorCode (*routine_create)(MatMFFD))
1568e884886fSBarry Smith 
1569e884886fSBarry Smith    Not Collective
1570e884886fSBarry Smith 
1571e884886fSBarry Smith    Input Parameters:
1572e884886fSBarry Smith +  name_solver - name of a new user-defined compute-h module
1573e884886fSBarry Smith .  path - path (either absolute or relative) the library containing this solver
1574e884886fSBarry Smith .  name_create - name of routine to create method context
1575e884886fSBarry Smith -  routine_create - routine to create method context
1576e884886fSBarry Smith 
1577e884886fSBarry Smith    Level: developer
1578e884886fSBarry Smith 
1579e884886fSBarry Smith    Notes:
1580e884886fSBarry Smith    MatMFFDRegisterDynamic() may be called multiple times to add several user-defined solvers.
1581e884886fSBarry Smith 
1582e884886fSBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
1583e884886fSBarry Smith    is ignored.
1584e884886fSBarry Smith 
1585e884886fSBarry Smith    Sample usage:
1586e884886fSBarry Smith .vb
1587e884886fSBarry Smith    MatMFFDRegisterDynamic("my_h",/home/username/my_lib/lib/libO/solaris/mylib.a,
1588e884886fSBarry Smith                "MyHCreate",MyHCreate);
1589e884886fSBarry Smith .ve
1590e884886fSBarry Smith 
1591e884886fSBarry Smith    Then, your solver can be chosen with the procedural interface via
1592e884886fSBarry Smith $     MatMFFDSetType(mfctx,"my_h")
1593e884886fSBarry Smith    or at runtime via the option
1594e884886fSBarry Smith $     -snes_mf_type my_h
1595e884886fSBarry Smith 
1596e884886fSBarry Smith .keywords: MatMFFD, register
1597e884886fSBarry Smith 
1598e884886fSBarry Smith .seealso: MatMFFDRegisterAll(), MatMFFDRegisterDestroy()
1599e884886fSBarry Smith M*/
1600e884886fSBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1601e884886fSBarry Smith #define MatMFFDRegisterDynamic(a,b,c,d) MatMFFDRegister(a,b,c,0)
1602e884886fSBarry Smith #else
1603e884886fSBarry Smith #define MatMFFDRegisterDynamic(a,b,c,d) MatMFFDRegister(a,b,c,d)
1604e884886fSBarry Smith #endif
1605e884886fSBarry Smith 
1606e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDRegisterAll(const char[]);
1607e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDRegisterDestroy(void);
1608e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDDefaultSetUmin(Mat,PetscReal);
1609e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDWPSetComputeNormU(Mat,PetscTruth);
1610e884886fSBarry Smith 
1611e884886fSBarry Smith 
1612be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *);
1613be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *);
16147dbadf16SMatthew Knepley 
1615e9fa29b7SSatish Balay PETSC_EXTERN_CXX_END
16162eac72dbSBarry Smith #endif
1617