xref: /petsc/include/petscmat.h (revision eeffb40d691afbdd57a8091619e7ddd44ac5fdca)
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 
70773941d6SBarry Smith 
719989ab13SBarry Smith /*E
72c7393fdbSBarry Smith     MatSolverPackage - String with the name of a PETSc matrix solver type.
739989ab13SBarry Smith 
749989ab13SBarry Smith     For example: "petsc" indicates what PETSc provides, "superlu" indicates either
759989ab13SBarry Smith        SuperLU or SuperLU_Dist etc.
769989ab13SBarry Smith 
779989ab13SBarry Smith 
789989ab13SBarry Smith    Level: beginner
799989ab13SBarry Smith 
805c9eb25fSBarry Smith .seealso: MatGetFactor(), Mat, MatSetType(), MatType
819989ab13SBarry Smith E*/
82c7393fdbSBarry Smith #define MatSolverPackage char*
83b24902e0SBarry Smith #define MAT_SOLVER_SPOOLES      "spooles"
84b24902e0SBarry Smith #define MAT_SOLVER_SUPERLU      "superlu"
85b24902e0SBarry Smith #define MAT_SOLVER_SUPERLU_DIST "superlu_dist"
86b24902e0SBarry Smith #define MAT_SOLVER_UMFPACK      "umfpack"
87b24902e0SBarry Smith #define MAT_SOLVER_ESSL         "essl"
88b24902e0SBarry Smith #define MAT_SOLVER_LUSOL        "lusol"
89b24902e0SBarry Smith #define MAT_SOLVER_MUMPS        "mumps"
903bf14a46SMatthew Knepley #define MAT_SOLVER_PASTIX       "pastix"
91b24902e0SBarry Smith #define MAT_SOLVER_DSCPACK      "dscpack"
92b24902e0SBarry Smith #define MAT_SOLVER_MATLAB       "matlab"
9343244d56SBarry Smith #define MAT_SOLVER_PETSC        "petsc"
94b6806ab0SHong Zhang #define MAT_SOLVER_PLAPACK      "plapack"
95b24902e0SBarry Smith 
96773941d6SBarry Smith 
97b24902e0SBarry Smith /*E
98b24902e0SBarry Smith     MatFactorType - indicates what type of factorization is requested
99b24902e0SBarry Smith 
100b24902e0SBarry Smith     Level: beginner
101b24902e0SBarry Smith 
102b24902e0SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
103b24902e0SBarry Smith 
104c7393fdbSBarry Smith .seealso: MatSolverPackage, MatGetFactor()
105b24902e0SBarry Smith E*/
106599ef60dSHong Zhang typedef enum {MAT_FACTOR_NONE, MAT_FACTOR_LU, MAT_FACTOR_CHOLESKY, MAT_FACTOR_ILU, MAT_FACTOR_ICC,MAT_FACTOR_ILUDT} MatFactorType;
107c7393fdbSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetFactor(Mat,const MatSolverPackage,MatFactorType,Mat*);
108db4efbfdSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetFactorAvailable(Mat,const MatSolverPackage,MatFactorType,PetscTruth*);
10935bd34faSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFactorGetSolverPackage(Mat,const MatSolverPackage*);
110b24902e0SBarry Smith 
1119989ab13SBarry Smith 
112c06d978dSMatthew Knepley /* Logging support */
113552e946dSBarry Smith #define    MAT_FILE_COOKIE 1211216    /* used to indicate matrices in binary files */
114be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_COOKIE;
115be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_FDCOLORING_COOKIE;
116be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_PARTITIONING_COOKIE;
117be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_NULLSPACE_COOKIE;
118e884886fSBarry Smith extern PetscCookie PETSCMAT_DLLEXPORT MATMFFD_COOKIE;
119c06d978dSMatthew Knepley 
120ceb03754SKris Buschelman /*E
121ceb03754SKris Buschelman     MatReuse - Indicates if matrices obtained from a previous call to MatGetSubMatrices()
122ceb03754SKris Buschelman      or MatGetSubMatrix() are to be reused to store the new matrix values.
123ceb03754SKris Buschelman 
124ceb03754SKris Buschelman     Level: beginner
125ceb03754SKris Buschelman 
126ceb03754SKris Buschelman    Any additions/changes here MUST also be made in include/finclude/petscmat.h
127ceb03754SKris Buschelman 
1280c451bc4SBarry Smith .seealso: MatGetSubMatrices(), MatGetSubMatrix(), MatDestroyMatrices(), MatConvert()
129ceb03754SKris Buschelman E*/
130dfe085dbSJed Brown typedef enum {MAT_INITIAL_MATRIX,MAT_REUSE_MATRIX,MAT_IGNORE_MATRIX} MatReuse;
131ceb03754SKris Buschelman 
1325494a064SHong Zhang /*E
1335494a064SHong Zhang     MatGetSubMatrixOption - Indicates if matrices obtained from a call to MatGetSubMatrices()
134829201f2SHong Zhang      include the matrix values. Currently it is only used by MatGetSeqNonzerostructure().
1355494a064SHong Zhang 
1365494a064SHong Zhang     Level: beginner
1375494a064SHong Zhang 
138829201f2SHong Zhang .seealso: MatGetSeqNonzerostructure()
1395494a064SHong Zhang E*/
1405494a064SHong Zhang typedef enum {MAT_DO_NOT_GET_VALUES,MAT_GET_VALUES} MatGetSubMatrixOption;
1415494a064SHong Zhang 
142e5bd5246SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInitializePackage(const char[]);
143c06d978dSMatthew Knepley 
144f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreate(MPI_Comm,Mat*);
145a69a960dSMatthew Knepley PetscPolymorphicFunction(MatCreate,(MPI_Comm comm),(comm,&A),Mat,A)
146a69a960dSMatthew Knepley PetscPolymorphicFunction(MatCreate,(),(PETSC_COMM_WORLD,&A),Mat,A)
147f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetSizes(Mat,PetscInt,PetscInt,PetscInt,PetscInt);
148a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetType(Mat,const MatType);
149be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetFromOptions(Mat);
150be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUpPreallocation(Mat);
151be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRegisterAll(const char[]);
152be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat));
15330de9b25SBarry Smith 
154f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetOptionsPrefix(Mat,const char[]);
155f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAppendOptionsPrefix(Mat,const char[]);
156f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOptionsPrefix(Mat,const char*[]);
157f69a0ea3SMatthew Knepley 
15830de9b25SBarry Smith /*MC
15930de9b25SBarry Smith    MatRegisterDynamic - Adds a new matrix type
16030de9b25SBarry Smith 
16130de9b25SBarry Smith    Synopsis:
162c1ac3661SBarry Smith    PetscErrorCode MatRegisterDynamic(char *name,char *path,char *name_create,PetscErrorCode (*routine_create)(Mat))
16330de9b25SBarry Smith 
16430de9b25SBarry Smith    Not Collective
16530de9b25SBarry Smith 
16630de9b25SBarry Smith    Input Parameters:
16730de9b25SBarry Smith +  name - name of a new user-defined matrix type
16830de9b25SBarry Smith .  path - path (either absolute or relative) the library containing this solver
16930de9b25SBarry Smith .  name_create - name of routine to create method context
17030de9b25SBarry Smith -  routine_create - routine to create method context
17130de9b25SBarry Smith 
17230de9b25SBarry Smith    Notes:
17330de9b25SBarry Smith    MatRegisterDynamic() may be called multiple times to add several user-defined solvers.
17430de9b25SBarry Smith 
17530de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
17630de9b25SBarry Smith    is ignored.
17730de9b25SBarry Smith 
17830de9b25SBarry Smith    Sample usage:
17930de9b25SBarry Smith .vb
18030de9b25SBarry Smith    MatRegisterDynamic("my_mat",/home/username/my_lib/lib/libO/solaris/mylib.a,
18130de9b25SBarry Smith                "MyMatCreate",MyMatCreate);
18230de9b25SBarry Smith .ve
18330de9b25SBarry Smith 
18430de9b25SBarry Smith    Then, your solver can be chosen with the procedural interface via
18530de9b25SBarry Smith $     MatSetType(Mat,"my_mat")
18630de9b25SBarry Smith    or at runtime via the option
18730de9b25SBarry Smith $     -mat_type my_mat
18830de9b25SBarry Smith 
18930de9b25SBarry Smith    Level: advanced
19030de9b25SBarry Smith 
191ab901514SBarry Smith    Notes: ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values.
19230de9b25SBarry Smith          If your function is not being put into a shared library then use VecRegister() instead
19330de9b25SBarry Smith 
19430de9b25SBarry Smith .keywords: Mat, register
19530de9b25SBarry Smith 
19630de9b25SBarry Smith .seealso: MatRegisterAll(), MatRegisterDestroy()
19730de9b25SBarry Smith 
19830de9b25SBarry Smith M*/
199273d9f13SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
200273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,0)
201273d9f13SBarry Smith #else
202273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,d)
20330de9b25SBarry Smith #endif
20430de9b25SBarry Smith 
205273d9f13SBarry Smith extern PetscTruth MatRegisterAllCalled;
206b0a32e0cSBarry Smith extern PetscFList MatList;
207b022a5c1SBarry Smith extern PetscFList MatColoringList;
208b022a5c1SBarry Smith extern PetscFList MatPartitioningList;
20928988994SBarry Smith 
2103b224e63SBarry Smith /*E
2113b224e63SBarry Smith     MatStructure - Indicates if the matrix has the same nonzero structure
2123b224e63SBarry Smith 
2133b224e63SBarry Smith     Level: beginner
2143b224e63SBarry Smith 
2153b224e63SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
2163b224e63SBarry Smith 
2173b224e63SBarry Smith .seealso: MatCopy(), KSPSetOperators(), PCSetOperators()
2183b224e63SBarry Smith E*/
2193b224e63SBarry Smith typedef enum {SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER,SUBSET_NONZERO_PATTERN} MatStructure;
2203b224e63SBarry Smith 
221be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqDense(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*);
222be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIDense(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*);
223be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
224ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,nz,nnz,&A),Mat,A)
225ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,&A),Mat,A)
226ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,0,nnz,&A),Mat,A)
227ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,&A),Mat,A)
228ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,A))
229ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,m,n,0,nnz,A))
230ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,A))
231be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
232ba966639SSatish 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)
233ba966639SSatish 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)
234ba966639SSatish 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)
235ba966639SSatish 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)
236ba966639SSatish 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))
237ba966639SSatish 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))
238ba966639SSatish 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))
239ba966639SSatish 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)
240ba966639SSatish 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)
241ba966639SSatish 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)
242ba966639SSatish 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)
243ba966639SSatish 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))
244ba966639SSatish 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))
245ba966639SSatish 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))
24663ab93bfSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *);
2478d7a6e47SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJWithSplitArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],PetscInt[],PetscInt[],PetscScalar[],Mat*);
2488d7a6e47SBarry Smith 
249be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIRowbs(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
250be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
251ba966639SSatish 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)
252ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A)
253ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A)
254ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A)
255ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A))
256ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A))
257ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A))
258be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
259ba966639SSatish 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)
260ba966639SSatish 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)
261ba966639SSatish 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)
262ba966639SSatish 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)
263ba966639SSatish 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))
264ba966639SSatish 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))
265ba966639SSatish 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))
266ba966639SSatish 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)
267ba966639SSatish 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)
268ba966639SSatish 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)
269ba966639SSatish 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)
270ba966639SSatish 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))
271ba966639SSatish 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))
272ba966639SSatish 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))
273be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAdj(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscInt[],Mat*);
274be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
275ba966639SSatish 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)
276ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A)
277ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A)
278ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A)
279ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A))
280ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A))
281ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A))
282be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPISBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
283ba966639SSatish 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)
284ba966639SSatish 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)
285ba966639SSatish 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)
286ba966639SSatish 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)
287ba966639SSatish 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))
288ba966639SSatish 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))
289ba966639SSatish 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))
290ba966639SSatish 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)
291ba966639SSatish 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)
292ba966639SSatish 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)
293ba966639SSatish 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)
294ba966639SSatish 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))
295ba966639SSatish 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))
296ba966639SSatish 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))
297be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateShell(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,void *,Mat*);
298ba966639SSatish 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)
299ba966639SSatish Balay PetscPolymorphicFunction(MatCreateShell,(PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx),(PETSC_COMM_WORLD,m,n,M,N,ctx,&A),Mat,A)
300be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateAdic(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,void (*)(void),Mat*);
301be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateNormal(Mat,Mat*);
302ba966639SSatish Balay PetscPolymorphicFunction(MatCreateNormal,(Mat mat),(mat,&A),Mat,A)
303ab92ecdeSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateLRC(Mat,Mat,Mat,Mat*);
304284134d9SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateIS(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,ISLocalToGlobalMapping,Mat*);
3051472f72bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
3061472f72bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPICRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
3072a6744ebSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateScatter(MPI_Comm,VecScatter,Mat*);
3082a6744ebSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScatterSetVecScatter(Mat,VecScatter);
3092a6744ebSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScatterGetVecScatter(Mat,VecScatter*);
3108cdf2d9bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateBlockMat(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,Mat*);
311793850ffSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompositeAddMat(Mat,Mat);
312b52f573bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompositeMerge(Mat);
313793850ffSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateComposite(MPI_Comm,PetscInt,const Mat*,Mat*);
3146d7c1e57SBarry Smith typedef enum {MAT_COMPOSITE_ADDITIVE,MAT_COMPOSITE_MULTIPLICATIVE} MatCompositeType;
3156d7c1e57SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompositeSetType(Mat,MatCompositeType);
3166d7c1e57SBarry Smith 
3175a7f1df3SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqFFTW(MPI_Comm,PetscInt,const PetscInt[],Mat*);
318f20085c4SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateTranspose(Mat,Mat*);
319ee1cef2cSJed Brown EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSubMatrix(Mat,IS,IS,Mat*);
320ee1cef2cSJed Brown EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSubMatrixUpdate(Mat,Mat,IS,IS);
3211472f72bSBarry Smith 
3221d6018f0SLisandro Dalcin EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreatePython(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const char[],Mat*);
3231d6018f0SLisandro Dalcin EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPythonSetType(Mat,const char[]);
3241d6018f0SLisandro Dalcin 
3251d6018f0SLisandro Dalcin 
326f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUp(Mat);
327be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroy(Mat);
32821c89e3eSBarry Smith 
3290c451bc4SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConjugate(Mat);
33099cafbc1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRealPart(Mat);
33199cafbc1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatImaginaryPart(Mat);
332bfc7d00aSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetDiagonalBlock(Mat,PetscTruth*,MatReuse,Mat*);
33399cafbc1SBarry Smith 
3348ed539a5SBarry Smith /* ------------------------------------------------------------*/
335be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
336be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
33787d4246cSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesRow(Mat,PetscInt,const PetscScalar[]);
338f259bd47SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesRowLocal(Mat,PetscInt,const PetscScalar[]);
33984cb2905SBarry Smith 
3402ef4de8bSBarry Smith /*S
3412ef4de8bSBarry Smith      MatStencil - Data structure (C struct) for storing information about a single row or
3422ef4de8bSBarry Smith         column of a matrix as index on an associated grid.
3432ef4de8bSBarry Smith 
3442ef4de8bSBarry Smith    Level: beginner
3452ef4de8bSBarry Smith 
3462ef4de8bSBarry Smith   Concepts: matrix; linear operator
3472ef4de8bSBarry Smith 
348d7bd93baSBarry Smith .seealso:  MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockStencil()
3492ef4de8bSBarry Smith S*/
350435da068SBarry Smith typedef struct {
351c1ac3661SBarry Smith   PetscInt k,j,i,c;
352435da068SBarry Smith } MatStencil;
3532ef4de8bSBarry Smith 
354be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
355be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
356be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt);
357435da068SBarry Smith 
358be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetColoring(Mat,ISColoring);
359be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesAdic(Mat,void*);
360be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesAdifor(Mat,PetscInt,void*);
3613a7fca6bSBarry Smith 
362d91e6319SBarry Smith /*E
363d91e6319SBarry Smith     MatAssemblyType - Indicates if the matrix is now to be used, or if you plan
364d91e6319SBarry Smith      to continue to add values to it
365d91e6319SBarry Smith 
366d91e6319SBarry Smith     Level: beginner
367d91e6319SBarry Smith 
368d91e6319SBarry Smith .seealso: MatAssemblyBegin(), MatAssemblyEnd()
369d91e6319SBarry Smith E*/
3706d4a8577SBarry Smith typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType;
371be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssemblyBegin(Mat,MatAssemblyType);
372be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssemblyEnd(Mat,MatAssemblyType);
373be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssembled(Mat,PetscTruth*);
3744f9c727eSBarry Smith 
3751d73ed98SBarry Smith 
37630de9b25SBarry Smith 
377d91e6319SBarry Smith /*E
378d91e6319SBarry Smith     MatOption - Options that may be set for a matrix and its behavior or storage
379d91e6319SBarry Smith 
380d91e6319SBarry Smith     Level: beginner
381d91e6319SBarry Smith 
3820a835dfdSSatish Balay    Any additions/changes here MUST also be made in include/finclude/petscmat.h
383d91e6319SBarry Smith 
384d91e6319SBarry Smith .seealso: MatSetOption()
385d91e6319SBarry Smith E*/
386512a5fc5SBarry Smith typedef enum {MAT_ROW_ORIENTED,MAT_NEW_NONZERO_LOCATIONS,
3874e0d8c25SBarry Smith               MAT_SYMMETRIC,
3884e0d8c25SBarry Smith               MAT_STRUCTURALLY_SYMMETRIC,
3898047e77bSBarry Smith               MAT_NEW_DIAGONALS,MAT_IGNORE_OFF_PROC_ENTRIES,
3904e0d8c25SBarry Smith               MAT_NEW_NONZERO_LOCATION_ERR,
3914e0d8c25SBarry Smith               MAT_NEW_NONZERO_ALLOCATION_ERR,MAT_USE_HASH_TABLE,
392a9817697SBarry Smith               MAT_KEEP_NONZERO_PATTERN,MAT_IGNORE_ZERO_ENTRIES,
3934e0d8c25SBarry Smith               MAT_USE_INODES,
3944e0d8c25SBarry Smith               MAT_HERMITIAN,
3954e0d8c25SBarry Smith               MAT_SYMMETRY_ETERNAL,
3964e0d8c25SBarry Smith               MAT_USE_COMPRESSEDROW,
3974e0d8c25SBarry Smith               MAT_IGNORE_LOWER_TRIANGULAR,MAT_ERROR_LOWER_TRIANGULAR,
39828b2fa4aSMatthew Knepley               MAT_GETROW_UPPERTRIANGULAR,MAT_UNUSED_NONZERO_LOCATION_ERR,
39928b2fa4aSMatthew Knepley               NUM_MAT_OPTIONS} MatOption;
400290bbb0aSBarry Smith extern const char *MatOptions[];
4014e0d8c25SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetOption(Mat,MatOption,PetscTruth);
402a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetType(Mat,const MatType*);
403a313700dSBarry Smith PetscPolymorphicFunction(MatGetType,(Mat mat),(mat,&t),const MatType,t)
40484cb2905SBarry Smith 
405be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]);
406be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
407be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
408f5edf698SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowUpperTriangular(Mat);
409f5edf698SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRowUpperTriangular(Mat);
410be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
411be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
412be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumnVector(Mat,Vec,PetscInt);
413be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetArray(Mat,PetscScalar *[]);
414ba966639SSatish Balay PetscPolymorphicFunction(MatGetArray,(Mat mat),(mat,&a),PetscScalar*,a)
415be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreArray(Mat,PetscScalar *[]);
416be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBlockSize(Mat,PetscInt *);
417ba966639SSatish Balay PetscPolymorphicFunction(MatGetBlockSize,(Mat mat),(mat,&a),PetscInt,a)
418be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetBlockSize(Mat,PetscInt);
4197b80b807SBarry Smith 
4201620fd73SBarry Smith 
421be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMult(Mat,Vec,Vec);
42289c6957cSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultDiagonalBlock(Mat,Vec,Vec);
423be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultAdd(Mat,Vec,Vec,Vec);
424ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
425be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTranspose(Mat,Vec,Vec);
426be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsTranspose(Mat,Mat,PetscReal,PetscTruth*);
427ba966639SSatish Balay PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B,PetscReal tol),(A,B,tol,&t),PetscTruth,t)
428e37f64cdSBarry Smith PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B),(A,B,0,&t),PetscTruth,t)
4291cbb95d3SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitianTranspose(Mat,Mat,PetscReal,PetscTruth*);
430be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeAdd(Mat,Vec,Vec,Vec);
431ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultTransposeAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
432be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultConstrained(Mat,Vec,Vec);
433be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeConstrained(Mat,Vec,Vec);
4342bebee5dSHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatSolve(Mat,Mat,Mat);
435f5edf698SHong Zhang 
436d91e6319SBarry Smith /*E
437d91e6319SBarry Smith     MatDuplicateOption - Indicates if a duplicated sparse matrix should have
438d91e6319SBarry Smith   its numerical values copied over or just its nonzero structure.
439d91e6319SBarry Smith 
440d91e6319SBarry Smith     Level: beginner
441d91e6319SBarry Smith 
442d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
443d91e6319SBarry Smith 
44470dcbbb9SBarry Smith $   MAT_SHARE_NONZERO_PATTERN - the i and j arrays in the new matrix will be shared with the original matrix
44570dcbbb9SBarry Smith $                               this also triggers the MAT_DO_NOT_COPY_VALUES option. This is used when you
44670dcbbb9SBarry Smith $                               have several matrices with the same nonzero pattern.
44770dcbbb9SBarry Smith 
448d91e6319SBarry Smith .seealso: MatDuplicate()
449d91e6319SBarry Smith E*/
45070dcbbb9SBarry Smith typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES,MAT_SHARE_NONZERO_PATTERN} MatDuplicateOption;
4512e8a6d31SBarry Smith 
452a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert(Mat,const MatType,MatReuse,Mat*);
453a313700dSBarry Smith PetscPolymorphicFunction(MatConvert,(Mat A,const MatType t),(A,t,MAT_INITIAL_MATRIX,&a),Mat,a)
454be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDuplicate(Mat,MatDuplicateOption,Mat*);
455ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A,MatDuplicateOption o),(A,o,&a),Mat,a)
456ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A),(A,MAT_COPY_VALUES,&a),Mat,a)
45794a9d846SBarry Smith 
458cb5b572fSBarry Smith 
459be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCopy(Mat,Mat,MatStructure);
460be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatView(Mat,PetscViewer);
461be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsSymmetric(Mat,PetscReal,PetscTruth*);
462ba966639SSatish Balay PetscPolymorphicFunction(MatIsSymmetric,(Mat A,PetscReal tol),(A,tol,&t),PetscTruth,t)
463e37f64cdSBarry Smith PetscPolymorphicFunction(MatIsSymmetric,(Mat A),(A,0,&t),PetscTruth,t)
464be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsStructurallySymmetric(Mat,PetscTruth*);
465ba966639SSatish Balay PetscPolymorphicFunction(MatIsStructurallySymmetric,(Mat A),(A,&t),PetscTruth,t)
4661cbb95d3SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitian(Mat,PetscReal,PetscTruth*);
4671cbb95d3SBarry Smith PetscPolymorphicFunction(MatIsHermitian,(Mat A),(A,0,&t),PetscTruth,t)
468be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsSymmetricKnown(Mat,PetscTruth*,PetscTruth*);
469be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitianKnown(Mat,PetscTruth*,PetscTruth*);
47064c62002SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMissingDiagonal(Mat,PetscTruth *,PetscInt *);
471a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLoad(PetscViewer,const MatType,Mat*);
472a313700dSBarry Smith PetscPolymorphicFunction(MatLoad,(PetscViewer v,const MatType t),(v,t,&a),Mat,a)
4737b80b807SBarry Smith 
4748f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *);
4758f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRowIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *);
4768f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumnIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *);
4778f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreColumnIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *);
478d4fbbf0eSBarry Smith 
479d91e6319SBarry Smith /*S
480d91e6319SBarry Smith      MatInfo - Context of matrix information, used with MatGetInfo()
481d91e6319SBarry Smith 
482d91e6319SBarry Smith    In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE
483d91e6319SBarry Smith 
484d91e6319SBarry Smith    Level: intermediate
485d91e6319SBarry Smith 
486d91e6319SBarry Smith   Concepts: matrix^nonzero information
487d91e6319SBarry Smith 
488d9274352SBarry Smith .seealso:  MatGetInfo(), MatInfoType
489d91e6319SBarry Smith S*/
4904e220ebcSLois Curfman McInnes typedef struct {
491b0a32e0cSBarry Smith   PetscLogDouble block_size;                         /* block size */
492b0a32e0cSBarry Smith   PetscLogDouble nz_allocated,nz_used,nz_unneeded;   /* number of nonzeros */
493b0a32e0cSBarry Smith   PetscLogDouble memory;                             /* memory allocated */
494b0a32e0cSBarry Smith   PetscLogDouble assemblies;                         /* number of matrix assemblies called */
495b0a32e0cSBarry Smith   PetscLogDouble mallocs;                            /* number of mallocs during MatSetValues() */
496b0a32e0cSBarry Smith   PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */
497b0a32e0cSBarry Smith   PetscLogDouble factor_mallocs;                     /* number of mallocs during factorization */
4984e220ebcSLois Curfman McInnes } MatInfo;
4994e220ebcSLois Curfman McInnes 
500d9274352SBarry Smith /*E
501d9274352SBarry Smith     MatInfoType - Indicates if you want information about the local part of the matrix,
502d9274352SBarry Smith      the entire parallel matrix or the maximum over all the local parts.
503d9274352SBarry Smith 
504d9274352SBarry Smith     Level: beginner
505d9274352SBarry Smith 
506d9274352SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
507d9274352SBarry Smith 
508d9274352SBarry Smith .seealso: MatGetInfo(), MatInfo
509d9274352SBarry Smith E*/
5107b80b807SBarry Smith typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType;
511be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetInfo(Mat,MatInfoType,MatInfo*);
512be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatValid(Mat,PetscTruth*);
513be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetDiagonal(Mat,Vec);
514985db425SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMax(Mat,Vec,PetscInt[]);
515985db425SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMin(Mat,Vec,PetscInt[]);
516985db425SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMaxAbs(Mat,Vec,PetscInt[]);
517c87e5d42SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMinAbs(Mat,Vec,PetscInt[]);
51886697f06SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowSum(Mat,Vec);
519fc4dec0aSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatTranspose(Mat,MatReuse,Mat*);
520fc4dec0aSBarry Smith PetscPolymorphicFunction(MatTranspose,(Mat A),(A,MAT_INITIAL_MATRIX,&t),Mat,t)
521*eeffb40dSHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatHermitianTranspose(Mat,MatReuse,Mat*);
522be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPermute(Mat,IS,IS,Mat *);
523ba966639SSatish Balay PetscPolymorphicFunction(MatPermute,(Mat A,IS is1,IS is2),(A,is1,is2,&t),Mat,t)
524be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPermuteSparsify(Mat,PetscInt,PetscReal,PetscReal,IS,IS,Mat *);
525be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalScale(Mat,Vec,Vec);
526be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalSet(Mat,Vec,InsertMode);
527be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatEqual(Mat,Mat,PetscTruth*);
528ba966639SSatish Balay PetscPolymorphicFunction(MatEqual,(Mat A,Mat B),(A,B,&t),PetscTruth,t)
529be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultEqual(Mat,Mat,PetscInt,PetscTruth*);
530be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultAddEqual(Mat,Mat,PetscInt,PetscTruth*);
531be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeEqual(Mat,Mat,PetscInt,PetscTruth*);
532be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscTruth*);
5337b80b807SBarry Smith 
534be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNorm(Mat,NormType,PetscReal *);
535ba966639SSatish Balay PetscPolymorphicFunction(MatNorm,(Mat A,NormType t),(A,t,&n),PetscReal,n)
536242f1d38SBarry Smith EXTERN PetscErrorCode MatGetColumnNorms(Mat,NormType,PetscReal *);
537be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroEntries(Mat);
538f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar);
539f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsIS(Mat,IS,PetscScalar);
540f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroColumns(Mat,PetscInt,const PetscInt [],const PetscScalar*);
541f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroColumnsIS(Mat,IS,const PetscScalar*);
5427b80b807SBarry Smith 
543be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatUseScaledForm(Mat,PetscTruth);
544be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScaleSystem(Mat,Vec,Vec);
545be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatUnScaleSystem(Mat,Vec,Vec);
5465ef9f2a5SBarry Smith 
547be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSize(Mat,PetscInt*,PetscInt*);
548be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalSize(Mat,PetscInt*,PetscInt*);
549be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRange(Mat,PetscInt*,PetscInt*);
5508ee2e534SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRanges(Mat,const PetscInt**);
551f4b8409dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRangeColumn(Mat,PetscInt*,PetscInt*);
55272ca8751SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRangesColumn(Mat,const PetscInt**);
5537b80b807SBarry Smith 
554be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
555be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroyMatrices(PetscInt,Mat *[]);
5564aa3045dSJed Brown EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrix(Mat,IS,IS,MatReuse,Mat *);
557f6d58c54SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSeqNonzeroStructure(Mat,Mat*);
558f6d58c54SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroySeqNonzeroStructure(Mat*);
5595494a064SHong Zhang 
560be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
561be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPI(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*);
562be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPISymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*);
563be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPINumeric(Mat,Mat);
564be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroy_MPIAIJ_SeqsToMPI(Mat);
565be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalMat(Mat,MatReuse,Mat*);
566be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*);
567be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,PetscInt*,Mat*);
568dd6ea824SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBrowsOfAoCols(Mat,Mat,MatReuse,PetscInt**,MatScalar**,Mat*);
56943eb5e2fSMatthew Knepley #if defined (PETSC_USE_CTABLE)
570cd116e26SSatish Balay #include "petscctable.h"
57143eb5e2fSMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetCommunicationStructs(Mat, Vec *, PetscTable *, VecScatter *);
57243eb5e2fSMatthew Knepley #else
57343eb5e2fSMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetCommunicationStructs(Mat, Vec *, PetscInt *[], VecScatter *);
57443eb5e2fSMatthew Knepley #endif
5758c7482ecSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetGhosts(Mat, PetscInt *,const PetscInt *[]);
5768efafbd8SBarry Smith 
577be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt);
5787b80b807SBarry Smith 
579be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*);
580be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultSymbolic(Mat,Mat,PetscReal,Mat*);
581be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultNumeric(Mat,Mat,Mat);
58222440eb1SKris Buschelman 
583be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*);
584be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAPSymbolic(Mat,Mat,PetscReal,Mat*);
585be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAPNumeric(Mat,Mat,Mat);
58622440eb1SKris Buschelman 
587be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTranspose(Mat,Mat,MatReuse,PetscReal,Mat*);
588be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTransposeSymbolic(Mat,Mat,PetscReal,Mat*);
589be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTransposeNumeric(Mat,Mat,Mat);
590bc011b1eSHong Zhang 
591f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAXPY(Mat,PetscScalar,Mat,MatStructure);
59204aac2b0SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAYPX(Mat,PetscScalar,Mat,MatStructure);
5931c741599SHong Zhang 
594f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScale(Mat,PetscScalar);
595f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShift(Mat,PetscScalar);
5967b80b807SBarry Smith 
597be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping);
598be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetLocalToGlobalMappingBlock(Mat,ISLocalToGlobalMapping);
599f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar);
600f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsLocalIS(Mat,IS,PetscScalar);
601be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
602be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
603052efed2SBarry Smith 
604be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStashSetInitialSize(Mat,PetscInt,PetscInt);
605be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
60690f02eecSBarry Smith 
607be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInterpolate(Mat,Vec,Vec);
6080c451bc4SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInterpolateAdd(Mat,Vec,Vec,Vec);
609ba966639SSatish Balay PetscPolymorphicSubroutine(MatInterpolateAdd,(Mat A,Vec x,Vec y),(A,x,y,y))
610be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestrict(Mat,Vec,Vec);
611be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetVecs(Mat,Vec*,Vec*);
61269db28dcSHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRedundantMatrix(Mat,PetscInt,MPI_Comm,PetscInt,MatReuse,Mat*);
613bd481603SBarry Smith 
614bd481603SBarry Smith /*MC
6151620fd73SBarry Smith    MatSetValue - Set a single entry into a matrix.
6161620fd73SBarry Smith 
6171620fd73SBarry Smith    Not collective
6181620fd73SBarry Smith 
6191620fd73SBarry Smith    Input Parameters:
6201620fd73SBarry Smith +  m - the matrix
6211620fd73SBarry Smith .  row - the row location of the entry
6221620fd73SBarry Smith .  col - the column location of the entry
6231620fd73SBarry Smith .  value - the value to insert
6241620fd73SBarry Smith -  mode - either INSERT_VALUES or ADD_VALUES
6251620fd73SBarry Smith 
6261620fd73SBarry Smith    Notes:
6271620fd73SBarry Smith    For efficiency one should use MatSetValues() and set several or many
6281620fd73SBarry Smith    values simultaneously if possible.
6291620fd73SBarry Smith 
6301620fd73SBarry Smith    Level: beginner
6311620fd73SBarry Smith 
6321620fd73SBarry Smith .seealso: MatSetValues(), MatSetValueLocal()
6331620fd73SBarry Smith M*/
6341620fd73SBarry 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);}
6351620fd73SBarry Smith 
6361620fd73SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatGetValue(Mat v,PetscInt i,PetscInt j,PetscScalar *va) {return MatGetValues(v,1,&i,1,&j,va);}
6371620fd73SBarry Smith 
6381620fd73SBarry 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);}
6391620fd73SBarry Smith 
6401620fd73SBarry Smith /*MC
6410d2a7ff2SSatish Balay    MatPreallocateInitialize - Begins the block of code that will count the number of nonzeros per
642bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
643bd481603SBarry Smith 
644bd481603SBarry Smith    Synopsis:
645c1ac3661SBarry Smith    PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)
646bd481603SBarry Smith 
647bd481603SBarry Smith    Collective on MPI_Comm
648bd481603SBarry Smith 
649bd481603SBarry Smith    Input Parameters:
650bd481603SBarry Smith +  comm - the communicator that will share the eventually allocated matrix
651859fcb39SBarry Smith .  nrows - the number of LOCAL rows in the matrix
652859fcb39SBarry Smith -  ncols - the number of LOCAL columns in the matrix
653bd481603SBarry Smith 
654bd481603SBarry Smith    Output Parameters:
655bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
656bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
657bd481603SBarry Smith 
658bd481603SBarry Smith 
659bd481603SBarry Smith    Level: intermediate
660bd481603SBarry Smith 
661bd481603SBarry Smith    Notes:
662bd481603SBarry Smith    See the chapter in the users manual on performance for more details
663bd481603SBarry Smith 
6641d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
665bd481603SBarry Smith 
666bd481603SBarry Smith    Use MatPreallocateInitializeSymmetric() for symmetric matrices (MPISBAIJ matrices)
667bd481603SBarry Smith 
6681620fd73SBarry Smith    This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize().
6691620fd73SBarry Smith 
670bd481603SBarry Smith   Concepts: preallocation^Matrix
671bd481603SBarry Smith 
672bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
673bd481603SBarry Smith           MatPreallocateInitializeSymmetric(), MatPreallocateSymmetricSetLocal()
674bd481603SBarry Smith M*/
675c4f061fbSSatish Balay #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \
6767c922b88SBarry Smith { \
677a89bc211SBarry Smith   PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__start,__end; \
678a89bc211SBarry Smith   _4_ierr = PetscMalloc2(__nrows,PetscInt,&dnz,__nrows,PetscInt,&onz);CHKERRQ(_4_ierr); \
679a89bc211SBarry Smith   _4_ierr = PetscMemzero(dnz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
680a89bc211SBarry Smith   _4_ierr = PetscMemzero(onz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
681c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ctmp;\
682a89bc211SBarry Smith   _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows;
6837c922b88SBarry Smith 
684bd481603SBarry Smith /*MC
6850d2a7ff2SSatish Balay    MatPreallocateSymmetricInitialize - Begins the block of code that will count the number of nonzeros per
686bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
687bd481603SBarry Smith 
688bd481603SBarry Smith    Synopsis:
689c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)
690bd481603SBarry Smith 
691bd481603SBarry Smith    Collective on MPI_Comm
692bd481603SBarry Smith 
693bd481603SBarry Smith    Input Parameters:
694bd481603SBarry Smith +  comm - the communicator that will share the eventually allocated matrix
695859fcb39SBarry Smith .  nrows - the number of LOCAL rows in the matrix
696859fcb39SBarry Smith -  ncols - the number of LOCAL columns in the matrix
697bd481603SBarry Smith 
698bd481603SBarry Smith    Output Parameters:
699bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
700bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
701bd481603SBarry Smith 
702bd481603SBarry Smith 
703bd481603SBarry Smith    Level: intermediate
704bd481603SBarry Smith 
705bd481603SBarry Smith    Notes:
706bd481603SBarry Smith    See the chapter in the users manual on performance for more details
707bd481603SBarry Smith 
7081d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
709bd481603SBarry Smith 
7101620fd73SBarry Smith    This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize().
7111620fd73SBarry Smith 
712bd481603SBarry Smith   Concepts: preallocation^Matrix
713bd481603SBarry Smith 
714bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
715bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal()
716bd481603SBarry Smith M*/
717222b16d4SBarry Smith #define MatPreallocateSymmetricInitialize(comm,nrows,ncols,dnz,onz) 0; \
718222b16d4SBarry Smith { \
719a89bc211SBarry Smith   PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__end; \
720a89bc211SBarry Smith   _4_ierr = PetscMalloc2(__nrows,PetscInt,&dnz,__nrows,PetscInt,&onz);CHKERRQ(_4_ierr); \
721a89bc211SBarry Smith   _4_ierr = PetscMemzero(dnz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
722a89bc211SBarry Smith   _4_ierr = PetscMemzero(onz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
723c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr);\
724a89bc211SBarry Smith   _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows;
725222b16d4SBarry Smith 
726bd481603SBarry Smith /*MC
727bd481603SBarry Smith    MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
728bd481603SBarry Smith        inserted using a local number of the rows and columns
729bd481603SBarry Smith 
730bd481603SBarry Smith    Synopsis:
731c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
732bd481603SBarry Smith 
733bd481603SBarry Smith    Not Collective
734bd481603SBarry Smith 
735bd481603SBarry Smith    Input Parameters:
736bd481603SBarry Smith +  map - the mapping between local numbering and global numbering
737bd481603SBarry Smith .  nrows - the number of rows indicated
7381d73ed98SBarry Smith .  rows - the indices of the rows
739bd481603SBarry Smith .  ncols - the number of columns in the matrix
740bd481603SBarry Smith .  cols - the columns indicated
741bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
742bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
743bd481603SBarry Smith 
744bd481603SBarry Smith 
745bd481603SBarry Smith    Level: intermediate
746bd481603SBarry Smith 
747bd481603SBarry Smith    Notes:
748bd481603SBarry Smith    See the chapter in the users manual on performance for more details
749bd481603SBarry Smith 
7501d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
751bd481603SBarry Smith 
752bd481603SBarry Smith   Concepts: preallocation^Matrix
753bd481603SBarry Smith 
754bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
755bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal()
756bd481603SBarry Smith M*/
757c4f061fbSSatish Balay #define MatPreallocateSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\
758c4f061fbSSatish Balay {\
759c1ac3661SBarry Smith   PetscInt __l;\
760ef66eb69SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\
761ef66eb69SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\
762c4f061fbSSatish Balay   for (__l=0;__l<nrows;__l++) {\
763ef66eb69SBarry Smith     _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
764c4f061fbSSatish Balay   }\
765c4f061fbSSatish Balay }
766c4f061fbSSatish Balay 
767bd481603SBarry Smith /*MC
768bd481603SBarry Smith    MatPreallocateSymmetricSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
769bd481603SBarry Smith        inserted using a local number of the rows and columns
770bd481603SBarry Smith 
771bd481603SBarry Smith    Synopsis:
772c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
773bd481603SBarry Smith 
774bd481603SBarry Smith    Not Collective
775bd481603SBarry Smith 
776bd481603SBarry Smith    Input Parameters:
777bd481603SBarry Smith +  map - the mapping between local numbering and global numbering
778bd481603SBarry Smith .  nrows - the number of rows indicated
7791d73ed98SBarry Smith .  rows - the indices of the rows
780bd481603SBarry Smith .  ncols - the number of columns in the matrix
781bd481603SBarry Smith .  cols - the columns indicated
782bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
783bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
784bd481603SBarry Smith 
785bd481603SBarry Smith 
786bd481603SBarry Smith    Level: intermediate
787bd481603SBarry Smith 
788bd481603SBarry Smith    Notes:
789bd481603SBarry Smith    See the chapter in the users manual on performance for more details
790bd481603SBarry Smith 
791bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
792bd481603SBarry Smith 
793bd481603SBarry Smith   Concepts: preallocation^Matrix
794bd481603SBarry Smith 
795bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
796bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
797bd481603SBarry Smith M*/
798d3d32019SBarry Smith #define MatPreallocateSymmetricSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\
799d3d32019SBarry Smith {\
800c1ac3661SBarry Smith   PetscInt __l;\
801d3d32019SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\
802d3d32019SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\
803d3d32019SBarry Smith   for (__l=0;__l<nrows;__l++) {\
804d3d32019SBarry Smith     _4_ierr = MatPreallocateSymmetricSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
805d3d32019SBarry Smith   }\
806d3d32019SBarry Smith }
807d3d32019SBarry Smith 
808bd481603SBarry Smith /*MC
809bd481603SBarry Smith    MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
810bd481603SBarry Smith        inserted using a local number of the rows and columns
811bd481603SBarry Smith 
812bd481603SBarry Smith    Synopsis:
813c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
814bd481603SBarry Smith 
815bd481603SBarry Smith    Not Collective
816bd481603SBarry Smith 
817bd481603SBarry Smith    Input Parameters:
81864075487SBarry Smith +  row - the row
819bd481603SBarry Smith .  ncols - the number of columns in the matrix
820a50f8bf6SHong Zhang -  cols - the columns indicated
821a50f8bf6SHong Zhang 
822a50f8bf6SHong Zhang    Output Parameters:
823a50f8bf6SHong Zhang +  dnz - the array that will be passed to the matrix preallocation routines
824bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
825bd481603SBarry Smith 
826bd481603SBarry Smith 
827bd481603SBarry Smith    Level: intermediate
828bd481603SBarry Smith 
829bd481603SBarry Smith    Notes:
830bd481603SBarry Smith    See the chapter in the users manual on performance for more details
831bd481603SBarry Smith 
832bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
833bd481603SBarry Smith 
8341620fd73SBarry Smith    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
8351620fd73SBarry Smith 
836bd481603SBarry Smith   Concepts: preallocation^Matrix
837bd481603SBarry Smith 
838bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
839bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
840bd481603SBarry Smith M*/
841c4f061fbSSatish Balay #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\
842c1ac3661SBarry Smith { PetscInt __i; \
8438ee2e534SBarry Smith   if (row < __rstart) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Trying to set preallocation for row %D less than first local row %D",row,__rstart);\
844a89bc211SBarry 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);\
8457c922b88SBarry Smith   for (__i=0; __i<nc; __i++) {\
84664075487SBarry Smith     if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \
8477cc688d7SBarry Smith     else dnz[row - __rstart]++;\
8487c922b88SBarry Smith   }\
8497c922b88SBarry Smith }
8507c922b88SBarry Smith 
851bd481603SBarry Smith /*MC
852bd481603SBarry Smith    MatPreallocateSymmetricSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
853bd481603SBarry Smith        inserted using a local number of the rows and columns
854bd481603SBarry Smith 
855bd481603SBarry Smith    Synopsis:
856c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
857bd481603SBarry Smith 
858bd481603SBarry Smith    Not Collective
859bd481603SBarry Smith 
860bd481603SBarry Smith    Input Parameters:
861bd481603SBarry Smith +  nrows - the number of rows indicated
8621d73ed98SBarry Smith .  rows - the indices of the rows
863bd481603SBarry Smith .  ncols - the number of columns in the matrix
864bd481603SBarry Smith .  cols - the columns indicated
865bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
866bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
867bd481603SBarry Smith 
868bd481603SBarry Smith 
869bd481603SBarry Smith    Level: intermediate
870bd481603SBarry Smith 
871bd481603SBarry Smith    Notes:
872bd481603SBarry Smith    See the chapter in the users manual on performance for more details
873bd481603SBarry Smith 
874bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
875bd481603SBarry Smith 
8761620fd73SBarry Smith    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
8771620fd73SBarry Smith 
878bd481603SBarry Smith   Concepts: preallocation^Matrix
879bd481603SBarry Smith 
880bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
881bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
882bd481603SBarry Smith M*/
883d3d32019SBarry Smith #define MatPreallocateSymmetricSet(row,nc,cols,dnz,onz) 0;\
884c1ac3661SBarry Smith { PetscInt __i; \
885d3d32019SBarry Smith   for (__i=0; __i<nc; __i++) {\
886d3d32019SBarry Smith     if (cols[__i] >= __end) onz[row - __rstart]++; \
887d3d32019SBarry Smith     else if (cols[__i] >= row) dnz[row - __rstart]++;\
888d3d32019SBarry Smith   }\
889d3d32019SBarry Smith }
890d3d32019SBarry Smith 
891bd481603SBarry Smith /*MC
89216371a99SBarry Smith    MatPreallocateLocation -  An alternative to MatPreallocationSet() that puts the nonzero locations into the matrix if it exists
89316371a99SBarry Smith 
89416371a99SBarry Smith    Synopsis:
89516371a99SBarry Smith    PetscErrorCode MatPreallocateLocations(Mat A,PetscInt row,PetscInt ncols,PetscInt *cols,PetscInt *dnz,PetscInt *onz)
89616371a99SBarry Smith 
89716371a99SBarry Smith    Not Collective
89816371a99SBarry Smith 
89916371a99SBarry Smith    Input Parameters:
90016371a99SBarry Smith .  A - matrix
90116371a99SBarry Smith .  row - row where values exist (must be local to this process)
90216371a99SBarry Smith .  ncols - number of columns
90316371a99SBarry Smith .  cols - columns with nonzeros
90416371a99SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
90516371a99SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
90616371a99SBarry Smith 
90716371a99SBarry Smith 
90816371a99SBarry Smith    Level: intermediate
90916371a99SBarry Smith 
91016371a99SBarry Smith    Notes:
91116371a99SBarry Smith    See the chapter in the users manual on performance for more details
91216371a99SBarry Smith 
91316371a99SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
91416371a99SBarry Smith 
91516371a99SBarry Smith    This is a MACRO not a function because it uses a bunch of variables private to the MatPreallocation.... routines.
91616371a99SBarry Smith 
91716371a99SBarry Smith   Concepts: preallocation^Matrix
91816371a99SBarry Smith 
91916371a99SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
92016371a99SBarry Smith           MatPreallocateSymmetricInitialize(), MatPreallocateSymmetricSetLocal()
92116371a99SBarry Smith M*/
92216371a99SBarry 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);}
92316371a99SBarry Smith 
92416371a99SBarry Smith 
92516371a99SBarry Smith /*MC
9260d2a7ff2SSatish Balay    MatPreallocateFinalize - Ends the block of code that will count the number of nonzeros per
927bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
928bd481603SBarry Smith 
929bd481603SBarry Smith    Synopsis:
930c1ac3661SBarry Smith    PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz)
931bd481603SBarry Smith 
932bd481603SBarry Smith    Collective on MPI_Comm
933bd481603SBarry Smith 
934bd481603SBarry Smith    Input Parameters:
93516371a99SBarry Smith +  dnz - the array that was be passed to the matrix preallocation routines
936bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
937bd481603SBarry Smith 
938bd481603SBarry Smith 
939bd481603SBarry Smith    Level: intermediate
940bd481603SBarry Smith 
941bd481603SBarry Smith    Notes:
942bd481603SBarry Smith    See the chapter in the users manual on performance for more details
943bd481603SBarry Smith 
944bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
945bd481603SBarry Smith 
9461620fd73SBarry Smith    This is a MACRO not a function because it closes the { started in MatPreallocateInitialize().
9471620fd73SBarry Smith 
948bd481603SBarry Smith   Concepts: preallocation^Matrix
949bd481603SBarry Smith 
950bd481603SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
951bd481603SBarry Smith           MatPreallocateSymmetricInitialize(), MatPreallocateSymmetricSetLocal()
952bd481603SBarry Smith M*/
953a89bc211SBarry Smith #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree2(dnz,onz);CHKERRQ(_4_ierr);}
9547c922b88SBarry Smith 
955bd481603SBarry Smith 
956bd481603SBarry Smith 
9577b80b807SBarry Smith /* Routines unique to particular data structures */
958be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellGetContext(Mat,void **);
959ba966639SSatish Balay PetscPolymorphicFunction(MatShellGetContext,(Mat A),(A,&t),void*,t)
960698d4c6aSKris Buschelman 
961be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInodeAdjustForInodes(Mat,IS*,IS*);
962be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *);
963698d4c6aSKris Buschelman 
964be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetColumnIndices(Mat,PetscInt[]);
965be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetColumnIndices(Mat,PetscInt[]);
966be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
967c75a6043SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
968c75a6043SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
9697b80b807SBarry Smith 
970a96a251dSBarry Smith #define MAT_SKIP_ALLOCATION -4
971a96a251dSBarry Smith 
972be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
973ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz))
974be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
975ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqSBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz))
976be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]);
977ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqAIJSetPreallocation,(Mat A,const PetscInt nnz[]),(A,0,nnz))
978273d9f13SBarry Smith 
979be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
980ba966639SSatish Balay PetscPolymorphicSubroutine(MatMPIBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[],const PetscInt onz[]),(A,bs,0,nnz,0,onz))
981be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
982be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
98353803f26SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []);
984725b52f3SLisandro Dalcin EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
985be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]);
986be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
987be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]);
988be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIDenseSetPreallocation(Mat,PetscScalar[]);
989284134d9SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation(Mat,PetscScalar[]);
990be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIRowbsSetPreallocation(Mat,PetscInt,const PetscInt[]);
991be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,PetscInt*[]);
992be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,PetscInt*[]);
993be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAdicSetLocalFunction(Mat,void (*)(void));
994273d9f13SBarry Smith 
995be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetLDA(Mat,PetscInt);
996ab92ecdeSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDenseGetLocalMatrix(Mat,Mat*);
9971b807ce4Svictorle 
998be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStoreValues(Mat);
999be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRetrieveValues(Mat);
10002e8a6d31SBarry Smith 
1001be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDAADSetCtx(Mat,void*);
10023a7fca6bSBarry Smith 
10037b80b807SBarry Smith /*
10047b80b807SBarry Smith   These routines are not usually accessed directly, rather solving is
100594b7f48cSBarry Smith   done through the KSP and PC interfaces.
10067b80b807SBarry Smith */
10077b80b807SBarry Smith 
1008d9274352SBarry Smith /*E
1009d9274352SBarry Smith     MatOrderingType - String with the name of a PETSc matrix ordering or the creation function
1010d9274352SBarry Smith        with an optional dynamic library name, for example
1011d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:orderingcreate()
1012d9274352SBarry Smith 
1013d9274352SBarry Smith    Level: beginner
1014d9274352SBarry Smith 
1015e5a9bf91SBarry Smith    Cannot use const because the PC objects manipulate the string
1016e5a9bf91SBarry Smith 
1017d9274352SBarry Smith .seealso: MatGetOrdering()
1018d9274352SBarry Smith E*/
10193eaad2caSSatish Balay #define MatOrderingType char*
1020b12f92e5SBarry Smith #define MATORDERING_NATURAL     "natural"
1021b12f92e5SBarry Smith #define MATORDERING_ND          "nd"
1022b12f92e5SBarry Smith #define MATORDERING_1WD         "1wd"
1023b12f92e5SBarry Smith #define MATORDERING_RCM         "rcm"
1024b12f92e5SBarry Smith #define MATORDERING_QMD         "qmd"
1025b12f92e5SBarry Smith #define MATORDERING_ROWLENGTH   "rowlength"
1026e106eecfSBarry Smith #define MATORDERING_DSC_ND      "dsc_nd"         /* these three are only for DSCPACK, see its documentation for details */
102762152c8bSBarry Smith #define MATORDERING_DSC_MMD     "dsc_mmd"
102862152c8bSBarry Smith #define MATORDERING_DSC_MDF     "dsc_mdf"
1029c06d978dSMatthew Knepley #define MATORDERING_CONSTRAINED "constrained"
1030c06d978dSMatthew Knepley #define MATORDERING_IDENTITY    "identity"
1031c06d978dSMatthew Knepley #define MATORDERING_REVERSE     "reverse"
10328fa24674SBarry Smith #define MATORDERING_FLOW        "flow"
1033b12f92e5SBarry Smith 
1034f2416520SSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOrdering(Mat,const MatOrderingType,IS*,IS*);
1035be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOrderingList(PetscFList *list);
1036f2416520SSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,const MatOrderingType,IS*,IS*));
103730de9b25SBarry Smith 
103830de9b25SBarry Smith /*MC
103930de9b25SBarry Smith    MatOrderingRegisterDynamic - Adds a new sparse matrix ordering to the
104030de9b25SBarry Smith                                matrix package.
104130de9b25SBarry Smith 
104230de9b25SBarry Smith    Synopsis:
1043c1ac3661SBarry Smith    PetscErrorCode MatOrderingRegisterDynamic(char *name_ordering,char *path,char *name_create,PetscErrorCode (*routine_create)(MatOrdering))
104430de9b25SBarry Smith 
104530de9b25SBarry Smith    Not Collective
104630de9b25SBarry Smith 
104730de9b25SBarry Smith    Input Parameters:
104830de9b25SBarry Smith +  sname - name of ordering (for example MATORDERING_ND)
104930de9b25SBarry Smith .  path - location of library where creation routine is
105030de9b25SBarry Smith .  name - name of function that creates the ordering type,a string
105130de9b25SBarry Smith -  function - function pointer that creates the ordering
105230de9b25SBarry Smith 
105330de9b25SBarry Smith    Level: developer
105430de9b25SBarry Smith 
105530de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
105630de9b25SBarry Smith    is ignored.
105730de9b25SBarry Smith 
105830de9b25SBarry Smith    Sample usage:
105930de9b25SBarry Smith .vb
106030de9b25SBarry Smith    MatOrderingRegisterDynamic("my_order",/home/username/my_lib/lib/libO/solaris/mylib.a,
106130de9b25SBarry Smith                "MyOrder",MyOrder);
106230de9b25SBarry Smith .ve
106330de9b25SBarry Smith 
106430de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
106530de9b25SBarry Smith $     MatOrderingSetType(part,"my_order)
106630de9b25SBarry Smith    or at runtime via the option
10672401956bSBarry Smith $     -pc_factor_mat_ordering_type my_order
106830de9b25SBarry Smith 
1069ab901514SBarry Smith    ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values.
107030de9b25SBarry Smith 
107130de9b25SBarry Smith .keywords: matrix, ordering, register
107230de9b25SBarry Smith 
107330de9b25SBarry Smith .seealso: MatOrderingRegisterDestroy(), MatOrderingRegisterAll()
107430de9b25SBarry Smith M*/
1075aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1076f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,0)
1077b12f92e5SBarry Smith #else
1078f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,d)
1079b12f92e5SBarry Smith #endif
108030de9b25SBarry Smith 
1081be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegisterDestroy(void);
1082be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegisterAll(const char[]);
10832bad1931SBarry Smith extern PetscTruth MatOrderingRegisterAllCalled;
1084b0a32e0cSBarry Smith extern PetscFList MatOrderingList;
1085d4fbbf0eSBarry Smith 
1086be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS);
1087a2ce50c7SBarry Smith 
1088d91e6319SBarry Smith /*S
10892401956bSBarry Smith    MatFactorInfo - Data passed into the matrix factorization routines
1090b00f7748SHong Zhang 
109161cad860SBarry Smith    In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE, that is use
109261cad860SBarry Smith $     MatFactorInfo  info(MAT_FACTORINFO_SIZE)
1093b00f7748SHong Zhang 
109415e8a5b3SHong Zhang    Notes: These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC.
1095b00f7748SHong Zhang 
109661cad860SBarry Smith       You can use MatFactorInfoInitialize() to set default values.
109761cad860SBarry Smith 
1098b00f7748SHong Zhang    Level: developer
1099b00f7748SHong Zhang 
1100d7d82daaSBarry Smith .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(),
1101d7d82daaSBarry Smith           MatFactorInfoInitialize()
1102b00f7748SHong Zhang 
1103b00f7748SHong Zhang S*/
1104b00f7748SHong Zhang typedef struct {
11050a29876aSHong Zhang   PetscReal     shiftnz;        /* scaling of identity added to matrix to prevent zero pivots */
1106fbf22428SSatish Balay   PetscReal     shiftpd;        /* if true, shift until positive pivots */
110715e8a5b3SHong Zhang   PetscReal     diagonal_fill;  /* force diagonal to fill in if initially not filled */
110885317021SBarry Smith   PetscReal     usedt;
110915e8a5b3SHong Zhang   PetscReal     dt;             /* drop tolerance */
1110b00f7748SHong Zhang   PetscReal     dtcol;          /* tolerance for pivoting */
111115e8a5b3SHong Zhang   PetscReal     dtcount;        /* maximum nonzeros to be allowed per row */
111267e28bfeSBarry Smith   PetscReal     fill;           /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */
1113348344bbSBarry Smith   PetscReal     levels;         /* ICC/ILU(levels) */
1114bcd9e38bSBarry Smith   PetscReal     pivotinblocks;  /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0
1115bcd9e38bSBarry Smith                                    factorization may be faster if do not pivot */
111662bba022SBarry Smith   PetscReal     shiftinblocks;  /* if block in block factorization has zero pivot then shift diagonal until non-singular */
111715e8a5b3SHong Zhang   PetscReal     zeropivot;      /* pivot is called zero if less than this */
111815e8a5b3SHong Zhang } MatFactorInfo;
1119ffa6d0a5SLois Curfman McInnes 
1120be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFactorInfoInitialize(MatFactorInfo*);
11210481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactor(Mat,IS,const MatFactorInfo*);
11220481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
11230481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactorNumeric(Mat,Mat,const MatFactorInfo*);
11240481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactor(Mat,IS,IS,const MatFactorInfo*);
11250481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUFactor(Mat,IS,IS,const MatFactorInfo*);
11260481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
11270481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
11280481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatICCFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
11290481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatICCFactor(Mat,IS,const MatFactorInfo*);
11300481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactorNumeric(Mat,Mat,const MatFactorInfo*);
11310481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUDTFactor(Mat,IS,IS,const MatFactorInfo*,Mat *);
11323c2a7987SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUDTFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
11333c2a7987SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUDTFactorNumeric(Mat,Mat,const MatFactorInfo*);
1134be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*);
1135be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolve(Mat,Vec,Vec);
1136be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatForwardSolve(Mat,Vec,Vec);
1137be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatBackwardSolve(Mat,Vec,Vec);
1138be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveAdd(Mat,Vec,Vec,Vec);
1139be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveTranspose(Mat,Vec,Vec);
1140be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveTransposeAdd(Mat,Vec,Vec,Vec);
1141be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolves(Mat,Vecs,Vecs);
11428ed539a5SBarry Smith 
1143be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUnfactored(Mat);
1144bb5a7306SBarry Smith 
1145d91e6319SBarry Smith /*E
1146d91e6319SBarry Smith     MatSORType - What type of (S)SOR to perform
1147bb1eb677SSatish Balay 
1148d91e6319SBarry Smith     Level: beginner
1149d91e6319SBarry Smith 
1150d9274352SBarry Smith    May be bitwise ORd together
1151d9274352SBarry Smith 
1152d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
1153d91e6319SBarry Smith 
11544e7234bfSBarry Smith    MatSORType may be bitwise ORd together, so do not change the numbers
11554e7234bfSBarry Smith 
115641f059aeSBarry Smith .seealso: MatSOR()
1157d91e6319SBarry Smith E*/
1158ee50ffe9SBarry Smith typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3,
1159ee50ffe9SBarry Smith               SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8,
1160ee50ffe9SBarry Smith               SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16,
116184cb2905SBarry Smith               SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType;
116241f059aeSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSOR(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
11638ed539a5SBarry Smith 
1164d4fbbf0eSBarry Smith /*
1165639f9d9dSBarry Smith     These routines are for efficiently computing Jacobians via finite differences.
1166639f9d9dSBarry Smith */
1167b12f92e5SBarry Smith 
1168d9274352SBarry Smith /*E
1169d9274352SBarry Smith     MatColoringType - String with the name of a PETSc matrix coloring or the creation function
1170d9274352SBarry Smith        with an optional dynamic library name, for example
1171d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:coloringcreate()
1172d9274352SBarry Smith 
1173d9274352SBarry Smith    Level: beginner
1174d9274352SBarry Smith 
1175d9274352SBarry Smith .seealso: MatGetColoring()
1176d9274352SBarry Smith E*/
1177a313700dSBarry Smith #define MatColoringType char*
1178b12f92e5SBarry Smith #define MATCOLORING_NATURAL "natural"
1179b12f92e5SBarry Smith #define MATCOLORING_SL      "sl"
1180b12f92e5SBarry Smith #define MATCOLORING_LF      "lf"
1181b12f92e5SBarry Smith #define MATCOLORING_ID      "id"
1182b12f92e5SBarry Smith 
1183ab970d9bSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColoring(Mat,const MatColoringType,ISColoring*);
11842e90c967SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,MatColoringType,ISColoring *));
118530de9b25SBarry Smith 
118630de9b25SBarry Smith /*MC
118730de9b25SBarry Smith    MatColoringRegisterDynamic - Adds a new sparse matrix coloring to the
118830de9b25SBarry Smith                                matrix package.
118930de9b25SBarry Smith 
119030de9b25SBarry Smith    Synopsis:
1191c1ac3661SBarry Smith    PetscErrorCode MatColoringRegisterDynamic(char *name_coloring,char *path,char *name_create,PetscErrorCode (*routine_create)(MatColoring))
119230de9b25SBarry Smith 
119330de9b25SBarry Smith    Not Collective
119430de9b25SBarry Smith 
119530de9b25SBarry Smith    Input Parameters:
119630de9b25SBarry Smith +  sname - name of Coloring (for example MATCOLORING_SL)
119730de9b25SBarry Smith .  path - location of library where creation routine is
119830de9b25SBarry Smith .  name - name of function that creates the Coloring type, a string
119930de9b25SBarry Smith -  function - function pointer that creates the coloring
120030de9b25SBarry Smith 
120130de9b25SBarry Smith    Level: developer
120230de9b25SBarry Smith 
120330de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
120430de9b25SBarry Smith    is ignored.
120530de9b25SBarry Smith 
120630de9b25SBarry Smith    Sample usage:
120730de9b25SBarry Smith .vb
120830de9b25SBarry Smith    MatColoringRegisterDynamic("my_color",/home/username/my_lib/lib/libO/solaris/mylib.a,
120930de9b25SBarry Smith                "MyColor",MyColor);
121030de9b25SBarry Smith .ve
121130de9b25SBarry Smith 
121230de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
121330de9b25SBarry Smith $     MatColoringSetType(part,"my_color")
121430de9b25SBarry Smith    or at runtime via the option
121530de9b25SBarry Smith $     -mat_coloring_type my_color
121630de9b25SBarry Smith 
1217ab901514SBarry Smith    $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
121830de9b25SBarry Smith 
121930de9b25SBarry Smith .keywords: matrix, Coloring, register
122030de9b25SBarry Smith 
122130de9b25SBarry Smith .seealso: MatColoringRegisterDestroy(), MatColoringRegisterAll()
122230de9b25SBarry Smith M*/
1223aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1224f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,0)
1225b12f92e5SBarry Smith #else
1226f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,d)
1227b12f92e5SBarry Smith #endif
122830de9b25SBarry Smith 
12292bad1931SBarry Smith extern PetscTruth MatColoringRegisterAllCalled;
1230f1144a30SSatish Balay 
1231be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegisterAll(const char[]);
1232be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegisterDestroy(void);
1233be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*);
1234639f9d9dSBarry Smith 
1235d9274352SBarry Smith /*S
1236d9274352SBarry Smith      MatFDColoring - Object for computing a sparse Jacobian via finite differences
1237d9274352SBarry Smith         and coloring
1238639f9d9dSBarry Smith 
1239d9274352SBarry Smith    Level: beginner
1240d9274352SBarry Smith 
1241d9274352SBarry Smith   Concepts: coloring, sparse Jacobian, finite differences
1242d9274352SBarry Smith 
1243d9274352SBarry Smith .seealso:  MatFDColoringCreate()
1244d9274352SBarry Smith S*/
1245e2a1c21fSSatish Balay typedef struct _p_MatFDColoring* MatFDColoring;
1246639f9d9dSBarry Smith 
1247be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringCreate(Mat,ISColoring,MatFDColoring *);
1248be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringDestroy(MatFDColoring);
1249be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringView(MatFDColoring,PetscViewer);
1250be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*);
12514a9d489dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**);
1252be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal);
1253be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFromOptions(MatFDColoring);
1254be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApply(Mat,MatFDColoring,Vec,MatStructure*,void *);
1255be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApplyTS(Mat,MatFDColoring,PetscReal,Vec,MatStructure*,void *);
1256be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetF(MatFDColoring,Vec);
1257be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,PetscInt*[]);
1258639f9d9dSBarry Smith /*
12590752156aSBarry Smith     These routines are for partitioning matrices: currently used only
12603eda8832SBarry Smith   for adjacency matrix, MatCreateMPIAdj().
12610752156aSBarry Smith */
1262ca161407SBarry Smith 
1263d9274352SBarry Smith /*S
1264d9274352SBarry Smith      MatPartitioning - Object for managing the partitioning of a matrix or graph
1265d9274352SBarry Smith 
1266d9274352SBarry Smith    Level: beginner
1267d9274352SBarry Smith 
1268d9274352SBarry Smith   Concepts: partitioning
1269d9274352SBarry Smith 
1270743a1026Svictorle .seealso:  MatPartitioningCreate(), MatPartitioningType
1271d9274352SBarry Smith S*/
127291e9ee9fSBarry Smith typedef struct _p_MatPartitioning* MatPartitioning;
1273d9274352SBarry Smith 
1274d9274352SBarry Smith /*E
12755bc6e7ccSvictorle     MatPartitioningType - String with the name of a PETSc matrix partitioning or the creation function
1276d9274352SBarry Smith        with an optional dynamic library name, for example
1277d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:partitioningcreate()
1278d9274352SBarry Smith 
1279d9274352SBarry Smith    Level: beginner
1280d9274352SBarry Smith 
1281b547a13bSvictorle .seealso: MatPartitioningCreate(), MatPartitioning
1282d9274352SBarry Smith E*/
1283a313700dSBarry Smith #define MatPartitioningType char*
12848ba1e511SMatthew Knepley #define MAT_PARTITIONING_CURRENT  "current"
1285c2ec2a73SHong Zhang #define MAT_PARTITIONING_SQUARE   "square"
12868ba1e511SMatthew Knepley #define MAT_PARTITIONING_PARMETIS "parmetis"
128771306c60Sjroman #define MAT_PARTITIONING_CHACO    "chaco"
128871306c60Sjroman #define MAT_PARTITIONING_JOSTLE   "jostle"
128971306c60Sjroman #define MAT_PARTITIONING_PARTY    "party"
129071306c60Sjroman #define MAT_PARTITIONING_SCOTCH   "scotch"
129171306c60Sjroman 
1292ca161407SBarry Smith 
1293be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningCreate(MPI_Comm,MatPartitioning*);
1294a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetType(MatPartitioning,const MatPartitioningType);
1295be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetNParts(MatPartitioning,PetscInt);
1296be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetAdjacency(MatPartitioning,Mat);
1297be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]);
1298be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []);
1299be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningApply(MatPartitioning,IS*);
1300be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningDestroy(MatPartitioning);
13012aabb6bbSBarry Smith 
1302be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatPartitioning));
130330de9b25SBarry Smith 
130430de9b25SBarry Smith /*MC
130530de9b25SBarry Smith    MatPartitioningRegisterDynamic - Adds a new sparse matrix partitioning to the
130630de9b25SBarry Smith    matrix package.
130730de9b25SBarry Smith 
130830de9b25SBarry Smith    Synopsis:
1309c1ac3661SBarry Smith    PetscErrorCode MatPartitioningRegisterDynamic(char *name_partitioning,char *path,char *name_create,PetscErrorCode (*routine_create)(MatPartitioning))
131030de9b25SBarry Smith 
131130de9b25SBarry Smith    Not Collective
131230de9b25SBarry Smith 
131330de9b25SBarry Smith    Input Parameters:
131430de9b25SBarry Smith +  sname - name of partitioning (for example MAT_PARTITIONING_CURRENT) or parmetis
131530de9b25SBarry Smith .  path - location of library where creation routine is
131630de9b25SBarry Smith .  name - name of function that creates the partitioning type, a string
131730de9b25SBarry Smith -  function - function pointer that creates the partitioning type
131830de9b25SBarry Smith 
131930de9b25SBarry Smith    Level: developer
132030de9b25SBarry Smith 
132130de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
132230de9b25SBarry Smith    is ignored.
132330de9b25SBarry Smith 
132430de9b25SBarry Smith    Sample usage:
132530de9b25SBarry Smith .vb
132630de9b25SBarry Smith    MatPartitioningRegisterDynamic("my_part",/home/username/my_lib/lib/libO/solaris/mylib.a,
132730de9b25SBarry Smith                "MyPartCreate",MyPartCreate);
132830de9b25SBarry Smith .ve
132930de9b25SBarry Smith 
133030de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
133130de9b25SBarry Smith $     MatPartitioningSetType(part,"my_part")
133230de9b25SBarry Smith    or at runtime via the option
133330de9b25SBarry Smith $     -mat_partitioning_type my_part
133430de9b25SBarry Smith 
1335ab901514SBarry Smith    $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
133630de9b25SBarry Smith 
133730de9b25SBarry Smith .keywords: matrix, partitioning, register
133830de9b25SBarry Smith 
133930de9b25SBarry Smith .seealso: MatPartitioningRegisterDestroy(), MatPartitioningRegisterAll()
134030de9b25SBarry Smith M*/
1341aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1342f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,0)
13432aabb6bbSBarry Smith #else
1344f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,d)
13452aabb6bbSBarry Smith #endif
13462aabb6bbSBarry Smith 
13472bad1931SBarry Smith extern PetscTruth MatPartitioningRegisterAllCalled;
1348f1144a30SSatish Balay 
1349be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegisterAll(const char[]);
1350be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegisterDestroy(void);
13512bad1931SBarry Smith 
1352be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningView(MatPartitioning,PetscViewer);
1353be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetFromOptions(MatPartitioning);
1354a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningGetType(MatPartitioning,const MatPartitioningType*);
1355ca161407SBarry Smith 
1356be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningParmetisSetCoarseSequential(MatPartitioning);
13570e8af3caSMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningParmetisGetEdgeCut(MatPartitioning, PetscInt *);
13580752156aSBarry Smith 
1359be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningJostleSetCoarseLevel(MatPartitioning,PetscReal);
1360be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningJostleSetCoarseSequential(MatPartitioning);
136171306c60Sjroman 
1362954b3ec6SBarry Smith typedef enum {MP_CHACO_MULTILEVEL_KL,MP_CHACO_SPECTRAL,MP_CHACO_LINEAR,MP_CHACO_RANDOM, MP_CHACO_SCATTERED} MPChacoGlobalType;
1363be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetGlobal(MatPartitioning, MPChacoGlobalType);
136471306c60Sjroman typedef enum { MP_CHACO_KERNIGHAN_LIN, MP_CHACO_NONE } MPChacoLocalType;
1365be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetLocal(MatPartitioning, MPChacoLocalType);
1366be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal);
136771306c60Sjroman typedef enum { MP_CHACO_LANCZOS, MP_CHACO_RQI_SYMMLQ } MPChacoEigenType;
1368be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType);
1369be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenTol(MatPartitioning, PetscReal);
1370be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenNumber(MatPartitioning, PetscInt);
137171306c60Sjroman 
137271306c60Sjroman #define MP_PARTY_OPT "opt"
137371306c60Sjroman #define MP_PARTY_LIN "lin"
137471306c60Sjroman #define MP_PARTY_SCA "sca"
137571306c60Sjroman #define MP_PARTY_RAN "ran"
137671306c60Sjroman #define MP_PARTY_GBF "gbf"
137771306c60Sjroman #define MP_PARTY_GCF "gcf"
137871306c60Sjroman #define MP_PARTY_BUB "bub"
137971306c60Sjroman #define MP_PARTY_DEF "def"
1380be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetGlobal(MatPartitioning, const char*);
138171306c60Sjroman #define MP_PARTY_HELPFUL_SETS "hs"
138271306c60Sjroman #define MP_PARTY_KERNIGHAN_LIN "kl"
138371306c60Sjroman #define MP_PARTY_NONE "no"
1384be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetLocal(MatPartitioning, const char*);
1385be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal);
1386be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetBipart(MatPartitioning,PetscTruth);
1387be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscTruth);
138871306c60Sjroman 
138971306c60Sjroman typedef enum { MP_SCOTCH_GREEDY, MP_SCOTCH_GPS, MP_SCOTCH_GR_GPS } MPScotchGlobalType;
1390be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetArch(MatPartitioning,const char*);
1391be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetMultilevel(MatPartitioning);
1392be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetGlobal(MatPartitioning,MPScotchGlobalType);
1393be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetCoarseLevel(MatPartitioning,PetscReal);
1394be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetHostList(MatPartitioning,const char*);
139571306c60Sjroman typedef enum { MP_SCOTCH_KERNIGHAN_LIN, MP_SCOTCH_NONE } MPScotchLocalType;
1396be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetLocal(MatPartitioning,MPScotchLocalType);
1397be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetMapping(MatPartitioning);
1398be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetStrategy(MatPartitioning,char*);
139971306c60Sjroman 
1400591294e4SBarry Smith EXTERN PetscErrorCode MatMeshToVertexGraph(Mat,PetscInt,Mat*);
140114196391SBarry Smith EXTERN PetscErrorCode MatMeshToCellGraph(Mat,PetscInt,Mat*);
1402591294e4SBarry Smith 
14030752156aSBarry Smith /*
14040a835dfdSSatish Balay     If you add entries here you must also add them to finclude/petscmat.h
1405d4fbbf0eSBarry Smith */
14061c1c02c0SLois Curfman McInnes typedef enum { MATOP_SET_VALUES=0,
14071c1c02c0SLois Curfman McInnes                MATOP_GET_ROW=1,
14081c1c02c0SLois Curfman McInnes                MATOP_RESTORE_ROW=2,
14091c1c02c0SLois Curfman McInnes                MATOP_MULT=3,
14101c1c02c0SLois Curfman McInnes                MATOP_MULT_ADD=4,
14117c922b88SBarry Smith                MATOP_MULT_TRANSPOSE=5,
14127c922b88SBarry Smith                MATOP_MULT_TRANSPOSE_ADD=6,
14131c1c02c0SLois Curfman McInnes                MATOP_SOLVE=7,
14141c1c02c0SLois Curfman McInnes                MATOP_SOLVE_ADD=8,
14157c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE=9,
14167c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE_ADD=10,
14171c1c02c0SLois Curfman McInnes                MATOP_LUFACTOR=11,
14181c1c02c0SLois Curfman McInnes                MATOP_CHOLESKYFACTOR=12,
14191c1c02c0SLois Curfman McInnes                MATOP_RELAX=13,
14201c1c02c0SLois Curfman McInnes                MATOP_TRANSPOSE=14,
14211c1c02c0SLois Curfman McInnes                MATOP_GETINFO=15,
14221c1c02c0SLois Curfman McInnes                MATOP_EQUAL=16,
14231c1c02c0SLois Curfman McInnes                MATOP_GET_DIAGONAL=17,
14241c1c02c0SLois Curfman McInnes                MATOP_DIAGONAL_SCALE=18,
14251c1c02c0SLois Curfman McInnes                MATOP_NORM=19,
14261c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_BEGIN=20,
14271c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_END=21,
1428d519adbfSMatthew Knepley                MATOP_SET_OPTION=22,
1429d519adbfSMatthew Knepley                MATOP_ZERO_ENTRIES=23,
1430d519adbfSMatthew Knepley                MATOP_ZERO_ROWS=24,
1431d519adbfSMatthew Knepley                MATOP_LUFACTOR_SYMBOLIC=25,
1432d519adbfSMatthew Knepley                MATOP_LUFACTOR_NUMERIC=26,
1433d519adbfSMatthew Knepley                MATOP_CHOLESKY_FACTOR_SYMBOLIC=27,
1434d519adbfSMatthew Knepley                MATOP_CHOLESKY_FACTOR_NUMERIC=28,
1435d519adbfSMatthew Knepley                MATOP_SETUP_PREALLOCATION=29,
1436d519adbfSMatthew Knepley                MATOP_ILUFACTOR_SYMBOLIC=30,
1437d519adbfSMatthew Knepley                MATOP_ICCFACTOR_SYMBOLIC=31,
1438d519adbfSMatthew Knepley                MATOP_GET_ARRAY=32,
1439d519adbfSMatthew Knepley                MATOP_RESTORE_ARRAY=33,
1440d519adbfSMatthew Knepley                MATOP_DUPLICATE=34,
1441d519adbfSMatthew Knepley                MATOP_FORWARD_SOLVE=35,
1442d519adbfSMatthew Knepley                MATOP_BACKWARD_SOLVE=36,
1443d519adbfSMatthew Knepley                MATOP_ILUFACTOR=37,
1444d519adbfSMatthew Knepley                MATOP_ICCFACTOR=38,
1445d519adbfSMatthew Knepley                MATOP_AXPY=39,
1446d519adbfSMatthew Knepley                MATOP_GET_SUBMATRICES=40,
1447d519adbfSMatthew Knepley                MATOP_INCREASE_OVERLAP=41,
1448d519adbfSMatthew Knepley                MATOP_GET_VALUES=42,
1449d519adbfSMatthew Knepley                MATOP_COPY=43,
1450d519adbfSMatthew Knepley                MATOP_GET_ROW_MAX=44,
1451d519adbfSMatthew Knepley                MATOP_SCALE=45,
1452d519adbfSMatthew Knepley                MATOP_SHIFT=46,
1453d519adbfSMatthew Knepley                MATOP_DIAGONAL_SHIFT=47,
1454d519adbfSMatthew Knepley                MATOP_ILUDT_FACTOR=48,
1455d519adbfSMatthew Knepley                MATOP_SET_BLOCK_SIZE=49,
1456d519adbfSMatthew Knepley                MATOP_GET_ROW_IJ=50,
1457d519adbfSMatthew Knepley                MATOP_RESTORE_ROW_IJ=51,
1458d519adbfSMatthew Knepley                MATOP_GET_COLUMN_IJ=52,
1459d519adbfSMatthew Knepley                MATOP_RESTORE_COLUMN_IJ=53,
1460d519adbfSMatthew Knepley                MATOP_FDCOLORING_CREATE=54,
1461d519adbfSMatthew Knepley                MATOP_COLORING_PATCH=55,
1462d519adbfSMatthew Knepley                MATOP_SET_UNFACTORED=56,
1463d519adbfSMatthew Knepley                MATOP_PERMUTE=57,
1464d519adbfSMatthew Knepley                MATOP_SET_VALUES_BLOCKED=58,
1465d519adbfSMatthew Knepley                MATOP_GET_SUBMATRIX=59,
1466d519adbfSMatthew Knepley                MATOP_DESTROY=60,
1467d519adbfSMatthew Knepley                MATOP_VIEW=61,
1468d519adbfSMatthew Knepley                MATOP_CONVERT_FROM=62,
1469d519adbfSMatthew Knepley                MATOP_USE_SCALED_FORM=63,
1470d519adbfSMatthew Knepley                MATOP_SCALE_SYSTEM=64,
1471d519adbfSMatthew Knepley                MATOP_UNSCALE_SYSTEM=65,
1472d519adbfSMatthew Knepley                MATOP_SET_LOCAL_TO_GLOBAL_MAP=66,
1473d519adbfSMatthew Knepley                MATOP_SET_VALUES_LOCAL=67,
1474d519adbfSMatthew Knepley                MATOP_ZERO_ROWS_LOCAL=68,
1475d519adbfSMatthew Knepley                MATOP_GET_ROW_MAX_ABS=69,
1476d519adbfSMatthew Knepley                MATOP_GET_ROW_MIN_ABS=70,
1477d519adbfSMatthew Knepley                MATOP_CONVERT=71,
1478d519adbfSMatthew Knepley                MATOP_SET_COLORING=72,
1479d519adbfSMatthew Knepley                MATOP_SET_VALUES_ADIC=73,
1480d519adbfSMatthew Knepley                MATOP_SET_VALUES_ADIFOR=74,
1481d519adbfSMatthew Knepley                MATOP_FD_COLORING_APPLY=75,
1482d519adbfSMatthew Knepley                MATOP_SET_FROM_OPTIONS=76,
1483d519adbfSMatthew Knepley                MATOP_MULT_CON=77,
1484d519adbfSMatthew Knepley                MATOP_MULT_TRANSPOSE_CON=78,
1485d519adbfSMatthew Knepley                MATOP_PERMUTE_SPARSIFY=79,
1486d519adbfSMatthew Knepley                MATOP_MULT_MULTIPLE=80,
1487d519adbfSMatthew Knepley                MATOP_SOLVE_MULTIPLE=81,
1488d519adbfSMatthew Knepley                MATOP_GET_INERTIA=82,
1489d519adbfSMatthew Knepley                MATOP_LOAD=83,
1490d519adbfSMatthew Knepley                MATOP_IS_SYMMETRIC=84,
1491d519adbfSMatthew Knepley                MATOP_IS_HERMITIAN=85,
1492d519adbfSMatthew Knepley                MATOP_IS_STRUCTURALLY_SYMMETRIC=86,
149341f059aeSBarry Smith                MATOP_DUMMY=87,
1494d519adbfSMatthew Knepley                MATOP_GET_VECS=88,
1495d519adbfSMatthew Knepley                MATOP_MAT_MULT=89,
1496d519adbfSMatthew Knepley                MATOP_MAT_MULT_SYMBOLIC=90,
1497d519adbfSMatthew Knepley                MATOP_MAT_MULT_NUMERIC=91,
1498d519adbfSMatthew Knepley                MATOP_PTAP=92,
1499d519adbfSMatthew Knepley                MATOP_PTAP_SYMBOLIC=93,
1500d519adbfSMatthew Knepley                MATOP_PTAP_NUMERIC=94,
1501d519adbfSMatthew Knepley                MATOP_MAT_MULTTRANSPOSE=95,
1502d519adbfSMatthew Knepley                MATOP_MAT_MULTTRANSPOSE_SYMBOLIC=96,
1503d519adbfSMatthew Knepley                MATOP_MAT_MULTTRANSPOSE_NUMERIC=97,
1504d519adbfSMatthew Knepley                MATOP_PTAP_SYMBOLIC_SEQAIJ=98,
1505d519adbfSMatthew Knepley                MATOP_PTAP_NUMERIC_SEQAIJ=99,
1506d519adbfSMatthew Knepley                MATOP_PTAP_SYMBOLIC_MPIAIJ=100,
1507d519adbfSMatthew Knepley                MATOP_PTAP_NUMERIC_MPIAIJ=101,
1508d519adbfSMatthew Knepley                MATOP_CONJUGATE=102,
1509d519adbfSMatthew Knepley                MATOP_SET_SIZES=103,
1510d519adbfSMatthew Knepley                MATOP_SET_VALUES_ROW = 104,
1511d519adbfSMatthew Knepley                MATOP_REAL_PART=105,
1512d519adbfSMatthew Knepley                MATOP_IMAG_PART=106,
1513d519adbfSMatthew Knepley                MATOP_GET_ROW_UTRIANGULAR=107,
1514d519adbfSMatthew Knepley                MATOP_RESTORE_ROW_UTRIANGULAR=108,
1515d519adbfSMatthew Knepley                MATOP_MATSOLVE=109,
1516d519adbfSMatthew Knepley                MATOP_GET_REDUNDANTMATRIX=110,
1517d519adbfSMatthew Knepley                MATOP_GET_ROW_MIN=111,
1518d519adbfSMatthew Knepley                MATOP_GET_COLUMN_VEC=112,
1519d519adbfSMatthew Knepley                MATOP_MISSING_DIAGONAL=113,
1520d519adbfSMatthew Knepley                MATOP_MATGETSEQNONZEROSTRUCTURE=114,
152189c6957cSBarry Smith                MATOP_CREATE=115,
152289c6957cSBarry Smith                MATOP_GET_GHOSTS=116,
152389c6957cSBarry Smith                MATOP_ILUDTFACTORSYMBOLIC=117,
152489c6957cSBarry Smith                MATOP_ILUDTFACTORNUMERIC=118,
1525*eeffb40dSHong Zhang                MATOP_MULT_DIAGONAL_BLOCK=119,
1526*eeffb40dSHong Zhang                MATOP_HERMITIANTRANSPOSE=120
1527fae171e0SBarry Smith              } MatOperation;
1528be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatHasOperation(Mat,MatOperation,PetscTruth*);
1529be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellSetOperation(Mat,MatOperation,void(*)(void));
1530be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellGetOperation(Mat,MatOperation,void(**)(void));
1531be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellSetContext(Mat,void*);
1532112a2221SBarry Smith 
153390ace30eSBarry Smith /*
153490ace30eSBarry Smith    Codes for matrices stored on disk. By default they are
153590ace30eSBarry Smith  stored in a universal format. By changing the format with
15367973ad04SBarry Smith  PetscViewerSetFormat(viewer,PETSC_VIEWER_NATIVE); the matrices will
153790ace30eSBarry Smith  be stored in a way natural for the matrix, for example dense matrices
153890ace30eSBarry Smith  would be stored as dense. Matrices stored this way may only be
153990ace30eSBarry Smith  read into matrices of the same time.
154090ace30eSBarry Smith */
154190ace30eSBarry Smith #define MATRIX_BINARY_FORMAT_DENSE -1
154290ace30eSBarry Smith 
1543be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetHashTableFactor(Mat,PetscReal);
1544be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatISGetLocalMat(Mat,Mat*);
15451f4e1ec7SBarry Smith 
1546d9274352SBarry Smith /*S
1547d9274352SBarry Smith      MatNullSpace - Object that removes a null space from a vector, i.e.
1548d9274352SBarry Smith          orthogonalizes the vector to a subsapce
1549d9274352SBarry Smith 
1550f7a9e4ceSBarry Smith    Level: advanced
1551d9274352SBarry Smith 
1552d9274352SBarry Smith   Concepts: matrix; linear operator, null space
1553d9274352SBarry Smith 
15546e1639daSBarry Smith   Users manual sections:
15556e1639daSBarry Smith .   sec_singular
15566e1639daSBarry Smith 
1557d9274352SBarry Smith .seealso:  MatNullSpaceCreate()
1558d9274352SBarry Smith S*/
155974637425SBarry Smith typedef struct _p_MatNullSpace* MatNullSpace;
1560d9274352SBarry Smith 
1561be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceCreate(MPI_Comm,PetscTruth,PetscInt,const Vec[],MatNullSpace*);
1562b22b330cSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(MatNullSpace,Vec,void*),void*);
1563be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceDestroy(MatNullSpace);
1564be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceRemove(MatNullSpace,Vec,Vec*);
1565be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceAttach(Mat,MatNullSpace);
156695902228SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceTest(MatNullSpace,Mat,PetscTruth *);
156774637425SBarry Smith 
1568be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatReorderingSeqSBAIJ(Mat,IS);
1569be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPISBAIJSetHashTableFactor(Mat,PetscReal);
1570be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetColumnIndices(Mat,PetscInt *);
157146129b97SKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJInvertBlockDiagonal(Mat);
15723f1d51d7SBarry Smith 
1573be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMAIJ(Mat,PetscInt,Mat*);
1574be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMAIJRedimension(Mat,PetscInt,Mat*);
1575be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMAIJGetAIJ(Mat,Mat*);
1576c4f061fbSSatish Balay 
1577be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatComputeExplicitOperator(Mat,Mat*);
1578b0a32e0cSBarry Smith 
1579be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalScaleLocal(Mat,Vec);
158004f1ad80SBarry Smith 
1581fef1beadSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*);
15823ec795f1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetBase(Mat,Vec,Vec);
15833ec795f1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*);
1584e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*));
1585e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec));
1586e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDAddNullSpace(Mat,MatNullSpace);
1587e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt);
1588e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDResetHHistory(Mat);
1589e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctionError(Mat,PetscReal);
1590e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetPeriod(Mat,PetscInt);
1591e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDGetH(Mat,PetscScalar *);
15926aa9148fSLisandro Dalcin EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetOptionsPrefix(Mat,const char[]);
1593e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFromOptions(Mat);
1594e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*);
1595e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*);
1596e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDDSSetUmin(Mat,PetscReal);
1597e884886fSBarry Smith 
1598e884886fSBarry Smith 
1599e884886fSBarry Smith typedef struct _p_MatMFFD* MatMFFD;
1600e884886fSBarry Smith 
1601e884886fSBarry Smith /*E
1602e884886fSBarry Smith     MatMFFDType - algorithm used to compute the h used in computing matrix-vector products via differencing of the function
1603e884886fSBarry Smith 
1604e884886fSBarry Smith    Level: beginner
1605e884886fSBarry Smith 
1606e884886fSBarry Smith .seealso: MatMFFDSetType(), MatMFFDRegister()
1607e884886fSBarry Smith E*/
1608a313700dSBarry Smith #define MatMFFDType char*
1609e884886fSBarry Smith #define MATMFFD_DS  "ds"
1610e884886fSBarry Smith #define MATMFFD_WP  "wp"
1611e884886fSBarry Smith 
1612a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetType(Mat,const MatMFFDType);
1613e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatMFFD));
1614e884886fSBarry Smith 
1615e884886fSBarry Smith /*MC
1616e884886fSBarry Smith    MatMFFDRegisterDynamic - Adds a method to the MatMFFD registry.
1617e884886fSBarry Smith 
1618e884886fSBarry Smith    Synopsis:
1619e884886fSBarry Smith    PetscErrorCode MatMFFDRegisterDynamic(char *name_solver,char *path,char *name_create,PetscErrorCode (*routine_create)(MatMFFD))
1620e884886fSBarry Smith 
1621e884886fSBarry Smith    Not Collective
1622e884886fSBarry Smith 
1623e884886fSBarry Smith    Input Parameters:
1624e884886fSBarry Smith +  name_solver - name of a new user-defined compute-h module
1625e884886fSBarry Smith .  path - path (either absolute or relative) the library containing this solver
1626e884886fSBarry Smith .  name_create - name of routine to create method context
1627e884886fSBarry Smith -  routine_create - routine to create method context
1628e884886fSBarry Smith 
1629e884886fSBarry Smith    Level: developer
1630e884886fSBarry Smith 
1631e884886fSBarry Smith    Notes:
1632e884886fSBarry Smith    MatMFFDRegisterDynamic() may be called multiple times to add several user-defined solvers.
1633e884886fSBarry Smith 
1634e884886fSBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
1635e884886fSBarry Smith    is ignored.
1636e884886fSBarry Smith 
1637e884886fSBarry Smith    Sample usage:
1638e884886fSBarry Smith .vb
1639e884886fSBarry Smith    MatMFFDRegisterDynamic("my_h",/home/username/my_lib/lib/libO/solaris/mylib.a,
1640e884886fSBarry Smith                "MyHCreate",MyHCreate);
1641e884886fSBarry Smith .ve
1642e884886fSBarry Smith 
1643e884886fSBarry Smith    Then, your solver can be chosen with the procedural interface via
1644e884886fSBarry Smith $     MatMFFDSetType(mfctx,"my_h")
1645e884886fSBarry Smith    or at runtime via the option
1646e884886fSBarry Smith $     -snes_mf_type my_h
1647e884886fSBarry Smith 
1648e884886fSBarry Smith .keywords: MatMFFD, register
1649e884886fSBarry Smith 
1650e884886fSBarry Smith .seealso: MatMFFDRegisterAll(), MatMFFDRegisterDestroy()
1651e884886fSBarry Smith M*/
1652e884886fSBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1653e884886fSBarry Smith #define MatMFFDRegisterDynamic(a,b,c,d) MatMFFDRegister(a,b,c,0)
1654e884886fSBarry Smith #else
1655e884886fSBarry Smith #define MatMFFDRegisterDynamic(a,b,c,d) MatMFFDRegister(a,b,c,d)
1656e884886fSBarry Smith #endif
1657e884886fSBarry Smith 
1658e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDRegisterAll(const char[]);
1659e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDRegisterDestroy(void);
1660e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDDefaultSetUmin(Mat,PetscReal);
1661e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDWPSetComputeNormU(Mat,PetscTruth);
1662e884886fSBarry Smith 
1663e884886fSBarry Smith 
1664be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *);
1665be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *);
16667dbadf16SMatthew Knepley 
1667e9fa29b7SSatish Balay PETSC_EXTERN_CXX_END
16682eac72dbSBarry Smith #endif
1669