xref: /petsc/include/petscmat.h (revision 014dd563d73e9fc78d056590fa6cf997782bf92d)
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"
72eac72dbSBarry Smith 
8d9274352SBarry Smith /*S
9d9274352SBarry Smith      Mat - Abstract PETSc matrix object
102eac72dbSBarry Smith 
11d91e6319SBarry Smith    Level: beginner
12d91e6319SBarry Smith 
13d9274352SBarry Smith   Concepts: matrix; linear operator
14d9274352SBarry Smith 
15d9274352SBarry Smith .seealso:  MatCreate(), MatType, MatSetType()
16d9274352SBarry Smith S*/
17d9274352SBarry Smith typedef struct _p_Mat*           Mat;
18d9274352SBarry Smith 
1976bdecfbSBarry Smith /*J
20d9274352SBarry Smith     MatType - String with the name of a PETSc matrix or the creation function
21d9274352SBarry Smith        with an optional dynamic library name, for example
22d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:mymatcreate()
23d9274352SBarry Smith 
24d9274352SBarry Smith    Level: beginner
25d9274352SBarry Smith 
26c7393fdbSBarry Smith .seealso: MatSetType(), Mat, MatSolverPackage
2776bdecfbSBarry Smith J*/
28a313700dSBarry Smith #define MatType char*
29273d9f13SBarry Smith #define MATSAME            "same"
305a11e1b2SBarry Smith #define MATMAIJ            "maij"
31273d9f13SBarry Smith #define MATSEQMAIJ         "seqmaij"
32273d9f13SBarry Smith #define MATMPIMAIJ         "mpimaij"
33273d9f13SBarry Smith #define MATIS              "is"
345a11e1b2SBarry Smith #define MATAIJ             "aij"
35273d9f13SBarry Smith #define MATSEQAIJ          "seqaij"
367d6a0e61SBarry Smith #define MATSEQAIJPTHREAD   "seqaijpthread"
37bf2c1783SBarry Smith #define MATAIJPTHREAD      "aijpthread"
38273d9f13SBarry Smith #define MATMPIAIJ          "mpiaij"
395a11e1b2SBarry Smith #define MATAIJCRL          "aijcrl"
405a11e1b2SBarry Smith #define MATSEQAIJCRL       "seqaijcrl"
415a11e1b2SBarry Smith #define MATMPIAIJCRL       "mpiaijcrl"
428154be41SBarry Smith #define MATAIJCUSP         "aijcusp"
438154be41SBarry Smith #define MATSEQAIJCUSP      "seqaijcusp"
448154be41SBarry Smith #define MATMPIAIJCUSP      "mpiaijcusp"
455a11e1b2SBarry Smith #define MATAIJPERM         "aijperm"
465a11e1b2SBarry Smith #define MATSEQAIJPERM      "seqaijperm"
475a11e1b2SBarry Smith #define MATMPIAIJPERM      "mpiaijperm"
48273d9f13SBarry Smith #define MATSHELL           "shell"
495a11e1b2SBarry Smith #define MATDENSE           "dense"
50209238afSKris Buschelman #define MATSEQDENSE        "seqdense"
51273d9f13SBarry Smith #define MATMPIDENSE        "mpidense"
525a11e1b2SBarry Smith #define MATBAIJ            "baij"
53273d9f13SBarry Smith #define MATSEQBAIJ         "seqbaij"
54273d9f13SBarry Smith #define MATMPIBAIJ         "mpibaij"
55273d9f13SBarry Smith #define MATMPIADJ          "mpiadj"
565a11e1b2SBarry Smith #define MATSBAIJ           "sbaij"
57273d9f13SBarry Smith #define MATSEQSBAIJ        "seqsbaij"
58273d9f13SBarry Smith #define MATMPISBAIJ        "mpisbaij"
59c0cdd4a1SDahai Guo #define MATSEQBSTRM        "seqbstrm"
60c0cdd4a1SDahai Guo #define MATMPIBSTRM        "mpibstrm"
61c0cdd4a1SDahai Guo #define MATBSTRM           "bstrm"
62aa5a9175SDahai Guo #define MATSEQSBSTRM       "seqsbstrm"
63aa5a9175SDahai Guo #define MATMPISBSTRM       "mpisbstrm"
64aa5a9175SDahai Guo #define MATSBSTRM          "sbstrm"
65cebc7f6cSBarry Smith #define MATDAAD            "daad"
66cebc7f6cSBarry Smith #define MATMFFD            "mffd"
67c8a8475eSBarry Smith #define MATNORMAL          "normal"
68ab92ecdeSBarry Smith #define MATLRC             "lrc"
692a6744ebSBarry Smith #define MATSCATTER         "scatter"
70421e10b8SBarry Smith #define MATBLOCKMAT        "blockmat"
71793850ffSBarry Smith #define MATCOMPOSITE       "composite"
721f8c7532SHong Zhang #define MATFFT             "fft"
731f8c7532SHong Zhang #define MATFFTW            "fftw"
74e133240eSMatthew G Knepley #define MATSEQCUFFT        "seqcufft"
75557cca28SSatish Balay #define MATTRANSPOSEMAT    "transpose"
7672ca8751SBarry Smith #define MATSCHURCOMPLEMENT "schurcomplement"
771d6018f0SLisandro Dalcin #define MATPYTHON          "python"
78f91d8e95SBarry Smith #define MATHYPRESTRUCT     "hyprestruct"
79a9e6138eSGlenn Hammond #define MATHYPRESSTRUCT    "hypresstruct"
80ee1cef2cSJed Brown #define MATSUBMATRIX       "submatrix"
812c0dbf93SJed Brown #define MATLOCALREF        "localref"
82d8588912SDave May #define MATNEST            "nest"
8315c1f1baSDmitry Karpeev #define MATIJ              "ij"
84773941d6SBarry Smith 
8576bdecfbSBarry Smith /*J
86c7393fdbSBarry Smith     MatSolverPackage - String with the name of a PETSc matrix solver type.
879989ab13SBarry Smith 
889989ab13SBarry Smith     For example: "petsc" indicates what PETSc provides, "superlu" indicates either
899989ab13SBarry Smith        SuperLU or SuperLU_Dist etc.
909989ab13SBarry Smith 
919989ab13SBarry Smith 
929989ab13SBarry Smith    Level: beginner
939989ab13SBarry Smith 
945c9eb25fSBarry Smith .seealso: MatGetFactor(), Mat, MatSetType(), MatType
9576bdecfbSBarry Smith J*/
96c7393fdbSBarry Smith #define MatSolverPackage char*
972692d6eeSBarry Smith #define MATSOLVERSPOOLES      "spooles"
982692d6eeSBarry Smith #define MATSOLVERSUPERLU      "superlu"
992692d6eeSBarry Smith #define MATSOLVERSUPERLU_DIST "superlu_dist"
1002692d6eeSBarry Smith #define MATSOLVERUMFPACK      "umfpack"
10120db9a53SJed Brown #define MATSOLVERCHOLMOD      "cholmod"
1022692d6eeSBarry Smith #define MATSOLVERESSL         "essl"
1032692d6eeSBarry Smith #define MATSOLVERLUSOL        "lusol"
1042692d6eeSBarry Smith #define MATSOLVERMUMPS        "mumps"
1052692d6eeSBarry Smith #define MATSOLVERPASTIX       "pastix"
1062692d6eeSBarry Smith #define MATSOLVERMATLAB       "matlab"
1072692d6eeSBarry Smith #define MATSOLVERPETSC        "petsc"
1082692d6eeSBarry Smith #define MATSOLVERPLAPACK      "plapack"
1092692d6eeSBarry Smith #define MATSOLVERBAS          "bas"
110773941d6SBarry Smith 
111aa5a9175SDahai Guo #define MATSOLVERBSTRM        "bstrm"
112aa5a9175SDahai Guo #define MATSOLVERSBSTRM       "sbstrm"
113c0cdd4a1SDahai Guo 
114b24902e0SBarry Smith /*E
115b24902e0SBarry Smith     MatFactorType - indicates what type of factorization is requested
116b24902e0SBarry Smith 
117b24902e0SBarry Smith     Level: beginner
118b24902e0SBarry Smith 
119b24902e0SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
120b24902e0SBarry Smith 
121c7393fdbSBarry Smith .seealso: MatSolverPackage, MatGetFactor()
122b24902e0SBarry Smith E*/
123599ef60dSHong Zhang typedef enum {MAT_FACTOR_NONE, MAT_FACTOR_LU, MAT_FACTOR_CHOLESKY, MAT_FACTOR_ILU, MAT_FACTOR_ICC,MAT_FACTOR_ILUDT} MatFactorType;
124*014dd563SJed Brown PETSC_EXTERN const char *const MatFactorTypes[];
125e92e720dSBarry Smith 
126*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactor(Mat,const MatSolverPackage,MatFactorType,Mat*);
127*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactorAvailable(Mat,const MatSolverPackage,MatFactorType,PetscBool *);
128*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFactorGetSolverPackage(Mat,const MatSolverPackage*);
129*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactorType(Mat,MatFactorType*);
1309989ab13SBarry Smith 
131c06d978dSMatthew Knepley /* Logging support */
1320700a824SBarry Smith #define    MAT_FILE_CLASSID 1211216    /* used to indicate matrices in binary files */
133*014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_CLASSID;
134*014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_FDCOLORING_CLASSID;
135*014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_TRANSPOSECOLORING_CLASSID;
136*014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_PARTITIONING_CLASSID;
137*014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_COARSEN_CLASSID;
138*014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_NULLSPACE_CLASSID;
139*014dd563SJed Brown PETSC_EXTERN PetscClassId MATMFFD_CLASSID;
140c06d978dSMatthew Knepley 
141ceb03754SKris Buschelman /*E
142ceb03754SKris Buschelman     MatReuse - Indicates if matrices obtained from a previous call to MatGetSubMatrices()
143d6eff37eSBarry Smith      or MatGetSubMatrix() are to be reused to store the new matrix values. For MatConvert() is used to indicate
144d6eff37eSBarry Smith      that the input matrix is to be replaced with the converted matrix.
145ceb03754SKris Buschelman 
146ceb03754SKris Buschelman     Level: beginner
147ceb03754SKris Buschelman 
148ceb03754SKris Buschelman    Any additions/changes here MUST also be made in include/finclude/petscmat.h
149ceb03754SKris Buschelman 
1500c451bc4SBarry Smith .seealso: MatGetSubMatrices(), MatGetSubMatrix(), MatDestroyMatrices(), MatConvert()
151ceb03754SKris Buschelman E*/
152dfe085dbSJed Brown typedef enum {MAT_INITIAL_MATRIX,MAT_REUSE_MATRIX,MAT_IGNORE_MATRIX} MatReuse;
153ceb03754SKris Buschelman 
1545494a064SHong Zhang /*E
1555494a064SHong Zhang     MatGetSubMatrixOption - Indicates if matrices obtained from a call to MatGetSubMatrices()
156829201f2SHong Zhang      include the matrix values. Currently it is only used by MatGetSeqNonzerostructure().
1575494a064SHong Zhang 
1585494a064SHong Zhang     Level: beginner
1595494a064SHong Zhang 
160829201f2SHong Zhang .seealso: MatGetSeqNonzerostructure()
1615494a064SHong Zhang E*/
1625494a064SHong Zhang typedef enum {MAT_DO_NOT_GET_VALUES,MAT_GET_VALUES} MatGetSubMatrixOption;
1635494a064SHong Zhang 
164*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInitializePackage(const char[]);
165c06d978dSMatthew Knepley 
166*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreate(MPI_Comm,Mat*);
167*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetSizes(Mat,PetscInt,PetscInt,PetscInt,PetscInt);
168*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetType(Mat,const MatType);
169*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetFromOptions(Mat);
170*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRegisterAll(const char[]);
171*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat));
172*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRegisterBaseName(const char[],const char[],const char[]);
173*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetOptionsPrefix(Mat,const char[]);
174*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAppendOptionsPrefix(Mat,const char[]);
175*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOptionsPrefix(Mat,const char*[]);
176f69a0ea3SMatthew Knepley 
17730de9b25SBarry Smith /*MC
17830de9b25SBarry Smith    MatRegisterDynamic - Adds a new matrix type
17930de9b25SBarry Smith 
18030de9b25SBarry Smith    Synopsis:
1811890ba74SBarry Smith    PetscErrorCode MatRegisterDynamic(const char *name,const char *path,const char *name_create,PetscErrorCode (*routine_create)(Mat))
18230de9b25SBarry Smith 
18330de9b25SBarry Smith    Not Collective
18430de9b25SBarry Smith 
18530de9b25SBarry Smith    Input Parameters:
18630de9b25SBarry Smith +  name - name of a new user-defined matrix type
18730de9b25SBarry Smith .  path - path (either absolute or relative) the library containing this solver
18830de9b25SBarry Smith .  name_create - name of routine to create method context
18930de9b25SBarry Smith -  routine_create - routine to create method context
19030de9b25SBarry Smith 
19130de9b25SBarry Smith    Notes:
19230de9b25SBarry Smith    MatRegisterDynamic() may be called multiple times to add several user-defined solvers.
19330de9b25SBarry Smith 
19430de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
19530de9b25SBarry Smith    is ignored.
19630de9b25SBarry Smith 
19730de9b25SBarry Smith    Sample usage:
19830de9b25SBarry Smith .vb
19930de9b25SBarry Smith    MatRegisterDynamic("my_mat",/home/username/my_lib/lib/libO/solaris/mylib.a,
20030de9b25SBarry Smith                "MyMatCreate",MyMatCreate);
20130de9b25SBarry Smith .ve
20230de9b25SBarry Smith 
20330de9b25SBarry Smith    Then, your solver can be chosen with the procedural interface via
20430de9b25SBarry Smith $     MatSetType(Mat,"my_mat")
20530de9b25SBarry Smith    or at runtime via the option
20630de9b25SBarry Smith $     -mat_type my_mat
20730de9b25SBarry Smith 
20830de9b25SBarry Smith    Level: advanced
20930de9b25SBarry Smith 
210ab901514SBarry Smith    Notes: ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values.
21130de9b25SBarry Smith          If your function is not being put into a shared library then use VecRegister() instead
21230de9b25SBarry Smith 
21330de9b25SBarry Smith .keywords: Mat, register
21430de9b25SBarry Smith 
21530de9b25SBarry Smith .seealso: MatRegisterAll(), MatRegisterDestroy()
21630de9b25SBarry Smith 
21730de9b25SBarry Smith M*/
218273d9f13SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
219273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,0)
220273d9f13SBarry Smith #else
221273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,d)
22230de9b25SBarry Smith #endif
22330de9b25SBarry Smith 
224*014dd563SJed Brown PETSC_EXTERN PetscBool MatRegisterAllCalled;
225*014dd563SJed Brown PETSC_EXTERN PetscFList MatList;
226*014dd563SJed Brown PETSC_EXTERN PetscFList MatColoringList;
227*014dd563SJed Brown PETSC_EXTERN PetscFList MatPartitioningList;
228*014dd563SJed Brown PETSC_EXTERN PetscFList MatCoarsenList;
22928988994SBarry Smith 
2303b224e63SBarry Smith /*E
2313b224e63SBarry Smith     MatStructure - Indicates if the matrix has the same nonzero structure
2323b224e63SBarry Smith 
2333b224e63SBarry Smith     Level: beginner
2343b224e63SBarry Smith 
2353b224e63SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
2363b224e63SBarry Smith 
2373b224e63SBarry Smith .seealso: MatCopy(), KSPSetOperators(), PCSetOperators()
2383b224e63SBarry Smith E*/
239214bc6a2SJed Brown typedef enum {DIFFERENT_NONZERO_PATTERN,SUBSET_NONZERO_PATTERN,SAME_NONZERO_PATTERN,SAME_PRECONDITIONER} MatStructure;
2403b224e63SBarry Smith 
241*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqDense(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*);
242*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateDense(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*);
243*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
244*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
245*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *);
246*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithSplitArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],PetscInt[],PetscInt[],PetscScalar[],Mat*);
2478d7a6e47SBarry Smith 
248*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
249*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
250*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat*);
251d21a29f3SJed Brown 
252*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAdj(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscInt[],Mat*);
253*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
254d21a29f3SJed Brown 
255*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
256*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPISBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *);
257*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
258*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatXAIJSetPreallocation(Mat,PetscInt,const PetscInt*,const PetscInt*,const PetscInt*,const PetscInt*);
259dfb205c3SBarry Smith 
260*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateShell(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,void *,Mat*);
261*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateAdic(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,void (*)(void),Mat*);
262*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateNormal(Mat,Mat*);
263*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateLRC(Mat,Mat,Mat,Mat*);
264*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateIS(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,ISLocalToGlobalMapping,Mat*);
265*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
266*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
267c0cdd4a1SDahai Guo 
268*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
269*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
270*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqSBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
271*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPISBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
272c0cdd4a1SDahai Guo 
273*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateScatter(MPI_Comm,VecScatter,Mat*);
274*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScatterSetVecScatter(Mat,VecScatter);
275*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScatterGetVecScatter(Mat,VecScatter*);
276*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateBlockMat(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,Mat*);
277*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeAddMat(Mat,Mat);
278*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeMerge(Mat);
279*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateComposite(MPI_Comm,PetscInt,const Mat*,Mat*);
2806d7c1e57SBarry Smith typedef enum {MAT_COMPOSITE_ADDITIVE,MAT_COMPOSITE_MULTIPLICATIVE} MatCompositeType;
281*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeSetType(Mat,MatCompositeType);
2826d7c1e57SBarry Smith 
283*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateFFT(MPI_Comm,PetscInt,const PetscInt[],const MatType,Mat*);
284*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqCUFFT(MPI_Comm,PetscInt,const PetscInt[],Mat*);
285dedccee8SHong Zhang 
286*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateTranspose(Mat,Mat*);
287*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSubMatrix(Mat,IS,IS,Mat*);
288*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSubMatrixUpdate(Mat,Mat,IS,IS);
289*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateLocalRef(Mat,IS,IS,Mat*);
2901472f72bSBarry Smith 
291*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPythonSetType(Mat,const char[]);
2921d6018f0SLisandro Dalcin 
293*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetUp(Mat);
294*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroy(Mat*);
29521c89e3eSBarry Smith 
296*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatConjugate(Mat);
297*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRealPart(Mat);
298*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatImaginaryPart(Mat);
299*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetDiagonalBlock(Mat,Mat*);
300*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetTrace(Mat,PetscScalar*);
301*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInvertBlockDiagonal(Mat,PetscScalar **);
30299cafbc1SBarry Smith 
3038ed539a5SBarry Smith /* ------------------------------------------------------------*/
304*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
305*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
306*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesRow(Mat,PetscInt,const PetscScalar[]);
307*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesRowLocal(Mat,PetscInt,const PetscScalar[]);
308*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBatch(Mat,PetscInt,PetscInt,PetscInt[],const PetscScalar[]);
30984cb2905SBarry Smith 
3102ef4de8bSBarry Smith /*S
3112ef4de8bSBarry Smith      MatStencil - Data structure (C struct) for storing information about a single row or
312be479b30SJed Brown         column of a matrix as indexed on an associated grid.
313be479b30SJed Brown 
314be479b30SJed Brown    Fortran usage is different, see MatSetValuesStencil() for details.
3152ef4de8bSBarry Smith 
3162ef4de8bSBarry Smith    Level: beginner
3172ef4de8bSBarry Smith 
3182ef4de8bSBarry Smith   Concepts: matrix; linear operator
3192ef4de8bSBarry Smith 
320be479b30SJed Brown .seealso:  MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockedStencil()
3212ef4de8bSBarry Smith S*/
322435da068SBarry Smith typedef struct {
323c1ac3661SBarry Smith   PetscInt k,j,i,c;
324435da068SBarry Smith } MatStencil;
3252ef4de8bSBarry Smith 
326*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
327*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode);
328*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt);
329435da068SBarry Smith 
330*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetColoring(Mat,ISColoring);
331*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesAdic(Mat,void*);
332*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesAdifor(Mat,PetscInt,void*);
3333a7fca6bSBarry Smith 
334d91e6319SBarry Smith /*E
335d91e6319SBarry Smith     MatAssemblyType - Indicates if the matrix is now to be used, or if you plan
336d91e6319SBarry Smith      to continue to add values to it
337d91e6319SBarry Smith 
338d91e6319SBarry Smith     Level: beginner
339d91e6319SBarry Smith 
340d91e6319SBarry Smith .seealso: MatAssemblyBegin(), MatAssemblyEnd()
341d91e6319SBarry Smith E*/
3426d4a8577SBarry Smith typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType;
343*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssemblyBegin(Mat,MatAssemblyType);
344*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssemblyEnd(Mat,MatAssemblyType);
345*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssembled(Mat,PetscBool *);
3464f9c727eSBarry Smith 
3471d73ed98SBarry Smith 
34830de9b25SBarry Smith 
349d91e6319SBarry Smith /*E
350d91e6319SBarry Smith     MatOption - Options that may be set for a matrix and its behavior or storage
351d91e6319SBarry Smith 
352d91e6319SBarry Smith     Level: beginner
353d91e6319SBarry Smith 
3540a835dfdSSatish Balay    Any additions/changes here MUST also be made in include/finclude/petscmat.h
355d91e6319SBarry Smith 
356d91e6319SBarry Smith .seealso: MatSetOption()
357d91e6319SBarry Smith E*/
358512a5fc5SBarry Smith typedef enum {MAT_ROW_ORIENTED,MAT_NEW_NONZERO_LOCATIONS,
3594e0d8c25SBarry Smith               MAT_SYMMETRIC,
3604e0d8c25SBarry Smith               MAT_STRUCTURALLY_SYMMETRIC,
3618047e77bSBarry Smith               MAT_NEW_DIAGONALS,MAT_IGNORE_OFF_PROC_ENTRIES,
3624e0d8c25SBarry Smith               MAT_NEW_NONZERO_LOCATION_ERR,
3634e0d8c25SBarry Smith               MAT_NEW_NONZERO_ALLOCATION_ERR,MAT_USE_HASH_TABLE,
364a9817697SBarry Smith               MAT_KEEP_NONZERO_PATTERN,MAT_IGNORE_ZERO_ENTRIES,
3654e0d8c25SBarry Smith               MAT_USE_INODES,
3664e0d8c25SBarry Smith               MAT_HERMITIAN,
3674e0d8c25SBarry Smith               MAT_SYMMETRY_ETERNAL,
368cd6b891eSBarry Smith               MAT_CHECK_COMPRESSED_ROW,
3694e0d8c25SBarry Smith               MAT_IGNORE_LOWER_TRIANGULAR,MAT_ERROR_LOWER_TRIANGULAR,
37028b2fa4aSMatthew Knepley               MAT_GETROW_UPPERTRIANGULAR,MAT_UNUSED_NONZERO_LOCATION_ERR,
3714cb17eb5SBarry Smith               MAT_SPD,MAT_NO_OFF_PROC_ENTRIES,MAT_NO_OFF_PROC_ZERO_ROWS,
37228b2fa4aSMatthew Knepley               NUM_MAT_OPTIONS} MatOption;
373*014dd563SJed Brown PETSC_EXTERN const char *MatOptions[];
374*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetOption(Mat,MatOption,PetscBool );
375*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetType(Mat,const MatType*);
37684cb2905SBarry Smith 
377*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]);
378*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
379*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
380*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowUpperTriangular(Mat);
381*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRowUpperTriangular(Mat);
382*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
383*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]);
384*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnVector(Mat,Vec,PetscInt);
385*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetArray(Mat,PetscScalar *[]);
386*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreArray(Mat,PetscScalar *[]);
387*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBlockSize(Mat,PetscInt *);
388*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSize(Mat,PetscInt);
389*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBlockSizes(Mat,PetscInt *,PetscInt *);
390*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSizes(Mat,PetscInt,PetscInt);
391*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNThreads(Mat,PetscInt);
392*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNThreads(Mat,PetscInt*);
3931620fd73SBarry Smith 
394*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMult(Mat,Vec,Vec);
395*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultDiagonalBlock(Mat,Vec,Vec);
396*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultAdd(Mat,Vec,Vec,Vec);
397*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTranspose(Mat,Vec,Vec);
398*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultHermitianTranspose(Mat,Vec,Vec);
399*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsTranspose(Mat,Mat,PetscReal,PetscBool *);
400*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitianTranspose(Mat,Mat,PetscReal,PetscBool *);
401*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeAdd(Mat,Vec,Vec,Vec);
402*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultHermitianTransposeAdd(Mat,Vec,Vec,Vec);
403*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultConstrained(Mat,Vec,Vec);
404*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeConstrained(Mat,Vec,Vec);
405*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatSolve(Mat,Mat,Mat);
406f5edf698SHong Zhang 
407d91e6319SBarry Smith /*E
408d91e6319SBarry Smith     MatDuplicateOption - Indicates if a duplicated sparse matrix should have
409d91e6319SBarry Smith   its numerical values copied over or just its nonzero structure.
410d91e6319SBarry Smith 
411d91e6319SBarry Smith     Level: beginner
412d91e6319SBarry Smith 
413d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
414d91e6319SBarry Smith 
41570dcbbb9SBarry Smith $   MAT_SHARE_NONZERO_PATTERN - the i and j arrays in the new matrix will be shared with the original matrix
41670dcbbb9SBarry Smith $                               this also triggers the MAT_DO_NOT_COPY_VALUES option. This is used when you
41770dcbbb9SBarry Smith $                               have several matrices with the same nonzero pattern.
41870dcbbb9SBarry Smith 
419d91e6319SBarry Smith .seealso: MatDuplicate()
420d91e6319SBarry Smith E*/
42170dcbbb9SBarry Smith typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES,MAT_SHARE_NONZERO_PATTERN} MatDuplicateOption;
4222e8a6d31SBarry Smith 
423*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatConvert(Mat,const MatType,MatReuse,Mat*);
424*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDuplicate(Mat,MatDuplicateOption,Mat*);
42594a9d846SBarry Smith 
426cb5b572fSBarry Smith 
427*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCopy(Mat,Mat,MatStructure);
428*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatView(Mat,PetscViewer);
429*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsSymmetric(Mat,PetscReal,PetscBool *);
430*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsStructurallySymmetric(Mat,PetscBool *);
431*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitian(Mat,PetscReal,PetscBool *);
432*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsSymmetricKnown(Mat,PetscBool *,PetscBool *);
433*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitianKnown(Mat,PetscBool *,PetscBool *);
434*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMissingDiagonal(Mat,PetscBool  *,PetscInt *);
435*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLoad(Mat, PetscViewer);
4367b80b807SBarry Smith 
437*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,PetscInt *[],PetscInt *[],PetscBool  *);
438*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,PetscInt *[],PetscInt *[],PetscBool  *);
439*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,PetscInt *[],PetscInt *[],PetscBool  *);
440*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,PetscInt *[],PetscInt *[],PetscBool  *);
441d4fbbf0eSBarry Smith 
442d91e6319SBarry Smith /*S
443d91e6319SBarry Smith      MatInfo - Context of matrix information, used with MatGetInfo()
444d91e6319SBarry Smith 
445d91e6319SBarry Smith    In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE
446d91e6319SBarry Smith 
447d91e6319SBarry Smith    Level: intermediate
448d91e6319SBarry Smith 
449d91e6319SBarry Smith   Concepts: matrix^nonzero information
450d91e6319SBarry Smith 
451d9274352SBarry Smith .seealso:  MatGetInfo(), MatInfoType
452d91e6319SBarry Smith S*/
4534e220ebcSLois Curfman McInnes typedef struct {
454b0a32e0cSBarry Smith   PetscLogDouble block_size;                         /* block size */
455b0a32e0cSBarry Smith   PetscLogDouble nz_allocated,nz_used,nz_unneeded;   /* number of nonzeros */
456b0a32e0cSBarry Smith   PetscLogDouble memory;                             /* memory allocated */
457b0a32e0cSBarry Smith   PetscLogDouble assemblies;                         /* number of matrix assemblies called */
458b0a32e0cSBarry Smith   PetscLogDouble mallocs;                            /* number of mallocs during MatSetValues() */
459b0a32e0cSBarry Smith   PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */
460b0a32e0cSBarry Smith   PetscLogDouble factor_mallocs;                     /* number of mallocs during factorization */
4614e220ebcSLois Curfman McInnes } MatInfo;
4624e220ebcSLois Curfman McInnes 
463d9274352SBarry Smith /*E
464d9274352SBarry Smith     MatInfoType - Indicates if you want information about the local part of the matrix,
465d9274352SBarry Smith      the entire parallel matrix or the maximum over all the local parts.
466d9274352SBarry Smith 
467d9274352SBarry Smith     Level: beginner
468d9274352SBarry Smith 
469d9274352SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
470d9274352SBarry Smith 
471d9274352SBarry Smith .seealso: MatGetInfo(), MatInfo
472d9274352SBarry Smith E*/
4737b80b807SBarry Smith typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType;
474*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetInfo(Mat,MatInfoType,MatInfo*);
475*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetDiagonal(Mat,Vec);
476*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMax(Mat,Vec,PetscInt[]);
477*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMin(Mat,Vec,PetscInt[]);
478*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMaxAbs(Mat,Vec,PetscInt[]);
479*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMinAbs(Mat,Vec,PetscInt[]);
480*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowSum(Mat,Vec);
481*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTranspose(Mat,MatReuse,Mat*);
482*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatHermitianTranspose(Mat,MatReuse,Mat*);
483*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPermute(Mat,IS,IS,Mat *);
484*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalScale(Mat,Vec,Vec);
485*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalSet(Mat,Vec,InsertMode);
486*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatEqual(Mat,Mat,PetscBool *);
487*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultEqual(Mat,Mat,PetscInt,PetscBool *);
488*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultAddEqual(Mat,Mat,PetscInt,PetscBool *);
489*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeEqual(Mat,Mat,PetscInt,PetscBool *);
490*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscBool *);
4917b80b807SBarry Smith 
492*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNorm(Mat,NormType,PetscReal *);
493*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnNorms(Mat,NormType,PetscReal *);
494*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroEntries(Mat);
495*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
496*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsIS(Mat,IS,PetscScalar,Vec,Vec);
497*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsStencil(Mat,PetscInt,const MatStencil [],PetscScalar,Vec,Vec);
498*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsStencil(Mat,PetscInt,const MatStencil[],PetscScalar,Vec,Vec);
499*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumns(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
500*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsIS(Mat,IS,PetscScalar,Vec,Vec);
5017b80b807SBarry Smith 
502*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatUseScaledForm(Mat,PetscBool );
503*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScaleSystem(Mat,Vec,Vec);
504*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatUnScaleSystem(Mat,Vec,Vec);
5055ef9f2a5SBarry Smith 
506*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSize(Mat,PetscInt*,PetscInt*);
507*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalSize(Mat,PetscInt*,PetscInt*);
508*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRange(Mat,PetscInt*,PetscInt*);
509*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRanges(Mat,const PetscInt**);
510*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRangeColumn(Mat,PetscInt*,PetscInt*);
511*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRangesColumn(Mat,const PetscInt**);
5127b80b807SBarry Smith 
513*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
514*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSubMatricesParallel(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
515*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroyMatrices(PetscInt,Mat *[]);
516*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSubMatrix(Mat,IS,IS,MatReuse,Mat *);
517*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalSubMatrix(Mat,IS,IS,Mat*);
518*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreLocalSubMatrix(Mat,IS,IS,Mat*);
519*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSeqNonzeroStructure(Mat,Mat*);
520*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroySeqNonzeroStructure(Mat*);
5215494a064SHong Zhang 
522*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJConcatenateSeqAIJ(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
523*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJConcatenateSeqAIJSymbolic(MPI_Comm,Mat,PetscInt,Mat*);
524*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJConcatenateSeqAIJNumeric(MPI_Comm,Mat,PetscInt,Mat);
525*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJ(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*);
526*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJSymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*);
527*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJNumeric(Mat,Mat);
528*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMat(Mat,MatReuse,Mat*);
529*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*);
530*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,Mat*);
53143eb5e2fSMatthew Knepley #if defined (PETSC_USE_CTABLE)
532*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetCommunicationStructs(Mat, Vec *, PetscTable *, VecScatter *);
53343eb5e2fSMatthew Knepley #else
534*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetCommunicationStructs(Mat, Vec *, PetscInt *[], VecScatter *);
53543eb5e2fSMatthew Knepley #endif
536*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetGhosts(Mat, PetscInt *,const PetscInt *[]);
5378efafbd8SBarry Smith 
538*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt);
5397b80b807SBarry Smith 
540*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*);
541*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMultSymbolic(Mat,Mat,PetscReal,Mat*);
542*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMultNumeric(Mat,Mat,Mat);
54322440eb1SKris Buschelman 
544*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*);
545*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAPSymbolic(Mat,Mat,PetscReal,Mat*);
546*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAPNumeric(Mat,Mat,Mat);
547*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARt(Mat,Mat,MatReuse,PetscReal,Mat*);
548*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARtSymbolic(Mat,Mat,PetscReal,Mat*);
549*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARtNumeric(Mat,Mat,Mat);
55022440eb1SKris Buschelman 
551*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeMatMult(Mat,Mat,MatReuse,PetscReal,Mat*);
552*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposetMatMultSymbolic(Mat,Mat,PetscReal,Mat*);
553*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposetMatMultNumeric(Mat,Mat,Mat);
554*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMult(Mat,Mat,MatReuse,PetscReal,Mat*);
555*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMultSymbolic(Mat,Mat,PetscReal,Mat*);
556*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMultNumeric(Mat,Mat,Mat);
557bc011b1eSHong Zhang 
558*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAXPY(Mat,PetscScalar,Mat,MatStructure);
559*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAYPX(Mat,PetscScalar,Mat,MatStructure);
5601c741599SHong Zhang 
561*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScale(Mat,PetscScalar);
562*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShift(Mat,PetscScalar);
5637b80b807SBarry Smith 
564*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping);
565*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetLocalToGlobalMappingBlock(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping);
566*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*);
567*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalToGlobalMappingBlock(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*);
568*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
569*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsLocalIS(Mat,IS,PetscScalar,Vec,Vec);
570*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec);
571*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocalIS(Mat,IS,PetscScalar,Vec,Vec);
572*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
573*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
574052efed2SBarry Smith 
575*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStashSetInitialSize(Mat,PetscInt,PetscInt);
576*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
57790f02eecSBarry Smith 
578*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInterpolate(Mat,Vec,Vec);
579*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInterpolateAdd(Mat,Vec,Vec,Vec);
580*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestrict(Mat,Vec,Vec);
581*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetVecs(Mat,Vec*,Vec*);
582*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRedundantMatrix(Mat,PetscInt,MPI_Comm,PetscInt,MatReuse,Mat*);
583*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetMultiProcBlock(Mat,MPI_Comm,MatReuse,Mat*);
584*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFindZeroDiagonals(Mat,IS*);
585bd481603SBarry Smith 
586bd481603SBarry Smith /*MC
5871620fd73SBarry Smith    MatSetValue - Set a single entry into a matrix.
5881620fd73SBarry Smith 
5891620fd73SBarry Smith    Not collective
5901620fd73SBarry Smith 
5911620fd73SBarry Smith    Input Parameters:
5921620fd73SBarry Smith +  m - the matrix
5931620fd73SBarry Smith .  row - the row location of the entry
5941620fd73SBarry Smith .  col - the column location of the entry
5951620fd73SBarry Smith .  value - the value to insert
5961620fd73SBarry Smith -  mode - either INSERT_VALUES or ADD_VALUES
5971620fd73SBarry Smith 
5981620fd73SBarry Smith    Notes:
5991620fd73SBarry Smith    For efficiency one should use MatSetValues() and set several or many
6001620fd73SBarry Smith    values simultaneously if possible.
6011620fd73SBarry Smith 
6021620fd73SBarry Smith    Level: beginner
6031620fd73SBarry Smith 
6041620fd73SBarry Smith .seealso: MatSetValues(), MatSetValueLocal()
6051620fd73SBarry Smith M*/
6061620fd73SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatSetValue(Mat v,PetscInt i,PetscInt j,PetscScalar va,InsertMode mode) {return MatSetValues(v,1,&i,1,&j,&va,mode);}
6071620fd73SBarry Smith 
6081620fd73SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatGetValue(Mat v,PetscInt i,PetscInt j,PetscScalar *va) {return MatGetValues(v,1,&i,1,&j,va);}
6091620fd73SBarry Smith 
6101620fd73SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatSetValueLocal(Mat v,PetscInt i,PetscInt j,PetscScalar va,InsertMode mode) {return MatSetValuesLocal(v,1,&i,1,&j,&va,mode);}
6111620fd73SBarry Smith 
6121620fd73SBarry Smith /*MC
6130d2a7ff2SSatish Balay    MatPreallocateInitialize - Begins the block of code that will count the number of nonzeros per
614bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
615bd481603SBarry Smith 
616bd481603SBarry Smith    Synopsis:
617c1ac3661SBarry Smith    PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz)
618bd481603SBarry Smith 
619bd481603SBarry Smith    Collective on MPI_Comm
620bd481603SBarry Smith 
621bd481603SBarry Smith    Input Parameters:
622bd481603SBarry Smith +  comm - the communicator that will share the eventually allocated matrix
623859fcb39SBarry Smith .  nrows - the number of LOCAL rows in the matrix
624859fcb39SBarry Smith -  ncols - the number of LOCAL columns in the matrix
625bd481603SBarry Smith 
626bd481603SBarry Smith    Output Parameters:
627bd481603SBarry Smith +  dnz - the array that will be passed to the matrix preallocation routines
628bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
629bd481603SBarry Smith 
630bd481603SBarry Smith 
631bd481603SBarry Smith    Level: intermediate
632bd481603SBarry Smith 
633bd481603SBarry Smith    Notes:
6340598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
635bd481603SBarry Smith 
6361d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
637bd481603SBarry Smith 
638bd481603SBarry Smith    Use MatPreallocateInitializeSymmetric() for symmetric matrices (MPISBAIJ matrices)
639bd481603SBarry Smith 
6401620fd73SBarry Smith    This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize().
6411620fd73SBarry Smith 
642bd481603SBarry Smith   Concepts: preallocation^Matrix
643bd481603SBarry Smith 
644bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
645bd481603SBarry Smith           MatPreallocateInitializeSymmetric(), MatPreallocateSymmetricSetLocal()
646bd481603SBarry Smith M*/
647c4f061fbSSatish Balay #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \
6487c922b88SBarry Smith { \
649a89bc211SBarry Smith   PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__start,__end; \
650a89bc211SBarry Smith   _4_ierr = PetscMalloc2(__nrows,PetscInt,&dnz,__nrows,PetscInt,&onz);CHKERRQ(_4_ierr); \
651a89bc211SBarry Smith   _4_ierr = PetscMemzero(dnz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\
652eabe889fSLisandro Dalcin   _4_ierr = PetscMemzero(onz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr); __start = 0; __end = __start; \
653c1ac3661SBarry Smith   _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ctmp;\
654a89bc211SBarry Smith   _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows;
6557c922b88SBarry Smith 
656bd481603SBarry Smith /*MC
657bd481603SBarry Smith    MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
658bd481603SBarry Smith        inserted using a local number of the rows and columns
659bd481603SBarry Smith 
660bd481603SBarry Smith    Synopsis:
661c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
662bd481603SBarry Smith 
663bd481603SBarry Smith    Not Collective
664bd481603SBarry Smith 
665bd481603SBarry Smith    Input Parameters:
666784ac674SJed Brown +  map - the row mapping from local numbering to global numbering
667bd481603SBarry Smith .  nrows - the number of rows indicated
6681d73ed98SBarry Smith .  rows - the indices of the rows
669784ac674SJed Brown .  cmap - the column mapping from local to global numbering
670bd481603SBarry Smith .  ncols - the number of columns in the matrix
671bd481603SBarry Smith .  cols - the columns indicated
672bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
673bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
674bd481603SBarry Smith 
675bd481603SBarry Smith 
676bd481603SBarry Smith    Level: intermediate
677bd481603SBarry Smith 
678bd481603SBarry Smith    Notes:
6790598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
680bd481603SBarry Smith 
6811d73ed98SBarry Smith    Do not malloc or free dnz and onz, that is handled internally by these routines
682bd481603SBarry Smith 
683bd481603SBarry Smith   Concepts: preallocation^Matrix
684bd481603SBarry Smith 
685bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
686bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal()
687bd481603SBarry Smith M*/
688784ac674SJed Brown #define MatPreallocateSetLocal(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \
689c4f061fbSSatish Balay {\
690c1ac3661SBarry Smith   PetscInt __l;\
691784ac674SJed Brown   _4_ierr = ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\
692784ac674SJed Brown   _4_ierr = ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\
693c4f061fbSSatish Balay   for (__l=0;__l<nrows;__l++) {\
694ef66eb69SBarry Smith     _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
695c4f061fbSSatish Balay   }\
696c4f061fbSSatish Balay }
697c4f061fbSSatish Balay 
698bd481603SBarry Smith /*MC
699bd481603SBarry Smith    MatPreallocateSymmetricSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be
700bd481603SBarry Smith        inserted using a local number of the rows and columns
701bd481603SBarry Smith 
702bd481603SBarry Smith    Synopsis:
703c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
704bd481603SBarry Smith 
705bd481603SBarry Smith    Not Collective
706bd481603SBarry Smith 
707bd481603SBarry Smith    Input Parameters:
708bd481603SBarry Smith +  map - the mapping between local numbering and global numbering
709bd481603SBarry Smith .  nrows - the number of rows indicated
7101d73ed98SBarry Smith .  rows - the indices of the rows
711bd481603SBarry Smith .  ncols - the number of columns in the matrix
712bd481603SBarry Smith .  cols - the columns indicated
713bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
714bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
715bd481603SBarry Smith 
716bd481603SBarry Smith 
717bd481603SBarry Smith    Level: intermediate
718bd481603SBarry Smith 
719bd481603SBarry Smith    Notes:
7200598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
721bd481603SBarry Smith 
722bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
723bd481603SBarry Smith 
724bd481603SBarry Smith   Concepts: preallocation^Matrix
725bd481603SBarry Smith 
726bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
727bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
728bd481603SBarry Smith M*/
729d3d32019SBarry Smith #define MatPreallocateSymmetricSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\
730d3d32019SBarry Smith {\
731c1ac3661SBarry Smith   PetscInt __l;\
732d3d32019SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\
733d3d32019SBarry Smith   _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\
734d3d32019SBarry Smith   for (__l=0;__l<nrows;__l++) {\
735d3d32019SBarry Smith     _4_ierr = MatPreallocateSymmetricSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\
736d3d32019SBarry Smith   }\
737d3d32019SBarry Smith }
738d3d32019SBarry Smith 
739bd481603SBarry Smith /*MC
740bd481603SBarry Smith    MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
741bd481603SBarry Smith        inserted using a local number of the rows and columns
742bd481603SBarry Smith 
743bd481603SBarry Smith    Synopsis:
744c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
745bd481603SBarry Smith 
746bd481603SBarry Smith    Not Collective
747bd481603SBarry Smith 
748bd481603SBarry Smith    Input Parameters:
74964075487SBarry Smith +  row - the row
750bd481603SBarry Smith .  ncols - the number of columns in the matrix
751a50f8bf6SHong Zhang -  cols - the columns indicated
752a50f8bf6SHong Zhang 
753a50f8bf6SHong Zhang    Output Parameters:
754a50f8bf6SHong Zhang +  dnz - the array that will be passed to the matrix preallocation routines
755bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
756bd481603SBarry Smith 
757bd481603SBarry Smith 
758bd481603SBarry Smith    Level: intermediate
759bd481603SBarry Smith 
760bd481603SBarry Smith    Notes:
7610598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
762bd481603SBarry Smith 
763bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
764bd481603SBarry Smith 
7651620fd73SBarry Smith    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
7661620fd73SBarry Smith 
767bd481603SBarry Smith   Concepts: preallocation^Matrix
768bd481603SBarry Smith 
769bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
770bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
771bd481603SBarry Smith M*/
772c4f061fbSSatish Balay #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\
773c1ac3661SBarry Smith { PetscInt __i; \
774e32f2f54SBarry Smith   if (row < __rstart) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Trying to set preallocation for row %D less than first local row %D",row,__rstart);\
775e32f2f54SBarry Smith   if (row >= __rstart+__nrows) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Trying to set preallocation for row %D greater than last local row %D",row,__rstart+__nrows-1);\
7767c922b88SBarry Smith   for (__i=0; __i<nc; __i++) {\
77764075487SBarry Smith     if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \
7787cc688d7SBarry Smith     else dnz[row - __rstart]++;\
7797c922b88SBarry Smith   }\
7807c922b88SBarry Smith }
7817c922b88SBarry Smith 
782bd481603SBarry Smith /*MC
783bd481603SBarry Smith    MatPreallocateSymmetricSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be
784bd481603SBarry Smith        inserted using a local number of the rows and columns
785bd481603SBarry Smith 
786bd481603SBarry Smith    Synopsis:
787c1ac3661SBarry Smith    PetscErrorCode MatPreallocateSymmetricSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz)
788bd481603SBarry Smith 
789bd481603SBarry Smith    Not Collective
790bd481603SBarry Smith 
791bd481603SBarry Smith    Input Parameters:
792bd481603SBarry Smith +  nrows - the number of rows indicated
7931d73ed98SBarry Smith .  rows - the indices of the rows
794bd481603SBarry Smith .  ncols - the number of columns in the matrix
795bd481603SBarry Smith .  cols - the columns indicated
796bd481603SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
797bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
798bd481603SBarry Smith 
799bd481603SBarry Smith 
800bd481603SBarry Smith    Level: intermediate
801bd481603SBarry Smith 
802bd481603SBarry Smith    Notes:
8030598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
804bd481603SBarry Smith 
805bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
806bd481603SBarry Smith 
8071620fd73SBarry Smith    This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize().
8081620fd73SBarry Smith 
809bd481603SBarry Smith   Concepts: preallocation^Matrix
810bd481603SBarry Smith 
811bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(),
812bd481603SBarry Smith           MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal()
813bd481603SBarry Smith M*/
814d3d32019SBarry Smith #define MatPreallocateSymmetricSet(row,nc,cols,dnz,onz) 0;\
815c1ac3661SBarry Smith { PetscInt __i; \
816d3d32019SBarry Smith   for (__i=0; __i<nc; __i++) {\
817d3d32019SBarry Smith     if (cols[__i] >= __end) onz[row - __rstart]++; \
818d3d32019SBarry Smith     else if (cols[__i] >= row) dnz[row - __rstart]++;\
819d3d32019SBarry Smith   }\
820d3d32019SBarry Smith }
821d3d32019SBarry Smith 
822bd481603SBarry Smith /*MC
82316371a99SBarry Smith    MatPreallocateLocation -  An alternative to MatPreallocationSet() that puts the nonzero locations into the matrix if it exists
82416371a99SBarry Smith 
82516371a99SBarry Smith    Synopsis:
82616371a99SBarry Smith    PetscErrorCode MatPreallocateLocations(Mat A,PetscInt row,PetscInt ncols,PetscInt *cols,PetscInt *dnz,PetscInt *onz)
82716371a99SBarry Smith 
82816371a99SBarry Smith    Not Collective
82916371a99SBarry Smith 
83016371a99SBarry Smith    Input Parameters:
83116371a99SBarry Smith .  A - matrix
83216371a99SBarry Smith .  row - row where values exist (must be local to this process)
83316371a99SBarry Smith .  ncols - number of columns
83416371a99SBarry Smith .  cols - columns with nonzeros
83516371a99SBarry Smith .  dnz - the array that will be passed to the matrix preallocation routines
83616371a99SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
83716371a99SBarry Smith 
83816371a99SBarry Smith 
83916371a99SBarry Smith    Level: intermediate
84016371a99SBarry Smith 
84116371a99SBarry Smith    Notes:
8420598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
84316371a99SBarry Smith 
84416371a99SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
84516371a99SBarry Smith 
84616371a99SBarry Smith    This is a MACRO not a function because it uses a bunch of variables private to the MatPreallocation.... routines.
84716371a99SBarry Smith 
84816371a99SBarry Smith   Concepts: preallocation^Matrix
84916371a99SBarry Smith 
85016371a99SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
851eabe889fSLisandro Dalcin           MatPreallocateSymmetricSetLocal()
85216371a99SBarry Smith M*/
85316371a99SBarry 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);}
85416371a99SBarry Smith 
85516371a99SBarry Smith 
85616371a99SBarry Smith /*MC
8570d2a7ff2SSatish Balay    MatPreallocateFinalize - Ends the block of code that will count the number of nonzeros per
858bd481603SBarry Smith        row in a matrix providing the data that one can use to correctly preallocate the matrix.
859bd481603SBarry Smith 
860bd481603SBarry Smith    Synopsis:
861c1ac3661SBarry Smith    PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz)
862bd481603SBarry Smith 
863bd481603SBarry Smith    Collective on MPI_Comm
864bd481603SBarry Smith 
865bd481603SBarry Smith    Input Parameters:
86616371a99SBarry Smith +  dnz - the array that was be passed to the matrix preallocation routines
867bd481603SBarry Smith -  ozn - the other array passed to the matrix preallocation routines
868bd481603SBarry Smith 
869bd481603SBarry Smith 
870bd481603SBarry Smith    Level: intermediate
871bd481603SBarry Smith 
872bd481603SBarry Smith    Notes:
8730598bfebSBarry Smith     See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details.
874bd481603SBarry Smith 
875bd481603SBarry Smith    Do not malloc or free dnz and onz that is handled internally by these routines
876bd481603SBarry Smith 
8771620fd73SBarry Smith    This is a MACRO not a function because it closes the { started in MatPreallocateInitialize().
8781620fd73SBarry Smith 
879bd481603SBarry Smith   Concepts: preallocation^Matrix
880bd481603SBarry Smith 
881bd481603SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(),
882eabe889fSLisandro Dalcin           MatPreallocateSymmetricSetLocal()
883bd481603SBarry Smith M*/
884a89bc211SBarry Smith #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree2(dnz,onz);CHKERRQ(_4_ierr);}
8857c922b88SBarry Smith 
886bd481603SBarry Smith 
887bd481603SBarry Smith 
8887b80b807SBarry Smith /* Routines unique to particular data structures */
889*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellGetContext(Mat,void *);
890698d4c6aSKris Buschelman 
891*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInodeAdjustForInodes(Mat,IS*,IS*);
892*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *);
893698d4c6aSKris Buschelman 
894*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetColumnIndices(Mat,PetscInt[]);
895*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetColumnIndices(Mat,PetscInt[]);
896*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
897*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
898*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*);
899*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJFromTriple(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*,PetscInt,PetscBool);
9007b80b807SBarry Smith 
901a96a251dSBarry Smith #define MAT_SKIP_ALLOCATION -4
902a96a251dSBarry Smith 
903*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
904*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]);
905*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]);
906273d9f13SBarry Smith 
907*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
908*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
909*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
910*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []);
911*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
912*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]);
913*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]);
914*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]);
915*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIDenseSetPreallocation(Mat,PetscScalar[]);
916*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqDenseSetPreallocation(Mat,PetscScalar[]);
917*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,PetscInt*[]);
918*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,PetscInt*[]);
919*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAdicSetLocalFunction(Mat,void (*)(void));
920*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAdjCreateNonemptySubcommMat(Mat,Mat*);
921273d9f13SBarry Smith 
922*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqDenseSetLDA(Mat,PetscInt);
923*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDenseGetLocalMatrix(Mat,Mat*);
9241b807ce4Svictorle 
925*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStoreValues(Mat);
926*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRetrieveValues(Mat);
9272e8a6d31SBarry Smith 
928*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDAADSetCtx(Mat,void*);
9293a7fca6bSBarry Smith 
930*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFindNonzeroRows(Mat,IS*);
9317b80b807SBarry Smith /*
9327b80b807SBarry Smith   These routines are not usually accessed directly, rather solving is
93394b7f48cSBarry Smith   done through the KSP and PC interfaces.
9347b80b807SBarry Smith */
9357b80b807SBarry Smith 
93676bdecfbSBarry Smith /*J
937d9274352SBarry Smith     MatOrderingType - String with the name of a PETSc matrix ordering or the creation function
938d9274352SBarry Smith        with an optional dynamic library name, for example
939d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:orderingcreate()
940d9274352SBarry Smith 
941d9274352SBarry Smith    Level: beginner
942d9274352SBarry Smith 
943e5a9bf91SBarry Smith    Cannot use const because the PC objects manipulate the string
944e5a9bf91SBarry Smith 
945d9274352SBarry Smith .seealso: MatGetOrdering()
94676bdecfbSBarry Smith J*/
9473eaad2caSSatish Balay #define MatOrderingType char*
9482692d6eeSBarry Smith #define MATORDERINGNATURAL     "natural"
9492692d6eeSBarry Smith #define MATORDERINGND          "nd"
9502692d6eeSBarry Smith #define MATORDERING1WD         "1wd"
9512692d6eeSBarry Smith #define MATORDERINGRCM         "rcm"
9522692d6eeSBarry Smith #define MATORDERINGQMD         "qmd"
9532692d6eeSBarry Smith #define MATORDERINGROWLENGTH   "rowlength"
954312e372aSBarry Smith #define MATORDERINGAMD         "amd"            /* only works if UMFPACK is installed with PETSc */
955b12f92e5SBarry Smith 
956*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOrdering(Mat,const MatOrderingType,IS*,IS*);
957*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOrderingList(PetscFList *list);
958*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatOrderingRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,const MatOrderingType,IS*,IS*));
95930de9b25SBarry Smith 
96030de9b25SBarry Smith /*MC
9611890ba74SBarry Smith    MatOrderingRegisterDynamic - Adds a new sparse matrix ordering to the matrix package.
96230de9b25SBarry Smith 
96330de9b25SBarry Smith    Synopsis:
9641890ba74SBarry Smith    PetscErrorCode MatOrderingRegisterDynamic(const char *name_ordering,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatOrdering))
96530de9b25SBarry Smith 
96630de9b25SBarry Smith    Not Collective
96730de9b25SBarry Smith 
96830de9b25SBarry Smith    Input Parameters:
9692692d6eeSBarry Smith +  sname - name of ordering (for example MATORDERINGND)
97030de9b25SBarry Smith .  path - location of library where creation routine is
97130de9b25SBarry Smith .  name - name of function that creates the ordering type,a string
97230de9b25SBarry Smith -  function - function pointer that creates the ordering
97330de9b25SBarry Smith 
97430de9b25SBarry Smith    Level: developer
97530de9b25SBarry Smith 
97630de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
97730de9b25SBarry Smith    is ignored.
97830de9b25SBarry Smith 
97930de9b25SBarry Smith    Sample usage:
98030de9b25SBarry Smith .vb
98130de9b25SBarry Smith    MatOrderingRegisterDynamic("my_order",/home/username/my_lib/lib/libO/solaris/mylib.a,
98230de9b25SBarry Smith                "MyOrder",MyOrder);
98330de9b25SBarry Smith .ve
98430de9b25SBarry Smith 
98530de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
98630de9b25SBarry Smith $     MatOrderingSetType(part,"my_order)
98730de9b25SBarry Smith    or at runtime via the option
9882401956bSBarry Smith $     -pc_factor_mat_ordering_type my_order
98930de9b25SBarry Smith 
990ab901514SBarry Smith    ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values.
99130de9b25SBarry Smith 
99230de9b25SBarry Smith .keywords: matrix, ordering, register
99330de9b25SBarry Smith 
99430de9b25SBarry Smith .seealso: MatOrderingRegisterDestroy(), MatOrderingRegisterAll()
99530de9b25SBarry Smith M*/
996aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
997f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,0)
998b12f92e5SBarry Smith #else
999f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,d)
1000b12f92e5SBarry Smith #endif
100130de9b25SBarry Smith 
1002*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatOrderingRegisterDestroy(void);
1003*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatOrderingRegisterAll(const char[]);
1004*014dd563SJed Brown PETSC_EXTERN PetscBool MatOrderingRegisterAllCalled;
1005*014dd563SJed Brown PETSC_EXTERN PetscFList MatOrderingList;
1006d4fbbf0eSBarry Smith 
1007*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS);
1008a2ce50c7SBarry Smith 
1009d91e6319SBarry Smith /*S
1010d90ac83dSHong Zhang     MatFactorShiftType - Numeric Shift.
1011d90ac83dSHong Zhang 
1012d90ac83dSHong Zhang    Level: beginner
1013d90ac83dSHong Zhang 
1014d90ac83dSHong Zhang S*/
1015d90ac83dSHong Zhang typedef enum {MAT_SHIFT_NONE,MAT_SHIFT_NONZERO,MAT_SHIFT_POSITIVE_DEFINITE,MAT_SHIFT_INBLOCKS} MatFactorShiftType;
1016*014dd563SJed Brown PETSC_EXTERN const char *MatFactorShiftTypes[];
1017d90ac83dSHong Zhang 
1018d90ac83dSHong Zhang /*S
10192401956bSBarry Smith    MatFactorInfo - Data passed into the matrix factorization routines
1020b00f7748SHong Zhang 
102161cad860SBarry Smith    In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE, that is use
102261cad860SBarry Smith $     MatFactorInfo  info(MAT_FACTORINFO_SIZE)
1023b00f7748SHong Zhang 
102415e8a5b3SHong Zhang    Notes: These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC.
1025b00f7748SHong Zhang 
102661cad860SBarry Smith       You can use MatFactorInfoInitialize() to set default values.
102761cad860SBarry Smith 
1028b00f7748SHong Zhang    Level: developer
1029b00f7748SHong Zhang 
1030d7d82daaSBarry Smith .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(),
1031d7d82daaSBarry Smith           MatFactorInfoInitialize()
1032b00f7748SHong Zhang 
1033b00f7748SHong Zhang S*/
1034b00f7748SHong Zhang typedef struct {
103515e8a5b3SHong Zhang   PetscReal     diagonal_fill;  /* force diagonal to fill in if initially not filled */
103685317021SBarry Smith   PetscReal     usedt;
103715e8a5b3SHong Zhang   PetscReal     dt;             /* drop tolerance */
1038b00f7748SHong Zhang   PetscReal     dtcol;          /* tolerance for pivoting */
103915e8a5b3SHong Zhang   PetscReal     dtcount;        /* maximum nonzeros to be allowed per row */
104067e28bfeSBarry Smith   PetscReal     fill;           /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */
1041348344bbSBarry Smith   PetscReal     levels;         /* ICC/ILU(levels) */
1042bcd9e38bSBarry Smith   PetscReal     pivotinblocks;  /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0
1043bcd9e38bSBarry Smith                                    factorization may be faster if do not pivot */
104415e8a5b3SHong Zhang   PetscReal     zeropivot;      /* pivot is called zero if less than this */
1045f4db908eSBarry Smith   PetscReal     shifttype;      /* type of shift added to matrix factor to prevent zero pivots */
1046f4db908eSBarry Smith   PetscReal     shiftamount;     /* how large the shift is */
104715e8a5b3SHong Zhang } MatFactorInfo;
1048ffa6d0a5SLois Curfman McInnes 
1049*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFactorInfoInitialize(MatFactorInfo*);
1050*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCholeskyFactor(Mat,IS,const MatFactorInfo*);
1051*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCholeskyFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
1052*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCholeskyFactorNumeric(Mat,Mat,const MatFactorInfo*);
1053*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLUFactor(Mat,IS,IS,const MatFactorInfo*);
1054*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatILUFactor(Mat,IS,IS,const MatFactorInfo*);
1055*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
1056*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatILUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*);
1057*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatICCFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*);
1058*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatICCFactor(Mat,IS,const MatFactorInfo*);
1059*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLUFactorNumeric(Mat,Mat,const MatFactorInfo*);
1060*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*);
1061*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolve(Mat,Vec,Vec);
1062*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatForwardSolve(Mat,Vec,Vec);
1063*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatBackwardSolve(Mat,Vec,Vec);
1064*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveAdd(Mat,Vec,Vec,Vec);
1065*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveTranspose(Mat,Vec,Vec);
1066*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveTransposeAdd(Mat,Vec,Vec,Vec);
1067*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolves(Mat,Vecs,Vecs);
10688ed539a5SBarry Smith 
1069*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetUnfactored(Mat);
1070bb5a7306SBarry Smith 
1071d91e6319SBarry Smith /*E
1072d91e6319SBarry Smith     MatSORType - What type of (S)SOR to perform
1073bb1eb677SSatish Balay 
1074d91e6319SBarry Smith     Level: beginner
1075d91e6319SBarry Smith 
1076d9274352SBarry Smith    May be bitwise ORd together
1077d9274352SBarry Smith 
1078d91e6319SBarry Smith    Any additions/changes here MUST also be made in include/finclude/petscmat.h
1079d91e6319SBarry Smith 
10804e7234bfSBarry Smith    MatSORType may be bitwise ORd together, so do not change the numbers
10814e7234bfSBarry Smith 
108241f059aeSBarry Smith .seealso: MatSOR()
1083d91e6319SBarry Smith E*/
1084ee50ffe9SBarry Smith typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3,
1085ee50ffe9SBarry Smith               SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8,
1086ee50ffe9SBarry Smith               SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16,
108784cb2905SBarry Smith               SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType;
1088*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSOR(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
10898ed539a5SBarry Smith 
1090d4fbbf0eSBarry Smith /*
1091639f9d9dSBarry Smith     These routines are for efficiently computing Jacobians via finite differences.
1092639f9d9dSBarry Smith */
1093b12f92e5SBarry Smith 
109476bdecfbSBarry Smith /*J
1095d9274352SBarry Smith     MatColoringType - String with the name of a PETSc matrix coloring or the creation function
1096d9274352SBarry Smith        with an optional dynamic library name, for example
1097d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:coloringcreate()
1098d9274352SBarry Smith 
1099d9274352SBarry Smith    Level: beginner
1100d9274352SBarry Smith 
1101d9274352SBarry Smith .seealso: MatGetColoring()
110276bdecfbSBarry Smith J*/
1103a313700dSBarry Smith #define MatColoringType char*
11042692d6eeSBarry Smith #define MATCOLORINGNATURAL "natural"
11052692d6eeSBarry Smith #define MATCOLORINGSL      "sl"
11062692d6eeSBarry Smith #define MATCOLORINGLF      "lf"
11072692d6eeSBarry Smith #define MATCOLORINGID      "id"
1108b12f92e5SBarry Smith 
1109*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColoring(Mat,const MatColoringType,ISColoring*);
1110*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatColoringRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,MatColoringType,ISColoring *));
111130de9b25SBarry Smith 
111230de9b25SBarry Smith /*MC
111330de9b25SBarry Smith    MatColoringRegisterDynamic - Adds a new sparse matrix coloring to the
111430de9b25SBarry Smith                                matrix package.
111530de9b25SBarry Smith 
111630de9b25SBarry Smith    Synopsis:
11171890ba74SBarry Smith    PetscErrorCode MatColoringRegisterDynamic(const char *name_coloring,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatColoring))
111830de9b25SBarry Smith 
111930de9b25SBarry Smith    Not Collective
112030de9b25SBarry Smith 
112130de9b25SBarry Smith    Input Parameters:
11222692d6eeSBarry Smith +  sname - name of Coloring (for example MATCOLORINGSL)
112330de9b25SBarry Smith .  path - location of library where creation routine is
112430de9b25SBarry Smith .  name - name of function that creates the Coloring type, a string
112530de9b25SBarry Smith -  function - function pointer that creates the coloring
112630de9b25SBarry Smith 
112730de9b25SBarry Smith    Level: developer
112830de9b25SBarry Smith 
112930de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
113030de9b25SBarry Smith    is ignored.
113130de9b25SBarry Smith 
113230de9b25SBarry Smith    Sample usage:
113330de9b25SBarry Smith .vb
113430de9b25SBarry Smith    MatColoringRegisterDynamic("my_color",/home/username/my_lib/lib/libO/solaris/mylib.a,
113530de9b25SBarry Smith                "MyColor",MyColor);
113630de9b25SBarry Smith .ve
113730de9b25SBarry Smith 
113830de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
113930de9b25SBarry Smith $     MatColoringSetType(part,"my_color")
114030de9b25SBarry Smith    or at runtime via the option
114130de9b25SBarry Smith $     -mat_coloring_type my_color
114230de9b25SBarry Smith 
1143ab901514SBarry Smith    $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
114430de9b25SBarry Smith 
114530de9b25SBarry Smith .keywords: matrix, Coloring, register
114630de9b25SBarry Smith 
114730de9b25SBarry Smith .seealso: MatColoringRegisterDestroy(), MatColoringRegisterAll()
114830de9b25SBarry Smith M*/
1149aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1150f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,0)
1151b12f92e5SBarry Smith #else
1152f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,d)
1153b12f92e5SBarry Smith #endif
115430de9b25SBarry Smith 
1155*014dd563SJed Brown PETSC_EXTERN PetscBool MatColoringRegisterAllCalled;
1156f1144a30SSatish Balay 
1157*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatColoringRegisterAll(const char[]);
1158*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatColoringRegisterDestroy(void);
1159*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*);
1160639f9d9dSBarry Smith 
1161d9274352SBarry Smith /*S
1162d9274352SBarry Smith      MatFDColoring - Object for computing a sparse Jacobian via finite differences
1163d9274352SBarry Smith         and coloring
1164639f9d9dSBarry Smith 
1165d9274352SBarry Smith    Level: beginner
1166d9274352SBarry Smith 
1167d9274352SBarry Smith   Concepts: coloring, sparse Jacobian, finite differences
1168d9274352SBarry Smith 
1169d9274352SBarry Smith .seealso:  MatFDColoringCreate()
1170d9274352SBarry Smith S*/
1171e2a1c21fSSatish Balay typedef struct _p_MatFDColoring* MatFDColoring;
1172639f9d9dSBarry Smith 
1173*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringCreate(Mat,ISColoring,MatFDColoring *);
1174*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringDestroy(MatFDColoring*);
1175*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringView(MatFDColoring,PetscViewer);
1176*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*);
1177*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**);
1178*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal);
1179*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetFromOptions(MatFDColoring);
1180*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringApply(Mat,MatFDColoring,Vec,MatStructure*,void *);
1181*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetF(MatFDColoring,Vec);
1182*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,PetscInt*[]);
1183b1683b59SHong Zhang 
1184b1683b59SHong Zhang /*S
1185b9af6bddSHong Zhang      MatTransposeColoring - Object for computing a sparse matrix product C=A*B^T via coloring
1186b1683b59SHong Zhang 
1187b1683b59SHong Zhang    Level: beginner
1188b1683b59SHong Zhang 
1189b1683b59SHong Zhang   Concepts: coloring, sparse matrix product
1190b1683b59SHong Zhang 
1191b9af6bddSHong Zhang .seealso:  MatTransposeColoringCreate()
1192b1683b59SHong Zhang S*/
1193b9af6bddSHong Zhang typedef struct _p_MatTransposeColoring* MatTransposeColoring;
1194b1683b59SHong Zhang 
1195*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeColoringCreate(Mat,ISColoring,MatTransposeColoring *);
1196*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransColoringApplySpToDen(MatTransposeColoring,Mat,Mat);
1197*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransColoringApplyDenToSp(MatTransposeColoring,Mat,Mat);
1198*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeColoringDestroy(MatTransposeColoring*);
1199b1683b59SHong Zhang 
1200639f9d9dSBarry Smith /*
12010752156aSBarry Smith     These routines are for partitioning matrices: currently used only
12023eda8832SBarry Smith   for adjacency matrix, MatCreateMPIAdj().
12030752156aSBarry Smith */
1204ca161407SBarry Smith 
1205d9274352SBarry Smith /*S
1206d9274352SBarry Smith      MatPartitioning - Object for managing the partitioning of a matrix or graph
1207d9274352SBarry Smith 
1208d9274352SBarry Smith    Level: beginner
1209d9274352SBarry Smith 
1210d9274352SBarry Smith   Concepts: partitioning
1211d9274352SBarry Smith 
1212743a1026Svictorle .seealso:  MatPartitioningCreate(), MatPartitioningType
1213d9274352SBarry Smith S*/
121491e9ee9fSBarry Smith typedef struct _p_MatPartitioning* MatPartitioning;
1215d9274352SBarry Smith 
121676bdecfbSBarry Smith /*J
12175bc6e7ccSvictorle     MatPartitioningType - String with the name of a PETSc matrix partitioning or the creation function
1218d9274352SBarry Smith        with an optional dynamic library name, for example
1219d9274352SBarry Smith        http://www.mcs.anl.gov/petsc/lib.a:partitioningcreate()
1220d9274352SBarry Smith 
1221d9274352SBarry Smith    Level: beginner
12220d04baf8SBarry Smith dm
1223b547a13bSvictorle .seealso: MatPartitioningCreate(), MatPartitioning
122476bdecfbSBarry Smith J*/
1225a313700dSBarry Smith #define MatPartitioningType char*
12262692d6eeSBarry Smith #define MATPARTITIONINGCURRENT  "current"
12272692d6eeSBarry Smith #define MATPARTITIONINGSQUARE   "square"
12282692d6eeSBarry Smith #define MATPARTITIONINGPARMETIS "parmetis"
12292692d6eeSBarry Smith #define MATPARTITIONINGCHACO    "chaco"
12302692d6eeSBarry Smith #define MATPARTITIONINGPARTY    "party"
123150d91057SBarry Smith #define MATPARTITIONINGPTSCOTCH "ptscotch"
123271306c60Sjroman 
1233ca161407SBarry Smith 
1234*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningCreate(MPI_Comm,MatPartitioning*);
1235*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetType(MatPartitioning,const MatPartitioningType);
1236*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetNParts(MatPartitioning,PetscInt);
1237*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetAdjacency(MatPartitioning,Mat);
1238*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]);
1239*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []);
1240*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningApply(MatPartitioning,IS*);
1241*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningDestroy(MatPartitioning*);
12422aabb6bbSBarry Smith 
1243*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatPartitioning));
124430de9b25SBarry Smith 
124530de9b25SBarry Smith /*MC
124630de9b25SBarry Smith    MatPartitioningRegisterDynamic - Adds a new sparse matrix partitioning to the
124730de9b25SBarry Smith    matrix package.
124830de9b25SBarry Smith 
124930de9b25SBarry Smith    Synopsis:
12501890ba74SBarry Smith    PetscErrorCode MatPartitioningRegisterDynamic(const char *name_partitioning,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatPartitioning))
125130de9b25SBarry Smith 
125230de9b25SBarry Smith    Not Collective
125330de9b25SBarry Smith 
125430de9b25SBarry Smith    Input Parameters:
12552692d6eeSBarry Smith +  sname - name of partitioning (for example MATPARTITIONINGCURRENT) or parmetis
125630de9b25SBarry Smith .  path - location of library where creation routine is
125730de9b25SBarry Smith .  name - name of function that creates the partitioning type, a string
125830de9b25SBarry Smith -  function - function pointer that creates the partitioning type
125930de9b25SBarry Smith 
126030de9b25SBarry Smith    Level: developer
126130de9b25SBarry Smith 
126230de9b25SBarry Smith    If dynamic libraries are used, then the fourth input argument (function)
126330de9b25SBarry Smith    is ignored.
126430de9b25SBarry Smith 
126530de9b25SBarry Smith    Sample usage:
126630de9b25SBarry Smith .vb
126730de9b25SBarry Smith    MatPartitioningRegisterDynamic("my_part",/home/username/my_lib/lib/libO/solaris/mylib.a,
126830de9b25SBarry Smith                "MyPartCreate",MyPartCreate);
126930de9b25SBarry Smith .ve
127030de9b25SBarry Smith 
127130de9b25SBarry Smith    Then, your partitioner can be chosen with the procedural interface via
127230de9b25SBarry Smith $     MatPartitioningSetType(part,"my_part")
127330de9b25SBarry Smith    or at runtime via the option
127430de9b25SBarry Smith $     -mat_partitioning_type my_part
127530de9b25SBarry Smith 
1276ab901514SBarry Smith    $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
127730de9b25SBarry Smith 
127830de9b25SBarry Smith .keywords: matrix, partitioning, register
127930de9b25SBarry Smith 
128030de9b25SBarry Smith .seealso: MatPartitioningRegisterDestroy(), MatPartitioningRegisterAll()
128130de9b25SBarry Smith M*/
1282aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1283f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,0)
12842aabb6bbSBarry Smith #else
1285f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,d)
12862aabb6bbSBarry Smith #endif
12872aabb6bbSBarry Smith 
1288*014dd563SJed Brown PETSC_EXTERN PetscBool MatPartitioningRegisterAllCalled;
1289f1144a30SSatish Balay 
1290*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningRegisterAll(const char[]);
1291*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningRegisterDestroy(void);
12922bad1931SBarry Smith 
1293*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningView(MatPartitioning,PetscViewer);
1294*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetFromOptions(MatPartitioning);
1295*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningGetType(MatPartitioning,const MatPartitioningType*);
1296ca161407SBarry Smith 
1297*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetCoarseSequential(MatPartitioning);
1298*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningParmetisGetEdgeCut(MatPartitioning, PetscInt *);
12990752156aSBarry Smith 
1300b6956c12SJose Roman typedef enum { MP_CHACO_MULTILEVEL=1,MP_CHACO_SPECTRAL=2,MP_CHACO_LINEAR=4,MP_CHACO_RANDOM=5,MP_CHACO_SCATTERED=6 } MPChacoGlobalType;
1301*014dd563SJed Brown PETSC_EXTERN const char *MPChacoGlobalTypes[];
1302b6956c12SJose Roman typedef enum { MP_CHACO_KERNIGHAN=1,MP_CHACO_NONE=2 } MPChacoLocalType;
1303*014dd563SJed Brown PETSC_EXTERN const char *MPChacoLocalTypes[];
1304b6956c12SJose Roman typedef enum { MP_CHACO_LANCZOS=0,MP_CHACO_RQI=1 } MPChacoEigenType;
1305*014dd563SJed Brown PETSC_EXTERN const char *MPChacoEigenTypes[];
1306b6956c12SJose Roman 
1307*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetGlobal(MatPartitioning,MPChacoGlobalType);
1308*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetGlobal(MatPartitioning,MPChacoGlobalType*);
1309*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetLocal(MatPartitioning,MPChacoLocalType);
1310*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetLocal(MatPartitioning,MPChacoLocalType*);
1311*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal);
1312*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType);
1313*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenSolver(MatPartitioning,MPChacoEigenType*);
1314*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenTol(MatPartitioning,PetscReal);
1315*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenTol(MatPartitioning,PetscReal*);
1316*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenNumber(MatPartitioning,PetscInt);
1317*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenNumber(MatPartitioning,PetscInt*);
131871306c60Sjroman 
131971306c60Sjroman #define MP_PARTY_OPT "opt"
132071306c60Sjroman #define MP_PARTY_LIN "lin"
132171306c60Sjroman #define MP_PARTY_SCA "sca"
132271306c60Sjroman #define MP_PARTY_RAN "ran"
132371306c60Sjroman #define MP_PARTY_GBF "gbf"
132471306c60Sjroman #define MP_PARTY_GCF "gcf"
132571306c60Sjroman #define MP_PARTY_BUB "bub"
132671306c60Sjroman #define MP_PARTY_DEF "def"
1327*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetGlobal(MatPartitioning,const char*);
132871306c60Sjroman #define MP_PARTY_HELPFUL_SETS "hs"
132971306c60Sjroman #define MP_PARTY_KERNIGHAN_LIN "kl"
133071306c60Sjroman #define MP_PARTY_NONE "no"
1331*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetLocal(MatPartitioning,const char*);
1332*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal);
1333*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetBipart(MatPartitioning,PetscBool);
1334*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscBool);
133571306c60Sjroman 
133650d91057SBarry Smith typedef enum { MP_PTSCOTCH_QUALITY,MP_PTSCOTCH_SPEED,MP_PTSCOTCH_BALANCE,MP_PTSCOTCH_SAFETY,MP_PTSCOTCH_SCALABILITY } MPPTScotchStrategyType;
1337*014dd563SJed Brown PETSC_EXTERN const char *MPPTScotchStrategyTypes[];
1338e0f1cffaSJose Roman 
1339*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetImbalance(MatPartitioning,PetscReal);
1340*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetImbalance(MatPartitioning,PetscReal*);
1341*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetStrategy(MatPartitioning,MPPTScotchStrategyType);
1342*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetStrategy(MatPartitioning,MPPTScotchStrategyType*);
134371306c60Sjroman 
1344b43b03e9SMark F. Adams /*
1345b43b03e9SMark F. Adams     These routines are for coarsening matrices:
1346b43b03e9SMark F. Adams */
1347b43b03e9SMark F. Adams 
1348b43b03e9SMark F. Adams /*S
1349b43b03e9SMark F. Adams      MatCoarsen - Object for managing the coarsening of a graph (symmetric matrix)
1350b43b03e9SMark F. Adams 
1351b43b03e9SMark F. Adams    Level: beginner
1352b43b03e9SMark F. Adams 
1353b43b03e9SMark F. Adams   Concepts: coarsen
1354b43b03e9SMark F. Adams 
1355b43b03e9SMark F. Adams .seealso:  MatCoarsenCreate), MatCoarsenType
1356b43b03e9SMark F. Adams S*/
1357b43b03e9SMark F. Adams typedef struct _p_MatCoarsen* MatCoarsen;
1358b43b03e9SMark F. Adams 
1359b43b03e9SMark F. Adams /*J
1360b43b03e9SMark F. Adams     MatCoarsenType - String with the name of a PETSc matrix coarsen or the creation function
1361b43b03e9SMark F. Adams        with an optional dynamic library name, for example
1362b43b03e9SMark F. Adams        http://www.mcs.anl.gov/petsc/lib.a:coarsencreate()
1363b43b03e9SMark F. Adams 
1364b43b03e9SMark F. Adams    Level: beginner
1365b43b03e9SMark F. Adams dm
1366b43b03e9SMark F. Adams .seealso: MatCoarsenCreate(), MatCoarsen
1367b43b03e9SMark F. Adams J*/
1368b43b03e9SMark F. Adams #define MatCoarsenType char*
1369b43b03e9SMark F. Adams #define MATCOARSENMIS  "mis"
13709057884aSMark F. Adams #define MATCOARSENHEM  "hem"
1371b43b03e9SMark F. Adams 
13720cbbd2e1SMark F. Adams /* linked list for aggregates */
137341b27cdeSMark F. Adams typedef struct _PetscCDIntNd{
137441b27cdeSMark F. Adams   struct _PetscCDIntNd *next;
13750cbbd2e1SMark F. Adams   PetscInt   gid;
137641b27cdeSMark F. Adams }PetscCDIntNd;
13770cbbd2e1SMark F. Adams 
13780cbbd2e1SMark F. Adams /* only used by node pool */
137941b27cdeSMark F. Adams typedef struct _PetscCDArrNd{
138041b27cdeSMark F. Adams   struct _PetscCDArrNd *next;
138141b27cdeSMark F. Adams   struct _PetscCDIntNd *array;
138241b27cdeSMark F. Adams }PetscCDArrNd;
13830cbbd2e1SMark F. Adams 
13840cbbd2e1SMark F. Adams typedef struct _PetscCoarsenData{
13850cbbd2e1SMark F. Adams   /* node pool */
138641b27cdeSMark F. Adams   PetscCDArrNd  pool_list;
138741b27cdeSMark F. Adams   PetscCDIntNd *new_node;
13880cbbd2e1SMark F. Adams   PetscInt   new_left;
13890cbbd2e1SMark F. Adams   PetscInt   chk_sz;
139041b27cdeSMark F. Adams   PetscCDIntNd *extra_nodes;
13910cbbd2e1SMark F. Adams   /* Array of lists */
139241b27cdeSMark F. Adams   PetscCDIntNd **array;
13930cbbd2e1SMark F. Adams   PetscInt size;
13940cbbd2e1SMark F. Adams   /* cache a Mat for communication data */
13950cbbd2e1SMark F. Adams   Mat mat;
1396a94c3b12SMark F. Adams   /* cache IS of removed equations */
1397a94c3b12SMark F. Adams   IS removedIS;
13980cbbd2e1SMark F. Adams }PetscCoarsenData;
13990cbbd2e1SMark F. Adams 
1400*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenCreate(MPI_Comm,MatCoarsen*);
1401*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetType(MatCoarsen,const MatCoarsenType);
1402*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetAdjacency(MatCoarsen,Mat);
1403*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetGreedyOrdering(MatCoarsen,const IS);
1404*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetStrictAggs(MatCoarsen,PetscBool);
1405*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetVerbose(MatCoarsen,PetscInt);
1406*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenGetData( MatCoarsen, PetscCoarsenData ** );
1407*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenApply(MatCoarsen);
1408*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenDestroy(MatCoarsen*);
1409b43b03e9SMark F. Adams 
1410*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatCoarsen));
1411b43b03e9SMark F. Adams 
1412b43b03e9SMark F. Adams /*MC
1413b43b03e9SMark F. Adams    MatCoarsenRegisterDynamic - Adds a new sparse matrix coarsen to the
1414b43b03e9SMark F. Adams    matrix package.
1415b43b03e9SMark F. Adams 
1416b43b03e9SMark F. Adams    Synopsis:
1417b43b03e9SMark F. Adams    PetscErrorCode MatCoarsenRegisterDynamic(const char *name_coarsen,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatCoarsen))
1418b43b03e9SMark F. Adams 
1419b43b03e9SMark F. Adams    Not Collective
1420b43b03e9SMark F. Adams 
1421b43b03e9SMark F. Adams    Input Parameters:
1422b43b03e9SMark F. Adams +  sname - name of coarsen (for example MATCOARSENMIS)
1423b43b03e9SMark F. Adams .  path - location of library where creation routine is
1424b43b03e9SMark F. Adams .  name - name of function that creates the coarsen type, a string
1425b43b03e9SMark F. Adams -  function - function pointer that creates the coarsen type
1426b43b03e9SMark F. Adams 
1427b43b03e9SMark F. Adams    Level: developer
1428b43b03e9SMark F. Adams 
1429b43b03e9SMark F. Adams    If dynamic libraries are used, then the fourth input argument (function)
1430b43b03e9SMark F. Adams    is ignored.
1431b43b03e9SMark F. Adams 
1432b43b03e9SMark F. Adams    Sample usage:
1433b43b03e9SMark F. Adams .vb
1434b43b03e9SMark F. Adams    MatCoarsenRegisterDynamic("my_agg",/home/username/my_lib/lib/libO/solaris/mylib.a,
1435b43b03e9SMark F. Adams                "MyAggCreate",MyAggCreate);
1436b43b03e9SMark F. Adams .ve
1437b43b03e9SMark F. Adams 
1438b43b03e9SMark F. Adams    Then, your aggregator can be chosen with the procedural interface via
1439b43b03e9SMark F. Adams $     MatCoarsenSetType(agg,"my_agg")
1440b43b03e9SMark F. Adams    or at runtime via the option
1441b43b03e9SMark F. Adams $     -mat_coarsen_type my_agg
1442b43b03e9SMark F. Adams 
1443b43b03e9SMark F. Adams    $PETSC_ARCH occuring in pathname will be replaced with appropriate values.
1444b43b03e9SMark F. Adams 
1445b43b03e9SMark F. Adams .keywords: matrix, coarsen, register
1446b43b03e9SMark F. Adams 
1447b43b03e9SMark F. Adams .seealso: MatCoarsenRegisterDestroy(), MatCoarsenRegisterAll()
1448b43b03e9SMark F. Adams M*/
1449b43b03e9SMark F. Adams #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1450b43b03e9SMark F. Adams #define MatCoarsenRegisterDynamic(a,b,c,d) MatCoarsenRegister(a,b,c,0)
1451b43b03e9SMark F. Adams #else
1452b43b03e9SMark F. Adams #define MatCoarsenRegisterDynamic(a,b,c,d) MatCoarsenRegister(a,b,c,d)
1453b43b03e9SMark F. Adams #endif
1454b43b03e9SMark F. Adams 
1455*014dd563SJed Brown PETSC_EXTERN PetscBool MatCoarsenRegisterAllCalled;
1456b43b03e9SMark F. Adams 
1457*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenRegisterAll(const char[]);
1458*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenRegisterDestroy(void);
1459b43b03e9SMark F. Adams 
1460*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenView(MatCoarsen,PetscViewer);
1461*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetFromOptions(MatCoarsen);
1462*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenGetType(MatCoarsen,const MatCoarsenType*);
1463b43b03e9SMark F. Adams 
1464b43b03e9SMark F. Adams 
1465*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMeshToVertexGraph(Mat,PetscInt,Mat*);
1466*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMeshToCellGraph(Mat,PetscInt,Mat*);
1467591294e4SBarry Smith 
14680752156aSBarry Smith /*
14690a835dfdSSatish Balay     If you add entries here you must also add them to finclude/petscmat.h
1470d4fbbf0eSBarry Smith */
14711c1c02c0SLois Curfman McInnes typedef enum { MATOP_SET_VALUES=0,
14721c1c02c0SLois Curfman McInnes                MATOP_GET_ROW=1,
14731c1c02c0SLois Curfman McInnes                MATOP_RESTORE_ROW=2,
14741c1c02c0SLois Curfman McInnes                MATOP_MULT=3,
14751c1c02c0SLois Curfman McInnes                MATOP_MULT_ADD=4,
14767c922b88SBarry Smith                MATOP_MULT_TRANSPOSE=5,
14777c922b88SBarry Smith                MATOP_MULT_TRANSPOSE_ADD=6,
14781c1c02c0SLois Curfman McInnes                MATOP_SOLVE=7,
14791c1c02c0SLois Curfman McInnes                MATOP_SOLVE_ADD=8,
14807c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE=9,
14817c922b88SBarry Smith                MATOP_SOLVE_TRANSPOSE_ADD=10,
14821c1c02c0SLois Curfman McInnes                MATOP_LUFACTOR=11,
14831c1c02c0SLois Curfman McInnes                MATOP_CHOLESKYFACTOR=12,
1484a714c06dSBarry Smith                MATOP_SOR=13,
14851c1c02c0SLois Curfman McInnes                MATOP_TRANSPOSE=14,
14861c1c02c0SLois Curfman McInnes                MATOP_GETINFO=15,
14871c1c02c0SLois Curfman McInnes                MATOP_EQUAL=16,
14881c1c02c0SLois Curfman McInnes                MATOP_GET_DIAGONAL=17,
14891c1c02c0SLois Curfman McInnes                MATOP_DIAGONAL_SCALE=18,
14901c1c02c0SLois Curfman McInnes                MATOP_NORM=19,
14911c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_BEGIN=20,
14921c1c02c0SLois Curfman McInnes                MATOP_ASSEMBLY_END=21,
1493d519adbfSMatthew Knepley                MATOP_SET_OPTION=22,
1494d519adbfSMatthew Knepley                MATOP_ZERO_ENTRIES=23,
1495d519adbfSMatthew Knepley                MATOP_ZERO_ROWS=24,
1496d519adbfSMatthew Knepley                MATOP_LUFACTOR_SYMBOLIC=25,
1497d519adbfSMatthew Knepley                MATOP_LUFACTOR_NUMERIC=26,
1498d519adbfSMatthew Knepley                MATOP_CHOLESKY_FACTOR_SYMBOLIC=27,
1499d519adbfSMatthew Knepley                MATOP_CHOLESKY_FACTOR_NUMERIC=28,
1500d519adbfSMatthew Knepley                MATOP_SETUP_PREALLOCATION=29,
1501d519adbfSMatthew Knepley                MATOP_ILUFACTOR_SYMBOLIC=30,
1502d519adbfSMatthew Knepley                MATOP_ICCFACTOR_SYMBOLIC=31,
1503d519adbfSMatthew Knepley                MATOP_GET_ARRAY=32,
1504d519adbfSMatthew Knepley                MATOP_RESTORE_ARRAY=33,
1505d519adbfSMatthew Knepley                MATOP_DUPLICATE=34,
1506d519adbfSMatthew Knepley                MATOP_FORWARD_SOLVE=35,
1507d519adbfSMatthew Knepley                MATOP_BACKWARD_SOLVE=36,
1508d519adbfSMatthew Knepley                MATOP_ILUFACTOR=37,
1509d519adbfSMatthew Knepley                MATOP_ICCFACTOR=38,
1510d519adbfSMatthew Knepley                MATOP_AXPY=39,
1511d519adbfSMatthew Knepley                MATOP_GET_SUBMATRICES=40,
1512d519adbfSMatthew Knepley                MATOP_INCREASE_OVERLAP=41,
1513d519adbfSMatthew Knepley                MATOP_GET_VALUES=42,
1514d519adbfSMatthew Knepley                MATOP_COPY=43,
1515d519adbfSMatthew Knepley                MATOP_GET_ROW_MAX=44,
1516d519adbfSMatthew Knepley                MATOP_SCALE=45,
1517d519adbfSMatthew Knepley                MATOP_SHIFT=46,
151835153367SBarry Smith                MATOP_DIAGONAL_SET=47,
1519d519adbfSMatthew Knepley                MATOP_ILUDT_FACTOR=48,
1520d519adbfSMatthew Knepley                MATOP_SET_BLOCK_SIZE=49,
1521d519adbfSMatthew Knepley                MATOP_GET_ROW_IJ=50,
1522d519adbfSMatthew Knepley                MATOP_RESTORE_ROW_IJ=51,
1523d519adbfSMatthew Knepley                MATOP_GET_COLUMN_IJ=52,
1524d519adbfSMatthew Knepley                MATOP_RESTORE_COLUMN_IJ=53,
1525d519adbfSMatthew Knepley                MATOP_FDCOLORING_CREATE=54,
1526d519adbfSMatthew Knepley                MATOP_COLORING_PATCH=55,
1527d519adbfSMatthew Knepley                MATOP_SET_UNFACTORED=56,
1528d519adbfSMatthew Knepley                MATOP_PERMUTE=57,
1529d519adbfSMatthew Knepley                MATOP_SET_VALUES_BLOCKED=58,
1530d519adbfSMatthew Knepley                MATOP_GET_SUBMATRIX=59,
1531d519adbfSMatthew Knepley                MATOP_DESTROY=60,
1532d519adbfSMatthew Knepley                MATOP_VIEW=61,
1533d519adbfSMatthew Knepley                MATOP_CONVERT_FROM=62,
1534d519adbfSMatthew Knepley                MATOP_USE_SCALED_FORM=63,
1535d519adbfSMatthew Knepley                MATOP_SCALE_SYSTEM=64,
1536d519adbfSMatthew Knepley                MATOP_UNSCALE_SYSTEM=65,
1537d519adbfSMatthew Knepley                MATOP_SET_LOCAL_TO_GLOBAL_MAP=66,
1538d519adbfSMatthew Knepley                MATOP_SET_VALUES_LOCAL=67,
1539d519adbfSMatthew Knepley                MATOP_ZERO_ROWS_LOCAL=68,
1540d519adbfSMatthew Knepley                MATOP_GET_ROW_MAX_ABS=69,
1541d519adbfSMatthew Knepley                MATOP_GET_ROW_MIN_ABS=70,
1542d519adbfSMatthew Knepley                MATOP_CONVERT=71,
1543d519adbfSMatthew Knepley                MATOP_SET_COLORING=72,
1544d519adbfSMatthew Knepley                MATOP_SET_VALUES_ADIC=73,
1545d519adbfSMatthew Knepley                MATOP_SET_VALUES_ADIFOR=74,
1546d519adbfSMatthew Knepley                MATOP_FD_COLORING_APPLY=75,
1547d519adbfSMatthew Knepley                MATOP_SET_FROM_OPTIONS=76,
1548d519adbfSMatthew Knepley                MATOP_MULT_CON=77,
1549d519adbfSMatthew Knepley                MATOP_MULT_TRANSPOSE_CON=78,
1550d519adbfSMatthew Knepley                MATOP_PERMUTE_SPARSIFY=79,
1551d519adbfSMatthew Knepley                MATOP_MULT_MULTIPLE=80,
1552d519adbfSMatthew Knepley                MATOP_SOLVE_MULTIPLE=81,
1553d519adbfSMatthew Knepley                MATOP_GET_INERTIA=82,
1554d519adbfSMatthew Knepley                MATOP_LOAD=83,
1555d519adbfSMatthew Knepley                MATOP_IS_SYMMETRIC=84,
1556d519adbfSMatthew Knepley                MATOP_IS_HERMITIAN=85,
1557d519adbfSMatthew Knepley                MATOP_IS_STRUCTURALLY_SYMMETRIC=86,
155841f059aeSBarry Smith                MATOP_DUMMY=87,
1559d519adbfSMatthew Knepley                MATOP_GET_VECS=88,
1560d519adbfSMatthew Knepley                MATOP_MAT_MULT=89,
1561d519adbfSMatthew Knepley                MATOP_MAT_MULT_SYMBOLIC=90,
1562d519adbfSMatthew Knepley                MATOP_MAT_MULT_NUMERIC=91,
1563d519adbfSMatthew Knepley                MATOP_PTAP=92,
1564d519adbfSMatthew Knepley                MATOP_PTAP_SYMBOLIC=93,
1565d519adbfSMatthew Knepley                MATOP_PTAP_NUMERIC=94,
1566d519adbfSMatthew Knepley                MATOP_MAT_MULTTRANSPOSE=95,
15671763ac4eSSatish Balay                MATOP_MAT_MULTTRANSPOSE_SYM=96,
15681763ac4eSSatish Balay                MATOP_MAT_MULTTRANSPOSE_NUM=97,
1569d519adbfSMatthew Knepley                MATOP_PTAP_SYMBOLIC_SEQAIJ=98,
1570d519adbfSMatthew Knepley                MATOP_PTAP_NUMERIC_SEQAIJ=99,
1571d519adbfSMatthew Knepley                MATOP_PTAP_SYMBOLIC_MPIAIJ=100,
1572d519adbfSMatthew Knepley                MATOP_PTAP_NUMERIC_MPIAIJ=101,
1573d519adbfSMatthew Knepley                MATOP_CONJUGATE=102,
1574d519adbfSMatthew Knepley                MATOP_SET_SIZES=103,
1575d519adbfSMatthew Knepley                MATOP_SET_VALUES_ROW=104,
1576d519adbfSMatthew Knepley                MATOP_REAL_PART=105,
1577d519adbfSMatthew Knepley                MATOP_IMAG_PART=106,
1578d519adbfSMatthew Knepley                MATOP_GET_ROW_UTRIANGULAR=107,
1579d519adbfSMatthew Knepley                MATOP_RESTORE_ROW_UTRIANGULAR=108,
1580d519adbfSMatthew Knepley                MATOP_MATSOLVE=109,
1581d519adbfSMatthew Knepley                MATOP_GET_REDUNDANTMATRIX=110,
1582d519adbfSMatthew Knepley                MATOP_GET_ROW_MIN=111,
1583d519adbfSMatthew Knepley                MATOP_GET_COLUMN_VEC=112,
1584d519adbfSMatthew Knepley                MATOP_MISSING_DIAGONAL=113,
1585d519adbfSMatthew Knepley                MATOP_MATGETSEQNONZEROSTRUCTURE=114,
158689c6957cSBarry Smith                MATOP_CREATE=115,
158789c6957cSBarry Smith                MATOP_GET_GHOSTS=116,
1588ea38565cSJed Brown                MATOP_GET_LOCALSUBMATRIX=117,
1589ea38565cSJed Brown                MATOP_RESTORE_LOCALSUBMATRIX=118,
1590eeffb40dSHong Zhang                MATOP_MULT_DIAGONAL_BLOCK=119,
1591547795f9SHong Zhang                MATOP_HERMITIANTRANSPOSE=120,
1592547795f9SHong Zhang                MATOP_MULTHERMITIANTRANSPOSE=121,
1593fbdbba38SShri Abhyankar                MATOP_MULTHERMITIANTRANSPOSEADD=122,
1594b9614d88SDmitry Karpeev                MATOP_GETMULTIPROCBLOCK=123,
15950716a85fSBarry Smith                MATOP_GETCOLUMNNORMS=125,
159637868618SMatthew G Knepley 	       MATOP_GET_SUBMATRICES_PARALLEL=128,
15972b8ad9a3SHong Zhang                MATOP_SET_VALUES_BATCH=129,
15982b8ad9a3SHong Zhang                MATOP_TRANSPOSEMATMULT=130,
15992b8ad9a3SHong Zhang                MATOP_TRANSPOSEMATMULT_SYMBOLIC=131,
16002b8ad9a3SHong Zhang                MATOP_TRANSPOSEMATMULT_NUMERIC=132,
16012b8ad9a3SHong Zhang                MATOP_TRANSPOSECOLORING_CREATE=133,
16022b8ad9a3SHong Zhang                MATOP_TRANSCOLORING_APPLY_SPTODEN=134,
16032b8ad9a3SHong Zhang                MATOP_TRANSCOLORING_APPLY_DENTOSP=135,
16042b8ad9a3SHong Zhang                MATOP_RARt=136,
16052b8ad9a3SHong Zhang                MATOP_RARt_SYMBOLIC=137,
16062274620bSMark F. Adams                MATOP_RARt_NUMERIC=138,
16072274620bSMark F. Adams                MATOP_SET_BLOCK_SIZES=139
1608fae171e0SBarry Smith              } MatOperation;
1609*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatHasOperation(Mat,MatOperation,PetscBool *);
1610*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellSetOperation(Mat,MatOperation,void(*)(void));
1611*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellGetOperation(Mat,MatOperation,void(**)(void));
1612*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellSetContext(Mat,void*);
1613112a2221SBarry Smith 
161490ace30eSBarry Smith /*
161590ace30eSBarry Smith    Codes for matrices stored on disk. By default they are
161690ace30eSBarry Smith    stored in a universal format. By changing the format with
16177973ad04SBarry Smith    PetscViewerSetFormat(viewer,PETSC_VIEWER_NATIVE); the matrices will
161890ace30eSBarry Smith    be stored in a way natural for the matrix, for example dense matrices
161990ace30eSBarry Smith    would be stored as dense. Matrices stored this way may only be
1620f4403165SShri Abhyankar    read into matrices of the same type.
162190ace30eSBarry Smith */
162290ace30eSBarry Smith #define MATRIX_BINARY_FORMAT_DENSE -1
162390ace30eSBarry Smith 
1624*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetHashTableFactor(Mat,PetscReal);
1625*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatISGetLocalMat(Mat,Mat*);
1626*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatISSetLocalMat(Mat,Mat);
16271f4e1ec7SBarry Smith 
1628d9274352SBarry Smith /*S
1629d9274352SBarry Smith      MatNullSpace - Object that removes a null space from a vector, i.e.
1630d9274352SBarry Smith          orthogonalizes the vector to a subsapce
1631d9274352SBarry Smith 
1632f7a9e4ceSBarry Smith    Level: advanced
1633d9274352SBarry Smith 
1634d9274352SBarry Smith   Concepts: matrix; linear operator, null space
1635d9274352SBarry Smith 
16366e1639daSBarry Smith   Users manual sections:
16376e1639daSBarry Smith .   sec_singular
16386e1639daSBarry Smith 
1639d9274352SBarry Smith .seealso:  MatNullSpaceCreate()
1640d9274352SBarry Smith S*/
164174637425SBarry Smith typedef struct _p_MatNullSpace* MatNullSpace;
1642d9274352SBarry Smith 
1643*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceCreate(MPI_Comm,PetscBool ,PetscInt,const Vec[],MatNullSpace*);
1644*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(MatNullSpace,Vec,void*),void*);
1645*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceDestroy(MatNullSpace*);
1646*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceRemove(MatNullSpace,Vec,Vec*);
1647*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNullSpace(Mat, MatNullSpace *);
1648*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNullSpace(Mat,MatNullSpace);
1649*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNearNullSpace(Mat,MatNullSpace);
1650*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNearNullSpace(Mat,MatNullSpace*);
1651*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceTest(MatNullSpace,Mat,PetscBool  *);
1652*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceView(MatNullSpace,PetscViewer);
1653*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceGetVecs(MatNullSpace,PetscBool*,PetscInt*,const Vec**);
1654*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceCreateRigidBody(Vec,MatNullSpace*);
165574637425SBarry Smith 
1656*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatReorderingSeqSBAIJ(Mat,IS);
1657*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetHashTableFactor(Mat,PetscReal);
1658*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetColumnIndices(Mat,PetscInt *);
1659*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJInvertBlockDiagonal(Mat);
16603f1d51d7SBarry Smith 
1661*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMAIJ(Mat,PetscInt,Mat*);
1662*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMAIJRedimension(Mat,PetscInt,Mat*);
1663*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMAIJGetAIJ(Mat,Mat*);
1664c4f061fbSSatish Balay 
1665*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatComputeExplicitOperator(Mat,Mat*);
1666b0a32e0cSBarry Smith 
1667*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalScaleLocal(Mat,Vec);
166804f1ad80SBarry Smith 
1669*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*);
1670*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetBase(Mat,Vec,Vec);
1671*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*);
1672*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*));
1673*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec));
1674*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDAddNullSpace(Mat,MatNullSpace);
1675*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt);
1676*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDResetHHistory(Mat);
1677*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctionError(Mat,PetscReal);
1678*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetPeriod(Mat,PetscInt);
1679*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDGetH(Mat,PetscScalar *);
1680*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetOptionsPrefix(Mat,const char[]);
1681*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*);
1682*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*);
1683e884886fSBarry Smith 
16846370053bSBarry Smith /*S
16856370053bSBarry Smith     MatMFFD - A data structured used to manage the computation of the h differencing parameter for matrix-free
16866370053bSBarry Smith               Jacobian vector products
1687e884886fSBarry Smith 
16886370053bSBarry Smith     Notes: MATMFFD is a specific MatType which uses the MatMFFD data structure
16896370053bSBarry Smith 
16906370053bSBarry Smith            MatMFFD*() methods actually take the Mat as their first argument. Not a MatMFFD data structure
16916370053bSBarry Smith 
16926370053bSBarry Smith     Level: developer
16936370053bSBarry Smith 
16946370053bSBarry Smith .seealso: MATMFFD, MatCreateMFFD(), MatMFFDSetFuction(), MatMFFDSetType(), MatMFFDRegister()
16956370053bSBarry Smith S*/
1696e884886fSBarry Smith typedef struct _p_MatMFFD* MatMFFD;
1697e884886fSBarry Smith 
169876bdecfbSBarry Smith /*J
1699e884886fSBarry Smith     MatMFFDType - algorithm used to compute the h used in computing matrix-vector products via differencing of the function
1700e884886fSBarry Smith 
1701e884886fSBarry Smith    Level: beginner
1702e884886fSBarry Smith 
1703e884886fSBarry Smith .seealso: MatMFFDSetType(), MatMFFDRegister()
170476bdecfbSBarry Smith J*/
1705a313700dSBarry Smith #define MatMFFDType char*
1706e884886fSBarry Smith #define MATMFFD_DS  "ds"
1707e884886fSBarry Smith #define MATMFFD_WP  "wp"
1708e884886fSBarry Smith 
1709*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetType(Mat,const MatMFFDType);
1710*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatMFFD));
1711e884886fSBarry Smith 
1712e884886fSBarry Smith /*MC
1713e884886fSBarry Smith    MatMFFDRegisterDynamic - Adds a method to the MatMFFD registry.
1714e884886fSBarry Smith 
1715e884886fSBarry Smith    Synopsis:
17161890ba74SBarry Smith    PetscErrorCode MatMFFDRegisterDynamic(const char *name_solver,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatMFFD))
1717e884886fSBarry Smith 
1718e884886fSBarry Smith    Not Collective
1719e884886fSBarry Smith 
1720e884886fSBarry Smith    Input Parameters:
1721e884886fSBarry Smith +  name_solver - name of a new user-defined compute-h module
1722e884886fSBarry Smith .  path - path (either absolute or relative) the library containing this solver
1723e884886fSBarry Smith .  name_create - name of routine to create method context
1724e884886fSBarry Smith -  routine_create - routine to create method context
1725e884886fSBarry Smith 
1726e884886fSBarry Smith    Level: developer
1727e884886fSBarry Smith 
1728e884886fSBarry Smith    Notes:
1729e884886fSBarry Smith    MatMFFDRegisterDynamic() may be called multiple times to add several user-defined solvers.
1730e884886fSBarry Smith 
1731e884886fSBarry Smith    If dynamic libraries are used, then the fourth input argument (routine_create)
1732e884886fSBarry Smith    is ignored.
1733e884886fSBarry Smith 
1734e884886fSBarry Smith    Sample usage:
1735e884886fSBarry Smith .vb
1736e884886fSBarry Smith    MatMFFDRegisterDynamic("my_h",/home/username/my_lib/lib/libO/solaris/mylib.a,
1737e884886fSBarry Smith                "MyHCreate",MyHCreate);
1738e884886fSBarry Smith .ve
1739e884886fSBarry Smith 
1740e884886fSBarry Smith    Then, your solver can be chosen with the procedural interface via
1741e884886fSBarry Smith $     MatMFFDSetType(mfctx,"my_h")
1742e884886fSBarry Smith    or at runtime via the option
1743e884886fSBarry Smith $     -snes_mf_type my_h
1744e884886fSBarry Smith 
1745e884886fSBarry Smith .keywords: MatMFFD, register
1746e884886fSBarry Smith 
1747e884886fSBarry Smith .seealso: MatMFFDRegisterAll(), MatMFFDRegisterDestroy()
1748e884886fSBarry Smith M*/
1749e884886fSBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
1750e884886fSBarry Smith #define MatMFFDRegisterDynamic(a,b,c,d) MatMFFDRegister(a,b,c,0)
1751e884886fSBarry Smith #else
1752e884886fSBarry Smith #define MatMFFDRegisterDynamic(a,b,c,d) MatMFFDRegister(a,b,c,d)
1753e884886fSBarry Smith #endif
1754e884886fSBarry Smith 
1755*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDRegisterAll(const char[]);
1756*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDRegisterDestroy(void);
1757*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDDSSetUmin(Mat,PetscReal);
1758*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDWPSetComputeNormU(Mat,PetscBool );
1759e884886fSBarry Smith 
1760e884886fSBarry Smith 
1761*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *);
1762*014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *);
17637dbadf16SMatthew Knepley 
176497969023SHong Zhang /*
176597969023SHong Zhang    PETSc interface to MUMPS
176697969023SHong Zhang */
176797969023SHong Zhang #ifdef PETSC_HAVE_MUMPS
1768*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMumpsSetIcntl(Mat,PetscInt,PetscInt);
176997969023SHong Zhang #endif
177023a5497aSJed Brown 
1771d954db57SHong Zhang /*
1772d954db57SHong Zhang    PETSc interface to SUPERLU
1773d954db57SHong Zhang */
1774d954db57SHong Zhang #ifdef PETSC_HAVE_SUPERLU
1775*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSuperluSetILUDropTol(Mat,PetscReal);
1776d954db57SHong Zhang #endif
1777d954db57SHong Zhang 
177890c53902SBarry Smith #if defined(PETSC_HAVE_CUSP)
1779*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCUSP(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*);
1780*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateAIJCUSP(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*);
178190c53902SBarry Smith #endif
178290c53902SBarry Smith 
178354efbe56SHong Zhang /*
178454efbe56SHong Zhang    PETSc interface to FFTW
178554efbe56SHong Zhang */
178654efbe56SHong Zhang #if defined(PETSC_HAVE_FFTW)
1787*014dd563SJed Brown PETSC_EXTERN PetscErrorCode VecScatterPetscToFFTW(Mat,Vec,Vec);
1788*014dd563SJed Brown PETSC_EXTERN PetscErrorCode VecScatterFFTWToPetsc(Mat,Vec,Vec);
1789*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetVecsFFTW(Mat,Vec*,Vec*,Vec*);
179054efbe56SHong Zhang #endif
179154efbe56SHong Zhang 
1792*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateNest(MPI_Comm,PetscInt,const IS[],PetscInt,const IS[],const Mat[],Mat*);
1793*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSize(Mat,PetscInt*,PetscInt*);
1794*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetISs(Mat,IS[],IS[]);
1795*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetLocalISs(Mat,IS[],IS[]);
1796*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSubMats(Mat,PetscInt*,PetscInt*,Mat***);
1797*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSubMat(Mat,PetscInt,PetscInt,Mat*);
1798*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestSetVecType(Mat,const VecType);
1799*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestSetSubMats(Mat,PetscInt,const IS[],PetscInt,const IS[],const Mat[]);
1800*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestSetSubMat(Mat,PetscInt,PetscInt,Mat);
1801d8588912SDave May 
180215c1f1baSDmitry Karpeev /*
180315c1f1baSDmitry Karpeev  MatIJ:
180415c1f1baSDmitry Karpeev  An unweighted directed pseudograph
180515c1f1baSDmitry Karpeev  An interpretation of this matrix as a (pseudo)graph allows us to define additional operations on it:
180615c1f1baSDmitry Karpeev  A MatIJ can act on sparse arrays: arrays of indices, or index arrays of integers, scalars, or integer-scalar pairs
180715c1f1baSDmitry Karpeev  by mapping the indices to the indices connected to them by the (pseudo)graph ed
180815c1f1baSDmitry Karpeev  */
180915c1f1baSDmitry Karpeev typedef enum {MATIJ_LOCAL, MATIJ_GLOBAL} MatIJIndexType;
1810*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIJSetMultivalued(Mat, PetscBool);
1811*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIJGetMultivalued(Mat, PetscBool*);
1812*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIJSetEdges(Mat, PetscInt, const PetscInt*, const PetscInt*);
1813*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIJGetEdges(Mat, PetscInt *, PetscInt **, PetscInt **);
1814*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIJSetEdgesIS(Mat, IS, IS);
1815*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIJGetEdgesIS(Mat, IS*, IS*);
1816*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIJGetRowSizes(Mat, MatIJIndexType, PetscInt, const PetscInt *, PetscInt **);
1817*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIJGetMinRowSize(Mat, PetscInt *);
1818*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIJGetMaxRowSize(Mat, PetscInt *);
1819*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIJGetSupport(Mat,  PetscInt *, PetscInt **);
1820*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIJGetSupportIS(Mat, IS *);
1821*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIJGetImage(Mat, PetscInt*, PetscInt**);
1822*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIJGetImageIS(Mat, IS *);
1823*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIJGetSupportSize(Mat, PetscInt *);
1824*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIJGetImageSize(Mat, PetscInt *);
182515c1f1baSDmitry Karpeev 
1826*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIJBinRenumber(Mat, Mat*);
182715c1f1baSDmitry Karpeev 
1828*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIJMap(Mat, MatIJIndexType, PetscInt,const PetscInt*,const PetscInt*,const PetscScalar*, MatIJIndexType,PetscInt*,PetscInt**,PetscInt**,PetscScalar**,PetscInt**);
1829*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIJBin(Mat, MatIJIndexType, PetscInt,const PetscInt*,const PetscInt*,const PetscScalar*,PetscInt*,PetscInt**,PetscInt**,PetscScalar**,PetscInt**);
1830*014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIJBinMap(Mat,Mat, MatIJIndexType,PetscInt,const PetscInt*,const PetscInt*,const PetscScalar*,MatIJIndexType,PetscInt*,PetscInt**,PetscInt**,PetscScalar**,PetscInt**);
183115c1f1baSDmitry Karpeev 
183223a5497aSJed Brown #endif
1833