12eac72dbSBarry Smith /* 22eac72dbSBarry Smith Include file for the matrix component of PETSc 32eac72dbSBarry Smith */ 40a835dfdSSatish Balay #ifndef __PETSCMAT_H 50a835dfdSSatish Balay #define __PETSCMAT_H 62c8e378dSBarry Smith #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*/ 2819fd82e9SBarry Smith typedef const char* MatType; 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" 459ae82921SPaul Mullowney #define MATAIJCUSPARSE "aijcusparse" 469ae82921SPaul Mullowney #define MATSEQAIJCUSPARSE "seqaijcusparse" 479ae82921SPaul Mullowney #define MATMPIAIJCUSPARSE "mpiaijcusparse" 485a11e1b2SBarry Smith #define MATAIJPERM "aijperm" 495a11e1b2SBarry Smith #define MATSEQAIJPERM "seqaijperm" 505a11e1b2SBarry Smith #define MATMPIAIJPERM "mpiaijperm" 51273d9f13SBarry Smith #define MATSHELL "shell" 525a11e1b2SBarry Smith #define MATDENSE "dense" 53209238afSKris Buschelman #define MATSEQDENSE "seqdense" 54273d9f13SBarry Smith #define MATMPIDENSE "mpidense" 55db31f6deSJed Brown #define MATELEMENTAL "elemental" 565a11e1b2SBarry Smith #define MATBAIJ "baij" 57273d9f13SBarry Smith #define MATSEQBAIJ "seqbaij" 58273d9f13SBarry Smith #define MATMPIBAIJ "mpibaij" 59273d9f13SBarry Smith #define MATMPIADJ "mpiadj" 605a11e1b2SBarry Smith #define MATSBAIJ "sbaij" 61273d9f13SBarry Smith #define MATSEQSBAIJ "seqsbaij" 62273d9f13SBarry Smith #define MATMPISBAIJ "mpisbaij" 63c0cdd4a1SDahai Guo #define MATSEQBSTRM "seqbstrm" 64c0cdd4a1SDahai Guo #define MATMPIBSTRM "mpibstrm" 65c0cdd4a1SDahai Guo #define MATBSTRM "bstrm" 66aa5a9175SDahai Guo #define MATSEQSBSTRM "seqsbstrm" 67aa5a9175SDahai Guo #define MATMPISBSTRM "mpisbstrm" 68aa5a9175SDahai Guo #define MATSBSTRM "sbstrm" 69cebc7f6cSBarry Smith #define MATDAAD "daad" 70cebc7f6cSBarry Smith #define MATMFFD "mffd" 71c8a8475eSBarry Smith #define MATNORMAL "normal" 72ab92ecdeSBarry Smith #define MATLRC "lrc" 732a6744ebSBarry Smith #define MATSCATTER "scatter" 74421e10b8SBarry Smith #define MATBLOCKMAT "blockmat" 75793850ffSBarry Smith #define MATCOMPOSITE "composite" 761f8c7532SHong Zhang #define MATFFT "fft" 771f8c7532SHong Zhang #define MATFFTW "fftw" 78e133240eSMatthew G Knepley #define MATSEQCUFFT "seqcufft" 79557cca28SSatish Balay #define MATTRANSPOSEMAT "transpose" 8072ca8751SBarry Smith #define MATSCHURCOMPLEMENT "schurcomplement" 811d6018f0SLisandro Dalcin #define MATPYTHON "python" 82f91d8e95SBarry Smith #define MATHYPRESTRUCT "hyprestruct" 83a9e6138eSGlenn Hammond #define MATHYPRESSTRUCT "hypresstruct" 84ee1cef2cSJed Brown #define MATSUBMATRIX "submatrix" 852c0dbf93SJed Brown #define MATLOCALREF "localref" 86d8588912SDave May #define MATNEST "nest" 8715c1f1baSDmitry Karpeev #define MATIJ "ij" 88773941d6SBarry Smith 8976bdecfbSBarry Smith /*J 90c7393fdbSBarry Smith MatSolverPackage - String with the name of a PETSc matrix solver type. 919989ab13SBarry Smith 929989ab13SBarry Smith For example: "petsc" indicates what PETSc provides, "superlu" indicates either 939989ab13SBarry Smith SuperLU or SuperLU_Dist etc. 949989ab13SBarry Smith 959989ab13SBarry Smith 969989ab13SBarry Smith Level: beginner 979989ab13SBarry Smith 985c9eb25fSBarry Smith .seealso: MatGetFactor(), Mat, MatSetType(), MatType 9976bdecfbSBarry Smith J*/ 100c7393fdbSBarry Smith #define MatSolverPackage char* 1012692d6eeSBarry Smith #define MATSOLVERSUPERLU "superlu" 1022692d6eeSBarry Smith #define MATSOLVERSUPERLU_DIST "superlu_dist" 1032692d6eeSBarry Smith #define MATSOLVERUMFPACK "umfpack" 10420db9a53SJed Brown #define MATSOLVERCHOLMOD "cholmod" 1052692d6eeSBarry Smith #define MATSOLVERESSL "essl" 1062692d6eeSBarry Smith #define MATSOLVERLUSOL "lusol" 1072692d6eeSBarry Smith #define MATSOLVERMUMPS "mumps" 1082692d6eeSBarry Smith #define MATSOLVERPASTIX "pastix" 1092692d6eeSBarry Smith #define MATSOLVERMATLAB "matlab" 1102692d6eeSBarry Smith #define MATSOLVERPETSC "petsc" 1112692d6eeSBarry Smith #define MATSOLVERBAS "bas" 1129ae82921SPaul Mullowney #define MATSOLVERCUSPARSE "cusparse" 113aa5a9175SDahai Guo #define MATSOLVERBSTRM "bstrm" 114aa5a9175SDahai Guo #define MATSOLVERSBSTRM "sbstrm" 11515767789SHong Zhang #define MATSOLVERELEMENTAL "elemental" 11617f1a0eaSHong Zhang #define MATSOLVERCLIQUE "clique" 117c0cdd4a1SDahai Guo 118b24902e0SBarry Smith /*E 119b24902e0SBarry Smith MatFactorType - indicates what type of factorization is requested 120b24902e0SBarry Smith 121b24902e0SBarry Smith Level: beginner 122b24902e0SBarry Smith 123b24902e0SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 124b24902e0SBarry Smith 125c7393fdbSBarry Smith .seealso: MatSolverPackage, MatGetFactor() 126b24902e0SBarry Smith E*/ 127599ef60dSHong Zhang typedef enum {MAT_FACTOR_NONE, MAT_FACTOR_LU, MAT_FACTOR_CHOLESKY, MAT_FACTOR_ILU, MAT_FACTOR_ICC,MAT_FACTOR_ILUDT} MatFactorType; 128014dd563SJed Brown PETSC_EXTERN const char *const MatFactorTypes[]; 129e92e720dSBarry Smith 130014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactor(Mat,const MatSolverPackage,MatFactorType,Mat*); 131014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactorAvailable(Mat,const MatSolverPackage,MatFactorType,PetscBool *); 132014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFactorGetSolverPackage(Mat,const MatSolverPackage*); 133014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactorType(Mat,MatFactorType*); 1349989ab13SBarry Smith 135c06d978dSMatthew Knepley /* Logging support */ 1360700a824SBarry Smith #define MAT_FILE_CLASSID 1211216 /* used to indicate matrices in binary files */ 137014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_CLASSID; 138014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_FDCOLORING_CLASSID; 139014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_TRANSPOSECOLORING_CLASSID; 140014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_PARTITIONING_CLASSID; 141014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_COARSEN_CLASSID; 142014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_NULLSPACE_CLASSID; 143014dd563SJed Brown PETSC_EXTERN PetscClassId MATMFFD_CLASSID; 144c06d978dSMatthew Knepley 145ceb03754SKris Buschelman /*E 146ceb03754SKris Buschelman MatReuse - Indicates if matrices obtained from a previous call to MatGetSubMatrices() 147d6eff37eSBarry Smith or MatGetSubMatrix() are to be reused to store the new matrix values. For MatConvert() is used to indicate 148d6eff37eSBarry Smith that the input matrix is to be replaced with the converted matrix. 149ceb03754SKris Buschelman 150ceb03754SKris Buschelman Level: beginner 151ceb03754SKris Buschelman 152ceb03754SKris Buschelman Any additions/changes here MUST also be made in include/finclude/petscmat.h 153ceb03754SKris Buschelman 1540c451bc4SBarry Smith .seealso: MatGetSubMatrices(), MatGetSubMatrix(), MatDestroyMatrices(), MatConvert() 155ceb03754SKris Buschelman E*/ 156dfe085dbSJed Brown typedef enum {MAT_INITIAL_MATRIX,MAT_REUSE_MATRIX,MAT_IGNORE_MATRIX} MatReuse; 157ceb03754SKris Buschelman 1585494a064SHong Zhang /*E 1595494a064SHong Zhang MatGetSubMatrixOption - Indicates if matrices obtained from a call to MatGetSubMatrices() 160829201f2SHong Zhang include the matrix values. Currently it is only used by MatGetSeqNonzerostructure(). 1615494a064SHong Zhang 1625494a064SHong Zhang Level: beginner 1635494a064SHong Zhang 164829201f2SHong Zhang .seealso: MatGetSeqNonzerostructure() 1655494a064SHong Zhang E*/ 1665494a064SHong Zhang typedef enum {MAT_DO_NOT_GET_VALUES,MAT_GET_VALUES} MatGetSubMatrixOption; 1675494a064SHong Zhang 168014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInitializePackage(const char[]); 169c06d978dSMatthew Knepley 170014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreate(MPI_Comm,Mat*); 171014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetSizes(Mat,PetscInt,PetscInt,PetscInt,PetscInt); 17219fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatSetType(Mat,MatType); 173014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetFromOptions(Mat); 1740d2bece7SBarry Smith PETSC_EXTERN PetscErrorCode MatViewFromOptions(Mat,const char[]); 175014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRegisterAll(const char[]); 176014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat)); 177014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRegisterBaseName(const char[],const char[],const char[]); 178014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetOptionsPrefix(Mat,const char[]); 179014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAppendOptionsPrefix(Mat,const char[]); 180014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOptionsPrefix(Mat,const char*[]); 181f69a0ea3SMatthew Knepley 18230de9b25SBarry Smith /*MC 18330de9b25SBarry Smith MatRegisterDynamic - Adds a new matrix type 18430de9b25SBarry Smith 18530de9b25SBarry Smith Synopsis: 1861890ba74SBarry Smith PetscErrorCode MatRegisterDynamic(const char *name,const char *path,const char *name_create,PetscErrorCode (*routine_create)(Mat)) 18730de9b25SBarry Smith 18830de9b25SBarry Smith Not Collective 18930de9b25SBarry Smith 19030de9b25SBarry Smith Input Parameters: 19130de9b25SBarry Smith + name - name of a new user-defined matrix type 19230de9b25SBarry Smith . path - path (either absolute or relative) the library containing this solver 19330de9b25SBarry Smith . name_create - name of routine to create method context 19430de9b25SBarry Smith - routine_create - routine to create method context 19530de9b25SBarry Smith 19630de9b25SBarry Smith Notes: 19730de9b25SBarry Smith MatRegisterDynamic() may be called multiple times to add several user-defined solvers. 19830de9b25SBarry Smith 19930de9b25SBarry Smith If dynamic libraries are used, then the fourth input argument (routine_create) 20030de9b25SBarry Smith is ignored. 20130de9b25SBarry Smith 20230de9b25SBarry Smith Sample usage: 20330de9b25SBarry Smith .vb 20430de9b25SBarry Smith MatRegisterDynamic("my_mat",/home/username/my_lib/lib/libO/solaris/mylib.a, 20530de9b25SBarry Smith "MyMatCreate",MyMatCreate); 20630de9b25SBarry Smith .ve 20730de9b25SBarry Smith 20830de9b25SBarry Smith Then, your solver can be chosen with the procedural interface via 20930de9b25SBarry Smith $ MatSetType(Mat,"my_mat") 21030de9b25SBarry Smith or at runtime via the option 21130de9b25SBarry Smith $ -mat_type my_mat 21230de9b25SBarry Smith 21330de9b25SBarry Smith Level: advanced 21430de9b25SBarry Smith 215ab901514SBarry Smith Notes: ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values. 21630de9b25SBarry Smith If your function is not being put into a shared library then use VecRegister() instead 21730de9b25SBarry Smith 21830de9b25SBarry Smith .keywords: Mat, register 21930de9b25SBarry Smith 22030de9b25SBarry Smith .seealso: MatRegisterAll(), MatRegisterDestroy() 22130de9b25SBarry Smith 22230de9b25SBarry Smith M*/ 223273d9f13SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 224273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,0) 225273d9f13SBarry Smith #else 226273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,d) 22730de9b25SBarry Smith #endif 22830de9b25SBarry Smith 229014dd563SJed Brown PETSC_EXTERN PetscBool MatRegisterAllCalled; 230014dd563SJed Brown PETSC_EXTERN PetscFList MatList; 231014dd563SJed Brown PETSC_EXTERN PetscFList MatColoringList; 232014dd563SJed Brown PETSC_EXTERN PetscFList MatPartitioningList; 233014dd563SJed Brown PETSC_EXTERN PetscFList MatCoarsenList; 23428988994SBarry Smith 2353b224e63SBarry Smith /*E 2363b224e63SBarry Smith MatStructure - Indicates if the matrix has the same nonzero structure 2373b224e63SBarry Smith 2383b224e63SBarry Smith Level: beginner 2393b224e63SBarry Smith 2403b224e63SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 2413b224e63SBarry Smith 2423b224e63SBarry Smith .seealso: MatCopy(), KSPSetOperators(), PCSetOperators() 2433b224e63SBarry Smith E*/ 244214bc6a2SJed Brown typedef enum {DIFFERENT_NONZERO_PATTERN,SUBSET_NONZERO_PATTERN,SAME_NONZERO_PATTERN,SAME_PRECONDITIONER} MatStructure; 2453b224e63SBarry Smith 246014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqDense(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*); 247014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateDense(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*); 248014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 249014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 250014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *); 251014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithSplitArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],PetscInt[],PetscInt[],PetscScalar[],Mat*); 2528d7a6e47SBarry Smith 253014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 254014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 255014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat*); 256d21a29f3SJed Brown 257014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAdj(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscInt[],Mat*); 258014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 259d21a29f3SJed Brown 260014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 261014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPISBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *); 262014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 263014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatXAIJSetPreallocation(Mat,PetscInt,const PetscInt*,const PetscInt*,const PetscInt*,const PetscInt*); 264dfb205c3SBarry Smith 265014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateShell(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,void *,Mat*); 266014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateNormal(Mat,Mat*); 267014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateLRC(Mat,Mat,Mat,Mat*); 268014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateIS(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,ISLocalToGlobalMapping,Mat*); 269014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 270014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 271c0cdd4a1SDahai Guo 272014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 273014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 274014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqSBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 275014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPISBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 276c0cdd4a1SDahai Guo 277014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateScatter(MPI_Comm,VecScatter,Mat*); 278014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScatterSetVecScatter(Mat,VecScatter); 279014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScatterGetVecScatter(Mat,VecScatter*); 280014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateBlockMat(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,Mat*); 281014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeAddMat(Mat,Mat); 282014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeMerge(Mat); 283014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateComposite(MPI_Comm,PetscInt,const Mat*,Mat*); 2846d7c1e57SBarry Smith typedef enum {MAT_COMPOSITE_ADDITIVE,MAT_COMPOSITE_MULTIPLICATIVE} MatCompositeType; 285014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeSetType(Mat,MatCompositeType); 2866d7c1e57SBarry Smith 28719fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatCreateFFT(MPI_Comm,PetscInt,const PetscInt[],MatType,Mat*); 288014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqCUFFT(MPI_Comm,PetscInt,const PetscInt[],Mat*); 289dedccee8SHong Zhang 290014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateTranspose(Mat,Mat*); 291014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSubMatrix(Mat,IS,IS,Mat*); 292014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSubMatrixUpdate(Mat,Mat,IS,IS); 293014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateLocalRef(Mat,IS,IS,Mat*); 2941472f72bSBarry Smith 295014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPythonSetType(Mat,const char[]); 2961d6018f0SLisandro Dalcin 297014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetUp(Mat); 298014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroy(Mat*); 29921c89e3eSBarry Smith 300014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatConjugate(Mat); 301014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRealPart(Mat); 302014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatImaginaryPart(Mat); 303014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetDiagonalBlock(Mat,Mat*); 304014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetTrace(Mat,PetscScalar*); 305713ccfa9SJed Brown PETSC_EXTERN PetscErrorCode MatInvertBlockDiagonal(Mat,const PetscScalar **); 30699cafbc1SBarry Smith 3078ed539a5SBarry Smith /* ------------------------------------------------------------*/ 308014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 309014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 310014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesRow(Mat,PetscInt,const PetscScalar[]); 311014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesRowLocal(Mat,PetscInt,const PetscScalar[]); 312014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBatch(Mat,PetscInt,PetscInt,PetscInt[],const PetscScalar[]); 31373a71a0fSBarry Smith PETSC_EXTERN PetscErrorCode MatSetRandom(Mat,PetscRandom); 31484cb2905SBarry Smith 3152ef4de8bSBarry Smith /*S 3162ef4de8bSBarry Smith MatStencil - Data structure (C struct) for storing information about a single row or 317be479b30SJed Brown column of a matrix as indexed on an associated grid. 318be479b30SJed Brown 319be479b30SJed Brown Fortran usage is different, see MatSetValuesStencil() for details. 3202ef4de8bSBarry Smith 3212ef4de8bSBarry Smith Level: beginner 3222ef4de8bSBarry Smith 3232ef4de8bSBarry Smith Concepts: matrix; linear operator 3242ef4de8bSBarry Smith 325be479b30SJed Brown .seealso: MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockedStencil() 3262ef4de8bSBarry Smith S*/ 327435da068SBarry Smith typedef struct { 328c1ac3661SBarry Smith PetscInt k,j,i,c; 329435da068SBarry Smith } MatStencil; 3302ef4de8bSBarry Smith 331014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode); 332014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode); 333014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt); 334435da068SBarry Smith 335014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetColoring(Mat,ISColoring); 336014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesAdifor(Mat,PetscInt,void*); 3373a7fca6bSBarry Smith 338d91e6319SBarry Smith /*E 339d91e6319SBarry Smith MatAssemblyType - Indicates if the matrix is now to be used, or if you plan 340d91e6319SBarry Smith to continue to add values to it 341d91e6319SBarry Smith 342d91e6319SBarry Smith Level: beginner 343d91e6319SBarry Smith 344d91e6319SBarry Smith .seealso: MatAssemblyBegin(), MatAssemblyEnd() 345d91e6319SBarry Smith E*/ 3466d4a8577SBarry Smith typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType; 347014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssemblyBegin(Mat,MatAssemblyType); 348014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssemblyEnd(Mat,MatAssemblyType); 349014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssembled(Mat,PetscBool *); 3504f9c727eSBarry Smith 3511d73ed98SBarry Smith 35230de9b25SBarry Smith 353d91e6319SBarry Smith /*E 354d91e6319SBarry Smith MatOption - Options that may be set for a matrix and its behavior or storage 355d91e6319SBarry Smith 356d91e6319SBarry Smith Level: beginner 357d91e6319SBarry Smith 3580a835dfdSSatish Balay Any additions/changes here MUST also be made in include/finclude/petscmat.h 359d91e6319SBarry Smith 360382164b0SBarry Smith Developer Notes: Entries that are negative need not be called collectively by all processes. 361382164b0SBarry Smith 362d91e6319SBarry Smith .seealso: MatSetOption() 363d91e6319SBarry Smith E*/ 364382164b0SBarry Smith typedef enum {MAT_OPTION_MIN = -8, 365382164b0SBarry Smith MAT_NEW_NONZERO_LOCATION_ERR = -7, 366382164b0SBarry Smith MAT_NO_OFF_PROC_ZERO_ROWS = -6, 367382164b0SBarry Smith MAT_NO_OFF_PROC_ENTRIES = -5, 368382164b0SBarry Smith MAT_UNUSED_NONZERO_LOCATION_ERR = -4, 369382164b0SBarry Smith MAT_NEW_NONZERO_ALLOCATION_ERR = -3, 370382164b0SBarry Smith MAT_ROW_ORIENTED = -2, 371382164b0SBarry Smith MAT_NEW_NONZERO_LOCATIONS = -1, 372382164b0SBarry Smith MAT_SYMMETRIC = 1, 373382164b0SBarry Smith MAT_STRUCTURALLY_SYMMETRIC = 2, 374382164b0SBarry Smith MAT_NEW_DIAGONALS = 3, 375382164b0SBarry Smith MAT_IGNORE_OFF_PROC_ENTRIES = 4, 376382164b0SBarry Smith MAT_USE_HASH_TABLE = 5, 377382164b0SBarry Smith MAT_KEEP_NONZERO_PATTERN = 6, 378382164b0SBarry Smith MAT_IGNORE_ZERO_ENTRIES = 7, 379382164b0SBarry Smith MAT_USE_INODES = 8, 380382164b0SBarry Smith MAT_HERMITIAN = 9, 381382164b0SBarry Smith MAT_SYMMETRY_ETERNAL = 10, 382382164b0SBarry Smith MAT_CHECK_COMPRESSED_ROW = 11, 383382164b0SBarry Smith MAT_IGNORE_LOWER_TRIANGULAR = 12, 384382164b0SBarry Smith MAT_ERROR_LOWER_TRIANGULAR = 13, 385382164b0SBarry Smith MAT_GETROW_UPPERTRIANGULAR = 14, 386382164b0SBarry Smith MAT_SPD = 15, 387382164b0SBarry Smith MAT_OPTION_MAX = 16} MatOption; 388e057df02SPaul Mullowney 389014dd563SJed Brown PETSC_EXTERN const char *MatOptions[]; 390014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetOption(Mat,MatOption,PetscBool ); 39119fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatGetType(Mat,MatType*); 39284cb2905SBarry Smith 393014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]); 394014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 395014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 396014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowUpperTriangular(Mat); 397014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRowUpperTriangular(Mat); 398014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 399014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 400014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnVector(Mat,Vec,PetscInt); 4018c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJGetArray(Mat,PetscScalar *[]); 4028c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJRestoreArray(Mat,PetscScalar *[]); 4038c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatDenseGetArray(Mat,PetscScalar *[]); 4048c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatDenseRestoreArray(Mat,PetscScalar *[]); 405014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBlockSize(Mat,PetscInt *); 406014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSize(Mat,PetscInt); 407014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBlockSizes(Mat,PetscInt *,PetscInt *); 408014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSizes(Mat,PetscInt,PetscInt); 409014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNThreads(Mat,PetscInt); 410014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNThreads(Mat,PetscInt*); 4111620fd73SBarry Smith 412014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMult(Mat,Vec,Vec); 413014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultDiagonalBlock(Mat,Vec,Vec); 414014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultAdd(Mat,Vec,Vec,Vec); 415014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTranspose(Mat,Vec,Vec); 416014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultHermitianTranspose(Mat,Vec,Vec); 417014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsTranspose(Mat,Mat,PetscReal,PetscBool *); 418014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitianTranspose(Mat,Mat,PetscReal,PetscBool *); 419014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeAdd(Mat,Vec,Vec,Vec); 420014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultHermitianTransposeAdd(Mat,Vec,Vec,Vec); 421014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultConstrained(Mat,Vec,Vec); 422014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeConstrained(Mat,Vec,Vec); 423014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatSolve(Mat,Mat,Mat); 424f5edf698SHong Zhang 425d91e6319SBarry Smith /*E 426d91e6319SBarry Smith MatDuplicateOption - Indicates if a duplicated sparse matrix should have 427d91e6319SBarry Smith its numerical values copied over or just its nonzero structure. 428d91e6319SBarry Smith 429d91e6319SBarry Smith Level: beginner 430d91e6319SBarry Smith 431d91e6319SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 432d91e6319SBarry Smith 43370dcbbb9SBarry Smith $ MAT_SHARE_NONZERO_PATTERN - the i and j arrays in the new matrix will be shared with the original matrix 43470dcbbb9SBarry Smith $ this also triggers the MAT_DO_NOT_COPY_VALUES option. This is used when you 43570dcbbb9SBarry Smith $ have several matrices with the same nonzero pattern. 43670dcbbb9SBarry Smith 437d91e6319SBarry Smith .seealso: MatDuplicate() 438d91e6319SBarry Smith E*/ 43970dcbbb9SBarry Smith typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES,MAT_SHARE_NONZERO_PATTERN} MatDuplicateOption; 4402e8a6d31SBarry Smith 44119fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatConvert(Mat,MatType,MatReuse,Mat*); 442014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDuplicate(Mat,MatDuplicateOption,Mat*); 44394a9d846SBarry Smith 444cb5b572fSBarry Smith 445014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCopy(Mat,Mat,MatStructure); 446014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatView(Mat,PetscViewer); 447014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsSymmetric(Mat,PetscReal,PetscBool *); 448014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsStructurallySymmetric(Mat,PetscBool *); 449014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitian(Mat,PetscReal,PetscBool *); 450014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsSymmetricKnown(Mat,PetscBool *,PetscBool *); 451014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitianKnown(Mat,PetscBool *,PetscBool *); 452014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMissingDiagonal(Mat,PetscBool *,PetscInt *); 453014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLoad(Mat, PetscViewer); 4547b80b807SBarry Smith 4551a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool *); 4561a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool *); 4571a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool *); 4581a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool *); 459d4fbbf0eSBarry Smith 460d91e6319SBarry Smith /*S 461d91e6319SBarry Smith MatInfo - Context of matrix information, used with MatGetInfo() 462d91e6319SBarry Smith 463d91e6319SBarry Smith In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE 464d91e6319SBarry Smith 465d91e6319SBarry Smith Level: intermediate 466d91e6319SBarry Smith 467d91e6319SBarry Smith Concepts: matrix^nonzero information 468d91e6319SBarry Smith 469d9274352SBarry Smith .seealso: MatGetInfo(), MatInfoType 470d91e6319SBarry Smith S*/ 4714e220ebcSLois Curfman McInnes typedef struct { 472b0a32e0cSBarry Smith PetscLogDouble block_size; /* block size */ 473b0a32e0cSBarry Smith PetscLogDouble nz_allocated,nz_used,nz_unneeded; /* number of nonzeros */ 474b0a32e0cSBarry Smith PetscLogDouble memory; /* memory allocated */ 475b0a32e0cSBarry Smith PetscLogDouble assemblies; /* number of matrix assemblies called */ 476b0a32e0cSBarry Smith PetscLogDouble mallocs; /* number of mallocs during MatSetValues() */ 477b0a32e0cSBarry Smith PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */ 478b0a32e0cSBarry Smith PetscLogDouble factor_mallocs; /* number of mallocs during factorization */ 4794e220ebcSLois Curfman McInnes } MatInfo; 4804e220ebcSLois Curfman McInnes 481d9274352SBarry Smith /*E 482d9274352SBarry Smith MatInfoType - Indicates if you want information about the local part of the matrix, 483d9274352SBarry Smith the entire parallel matrix or the maximum over all the local parts. 484d9274352SBarry Smith 485d9274352SBarry Smith Level: beginner 486d9274352SBarry Smith 487d9274352SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 488d9274352SBarry Smith 489d9274352SBarry Smith .seealso: MatGetInfo(), MatInfo 490d9274352SBarry Smith E*/ 4917b80b807SBarry Smith typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType; 492014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetInfo(Mat,MatInfoType,MatInfo*); 493014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetDiagonal(Mat,Vec); 494014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMax(Mat,Vec,PetscInt[]); 495014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMin(Mat,Vec,PetscInt[]); 496014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMaxAbs(Mat,Vec,PetscInt[]); 497014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMinAbs(Mat,Vec,PetscInt[]); 498014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowSum(Mat,Vec); 499014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTranspose(Mat,MatReuse,Mat*); 500014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatHermitianTranspose(Mat,MatReuse,Mat*); 501014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPermute(Mat,IS,IS,Mat *); 502014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalScale(Mat,Vec,Vec); 503014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalSet(Mat,Vec,InsertMode); 504014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatEqual(Mat,Mat,PetscBool *); 505014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultEqual(Mat,Mat,PetscInt,PetscBool *); 506014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultAddEqual(Mat,Mat,PetscInt,PetscBool *); 507014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeEqual(Mat,Mat,PetscInt,PetscBool *); 508014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscBool *); 5097b80b807SBarry Smith 510014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNorm(Mat,NormType,PetscReal *); 511014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnNorms(Mat,NormType,PetscReal *); 512014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroEntries(Mat); 513014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 514014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsIS(Mat,IS,PetscScalar,Vec,Vec); 515014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsStencil(Mat,PetscInt,const MatStencil [],PetscScalar,Vec,Vec); 516014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsStencil(Mat,PetscInt,const MatStencil[],PetscScalar,Vec,Vec); 517014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumns(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 518014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsIS(Mat,IS,PetscScalar,Vec,Vec); 5197b80b807SBarry Smith 520014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSize(Mat,PetscInt*,PetscInt*); 521014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalSize(Mat,PetscInt*,PetscInt*); 522014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRange(Mat,PetscInt*,PetscInt*); 523014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRanges(Mat,const PetscInt**); 524014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRangeColumn(Mat,PetscInt*,PetscInt*); 525014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRangesColumn(Mat,const PetscInt**); 526566876eaSJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipIS(Mat,IS*,IS*); 5275ef9f2a5SBarry Smith 528014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]); 529014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSubMatricesParallel(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]); 530014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroyMatrices(PetscInt,Mat *[]); 531014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSubMatrix(Mat,IS,IS,MatReuse,Mat *); 532014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalSubMatrix(Mat,IS,IS,Mat*); 533014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreLocalSubMatrix(Mat,IS,IS,Mat*); 534014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSeqNonzeroStructure(Mat,Mat*); 535014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroySeqNonzeroStructure(Mat*); 5367b80b807SBarry Smith 537014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJConcatenateSeqAIJ(MPI_Comm,Mat,PetscInt,MatReuse,Mat*); 538014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJConcatenateSeqAIJSymbolic(MPI_Comm,Mat,PetscInt,Mat*); 539014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJConcatenateSeqAIJNumeric(MPI_Comm,Mat,PetscInt,Mat); 540014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJ(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*); 541014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJSymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*); 542014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJNumeric(Mat,Mat); 543014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMat(Mat,MatReuse,Mat*); 544014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*); 545014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,Mat*); 54643eb5e2fSMatthew Knepley #if defined (PETSC_USE_CTABLE) 547014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetCommunicationStructs(Mat, Vec *, PetscTable *, VecScatter *); 54843eb5e2fSMatthew Knepley #else 549014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetCommunicationStructs(Mat, Vec *, PetscInt *[], VecScatter *); 55043eb5e2fSMatthew Knepley #endif 551014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetGhosts(Mat, PetscInt *,const PetscInt *[]); 5528efafbd8SBarry Smith 553014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt); 5547b80b807SBarry Smith 555014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*); 556014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMultSymbolic(Mat,Mat,PetscReal,Mat*); 557014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMultNumeric(Mat,Mat,Mat); 55822440eb1SKris Buschelman 5597bab7c10SHong Zhang PETSC_EXTERN PetscErrorCode MatMatMatMult(Mat,Mat,Mat,MatReuse,PetscReal,Mat*); 5607bab7c10SHong Zhang PETSC_EXTERN PetscErrorCode MatMatMatMultSymbolic(Mat,Mat,Mat,PetscReal,Mat*); 5617bab7c10SHong Zhang PETSC_EXTERN PetscErrorCode MatMatMatMultNumeric(Mat,Mat,Mat,Mat); 5627bab7c10SHong Zhang 563014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*); 564014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAPSymbolic(Mat,Mat,PetscReal,Mat*); 565014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAPNumeric(Mat,Mat,Mat); 566014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARt(Mat,Mat,MatReuse,PetscReal,Mat*); 567014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARtSymbolic(Mat,Mat,PetscReal,Mat*); 568014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARtNumeric(Mat,Mat,Mat); 56922440eb1SKris Buschelman 570014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeMatMult(Mat,Mat,MatReuse,PetscReal,Mat*); 571014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposetMatMultSymbolic(Mat,Mat,PetscReal,Mat*); 572014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposetMatMultNumeric(Mat,Mat,Mat); 573014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMult(Mat,Mat,MatReuse,PetscReal,Mat*); 574014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMultSymbolic(Mat,Mat,PetscReal,Mat*); 575014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMultNumeric(Mat,Mat,Mat); 576bc011b1eSHong Zhang 577014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAXPY(Mat,PetscScalar,Mat,MatStructure); 578014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAYPX(Mat,PetscScalar,Mat,MatStructure); 5791c741599SHong Zhang 580014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScale(Mat,PetscScalar); 581014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShift(Mat,PetscScalar); 5827b80b807SBarry Smith 583014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping); 584014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetLocalToGlobalMappingBlock(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping); 585014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*); 586014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalToGlobalMappingBlock(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*); 587014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 588014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsLocalIS(Mat,IS,PetscScalar,Vec,Vec); 589014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 590014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocalIS(Mat,IS,PetscScalar,Vec,Vec); 591014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 592014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 593052efed2SBarry Smith 594014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStashSetInitialSize(Mat,PetscInt,PetscInt); 595014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*); 59690f02eecSBarry Smith 597014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInterpolate(Mat,Vec,Vec); 598014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInterpolateAdd(Mat,Vec,Vec,Vec); 599014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestrict(Mat,Vec,Vec); 600014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetVecs(Mat,Vec*,Vec*); 601014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRedundantMatrix(Mat,PetscInt,MPI_Comm,PetscInt,MatReuse,Mat*); 602014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetMultiProcBlock(Mat,MPI_Comm,MatReuse,Mat*); 603014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFindZeroDiagonals(Mat,IS*); 604bd481603SBarry Smith 605bd481603SBarry Smith /*MC 6061620fd73SBarry Smith MatSetValue - Set a single entry into a matrix. 6071620fd73SBarry Smith 6081620fd73SBarry Smith Not collective 6091620fd73SBarry Smith 6101620fd73SBarry Smith Input Parameters: 6111620fd73SBarry Smith + m - the matrix 6121620fd73SBarry Smith . row - the row location of the entry 6131620fd73SBarry Smith . col - the column location of the entry 6141620fd73SBarry Smith . value - the value to insert 6151620fd73SBarry Smith - mode - either INSERT_VALUES or ADD_VALUES 6161620fd73SBarry Smith 6171620fd73SBarry Smith Notes: 6181620fd73SBarry Smith For efficiency one should use MatSetValues() and set several or many 6191620fd73SBarry Smith values simultaneously if possible. 6201620fd73SBarry Smith 6211620fd73SBarry Smith Level: beginner 6221620fd73SBarry Smith 6231620fd73SBarry Smith .seealso: MatSetValues(), MatSetValueLocal() 6241620fd73SBarry Smith M*/ 6251620fd73SBarry 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);} 6261620fd73SBarry Smith 6271620fd73SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatGetValue(Mat v,PetscInt i,PetscInt j,PetscScalar *va) {return MatGetValues(v,1,&i,1,&j,va);} 6281620fd73SBarry Smith 6291620fd73SBarry 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);} 6301620fd73SBarry Smith 6311620fd73SBarry Smith /*MC 6320d2a7ff2SSatish Balay MatPreallocateInitialize - Begins the block of code that will count the number of nonzeros per 633bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 634bd481603SBarry Smith 635bd481603SBarry Smith Synopsis: 636c1ac3661SBarry Smith PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz) 637bd481603SBarry Smith 638bd481603SBarry Smith Collective on MPI_Comm 639bd481603SBarry Smith 640bd481603SBarry Smith Input Parameters: 641bd481603SBarry Smith + comm - the communicator that will share the eventually allocated matrix 642859fcb39SBarry Smith . nrows - the number of LOCAL rows in the matrix 643859fcb39SBarry Smith - ncols - the number of LOCAL columns in the matrix 644bd481603SBarry Smith 645bd481603SBarry Smith Output Parameters: 646bd481603SBarry Smith + dnz - the array that will be passed to the matrix preallocation routines 647bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 648bd481603SBarry Smith 649bd481603SBarry Smith 650bd481603SBarry Smith Level: intermediate 651bd481603SBarry Smith 652bd481603SBarry Smith Notes: 6530598bfebSBarry Smith See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details. 654bd481603SBarry Smith 6551d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 656bd481603SBarry Smith 657bd481603SBarry Smith Use MatPreallocateInitializeSymmetric() for symmetric matrices (MPISBAIJ matrices) 658bd481603SBarry Smith 6591620fd73SBarry Smith This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize(). 6601620fd73SBarry Smith 661bd481603SBarry Smith Concepts: preallocation^Matrix 662bd481603SBarry Smith 663bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(), 664bd481603SBarry Smith MatPreallocateInitializeSymmetric(), MatPreallocateSymmetricSetLocal() 665bd481603SBarry Smith M*/ 666c4f061fbSSatish Balay #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \ 6677c922b88SBarry Smith { \ 668a89bc211SBarry Smith PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__start,__end; \ 669a89bc211SBarry Smith _4_ierr = PetscMalloc2(__nrows,PetscInt,&dnz,__nrows,PetscInt,&onz);CHKERRQ(_4_ierr); \ 670a89bc211SBarry Smith _4_ierr = PetscMemzero(dnz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\ 671eabe889fSLisandro Dalcin _4_ierr = PetscMemzero(onz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr); __start = 0; __end = __start; \ 672c1ac3661SBarry Smith _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ctmp;\ 673a89bc211SBarry Smith _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows; 6747c922b88SBarry Smith 675bd481603SBarry Smith /*MC 676bd481603SBarry Smith MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be 677bd481603SBarry Smith inserted using a local number of the rows and columns 678bd481603SBarry Smith 679bd481603SBarry Smith Synopsis: 680c1ac3661SBarry Smith PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 681bd481603SBarry Smith 682bd481603SBarry Smith Not Collective 683bd481603SBarry Smith 684bd481603SBarry Smith Input Parameters: 685784ac674SJed Brown + map - the row mapping from local numbering to global numbering 686bd481603SBarry Smith . nrows - the number of rows indicated 6871d73ed98SBarry Smith . rows - the indices of the rows 688784ac674SJed Brown . cmap - the column mapping from local to global numbering 689bd481603SBarry Smith . ncols - the number of columns in the matrix 690bd481603SBarry Smith . cols - the columns indicated 691bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 692bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 693bd481603SBarry Smith 694bd481603SBarry Smith 695bd481603SBarry Smith Level: intermediate 696bd481603SBarry Smith 697bd481603SBarry Smith Notes: 6980598bfebSBarry Smith See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details. 699bd481603SBarry Smith 7001d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 701bd481603SBarry Smith 702bd481603SBarry Smith Concepts: preallocation^Matrix 703bd481603SBarry Smith 704bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(), 705bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal() 706bd481603SBarry Smith M*/ 707784ac674SJed Brown #define MatPreallocateSetLocal(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \ 708c4f061fbSSatish Balay {\ 709c1ac3661SBarry Smith PetscInt __l;\ 710784ac674SJed Brown _4_ierr = ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\ 711784ac674SJed Brown _4_ierr = ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\ 712c4f061fbSSatish Balay for (__l=0;__l<nrows;__l++) {\ 713ef66eb69SBarry Smith _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 714c4f061fbSSatish Balay }\ 715c4f061fbSSatish Balay } 716c4f061fbSSatish Balay 717bd481603SBarry Smith /*MC 718bd481603SBarry Smith MatPreallocateSymmetricSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be 719bd481603SBarry Smith inserted using a local number of the rows and columns 720bd481603SBarry Smith 721bd481603SBarry Smith Synopsis: 722c1ac3661SBarry Smith PetscErrorCode MatPreallocateSymmetricSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 723bd481603SBarry Smith 724bd481603SBarry Smith Not Collective 725bd481603SBarry Smith 726bd481603SBarry Smith Input Parameters: 727bd481603SBarry Smith + map - the mapping between local numbering and global numbering 728bd481603SBarry Smith . nrows - the number of rows indicated 7291d73ed98SBarry Smith . rows - the indices of the rows 730bd481603SBarry Smith . ncols - the number of columns in the matrix 731bd481603SBarry Smith . cols - the columns indicated 732bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 733bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 734bd481603SBarry Smith 735bd481603SBarry Smith 736bd481603SBarry Smith Level: intermediate 737bd481603SBarry Smith 738bd481603SBarry Smith Notes: 7390598bfebSBarry Smith See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details. 740bd481603SBarry Smith 741bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 742bd481603SBarry Smith 743bd481603SBarry Smith Concepts: preallocation^Matrix 744bd481603SBarry Smith 745bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(), 746bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal() 747bd481603SBarry Smith M*/ 748d3d32019SBarry Smith #define MatPreallocateSymmetricSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\ 749d3d32019SBarry Smith {\ 750c1ac3661SBarry Smith PetscInt __l;\ 751d3d32019SBarry Smith _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\ 752d3d32019SBarry Smith _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\ 753d3d32019SBarry Smith for (__l=0;__l<nrows;__l++) {\ 754d3d32019SBarry Smith _4_ierr = MatPreallocateSymmetricSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 755d3d32019SBarry Smith }\ 756d3d32019SBarry Smith } 757d3d32019SBarry Smith 758bd481603SBarry Smith /*MC 759bd481603SBarry Smith MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be 760bd481603SBarry Smith inserted using a local number of the rows and columns 761bd481603SBarry Smith 762bd481603SBarry Smith Synopsis: 763c1ac3661SBarry Smith PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 764bd481603SBarry Smith 765bd481603SBarry Smith Not Collective 766bd481603SBarry Smith 767bd481603SBarry Smith Input Parameters: 76864075487SBarry Smith + row - the row 769bd481603SBarry Smith . ncols - the number of columns in the matrix 770a50f8bf6SHong Zhang - cols - the columns indicated 771a50f8bf6SHong Zhang 772a50f8bf6SHong Zhang Output Parameters: 773a50f8bf6SHong Zhang + dnz - the array that will be passed to the matrix preallocation routines 774bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 775bd481603SBarry Smith 776bd481603SBarry Smith 777bd481603SBarry Smith Level: intermediate 778bd481603SBarry Smith 779bd481603SBarry Smith Notes: 7800598bfebSBarry Smith See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details. 781bd481603SBarry Smith 782bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 783bd481603SBarry Smith 7841620fd73SBarry Smith This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize(). 7851620fd73SBarry Smith 786bd481603SBarry Smith Concepts: preallocation^Matrix 787bd481603SBarry Smith 788bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(), 789bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal() 790bd481603SBarry Smith M*/ 791c4f061fbSSatish Balay #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\ 792c1ac3661SBarry Smith { PetscInt __i; \ 793e32f2f54SBarry 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);\ 794e32f2f54SBarry 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);\ 7957c922b88SBarry Smith for (__i=0; __i<nc; __i++) {\ 79664075487SBarry Smith if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \ 7977cc688d7SBarry Smith else dnz[row - __rstart]++;\ 7987c922b88SBarry Smith }\ 7997c922b88SBarry Smith } 8007c922b88SBarry Smith 801bd481603SBarry Smith /*MC 802bd481603SBarry Smith MatPreallocateSymmetricSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be 803bd481603SBarry Smith inserted using a local number of the rows and columns 804bd481603SBarry Smith 805bd481603SBarry Smith Synopsis: 806c1ac3661SBarry Smith PetscErrorCode MatPreallocateSymmetricSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 807bd481603SBarry Smith 808bd481603SBarry Smith Not Collective 809bd481603SBarry Smith 810bd481603SBarry Smith Input Parameters: 811bd481603SBarry Smith + nrows - the number of rows indicated 8121d73ed98SBarry Smith . rows - the indices of the rows 813bd481603SBarry Smith . ncols - the number of columns in the matrix 814bd481603SBarry Smith . cols - the columns indicated 815bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 816bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 817bd481603SBarry Smith 818bd481603SBarry Smith 819bd481603SBarry Smith Level: intermediate 820bd481603SBarry Smith 821bd481603SBarry Smith Notes: 8220598bfebSBarry Smith See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details. 823bd481603SBarry Smith 824bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 825bd481603SBarry Smith 8261620fd73SBarry Smith This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize(). 8271620fd73SBarry Smith 828bd481603SBarry Smith Concepts: preallocation^Matrix 829bd481603SBarry Smith 830bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(), 831bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal() 832bd481603SBarry Smith M*/ 833d3d32019SBarry Smith #define MatPreallocateSymmetricSet(row,nc,cols,dnz,onz) 0;\ 834c1ac3661SBarry Smith { PetscInt __i; \ 835d3d32019SBarry Smith for (__i=0; __i<nc; __i++) {\ 836d3d32019SBarry Smith if (cols[__i] >= __end) onz[row - __rstart]++; \ 837d3d32019SBarry Smith else if (cols[__i] >= row) dnz[row - __rstart]++;\ 838d3d32019SBarry Smith }\ 839d3d32019SBarry Smith } 840d3d32019SBarry Smith 841bd481603SBarry Smith /*MC 84216371a99SBarry Smith MatPreallocateLocation - An alternative to MatPreallocationSet() that puts the nonzero locations into the matrix if it exists 84316371a99SBarry Smith 84416371a99SBarry Smith Synopsis: 84516371a99SBarry Smith PetscErrorCode MatPreallocateLocations(Mat A,PetscInt row,PetscInt ncols,PetscInt *cols,PetscInt *dnz,PetscInt *onz) 84616371a99SBarry Smith 84716371a99SBarry Smith Not Collective 84816371a99SBarry Smith 84916371a99SBarry Smith Input Parameters: 85016371a99SBarry Smith . A - matrix 85116371a99SBarry Smith . row - row where values exist (must be local to this process) 85216371a99SBarry Smith . ncols - number of columns 85316371a99SBarry Smith . cols - columns with nonzeros 85416371a99SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 85516371a99SBarry Smith - ozn - the other array passed to the matrix preallocation routines 85616371a99SBarry Smith 85716371a99SBarry Smith 85816371a99SBarry Smith Level: intermediate 85916371a99SBarry Smith 86016371a99SBarry Smith Notes: 8610598bfebSBarry Smith See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details. 86216371a99SBarry Smith 86316371a99SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 86416371a99SBarry Smith 86516371a99SBarry Smith This is a MACRO not a function because it uses a bunch of variables private to the MatPreallocation.... routines. 86616371a99SBarry Smith 86716371a99SBarry Smith Concepts: preallocation^Matrix 86816371a99SBarry Smith 86916371a99SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(), 870eabe889fSLisandro Dalcin MatPreallocateSymmetricSetLocal() 87116371a99SBarry Smith M*/ 87216371a99SBarry 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);} 87316371a99SBarry Smith 87416371a99SBarry Smith 87516371a99SBarry Smith /*MC 8760d2a7ff2SSatish Balay MatPreallocateFinalize - Ends the block of code that will count the number of nonzeros per 877bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 878bd481603SBarry Smith 879bd481603SBarry Smith Synopsis: 880c1ac3661SBarry Smith PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz) 881bd481603SBarry Smith 882bd481603SBarry Smith Collective on MPI_Comm 883bd481603SBarry Smith 884bd481603SBarry Smith Input Parameters: 88516371a99SBarry Smith + dnz - the array that was be passed to the matrix preallocation routines 886bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 887bd481603SBarry Smith 888bd481603SBarry Smith 889bd481603SBarry Smith Level: intermediate 890bd481603SBarry Smith 891bd481603SBarry Smith Notes: 8920598bfebSBarry Smith See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details. 893bd481603SBarry Smith 894bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 895bd481603SBarry Smith 8961620fd73SBarry Smith This is a MACRO not a function because it closes the { started in MatPreallocateInitialize(). 8971620fd73SBarry Smith 898bd481603SBarry Smith Concepts: preallocation^Matrix 899bd481603SBarry Smith 900bd481603SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(), 901eabe889fSLisandro Dalcin MatPreallocateSymmetricSetLocal() 902bd481603SBarry Smith M*/ 903a89bc211SBarry Smith #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree2(dnz,onz);CHKERRQ(_4_ierr);} 9047c922b88SBarry Smith 905bd481603SBarry Smith 906bd481603SBarry Smith 9077b80b807SBarry Smith /* Routines unique to particular data structures */ 908014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellGetContext(Mat,void *); 909698d4c6aSKris Buschelman 910014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInodeAdjustForInodes(Mat,IS*,IS*); 911014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *); 912698d4c6aSKris Buschelman 913014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetColumnIndices(Mat,PetscInt[]); 914014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetColumnIndices(Mat,PetscInt[]); 915014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 916014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 917014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 918014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJFromTriple(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*,PetscInt,PetscBool); 9197b80b807SBarry Smith 920a96a251dSBarry Smith #define MAT_SKIP_ALLOCATION -4 921a96a251dSBarry Smith 922014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]); 923014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]); 924014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]); 925273d9f13SBarry Smith 926014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 927014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 928014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 929014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []); 930014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 931014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]); 932014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 933014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]); 934014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIDenseSetPreallocation(Mat,PetscScalar[]); 935014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqDenseSetPreallocation(Mat,PetscScalar[]); 9369230625dSJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,const PetscInt*[]); 9379230625dSJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,const PetscInt*[]); 938014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAdjCreateNonemptySubcommMat(Mat,Mat*); 939273d9f13SBarry Smith 940014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqDenseSetLDA(Mat,PetscInt); 941014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDenseGetLocalMatrix(Mat,Mat*); 9421b807ce4Svictorle 943014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStoreValues(Mat); 944014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRetrieveValues(Mat); 9452e8a6d31SBarry Smith 946014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDAADSetCtx(Mat,void*); 9473a7fca6bSBarry Smith 948014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFindNonzeroRows(Mat,IS*); 9497b80b807SBarry Smith /* 9507b80b807SBarry Smith These routines are not usually accessed directly, rather solving is 95194b7f48cSBarry Smith done through the KSP and PC interfaces. 9527b80b807SBarry Smith */ 9537b80b807SBarry Smith 95476bdecfbSBarry Smith /*J 955d9274352SBarry Smith MatOrderingType - String with the name of a PETSc matrix ordering or the creation function 956d9274352SBarry Smith with an optional dynamic library name, for example 957d9274352SBarry Smith http://www.mcs.anl.gov/petsc/lib.a:orderingcreate() 958d9274352SBarry Smith 959d9274352SBarry Smith Level: beginner 960d9274352SBarry Smith 961e5a9bf91SBarry Smith Cannot use const because the PC objects manipulate the string 962e5a9bf91SBarry Smith 963d9274352SBarry Smith .seealso: MatGetOrdering() 96476bdecfbSBarry Smith J*/ 96519fd82e9SBarry Smith typedef const char* MatOrderingType; 9662692d6eeSBarry Smith #define MATORDERINGNATURAL "natural" 9672692d6eeSBarry Smith #define MATORDERINGND "nd" 9682692d6eeSBarry Smith #define MATORDERING1WD "1wd" 9692692d6eeSBarry Smith #define MATORDERINGRCM "rcm" 9702692d6eeSBarry Smith #define MATORDERINGQMD "qmd" 9712692d6eeSBarry Smith #define MATORDERINGROWLENGTH "rowlength" 972312e372aSBarry Smith #define MATORDERINGAMD "amd" /* only works if UMFPACK is installed with PETSc */ 973b12f92e5SBarry Smith 97419fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatGetOrdering(Mat,MatOrderingType,IS*,IS*); 975014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOrderingList(PetscFList *list); 97619fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatOrderingRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,MatOrderingType,IS*,IS*)); 97730de9b25SBarry Smith 97830de9b25SBarry Smith /*MC 9791890ba74SBarry Smith MatOrderingRegisterDynamic - Adds a new sparse matrix ordering to the matrix package. 98030de9b25SBarry Smith 98130de9b25SBarry Smith Synopsis: 9821890ba74SBarry Smith PetscErrorCode MatOrderingRegisterDynamic(const char *name_ordering,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatOrdering)) 98330de9b25SBarry Smith 98430de9b25SBarry Smith Not Collective 98530de9b25SBarry Smith 98630de9b25SBarry Smith Input Parameters: 9872692d6eeSBarry Smith + sname - name of ordering (for example MATORDERINGND) 98830de9b25SBarry Smith . path - location of library where creation routine is 98930de9b25SBarry Smith . name - name of function that creates the ordering type,a string 99030de9b25SBarry Smith - function - function pointer that creates the ordering 99130de9b25SBarry Smith 99230de9b25SBarry Smith Level: developer 99330de9b25SBarry Smith 99430de9b25SBarry Smith If dynamic libraries are used, then the fourth input argument (function) 99530de9b25SBarry Smith is ignored. 99630de9b25SBarry Smith 99730de9b25SBarry Smith Sample usage: 99830de9b25SBarry Smith .vb 99930de9b25SBarry Smith MatOrderingRegisterDynamic("my_order",/home/username/my_lib/lib/libO/solaris/mylib.a, 100030de9b25SBarry Smith "MyOrder",MyOrder); 100130de9b25SBarry Smith .ve 100230de9b25SBarry Smith 100330de9b25SBarry Smith Then, your partitioner can be chosen with the procedural interface via 100430de9b25SBarry Smith $ MatOrderingSetType(part,"my_order) 100530de9b25SBarry Smith or at runtime via the option 10062401956bSBarry Smith $ -pc_factor_mat_ordering_type my_order 100730de9b25SBarry Smith 1008ab901514SBarry Smith ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values. 100930de9b25SBarry Smith 101030de9b25SBarry Smith .keywords: matrix, ordering, register 101130de9b25SBarry Smith 101230de9b25SBarry Smith .seealso: MatOrderingRegisterDestroy(), MatOrderingRegisterAll() 101330de9b25SBarry Smith M*/ 1014aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 1015f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,0) 1016b12f92e5SBarry Smith #else 1017f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,d) 1018b12f92e5SBarry Smith #endif 101930de9b25SBarry Smith 1020014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatOrderingRegisterDestroy(void); 1021014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatOrderingRegisterAll(const char[]); 1022014dd563SJed Brown PETSC_EXTERN PetscBool MatOrderingRegisterAllCalled; 1023014dd563SJed Brown PETSC_EXTERN PetscFList MatOrderingList; 1024d4fbbf0eSBarry Smith 1025014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS); 1026a2ce50c7SBarry Smith 1027d91e6319SBarry Smith /*S 1028d90ac83dSHong Zhang MatFactorShiftType - Numeric Shift. 1029d90ac83dSHong Zhang 1030d90ac83dSHong Zhang Level: beginner 1031d90ac83dSHong Zhang 1032d90ac83dSHong Zhang S*/ 1033d90ac83dSHong Zhang typedef enum {MAT_SHIFT_NONE,MAT_SHIFT_NONZERO,MAT_SHIFT_POSITIVE_DEFINITE,MAT_SHIFT_INBLOCKS} MatFactorShiftType; 10346a6fc655SJed Brown PETSC_EXTERN const char *const MatFactorShiftTypes[]; 1035d90ac83dSHong Zhang 1036d90ac83dSHong Zhang /*S 10372401956bSBarry Smith MatFactorInfo - Data passed into the matrix factorization routines 1038b00f7748SHong Zhang 103961cad860SBarry Smith In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE, that is use 104061cad860SBarry Smith $ MatFactorInfo info(MAT_FACTORINFO_SIZE) 1041b00f7748SHong Zhang 104215e8a5b3SHong Zhang Notes: These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC. 1043b00f7748SHong Zhang 104461cad860SBarry Smith You can use MatFactorInfoInitialize() to set default values. 104561cad860SBarry Smith 1046b00f7748SHong Zhang Level: developer 1047b00f7748SHong Zhang 1048d7d82daaSBarry Smith .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(), 1049d7d82daaSBarry Smith MatFactorInfoInitialize() 1050b00f7748SHong Zhang 1051b00f7748SHong Zhang S*/ 1052b00f7748SHong Zhang typedef struct { 105315e8a5b3SHong Zhang PetscReal diagonal_fill; /* force diagonal to fill in if initially not filled */ 105485317021SBarry Smith PetscReal usedt; 105515e8a5b3SHong Zhang PetscReal dt; /* drop tolerance */ 1056b00f7748SHong Zhang PetscReal dtcol; /* tolerance for pivoting */ 105715e8a5b3SHong Zhang PetscReal dtcount; /* maximum nonzeros to be allowed per row */ 105867e28bfeSBarry Smith PetscReal fill; /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */ 1059348344bbSBarry Smith PetscReal levels; /* ICC/ILU(levels) */ 1060bcd9e38bSBarry Smith PetscReal pivotinblocks; /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0 1061bcd9e38bSBarry Smith factorization may be faster if do not pivot */ 106215e8a5b3SHong Zhang PetscReal zeropivot; /* pivot is called zero if less than this */ 1063f4db908eSBarry Smith PetscReal shifttype; /* type of shift added to matrix factor to prevent zero pivots */ 1064f4db908eSBarry Smith PetscReal shiftamount; /* how large the shift is */ 106515e8a5b3SHong Zhang } MatFactorInfo; 1066ffa6d0a5SLois Curfman McInnes 1067014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFactorInfoInitialize(MatFactorInfo*); 1068014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCholeskyFactor(Mat,IS,const MatFactorInfo*); 1069014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCholeskyFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*); 1070014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCholeskyFactorNumeric(Mat,Mat,const MatFactorInfo*); 1071014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLUFactor(Mat,IS,IS,const MatFactorInfo*); 1072014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatILUFactor(Mat,IS,IS,const MatFactorInfo*); 1073014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*); 1074014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatILUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*); 1075014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatICCFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*); 1076014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatICCFactor(Mat,IS,const MatFactorInfo*); 1077014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLUFactorNumeric(Mat,Mat,const MatFactorInfo*); 1078014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*); 1079014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolve(Mat,Vec,Vec); 1080014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatForwardSolve(Mat,Vec,Vec); 1081014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatBackwardSolve(Mat,Vec,Vec); 1082014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveAdd(Mat,Vec,Vec,Vec); 1083014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveTranspose(Mat,Vec,Vec); 1084014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveTransposeAdd(Mat,Vec,Vec,Vec); 1085014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolves(Mat,Vecs,Vecs); 10868ed539a5SBarry Smith 1087014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetUnfactored(Mat); 1088bb5a7306SBarry Smith 1089d91e6319SBarry Smith /*E 1090d91e6319SBarry Smith MatSORType - What type of (S)SOR to perform 1091bb1eb677SSatish Balay 1092d91e6319SBarry Smith Level: beginner 1093d91e6319SBarry Smith 1094d9274352SBarry Smith May be bitwise ORd together 1095d9274352SBarry Smith 1096d91e6319SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 1097d91e6319SBarry Smith 10984e7234bfSBarry Smith MatSORType may be bitwise ORd together, so do not change the numbers 10994e7234bfSBarry Smith 110041f059aeSBarry Smith .seealso: MatSOR() 1101d91e6319SBarry Smith E*/ 1102ee50ffe9SBarry Smith typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3, 1103ee50ffe9SBarry Smith SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8, 1104ee50ffe9SBarry Smith SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16, 110584cb2905SBarry Smith SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType; 1106014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSOR(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec); 11078ed539a5SBarry Smith 1108d4fbbf0eSBarry Smith /* 1109639f9d9dSBarry Smith These routines are for efficiently computing Jacobians via finite differences. 1110639f9d9dSBarry Smith */ 1111b12f92e5SBarry Smith 111276bdecfbSBarry Smith /*J 1113d9274352SBarry Smith MatColoringType - String with the name of a PETSc matrix coloring or the creation function 1114d9274352SBarry Smith with an optional dynamic library name, for example 1115d9274352SBarry Smith http://www.mcs.anl.gov/petsc/lib.a:coloringcreate() 1116d9274352SBarry Smith 1117d9274352SBarry Smith Level: beginner 1118d9274352SBarry Smith 1119d9274352SBarry Smith .seealso: MatGetColoring() 112076bdecfbSBarry Smith J*/ 112119fd82e9SBarry Smith typedef const char* MatColoringType; 11222692d6eeSBarry Smith #define MATCOLORINGNATURAL "natural" 11232692d6eeSBarry Smith #define MATCOLORINGSL "sl" 11242692d6eeSBarry Smith #define MATCOLORINGLF "lf" 11252692d6eeSBarry Smith #define MATCOLORINGID "id" 1126b12f92e5SBarry Smith 112719fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatGetColoring(Mat,MatColoringType,ISColoring*); 1128014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatColoringRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,MatColoringType,ISColoring *)); 112930de9b25SBarry Smith 113030de9b25SBarry Smith /*MC 113130de9b25SBarry Smith MatColoringRegisterDynamic - Adds a new sparse matrix coloring to the 113230de9b25SBarry Smith matrix package. 113330de9b25SBarry Smith 113430de9b25SBarry Smith Synopsis: 11351890ba74SBarry Smith PetscErrorCode MatColoringRegisterDynamic(const char *name_coloring,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatColoring)) 113630de9b25SBarry Smith 113730de9b25SBarry Smith Not Collective 113830de9b25SBarry Smith 113930de9b25SBarry Smith Input Parameters: 11402692d6eeSBarry Smith + sname - name of Coloring (for example MATCOLORINGSL) 114130de9b25SBarry Smith . path - location of library where creation routine is 114230de9b25SBarry Smith . name - name of function that creates the Coloring type, a string 114330de9b25SBarry Smith - function - function pointer that creates the coloring 114430de9b25SBarry Smith 114530de9b25SBarry Smith Level: developer 114630de9b25SBarry Smith 114730de9b25SBarry Smith If dynamic libraries are used, then the fourth input argument (function) 114830de9b25SBarry Smith is ignored. 114930de9b25SBarry Smith 115030de9b25SBarry Smith Sample usage: 115130de9b25SBarry Smith .vb 115230de9b25SBarry Smith MatColoringRegisterDynamic("my_color",/home/username/my_lib/lib/libO/solaris/mylib.a, 115330de9b25SBarry Smith "MyColor",MyColor); 115430de9b25SBarry Smith .ve 115530de9b25SBarry Smith 115630de9b25SBarry Smith Then, your partitioner can be chosen with the procedural interface via 115730de9b25SBarry Smith $ MatColoringSetType(part,"my_color") 115830de9b25SBarry Smith or at runtime via the option 115930de9b25SBarry Smith $ -mat_coloring_type my_color 116030de9b25SBarry Smith 1161ab901514SBarry Smith $PETSC_ARCH occuring in pathname will be replaced with appropriate values. 116230de9b25SBarry Smith 116330de9b25SBarry Smith .keywords: matrix, Coloring, register 116430de9b25SBarry Smith 116530de9b25SBarry Smith .seealso: MatColoringRegisterDestroy(), MatColoringRegisterAll() 116630de9b25SBarry Smith M*/ 1167aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 1168f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,0) 1169b12f92e5SBarry Smith #else 1170f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,d) 1171b12f92e5SBarry Smith #endif 117230de9b25SBarry Smith 1173014dd563SJed Brown PETSC_EXTERN PetscBool MatColoringRegisterAllCalled; 1174f1144a30SSatish Balay 1175014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatColoringRegisterAll(const char[]); 1176014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatColoringRegisterDestroy(void); 1177014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*); 1178639f9d9dSBarry Smith 1179d9274352SBarry Smith /*S 1180d9274352SBarry Smith MatFDColoring - Object for computing a sparse Jacobian via finite differences 1181d9274352SBarry Smith and coloring 1182639f9d9dSBarry Smith 1183d9274352SBarry Smith Level: beginner 1184d9274352SBarry Smith 1185d9274352SBarry Smith Concepts: coloring, sparse Jacobian, finite differences 1186d9274352SBarry Smith 1187d9274352SBarry Smith .seealso: MatFDColoringCreate() 1188d9274352SBarry Smith S*/ 1189e2a1c21fSSatish Balay typedef struct _p_MatFDColoring* MatFDColoring; 1190639f9d9dSBarry Smith 1191014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringCreate(Mat,ISColoring,MatFDColoring *); 1192014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringDestroy(MatFDColoring*); 1193014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringView(MatFDColoring,PetscViewer); 1194014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*); 1195014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**); 1196014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal); 1197014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetFromOptions(MatFDColoring); 1198014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringApply(Mat,MatFDColoring,Vec,MatStructure*,void *); 1199014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetF(MatFDColoring,Vec); 1200014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,PetscInt*[]); 1201b1683b59SHong Zhang 1202b1683b59SHong Zhang /*S 1203b9af6bddSHong Zhang MatTransposeColoring - Object for computing a sparse matrix product C=A*B^T via coloring 1204b1683b59SHong Zhang 1205b1683b59SHong Zhang Level: beginner 1206b1683b59SHong Zhang 1207b1683b59SHong Zhang Concepts: coloring, sparse matrix product 1208b1683b59SHong Zhang 1209b9af6bddSHong Zhang .seealso: MatTransposeColoringCreate() 1210b1683b59SHong Zhang S*/ 1211b9af6bddSHong Zhang typedef struct _p_MatTransposeColoring* MatTransposeColoring; 1212b1683b59SHong Zhang 1213014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeColoringCreate(Mat,ISColoring,MatTransposeColoring *); 1214014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransColoringApplySpToDen(MatTransposeColoring,Mat,Mat); 1215014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransColoringApplyDenToSp(MatTransposeColoring,Mat,Mat); 1216014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeColoringDestroy(MatTransposeColoring*); 1217b1683b59SHong Zhang 1218639f9d9dSBarry Smith /* 12190752156aSBarry Smith These routines are for partitioning matrices: currently used only 12203eda8832SBarry Smith for adjacency matrix, MatCreateMPIAdj(). 12210752156aSBarry Smith */ 1222ca161407SBarry Smith 1223d9274352SBarry Smith /*S 1224d9274352SBarry Smith MatPartitioning - Object for managing the partitioning of a matrix or graph 1225d9274352SBarry Smith 1226d9274352SBarry Smith Level: beginner 1227d9274352SBarry Smith 1228d9274352SBarry Smith Concepts: partitioning 1229d9274352SBarry Smith 1230743a1026Svictorle .seealso: MatPartitioningCreate(), MatPartitioningType 1231d9274352SBarry Smith S*/ 123291e9ee9fSBarry Smith typedef struct _p_MatPartitioning* MatPartitioning; 1233d9274352SBarry Smith 123476bdecfbSBarry Smith /*J 12355bc6e7ccSvictorle MatPartitioningType - String with the name of a PETSc matrix partitioning or the creation function 1236d9274352SBarry Smith with an optional dynamic library name, for example 1237d9274352SBarry Smith http://www.mcs.anl.gov/petsc/lib.a:partitioningcreate() 1238d9274352SBarry Smith 1239d9274352SBarry Smith Level: beginner 12400d04baf8SBarry Smith dm 1241b547a13bSvictorle .seealso: MatPartitioningCreate(), MatPartitioning 124276bdecfbSBarry Smith J*/ 124319fd82e9SBarry Smith typedef const char* MatPartitioningType; 12442692d6eeSBarry Smith #define MATPARTITIONINGCURRENT "current" 12452692d6eeSBarry Smith #define MATPARTITIONINGSQUARE "square" 12462692d6eeSBarry Smith #define MATPARTITIONINGPARMETIS "parmetis" 12472692d6eeSBarry Smith #define MATPARTITIONINGCHACO "chaco" 12482692d6eeSBarry Smith #define MATPARTITIONINGPARTY "party" 124950d91057SBarry Smith #define MATPARTITIONINGPTSCOTCH "ptscotch" 125071306c60Sjroman 1251ca161407SBarry Smith 1252014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningCreate(MPI_Comm,MatPartitioning*); 125319fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningSetType(MatPartitioning,MatPartitioningType); 1254014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetNParts(MatPartitioning,PetscInt); 1255014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetAdjacency(MatPartitioning,Mat); 1256014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]); 1257014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []); 1258014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningApply(MatPartitioning,IS*); 1259014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningDestroy(MatPartitioning*); 12602aabb6bbSBarry Smith 1261014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatPartitioning)); 126230de9b25SBarry Smith 126330de9b25SBarry Smith /*MC 126430de9b25SBarry Smith MatPartitioningRegisterDynamic - Adds a new sparse matrix partitioning to the 126530de9b25SBarry Smith matrix package. 126630de9b25SBarry Smith 126730de9b25SBarry Smith Synopsis: 12681890ba74SBarry Smith PetscErrorCode MatPartitioningRegisterDynamic(const char *name_partitioning,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatPartitioning)) 126930de9b25SBarry Smith 127030de9b25SBarry Smith Not Collective 127130de9b25SBarry Smith 127230de9b25SBarry Smith Input Parameters: 12732692d6eeSBarry Smith + sname - name of partitioning (for example MATPARTITIONINGCURRENT) or parmetis 127430de9b25SBarry Smith . path - location of library where creation routine is 127530de9b25SBarry Smith . name - name of function that creates the partitioning type, a string 127630de9b25SBarry Smith - function - function pointer that creates the partitioning type 127730de9b25SBarry Smith 127830de9b25SBarry Smith Level: developer 127930de9b25SBarry Smith 128030de9b25SBarry Smith If dynamic libraries are used, then the fourth input argument (function) 128130de9b25SBarry Smith is ignored. 128230de9b25SBarry Smith 128330de9b25SBarry Smith Sample usage: 128430de9b25SBarry Smith .vb 128530de9b25SBarry Smith MatPartitioningRegisterDynamic("my_part",/home/username/my_lib/lib/libO/solaris/mylib.a, 128630de9b25SBarry Smith "MyPartCreate",MyPartCreate); 128730de9b25SBarry Smith .ve 128830de9b25SBarry Smith 128930de9b25SBarry Smith Then, your partitioner can be chosen with the procedural interface via 129030de9b25SBarry Smith $ MatPartitioningSetType(part,"my_part") 129130de9b25SBarry Smith or at runtime via the option 129230de9b25SBarry Smith $ -mat_partitioning_type my_part 129330de9b25SBarry Smith 1294ab901514SBarry Smith $PETSC_ARCH occuring in pathname will be replaced with appropriate values. 129530de9b25SBarry Smith 129630de9b25SBarry Smith .keywords: matrix, partitioning, register 129730de9b25SBarry Smith 129830de9b25SBarry Smith .seealso: MatPartitioningRegisterDestroy(), MatPartitioningRegisterAll() 129930de9b25SBarry Smith M*/ 1300aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 1301f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,0) 13022aabb6bbSBarry Smith #else 1303f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,d) 13042aabb6bbSBarry Smith #endif 13052aabb6bbSBarry Smith 1306014dd563SJed Brown PETSC_EXTERN PetscBool MatPartitioningRegisterAllCalled; 1307f1144a30SSatish Balay 1308014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningRegisterAll(const char[]); 1309014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningRegisterDestroy(void); 13102bad1931SBarry Smith 1311014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningView(MatPartitioning,PetscViewer); 1312014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetFromOptions(MatPartitioning); 131319fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningGetType(MatPartitioning,MatPartitioningType*); 1314ca161407SBarry Smith 1315014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetCoarseSequential(MatPartitioning); 1316014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningParmetisGetEdgeCut(MatPartitioning, PetscInt *); 13170752156aSBarry Smith 1318b6956c12SJose Roman typedef enum { MP_CHACO_MULTILEVEL=1,MP_CHACO_SPECTRAL=2,MP_CHACO_LINEAR=4,MP_CHACO_RANDOM=5,MP_CHACO_SCATTERED=6 } MPChacoGlobalType; 13196a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoGlobalTypes[]; 1320b6956c12SJose Roman typedef enum { MP_CHACO_KERNIGHAN=1,MP_CHACO_NONE=2 } MPChacoLocalType; 13216a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoLocalTypes[]; 1322b6956c12SJose Roman typedef enum { MP_CHACO_LANCZOS=0,MP_CHACO_RQI=1 } MPChacoEigenType; 13236a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoEigenTypes[]; 1324b6956c12SJose Roman 1325014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetGlobal(MatPartitioning,MPChacoGlobalType); 1326014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetGlobal(MatPartitioning,MPChacoGlobalType*); 1327014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetLocal(MatPartitioning,MPChacoLocalType); 1328014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetLocal(MatPartitioning,MPChacoLocalType*); 1329014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal); 1330014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType); 1331014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenSolver(MatPartitioning,MPChacoEigenType*); 1332014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenTol(MatPartitioning,PetscReal); 1333014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenTol(MatPartitioning,PetscReal*); 1334014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenNumber(MatPartitioning,PetscInt); 1335014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenNumber(MatPartitioning,PetscInt*); 133671306c60Sjroman 133771306c60Sjroman #define MP_PARTY_OPT "opt" 133871306c60Sjroman #define MP_PARTY_LIN "lin" 133971306c60Sjroman #define MP_PARTY_SCA "sca" 134071306c60Sjroman #define MP_PARTY_RAN "ran" 134171306c60Sjroman #define MP_PARTY_GBF "gbf" 134271306c60Sjroman #define MP_PARTY_GCF "gcf" 134371306c60Sjroman #define MP_PARTY_BUB "bub" 134471306c60Sjroman #define MP_PARTY_DEF "def" 1345014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetGlobal(MatPartitioning,const char*); 134671306c60Sjroman #define MP_PARTY_HELPFUL_SETS "hs" 134771306c60Sjroman #define MP_PARTY_KERNIGHAN_LIN "kl" 134871306c60Sjroman #define MP_PARTY_NONE "no" 1349014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetLocal(MatPartitioning,const char*); 1350014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal); 1351014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetBipart(MatPartitioning,PetscBool); 1352014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscBool); 135371306c60Sjroman 135450d91057SBarry Smith typedef enum { MP_PTSCOTCH_QUALITY,MP_PTSCOTCH_SPEED,MP_PTSCOTCH_BALANCE,MP_PTSCOTCH_SAFETY,MP_PTSCOTCH_SCALABILITY } MPPTScotchStrategyType; 13556a6fc655SJed Brown PETSC_EXTERN const char *const MPPTScotchStrategyTypes[]; 1356e0f1cffaSJose Roman 1357014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetImbalance(MatPartitioning,PetscReal); 1358014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetImbalance(MatPartitioning,PetscReal*); 1359014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetStrategy(MatPartitioning,MPPTScotchStrategyType); 1360014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetStrategy(MatPartitioning,MPPTScotchStrategyType*); 136171306c60Sjroman 1362b43b03e9SMark F. Adams /* 1363b43b03e9SMark F. Adams These routines are for coarsening matrices: 1364b43b03e9SMark F. Adams */ 1365b43b03e9SMark F. Adams 1366b43b03e9SMark F. Adams /*S 1367b43b03e9SMark F. Adams MatCoarsen - Object for managing the coarsening of a graph (symmetric matrix) 1368b43b03e9SMark F. Adams 1369b43b03e9SMark F. Adams Level: beginner 1370b43b03e9SMark F. Adams 1371b43b03e9SMark F. Adams Concepts: coarsen 1372b43b03e9SMark F. Adams 1373b43b03e9SMark F. Adams .seealso: MatCoarsenCreate), MatCoarsenType 1374b43b03e9SMark F. Adams S*/ 1375b43b03e9SMark F. Adams typedef struct _p_MatCoarsen* MatCoarsen; 1376b43b03e9SMark F. Adams 1377b43b03e9SMark F. Adams /*J 1378b43b03e9SMark F. Adams MatCoarsenType - String with the name of a PETSc matrix coarsen or the creation function 1379b43b03e9SMark F. Adams with an optional dynamic library name, for example 1380b43b03e9SMark F. Adams http://www.mcs.anl.gov/petsc/lib.a:coarsencreate() 1381b43b03e9SMark F. Adams 1382b43b03e9SMark F. Adams Level: beginner 1383b43b03e9SMark F. Adams dm 1384b43b03e9SMark F. Adams .seealso: MatCoarsenCreate(), MatCoarsen 1385b43b03e9SMark F. Adams J*/ 138619fd82e9SBarry Smith typedef const char* MatCoarsenType; 1387b43b03e9SMark F. Adams #define MATCOARSENMIS "mis" 13889057884aSMark F. Adams #define MATCOARSENHEM "hem" 1389b43b03e9SMark F. Adams 13900cbbd2e1SMark F. Adams /* linked list for aggregates */ 139141b27cdeSMark F. Adams typedef struct _PetscCDIntNd{ 139241b27cdeSMark F. Adams struct _PetscCDIntNd *next; 13930cbbd2e1SMark F. Adams PetscInt gid; 139441b27cdeSMark F. Adams }PetscCDIntNd; 13950cbbd2e1SMark F. Adams 13960cbbd2e1SMark F. Adams /* only used by node pool */ 139741b27cdeSMark F. Adams typedef struct _PetscCDArrNd{ 139841b27cdeSMark F. Adams struct _PetscCDArrNd *next; 139941b27cdeSMark F. Adams struct _PetscCDIntNd *array; 140041b27cdeSMark F. Adams }PetscCDArrNd; 14010cbbd2e1SMark F. Adams 14020cbbd2e1SMark F. Adams typedef struct _PetscCoarsenData{ 14030cbbd2e1SMark F. Adams /* node pool */ 140441b27cdeSMark F. Adams PetscCDArrNd pool_list; 140541b27cdeSMark F. Adams PetscCDIntNd *new_node; 14060cbbd2e1SMark F. Adams PetscInt new_left; 14070cbbd2e1SMark F. Adams PetscInt chk_sz; 140841b27cdeSMark F. Adams PetscCDIntNd *extra_nodes; 14090cbbd2e1SMark F. Adams /* Array of lists */ 141041b27cdeSMark F. Adams PetscCDIntNd **array; 14110cbbd2e1SMark F. Adams PetscInt size; 14120cbbd2e1SMark F. Adams /* cache a Mat for communication data */ 14130cbbd2e1SMark F. Adams Mat mat; 1414a94c3b12SMark F. Adams /* cache IS of removed equations */ 1415e696ed0bSMark F. Adams /* IS removedIS; */ 14160cbbd2e1SMark F. Adams }PetscCoarsenData; 14170cbbd2e1SMark F. Adams 1418014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenCreate(MPI_Comm,MatCoarsen*); 141919fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenSetType(MatCoarsen,MatCoarsenType); 1420014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetAdjacency(MatCoarsen,Mat); 1421014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetGreedyOrdering(MatCoarsen,const IS); 1422014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetStrictAggs(MatCoarsen,PetscBool); 1423014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetVerbose(MatCoarsen,PetscInt); 1424014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenGetData( MatCoarsen, PetscCoarsenData ** ); 1425014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenApply(MatCoarsen); 1426014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenDestroy(MatCoarsen*); 1427b43b03e9SMark F. Adams 1428014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatCoarsen)); 1429b43b03e9SMark F. Adams 1430b43b03e9SMark F. Adams /*MC 1431b43b03e9SMark F. Adams MatCoarsenRegisterDynamic - Adds a new sparse matrix coarsen to the 1432b43b03e9SMark F. Adams matrix package. 1433b43b03e9SMark F. Adams 1434b43b03e9SMark F. Adams Synopsis: 1435b43b03e9SMark F. Adams PetscErrorCode MatCoarsenRegisterDynamic(const char *name_coarsen,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatCoarsen)) 1436b43b03e9SMark F. Adams 1437b43b03e9SMark F. Adams Not Collective 1438b43b03e9SMark F. Adams 1439b43b03e9SMark F. Adams Input Parameters: 1440b43b03e9SMark F. Adams + sname - name of coarsen (for example MATCOARSENMIS) 1441b43b03e9SMark F. Adams . path - location of library where creation routine is 1442b43b03e9SMark F. Adams . name - name of function that creates the coarsen type, a string 1443b43b03e9SMark F. Adams - function - function pointer that creates the coarsen type 1444b43b03e9SMark F. Adams 1445b43b03e9SMark F. Adams Level: developer 1446b43b03e9SMark F. Adams 1447b43b03e9SMark F. Adams If dynamic libraries are used, then the fourth input argument (function) 1448b43b03e9SMark F. Adams is ignored. 1449b43b03e9SMark F. Adams 1450b43b03e9SMark F. Adams Sample usage: 1451b43b03e9SMark F. Adams .vb 1452b43b03e9SMark F. Adams MatCoarsenRegisterDynamic("my_agg",/home/username/my_lib/lib/libO/solaris/mylib.a, 1453b43b03e9SMark F. Adams "MyAggCreate",MyAggCreate); 1454b43b03e9SMark F. Adams .ve 1455b43b03e9SMark F. Adams 1456b43b03e9SMark F. Adams Then, your aggregator can be chosen with the procedural interface via 1457b43b03e9SMark F. Adams $ MatCoarsenSetType(agg,"my_agg") 1458b43b03e9SMark F. Adams or at runtime via the option 1459b43b03e9SMark F. Adams $ -mat_coarsen_type my_agg 1460b43b03e9SMark F. Adams 1461b43b03e9SMark F. Adams $PETSC_ARCH occuring in pathname will be replaced with appropriate values. 1462b43b03e9SMark F. Adams 1463b43b03e9SMark F. Adams .keywords: matrix, coarsen, register 1464b43b03e9SMark F. Adams 1465b43b03e9SMark F. Adams .seealso: MatCoarsenRegisterDestroy(), MatCoarsenRegisterAll() 1466b43b03e9SMark F. Adams M*/ 1467b43b03e9SMark F. Adams #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 1468b43b03e9SMark F. Adams #define MatCoarsenRegisterDynamic(a,b,c,d) MatCoarsenRegister(a,b,c,0) 1469b43b03e9SMark F. Adams #else 1470b43b03e9SMark F. Adams #define MatCoarsenRegisterDynamic(a,b,c,d) MatCoarsenRegister(a,b,c,d) 1471b43b03e9SMark F. Adams #endif 1472b43b03e9SMark F. Adams 1473014dd563SJed Brown PETSC_EXTERN PetscBool MatCoarsenRegisterAllCalled; 1474b43b03e9SMark F. Adams 1475014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenRegisterAll(const char[]); 1476014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenRegisterDestroy(void); 1477b43b03e9SMark F. Adams 1478014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenView(MatCoarsen,PetscViewer); 1479014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetFromOptions(MatCoarsen); 148019fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenGetType(MatCoarsen,MatCoarsenType*); 1481b43b03e9SMark F. Adams 1482b43b03e9SMark F. Adams 1483014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMeshToVertexGraph(Mat,PetscInt,Mat*); 1484014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMeshToCellGraph(Mat,PetscInt,Mat*); 1485591294e4SBarry Smith 14860752156aSBarry Smith /* 14870a835dfdSSatish Balay If you add entries here you must also add them to finclude/petscmat.h 1488d4fbbf0eSBarry Smith */ 14891c1c02c0SLois Curfman McInnes typedef enum { MATOP_SET_VALUES=0, 14901c1c02c0SLois Curfman McInnes MATOP_GET_ROW=1, 14911c1c02c0SLois Curfman McInnes MATOP_RESTORE_ROW=2, 14921c1c02c0SLois Curfman McInnes MATOP_MULT=3, 14931c1c02c0SLois Curfman McInnes MATOP_MULT_ADD=4, 14947c922b88SBarry Smith MATOP_MULT_TRANSPOSE=5, 14957c922b88SBarry Smith MATOP_MULT_TRANSPOSE_ADD=6, 14961c1c02c0SLois Curfman McInnes MATOP_SOLVE=7, 14971c1c02c0SLois Curfman McInnes MATOP_SOLVE_ADD=8, 14987c922b88SBarry Smith MATOP_SOLVE_TRANSPOSE=9, 14997c922b88SBarry Smith MATOP_SOLVE_TRANSPOSE_ADD=10, 15001c1c02c0SLois Curfman McInnes MATOP_LUFACTOR=11, 15011c1c02c0SLois Curfman McInnes MATOP_CHOLESKYFACTOR=12, 1502a714c06dSBarry Smith MATOP_SOR=13, 15031c1c02c0SLois Curfman McInnes MATOP_TRANSPOSE=14, 15041c1c02c0SLois Curfman McInnes MATOP_GETINFO=15, 15051c1c02c0SLois Curfman McInnes MATOP_EQUAL=16, 15061c1c02c0SLois Curfman McInnes MATOP_GET_DIAGONAL=17, 15071c1c02c0SLois Curfman McInnes MATOP_DIAGONAL_SCALE=18, 15081c1c02c0SLois Curfman McInnes MATOP_NORM=19, 15091c1c02c0SLois Curfman McInnes MATOP_ASSEMBLY_BEGIN=20, 15101c1c02c0SLois Curfman McInnes MATOP_ASSEMBLY_END=21, 1511d519adbfSMatthew Knepley MATOP_SET_OPTION=22, 1512d519adbfSMatthew Knepley MATOP_ZERO_ENTRIES=23, 1513d519adbfSMatthew Knepley MATOP_ZERO_ROWS=24, 1514d519adbfSMatthew Knepley MATOP_LUFACTOR_SYMBOLIC=25, 1515d519adbfSMatthew Knepley MATOP_LUFACTOR_NUMERIC=26, 1516d519adbfSMatthew Knepley MATOP_CHOLESKY_FACTOR_SYMBOLIC=27, 1517d519adbfSMatthew Knepley MATOP_CHOLESKY_FACTOR_NUMERIC=28, 1518d519adbfSMatthew Knepley MATOP_SETUP_PREALLOCATION=29, 1519d519adbfSMatthew Knepley MATOP_ILUFACTOR_SYMBOLIC=30, 1520d519adbfSMatthew Knepley MATOP_ICCFACTOR_SYMBOLIC=31, 1521d519adbfSMatthew Knepley MATOP_GET_ARRAY=32, 1522d519adbfSMatthew Knepley MATOP_RESTORE_ARRAY=33, 1523d519adbfSMatthew Knepley MATOP_DUPLICATE=34, 1524d519adbfSMatthew Knepley MATOP_FORWARD_SOLVE=35, 1525d519adbfSMatthew Knepley MATOP_BACKWARD_SOLVE=36, 1526d519adbfSMatthew Knepley MATOP_ILUFACTOR=37, 1527d519adbfSMatthew Knepley MATOP_ICCFACTOR=38, 1528d519adbfSMatthew Knepley MATOP_AXPY=39, 1529d519adbfSMatthew Knepley MATOP_GET_SUBMATRICES=40, 1530d519adbfSMatthew Knepley MATOP_INCREASE_OVERLAP=41, 1531d519adbfSMatthew Knepley MATOP_GET_VALUES=42, 1532d519adbfSMatthew Knepley MATOP_COPY=43, 1533d519adbfSMatthew Knepley MATOP_GET_ROW_MAX=44, 1534d519adbfSMatthew Knepley MATOP_SCALE=45, 1535d519adbfSMatthew Knepley MATOP_SHIFT=46, 153635153367SBarry Smith MATOP_DIAGONAL_SET=47, 1537d519adbfSMatthew Knepley MATOP_ILUDT_FACTOR=48, 1538d519adbfSMatthew Knepley MATOP_SET_BLOCK_SIZE=49, 1539d519adbfSMatthew Knepley MATOP_GET_ROW_IJ=50, 1540d519adbfSMatthew Knepley MATOP_RESTORE_ROW_IJ=51, 1541d519adbfSMatthew Knepley MATOP_GET_COLUMN_IJ=52, 1542d519adbfSMatthew Knepley MATOP_RESTORE_COLUMN_IJ=53, 1543d519adbfSMatthew Knepley MATOP_FDCOLORING_CREATE=54, 1544d519adbfSMatthew Knepley MATOP_COLORING_PATCH=55, 1545d519adbfSMatthew Knepley MATOP_SET_UNFACTORED=56, 1546d519adbfSMatthew Knepley MATOP_PERMUTE=57, 1547d519adbfSMatthew Knepley MATOP_SET_VALUES_BLOCKED=58, 1548d519adbfSMatthew Knepley MATOP_GET_SUBMATRIX=59, 1549d519adbfSMatthew Knepley MATOP_DESTROY=60, 1550d519adbfSMatthew Knepley MATOP_VIEW=61, 1551d519adbfSMatthew Knepley MATOP_CONVERT_FROM=62, 15527bab7c10SHong Zhang MATOP_MATMAT_MULT=63, 15537bab7c10SHong Zhang MATOP_MATMAT_MULT_SYMBOLIC=64, 15547bab7c10SHong Zhang MATOP_MATMAT_MULT_NUMERIC=65, 1555d519adbfSMatthew Knepley MATOP_SET_LOCAL_TO_GLOBAL_MAP=66, 1556d519adbfSMatthew Knepley MATOP_SET_VALUES_LOCAL=67, 1557d519adbfSMatthew Knepley MATOP_ZERO_ROWS_LOCAL=68, 1558d519adbfSMatthew Knepley MATOP_GET_ROW_MAX_ABS=69, 1559d519adbfSMatthew Knepley MATOP_GET_ROW_MIN_ABS=70, 1560d519adbfSMatthew Knepley MATOP_CONVERT=71, 1561d519adbfSMatthew Knepley MATOP_SET_COLORING=72, 1562*bcaeba4dSBarry Smith MATOP_PLACEHOLDER=73, 1563d519adbfSMatthew Knepley MATOP_SET_VALUES_ADIFOR=74, 1564d519adbfSMatthew Knepley MATOP_FD_COLORING_APPLY=75, 1565d519adbfSMatthew Knepley MATOP_SET_FROM_OPTIONS=76, 1566bda6c4cbSBarry Smith MATOP_MULT_CONSTRAINED=77, 1567bda6c4cbSBarry Smith MATOP_MULT_TRANSPOSE_CONSTRAIN=78, 1568d519adbfSMatthew Knepley MATOP_PERMUTE_SPARSIFY=79, 1569d519adbfSMatthew Knepley MATOP_MULT_MULTIPLE=80, 1570d519adbfSMatthew Knepley MATOP_SOLVE_MULTIPLE=81, 1571d519adbfSMatthew Knepley MATOP_GET_INERTIA=82, 1572d519adbfSMatthew Knepley MATOP_LOAD=83, 1573d519adbfSMatthew Knepley MATOP_IS_SYMMETRIC=84, 1574d519adbfSMatthew Knepley MATOP_IS_HERMITIAN=85, 1575d519adbfSMatthew Knepley MATOP_IS_STRUCTURALLY_SYMMETRIC=86, 1576820469bcSHong Zhang MATOP_SET_VALUES_BLOCKEDLOCAL=87, 1577d519adbfSMatthew Knepley MATOP_GET_VECS=88, 1578d519adbfSMatthew Knepley MATOP_MAT_MULT=89, 1579d519adbfSMatthew Knepley MATOP_MAT_MULT_SYMBOLIC=90, 1580d519adbfSMatthew Knepley MATOP_MAT_MULT_NUMERIC=91, 1581d519adbfSMatthew Knepley MATOP_PTAP=92, 1582d519adbfSMatthew Knepley MATOP_PTAP_SYMBOLIC=93, 1583d519adbfSMatthew Knepley MATOP_PTAP_NUMERIC=94, 1584bda6c4cbSBarry Smith MATOP_MAT_TRANSPOSE_MULT=95, 1585bda6c4cbSBarry Smith MATOP_MAT_TRANSPOSE_MULT_SYMBO=96, 1586bda6c4cbSBarry Smith MATOP_MAT_TRANSPOSE_MULT_NUMER=97, 1587820469bcSHong Zhang MATOP_DUMMY98=98, 1588820469bcSHong Zhang MATOP_DUMMY99=99, 1589820469bcSHong Zhang MATOP_DUMMY100=100, 1590820469bcSHong Zhang MATOP_DUMMY101=101, 1591d519adbfSMatthew Knepley MATOP_CONJUGATE=102, 1592d519adbfSMatthew Knepley MATOP_SET_SIZES=103, 1593d519adbfSMatthew Knepley MATOP_SET_VALUES_ROW=104, 1594d519adbfSMatthew Knepley MATOP_REAL_PART=105, 1595bda6c4cbSBarry Smith MATOP_IMAGINARY_PART=106, 1596bda6c4cbSBarry Smith MATOP_GET_ROW_UPPER_TRIANGULAR=107, 1597bda6c4cbSBarry Smith MATOP_RESTORE_ROW_UPPER_TRIANG=108, 1598bda6c4cbSBarry Smith MATOP_MAT_SOLVE=109, 15990e8d04f7SBarry Smith MATOP_GET_REDUNDANT_MATRIX=110, 1600d519adbfSMatthew Knepley MATOP_GET_ROW_MIN=111, 16010e8d04f7SBarry Smith MATOP_GET_COLUMN_VECTOR=112, 1602d519adbfSMatthew Knepley MATOP_MISSING_DIAGONAL=113, 16030e8d04f7SBarry Smith MATOP_GET_SEQ_NONZERO_STRUCTUR=114, 160489c6957cSBarry Smith MATOP_CREATE=115, 160589c6957cSBarry Smith MATOP_GET_GHOSTS=116, 16060e8d04f7SBarry Smith MATOP_GET_LOCAL_SUB_MATRIX=117, 16070e8d04f7SBarry Smith MATOP_RESTORE_LOCALSUB_MATRIX=118, 1608eeffb40dSHong Zhang MATOP_MULT_DIAGONAL_BLOCK=119, 16090e8d04f7SBarry Smith MATOP_HERMITIAN_TRANSPOSE=120, 16100e8d04f7SBarry Smith MATOP_MULT_HERMITIAN_TRANSPOSE=121, 16110e8d04f7SBarry Smith MATOP_MULT_HERMITIAN_TRANS_ADD=122, 16120e8d04f7SBarry Smith MATOP_GET_MULTI_PROC_BLOCK=123, 16130e8d04f7SBarry Smith MATOP_GET_COLUMN_NORMS=125, 16140e8d04f7SBarry Smith MATOP_GET_SUB_MATRICES_PARALLE=128, 16152b8ad9a3SHong Zhang MATOP_SET_VALUES_BATCH=129, 16160e8d04f7SBarry Smith MATOP_TRANSPOSE_MAT_MULT=130, 16170e8d04f7SBarry Smith MATOP_TRANSPOSE_MAT_MULT_SYMBO=131, 16180e8d04f7SBarry Smith MATOP_TRANSPOSE_MAT_MULT_NUMER=132, 16190e8d04f7SBarry Smith MATOP_TRANSPOSE_COLORING_CREAT=133, 16200e8d04f7SBarry Smith MATOP_TRANS_COLORING_APPLY_SPT=134, 16210e8d04f7SBarry Smith MATOP_TRANS_COLORING_APPLY_DEN=135, 16220e8d04f7SBarry Smith MATOP_RART=136, 16230e8d04f7SBarry Smith MATOP_RART_SYMBOLIC=137, 16240e8d04f7SBarry Smith MATOP_RART_NUMERIC=138, 1625e09a3074SHong Zhang MATOP_SET_BLOCK_SIZES=139, 1626e09a3074SHong Zhang MATOP_AYPX=140 1627fae171e0SBarry Smith } MatOperation; 1628014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatHasOperation(Mat,MatOperation,PetscBool *); 1629014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellSetOperation(Mat,MatOperation,void(*)(void)); 1630014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellGetOperation(Mat,MatOperation,void(**)(void)); 1631014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellSetContext(Mat,void*); 1632112a2221SBarry Smith 163390ace30eSBarry Smith /* 163490ace30eSBarry Smith Codes for matrices stored on disk. By default they are 163590ace30eSBarry Smith stored in a universal format. By changing the format with 16367973ad04SBarry Smith PetscViewerSetFormat(viewer,PETSC_VIEWER_NATIVE); the matrices will 163790ace30eSBarry Smith be stored in a way natural for the matrix, for example dense matrices 163890ace30eSBarry Smith would be stored as dense. Matrices stored this way may only be 1639f4403165SShri Abhyankar read into matrices of the same type. 164090ace30eSBarry Smith */ 164190ace30eSBarry Smith #define MATRIX_BINARY_FORMAT_DENSE -1 164290ace30eSBarry Smith 1643014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetHashTableFactor(Mat,PetscReal); 1644014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatISGetLocalMat(Mat,Mat*); 1645014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatISSetLocalMat(Mat,Mat); 16461f4e1ec7SBarry Smith 1647d9274352SBarry Smith /*S 1648d9274352SBarry Smith MatNullSpace - Object that removes a null space from a vector, i.e. 1649d9274352SBarry Smith orthogonalizes the vector to a subsapce 1650d9274352SBarry Smith 1651f7a9e4ceSBarry Smith Level: advanced 1652d9274352SBarry Smith 1653d9274352SBarry Smith Concepts: matrix; linear operator, null space 1654d9274352SBarry Smith 16556e1639daSBarry Smith Users manual sections: 16566e1639daSBarry Smith . sec_singular 16576e1639daSBarry Smith 1658d9274352SBarry Smith .seealso: MatNullSpaceCreate() 1659d9274352SBarry Smith S*/ 166074637425SBarry Smith typedef struct _p_MatNullSpace* MatNullSpace; 1661d9274352SBarry Smith 1662014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceCreate(MPI_Comm,PetscBool ,PetscInt,const Vec[],MatNullSpace*); 1663014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(MatNullSpace,Vec,void*),void*); 1664014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceDestroy(MatNullSpace*); 1665014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceRemove(MatNullSpace,Vec,Vec*); 1666014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNullSpace(Mat, MatNullSpace *); 1667014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNullSpace(Mat,MatNullSpace); 1668014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNearNullSpace(Mat,MatNullSpace); 1669014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNearNullSpace(Mat,MatNullSpace*); 1670014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceTest(MatNullSpace,Mat,PetscBool *); 1671014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceView(MatNullSpace,PetscViewer); 1672014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceGetVecs(MatNullSpace,PetscBool*,PetscInt*,const Vec**); 1673014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceCreateRigidBody(Vec,MatNullSpace*); 167474637425SBarry Smith 1675014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatReorderingSeqSBAIJ(Mat,IS); 1676014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetHashTableFactor(Mat,PetscReal); 1677014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetColumnIndices(Mat,PetscInt *); 1678014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJInvertBlockDiagonal(Mat); 16793f1d51d7SBarry Smith 1680014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMAIJ(Mat,PetscInt,Mat*); 1681014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMAIJRedimension(Mat,PetscInt,Mat*); 1682014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMAIJGetAIJ(Mat,Mat*); 1683c4f061fbSSatish Balay 1684014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatComputeExplicitOperator(Mat,Mat*); 1685b0a32e0cSBarry Smith 1686014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalScaleLocal(Mat,Vec); 168704f1ad80SBarry Smith 1688014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*); 1689014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetBase(Mat,Vec,Vec); 1690014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*); 1691014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*)); 1692014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec)); 1693014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDAddNullSpace(Mat,MatNullSpace); 1694014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt); 1695014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDResetHHistory(Mat); 1696014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctionError(Mat,PetscReal); 1697014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetPeriod(Mat,PetscInt); 1698014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDGetH(Mat,PetscScalar *); 1699014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetOptionsPrefix(Mat,const char[]); 1700014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*); 1701014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*); 1702e884886fSBarry Smith 17036370053bSBarry Smith /*S 17046370053bSBarry Smith MatMFFD - A data structured used to manage the computation of the h differencing parameter for matrix-free 17056370053bSBarry Smith Jacobian vector products 1706e884886fSBarry Smith 17076370053bSBarry Smith Notes: MATMFFD is a specific MatType which uses the MatMFFD data structure 17086370053bSBarry Smith 17096370053bSBarry Smith MatMFFD*() methods actually take the Mat as their first argument. Not a MatMFFD data structure 17106370053bSBarry Smith 17116370053bSBarry Smith Level: developer 17126370053bSBarry Smith 17136370053bSBarry Smith .seealso: MATMFFD, MatCreateMFFD(), MatMFFDSetFuction(), MatMFFDSetType(), MatMFFDRegister() 17146370053bSBarry Smith S*/ 1715e884886fSBarry Smith typedef struct _p_MatMFFD* MatMFFD; 1716e884886fSBarry Smith 171776bdecfbSBarry Smith /*J 1718e884886fSBarry Smith MatMFFDType - algorithm used to compute the h used in computing matrix-vector products via differencing of the function 1719e884886fSBarry Smith 1720e884886fSBarry Smith Level: beginner 1721e884886fSBarry Smith 1722e884886fSBarry Smith .seealso: MatMFFDSetType(), MatMFFDRegister() 172376bdecfbSBarry Smith J*/ 172419fd82e9SBarry Smith typedef const char* MatMFFDType; 1725e884886fSBarry Smith #define MATMFFD_DS "ds" 1726e884886fSBarry Smith #define MATMFFD_WP "wp" 1727e884886fSBarry Smith 172819fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDSetType(Mat,MatMFFDType); 1729014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatMFFD)); 1730e884886fSBarry Smith 1731e884886fSBarry Smith /*MC 1732e884886fSBarry Smith MatMFFDRegisterDynamic - Adds a method to the MatMFFD registry. 1733e884886fSBarry Smith 1734e884886fSBarry Smith Synopsis: 17351890ba74SBarry Smith PetscErrorCode MatMFFDRegisterDynamic(const char *name_solver,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatMFFD)) 1736e884886fSBarry Smith 1737e884886fSBarry Smith Not Collective 1738e884886fSBarry Smith 1739e884886fSBarry Smith Input Parameters: 1740e884886fSBarry Smith + name_solver - name of a new user-defined compute-h module 1741e884886fSBarry Smith . path - path (either absolute or relative) the library containing this solver 1742e884886fSBarry Smith . name_create - name of routine to create method context 1743e884886fSBarry Smith - routine_create - routine to create method context 1744e884886fSBarry Smith 1745e884886fSBarry Smith Level: developer 1746e884886fSBarry Smith 1747e884886fSBarry Smith Notes: 1748e884886fSBarry Smith MatMFFDRegisterDynamic() may be called multiple times to add several user-defined solvers. 1749e884886fSBarry Smith 1750e884886fSBarry Smith If dynamic libraries are used, then the fourth input argument (routine_create) 1751e884886fSBarry Smith is ignored. 1752e884886fSBarry Smith 1753e884886fSBarry Smith Sample usage: 1754e884886fSBarry Smith .vb 1755e884886fSBarry Smith MatMFFDRegisterDynamic("my_h",/home/username/my_lib/lib/libO/solaris/mylib.a, 1756e884886fSBarry Smith "MyHCreate",MyHCreate); 1757e884886fSBarry Smith .ve 1758e884886fSBarry Smith 1759e884886fSBarry Smith Then, your solver can be chosen with the procedural interface via 1760e884886fSBarry Smith $ MatMFFDSetType(mfctx,"my_h") 1761e884886fSBarry Smith or at runtime via the option 1762e884886fSBarry Smith $ -snes_mf_type my_h 1763e884886fSBarry Smith 1764e884886fSBarry Smith .keywords: MatMFFD, register 1765e884886fSBarry Smith 1766e884886fSBarry Smith .seealso: MatMFFDRegisterAll(), MatMFFDRegisterDestroy() 1767e884886fSBarry Smith M*/ 1768e884886fSBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 1769e884886fSBarry Smith #define MatMFFDRegisterDynamic(a,b,c,d) MatMFFDRegister(a,b,c,0) 1770e884886fSBarry Smith #else 1771e884886fSBarry Smith #define MatMFFDRegisterDynamic(a,b,c,d) MatMFFDRegister(a,b,c,d) 1772e884886fSBarry Smith #endif 1773e884886fSBarry Smith 1774014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDRegisterAll(const char[]); 1775014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDRegisterDestroy(void); 1776014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDDSSetUmin(Mat,PetscReal); 1777014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDWPSetComputeNormU(Mat,PetscBool ); 1778e884886fSBarry Smith 1779e884886fSBarry Smith 1780014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *); 1781014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *); 17827dbadf16SMatthew Knepley 178397969023SHong Zhang /* 178497969023SHong Zhang PETSc interface to MUMPS 178597969023SHong Zhang */ 178697969023SHong Zhang #ifdef PETSC_HAVE_MUMPS 1787014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMumpsSetIcntl(Mat,PetscInt,PetscInt); 178897969023SHong Zhang #endif 178923a5497aSJed Brown 1790d954db57SHong Zhang /* 1791d954db57SHong Zhang PETSc interface to SUPERLU 1792d954db57SHong Zhang */ 1793d954db57SHong Zhang #ifdef PETSC_HAVE_SUPERLU 1794014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSuperluSetILUDropTol(Mat,PetscReal); 1795d954db57SHong Zhang #endif 1796d954db57SHong Zhang 17979ae82921SPaul Mullowney #if defined PETSC_HAVE_TXPETSCGPU 1798e057df02SPaul Mullowney 17993f9c0db1SPaul Mullowney /*E 1800e057df02SPaul Mullowney MatCUSPARSEStorageFormat - indicates the storage format for CUSPARSE (GPU) 1801e057df02SPaul Mullowney matrices. 1802e057df02SPaul Mullowney 1803e057df02SPaul Mullowney Not Collective 1804e057df02SPaul Mullowney 1805e057df02SPaul Mullowney + MAT_CUSPARSE_CSR : Compressed Sparse Row 1806e057df02SPaul Mullowney . MAT_CUSPARSE_ELL : Ellpack 1807e057df02SPaul Mullowney - MAT_CUSPARSE_HYB : Hybrid, a combination of Ellpack and Coordinate format. 1808e057df02SPaul Mullowney 1809e057df02SPaul Mullowney Level: intermediate 1810e057df02SPaul Mullowney 1811e057df02SPaul Mullowney Any additions/changes here MUST also be made in include/finclude/petscmat.h 1812e057df02SPaul Mullowney 1813e057df02SPaul Mullowney .seealso: MatCUSPARSESetFormat(), MatCUSPARSEFormatOperation 1814e057df02SPaul Mullowney E*/ 1815e057df02SPaul Mullowney 18163f9c0db1SPaul Mullowney typedef enum {MAT_CUSPARSE_CSR, MAT_CUSPARSE_ELL, MAT_CUSPARSE_HYB} MatCUSPARSEStorageFormat; 1817e057df02SPaul Mullowney 1818e057df02SPaul Mullowney /* these will be strings associated with enumerated type defined above */ 1819e057df02SPaul Mullowney PETSC_EXTERN const char *const MatCUSPARSEStorageFormats[]; 1820e057df02SPaul Mullowney 18213f9c0db1SPaul Mullowney /*E 1822e057df02SPaul Mullowney MatCUSPARSEFormatOperation - indicates the operation of CUSPARSE (GPU) 1823e057df02SPaul Mullowney matrices whose operation should use a particular storage format. 1824e057df02SPaul Mullowney 1825e057df02SPaul Mullowney Not Collective 1826e057df02SPaul Mullowney 1827e057df02SPaul Mullowney + MAT_CUSPARSE_MULT_DIAG : sets the storage format for the diagonal matrix in the parallel MatMult 1828e057df02SPaul Mullowney . MAT_CUSPARSE_MULT_OFFDIAG : sets the storage format for the offdiagonal matrix in the parallel MatMult 18293f9c0db1SPaul Mullowney . MAT_CUSPARSE_MULT : sets the storage format for the entire matrix in the serial (single GPU) MatMult 18303f9c0db1SPaul Mullowney . MAT_CUSPARSE_SOLVE : sets the storage format for the triangular factors in the serial (single GPU) MatSolve 18313f9c0db1SPaul Mullowney - MAT_CUSPARSE_ALL : sets the storage format for all CUSPARSE (GPU) matrices 1832e057df02SPaul Mullowney 1833e057df02SPaul Mullowney Level: intermediate 1834e057df02SPaul Mullowney 1835e057df02SPaul Mullowney .seealso: MatCUSPARSESetFormat(), MatCUSPARSEStorageFormat 1836e057df02SPaul Mullowney E*/ 18373f9c0db1SPaul Mullowney typedef enum {MAT_CUSPARSE_MULT_DIAG, MAT_CUSPARSE_MULT_OFFDIAG, MAT_CUSPARSE_MULT, MAT_CUSPARSE_SOLVE, MAT_CUSPARSE_ALL} MatCUSPARSEFormatOperation; 1838e057df02SPaul Mullowney 183910e9db80SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 184010e9db80SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCreateAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 1841e057df02SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCUSPARSESetFormat(Mat,MatCUSPARSEFormatOperation,MatCUSPARSEStorageFormat); 1842e057df02SPaul Mullowney 18439ae82921SPaul Mullowney #endif 18449ae82921SPaul Mullowney 184590c53902SBarry Smith #if defined(PETSC_HAVE_CUSP) 1846014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCUSP(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 1847014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateAIJCUSP(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 1848e057df02SPaul Mullowney 18493f9c0db1SPaul Mullowney /*E 1850e057df02SPaul Mullowney MatCUSPStorageFormat - indicates the storage format for CUSP (GPU) 1851e057df02SPaul Mullowney matrices. 1852e057df02SPaul Mullowney 1853e057df02SPaul Mullowney Not Collective 1854e057df02SPaul Mullowney 18553f9c0db1SPaul Mullowney + MAT_CUSP_CSR : Compressed Sparse Row 18563f9c0db1SPaul Mullowney . MAT_CUSP_DIA : Diagonal 18573f9c0db1SPaul Mullowney - MAT_CUSP_ELL : Ellpack 1858e057df02SPaul Mullowney 1859e057df02SPaul Mullowney Level: intermediate 1860e057df02SPaul Mullowney 1861e057df02SPaul Mullowney Any additions/changes here MUST also be made in include/finclude/petscmat.h 1862e057df02SPaul Mullowney 1863e057df02SPaul Mullowney .seealso: MatCUSPSetFormat(), MatCUSPFormatOperation 1864e057df02SPaul Mullowney E*/ 18653f9c0db1SPaul Mullowney typedef enum {MAT_CUSP_CSR, MAT_CUSP_DIA, MAT_CUSP_ELL} MatCUSPStorageFormat; 1866e057df02SPaul Mullowney 1867e057df02SPaul Mullowney /* these will be strings associated with enumerated type defined above */ 1868e057df02SPaul Mullowney PETSC_EXTERN const char *const MatCUSPStorageFormats[]; 1869e057df02SPaul Mullowney 18703f9c0db1SPaul Mullowney /*E 1871e057df02SPaul Mullowney MatCUSPFormatOperation - indicates the operation of CUSP (GPU) 1872e057df02SPaul Mullowney matrices whose operation should use a particular storage format. 1873e057df02SPaul Mullowney 1874e057df02SPaul Mullowney Not Collective 1875e057df02SPaul Mullowney 1876e057df02SPaul Mullowney + MAT_CUSP_MULT_DIAG : sets the storage format for the diagonal matrix in the parallel MatMult 1877e057df02SPaul Mullowney . MAT_CUSP_MULT_OFFDIAG : sets the storage format for the offdiagonal matrix in the parallel MatMult 18783f9c0db1SPaul Mullowney . MAT_CUSP_MULT : sets the storage format for the entire matrix in the serial (single GPU) MatMult 18793f9c0db1SPaul Mullowney . MAT_CUSP_SOLVE : sets the storage format for the triangular factors in the serial (single GPU) MatSolve 18803f9c0db1SPaul Mullowney - MAT_CUSP_ALL : sets the storage format for all CUSP (GPU) matrices 1881e057df02SPaul Mullowney 1882e057df02SPaul Mullowney Level: intermediate 1883e057df02SPaul Mullowney 1884e057df02SPaul Mullowney Any additions/changes here MUST also be made in include/finclude/petscmat.h 1885e057df02SPaul Mullowney 1886e057df02SPaul Mullowney .seealso: MatCUSPSetFormat(), MatCUSPStorageFormat 1887e057df02SPaul Mullowney E*/ 18883f9c0db1SPaul Mullowney typedef enum {MAT_CUSP_MULT_DIAG, MAT_CUSP_MULT_OFFDIAG, MAT_CUSP_MULT, MAT_CUSP_SOLVE, MAT_CUSP_ALL} MatCUSPFormatOperation; 1889e057df02SPaul Mullowney 1890e057df02SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCUSPSetFormat(Mat,MatCUSPFormatOperation,MatCUSPStorageFormat); 1891e057df02SPaul Mullowney #endif 189290c53902SBarry Smith 189354efbe56SHong Zhang /* 189454efbe56SHong Zhang PETSc interface to FFTW 189554efbe56SHong Zhang */ 189654efbe56SHong Zhang #if defined(PETSC_HAVE_FFTW) 1897014dd563SJed Brown PETSC_EXTERN PetscErrorCode VecScatterPetscToFFTW(Mat,Vec,Vec); 1898014dd563SJed Brown PETSC_EXTERN PetscErrorCode VecScatterFFTWToPetsc(Mat,Vec,Vec); 1899014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetVecsFFTW(Mat,Vec*,Vec*,Vec*); 190054efbe56SHong Zhang #endif 190154efbe56SHong Zhang 1902db31f6deSJed Brown #if defined(PETSC_HAVE_ELEMENTAL) 1903db31f6deSJed Brown PETSC_EXTERN PetscErrorCode PetscElementalInitializePackage(const char*); 1904db31f6deSJed Brown PETSC_EXTERN PetscErrorCode PetscElementalFinalizePackage(void); 1905db31f6deSJed Brown #endif 1906db31f6deSJed Brown 1907014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateNest(MPI_Comm,PetscInt,const IS[],PetscInt,const IS[],const Mat[],Mat*); 1908014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSize(Mat,PetscInt*,PetscInt*); 1909014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetISs(Mat,IS[],IS[]); 1910014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetLocalISs(Mat,IS[],IS[]); 1911014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSubMats(Mat,PetscInt*,PetscInt*,Mat***); 1912014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSubMat(Mat,PetscInt,PetscInt,Mat*); 191319fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatNestSetVecType(Mat,VecType); 1914014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestSetSubMats(Mat,PetscInt,const IS[],PetscInt,const IS[],const Mat[]); 1915014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestSetSubMat(Mat,PetscInt,PetscInt,Mat); 1916d8588912SDave May 191715c1f1baSDmitry Karpeev /* 191815c1f1baSDmitry Karpeev MatIJ: 191915c1f1baSDmitry Karpeev An unweighted directed pseudograph 192015c1f1baSDmitry Karpeev An interpretation of this matrix as a (pseudo)graph allows us to define additional operations on it: 192115c1f1baSDmitry Karpeev A MatIJ can act on sparse arrays: arrays of indices, or index arrays of integers, scalars, or integer-scalar pairs 192215c1f1baSDmitry Karpeev by mapping the indices to the indices connected to them by the (pseudo)graph ed 192315c1f1baSDmitry Karpeev */ 192415c1f1baSDmitry Karpeev typedef enum {MATIJ_LOCAL, MATIJ_GLOBAL} MatIJIndexType; 1925014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIJSetMultivalued(Mat, PetscBool); 1926014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIJGetMultivalued(Mat, PetscBool*); 1927014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIJSetEdges(Mat, PetscInt, const PetscInt*, const PetscInt*); 1928014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIJGetEdges(Mat, PetscInt *, PetscInt **, PetscInt **); 1929014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIJSetEdgesIS(Mat, IS, IS); 1930014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIJGetEdgesIS(Mat, IS*, IS*); 1931014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIJGetRowSizes(Mat, MatIJIndexType, PetscInt, const PetscInt *, PetscInt **); 1932014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIJGetMinRowSize(Mat, PetscInt *); 1933014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIJGetMaxRowSize(Mat, PetscInt *); 1934014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIJGetSupport(Mat, PetscInt *, PetscInt **); 1935014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIJGetSupportIS(Mat, IS *); 1936014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIJGetImage(Mat, PetscInt*, PetscInt**); 1937014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIJGetImageIS(Mat, IS *); 1938014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIJGetSupportSize(Mat, PetscInt *); 1939014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIJGetImageSize(Mat, PetscInt *); 194015c1f1baSDmitry Karpeev 1941014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIJBinRenumber(Mat, Mat*); 194215c1f1baSDmitry Karpeev 1943014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIJMap(Mat, MatIJIndexType, PetscInt,const PetscInt*,const PetscInt*,const PetscScalar*, MatIJIndexType,PetscInt*,PetscInt**,PetscInt**,PetscScalar**,PetscInt**); 1944014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIJBin(Mat, MatIJIndexType, PetscInt,const PetscInt*,const PetscInt*,const PetscScalar*,PetscInt*,PetscInt**,PetscInt**,PetscScalar**,PetscInt**); 1945014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIJBinMap(Mat,Mat, MatIJIndexType,PetscInt,const PetscInt*,const PetscInt*,const PetscScalar*,MatIJIndexType,PetscInt*,PetscInt**,PetscInt**,PetscScalar**,PetscInt**); 194615c1f1baSDmitry Karpeev 19474325cce7SMatthew G Knepley PETSC_EXTERN PetscErrorCode MatChop(Mat,PetscReal); 19484325cce7SMatthew G Knepley 194923a5497aSJed Brown #endif 1950