xref: /petsc/include/petscmat.h (revision b022a5c136836438e4f06517ed69d7947f7243d0)
12eac72dbSBarry Smith /*
22eac72dbSBarry Smith      Include file for the matrix component of PETSc
32eac72dbSBarry Smith */
40a835dfdSSatish Balay #ifndef __PETSCMAT_H
50a835dfdSSatish Balay #define __PETSCMAT_H
60a835dfdSSatish Balay #include "petscvec.h"
7e9fa29b7SSatish Balay PETSC_EXTERN_CXX_BEGIN
82eac72dbSBarry Smith 
9d9274352SBarry Smith /*S
10d9274352SBarry Smith      Mat - Abstract PETSc matrix object
112eac72dbSBarry Smith 
12d91e6319SBarry Smith    Level: beginner
13d91e6319SBarry Smith 
14d9274352SBarry Smith   Concepts: matrix; linear operator
15d9274352SBarry Smith 
16d9274352SBarry Smith .seealso:  MatCreate(), MatType, MatSetType()
17d9274352SBarry Smith S*/
18d9274352SBarry Smith typedef struct _p_Mat*           Mat;
19d9274352SBarry Smith 
20d9274352SBarry Smith /*E
21d9274352SBarry Smith     MatType - String with the name of a PETSc matrix or the creation function
22d9274352SBarry Smith        with an optional dynamic library name, for example
23d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:mymatcreate()
24d9274352SBarry Smith 
25d9274352SBarry Smith    Level: beginner
26d9274352SBarry Smith 
27c7393fdbSBarry Smith .seealso: MatSetType(), Mat, MatSolverPackage
28d91e6319SBarry Smith E*/
29a313700dSBarry Smith #define MatType char*
30273d9f13SBarry Smith #define MATSAME            "same"
31273d9f13SBarry Smith #define MATSEQMAIJ         "seqmaij"
32273d9f13SBarry Smith #define MATMPIMAIJ         "mpimaij"
33209238afSKris Buschelman #define MATMAIJ            "maij"
34273d9f13SBarry Smith #define MATIS              "is"
35273d9f13SBarry Smith #define MATMPIROWBS        "mpirowbs"
36273d9f13SBarry Smith #define MATSEQAIJ          "seqaij"
37273d9f13SBarry Smith #define MATMPIAIJ          "mpiaij"
38209238afSKris Buschelman #define MATAIJ             "aij"
39273d9f13SBarry Smith #define MATSHELL           "shell"
40209238afSKris Buschelman #define MATSEQDENSE        "seqdense"
41273d9f13SBarry Smith #define MATMPIDENSE        "mpidense"
42209238afSKris Buschelman #define MATDENSE           "dense"
43273d9f13SBarry Smith #define MATSEQBAIJ         "seqbaij"
44273d9f13SBarry Smith #define MATMPIBAIJ         "mpibaij"
45209238afSKris Buschelman #define MATBAIJ            "baij"
46273d9f13SBarry Smith #define MATMPIADJ          "mpiadj"
47273d9f13SBarry Smith #define MATSEQSBAIJ        "seqsbaij"
48273d9f13SBarry Smith #define MATMPISBAIJ        "mpisbaij"
49209238afSKris Buschelman #define MATSBAIJ           "sbaij"
50cebc7f6cSBarry Smith #define MATDAAD            "daad"
51cebc7f6cSBarry Smith #define MATMFFD            "mffd"
52c8a8475eSBarry Smith #define MATNORMAL          "normal"
53ab92ecdeSBarry Smith #define MATLRC             "lrc"
544b8d99adSRichard Tran Mills #define MATSEQCSRPERM      "seqcsrperm"
5518def556SRichard Tran Mills #define MATMPICSRPERM      "mpicsrperm"
56814655b8SBarry Smith #define MATCSRPERM         "csrperm"
5781824310SBarry Smith #define MATSEQCRL          "seqcrl"
58c02b24c5SBarry Smith #define MATMPICRL          "mpicrl"
59c02b24c5SBarry Smith #define MATCRL             "crl"
602a6744ebSBarry Smith #define MATSCATTER         "scatter"
61421e10b8SBarry Smith #define MATBLOCKMAT        "blockmat"
62793850ffSBarry Smith #define MATCOMPOSITE       "composite"
635a7f1df3SHong Zhang #define MATSEQFFTW         "seqfftw"
64557cca28SSatish Balay #define MATTRANSPOSEMAT    "transpose"
6572ca8751SBarry Smith #define MATSCHURCOMPLEMENT "schurcomplement"
661d6018f0SLisandro Dalcin #define MATPYTHON          "python"
67f91d8e95SBarry Smith #define MATHYPRESTRUCT     "hyprestruct"
68ee1cef2cSJed Brown #define MATSUBMATRIX       "submatrix"
69d91e6319SBarry Smith 
709989ab13SBarry Smith /*E
71c7393fdbSBarry Smith     MatSolverPackage - String with the name of a PETSc matrix solver type.
729989ab13SBarry Smith 
739989ab13SBarry Smith     For example: "petsc" indicates what PETSc provides, "superlu" indicates either
749989ab13SBarry Smith        SuperLU or SuperLU_Dist etc.
759989ab13SBarry Smith 
769989ab13SBarry Smith 
779989ab13SBarry Smith    Level: beginner
789989ab13SBarry Smith 
795c9eb25fSBarry Smith .seealso: MatGetFactor(), Mat, MatSetType(), MatType
809989ab13SBarry Smith E*/
81c7393fdbSBarry Smith #define MatSolverPackage char*
82b24902e0SBarry Smith #define MAT_SOLVER_SPOOLES      "spooles"
83b24902e0SBarry Smith #define MAT_SOLVER_SUPERLU      "superlu"
84b24902e0SBarry Smith #define MAT_SOLVER_SUPERLU_DIST "superlu_dist"
85b24902e0SBarry Smith #define MAT_SOLVER_UMFPACK      "umfpack"
86b24902e0SBarry Smith #define MAT_SOLVER_ESSL         "essl"
87b24902e0SBarry Smith #define MAT_SOLVER_LUSOL        "lusol"
88b24902e0SBarry Smith #define MAT_SOLVER_MUMPS        "mumps"
893bf14a46SMatthew Knepley #define MAT_SOLVER_PASTIX       "pastix"
90b24902e0SBarry Smith #define MAT_SOLVER_DSCPACK      "dscpack"
91b24902e0SBarry Smith #define MAT_SOLVER_MATLAB       "matlab"
9243244d56SBarry Smith #define MAT_SOLVER_PETSC        "petsc"
93b6806ab0SHong Zhang #define MAT_SOLVER_PLAPACK      "plapack"
94b24902e0SBarry Smith 
95b24902e0SBarry Smith /*E
96b24902e0SBarry Smith     MatFactorType - indicates what type of factorization is requested
97b24902e0SBarry Smith 
98b24902e0SBarry Smith     Level: beginner
99b24902e0SBarry Smith 
100b24902e0SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
101b24902e0SBarry Smith 
102c7393fdbSBarry Smith .seealso: MatSolverPackage, MatGetFactor()
103b24902e0SBarry Smith E*/
104599ef60dSHong Zhang typedef enum {MAT_FACTOR_NONE, MAT_FACTOR_LU, MAT_FACTOR_CHOLESKY, MAT_FACTOR_ILU, MAT_FACTOR_ICC,MAT_FACTOR_ILUDT} MatFactorType;
105c7393fdbSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetFactor(Mat,const MatSolverPackage,MatFactorType,Mat*);
106db4efbfdSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetFactorAvailable(Mat,const MatSolverPackage,MatFactorType,PetscTruth*);
10735bd34faSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFactorGetSolverPackage(Mat,const MatSolverPackage*);
108b24902e0SBarry Smith 
1099989ab13SBarry Smith 
110c06d978dSMatthew Knepley /* Logging support */
111552e946dSBarry Smith #define    MAT_FILE_COOKIE 1211216    /* used to indicate matrices in binary files */
112be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_COOKIE;
113be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_FDCOLORING_COOKIE;
114be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_PARTITIONING_COOKIE;
115be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_NULLSPACE_COOKIE;
116e884886fSBarry Smith extern PetscCookie PETSCMAT_DLLEXPORT MATMFFD_COOKIE;
117c06d978dSMatthew Knepley 
118ceb03754SKris Buschelman /*E
119ceb03754SKris Buschelman     MatReuse - Indicates if matrices obtained from a previous call to MatGetSubMatrices()
120ceb03754SKris Buschelman      or MatGetSubMatrix() are to be reused to store the new matrix values.
121ceb03754SKris Buschelman 
122ceb03754SKris Buschelman     Level: beginner
123ceb03754SKris Buschelman 
124ceb03754SKris Buschelman    Any additions/changes here MUST also be made in include/finclude/petscmat.h
125ceb03754SKris Buschelman 
1260c451bc4SBarry Smith .seealso: MatGetSubMatrices(), MatGetSubMatrix(), MatDestroyMatrices(), MatConvert()
127ceb03754SKris Buschelman E*/
128dfe085dbSJed Brown typedef enum {MAT_INITIAL_MATRIX,MAT_REUSE_MATRIX,MAT_IGNORE_MATRIX} MatReuse;
129ceb03754SKris Buschelman 
1305494a064SHong Zhang /*E
1315494a064SHong Zhang     MatGetSubMatrixOption - Indicates if matrices obtained from a call to MatGetSubMatrices()
132829201f2SHong Zhang      include the matrix values. Currently it is only used by MatGetSeqNonzerostructure().
1335494a064SHong Zhang 
1345494a064SHong Zhang     Level: beginner
1355494a064SHong Zhang 
136829201f2SHong Zhang .seealso: MatGetSeqNonzerostructure()
1375494a064SHong Zhang E*/
1385494a064SHong Zhang typedef enum {MAT_DO_NOT_GET_VALUES,MAT_GET_VALUES} MatGetSubMatrixOption;
1395494a064SHong Zhang 
140e5bd5246SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInitializePackage(const char[]);
141c06d978dSMatthew Knepley 
142f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreate(MPI_Comm,Mat*);
143a69a960dSMatthew Knepley PetscPolymorphicFunction(MatCreate,(MPI_Comm comm),(comm,&A),Mat,A)
144a69a960dSMatthew Knepley PetscPolymorphicFunction(MatCreate,(),(PETSC_COMM_WORLD,&A),Mat,A)
145f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetSizes(Mat,PetscInt,PetscInt,PetscInt,PetscInt);
146a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetType(Mat,const MatType);
147be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetFromOptions(Mat);
148be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUpPreallocation(Mat);
149be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRegisterAll(const char[]);
150be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat));
15130de9b25SBarry Smith 
152f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetOptionsPrefix(Mat,const char[]);
153f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAppendOptionsPrefix(Mat,const char[]);
154f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOptionsPrefix(Mat,const char*[]);
155f69a0ea3SMatthew Knepley 
15630de9b25SBarry Smith /*MC
15730de9b25SBarry Smith    MatRegisterDynamic - Adds a new matrix type
15830de9b25SBarry Smith 
15930de9b25SBarry Smith    Synopsis:
160c1ac3661SBarry Smith    PetscErrorCode MatRegisterDynamic(char *name,char *path,char *name_create,PetscErrorCode (*routine_create)(Mat))
16130de9b25SBarry Smith 
16230de9b25SBarry Smith    Not Collective
16330de9b25SBarry Smith 
16430de9b25SBarry Smith    Input Parameters:
16530de9b25SBarry Smith +  name - name of a new user-defined matrix type
16630de9b25SBarry Smith .  path - path (either absolute or relative) the library containing this solver
16730de9b25SBarry Smith .  name_create - name of routine to create method context
16830de9b25SBarry Smith -  routine_create - routine to create method context
16930de9b25SBarry Smith 
17030de9b25SBarry Smith    Notes:
17130de9b25SBarry Smith    MatRegisterDynamic() may be called multiple times to add several user-defined solvers.
17230de9b25SBarry Smith 
17330de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
17430de9b25SBarry Smith    is ignored.
17530de9b25SBarry Smith 
17630de9b25SBarry Smith    Sample usage:
17730de9b25SBarry Smith .vb
17830de9b25SBarry Smith    MatRegisterDynamic("my_mat",/home/username/my_lib/lib/libO/solaris/mylib.a,
17930de9b25SBarry Smith                "MyMatCreate",MyMatCreate);
18030de9b25SBarry Smith .ve
18130de9b25SBarry Smith 
18230de9b25SBarry Smith    Then, your solver can be chosen with the procedural interface via
18330de9b25SBarry Smith $     MatSetType(Mat,"my_mat")
18430de9b25SBarry Smith    or at runtime via the option
18530de9b25SBarry Smith $     -mat_type my_mat
18630de9b25SBarry Smith 
18730de9b25SBarry Smith    Level: advanced
18830de9b25SBarry Smith 
189ab901514SBarry Smith    Notes: ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values.
19030de9b25SBarry Smith          If your function is not being put into a shared library then use VecRegister() instead
19130de9b25SBarry Smith 
19230de9b25SBarry Smith .keywords: Mat, register
19330de9b25SBarry Smith 
19430de9b25SBarry Smith .seealso: MatRegisterAll(), MatRegisterDestroy()
19530de9b25SBarry Smith 
19630de9b25SBarry Smith M*/
197273d9f13SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
198273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,0)
199273d9f13SBarry Smith #else
200273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,d)
20130de9b25SBarry Smith #endif
20230de9b25SBarry Smith 
203273d9f13SBarry Smith extern PetscTruth MatRegisterAllCalled;
204b0a32e0cSBarry Smith extern PetscFList MatList;
205*b022a5c1SBarry Smith extern PetscFList MatColoringList;
206*b022a5c1SBarry Smith extern PetscFList MatPartitioningList;
20728988994SBarry Smith 
2083b224e63SBarry Smith /*E
2093b224e63SBarry Smith     MatStructure - Indicates if the matrix has the same nonzero structure
2103b224e63SBarry Smith 
2113b224e63SBarry Smith     Level: beginner
2123b224e63SBarry Smith 
2133b224e63SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
2143b224e63SBarry Smith 
2153b224e63SBarry Smith .seealso: MatCopy(), KSPSetOperators(), PCSetOperators()
2163b224e63SBarry Smith E*/
2173b224e63SBarry Smith typedef enum {SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER,SUBSET_NONZERO_PATTERN} MatStructure;
2183b224e63SBarry Smith 
219be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqDense(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*);
220be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIDense(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*);
221be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
222ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,nz,nnz,&A),Mat,A)
223ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,&A),Mat,A)
224ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,0,nnz,&A),Mat,A)
225ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,&A),Mat,A)
226ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,A))
227ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,m,n,0,nnz,A))
228ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,A))
229be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
230ba966639SSatish 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)
231ba966639SSatish 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)
232ba966639SSatish 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)
233ba966639SSatish 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)
234ba966639SSatish 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))
235ba966639SSatish 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))
236ba966639SSatish 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))
237ba966639SSatish 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)
238ba966639SSatish 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)
239ba966639SSatish 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)
240ba966639SSatish 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)
241ba966639SSatish 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))
242ba966639SSatish 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))
243ba966639SSatish 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))
24463ab93bfSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *);
2458d7a6e47SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJWithSplitArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],PetscInt[],PetscInt[],PetscScalar[],Mat*);
2468d7a6e47SBarry Smith 
247be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIRowbs(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
248be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
249ba966639SSatish 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)
250ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A)
251ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A)
252ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A)
253ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A))
254ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A))
255ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A))
256be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
257ba966639SSatish 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)
258ba966639SSatish 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)
259ba966639SSatish 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)
260ba966639SSatish 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)
261ba966639SSatish 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))
262ba966639SSatish 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))
263ba966639SSatish 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))
264ba966639SSatish 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)
265ba966639SSatish 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)
266ba966639SSatish 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)
267ba966639SSatish 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)
268ba966639SSatish 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))
269ba966639SSatish 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))
270ba966639SSatish 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))
271be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAdj(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscInt[],Mat*);
272be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
273ba966639SSatish 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)
274ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A)
275ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A)
276ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A)
277ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A))
278ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A))
279ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A))
280be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPISBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
281ba966639SSatish 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)
282ba966639SSatish 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)
283ba966639SSatish 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)
284ba966639SSatish 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)
285ba966639SSatish 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))
286ba966639SSatish 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))
287ba966639SSatish 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))
288ba966639SSatish 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)
289ba966639SSatish 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)
290ba966639SSatish 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)
291ba966639SSatish 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)
292ba966639SSatish 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))
293ba966639SSatish 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))
294ba966639SSatish 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))
295be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateShell(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,void *,Mat*);
296ba966639SSatish 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)
297ba966639SSatish Balay PetscPolymorphicFunction(MatCreateShell,(PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx),(PETSC_COMM_WORLD,m,n,M,N,ctx,&A),Mat,A)
298be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateAdic(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,void (*)(void),Mat*);
299be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateNormal(Mat,Mat*);
300ba966639SSatish Balay PetscPolymorphicFunction(MatCreateNormal,(Mat mat),(mat,&A),Mat,A)
301ab92ecdeSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateLRC(Mat,Mat,Mat,Mat*);
302284134d9SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateIS(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,ISLocalToGlobalMapping,Mat*);
3031472f72bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
3041472f72bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPICRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
3052a6744ebSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateScatter(MPI_Comm,VecScatter,Mat*);
3062a6744ebSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScatterSetVecScatter(Mat,VecScatter);
3072a6744ebSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScatterGetVecScatter(Mat,VecScatter*);
3088cdf2d9bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateBlockMat(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,Mat*);
309793850ffSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompositeAddMat(Mat,Mat);
310b52f573bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompositeMerge(Mat);
311793850ffSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateComposite(MPI_Comm,PetscInt,const Mat*,Mat*);
3126d7c1e57SBarry Smith typedef enum {MAT_COMPOSITE_ADDITIVE,MAT_COMPOSITE_MULTIPLICATIVE} MatCompositeType;
3136d7c1e57SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompositeSetType(Mat,MatCompositeType);
3146d7c1e57SBarry Smith 
3155a7f1df3SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqFFTW(MPI_Comm,PetscInt,const PetscInt[],Mat*);
316f20085c4SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateTranspose(Mat,Mat*);
317ee1cef2cSJed Brown EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSubMatrix(Mat,IS,IS,Mat*);
318ee1cef2cSJed Brown EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSubMatrixUpdate(Mat,Mat,IS,IS);
3191472f72bSBarry Smith 
3201d6018f0SLisandro Dalcin EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreatePython(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const char[],Mat*);
3211d6018f0SLisandro Dalcin EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPythonSetType(Mat,const char[]);
3221d6018f0SLisandro Dalcin 
3231d6018f0SLisandro Dalcin 
324f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUp(Mat);
325be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroy(Mat);
32621c89e3eSBarry Smith 
3270c451bc4SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConjugate(Mat);
32899cafbc1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRealPart(Mat);
32999cafbc1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatImaginaryPart(Mat);
330bfc7d00aSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetDiagonalBlock(Mat,PetscTruth*,MatReuse,Mat*);
33199cafbc1SBarry Smith 
3328ed539a5SBarry Smith /* ------------------------------------------------------------*/
333be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
334be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
33587d4246cSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesRow(Mat,PetscInt,const PetscScalar[]);
336f259bd47SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesRowLocal(Mat,PetscInt,const PetscScalar[]);
33784cb2905SBarry Smith 
3382ef4de8bSBarry Smith /*S
3392ef4de8bSBarry Smith      MatStencil - Data structure (C struct) for storing information about a single row or
3402ef4de8bSBarry Smith         column of a matrix as index on an associated grid.
3412ef4de8bSBarry Smith 
3422ef4de8bSBarry Smith    Level: beginner
3432ef4de8bSBarry Smith 
3442ef4de8bSBarry Smith   Concepts: matrix; linear operator
3452ef4de8bSBarry Smith 
346d7bd93baSBarry Smith .seealso:  MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockStencil()
3472ef4de8bSBarry Smith S*/
348435da068SBarry Smith typedef struct {
349c1ac3661SBarry Smith   PetscInt k,j,i,c;
350435da068SBarry Smith } MatStencil;
3512ef4de8bSBarry Smith 
352be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
353be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
354be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt);
355435da068SBarry Smith 
356be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetColoring(Mat,ISColoring);
357be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesAdic(Mat,void*);
358be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesAdifor(Mat,PetscInt,void*);
3593a7fca6bSBarry Smith 
360d91e6319SBarry Smith /*E
361d91e6319SBarry Smith     MatAssemblyType - Indicates if the matrix is now to be used, or if you plan
362d91e6319SBarry Smith      to continue to add values to it
363d91e6319SBarry Smith 
364d91e6319SBarry Smith     Level: beginner
365d91e6319SBarry Smith 
366d91e6319SBarry Smith .seealso: MatAssemblyBegin(), MatAssemblyEnd()
367d91e6319SBarry Smith E*/
3686d4a8577SBarry Smith typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType;
369be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssemblyBegin(Mat,MatAssemblyType);
370be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssemblyEnd(Mat,MatAssemblyType);
371be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssembled(Mat,PetscTruth*);
3724f9c727eSBarry Smith 
3731d73ed98SBarry Smith 
37430de9b25SBarry Smith 
375d91e6319SBarry Smith /*E
376d91e6319SBarry Smith     MatOption - Options that may be set for a matrix and its behavior or storage
377d91e6319SBarry Smith 
378d91e6319SBarry Smith     Level: beginner
379d91e6319SBarry Smith 
3800a835dfdSSatish Balay    Any additions/changes here MUST also be made in include/finclude/petscmat.h
381d91e6319SBarry Smith 
382d91e6319SBarry Smith .seealso: MatSetOption()
383d91e6319SBarry Smith E*/
384512a5fc5SBarry Smith typedef enum {MAT_ROW_ORIENTED,MAT_NEW_NONZERO_LOCATIONS,
3854e0d8c25SBarry Smith               MAT_SYMMETRIC,
3864e0d8c25SBarry Smith               MAT_STRUCTURALLY_SYMMETRIC,
3878047e77bSBarry Smith               MAT_NEW_DIAGONALS,MAT_IGNORE_OFF_PROC_ENTRIES,
3884e0d8c25SBarry Smith               MAT_NEW_NONZERO_LOCATION_ERR,
3894e0d8c25SBarry Smith               MAT_NEW_NONZERO_ALLOCATION_ERR,MAT_USE_HASH_TABLE,
390a9817697SBarry Smith               MAT_KEEP_NONZERO_PATTERN,MAT_IGNORE_ZERO_ENTRIES,
3914e0d8c25SBarry Smith               MAT_USE_INODES,
3924e0d8c25SBarry Smith               MAT_HERMITIAN,
3934e0d8c25SBarry Smith               MAT_SYMMETRY_ETERNAL,
3944e0d8c25SBarry Smith               MAT_USE_COMPRESSEDROW,
3954e0d8c25SBarry Smith               MAT_IGNORE_LOWER_TRIANGULAR,MAT_ERROR_LOWER_TRIANGULAR,
39628b2fa4aSMatthew Knepley               MAT_GETROW_UPPERTRIANGULAR,MAT_UNUSED_NONZERO_LOCATION_ERR,
39728b2fa4aSMatthew Knepley               NUM_MAT_OPTIONS} MatOption;
398290bbb0aSBarry Smith extern const char *MatOptions[];
3994e0d8c25SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetOption(Mat,MatOption,PetscTruth);
400a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetType(Mat,const MatType*);
401a313700dSBarry Smith PetscPolymorphicFunction(MatGetType,(Mat mat),(mat,&t),const MatType,t)
40284cb2905SBarry Smith 
403be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]);
404be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
405be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
406f5edf698SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowUpperTriangular(Mat);
407f5edf698SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRowUpperTriangular(Mat);
408be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
409be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
410be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumnVector(Mat,Vec,PetscInt);
411be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetArray(Mat,PetscScalar *[]);
412ba966639SSatish Balay PetscPolymorphicFunction(MatGetArray,(Mat mat),(mat,&a),PetscScalar*,a)
413be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreArray(Mat,PetscScalar *[]);
414be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBlockSize(Mat,PetscInt *);
415ba966639SSatish Balay PetscPolymorphicFunction(MatGetBlockSize,(Mat mat),(mat,&a),PetscInt,a)
416be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetBlockSize(Mat,PetscInt);
4177b80b807SBarry Smith 
4181620fd73SBarry Smith 
419be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMult(Mat,Vec,Vec);
42089c6957cSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultDiagonalBlock(Mat,Vec,Vec);
421be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultAdd(Mat,Vec,Vec,Vec);
422ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
423be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTranspose(Mat,Vec,Vec);
424be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsTranspose(Mat,Mat,PetscReal,PetscTruth*);
425ba966639SSatish Balay PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B,PetscReal tol),(A,B,tol,&t),PetscTruth,t)
426e37f64cdSBarry Smith PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B),(A,B,0,&t),PetscTruth,t)
4271cbb95d3SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitianTranspose(Mat,Mat,PetscReal,PetscTruth*);
428be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeAdd(Mat,Vec,Vec,Vec);
429ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultTransposeAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
430be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultConstrained(Mat,Vec,Vec);
431be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeConstrained(Mat,Vec,Vec);
4322bebee5dSHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatSolve(Mat,Mat,Mat);
433f5edf698SHong Zhang 
434d91e6319SBarry Smith /*E
435d91e6319SBarry Smith     MatDuplicateOption - Indicates if a duplicated sparse matrix should have
436d91e6319SBarry Smith   its numerical values copied over or just its nonzero structure.
437d91e6319SBarry Smith 
438d91e6319SBarry Smith     Level: beginner
439d91e6319SBarry Smith 
440d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
441d91e6319SBarry Smith 
442d91e6319SBarry Smith .seealso: MatDuplicate()
443d91e6319SBarry Smith E*/
4442e8a6d31SBarry Smith typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES} MatDuplicateOption;
4452e8a6d31SBarry Smith 
446a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert(Mat,const MatType,MatReuse,Mat*);
447a313700dSBarry Smith PetscPolymorphicFunction(MatConvert,(Mat A,const MatType t),(A,t,MAT_INITIAL_MATRIX,&a),Mat,a)
448be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDuplicate(Mat,MatDuplicateOption,Mat*);
449ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A,MatDuplicateOption o),(A,o,&a),Mat,a)
450ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A),(A,MAT_COPY_VALUES,&a),Mat,a)
45194a9d846SBarry Smith 
452cb5b572fSBarry Smith 
453be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCopy(Mat,Mat,MatStructure);
454be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatView(Mat,PetscViewer);
455be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsSymmetric(Mat,PetscReal,PetscTruth*);
456ba966639SSatish Balay PetscPolymorphicFunction(MatIsSymmetric,(Mat A,PetscReal tol),(A,tol,&t),PetscTruth,t)
457e37f64cdSBarry Smith PetscPolymorphicFunction(MatIsSymmetric,(Mat A),(A,0,&t),PetscTruth,t)
458be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsStructurallySymmetric(Mat,PetscTruth*);
459ba966639SSatish Balay PetscPolymorphicFunction(MatIsStructurallySymmetric,(Mat A),(A,&t),PetscTruth,t)
4601cbb95d3SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitian(Mat,PetscReal,PetscTruth*);
4611cbb95d3SBarry Smith PetscPolymorphicFunction(MatIsHermitian,(Mat A),(A,0,&t),PetscTruth,t)
462be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsSymmetricKnown(Mat,PetscTruth*,PetscTruth*);
463be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitianKnown(Mat,PetscTruth*,PetscTruth*);
46464c62002SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMissingDiagonal(Mat,PetscTruth *,PetscInt *);
465a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLoad(PetscViewer,const MatType,Mat*);
466a313700dSBarry Smith PetscPolymorphicFunction(MatLoad,(PetscViewer v,const MatType t),(v,t,&a),Mat,a)
4677b80b807SBarry Smith 
4688f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *);
4698f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRowIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *);
4708f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumnIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *);
4718f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreColumnIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *);
472d4fbbf0eSBarry Smith 
473d91e6319SBarry Smith /*S
474d91e6319SBarry Smith      MatInfo - Context of matrix information, used with MatGetInfo()
475d91e6319SBarry Smith 
476d91e6319SBarry Smith    In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE
477d91e6319SBarry Smith 
478d91e6319SBarry Smith    Level: intermediate
479d91e6319SBarry Smith 
480d91e6319SBarry Smith   Concepts: matrix^nonzero information
481d91e6319SBarry Smith 
482d9274352SBarry Smith .seealso:  MatGetInfo(), MatInfoType
483d91e6319SBarry Smith S*/
4844e220ebcSLois Curfman McInnes typedef struct {
485b0a32e0cSBarry Smith   PetscLogDouble block_size;                         /* block size */
486b0a32e0cSBarry Smith   PetscLogDouble nz_allocated,nz_used,nz_unneeded;   /* number of nonzeros */
487b0a32e0cSBarry Smith   PetscLogDouble memory;                             /* memory allocated */
488b0a32e0cSBarry Smith   PetscLogDouble assemblies;                         /* number of matrix assemblies called */
489b0a32e0cSBarry Smith   PetscLogDouble mallocs;                            /* number of mallocs during MatSetValues() */
490b0a32e0cSBarry Smith   PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */
491b0a32e0cSBarry Smith   PetscLogDouble factor_mallocs;                     /* number of mallocs during factorization */
4924e220ebcSLois Curfman McInnes } MatInfo;
4934e220ebcSLois Curfman McInnes 
494d9274352SBarry Smith /*E
495d9274352SBarry Smith     MatInfoType - Indicates if you want information about the local part of the matrix,
496d9274352SBarry Smith      the entire parallel matrix or the maximum over all the local parts.
497d9274352SBarry Smith 
498d9274352SBarry Smith     Level: beginner
499d9274352SBarry Smith 
500d9274352SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
501d9274352SBarry Smith 
502d9274352SBarry Smith .seealso: MatGetInfo(), MatInfo
503d9274352SBarry Smith E*/
5047b80b807SBarry Smith typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType;
505be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetInfo(Mat,MatInfoType,MatInfo*);
506be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatValid(Mat,PetscTruth*);
507be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetDiagonal(Mat,Vec);
508985db425SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMax(Mat,Vec,PetscInt[]);
509985db425SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMin(Mat,Vec,PetscInt[]);
510985db425SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMaxAbs(Mat,Vec,PetscInt[]);
511c87e5d42SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMinAbs(Mat,Vec,PetscInt[]);
51286697f06SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowSum(Mat,Vec);
513fc4dec0aSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatTranspose(Mat,MatReuse,Mat*);
514fc4dec0aSBarry Smith PetscPolymorphicFunction(MatTranspose,(Mat A),(A,MAT_INITIAL_MATRIX,&t),Mat,t)
515be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPermute(Mat,IS,IS,Mat *);
516ba966639SSatish Balay PetscPolymorphicFunction(MatPermute,(Mat A,IS is1,IS is2),(A,is1,is2,&t),Mat,t)
517be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPermuteSparsify(Mat,PetscInt,PetscReal,PetscReal,IS,IS,Mat *);
518be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalScale(Mat,Vec,Vec);
519be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalSet(Mat,Vec,InsertMode);
520be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatEqual(Mat,Mat,PetscTruth*);
521ba966639SSatish Balay PetscPolymorphicFunction(MatEqual,(Mat A,Mat B),(A,B,&t),PetscTruth,t)
522be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultEqual(Mat,Mat,PetscInt,PetscTruth*);
523be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultAddEqual(Mat,Mat,PetscInt,PetscTruth*);
524be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeEqual(Mat,Mat,PetscInt,PetscTruth*);
525be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscTruth*);
5267b80b807SBarry Smith 
527be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNorm(Mat,NormType,PetscReal *);
528ba966639SSatish Balay PetscPolymorphicFunction(MatNorm,(Mat A,NormType t),(A,t,&n),PetscReal,n)
529242f1d38SBarry Smith EXTERN PetscErrorCode MatGetColumnNorms(Mat,NormType,PetscReal *);
530be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroEntries(Mat);
531f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar);
532f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsIS(Mat,IS,PetscScalar);
533f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroColumns(Mat,PetscInt,const PetscInt [],const PetscScalar*);
534f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroColumnsIS(Mat,IS,const PetscScalar*);
5357b80b807SBarry Smith 
536be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatUseScaledForm(Mat,PetscTruth);
537be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScaleSystem(Mat,Vec,Vec);
538be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatUnScaleSystem(Mat,Vec,Vec);
5395ef9f2a5SBarry Smith 
540be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSize(Mat,PetscInt*,PetscInt*);
541be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalSize(Mat,PetscInt*,PetscInt*);
542be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRange(Mat,PetscInt*,PetscInt*);
5438ee2e534SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRanges(Mat,const PetscInt**);
544f4b8409dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRangeColumn(Mat,PetscInt*,PetscInt*);
54572ca8751SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRangesColumn(Mat,const PetscInt**);
5467b80b807SBarry Smith 
547be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
548be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroyMatrices(PetscInt,Mat *[]);
5494aa3045dSJed Brown EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrix(Mat,IS,IS,MatReuse,Mat *);
550f6d58c54SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSeqNonzeroStructure(Mat,Mat*);
551f6d58c54SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroySeqNonzeroStructure(Mat*);
5525494a064SHong Zhang 
553be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
554be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPI(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*);
555be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPISymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*);
556be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPINumeric(Mat,Mat);
557be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroy_MPIAIJ_SeqsToMPI(Mat);
558be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalMat(Mat,MatReuse,Mat*);
559be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*);
560be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,PetscInt*,Mat*);
561dd6ea824SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBrowsOfAoCols(Mat,Mat,MatReuse,PetscInt**,MatScalar**,Mat*);
56243eb5e2fSMatthew Knepley #if defined (PETSC_USE_CTABLE)
563cd116e26SSatish Balay #include "petscctable.h"
56443eb5e2fSMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetCommunicationStructs(Mat, Vec *, PetscTable *, VecScatter *);
56543eb5e2fSMatthew Knepley #else
56643eb5e2fSMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetCommunicationStructs(Mat, Vec *, PetscInt *[], VecScatter *);
56743eb5e2fSMatthew Knepley #endif
5688c7482ecSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetGhosts(Mat, PetscInt *,const PetscInt *[]);
5698efafbd8SBarry Smith 
570be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt);
5717b80b807SBarry Smith 
572be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*);
573be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultSymbolic(Mat,Mat,PetscReal,Mat*);
574be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultNumeric(Mat,Mat,Mat);
57522440eb1SKris Buschelman 
576be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*);
577be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAPSymbolic(Mat,Mat,PetscReal,Mat*);
578be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAPNumeric(Mat,Mat,Mat);
57922440eb1SKris Buschelman 
580be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTranspose(Mat,Mat,MatReuse,PetscReal,Mat*);
581be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTransposeSymbolic(Mat,Mat,PetscReal,Mat*);
582be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTransposeNumeric(Mat,Mat,Mat);
583bc011b1eSHong Zhang 
584f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAXPY(Mat,PetscScalar,Mat,MatStructure);
58504aac2b0SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAYPX(Mat,PetscScalar,Mat,MatStructure);
5861c741599SHong Zhang 
587f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScale(Mat,PetscScalar);
588f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShift(Mat,PetscScalar);
5897b80b807SBarry Smith 
590be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping);
591be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetLocalToGlobalMappingBlock(Mat,ISLocalToGlobalMapping);
592f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar);
593f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsLocalIS(Mat,IS,PetscScalar);
594be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
595be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
596052efed2SBarry Smith 
597be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStashSetInitialSize(Mat,PetscInt,PetscInt);
598be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
59990f02eecSBarry Smith 
600be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInterpolate(Mat,Vec,Vec);
6010c451bc4SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInterpolateAdd(Mat,Vec,Vec,Vec);
602ba966639SSatish Balay PetscPolymorphicSubroutine(MatInterpolateAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
603be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestrict(Mat,Vec,Vec);
604be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetVecs(Mat,Vec*,Vec*);
60569db28dcSHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRedundantMatrix(Mat,PetscInt,MPI_Comm,PetscInt,MatReuse,Mat*);
606bd481603SBarry Smith 
607bd481603SBarry Smith /*MC
6081620fd73SBarry Smith    MatSetValue - Set a single entry into a matrix.
6091620fd73SBarry Smith 
6101620fd73SBarry Smith    Not collective
6111620fd73SBarry Smith 
6121620fd73SBarry Smith    Input Parameters:
6131620fd73SBarry Smith +  m - the matrix
6141620fd73SBarry Smith .  row - the row location of the entry
6151620fd73SBarry Smith .  col - the column location of the entry
6161620fd73SBarry Smith .  value - the value to insert
6171620fd73SBarry Smith -  mode - either INSERT_VALUES or ADD_VALUES
6181620fd73SBarry Smith 
6191620fd73SBarry Smith    Notes:
6201620fd73SBarry Smith    For efficiency one should use MatSetValues() and set several or many
6211620fd73SBarry Smith    values simultaneously if possible.
6221620fd73SBarry Smith 
6231620fd73SBarry Smith    Level: beginner
6241620fd73SBarry Smith 
6251620fd73SBarry Smith .seealso: MatSetValues(), MatSetValueLocal()
6261620fd73SBarry Smith M*/
6271620fd73SBarry 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);}
6281620fd73SBarry Smith 
6291620fd73SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatGetValue(Mat v,PetscInt i,PetscInt j,PetscScalar *va) {return MatGetValues(v,1,&i,1,&j,va);}
6301620fd73SBarry Smith 
6311620fd73SBarry 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);}
6321620fd73SBarry Smith 
6331620fd73SBarry Smith /*MC
6340d2a7ff2SSatish Balay    MatPreallocateInitialize - Begins the block of code that will count the number of nonzeros per
635bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
636bd481603SBarry Smith 
637bd481603SBarry Smith    Synopsis:
638c1ac3661SBarry Smith    PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)
639bd481603SBarry Smith 
640bd481603SBarry Smith    Collective on MPI_Comm
641bd481603SBarry Smith 
642bd481603SBarry Smith    Input Parameters:
643bd481603SBarry Smith +  comm - the communicator that will share the eventually allocated matrix
644859fcb39SBarry Smith .  nrows - the number of LOCAL rows in the matrix
645859fcb39SBarry Smith -  ncols - the number of LOCAL columns in the matrix
646bd481603SBarry Smith 
647bd481603SBarry Smith    Output Parameters:
648bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
649bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
650bd481603SBarry Smith 
651bd481603SBarry Smith 
652bd481603SBarry Smith    Level: intermediate
653bd481603SBarry Smith 
654bd481603SBarry Smith    Notes:
655bd481603SBarry Smith    See the chapter in the users manual on performance for more details
656bd481603SBarry Smith 
6571d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
658bd481603SBarry Smith 
659bd481603SBarry Smith    Use MatPreallocateInitializeSymmetric() for symmetric matrices (MPISBAIJ matrices)
660bd481603SBarry Smith 
6611620fd73SBarry Smith    This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize().
6621620fd73SBarry Smith 
663bd481603SBarry Smith   Concepts: preallocation^Matrix
664bd481603SBarry Smith 
665bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
666bd481603SBarry Smith           MatPreallocateInitializeSymmetric(), MatPreallocateSymmetricSetLocal()
667bd481603SBarry Smith M*/
668c4f061fbSSatish Balay #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \
6697c922b88SBarry Smith { \
670a89bc211SBarry Smith   PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__start,__end; \
671a89bc211SBarry Smith   _4_ierr = PetscMalloc2(__nrows,PetscInt,&dnz,__nrows,PetscInt,&onz);CHKERRQ(_4_ierr); \
672a89bc211SBarry Smith   _4_ierr = PetscMemzero(dnz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
673a89bc211SBarry Smith   _4_ierr = PetscMemzero(onz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
674c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ctmp;\
675a89bc211SBarry Smith   _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows;
6767c922b88SBarry Smith 
677bd481603SBarry Smith /*MC
6780d2a7ff2SSatish Balay    MatPreallocateSymmetricInitialize - Begins the block of code that will count the number of nonzeros per
679bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
680bd481603SBarry Smith 
681bd481603SBarry Smith    Synopsis:
682c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)
683bd481603SBarry Smith 
684bd481603SBarry Smith    Collective on MPI_Comm
685bd481603SBarry Smith 
686bd481603SBarry Smith    Input Parameters:
687bd481603SBarry Smith +  comm - the communicator that will share the eventually allocated matrix
688859fcb39SBarry Smith .  nrows - the number of LOCAL rows in the matrix
689859fcb39SBarry Smith -  ncols - the number of LOCAL columns in the matrix
690bd481603SBarry Smith 
691bd481603SBarry Smith    Output Parameters:
692bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
693bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
694bd481603SBarry Smith 
695bd481603SBarry Smith 
696bd481603SBarry Smith    Level: intermediate
697bd481603SBarry Smith 
698bd481603SBarry Smith    Notes:
699bd481603SBarry Smith    See the chapter in the users manual on performance for more details
700bd481603SBarry Smith 
7011d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
702bd481603SBarry Smith 
7031620fd73SBarry Smith    This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize().
7041620fd73SBarry Smith 
705bd481603SBarry Smith   Concepts: preallocation^Matrix
706bd481603SBarry Smith 
707bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
708bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal()
709bd481603SBarry Smith M*/
710222b16d4SBarry Smith #define MatPreallocateSymmetricInitialize(comm,nrows,ncols,dnz,onz) 0; \
711222b16d4SBarry Smith { \
712a89bc211SBarry Smith   PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__end; \
713a89bc211SBarry Smith   _4_ierr = PetscMalloc2(__nrows,PetscInt,&dnz,__nrows,PetscInt,&onz);CHKERRQ(_4_ierr); \
714a89bc211SBarry Smith   _4_ierr = PetscMemzero(dnz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
715a89bc211SBarry Smith   _4_ierr = PetscMemzero(onz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
716c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr);\
717a89bc211SBarry Smith   _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows;
718222b16d4SBarry Smith 
719bd481603SBarry Smith /*MC
720bd481603SBarry Smith    MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
721bd481603SBarry Smith        inserted using a local number of the rows and columns
722bd481603SBarry Smith 
723bd481603SBarry Smith    Synopsis:
724c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
725bd481603SBarry Smith 
726bd481603SBarry Smith    Not Collective
727bd481603SBarry Smith 
728bd481603SBarry Smith    Input Parameters:
729bd481603SBarry Smith +  map - the mapping between local numbering and global numbering
730bd481603SBarry Smith .  nrows - the number of rows indicated
7311d73ed98SBarry Smith .  rows - the indices of the rows
732bd481603SBarry Smith .  ncols - the number of columns in the matrix
733bd481603SBarry Smith .  cols - the columns indicated
734bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
735bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
736bd481603SBarry Smith 
737bd481603SBarry Smith 
738bd481603SBarry Smith    Level: intermediate
739bd481603SBarry Smith 
740bd481603SBarry Smith    Notes:
741bd481603SBarry Smith    See the chapter in the users manual on performance for more details
742bd481603SBarry Smith 
7431d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
744bd481603SBarry Smith 
745bd481603SBarry Smith   Concepts: preallocation^Matrix
746bd481603SBarry Smith 
747bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
748bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal()
749bd481603SBarry Smith M*/
750c4f061fbSSatish Balay #define MatPreallocateSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\
751c4f061fbSSatish Balay {\
752c1ac3661SBarry Smith   PetscInt __l;\
753ef66eb69SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\
754ef66eb69SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\
755c4f061fbSSatish Balay   for (__l=0;__l<nrows;__l++) {\
756ef66eb69SBarry Smith     _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
757c4f061fbSSatish Balay   }\
758c4f061fbSSatish Balay }
759c4f061fbSSatish Balay 
760bd481603SBarry Smith /*MC
761bd481603SBarry Smith    MatPreallocateSymmetricSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
762bd481603SBarry Smith        inserted using a local number of the rows and columns
763bd481603SBarry Smith 
764bd481603SBarry Smith    Synopsis:
765c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
766bd481603SBarry Smith 
767bd481603SBarry Smith    Not Collective
768bd481603SBarry Smith 
769bd481603SBarry Smith    Input Parameters:
770bd481603SBarry Smith +  map - the mapping between local numbering and global numbering
771bd481603SBarry Smith .  nrows - the number of rows indicated
7721d73ed98SBarry Smith .  rows - the indices of the rows
773bd481603SBarry Smith .  ncols - the number of columns in the matrix
774bd481603SBarry Smith .  cols - the columns indicated
775bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
776bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
777bd481603SBarry Smith 
778bd481603SBarry Smith 
779bd481603SBarry Smith    Level: intermediate
780bd481603SBarry Smith 
781bd481603SBarry Smith    Notes:
782bd481603SBarry Smith    See the chapter in the users manual on performance for more details
783bd481603SBarry Smith 
784bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
785bd481603SBarry Smith 
786bd481603SBarry Smith   Concepts: preallocation^Matrix
787bd481603SBarry Smith 
788bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
789bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
790bd481603SBarry Smith M*/
791d3d32019SBarry Smith #define MatPreallocateSymmetricSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\
792d3d32019SBarry Smith {\
793c1ac3661SBarry Smith   PetscInt __l;\
794d3d32019SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\
795d3d32019SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\
796d3d32019SBarry Smith   for (__l=0;__l<nrows;__l++) {\
797d3d32019SBarry Smith     _4_ierr = MatPreallocateSymmetricSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
798d3d32019SBarry Smith   }\
799d3d32019SBarry Smith }
800d3d32019SBarry Smith 
801bd481603SBarry Smith /*MC
802bd481603SBarry Smith    MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
803bd481603SBarry Smith        inserted using a local number of the rows and columns
804bd481603SBarry Smith 
805bd481603SBarry Smith    Synopsis:
806c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
807bd481603SBarry Smith 
808bd481603SBarry Smith    Not Collective
809bd481603SBarry Smith 
810bd481603SBarry Smith    Input Parameters:
81164075487SBarry Smith +  row - the row
812bd481603SBarry Smith .  ncols - the number of columns in the matrix
813a50f8bf6SHong Zhang -  cols - the columns indicated
814a50f8bf6SHong Zhang 
815a50f8bf6SHong Zhang    Output Parameters:
816a50f8bf6SHong Zhang +  dnz - the array that will be passed to the matrix preallocation routines
817bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
818bd481603SBarry Smith 
819bd481603SBarry Smith 
820bd481603SBarry Smith    Level: intermediate
821bd481603SBarry Smith 
822bd481603SBarry Smith    Notes:
823bd481603SBarry Smith    See the chapter in the users manual on performance for more details
824bd481603SBarry Smith 
825bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
826bd481603SBarry Smith 
8271620fd73SBarry Smith    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
8281620fd73SBarry Smith 
829bd481603SBarry Smith   Concepts: preallocation^Matrix
830bd481603SBarry Smith 
831bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
832bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
833bd481603SBarry Smith M*/
834c4f061fbSSatish Balay #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\
835c1ac3661SBarry Smith { PetscInt __i; \
8368ee2e534SBarry Smith   if (row < __rstart) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Trying to set preallocation for row %D less than first local row %D",row,__rstart);\
837a89bc211SBarry 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);\
8387c922b88SBarry Smith   for (__i=0; __i<nc; __i++) {\
83964075487SBarry Smith     if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \
8407cc688d7SBarry Smith     else dnz[row - __rstart]++;\
8417c922b88SBarry Smith   }\
8427c922b88SBarry Smith }
8437c922b88SBarry Smith 
844bd481603SBarry Smith /*MC
845bd481603SBarry Smith    MatPreallocateSymmetricSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
846bd481603SBarry Smith        inserted using a local number of the rows and columns
847bd481603SBarry Smith 
848bd481603SBarry Smith    Synopsis:
849c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
850bd481603SBarry Smith 
851bd481603SBarry Smith    Not Collective
852bd481603SBarry Smith 
853bd481603SBarry Smith    Input Parameters:
854bd481603SBarry Smith +  nrows - the number of rows indicated
8551d73ed98SBarry Smith .  rows - the indices of the rows
856bd481603SBarry Smith .  ncols - the number of columns in the matrix
857bd481603SBarry Smith .  cols - the columns indicated
858bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
859bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
860bd481603SBarry Smith 
861bd481603SBarry Smith 
862bd481603SBarry Smith    Level: intermediate
863bd481603SBarry Smith 
864bd481603SBarry Smith    Notes:
865bd481603SBarry Smith    See the chapter in the users manual on performance for more details
866bd481603SBarry Smith 
867bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
868bd481603SBarry Smith 
8691620fd73SBarry Smith    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
8701620fd73SBarry Smith 
871bd481603SBarry Smith   Concepts: preallocation^Matrix
872bd481603SBarry Smith 
873bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
874bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
875bd481603SBarry Smith M*/
876d3d32019SBarry Smith #define MatPreallocateSymmetricSet(row,nc,cols,dnz,onz) 0;\
877c1ac3661SBarry Smith { PetscInt __i; \
878d3d32019SBarry Smith   for (__i=0; __i<nc; __i++) {\
879d3d32019SBarry Smith     if (cols[__i] >= __end) onz[row - __rstart]++; \
880d3d32019SBarry Smith     else if (cols[__i] >= row) dnz[row - __rstart]++;\
881d3d32019SBarry Smith   }\
882d3d32019SBarry Smith }
883d3d32019SBarry Smith 
884bd481603SBarry Smith /*MC
88516371a99SBarry Smith    MatPreallocateLocation -  An alternative to MatPreallocationSet() that puts the nonzero locations into the matrix if it exists
88616371a99SBarry Smith 
88716371a99SBarry Smith    Synopsis:
88816371a99SBarry Smith    PetscErrorCode MatPreallocateLocations(Mat A,PetscInt row,PetscInt ncols,PetscInt *cols,PetscInt *dnz,PetscInt *onz)
88916371a99SBarry Smith 
89016371a99SBarry Smith    Not Collective
89116371a99SBarry Smith 
89216371a99SBarry Smith    Input Parameters:
89316371a99SBarry Smith .  A - matrix
89416371a99SBarry Smith .  row - row where values exist (must be local to this process)
89516371a99SBarry Smith .  ncols - number of columns
89616371a99SBarry Smith .  cols - columns with nonzeros
89716371a99SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
89816371a99SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
89916371a99SBarry Smith 
90016371a99SBarry Smith 
90116371a99SBarry Smith    Level: intermediate
90216371a99SBarry Smith 
90316371a99SBarry Smith    Notes:
90416371a99SBarry Smith    See the chapter in the users manual on performance for more details
90516371a99SBarry Smith 
90616371a99SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
90716371a99SBarry Smith 
90816371a99SBarry Smith    This is a MACRO not a function because it uses a bunch of variables private to the MatPreallocation.... routines.
90916371a99SBarry Smith 
91016371a99SBarry Smith   Concepts: preallocation^Matrix
91116371a99SBarry Smith 
91216371a99SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
91316371a99SBarry Smith           MatPreallocateSymmetricInitialize(), MatPreallocateSymmetricSetLocal()
91416371a99SBarry Smith M*/
91516371a99SBarry 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);}
91616371a99SBarry Smith 
91716371a99SBarry Smith 
91816371a99SBarry Smith /*MC
9190d2a7ff2SSatish Balay    MatPreallocateFinalize - Ends the block of code that will count the number of nonzeros per
920bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
921bd481603SBarry Smith 
922bd481603SBarry Smith    Synopsis:
923c1ac3661SBarry Smith    PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz)
924bd481603SBarry Smith 
925bd481603SBarry Smith    Collective on MPI_Comm
926bd481603SBarry Smith 
927bd481603SBarry Smith    Input Parameters:
92816371a99SBarry Smith +  dnz - the array that was be passed to the matrix preallocation routines
929bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
930bd481603SBarry Smith 
931bd481603SBarry Smith 
932bd481603SBarry Smith    Level: intermediate
933bd481603SBarry Smith 
934bd481603SBarry Smith    Notes:
935bd481603SBarry Smith    See the chapter in the users manual on performance for more details
936bd481603SBarry Smith 
937bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
938bd481603SBarry Smith 
9391620fd73SBarry Smith    This is a MACRO not a function because it closes the { started in MatPreallocateInitialize().
9401620fd73SBarry Smith 
941bd481603SBarry Smith   Concepts: preallocation^Matrix
942bd481603SBarry Smith 
943bd481603SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
944bd481603SBarry Smith           MatPreallocateSymmetricInitialize(), MatPreallocateSymmetricSetLocal()
945bd481603SBarry Smith M*/
946a89bc211SBarry Smith #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree2(dnz,onz);CHKERRQ(_4_ierr);}
9477c922b88SBarry Smith 
948bd481603SBarry Smith 
949bd481603SBarry Smith 
9507b80b807SBarry Smith /* Routines unique to particular data structures */
951be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellGetContext(Mat,void **);
952ba966639SSatish Balay PetscPolymorphicFunction(MatShellGetContext,(Mat A),(A,&t),void*,t)
953698d4c6aSKris Buschelman 
954be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInodeAdjustForInodes(Mat,IS*,IS*);
955be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *);
956698d4c6aSKris Buschelman 
957be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetColumnIndices(Mat,PetscInt[]);
958be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetColumnIndices(Mat,PetscInt[]);
959be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
960c75a6043SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
961c75a6043SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
9627b80b807SBarry Smith 
963a96a251dSBarry Smith #define MAT_SKIP_ALLOCATION -4
964a96a251dSBarry Smith 
965be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
966ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz))
967be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
968ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqSBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz))
969be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]);
970ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqAIJSetPreallocation,(Mat A,const PetscInt nnz[]),(A,0,nnz))
971273d9f13SBarry Smith 
972be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
973ba966639SSatish Balay PetscPolymorphicSubroutine(MatMPIBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[],const PetscInt onz[]),(A,bs,0,nnz,0,onz))
974be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
975be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
97653803f26SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []);
977725b52f3SLisandro Dalcin EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
978be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]);
979be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
980be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]);
981be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIDenseSetPreallocation(Mat,PetscScalar[]);
982284134d9SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation(Mat,PetscScalar[]);
983be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIRowbsSetPreallocation(Mat,PetscInt,const PetscInt[]);
984be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,PetscInt*[]);
985be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,PetscInt*[]);
986be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAdicSetLocalFunction(Mat,void (*)(void));
987273d9f13SBarry Smith 
988be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetLDA(Mat,PetscInt);
989ab92ecdeSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDenseGetLocalMatrix(Mat,Mat*);
9901b807ce4Svictorle 
991be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStoreValues(Mat);
992be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRetrieveValues(Mat);
9932e8a6d31SBarry Smith 
994be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDAADSetCtx(Mat,void*);
9953a7fca6bSBarry Smith 
9967b80b807SBarry Smith /*
9977b80b807SBarry Smith   These routines are not usually accessed directly, rather solving is
99894b7f48cSBarry Smith   done through the KSP and PC interfaces.
9997b80b807SBarry Smith */
10007b80b807SBarry Smith 
1001d9274352SBarry Smith /*E
1002d9274352SBarry Smith     MatOrderingType - String with the name of a PETSc matrix ordering or the creation function
1003d9274352SBarry Smith        with an optional dynamic library name, for example
1004d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:orderingcreate()
1005d9274352SBarry Smith 
1006d9274352SBarry Smith    Level: beginner
1007d9274352SBarry Smith 
1008e5a9bf91SBarry Smith    Cannot use const because the PC objects manipulate the string
1009e5a9bf91SBarry Smith 
1010d9274352SBarry Smith .seealso: MatGetOrdering()
1011d9274352SBarry Smith E*/
10123eaad2caSSatish Balay #define MatOrderingType char*
1013b12f92e5SBarry Smith #define MATORDERING_NATURAL     "natural"
1014b12f92e5SBarry Smith #define MATORDERING_ND          "nd"
1015b12f92e5SBarry Smith #define MATORDERING_1WD         "1wd"
1016b12f92e5SBarry Smith #define MATORDERING_RCM         "rcm"
1017b12f92e5SBarry Smith #define MATORDERING_QMD         "qmd"
1018b12f92e5SBarry Smith #define MATORDERING_ROWLENGTH   "rowlength"
1019e106eecfSBarry Smith #define MATORDERING_DSC_ND      "dsc_nd"         /* these three are only for DSCPACK, see its documentation for details */
102062152c8bSBarry Smith #define MATORDERING_DSC_MMD     "dsc_mmd"
102162152c8bSBarry Smith #define MATORDERING_DSC_MDF     "dsc_mdf"
1022c06d978dSMatthew Knepley #define MATORDERING_CONSTRAINED "constrained"
1023c06d978dSMatthew Knepley #define MATORDERING_IDENTITY    "identity"
1024c06d978dSMatthew Knepley #define MATORDERING_REVERSE     "reverse"
10258fa24674SBarry Smith #define MATORDERING_FLOW        "flow"
1026b12f92e5SBarry Smith 
1027f2416520SSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOrdering(Mat,const MatOrderingType,IS*,IS*);
1028be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOrderingList(PetscFList *list);
1029f2416520SSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,const MatOrderingType,IS*,IS*));
103030de9b25SBarry Smith 
103130de9b25SBarry Smith /*MC
103230de9b25SBarry Smith    MatOrderingRegisterDynamic - Adds a new sparse matrix ordering to the
103330de9b25SBarry Smith                                matrix package.
103430de9b25SBarry Smith 
103530de9b25SBarry Smith    Synopsis:
1036c1ac3661SBarry Smith    PetscErrorCode MatOrderingRegisterDynamic(char *name_ordering,char *path,char *name_create,PetscErrorCode (*routine_create)(MatOrdering))
103730de9b25SBarry Smith 
103830de9b25SBarry Smith    Not Collective
103930de9b25SBarry Smith 
104030de9b25SBarry Smith    Input Parameters:
104130de9b25SBarry Smith +  sname - name of ordering (for example MATORDERING_ND)
104230de9b25SBarry Smith .  path - location of library where creation routine is
104330de9b25SBarry Smith .  name - name of function that creates the ordering type,a string
104430de9b25SBarry Smith -  function - function pointer that creates the ordering
104530de9b25SBarry Smith 
104630de9b25SBarry Smith    Level: developer
104730de9b25SBarry Smith 
104830de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
104930de9b25SBarry Smith    is ignored.
105030de9b25SBarry Smith 
105130de9b25SBarry Smith    Sample usage:
105230de9b25SBarry Smith .vb
105330de9b25SBarry Smith    MatOrderingRegisterDynamic("my_order",/home/username/my_lib/lib/libO/solaris/mylib.a,
105430de9b25SBarry Smith                "MyOrder",MyOrder);
105530de9b25SBarry Smith .ve
105630de9b25SBarry Smith 
105730de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
105830de9b25SBarry Smith $     MatOrderingSetType(part,"my_order)
105930de9b25SBarry Smith    or at runtime via the option
10602401956bSBarry Smith $     -pc_factor_mat_ordering_type my_order
106130de9b25SBarry Smith 
1062ab901514SBarry Smith    ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values.
106330de9b25SBarry Smith 
106430de9b25SBarry Smith .keywords: matrix, ordering, register
106530de9b25SBarry Smith 
106630de9b25SBarry Smith .seealso: MatOrderingRegisterDestroy(), MatOrderingRegisterAll()
106730de9b25SBarry Smith M*/
1068aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1069f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,0)
1070b12f92e5SBarry Smith #else
1071f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,d)
1072b12f92e5SBarry Smith #endif
107330de9b25SBarry Smith 
1074be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegisterDestroy(void);
1075be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegisterAll(const char[]);
10762bad1931SBarry Smith extern PetscTruth MatOrderingRegisterAllCalled;
1077b0a32e0cSBarry Smith extern PetscFList MatOrderingList;
1078d4fbbf0eSBarry Smith 
1079be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS);
1080a2ce50c7SBarry Smith 
1081d91e6319SBarry Smith /*S
10822401956bSBarry Smith    MatFactorInfo - Data passed into the matrix factorization routines
1083b00f7748SHong Zhang 
108461cad860SBarry Smith    In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE, that is use
108561cad860SBarry Smith $     MatFactorInfo  info(MAT_FACTORINFO_SIZE)
1086b00f7748SHong Zhang 
108715e8a5b3SHong Zhang    Notes: These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC.
1088b00f7748SHong Zhang 
108961cad860SBarry Smith       You can use MatFactorInfoInitialize() to set default values.
109061cad860SBarry Smith 
1091b00f7748SHong Zhang    Level: developer
1092b00f7748SHong Zhang 
1093d7d82daaSBarry Smith .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(),
1094d7d82daaSBarry Smith           MatFactorInfoInitialize()
1095b00f7748SHong Zhang 
1096b00f7748SHong Zhang S*/
1097b00f7748SHong Zhang typedef struct {
10980a29876aSHong Zhang   PetscReal     shiftnz;        /* scaling of identity added to matrix to prevent zero pivots */
1099fbf22428SSatish Balay   PetscReal     shiftpd;        /* if true, shift until positive pivots */
110015e8a5b3SHong Zhang   PetscReal     diagonal_fill;  /* force diagonal to fill in if initially not filled */
110185317021SBarry Smith   PetscReal     usedt;
110215e8a5b3SHong Zhang   PetscReal     dt;             /* drop tolerance */
1103b00f7748SHong Zhang   PetscReal     dtcol;          /* tolerance for pivoting */
110415e8a5b3SHong Zhang   PetscReal     dtcount;        /* maximum nonzeros to be allowed per row */
110567e28bfeSBarry Smith   PetscReal     fill;           /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */
1106348344bbSBarry Smith   PetscReal     levels;         /* ICC/ILU(levels) */
1107bcd9e38bSBarry Smith   PetscReal     pivotinblocks;  /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0
1108bcd9e38bSBarry Smith                                    factorization may be faster if do not pivot */
110962bba022SBarry Smith   PetscReal     shiftinblocks;  /* if block in block factorization has zero pivot then shift diagonal until non-singular */
111015e8a5b3SHong Zhang   PetscReal     zeropivot;      /* pivot is called zero if less than this */
111115e8a5b3SHong Zhang } MatFactorInfo;
1112ffa6d0a5SLois Curfman McInnes 
1113be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFactorInfoInitialize(MatFactorInfo*);
11140481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactor(Mat,IS,const MatFactorInfo*);
11150481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
11160481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactorNumeric(Mat,Mat,const MatFactorInfo*);
11170481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactor(Mat,IS,IS,const MatFactorInfo*);
11180481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUFactor(Mat,IS,IS,const MatFactorInfo*);
11190481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
11200481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
11210481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatICCFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
11220481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatICCFactor(Mat,IS,const MatFactorInfo*);
11230481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactorNumeric(Mat,Mat,const MatFactorInfo*);
11240481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUDTFactor(Mat,IS,IS,const MatFactorInfo*,Mat *);
11253c2a7987SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUDTFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
11263c2a7987SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUDTFactorNumeric(Mat,Mat,const MatFactorInfo*);
1127be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*);
1128be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolve(Mat,Vec,Vec);
1129be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatForwardSolve(Mat,Vec,Vec);
1130be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatBackwardSolve(Mat,Vec,Vec);
1131be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveAdd(Mat,Vec,Vec,Vec);
1132be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveTranspose(Mat,Vec,Vec);
1133be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveTransposeAdd(Mat,Vec,Vec,Vec);
1134be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolves(Mat,Vecs,Vecs);
11358ed539a5SBarry Smith 
1136be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUnfactored(Mat);
1137bb5a7306SBarry Smith 
1138d91e6319SBarry Smith /*E
1139d91e6319SBarry Smith     MatSORType - What type of (S)SOR to perform
1140bb1eb677SSatish Balay 
1141d91e6319SBarry Smith     Level: beginner
1142d91e6319SBarry Smith 
1143d9274352SBarry Smith    May be bitwise ORd together
1144d9274352SBarry Smith 
1145d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
1146d91e6319SBarry Smith 
11474e7234bfSBarry Smith    MatSORType may be bitwise ORd together, so do not change the numbers
11484e7234bfSBarry Smith 
1149a1d92eedSBarry Smith .seealso: MatRelax(), MatPBRelax()
1150d91e6319SBarry Smith E*/
1151ee50ffe9SBarry Smith typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3,
1152ee50ffe9SBarry Smith               SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8,
1153ee50ffe9SBarry Smith               SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16,
115484cb2905SBarry Smith               SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType;
1155be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRelax(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
1156be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPBRelax(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
11578ed539a5SBarry Smith 
1158d4fbbf0eSBarry Smith /*
1159639f9d9dSBarry Smith     These routines are for efficiently computing Jacobians via finite differences.
1160639f9d9dSBarry Smith */
1161b12f92e5SBarry Smith 
1162d9274352SBarry Smith /*E
1163d9274352SBarry Smith     MatColoringType - String with the name of a PETSc matrix coloring or the creation function
1164d9274352SBarry Smith        with an optional dynamic library name, for example
1165d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:coloringcreate()
1166d9274352SBarry Smith 
1167d9274352SBarry Smith    Level: beginner
1168d9274352SBarry Smith 
1169d9274352SBarry Smith .seealso: MatGetColoring()
1170d9274352SBarry Smith E*/
1171a313700dSBarry Smith #define MatColoringType char*
1172b12f92e5SBarry Smith #define MATCOLORING_NATURAL "natural"
1173b12f92e5SBarry Smith #define MATCOLORING_SL      "sl"
1174b12f92e5SBarry Smith #define MATCOLORING_LF      "lf"
1175b12f92e5SBarry Smith #define MATCOLORING_ID      "id"
1176b12f92e5SBarry Smith 
1177ab970d9bSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColoring(Mat,const MatColoringType,ISColoring*);
11782e90c967SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,MatColoringType,ISColoring *));
117930de9b25SBarry Smith 
118030de9b25SBarry Smith /*MC
118130de9b25SBarry Smith    MatColoringRegisterDynamic - Adds a new sparse matrix coloring to the
118230de9b25SBarry Smith                                matrix package.
118330de9b25SBarry Smith 
118430de9b25SBarry Smith    Synopsis:
1185c1ac3661SBarry Smith    PetscErrorCode MatColoringRegisterDynamic(char *name_coloring,char *path,char *name_create,PetscErrorCode (*routine_create)(MatColoring))
118630de9b25SBarry Smith 
118730de9b25SBarry Smith    Not Collective
118830de9b25SBarry Smith 
118930de9b25SBarry Smith    Input Parameters:
119030de9b25SBarry Smith +  sname - name of Coloring (for example MATCOLORING_SL)
119130de9b25SBarry Smith .  path - location of library where creation routine is
119230de9b25SBarry Smith .  name - name of function that creates the Coloring type, a string
119330de9b25SBarry Smith -  function - function pointer that creates the coloring
119430de9b25SBarry Smith 
119530de9b25SBarry Smith    Level: developer
119630de9b25SBarry Smith 
119730de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
119830de9b25SBarry Smith    is ignored.
119930de9b25SBarry Smith 
120030de9b25SBarry Smith    Sample usage:
120130de9b25SBarry Smith .vb
120230de9b25SBarry Smith    MatColoringRegisterDynamic("my_color",/home/username/my_lib/lib/libO/solaris/mylib.a,
120330de9b25SBarry Smith                "MyColor",MyColor);
120430de9b25SBarry Smith .ve
120530de9b25SBarry Smith 
120630de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
120730de9b25SBarry Smith $     MatColoringSetType(part,"my_color")
120830de9b25SBarry Smith    or at runtime via the option
120930de9b25SBarry Smith $     -mat_coloring_type my_color
121030de9b25SBarry Smith 
1211ab901514SBarry Smith    $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
121230de9b25SBarry Smith 
121330de9b25SBarry Smith .keywords: matrix, Coloring, register
121430de9b25SBarry Smith 
121530de9b25SBarry Smith .seealso: MatColoringRegisterDestroy(), MatColoringRegisterAll()
121630de9b25SBarry Smith M*/
1217aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1218f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,0)
1219b12f92e5SBarry Smith #else
1220f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,d)
1221b12f92e5SBarry Smith #endif
122230de9b25SBarry Smith 
12232bad1931SBarry Smith extern PetscTruth MatColoringRegisterAllCalled;
1224f1144a30SSatish Balay 
1225be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegisterAll(const char[]);
1226be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegisterDestroy(void);
1227be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*);
1228639f9d9dSBarry Smith 
1229d9274352SBarry Smith /*S
1230d9274352SBarry Smith      MatFDColoring - Object for computing a sparse Jacobian via finite differences
1231d9274352SBarry Smith         and coloring
1232639f9d9dSBarry Smith 
1233d9274352SBarry Smith    Level: beginner
1234d9274352SBarry Smith 
1235d9274352SBarry Smith   Concepts: coloring, sparse Jacobian, finite differences
1236d9274352SBarry Smith 
1237d9274352SBarry Smith .seealso:  MatFDColoringCreate()
1238d9274352SBarry Smith S*/
1239e2a1c21fSSatish Balay typedef struct _p_MatFDColoring* MatFDColoring;
1240639f9d9dSBarry Smith 
1241be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringCreate(Mat,ISColoring,MatFDColoring *);
1242be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringDestroy(MatFDColoring);
1243be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringView(MatFDColoring,PetscViewer);
1244be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*);
12454a9d489dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**);
1246be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal);
1247be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFromOptions(MatFDColoring);
1248be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApply(Mat,MatFDColoring,Vec,MatStructure*,void *);
1249be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApplyTS(Mat,MatFDColoring,PetscReal,Vec,MatStructure*,void *);
1250be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetF(MatFDColoring,Vec);
1251be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,PetscInt*[]);
1252639f9d9dSBarry Smith /*
12530752156aSBarry Smith     These routines are for partitioning matrices: currently used only
12543eda8832SBarry Smith   for adjacency matrix, MatCreateMPIAdj().
12550752156aSBarry Smith */
1256ca161407SBarry Smith 
1257d9274352SBarry Smith /*S
1258d9274352SBarry Smith      MatPartitioning - Object for managing the partitioning of a matrix or graph
1259d9274352SBarry Smith 
1260d9274352SBarry Smith    Level: beginner
1261d9274352SBarry Smith 
1262d9274352SBarry Smith   Concepts: partitioning
1263d9274352SBarry Smith 
1264743a1026Svictorle .seealso:  MatPartitioningCreate(), MatPartitioningType
1265d9274352SBarry Smith S*/
126691e9ee9fSBarry Smith typedef struct _p_MatPartitioning* MatPartitioning;
1267d9274352SBarry Smith 
1268d9274352SBarry Smith /*E
12695bc6e7ccSvictorle     MatPartitioningType - String with the name of a PETSc matrix partitioning or the creation function
1270d9274352SBarry Smith        with an optional dynamic library name, for example
1271d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:partitioningcreate()
1272d9274352SBarry Smith 
1273d9274352SBarry Smith    Level: beginner
1274d9274352SBarry Smith 
1275b547a13bSvictorle .seealso: MatPartitioningCreate(), MatPartitioning
1276d9274352SBarry Smith E*/
1277a313700dSBarry Smith #define MatPartitioningType char*
12788ba1e511SMatthew Knepley #define MAT_PARTITIONING_CURRENT  "current"
1279c2ec2a73SHong Zhang #define MAT_PARTITIONING_SQUARE   "square"
12808ba1e511SMatthew Knepley #define MAT_PARTITIONING_PARMETIS "parmetis"
128171306c60Sjroman #define MAT_PARTITIONING_CHACO    "chaco"
128271306c60Sjroman #define MAT_PARTITIONING_JOSTLE   "jostle"
128371306c60Sjroman #define MAT_PARTITIONING_PARTY    "party"
128471306c60Sjroman #define MAT_PARTITIONING_SCOTCH   "scotch"
128571306c60Sjroman 
1286ca161407SBarry Smith 
1287be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningCreate(MPI_Comm,MatPartitioning*);
1288a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetType(MatPartitioning,const MatPartitioningType);
1289be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetNParts(MatPartitioning,PetscInt);
1290be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetAdjacency(MatPartitioning,Mat);
1291be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]);
1292be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []);
1293be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningApply(MatPartitioning,IS*);
1294be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningDestroy(MatPartitioning);
12952aabb6bbSBarry Smith 
1296be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatPartitioning));
129730de9b25SBarry Smith 
129830de9b25SBarry Smith /*MC
129930de9b25SBarry Smith    MatPartitioningRegisterDynamic - Adds a new sparse matrix partitioning to the
130030de9b25SBarry Smith    matrix package.
130130de9b25SBarry Smith 
130230de9b25SBarry Smith    Synopsis:
1303c1ac3661SBarry Smith    PetscErrorCode MatPartitioningRegisterDynamic(char *name_partitioning,char *path,char *name_create,PetscErrorCode (*routine_create)(MatPartitioning))
130430de9b25SBarry Smith 
130530de9b25SBarry Smith    Not Collective
130630de9b25SBarry Smith 
130730de9b25SBarry Smith    Input Parameters:
130830de9b25SBarry Smith +  sname - name of partitioning (for example MAT_PARTITIONING_CURRENT) or parmetis
130930de9b25SBarry Smith .  path - location of library where creation routine is
131030de9b25SBarry Smith .  name - name of function that creates the partitioning type, a string
131130de9b25SBarry Smith -  function - function pointer that creates the partitioning type
131230de9b25SBarry Smith 
131330de9b25SBarry Smith    Level: developer
131430de9b25SBarry Smith 
131530de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
131630de9b25SBarry Smith    is ignored.
131730de9b25SBarry Smith 
131830de9b25SBarry Smith    Sample usage:
131930de9b25SBarry Smith .vb
132030de9b25SBarry Smith    MatPartitioningRegisterDynamic("my_part",/home/username/my_lib/lib/libO/solaris/mylib.a,
132130de9b25SBarry Smith                "MyPartCreate",MyPartCreate);
132230de9b25SBarry Smith .ve
132330de9b25SBarry Smith 
132430de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
132530de9b25SBarry Smith $     MatPartitioningSetType(part,"my_part")
132630de9b25SBarry Smith    or at runtime via the option
132730de9b25SBarry Smith $     -mat_partitioning_type my_part
132830de9b25SBarry Smith 
1329ab901514SBarry Smith    $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
133030de9b25SBarry Smith 
133130de9b25SBarry Smith .keywords: matrix, partitioning, register
133230de9b25SBarry Smith 
133330de9b25SBarry Smith .seealso: MatPartitioningRegisterDestroy(), MatPartitioningRegisterAll()
133430de9b25SBarry Smith M*/
1335aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1336f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,0)
13372aabb6bbSBarry Smith #else
1338f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,d)
13392aabb6bbSBarry Smith #endif
13402aabb6bbSBarry Smith 
13412bad1931SBarry Smith extern PetscTruth MatPartitioningRegisterAllCalled;
1342f1144a30SSatish Balay 
1343be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegisterAll(const char[]);
1344be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegisterDestroy(void);
13452bad1931SBarry Smith 
1346be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningView(MatPartitioning,PetscViewer);
1347be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetFromOptions(MatPartitioning);
1348a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningGetType(MatPartitioning,const MatPartitioningType*);
1349ca161407SBarry Smith 
1350be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningParmetisSetCoarseSequential(MatPartitioning);
13510e8af3caSMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningParmetisGetEdgeCut(MatPartitioning, PetscInt *);
13520752156aSBarry Smith 
1353be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningJostleSetCoarseLevel(MatPartitioning,PetscReal);
1354be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningJostleSetCoarseSequential(MatPartitioning);
135571306c60Sjroman 
1356954b3ec6SBarry Smith typedef enum {MP_CHACO_MULTILEVEL_KL,MP_CHACO_SPECTRAL,MP_CHACO_LINEAR,MP_CHACO_RANDOM, MP_CHACO_SCATTERED} MPChacoGlobalType;
1357be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetGlobal(MatPartitioning, MPChacoGlobalType);
135871306c60Sjroman typedef enum { MP_CHACO_KERNIGHAN_LIN, MP_CHACO_NONE } MPChacoLocalType;
1359be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetLocal(MatPartitioning, MPChacoLocalType);
1360be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal);
136171306c60Sjroman typedef enum { MP_CHACO_LANCZOS, MP_CHACO_RQI_SYMMLQ } MPChacoEigenType;
1362be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType);
1363be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenTol(MatPartitioning, PetscReal);
1364be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenNumber(MatPartitioning, PetscInt);
136571306c60Sjroman 
136671306c60Sjroman #define MP_PARTY_OPT "opt"
136771306c60Sjroman #define MP_PARTY_LIN "lin"
136871306c60Sjroman #define MP_PARTY_SCA "sca"
136971306c60Sjroman #define MP_PARTY_RAN "ran"
137071306c60Sjroman #define MP_PARTY_GBF "gbf"
137171306c60Sjroman #define MP_PARTY_GCF "gcf"
137271306c60Sjroman #define MP_PARTY_BUB "bub"
137371306c60Sjroman #define MP_PARTY_DEF "def"
1374be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetGlobal(MatPartitioning, const char*);
137571306c60Sjroman #define MP_PARTY_HELPFUL_SETS "hs"
137671306c60Sjroman #define MP_PARTY_KERNIGHAN_LIN "kl"
137771306c60Sjroman #define MP_PARTY_NONE "no"
1378be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetLocal(MatPartitioning, const char*);
1379be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal);
1380be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetBipart(MatPartitioning,PetscTruth);
1381be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscTruth);
138271306c60Sjroman 
138371306c60Sjroman typedef enum { MP_SCOTCH_GREEDY, MP_SCOTCH_GPS, MP_SCOTCH_GR_GPS } MPScotchGlobalType;
1384be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetArch(MatPartitioning,const char*);
1385be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetMultilevel(MatPartitioning);
1386be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetGlobal(MatPartitioning,MPScotchGlobalType);
1387be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetCoarseLevel(MatPartitioning,PetscReal);
1388be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetHostList(MatPartitioning,const char*);
138971306c60Sjroman typedef enum { MP_SCOTCH_KERNIGHAN_LIN, MP_SCOTCH_NONE } MPScotchLocalType;
1390be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetLocal(MatPartitioning,MPScotchLocalType);
1391be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetMapping(MatPartitioning);
1392be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetStrategy(MatPartitioning,char*);
139371306c60Sjroman 
13940752156aSBarry Smith /*
13950a835dfdSSatish Balay     If you add entries here you must also add them to finclude/petscmat.h
1396d4fbbf0eSBarry Smith */
13971c1c02c0SLois Curfman McInnes typedef enum { MATOP_SET_VALUES=0,
13981c1c02c0SLois Curfman McInnes                MATOP_GET_ROW=1,
13991c1c02c0SLois Curfman McInnes                MATOP_RESTORE_ROW=2,
14001c1c02c0SLois Curfman McInnes                MATOP_MULT=3,
14011c1c02c0SLois Curfman McInnes                MATOP_MULT_ADD=4,
14027c922b88SBarry Smith                MATOP_MULT_TRANSPOSE=5,
14037c922b88SBarry Smith                MATOP_MULT_TRANSPOSE_ADD=6,
14041c1c02c0SLois Curfman McInnes                MATOP_SOLVE=7,
14051c1c02c0SLois Curfman McInnes                MATOP_SOLVE_ADD=8,
14067c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE=9,
14077c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE_ADD=10,
14081c1c02c0SLois Curfman McInnes                MATOP_LUFACTOR=11,
14091c1c02c0SLois Curfman McInnes                MATOP_CHOLESKYFACTOR=12,
14101c1c02c0SLois Curfman McInnes                MATOP_RELAX=13,
14111c1c02c0SLois Curfman McInnes                MATOP_TRANSPOSE=14,
14121c1c02c0SLois Curfman McInnes                MATOP_GETINFO=15,
14131c1c02c0SLois Curfman McInnes                MATOP_EQUAL=16,
14141c1c02c0SLois Curfman McInnes                MATOP_GET_DIAGONAL=17,
14151c1c02c0SLois Curfman McInnes                MATOP_DIAGONAL_SCALE=18,
14161c1c02c0SLois Curfman McInnes                MATOP_NORM=19,
14171c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_BEGIN=20,
14181c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_END=21,
1419d519adbfSMatthew Knepley                MATOP_SET_OPTION=22,
1420d519adbfSMatthew Knepley                MATOP_ZERO_ENTRIES=23,
1421d519adbfSMatthew Knepley                MATOP_ZERO_ROWS=24,
1422d519adbfSMatthew Knepley                MATOP_LUFACTOR_SYMBOLIC=25,
1423d519adbfSMatthew Knepley                MATOP_LUFACTOR_NUMERIC=26,
1424d519adbfSMatthew Knepley                MATOP_CHOLESKY_FACTOR_SYMBOLIC=27,
1425d519adbfSMatthew Knepley                MATOP_CHOLESKY_FACTOR_NUMERIC=28,
1426d519adbfSMatthew Knepley                MATOP_SETUP_PREALLOCATION=29,
1427d519adbfSMatthew Knepley                MATOP_ILUFACTOR_SYMBOLIC=30,
1428d519adbfSMatthew Knepley                MATOP_ICCFACTOR_SYMBOLIC=31,
1429d519adbfSMatthew Knepley                MATOP_GET_ARRAY=32,
1430d519adbfSMatthew Knepley                MATOP_RESTORE_ARRAY=33,
1431d519adbfSMatthew Knepley                MATOP_DUPLICATE=34,
1432d519adbfSMatthew Knepley                MATOP_FORWARD_SOLVE=35,
1433d519adbfSMatthew Knepley                MATOP_BACKWARD_SOLVE=36,
1434d519adbfSMatthew Knepley                MATOP_ILUFACTOR=37,
1435d519adbfSMatthew Knepley                MATOP_ICCFACTOR=38,
1436d519adbfSMatthew Knepley                MATOP_AXPY=39,
1437d519adbfSMatthew Knepley                MATOP_GET_SUBMATRICES=40,
1438d519adbfSMatthew Knepley                MATOP_INCREASE_OVERLAP=41,
1439d519adbfSMatthew Knepley                MATOP_GET_VALUES=42,
1440d519adbfSMatthew Knepley                MATOP_COPY=43,
1441d519adbfSMatthew Knepley                MATOP_GET_ROW_MAX=44,
1442d519adbfSMatthew Knepley                MATOP_SCALE=45,
1443d519adbfSMatthew Knepley                MATOP_SHIFT=46,
1444d519adbfSMatthew Knepley                MATOP_DIAGONAL_SHIFT=47,
1445d519adbfSMatthew Knepley                MATOP_ILUDT_FACTOR=48,
1446d519adbfSMatthew Knepley                MATOP_SET_BLOCK_SIZE=49,
1447d519adbfSMatthew Knepley                MATOP_GET_ROW_IJ=50,
1448d519adbfSMatthew Knepley                MATOP_RESTORE_ROW_IJ=51,
1449d519adbfSMatthew Knepley                MATOP_GET_COLUMN_IJ=52,
1450d519adbfSMatthew Knepley                MATOP_RESTORE_COLUMN_IJ=53,
1451d519adbfSMatthew Knepley                MATOP_FDCOLORING_CREATE=54,
1452d519adbfSMatthew Knepley                MATOP_COLORING_PATCH=55,
1453d519adbfSMatthew Knepley                MATOP_SET_UNFACTORED=56,
1454d519adbfSMatthew Knepley                MATOP_PERMUTE=57,
1455d519adbfSMatthew Knepley                MATOP_SET_VALUES_BLOCKED=58,
1456d519adbfSMatthew Knepley                MATOP_GET_SUBMATRIX=59,
1457d519adbfSMatthew Knepley                MATOP_DESTROY=60,
1458d519adbfSMatthew Knepley                MATOP_VIEW=61,
1459d519adbfSMatthew Knepley                MATOP_CONVERT_FROM=62,
1460d519adbfSMatthew Knepley                MATOP_USE_SCALED_FORM=63,
1461d519adbfSMatthew Knepley                MATOP_SCALE_SYSTEM=64,
1462d519adbfSMatthew Knepley                MATOP_UNSCALE_SYSTEM=65,
1463d519adbfSMatthew Knepley                MATOP_SET_LOCAL_TO_GLOBAL_MAP=66,
1464d519adbfSMatthew Knepley                MATOP_SET_VALUES_LOCAL=67,
1465d519adbfSMatthew Knepley                MATOP_ZERO_ROWS_LOCAL=68,
1466d519adbfSMatthew Knepley                MATOP_GET_ROW_MAX_ABS=69,
1467d519adbfSMatthew Knepley                MATOP_GET_ROW_MIN_ABS=70,
1468d519adbfSMatthew Knepley                MATOP_CONVERT=71,
1469d519adbfSMatthew Knepley                MATOP_SET_COLORING=72,
1470d519adbfSMatthew Knepley                MATOP_SET_VALUES_ADIC=73,
1471d519adbfSMatthew Knepley                MATOP_SET_VALUES_ADIFOR=74,
1472d519adbfSMatthew Knepley                MATOP_FD_COLORING_APPLY=75,
1473d519adbfSMatthew Knepley                MATOP_SET_FROM_OPTIONS=76,
1474d519adbfSMatthew Knepley                MATOP_MULT_CON=77,
1475d519adbfSMatthew Knepley                MATOP_MULT_TRANSPOSE_CON=78,
1476d519adbfSMatthew Knepley                MATOP_PERMUTE_SPARSIFY=79,
1477d519adbfSMatthew Knepley                MATOP_MULT_MULTIPLE=80,
1478d519adbfSMatthew Knepley                MATOP_SOLVE_MULTIPLE=81,
1479d519adbfSMatthew Knepley                MATOP_GET_INERTIA=82,
1480d519adbfSMatthew Knepley                MATOP_LOAD=83,
1481d519adbfSMatthew Knepley                MATOP_IS_SYMMETRIC=84,
1482d519adbfSMatthew Knepley                MATOP_IS_HERMITIAN=85,
1483d519adbfSMatthew Knepley                MATOP_IS_STRUCTURALLY_SYMMETRIC=86,
1484d519adbfSMatthew Knepley                MATOP_PB_RELAX=87,
1485d519adbfSMatthew Knepley                MATOP_GET_VECS=88,
1486d519adbfSMatthew Knepley                MATOP_MAT_MULT=89,
1487d519adbfSMatthew Knepley                MATOP_MAT_MULT_SYMBOLIC=90,
1488d519adbfSMatthew Knepley                MATOP_MAT_MULT_NUMERIC=91,
1489d519adbfSMatthew Knepley                MATOP_PTAP=92,
1490d519adbfSMatthew Knepley                MATOP_PTAP_SYMBOLIC=93,
1491d519adbfSMatthew Knepley                MATOP_PTAP_NUMERIC=94,
1492d519adbfSMatthew Knepley                MATOP_MAT_MULTTRANSPOSE=95,
1493d519adbfSMatthew Knepley                MATOP_MAT_MULTTRANSPOSE_SYMBOLIC=96,
1494d519adbfSMatthew Knepley                MATOP_MAT_MULTTRANSPOSE_NUMERIC=97,
1495d519adbfSMatthew Knepley                MATOP_PTAP_SYMBOLIC_SEQAIJ=98,
1496d519adbfSMatthew Knepley                MATOP_PTAP_NUMERIC_SEQAIJ=99,
1497d519adbfSMatthew Knepley                MATOP_PTAP_SYMBOLIC_MPIAIJ=100,
1498d519adbfSMatthew Knepley                MATOP_PTAP_NUMERIC_MPIAIJ=101,
1499d519adbfSMatthew Knepley                MATOP_CONJUGATE=102,
1500d519adbfSMatthew Knepley                MATOP_SET_SIZES=103,
1501d519adbfSMatthew Knepley                MATOP_SET_VALUES_ROW = 104,
1502d519adbfSMatthew Knepley                MATOP_REAL_PART=105,
1503d519adbfSMatthew Knepley                MATOP_IMAG_PART=106,
1504d519adbfSMatthew Knepley                MATOP_GET_ROW_UTRIANGULAR=107,
1505d519adbfSMatthew Knepley                MATOP_RESTORE_ROW_UTRIANGULAR=108,
1506d519adbfSMatthew Knepley                MATOP_MATSOLVE=109,
1507d519adbfSMatthew Knepley                MATOP_GET_REDUNDANTMATRIX=110,
1508d519adbfSMatthew Knepley                MATOP_GET_ROW_MIN=111,
1509d519adbfSMatthew Knepley                MATOP_GET_COLUMN_VEC=112,
1510d519adbfSMatthew Knepley                MATOP_MISSING_DIAGONAL=113,
1511d519adbfSMatthew Knepley                MATOP_MATGETSEQNONZEROSTRUCTURE=114,
151289c6957cSBarry Smith                MATOP_CREATE=115,
151389c6957cSBarry Smith                MATOP_GET_GHOSTS=116,
151489c6957cSBarry Smith                MATOP_ILUDTFACTORSYMBOLIC=117,
151589c6957cSBarry Smith                MATOP_ILUDTFACTORNUMERIC=118,
151689c6957cSBarry Smith                MATOP_MULT_DIAGONAL_BLOCK=119
1517fae171e0SBarry Smith              } MatOperation;
1518be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatHasOperation(Mat,MatOperation,PetscTruth*);
1519be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellSetOperation(Mat,MatOperation,void(*)(void));
1520be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellGetOperation(Mat,MatOperation,void(**)(void));
1521be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellSetContext(Mat,void*);
1522112a2221SBarry Smith 
152390ace30eSBarry Smith /*
152490ace30eSBarry Smith    Codes for matrices stored on disk. By default they are
152590ace30eSBarry Smith  stored in a universal format. By changing the format with
15267973ad04SBarry Smith  PetscViewerSetFormat(viewer,PETSC_VIEWER_NATIVE); the matrices will
152790ace30eSBarry Smith  be stored in a way natural for the matrix, for example dense matrices
152890ace30eSBarry Smith  would be stored as dense. Matrices stored this way may only be
152990ace30eSBarry Smith  read into matrices of the same time.
153090ace30eSBarry Smith */
153190ace30eSBarry Smith #define MATRIX_BINARY_FORMAT_DENSE -1
153290ace30eSBarry Smith 
1533be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetHashTableFactor(Mat,PetscReal);
1534be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatISGetLocalMat(Mat,Mat*);
15351f4e1ec7SBarry Smith 
1536d9274352SBarry Smith /*S
1537d9274352SBarry Smith      MatNullSpace - Object that removes a null space from a vector, i.e.
1538d9274352SBarry Smith          orthogonalizes the vector to a subsapce
1539d9274352SBarry Smith 
1540f7a9e4ceSBarry Smith    Level: advanced
1541d9274352SBarry Smith 
1542d9274352SBarry Smith   Concepts: matrix; linear operator, null space
1543d9274352SBarry Smith 
15446e1639daSBarry Smith   Users manual sections:
15456e1639daSBarry Smith .   sec_singular
15466e1639daSBarry Smith 
1547d9274352SBarry Smith .seealso:  MatNullSpaceCreate()
1548d9274352SBarry Smith S*/
154974637425SBarry Smith typedef struct _p_MatNullSpace* MatNullSpace;
1550d9274352SBarry Smith 
1551be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceCreate(MPI_Comm,PetscTruth,PetscInt,const Vec[],MatNullSpace*);
1552281d1b2eSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(Vec,void*),void*);
1553be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceDestroy(MatNullSpace);
1554be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceRemove(MatNullSpace,Vec,Vec*);
1555be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceAttach(Mat,MatNullSpace);
155695902228SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceTest(MatNullSpace,Mat,PetscTruth *);
155774637425SBarry Smith 
1558be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatReorderingSeqSBAIJ(Mat,IS);
1559be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPISBAIJSetHashTableFactor(Mat,PetscReal);
1560be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetColumnIndices(Mat,PetscInt *);
156146129b97SKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJInvertBlockDiagonal(Mat);
15623f1d51d7SBarry Smith 
1563be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMAIJ(Mat,PetscInt,Mat*);
1564be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMAIJRedimension(Mat,PetscInt,Mat*);
1565be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMAIJGetAIJ(Mat,Mat*);
1566c4f061fbSSatish Balay 
1567be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatComputeExplicitOperator(Mat,Mat*);
1568b0a32e0cSBarry Smith 
1569be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalScaleLocal(Mat,Vec);
157004f1ad80SBarry Smith 
1571fef1beadSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*);
15723ec795f1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetBase(Mat,Vec,Vec);
15733ec795f1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*);
1574e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*));
1575e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec));
1576e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDAddNullSpace(Mat,MatNullSpace);
1577e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt);
1578e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDResetHHistory(Mat);
1579e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctionError(Mat,PetscReal);
1580e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetPeriod(Mat,PetscInt);
1581e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDGetH(Mat,PetscScalar *);
15826aa9148fSLisandro Dalcin EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetOptionsPrefix(Mat,const char[]);
1583e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFromOptions(Mat);
1584e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*);
1585e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*);
1586e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDDSSetUmin(Mat,PetscReal);
1587e884886fSBarry Smith 
1588e884886fSBarry Smith 
1589e884886fSBarry Smith typedef struct _p_MatMFFD* MatMFFD;
1590e884886fSBarry Smith 
1591e884886fSBarry Smith /*E
1592e884886fSBarry Smith     MatMFFDType - algorithm used to compute the h used in computing matrix-vector products via differencing of the function
1593e884886fSBarry Smith 
1594e884886fSBarry Smith    Level: beginner
1595e884886fSBarry Smith 
1596e884886fSBarry Smith .seealso: MatMFFDSetType(), MatMFFDRegister()
1597e884886fSBarry Smith E*/
1598a313700dSBarry Smith #define MatMFFDType char*
1599e884886fSBarry Smith #define MATMFFD_DS  "ds"
1600e884886fSBarry Smith #define MATMFFD_WP  "wp"
1601e884886fSBarry Smith 
1602a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetType(Mat,const MatMFFDType);
1603e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatMFFD));
1604e884886fSBarry Smith 
1605e884886fSBarry Smith /*MC
1606e884886fSBarry Smith    MatMFFDRegisterDynamic - Adds a method to the MatMFFD registry.
1607e884886fSBarry Smith 
1608e884886fSBarry Smith    Synopsis:
1609e884886fSBarry Smith    PetscErrorCode MatMFFDRegisterDynamic(char *name_solver,char *path,char *name_create,PetscErrorCode (*routine_create)(MatMFFD))
1610e884886fSBarry Smith 
1611e884886fSBarry Smith    Not Collective
1612e884886fSBarry Smith 
1613e884886fSBarry Smith    Input Parameters:
1614e884886fSBarry Smith +  name_solver - name of a new user-defined compute-h module
1615e884886fSBarry Smith .  path - path (either absolute or relative) the library containing this solver
1616e884886fSBarry Smith .  name_create - name of routine to create method context
1617e884886fSBarry Smith -  routine_create - routine to create method context
1618e884886fSBarry Smith 
1619e884886fSBarry Smith    Level: developer
1620e884886fSBarry Smith 
1621e884886fSBarry Smith    Notes:
1622e884886fSBarry Smith    MatMFFDRegisterDynamic() may be called multiple times to add several user-defined solvers.
1623e884886fSBarry Smith 
1624e884886fSBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
1625e884886fSBarry Smith    is ignored.
1626e884886fSBarry Smith 
1627e884886fSBarry Smith    Sample usage:
1628e884886fSBarry Smith .vb
1629e884886fSBarry Smith    MatMFFDRegisterDynamic("my_h",/home/username/my_lib/lib/libO/solaris/mylib.a,
1630e884886fSBarry Smith                "MyHCreate",MyHCreate);
1631e884886fSBarry Smith .ve
1632e884886fSBarry Smith 
1633e884886fSBarry Smith    Then, your solver can be chosen with the procedural interface via
1634e884886fSBarry Smith $     MatMFFDSetType(mfctx,"my_h")
1635e884886fSBarry Smith    or at runtime via the option
1636e884886fSBarry Smith $     -snes_mf_type my_h
1637e884886fSBarry Smith 
1638e884886fSBarry Smith .keywords: MatMFFD, register
1639e884886fSBarry Smith 
1640e884886fSBarry Smith .seealso: MatMFFDRegisterAll(), MatMFFDRegisterDestroy()
1641e884886fSBarry Smith M*/
1642e884886fSBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1643e884886fSBarry Smith #define MatMFFDRegisterDynamic(a,b,c,d) MatMFFDRegister(a,b,c,0)
1644e884886fSBarry Smith #else
1645e884886fSBarry Smith #define MatMFFDRegisterDynamic(a,b,c,d) MatMFFDRegister(a,b,c,d)
1646e884886fSBarry Smith #endif
1647e884886fSBarry Smith 
1648e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDRegisterAll(const char[]);
1649e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDRegisterDestroy(void);
1650e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDDefaultSetUmin(Mat,PetscReal);
1651e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDWPSetComputeNormU(Mat,PetscTruth);
1652e884886fSBarry Smith 
1653e884886fSBarry Smith 
1654be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *);
1655be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *);
16567dbadf16SMatthew Knepley 
1657e9fa29b7SSatish Balay PETSC_EXTERN_CXX_END
16582eac72dbSBarry Smith #endif
1659