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 98f6c3df8SBarry Smith Mat - Abstract PETSc matrix object used to manage all linear operators in PETSc, even those without 108f6c3df8SBarry Smith an explicit sparse representation (such as matrix-free operators) 112eac72dbSBarry Smith 12d91e6319SBarry Smith Level: beginner 13d91e6319SBarry Smith 14d9274352SBarry Smith Concepts: matrix; linear operator 15d9274352SBarry Smith 168f6c3df8SBarry Smith .seealso: MatCreate(), MatType, MatSetType(), MatDestroy() 17d9274352SBarry Smith S*/ 18d9274352SBarry Smith typedef struct _p_Mat* Mat; 19d9274352SBarry Smith 2076bdecfbSBarry Smith /*J 218f6c3df8SBarry Smith MatType - String with the name of a PETSc matrix type 22d9274352SBarry Smith 23d9274352SBarry Smith Level: beginner 24d9274352SBarry Smith 258f6c3df8SBarry Smith .seealso: MatSetType(), Mat, MatSolverPackage, MatRegister() 2676bdecfbSBarry Smith J*/ 2719fd82e9SBarry Smith typedef const char* MatType; 28273d9f13SBarry Smith #define MATSAME "same" 295a11e1b2SBarry Smith #define MATMAIJ "maij" 30273d9f13SBarry Smith #define MATSEQMAIJ "seqmaij" 31273d9f13SBarry Smith #define MATMPIMAIJ "mpimaij" 32273d9f13SBarry Smith #define MATIS "is" 335a11e1b2SBarry Smith #define MATAIJ "aij" 34273d9f13SBarry Smith #define MATSEQAIJ "seqaij" 357d6a0e61SBarry Smith #define MATSEQAIJPTHREAD "seqaijpthread" 36bf2c1783SBarry Smith #define MATAIJPTHREAD "aijpthread" 37273d9f13SBarry Smith #define MATMPIAIJ "mpiaij" 385a11e1b2SBarry Smith #define MATAIJCRL "aijcrl" 395a11e1b2SBarry Smith #define MATSEQAIJCRL "seqaijcrl" 405a11e1b2SBarry Smith #define MATMPIAIJCRL "mpiaijcrl" 418154be41SBarry Smith #define MATAIJCUSP "aijcusp" 428154be41SBarry Smith #define MATSEQAIJCUSP "seqaijcusp" 438154be41SBarry Smith #define MATMPIAIJCUSP "mpiaijcusp" 449ae82921SPaul Mullowney #define MATAIJCUSPARSE "aijcusparse" 459ae82921SPaul Mullowney #define MATSEQAIJCUSPARSE "seqaijcusparse" 469ae82921SPaul Mullowney #define MATMPIAIJCUSPARSE "mpiaijcusparse" 47d67ff14aSKarl Rupp #define MATAIJVIENNACL "aijviennacl" 48d67ff14aSKarl Rupp #define MATSEQAIJVIENNACL "seqaijviennacl" 49d67ff14aSKarl Rupp #define MATMPIAIJVIENNACL "mpiaijviennacl" 505a11e1b2SBarry Smith #define MATAIJPERM "aijperm" 515a11e1b2SBarry Smith #define MATSEQAIJPERM "seqaijperm" 525a11e1b2SBarry Smith #define MATMPIAIJPERM "mpiaijperm" 53273d9f13SBarry Smith #define MATSHELL "shell" 545a11e1b2SBarry Smith #define MATDENSE "dense" 55209238afSKris Buschelman #define MATSEQDENSE "seqdense" 56273d9f13SBarry Smith #define MATMPIDENSE "mpidense" 57db31f6deSJed Brown #define MATELEMENTAL "elemental" 585a11e1b2SBarry Smith #define MATBAIJ "baij" 59273d9f13SBarry Smith #define MATSEQBAIJ "seqbaij" 60273d9f13SBarry Smith #define MATMPIBAIJ "mpibaij" 61273d9f13SBarry Smith #define MATMPIADJ "mpiadj" 625a11e1b2SBarry Smith #define MATSBAIJ "sbaij" 63273d9f13SBarry Smith #define MATSEQSBAIJ "seqsbaij" 64273d9f13SBarry Smith #define MATMPISBAIJ "mpisbaij" 65c0cdd4a1SDahai Guo #define MATSEQBSTRM "seqbstrm" 66c0cdd4a1SDahai Guo #define MATMPIBSTRM "mpibstrm" 67c0cdd4a1SDahai Guo #define MATBSTRM "bstrm" 68aa5a9175SDahai Guo #define MATSEQSBSTRM "seqsbstrm" 69aa5a9175SDahai Guo #define MATMPISBSTRM "mpisbstrm" 70aa5a9175SDahai Guo #define MATSBSTRM "sbstrm" 71cebc7f6cSBarry Smith #define MATDAAD "daad" 72cebc7f6cSBarry Smith #define MATMFFD "mffd" 73c8a8475eSBarry Smith #define MATNORMAL "normal" 74ab92ecdeSBarry Smith #define MATLRC "lrc" 752a6744ebSBarry Smith #define MATSCATTER "scatter" 76421e10b8SBarry Smith #define MATBLOCKMAT "blockmat" 77793850ffSBarry Smith #define MATCOMPOSITE "composite" 781f8c7532SHong Zhang #define MATFFT "fft" 791f8c7532SHong Zhang #define MATFFTW "fftw" 80e133240eSMatthew G Knepley #define MATSEQCUFFT "seqcufft" 81557cca28SSatish Balay #define MATTRANSPOSEMAT "transpose" 8272ca8751SBarry Smith #define MATSCHURCOMPLEMENT "schurcomplement" 831d6018f0SLisandro Dalcin #define MATPYTHON "python" 84f91d8e95SBarry Smith #define MATHYPRESTRUCT "hyprestruct" 85a9e6138eSGlenn Hammond #define MATHYPRESSTRUCT "hypresstruct" 86ee1cef2cSJed Brown #define MATSUBMATRIX "submatrix" 872c0dbf93SJed Brown #define MATLOCALREF "localref" 88d8588912SDave May #define MATNEST "nest" 89773941d6SBarry Smith 9076bdecfbSBarry Smith /*J 91c7393fdbSBarry Smith MatSolverPackage - String with the name of a PETSc matrix solver type. 929989ab13SBarry Smith 939989ab13SBarry Smith For example: "petsc" indicates what PETSc provides, "superlu" indicates either 949989ab13SBarry Smith SuperLU or SuperLU_Dist etc. 959989ab13SBarry Smith 969989ab13SBarry Smith 979989ab13SBarry Smith Level: beginner 989989ab13SBarry Smith 995c9eb25fSBarry Smith .seealso: MatGetFactor(), Mat, MatSetType(), MatType 10076bdecfbSBarry Smith J*/ 101c7393fdbSBarry Smith #define MatSolverPackage char* 1022692d6eeSBarry Smith #define MATSOLVERSUPERLU "superlu" 1032692d6eeSBarry Smith #define MATSOLVERSUPERLU_DIST "superlu_dist" 1042692d6eeSBarry Smith #define MATSOLVERUMFPACK "umfpack" 10520db9a53SJed Brown #define MATSOLVERCHOLMOD "cholmod" 1062692d6eeSBarry Smith #define MATSOLVERESSL "essl" 1072692d6eeSBarry Smith #define MATSOLVERLUSOL "lusol" 1082692d6eeSBarry Smith #define MATSOLVERMUMPS "mumps" 109d615f992SSatish Balay #define MATSOLVERMKL_PARDISO "mkl_pardiso" 1102692d6eeSBarry Smith #define MATSOLVERPASTIX "pastix" 1112692d6eeSBarry Smith #define MATSOLVERMATLAB "matlab" 1122692d6eeSBarry Smith #define MATSOLVERPETSC "petsc" 1132692d6eeSBarry Smith #define MATSOLVERBAS "bas" 1149ae82921SPaul Mullowney #define MATSOLVERCUSPARSE "cusparse" 115aa5a9175SDahai Guo #define MATSOLVERBSTRM "bstrm" 116aa5a9175SDahai Guo #define MATSOLVERSBSTRM "sbstrm" 11715767789SHong Zhang #define MATSOLVERELEMENTAL "elemental" 11817f1a0eaSHong Zhang #define MATSOLVERCLIQUE "clique" 11969e15a41SShri Abhyankar #define MATSOLVERKLU "klu" 120c0cdd4a1SDahai Guo 121b24902e0SBarry Smith /*E 122b24902e0SBarry Smith MatFactorType - indicates what type of factorization is requested 123b24902e0SBarry Smith 124b24902e0SBarry Smith Level: beginner 125b24902e0SBarry Smith 126b24902e0SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 127b24902e0SBarry Smith 128c7393fdbSBarry Smith .seealso: MatSolverPackage, MatGetFactor() 129b24902e0SBarry Smith E*/ 130599ef60dSHong Zhang typedef enum {MAT_FACTOR_NONE, MAT_FACTOR_LU, MAT_FACTOR_CHOLESKY, MAT_FACTOR_ILU, MAT_FACTOR_ICC,MAT_FACTOR_ILUDT} MatFactorType; 131014dd563SJed Brown PETSC_EXTERN const char *const MatFactorTypes[]; 132e92e720dSBarry Smith 133014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactor(Mat,const MatSolverPackage,MatFactorType,Mat*); 134014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactorAvailable(Mat,const MatSolverPackage,MatFactorType,PetscBool *); 135014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFactorGetSolverPackage(Mat,const MatSolverPackage*); 136014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactorType(Mat,MatFactorType*); 137*18713533SBarry Smith PETSC_EXTERN PetscErrorCode MatSolverPackageRegister(const MatSolverPackage,const MatType,MatFactorType,PetscErrorCode(*)(Mat,MatFactorType,Mat*)); 13842c9c57cSBarry Smith PETSC_EXTERN PetscErrorCode MatSolverPackageGet(const MatSolverPackage,const MatType,MatFactorType,PetscBool*,PetscBool*,PetscErrorCode (**)(Mat,MatFactorType,Mat*)); 1399989ab13SBarry Smith 140c06d978dSMatthew Knepley /* Logging support */ 1410700a824SBarry Smith #define MAT_FILE_CLASSID 1211216 /* used to indicate matrices in binary files */ 142014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_CLASSID; 143335efc43SPeter Brune PETSC_EXTERN PetscClassId MAT_COLORING_CLASSID; 144014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_FDCOLORING_CLASSID; 145014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_TRANSPOSECOLORING_CLASSID; 146014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_PARTITIONING_CLASSID; 147014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_COARSEN_CLASSID; 148014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_NULLSPACE_CLASSID; 149014dd563SJed Brown PETSC_EXTERN PetscClassId MATMFFD_CLASSID; 150c06d978dSMatthew Knepley 151ceb03754SKris Buschelman /*E 152ceb03754SKris Buschelman MatReuse - Indicates if matrices obtained from a previous call to MatGetSubMatrices() 153d6eff37eSBarry Smith or MatGetSubMatrix() are to be reused to store the new matrix values. For MatConvert() is used to indicate 154d6eff37eSBarry Smith that the input matrix is to be replaced with the converted matrix. 155ceb03754SKris Buschelman 156ceb03754SKris Buschelman Level: beginner 157ceb03754SKris Buschelman 158ceb03754SKris Buschelman Any additions/changes here MUST also be made in include/finclude/petscmat.h 159ceb03754SKris Buschelman 1600c451bc4SBarry Smith .seealso: MatGetSubMatrices(), MatGetSubMatrix(), MatDestroyMatrices(), MatConvert() 161ceb03754SKris Buschelman E*/ 162dfe085dbSJed Brown typedef enum {MAT_INITIAL_MATRIX,MAT_REUSE_MATRIX,MAT_IGNORE_MATRIX} MatReuse; 163ceb03754SKris Buschelman 1645494a064SHong Zhang /*E 1655494a064SHong Zhang MatGetSubMatrixOption - Indicates if matrices obtained from a call to MatGetSubMatrices() 166829201f2SHong Zhang include the matrix values. Currently it is only used by MatGetSeqNonzerostructure(). 1675494a064SHong Zhang 1685494a064SHong Zhang Level: beginner 1695494a064SHong Zhang 170829201f2SHong Zhang .seealso: MatGetSeqNonzerostructure() 1715494a064SHong Zhang E*/ 1725494a064SHong Zhang typedef enum {MAT_DO_NOT_GET_VALUES,MAT_GET_VALUES} MatGetSubMatrixOption; 1735494a064SHong Zhang 174607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatInitializePackage(void); 175c06d978dSMatthew Knepley 176014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreate(MPI_Comm,Mat*); 177014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetSizes(Mat,PetscInt,PetscInt,PetscInt,PetscInt); 17819fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatSetType(Mat,MatType); 179014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetFromOptions(Mat); 180ce1779c8SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatViewFromOptions(Mat A,const char prefix[],const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,prefix,name);} 181607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatRegisterAll(void); 182bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatRegister(const char[],PetscErrorCode(*)(Mat)); 183014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRegisterBaseName(const char[],const char[],const char[]); 184014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetOptionsPrefix(Mat,const char[]); 185014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAppendOptionsPrefix(Mat,const char[]); 186014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOptionsPrefix(Mat,const char*[]); 187f69a0ea3SMatthew Knepley 188014dd563SJed Brown PETSC_EXTERN PetscBool MatRegisterAllCalled; 189140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatList; 190140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatColoringList; 191140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatPartitioningList; 192140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatCoarsenList; 19328988994SBarry Smith 1943b224e63SBarry Smith /*E 195aa6c7ce3SBarry Smith MatStructure - Indicates if two matrices have the same nonzero structure 1963b224e63SBarry Smith 1973b224e63SBarry Smith Level: beginner 1983b224e63SBarry Smith 1993b224e63SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 2003b224e63SBarry Smith 201aa6c7ce3SBarry Smith .seealso: MatCopy(), MatAXPY() 2023b224e63SBarry Smith E*/ 203aa6c7ce3SBarry Smith typedef enum {DIFFERENT_NONZERO_PATTERN,SUBSET_NONZERO_PATTERN,SAME_NONZERO_PATTERN} MatStructure; 2043b224e63SBarry Smith 205014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqDense(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*); 206014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateDense(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*); 207014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 208014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 209014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *); 210014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithSplitArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],PetscInt[],PetscInt[],PetscScalar[],Mat*); 2118d7a6e47SBarry Smith 212014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 213014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 214014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat*); 215d21a29f3SJed Brown 216014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAdj(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscInt[],Mat*); 217014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 218d21a29f3SJed Brown 219014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 220014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPISBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *); 22138f409ebSLisandro Dalcin PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 222014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 2234cce697cSJed Brown PETSC_EXTERN PetscErrorCode MatXAIJSetPreallocation(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscInt[],const PetscInt[]); 224dfb205c3SBarry Smith 225014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateShell(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,void *,Mat*); 226014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateNormal(Mat,Mat*); 227014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateLRC(Mat,Mat,Mat,Mat*); 228014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateIS(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,ISLocalToGlobalMapping,Mat*); 229014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 230014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 231c0cdd4a1SDahai Guo 232014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 233014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 234014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqSBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 235014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPISBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 236c0cdd4a1SDahai Guo 237014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateScatter(MPI_Comm,VecScatter,Mat*); 238014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScatterSetVecScatter(Mat,VecScatter); 239014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScatterGetVecScatter(Mat,VecScatter*); 240014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateBlockMat(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,Mat*); 241014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeAddMat(Mat,Mat); 242014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeMerge(Mat); 243014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateComposite(MPI_Comm,PetscInt,const Mat*,Mat*); 2446d7c1e57SBarry Smith typedef enum {MAT_COMPOSITE_ADDITIVE,MAT_COMPOSITE_MULTIPLICATIVE} MatCompositeType; 245014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeSetType(Mat,MatCompositeType); 2466d7c1e57SBarry Smith 24719fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatCreateFFT(MPI_Comm,PetscInt,const PetscInt[],MatType,Mat*); 248014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqCUFFT(MPI_Comm,PetscInt,const PetscInt[],Mat*); 249dedccee8SHong Zhang 250014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateTranspose(Mat,Mat*); 251014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSubMatrix(Mat,IS,IS,Mat*); 252014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSubMatrixUpdate(Mat,Mat,IS,IS); 253014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateLocalRef(Mat,IS,IS,Mat*); 2541472f72bSBarry Smith 255014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPythonSetType(Mat,const char[]); 2561d6018f0SLisandro Dalcin 257014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetUp(Mat); 258014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroy(Mat*); 259e56f5c9eSBarry Smith PETSC_EXTERN PetscErrorCode MatGetNonzeroState(Mat,PetscObjectState*); 26021c89e3eSBarry Smith 261014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatConjugate(Mat); 262014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRealPart(Mat); 263014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatImaginaryPart(Mat); 264014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetDiagonalBlock(Mat,Mat*); 265014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetTrace(Mat,PetscScalar*); 266713ccfa9SJed Brown PETSC_EXTERN PetscErrorCode MatInvertBlockDiagonal(Mat,const PetscScalar **); 26799cafbc1SBarry Smith 2688ed539a5SBarry Smith /* ------------------------------------------------------------*/ 269014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 270014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 271014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesRow(Mat,PetscInt,const PetscScalar[]); 272014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesRowLocal(Mat,PetscInt,const PetscScalar[]); 273014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBatch(Mat,PetscInt,PetscInt,PetscInt[],const PetscScalar[]); 27473a71a0fSBarry Smith PETSC_EXTERN PetscErrorCode MatSetRandom(Mat,PetscRandom); 27584cb2905SBarry Smith 2762ef4de8bSBarry Smith /*S 2772ef4de8bSBarry Smith MatStencil - Data structure (C struct) for storing information about a single row or 278be479b30SJed Brown column of a matrix as indexed on an associated grid. 279be479b30SJed Brown 280be479b30SJed Brown Fortran usage is different, see MatSetValuesStencil() for details. 2812ef4de8bSBarry Smith 2822ef4de8bSBarry Smith Level: beginner 2832ef4de8bSBarry Smith 2842ef4de8bSBarry Smith Concepts: matrix; linear operator 2852ef4de8bSBarry Smith 286be479b30SJed Brown .seealso: MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockedStencil() 2872ef4de8bSBarry Smith S*/ 288435da068SBarry Smith typedef struct { 289c1ac3661SBarry Smith PetscInt k,j,i,c; 290435da068SBarry Smith } MatStencil; 2912ef4de8bSBarry Smith 292014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode); 293014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode); 294014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt); 295435da068SBarry Smith 296014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetColoring(Mat,ISColoring); 297014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesAdifor(Mat,PetscInt,void*); 2983a7fca6bSBarry Smith 299d91e6319SBarry Smith /*E 300d91e6319SBarry Smith MatAssemblyType - Indicates if the matrix is now to be used, or if you plan 301d91e6319SBarry Smith to continue to add values to it 302d91e6319SBarry Smith 303d91e6319SBarry Smith Level: beginner 304d91e6319SBarry Smith 305d91e6319SBarry Smith .seealso: MatAssemblyBegin(), MatAssemblyEnd() 306d91e6319SBarry Smith E*/ 3076d4a8577SBarry Smith typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType; 308014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssemblyBegin(Mat,MatAssemblyType); 309014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssemblyEnd(Mat,MatAssemblyType); 310014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssembled(Mat,PetscBool *); 3114f9c727eSBarry Smith 3121d73ed98SBarry Smith 31330de9b25SBarry Smith 314d91e6319SBarry Smith /*E 315d91e6319SBarry Smith MatOption - Options that may be set for a matrix and its behavior or storage 316d91e6319SBarry Smith 317d91e6319SBarry Smith Level: beginner 318d91e6319SBarry Smith 3190a835dfdSSatish Balay Any additions/changes here MUST also be made in include/finclude/petscmat.h 320d91e6319SBarry Smith 321382164b0SBarry Smith Developer Notes: Entries that are negative need not be called collectively by all processes. 322382164b0SBarry Smith 323d91e6319SBarry Smith .seealso: MatSetOption() 324d91e6319SBarry Smith E*/ 32592d4d306SBarry Smith typedef enum {MAT_OPTION_MIN = -5, 32692d4d306SBarry Smith MAT_NEW_NONZERO_LOCATION_ERR = -4, 32792d4d306SBarry Smith MAT_UNUSED_NONZERO_LOCATION_ERR = -3, 32892d4d306SBarry Smith MAT_NEW_NONZERO_ALLOCATION_ERR = -2, 32992d4d306SBarry Smith MAT_ROW_ORIENTED = -1, 330382164b0SBarry Smith MAT_SYMMETRIC = 1, 331382164b0SBarry Smith MAT_STRUCTURALLY_SYMMETRIC = 2, 332382164b0SBarry Smith MAT_NEW_DIAGONALS = 3, 333382164b0SBarry Smith MAT_IGNORE_OFF_PROC_ENTRIES = 4, 334382164b0SBarry Smith MAT_USE_HASH_TABLE = 5, 335382164b0SBarry Smith MAT_KEEP_NONZERO_PATTERN = 6, 336382164b0SBarry Smith MAT_IGNORE_ZERO_ENTRIES = 7, 337382164b0SBarry Smith MAT_USE_INODES = 8, 338382164b0SBarry Smith MAT_HERMITIAN = 9, 339382164b0SBarry Smith MAT_SYMMETRY_ETERNAL = 10, 34011e456e1SBarry Smith MAT_DUMMY = 11, 341382164b0SBarry Smith MAT_IGNORE_LOWER_TRIANGULAR = 12, 342382164b0SBarry Smith MAT_ERROR_LOWER_TRIANGULAR = 13, 343382164b0SBarry Smith MAT_GETROW_UPPERTRIANGULAR = 14, 344382164b0SBarry Smith MAT_SPD = 15, 34592d4d306SBarry Smith MAT_NO_OFF_PROC_ZERO_ROWS = 16, 34692d4d306SBarry Smith MAT_NO_OFF_PROC_ENTRIES = 17, 34792d4d306SBarry Smith MAT_NEW_NONZERO_LOCATIONS = 18, 34892d4d306SBarry Smith MAT_OPTION_MAX = 19} MatOption; 349e057df02SPaul Mullowney 350014dd563SJed Brown PETSC_EXTERN const char *MatOptions[]; 351014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetOption(Mat,MatOption,PetscBool ); 35219fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatGetType(Mat,MatType*); 35384cb2905SBarry Smith 354014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]); 355014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 356014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 357014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowUpperTriangular(Mat); 358014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRowUpperTriangular(Mat); 359014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 360014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 361014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnVector(Mat,Vec,PetscInt); 3628c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJGetArray(Mat,PetscScalar *[]); 3638c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJRestoreArray(Mat,PetscScalar *[]); 36421e72a00SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJGetMaxRowNonzeros(Mat,PetscInt*); 365bd04181cSBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJSetValuesLocalFast(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 3668c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatDenseGetArray(Mat,PetscScalar *[]); 3678c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatDenseRestoreArray(Mat,PetscScalar *[]); 368014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBlockSize(Mat,PetscInt *); 369014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSize(Mat,PetscInt); 370014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBlockSizes(Mat,PetscInt *,PetscInt *); 371014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSizes(Mat,PetscInt,PetscInt); 37233d57670SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSizesFromMats(Mat,Mat,Mat); 373014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNThreads(Mat,PetscInt); 374014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNThreads(Mat,PetscInt*); 3751620fd73SBarry Smith 376014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMult(Mat,Vec,Vec); 377014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultDiagonalBlock(Mat,Vec,Vec); 378014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultAdd(Mat,Vec,Vec,Vec); 379014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTranspose(Mat,Vec,Vec); 380014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultHermitianTranspose(Mat,Vec,Vec); 381014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsTranspose(Mat,Mat,PetscReal,PetscBool *); 382014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitianTranspose(Mat,Mat,PetscReal,PetscBool *); 383014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeAdd(Mat,Vec,Vec,Vec); 384014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultHermitianTransposeAdd(Mat,Vec,Vec,Vec); 385014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultConstrained(Mat,Vec,Vec); 386014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeConstrained(Mat,Vec,Vec); 387014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatSolve(Mat,Mat,Mat); 388f9426fe0SMark Adams PETSC_EXTERN PetscErrorCode MatResidual(Mat,Vec,Vec,Vec); 389f5edf698SHong Zhang 390d91e6319SBarry Smith /*E 391d91e6319SBarry Smith MatDuplicateOption - Indicates if a duplicated sparse matrix should have 392d91e6319SBarry Smith its numerical values copied over or just its nonzero structure. 393d91e6319SBarry Smith 394d91e6319SBarry Smith Level: beginner 395d91e6319SBarry Smith 396d91e6319SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 397d91e6319SBarry Smith 39870dcbbb9SBarry Smith $ MAT_SHARE_NONZERO_PATTERN - the i and j arrays in the new matrix will be shared with the original matrix 39970dcbbb9SBarry Smith $ this also triggers the MAT_DO_NOT_COPY_VALUES option. This is used when you 40070dcbbb9SBarry Smith $ have several matrices with the same nonzero pattern. 40170dcbbb9SBarry Smith 402d91e6319SBarry Smith .seealso: MatDuplicate() 403d91e6319SBarry Smith E*/ 40470dcbbb9SBarry Smith typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES,MAT_SHARE_NONZERO_PATTERN} MatDuplicateOption; 4052e8a6d31SBarry Smith 40619fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatConvert(Mat,MatType,MatReuse,Mat*); 407014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDuplicate(Mat,MatDuplicateOption,Mat*); 40894a9d846SBarry Smith 409cb5b572fSBarry Smith 410014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCopy(Mat,Mat,MatStructure); 411014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatView(Mat,PetscViewer); 412014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsSymmetric(Mat,PetscReal,PetscBool *); 413014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsStructurallySymmetric(Mat,PetscBool *); 414014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitian(Mat,PetscReal,PetscBool *); 415014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsSymmetricKnown(Mat,PetscBool *,PetscBool *); 416014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitianKnown(Mat,PetscBool *,PetscBool *); 417014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMissingDiagonal(Mat,PetscBool *,PetscInt *); 418014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLoad(Mat, PetscViewer); 4197b80b807SBarry Smith 4201a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool *); 4211a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool *); 4221a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool *); 4231a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool *); 424d4fbbf0eSBarry Smith 425d91e6319SBarry Smith /*S 426d91e6319SBarry Smith MatInfo - Context of matrix information, used with MatGetInfo() 427d91e6319SBarry Smith 428d91e6319SBarry Smith In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE 429d91e6319SBarry Smith 430d91e6319SBarry Smith Level: intermediate 431d91e6319SBarry Smith 432d91e6319SBarry Smith Concepts: matrix^nonzero information 433d91e6319SBarry Smith 434d9274352SBarry Smith .seealso: MatGetInfo(), MatInfoType 435d91e6319SBarry Smith S*/ 4364e220ebcSLois Curfman McInnes typedef struct { 437b0a32e0cSBarry Smith PetscLogDouble block_size; /* block size */ 438b0a32e0cSBarry Smith PetscLogDouble nz_allocated,nz_used,nz_unneeded; /* number of nonzeros */ 439b0a32e0cSBarry Smith PetscLogDouble memory; /* memory allocated */ 440b0a32e0cSBarry Smith PetscLogDouble assemblies; /* number of matrix assemblies called */ 441b0a32e0cSBarry Smith PetscLogDouble mallocs; /* number of mallocs during MatSetValues() */ 442b0a32e0cSBarry Smith PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */ 443b0a32e0cSBarry Smith PetscLogDouble factor_mallocs; /* number of mallocs during factorization */ 4444e220ebcSLois Curfman McInnes } MatInfo; 4454e220ebcSLois Curfman McInnes 446d9274352SBarry Smith /*E 447d9274352SBarry Smith MatInfoType - Indicates if you want information about the local part of the matrix, 448d9274352SBarry Smith the entire parallel matrix or the maximum over all the local parts. 449d9274352SBarry Smith 450d9274352SBarry Smith Level: beginner 451d9274352SBarry Smith 452d9274352SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 453d9274352SBarry Smith 454d9274352SBarry Smith .seealso: MatGetInfo(), MatInfo 455d9274352SBarry Smith E*/ 4567b80b807SBarry Smith typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType; 457014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetInfo(Mat,MatInfoType,MatInfo*); 458014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetDiagonal(Mat,Vec); 459014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMax(Mat,Vec,PetscInt[]); 460014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMin(Mat,Vec,PetscInt[]); 461014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMaxAbs(Mat,Vec,PetscInt[]); 462014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMinAbs(Mat,Vec,PetscInt[]); 463014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowSum(Mat,Vec); 464014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTranspose(Mat,MatReuse,Mat*); 465014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatHermitianTranspose(Mat,MatReuse,Mat*); 466014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPermute(Mat,IS,IS,Mat *); 467014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalScale(Mat,Vec,Vec); 468014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalSet(Mat,Vec,InsertMode); 469014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatEqual(Mat,Mat,PetscBool *); 470014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultEqual(Mat,Mat,PetscInt,PetscBool *); 471014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultAddEqual(Mat,Mat,PetscInt,PetscBool *); 472014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeEqual(Mat,Mat,PetscInt,PetscBool *); 473014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscBool *); 4747b80b807SBarry Smith 475014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNorm(Mat,NormType,PetscReal *); 476014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnNorms(Mat,NormType,PetscReal *); 477014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroEntries(Mat); 478014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 479014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsIS(Mat,IS,PetscScalar,Vec,Vec); 480014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsStencil(Mat,PetscInt,const MatStencil [],PetscScalar,Vec,Vec); 481014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsStencil(Mat,PetscInt,const MatStencil[],PetscScalar,Vec,Vec); 482014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumns(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 483014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsIS(Mat,IS,PetscScalar,Vec,Vec); 4847b80b807SBarry Smith 485014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSize(Mat,PetscInt*,PetscInt*); 486014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalSize(Mat,PetscInt*,PetscInt*); 487014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRange(Mat,PetscInt*,PetscInt*); 488014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRanges(Mat,const PetscInt**); 489014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRangeColumn(Mat,PetscInt*,PetscInt*); 490014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRangesColumn(Mat,const PetscInt**); 491566876eaSJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipIS(Mat,IS*,IS*); 4925ef9f2a5SBarry Smith 493014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]); 494014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSubMatricesParallel(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]); 495014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroyMatrices(PetscInt,Mat *[]); 496014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSubMatrix(Mat,IS,IS,MatReuse,Mat *); 497014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalSubMatrix(Mat,IS,IS,Mat*); 498014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreLocalSubMatrix(Mat,IS,IS,Mat*); 499014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSeqNonzeroStructure(Mat,Mat*); 500014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroySeqNonzeroStructure(Mat*); 5017b80b807SBarry Smith 502014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJ(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*); 503014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJSymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*); 504014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJNumeric(Mat,Mat); 505014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMat(Mat,MatReuse,Mat*); 506014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*); 507014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,Mat*); 508014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetGhosts(Mat, PetscInt *,const PetscInt *[]); 5098efafbd8SBarry Smith 510014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt); 5117b80b807SBarry Smith 512014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*); 513014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMultSymbolic(Mat,Mat,PetscReal,Mat*); 514014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMultNumeric(Mat,Mat,Mat); 51522440eb1SKris Buschelman 5167bab7c10SHong Zhang PETSC_EXTERN PetscErrorCode MatMatMatMult(Mat,Mat,Mat,MatReuse,PetscReal,Mat*); 5177bab7c10SHong Zhang PETSC_EXTERN PetscErrorCode MatMatMatMultSymbolic(Mat,Mat,Mat,PetscReal,Mat*); 5187bab7c10SHong Zhang PETSC_EXTERN PetscErrorCode MatMatMatMultNumeric(Mat,Mat,Mat,Mat); 5197bab7c10SHong Zhang 520014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*); 521014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAPSymbolic(Mat,Mat,PetscReal,Mat*); 522014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAPNumeric(Mat,Mat,Mat); 523014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARt(Mat,Mat,MatReuse,PetscReal,Mat*); 524014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARtSymbolic(Mat,Mat,PetscReal,Mat*); 525014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARtNumeric(Mat,Mat,Mat); 52622440eb1SKris Buschelman 527014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeMatMult(Mat,Mat,MatReuse,PetscReal,Mat*); 528014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposetMatMultSymbolic(Mat,Mat,PetscReal,Mat*); 529014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposetMatMultNumeric(Mat,Mat,Mat); 530014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMult(Mat,Mat,MatReuse,PetscReal,Mat*); 531014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMultSymbolic(Mat,Mat,PetscReal,Mat*); 532014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMultNumeric(Mat,Mat,Mat); 533bc011b1eSHong Zhang 534014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAXPY(Mat,PetscScalar,Mat,MatStructure); 535014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAYPX(Mat,PetscScalar,Mat,MatStructure); 5361c741599SHong Zhang 537014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScale(Mat,PetscScalar); 538014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShift(Mat,PetscScalar); 5397b80b807SBarry Smith 540014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping); 541014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*); 542a93ff8c4SPeter Brune PETSC_EXTERN PetscErrorCode MatGetLayouts(Mat,PetscLayout*,PetscLayout*); 543014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 544014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsLocalIS(Mat,IS,PetscScalar,Vec,Vec); 545014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 546014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocalIS(Mat,IS,PetscScalar,Vec,Vec); 547014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 548014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 549052efed2SBarry Smith 550014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStashSetInitialSize(Mat,PetscInt,PetscInt); 551014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*); 55290f02eecSBarry Smith 553014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInterpolate(Mat,Vec,Vec); 554014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInterpolateAdd(Mat,Vec,Vec,Vec); 555014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestrict(Mat,Vec,Vec); 5562a7a6963SBarry Smith PETSC_EXTERN PetscErrorCode MatCreateVecs(Mat,Vec*,Vec*); 55784cc2a5aSBarry Smith PETSC_DEPRECATED("Use MatCreateVecs()") PETSC_STATIC_INLINE PetscErrorCode MatGetVecs(Mat mat,Vec *x,Vec *y) {return MatCreateVecs(mat,x,y);} 55853cd1579SHong Zhang PETSC_EXTERN PetscErrorCode MatCreateRedundantMatrix(Mat,PetscInt,MPI_Comm,MatReuse,Mat*); 559014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetMultiProcBlock(Mat,MPI_Comm,MatReuse,Mat*); 560014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFindZeroDiagonals(Mat,IS*); 5613a062f41SBarry Smith PETSC_EXTERN PetscErrorCode MatFindOffBlockDiagonalEntries(Mat,IS*); 5629c8f2541SHong Zhang PETSC_EXTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat(MPI_Comm,Mat,PetscInt,MatReuse,Mat*); 563bd481603SBarry Smith 564bd481603SBarry Smith /*MC 5651620fd73SBarry Smith MatSetValue - Set a single entry into a matrix. 5661620fd73SBarry Smith 5671620fd73SBarry Smith Not collective 5681620fd73SBarry Smith 569a9834a7dSSatish Balay Synopsis: 570a9834a7dSSatish Balay #include <petscmat.h> 571a9834a7dSSatish Balay PetscErrorCode MatSetValue(Mat m,PetscInt row,PetscInt col,PetscScalar value,InsertMode mode) 572a9834a7dSSatish Balay 5731620fd73SBarry Smith Input Parameters: 5741620fd73SBarry Smith + m - the matrix 5751620fd73SBarry Smith . row - the row location of the entry 5761620fd73SBarry Smith . col - the column location of the entry 5771620fd73SBarry Smith . value - the value to insert 5781620fd73SBarry Smith - mode - either INSERT_VALUES or ADD_VALUES 5791620fd73SBarry Smith 5801620fd73SBarry Smith Notes: 5811620fd73SBarry Smith For efficiency one should use MatSetValues() and set several or many 5821620fd73SBarry Smith values simultaneously if possible. 5831620fd73SBarry Smith 5841620fd73SBarry Smith Level: beginner 5851620fd73SBarry Smith 5861620fd73SBarry Smith .seealso: MatSetValues(), MatSetValueLocal() 5871620fd73SBarry Smith M*/ 5881620fd73SBarry 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);} 5891620fd73SBarry Smith 5901620fd73SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatGetValue(Mat v,PetscInt i,PetscInt j,PetscScalar *va) {return MatGetValues(v,1,&i,1,&j,va);} 5911620fd73SBarry Smith 5921620fd73SBarry 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);} 5931620fd73SBarry Smith 5941620fd73SBarry Smith /*MC 5950d2a7ff2SSatish Balay MatPreallocateInitialize - Begins the block of code that will count the number of nonzeros per 596bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 597bd481603SBarry Smith 598bd481603SBarry Smith Synopsis: 599aaa7dc30SBarry Smith #include <petscmat.h> 600c1ac3661SBarry Smith PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz) 601bd481603SBarry Smith 602bd481603SBarry Smith Collective on MPI_Comm 603bd481603SBarry Smith 604bd481603SBarry Smith Input Parameters: 605bd481603SBarry Smith + comm - the communicator that will share the eventually allocated matrix 606859fcb39SBarry Smith . nrows - the number of LOCAL rows in the matrix 607859fcb39SBarry Smith - ncols - the number of LOCAL columns in the matrix 608bd481603SBarry Smith 609bd481603SBarry Smith Output Parameters: 610bd481603SBarry Smith + dnz - the array that will be passed to the matrix preallocation routines 611bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 612bd481603SBarry Smith 613bd481603SBarry Smith Level: intermediate 614bd481603SBarry Smith 615bd481603SBarry Smith Notes: 616a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 617bd481603SBarry Smith 6181d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 619bd481603SBarry Smith 620bd481603SBarry Smith Use MatPreallocateInitializeSymmetric() for symmetric matrices (MPISBAIJ matrices) 621bd481603SBarry Smith 6221620fd73SBarry Smith This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize(). 6231620fd73SBarry Smith 624bd481603SBarry Smith Concepts: preallocation^Matrix 625bd481603SBarry Smith 626d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(), 627d6e23781SBarry Smith MatPreallocateInitializeSymmetric(), MatPreallocateSymmetricSetLocalBlock() 628bd481603SBarry Smith M*/ 629c4f061fbSSatish Balay #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \ 6307c922b88SBarry Smith { \ 631a89bc211SBarry Smith PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__start,__end; \ 6321795a4d1SJed Brown _4_ierr = PetscCalloc2(__nrows,&dnz,__nrows,&onz);CHKERRQ(_4_ierr); \ 6331795a4d1SJed Brown __start = 0; __end = __start; \ 634c1ac3661SBarry Smith _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ctmp;\ 635a89bc211SBarry Smith _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows; 6367c922b88SBarry Smith 637bd481603SBarry Smith /*MC 638bd481603SBarry Smith MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be 639bd481603SBarry Smith inserted using a local number of the rows and columns 640bd481603SBarry Smith 641bd481603SBarry Smith Synopsis: 642aaa7dc30SBarry Smith #include <petscmat.h> 643c1ac3661SBarry Smith PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 644bd481603SBarry Smith 645bd481603SBarry Smith Not Collective 646bd481603SBarry Smith 647bd481603SBarry Smith Input Parameters: 648784ac674SJed Brown + map - the row mapping from local numbering to global numbering 649bd481603SBarry Smith . nrows - the number of rows indicated 6501d73ed98SBarry Smith . rows - the indices of the rows 651784ac674SJed Brown . cmap - the column mapping from local to global numbering 652bd481603SBarry Smith . ncols - the number of columns in the matrix 653bd481603SBarry Smith . cols - the columns indicated 654bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 655bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 656bd481603SBarry Smith 657bd481603SBarry Smith Level: intermediate 658bd481603SBarry Smith 659bd481603SBarry Smith Notes: 660a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 661bd481603SBarry Smith 6621d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 663bd481603SBarry Smith 664bd481603SBarry Smith Concepts: preallocation^Matrix 665bd481603SBarry Smith 666d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(), 667d6e23781SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock() 668bd481603SBarry Smith M*/ 669784ac674SJed Brown #define MatPreallocateSetLocal(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \ 670c4f061fbSSatish Balay {\ 671c1ac3661SBarry Smith PetscInt __l;\ 672784ac674SJed Brown _4_ierr = ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\ 673784ac674SJed Brown _4_ierr = ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\ 674c4f061fbSSatish Balay for (__l=0;__l<nrows;__l++) {\ 675ef66eb69SBarry Smith _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 676c4f061fbSSatish Balay }\ 677c4f061fbSSatish Balay } 678c4f061fbSSatish Balay 679bd481603SBarry Smith /*MC 680d6e23781SBarry Smith MatPreallocateSetLocalBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be 681bd481603SBarry Smith inserted using a local number of the rows and columns 682bd481603SBarry Smith 683bd481603SBarry Smith Synopsis: 684aaa7dc30SBarry Smith #include <petscmat.h> 685d6e23781SBarry Smith PetscErrorCode MatPreallocateSetLocalBlock(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 686d6e23781SBarry Smith 687d6e23781SBarry Smith Not Collective 688d6e23781SBarry Smith 689d6e23781SBarry Smith Input Parameters: 690d6e23781SBarry Smith + map - the row mapping from local numbering to global numbering 691d6e23781SBarry Smith . nrows - the number of rows indicated 692d6e23781SBarry Smith . rows - the indices of the rows 693d6e23781SBarry Smith . cmap - the column mapping from local to global numbering 694d6e23781SBarry Smith . ncols - the number of columns in the matrix 695d6e23781SBarry Smith . cols - the columns indicated 696d6e23781SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 697d6e23781SBarry Smith - ozn - the other array passed to the matrix preallocation routines 698d6e23781SBarry Smith 699d6e23781SBarry Smith Level: intermediate 700d6e23781SBarry Smith 701d6e23781SBarry Smith Notes: 702d6e23781SBarry Smith See Users-Manual: ch_performance for more details. 703d6e23781SBarry Smith 704d6e23781SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 705d6e23781SBarry Smith 706d6e23781SBarry Smith Concepts: preallocation^Matrix 707d6e23781SBarry Smith 708d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(), 709d6e23781SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock() 710d6e23781SBarry Smith M*/ 711d6e23781SBarry Smith #define MatPreallocateSetLocalBlock(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \ 712d6e23781SBarry Smith {\ 713d6e23781SBarry Smith PetscInt __l;\ 714d6e23781SBarry Smith _4_ierr = ISLocalToGlobalMappingApplyBlock(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\ 715d6e23781SBarry Smith _4_ierr = ISLocalToGlobalMappingApplyBlock(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\ 716d6e23781SBarry Smith for (__l=0;__l<nrows;__l++) {\ 717d6e23781SBarry Smith _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 718d6e23781SBarry Smith }\ 719d6e23781SBarry Smith } 720d6e23781SBarry Smith 721d6e23781SBarry Smith /*MC 722d6e23781SBarry Smith MatPreallocateSymmetricSetLocalBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be 723d6e23781SBarry Smith inserted using a local number of the rows and columns 724d6e23781SBarry Smith 725d6e23781SBarry Smith Synopsis: 726d6e23781SBarry Smith #include <petscmat.h> 727d6e23781SBarry Smith PetscErrorCode MatPreallocateSymmetricSetLocalBlock(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 728bd481603SBarry Smith 729bd481603SBarry Smith Not Collective 730bd481603SBarry Smith 731bd481603SBarry Smith Input Parameters: 732bd481603SBarry Smith + map - the mapping between local numbering and global numbering 733bd481603SBarry Smith . nrows - the number of rows indicated 7341d73ed98SBarry Smith . rows - the indices of the rows 735bd481603SBarry Smith . ncols - the number of columns in the matrix 736bd481603SBarry Smith . cols - the columns indicated 737bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 738bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 739bd481603SBarry Smith 740bd481603SBarry Smith Level: intermediate 741bd481603SBarry Smith 742bd481603SBarry Smith Notes: 743a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 744bd481603SBarry Smith 745bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 746bd481603SBarry Smith 747bd481603SBarry Smith Concepts: preallocation^Matrix 748bd481603SBarry Smith 749d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateInitialize(), 750d6e23781SBarry Smith MatPreallocateInitialize(), MatPreallocateSetLocal() 751bd481603SBarry Smith M*/ 752d6e23781SBarry Smith #define MatPreallocateSymmetricSetLocalBlock(map,nrows,rows,ncols,cols,dnz,onz) 0;\ 753d3d32019SBarry Smith {\ 754c1ac3661SBarry Smith PetscInt __l;\ 755d6e23781SBarry Smith _4_ierr = ISLocalToGlobalMappingApplyBlock(map,nrows,rows,rows);CHKERRQ(_4_ierr);\ 756d6e23781SBarry Smith _4_ierr = ISLocalToGlobalMappingApplyBlock(map,ncols,cols,cols);CHKERRQ(_4_ierr);\ 757d3d32019SBarry Smith for (__l=0;__l<nrows;__l++) {\ 758d6e23781SBarry Smith _4_ierr = MatPreallocateSymmetricSetBlock((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 759d3d32019SBarry Smith }\ 760d3d32019SBarry Smith } 761bd481603SBarry Smith /*MC 762bd481603SBarry Smith MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be 763bd481603SBarry Smith inserted using a local number of the rows and columns 764bd481603SBarry Smith 765bd481603SBarry Smith Synopsis: 766aaa7dc30SBarry Smith #include <petscmat.h> 767c1ac3661SBarry Smith PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 768bd481603SBarry Smith 769bd481603SBarry Smith Not Collective 770bd481603SBarry Smith 771bd481603SBarry Smith Input Parameters: 77264075487SBarry Smith + row - the row 773bd481603SBarry Smith . ncols - the number of columns in the matrix 774a50f8bf6SHong Zhang - cols - the columns indicated 775a50f8bf6SHong Zhang 776a50f8bf6SHong Zhang Output Parameters: 777a50f8bf6SHong Zhang + dnz - the array that will be passed to the matrix preallocation routines 778bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 779bd481603SBarry Smith 780bd481603SBarry Smith Level: intermediate 781bd481603SBarry Smith 782bd481603SBarry Smith Notes: 783a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 784bd481603SBarry Smith 785bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 786bd481603SBarry Smith 7871620fd73SBarry Smith This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize(). 7881620fd73SBarry Smith 789bd481603SBarry Smith Concepts: preallocation^Matrix 790bd481603SBarry Smith 791d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(), 792d6e23781SBarry Smith MatPreallocateInitialize(), MatPreallocateSetLocal() 793bd481603SBarry Smith M*/ 794c4f061fbSSatish Balay #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\ 795c1ac3661SBarry Smith { PetscInt __i; \ 796e32f2f54SBarry 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);\ 797e32f2f54SBarry 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);\ 7987c922b88SBarry Smith for (__i=0; __i<nc; __i++) {\ 79964075487SBarry Smith if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \ 8007cc688d7SBarry Smith else dnz[row - __rstart]++;\ 8017c922b88SBarry Smith }\ 8027c922b88SBarry Smith } 8037c922b88SBarry Smith 804bd481603SBarry Smith /*MC 805d6e23781SBarry Smith MatPreallocateSymmetricSetBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be 806bd481603SBarry Smith inserted using a local number of the rows and columns 807bd481603SBarry Smith 808bd481603SBarry Smith Synopsis: 809aaa7dc30SBarry Smith #include <petscmat.h> 810d6e23781SBarry Smith PetscErrorCode MatPreallocateSymmetricSetBlock(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 811bd481603SBarry Smith 812bd481603SBarry Smith Not Collective 813bd481603SBarry Smith 814bd481603SBarry Smith Input Parameters: 815bd481603SBarry Smith + nrows - the number of rows indicated 8161d73ed98SBarry Smith . rows - the indices of the rows 817bd481603SBarry Smith . ncols - the number of columns in the matrix 818bd481603SBarry Smith . cols - the columns indicated 819bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 820bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 821bd481603SBarry Smith 822bd481603SBarry Smith Level: intermediate 823bd481603SBarry Smith 824bd481603SBarry Smith Notes: 825a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 826bd481603SBarry Smith 827bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 828bd481603SBarry Smith 8291620fd73SBarry Smith This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize(). 8301620fd73SBarry Smith 831bd481603SBarry Smith Concepts: preallocation^Matrix 832bd481603SBarry Smith 833d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateInitialize(), 834d6e23781SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocal() 835bd481603SBarry Smith M*/ 836d6e23781SBarry Smith #define MatPreallocateSymmetricSetBlock(row,nc,cols,dnz,onz) 0;\ 837c1ac3661SBarry Smith { PetscInt __i; \ 838d3d32019SBarry Smith for (__i=0; __i<nc; __i++) {\ 839d3d32019SBarry Smith if (cols[__i] >= __end) onz[row - __rstart]++; \ 840d3d32019SBarry Smith else if (cols[__i] >= row) dnz[row - __rstart]++;\ 841d3d32019SBarry Smith }\ 842d3d32019SBarry Smith } 843d3d32019SBarry Smith 844bd481603SBarry Smith /*MC 84516371a99SBarry Smith MatPreallocateLocation - An alternative to MatPreallocationSet() that puts the nonzero locations into the matrix if it exists 84616371a99SBarry Smith 84716371a99SBarry Smith Synopsis: 848aaa7dc30SBarry Smith #include <petscmat.h> 84916371a99SBarry Smith PetscErrorCode MatPreallocateLocations(Mat A,PetscInt row,PetscInt ncols,PetscInt *cols,PetscInt *dnz,PetscInt *onz) 85016371a99SBarry Smith 85116371a99SBarry Smith Not Collective 85216371a99SBarry Smith 85316371a99SBarry Smith Input Parameters: 85416371a99SBarry Smith . A - matrix 85516371a99SBarry Smith . row - row where values exist (must be local to this process) 85616371a99SBarry Smith . ncols - number of columns 85716371a99SBarry Smith . cols - columns with nonzeros 85816371a99SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 85916371a99SBarry Smith - ozn - the other array passed to the matrix preallocation routines 86016371a99SBarry Smith 86116371a99SBarry Smith Level: intermediate 86216371a99SBarry Smith 86316371a99SBarry Smith Notes: 864a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 86516371a99SBarry Smith 86616371a99SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 86716371a99SBarry Smith 86816371a99SBarry Smith This is a MACRO not a function because it uses a bunch of variables private to the MatPreallocation.... routines. 86916371a99SBarry Smith 87016371a99SBarry Smith Concepts: preallocation^Matrix 87116371a99SBarry Smith 872d6e23781SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(), 873d6e23781SBarry Smith MatPreallocateSymmetricSetLocalBlock() 87416371a99SBarry Smith M*/ 8750298fd71SBarry Smith #define MatPreallocateLocation(A,row,ncols,cols,dnz,onz) 0;if (A) {ierr = MatSetValues(A,1,&row,ncols,cols,NULL,INSERT_VALUES);CHKERRQ(ierr);} else {ierr = MatPreallocateSet(row,ncols,cols,dnz,onz);CHKERRQ(ierr);} 87616371a99SBarry Smith 87716371a99SBarry Smith 87816371a99SBarry Smith /*MC 8790d2a7ff2SSatish Balay MatPreallocateFinalize - Ends the block of code that will count the number of nonzeros per 880bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 881bd481603SBarry Smith 882bd481603SBarry Smith Synopsis: 883aaa7dc30SBarry Smith #include <petscmat.h> 884c1ac3661SBarry Smith PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz) 885bd481603SBarry Smith 886bd481603SBarry Smith Collective on MPI_Comm 887bd481603SBarry Smith 888bd481603SBarry Smith Input Parameters: 88916371a99SBarry Smith + dnz - the array that was be passed to the matrix preallocation routines 890bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 891bd481603SBarry Smith 892bd481603SBarry Smith Level: intermediate 893bd481603SBarry Smith 894bd481603SBarry Smith Notes: 895a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 896bd481603SBarry Smith 897bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 898bd481603SBarry Smith 8991620fd73SBarry Smith This is a MACRO not a function because it closes the { started in MatPreallocateInitialize(). 9001620fd73SBarry Smith 901bd481603SBarry Smith Concepts: preallocation^Matrix 902bd481603SBarry Smith 903d6e23781SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(), 904d6e23781SBarry Smith MatPreallocateSymmetricSetLocalBlock() 905bd481603SBarry Smith M*/ 906a89bc211SBarry Smith #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree2(dnz,onz);CHKERRQ(_4_ierr);} 9077c922b88SBarry Smith 9087b80b807SBarry Smith /* Routines unique to particular data structures */ 909014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellGetContext(Mat,void *); 910698d4c6aSKris Buschelman 911014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInodeAdjustForInodes(Mat,IS*,IS*); 912014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *); 913698d4c6aSKris Buschelman 914014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetColumnIndices(Mat,PetscInt[]); 915014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetColumnIndices(Mat,PetscInt[]); 916014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 917014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 918014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 919014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJFromTriple(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*,PetscInt,PetscBool); 9207b80b807SBarry Smith 921a96a251dSBarry Smith #define MAT_SKIP_ALLOCATION -4 922a96a251dSBarry Smith 923014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]); 924014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]); 925014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]); 926273d9f13SBarry Smith 927014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 928014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 929014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 930014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []); 931014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 932014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]); 933014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 934014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]); 935014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIDenseSetPreallocation(Mat,PetscScalar[]); 936014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqDenseSetPreallocation(Mat,PetscScalar[]); 9379230625dSJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,const PetscInt*[]); 9389230625dSJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,const PetscInt*[]); 939014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAdjCreateNonemptySubcommMat(Mat,Mat*); 940273d9f13SBarry Smith 9412e1947a5SStefano Zampini PETSC_EXTERN PetscErrorCode MatISSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 942014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqDenseSetLDA(Mat,PetscInt); 943014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDenseGetLocalMatrix(Mat,Mat*); 9441b807ce4Svictorle 945014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStoreValues(Mat); 946014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRetrieveValues(Mat); 9472e8a6d31SBarry Smith 948014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDAADSetCtx(Mat,void*); 9493a7fca6bSBarry Smith 950014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFindNonzeroRows(Mat,IS*); 9517b80b807SBarry Smith /* 9527b80b807SBarry Smith These routines are not usually accessed directly, rather solving is 95394b7f48cSBarry Smith done through the KSP and PC interfaces. 9547b80b807SBarry Smith */ 9557b80b807SBarry Smith 95676bdecfbSBarry Smith /*J 9578f6c3df8SBarry Smith MatOrderingType - String with the name of a PETSc matrix ordering 958d9274352SBarry Smith 959d9274352SBarry Smith Level: beginner 960d9274352SBarry Smith 961d9274352SBarry Smith .seealso: MatGetOrdering() 96276bdecfbSBarry Smith J*/ 96319fd82e9SBarry Smith typedef const char* MatOrderingType; 9642692d6eeSBarry Smith #define MATORDERINGNATURAL "natural" 9652692d6eeSBarry Smith #define MATORDERINGND "nd" 9662692d6eeSBarry Smith #define MATORDERING1WD "1wd" 9672692d6eeSBarry Smith #define MATORDERINGRCM "rcm" 9682692d6eeSBarry Smith #define MATORDERINGQMD "qmd" 9692692d6eeSBarry Smith #define MATORDERINGROWLENGTH "rowlength" 970510184a7SJed Brown #define MATORDERINGWBM "wbm" 971fc1bef75SJed Brown #define MATORDERINGSPECTRAL "spectral" 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*); 975140e18c1SBarry Smith PETSC_EXTERN PetscErrorCode MatGetOrderingList(PetscFunctionList*); 976bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatOrderingRegister(const char[],PetscErrorCode(*)(Mat,MatOrderingType,IS*,IS*)); 977607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatOrderingRegisterAll(void); 978014dd563SJed Brown PETSC_EXTERN PetscBool MatOrderingRegisterAllCalled; 979140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatOrderingList; 980d4fbbf0eSBarry Smith 981014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS); 982fc1bef75SJed Brown PETSC_EXTERN PetscErrorCode MatCreateLaplacian(Mat,PetscReal,PetscBool,Mat*); 983a2ce50c7SBarry Smith 984d91e6319SBarry Smith /*S 985d90ac83dSHong Zhang MatFactorShiftType - Numeric Shift. 986d90ac83dSHong Zhang 987d90ac83dSHong Zhang Level: beginner 988d90ac83dSHong Zhang 989d90ac83dSHong Zhang S*/ 990d90ac83dSHong Zhang typedef enum {MAT_SHIFT_NONE,MAT_SHIFT_NONZERO,MAT_SHIFT_POSITIVE_DEFINITE,MAT_SHIFT_INBLOCKS} MatFactorShiftType; 9916a6fc655SJed Brown PETSC_EXTERN const char *const MatFactorShiftTypes[]; 9925e9742b9SJed Brown PETSC_EXTERN const char *const MatFactorShiftTypesDetail[]; 993d90ac83dSHong Zhang 994d90ac83dSHong Zhang /*S 9952401956bSBarry Smith MatFactorInfo - Data passed into the matrix factorization routines 996b00f7748SHong Zhang 99761cad860SBarry Smith In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE, that is use 99861cad860SBarry Smith $ MatFactorInfo info(MAT_FACTORINFO_SIZE) 999b00f7748SHong Zhang 100015e8a5b3SHong Zhang Notes: These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC. 1001b00f7748SHong Zhang 100261cad860SBarry Smith You can use MatFactorInfoInitialize() to set default values. 100361cad860SBarry Smith 1004b00f7748SHong Zhang Level: developer 1005b00f7748SHong Zhang 1006d7d82daaSBarry Smith .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(), 1007d7d82daaSBarry Smith MatFactorInfoInitialize() 1008b00f7748SHong Zhang 1009b00f7748SHong Zhang S*/ 1010b00f7748SHong Zhang typedef struct { 101115e8a5b3SHong Zhang PetscReal diagonal_fill; /* force diagonal to fill in if initially not filled */ 101285317021SBarry Smith PetscReal usedt; 101315e8a5b3SHong Zhang PetscReal dt; /* drop tolerance */ 1014b00f7748SHong Zhang PetscReal dtcol; /* tolerance for pivoting */ 101515e8a5b3SHong Zhang PetscReal dtcount; /* maximum nonzeros to be allowed per row */ 101667e28bfeSBarry Smith PetscReal fill; /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */ 1017348344bbSBarry Smith PetscReal levels; /* ICC/ILU(levels) */ 1018bcd9e38bSBarry Smith PetscReal pivotinblocks; /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0 1019bcd9e38bSBarry Smith factorization may be faster if do not pivot */ 102015e8a5b3SHong Zhang PetscReal zeropivot; /* pivot is called zero if less than this */ 1021f4db908eSBarry Smith PetscReal shifttype; /* type of shift added to matrix factor to prevent zero pivots */ 1022f4db908eSBarry Smith PetscReal shiftamount; /* how large the shift is */ 102315e8a5b3SHong Zhang } MatFactorInfo; 1024ffa6d0a5SLois Curfman McInnes 1025014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFactorInfoInitialize(MatFactorInfo*); 1026014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCholeskyFactor(Mat,IS,const MatFactorInfo*); 1027014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCholeskyFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*); 1028014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCholeskyFactorNumeric(Mat,Mat,const MatFactorInfo*); 1029014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLUFactor(Mat,IS,IS,const MatFactorInfo*); 1030014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatILUFactor(Mat,IS,IS,const MatFactorInfo*); 1031014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*); 1032014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatILUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*); 1033014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatICCFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*); 1034014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatICCFactor(Mat,IS,const MatFactorInfo*); 1035014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLUFactorNumeric(Mat,Mat,const MatFactorInfo*); 1036014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*); 1037014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolve(Mat,Vec,Vec); 1038014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatForwardSolve(Mat,Vec,Vec); 1039014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatBackwardSolve(Mat,Vec,Vec); 1040014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveAdd(Mat,Vec,Vec,Vec); 1041014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveTranspose(Mat,Vec,Vec); 1042014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveTransposeAdd(Mat,Vec,Vec,Vec); 1043014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolves(Mat,Vecs,Vecs); 10448ed539a5SBarry Smith 1045014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetUnfactored(Mat); 1046bb5a7306SBarry Smith 1047d91e6319SBarry Smith /*E 1048d91e6319SBarry Smith MatSORType - What type of (S)SOR to perform 1049bb1eb677SSatish Balay 1050d91e6319SBarry Smith Level: beginner 1051d91e6319SBarry Smith 1052d9274352SBarry Smith May be bitwise ORd together 1053d9274352SBarry Smith 1054d91e6319SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 1055d91e6319SBarry Smith 10564e7234bfSBarry Smith MatSORType may be bitwise ORd together, so do not change the numbers 10574e7234bfSBarry Smith 105841f059aeSBarry Smith .seealso: MatSOR() 1059d91e6319SBarry Smith E*/ 1060ee50ffe9SBarry Smith typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3, 1061ee50ffe9SBarry Smith SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8, 1062ee50ffe9SBarry Smith SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16, 106384cb2905SBarry Smith SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType; 1064014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSOR(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec); 10658ed539a5SBarry Smith 1066d4fbbf0eSBarry Smith /* 1067639f9d9dSBarry Smith These routines are for efficiently computing Jacobians via finite differences. 1068639f9d9dSBarry Smith */ 1069b12f92e5SBarry Smith 10707086a01eSPeter Brune /*S 10717086a01eSPeter Brune MatColoring - Object for managing the coloring of matrices. 10727086a01eSPeter Brune 10737086a01eSPeter Brune Level: beginner 10747086a01eSPeter Brune 10757086a01eSPeter Brune Concepts: matrix, coloring 10767086a01eSPeter Brune 10777086a01eSPeter Brune .seealso: MatFDColoringCreate() ISColoring MatFDColoring 10787086a01eSPeter Brune S*/ 10797086a01eSPeter Brune typedef struct _p_MatColoring* MatColoring; 108076bdecfbSBarry Smith /*J 10818f6c3df8SBarry Smith MatColoringType - String with the name of a PETSc matrix coloring 1082d9274352SBarry Smith 1083d9274352SBarry Smith Level: beginner 1084d9274352SBarry Smith 10857086a01eSPeter Brune .seealso: MatColoringSetType(), MatColoring 108676bdecfbSBarry Smith J*/ 10877086a01eSPeter Brune 108819fd82e9SBarry Smith typedef const char* MatColoringType; 10895567d71aSPeter Brune #define MATCOLORINGJP "jp" 1090b9acb617SPeter Brune #define MATCOLORINGPOWER "power" 10912692d6eeSBarry Smith #define MATCOLORINGNATURAL "natural" 10922692d6eeSBarry Smith #define MATCOLORINGSL "sl" 10932692d6eeSBarry Smith #define MATCOLORINGLF "lf" 10942692d6eeSBarry Smith #define MATCOLORINGID "id" 10958121bdceSPeter Brune #define MATCOLORINGGREEDY "greedy" 1096b12f92e5SBarry Smith 10978ac6da0aSPeter Brune /*E 10988ac6da0aSPeter Brune MatColoringWeightType - Type of weight scheme 10998ac6da0aSPeter Brune 11008ac6da0aSPeter Brune Not Collective 11018ac6da0aSPeter Brune 11028ac6da0aSPeter Brune + MAT_COLORING_RANDOM - Random weights 11038ac6da0aSPeter Brune . MAT_COLORING_LEXICAL - Lexical weighting based upon global numbering. 11048ac6da0aSPeter Brune - MAT_COLORING_LF - Last-first weighting. 11058ac6da0aSPeter Brune 11068ac6da0aSPeter Brune Level: intermediate 11078ac6da0aSPeter Brune 11088ac6da0aSPeter Brune Any additions/changes here MUST also be made in include/finclude/petscmat.h 11098ac6da0aSPeter Brune 11108ac6da0aSPeter Brune .seealso: MatCUSPSetFormat(), MatCUSPFormatOperation 11118ac6da0aSPeter Brune E*/ 1112301c140bSPeter Brune typedef enum {MAT_COLORING_WEIGHT_RANDOM,MAT_COLORING_WEIGHT_LEXICAL,MAT_COLORING_WEIGHT_LF,MAT_COLORING_WEIGHT_SL} MatColoringWeightType; 11138ac6da0aSPeter Brune 1114335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringCreate(Mat,MatColoring*); 1115d7f2a307SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringGetDegrees(Mat,PetscInt,PetscInt*); 1116335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringDestroy(MatColoring*); 1117335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringView(MatColoring,PetscViewer); 1118335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetType(MatColoring,MatColoringType); 1119335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetFromOptions(MatColoring); 1120335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetDistance(MatColoring,PetscInt); 1121335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringGetDistance(MatColoring,PetscInt*); 1122335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetMaxColors(MatColoring,PetscInt); 1123335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringGetMaxColors(MatColoring,PetscInt*); 1124335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringApply(MatColoring,ISColoring*); 1125607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatColoringRegisterAll(void); 1126335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringRegister(const char[],PetscErrorCode(*)(MatColoring)); 1127335efc43SPeter Brune PETSC_EXTERN PetscBool MatColoringRegisterAllCalled; 1128014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*); 11298ac6da0aSPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetWeightType(MatColoring,MatColoringWeightType); 1130c29971fcSPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetWeights(MatColoring,PetscReal*,PetscInt*); 11318ac6da0aSPeter Brune PETSC_EXTERN PetscErrorCode MatColoringCreateWeights(MatColoring,PetscReal **,PetscInt **lperm); 1132639f9d9dSBarry Smith 1133d9274352SBarry Smith /*S 1134d9274352SBarry Smith MatFDColoring - Object for computing a sparse Jacobian via finite differences 1135d9274352SBarry Smith and coloring 1136639f9d9dSBarry Smith 1137d9274352SBarry Smith Level: beginner 1138d9274352SBarry Smith 1139d9274352SBarry Smith Concepts: coloring, sparse Jacobian, finite differences 1140d9274352SBarry Smith 1141d9274352SBarry Smith .seealso: MatFDColoringCreate() 1142d9274352SBarry Smith S*/ 1143e2a1c21fSSatish Balay typedef struct _p_MatFDColoring* MatFDColoring; 1144639f9d9dSBarry Smith 1145014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringCreate(Mat,ISColoring,MatFDColoring *); 1146014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringDestroy(MatFDColoring*); 1147014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringView(MatFDColoring,PetscViewer); 1148014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*); 1149014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**); 1150014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal); 1151014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetFromOptions(MatFDColoring); 1152d1e9a80fSBarry Smith PETSC_EXTERN PetscErrorCode MatFDColoringApply(Mat,MatFDColoring,Vec,void *); 1153014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetF(MatFDColoring,Vec); 1154014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,PetscInt*[]); 1155f86b9fbaSHong Zhang PETSC_EXTERN PetscErrorCode MatFDColoringSetUp(Mat,ISColoring,MatFDColoring); 1156f86b9fbaSHong Zhang PETSC_EXTERN PetscErrorCode MatFDColoringSetBlockSize(MatFDColoring,PetscInt,PetscInt); 1157f86b9fbaSHong Zhang 1158b1683b59SHong Zhang 1159b1683b59SHong Zhang /*S 1160b9af6bddSHong Zhang MatTransposeColoring - Object for computing a sparse matrix product C=A*B^T via coloring 1161b1683b59SHong Zhang 1162b1683b59SHong Zhang Level: beginner 1163b1683b59SHong Zhang 1164b1683b59SHong Zhang Concepts: coloring, sparse matrix product 1165b1683b59SHong Zhang 1166b9af6bddSHong Zhang .seealso: MatTransposeColoringCreate() 1167b1683b59SHong Zhang S*/ 1168b9af6bddSHong Zhang typedef struct _p_MatTransposeColoring* MatTransposeColoring; 1169b1683b59SHong Zhang 1170014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeColoringCreate(Mat,ISColoring,MatTransposeColoring *); 1171014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransColoringApplySpToDen(MatTransposeColoring,Mat,Mat); 1172014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransColoringApplyDenToSp(MatTransposeColoring,Mat,Mat); 1173014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeColoringDestroy(MatTransposeColoring*); 1174b1683b59SHong Zhang 1175639f9d9dSBarry Smith /* 11760752156aSBarry Smith These routines are for partitioning matrices: currently used only 11773eda8832SBarry Smith for adjacency matrix, MatCreateMPIAdj(). 11780752156aSBarry Smith */ 1179ca161407SBarry Smith 1180d9274352SBarry Smith /*S 1181d9274352SBarry Smith MatPartitioning - Object for managing the partitioning of a matrix or graph 1182d9274352SBarry Smith 1183d9274352SBarry Smith Level: beginner 1184d9274352SBarry Smith 1185d9274352SBarry Smith Concepts: partitioning 1186d9274352SBarry Smith 1187743a1026Svictorle .seealso: MatPartitioningCreate(), MatPartitioningType 1188d9274352SBarry Smith S*/ 118991e9ee9fSBarry Smith typedef struct _p_MatPartitioning* MatPartitioning; 1190d9274352SBarry Smith 119176bdecfbSBarry Smith /*J 11928f6c3df8SBarry Smith MatPartitioningType - String with the name of a PETSc matrix partitioning 1193d9274352SBarry Smith 1194d9274352SBarry Smith Level: beginner 11950d04baf8SBarry Smith dm 1196b547a13bSvictorle .seealso: MatPartitioningCreate(), MatPartitioning 119776bdecfbSBarry Smith J*/ 119819fd82e9SBarry Smith typedef const char* MatPartitioningType; 11992692d6eeSBarry Smith #define MATPARTITIONINGCURRENT "current" 12002692d6eeSBarry Smith #define MATPARTITIONINGSQUARE "square" 12012692d6eeSBarry Smith #define MATPARTITIONINGPARMETIS "parmetis" 12022692d6eeSBarry Smith #define MATPARTITIONINGCHACO "chaco" 12032692d6eeSBarry Smith #define MATPARTITIONINGPARTY "party" 120450d91057SBarry Smith #define MATPARTITIONINGPTSCOTCH "ptscotch" 120571306c60Sjroman 1206ca161407SBarry Smith 1207014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningCreate(MPI_Comm,MatPartitioning*); 120819fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningSetType(MatPartitioning,MatPartitioningType); 1209014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetNParts(MatPartitioning,PetscInt); 1210014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetAdjacency(MatPartitioning,Mat); 1211014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]); 1212014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []); 1213014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningApply(MatPartitioning,IS*); 1214014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningDestroy(MatPartitioning*); 12152aabb6bbSBarry Smith 1216bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningRegister(const char[],PetscErrorCode (*)(MatPartitioning)); 121730de9b25SBarry Smith 1218014dd563SJed Brown PETSC_EXTERN PetscBool MatPartitioningRegisterAllCalled; 1219f1144a30SSatish Balay 1220607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningRegisterAll(void); 12212bad1931SBarry Smith 1222014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningView(MatPartitioning,PetscViewer); 1223014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetFromOptions(MatPartitioning); 122419fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningGetType(MatPartitioning,MatPartitioningType*); 1225ca161407SBarry Smith 1226014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetCoarseSequential(MatPartitioning); 1227014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningParmetisGetEdgeCut(MatPartitioning, PetscInt *); 12280752156aSBarry Smith 1229b6956c12SJose Roman typedef enum { MP_CHACO_MULTILEVEL=1,MP_CHACO_SPECTRAL=2,MP_CHACO_LINEAR=4,MP_CHACO_RANDOM=5,MP_CHACO_SCATTERED=6 } MPChacoGlobalType; 12306a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoGlobalTypes[]; 1231b6956c12SJose Roman typedef enum { MP_CHACO_KERNIGHAN=1,MP_CHACO_NONE=2 } MPChacoLocalType; 12326a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoLocalTypes[]; 1233b6956c12SJose Roman typedef enum { MP_CHACO_LANCZOS=0,MP_CHACO_RQI=1 } MPChacoEigenType; 12346a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoEigenTypes[]; 1235b6956c12SJose Roman 1236014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetGlobal(MatPartitioning,MPChacoGlobalType); 1237014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetGlobal(MatPartitioning,MPChacoGlobalType*); 1238014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetLocal(MatPartitioning,MPChacoLocalType); 1239014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetLocal(MatPartitioning,MPChacoLocalType*); 1240014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal); 1241014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType); 1242014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenSolver(MatPartitioning,MPChacoEigenType*); 1243014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenTol(MatPartitioning,PetscReal); 1244014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenTol(MatPartitioning,PetscReal*); 1245014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenNumber(MatPartitioning,PetscInt); 1246014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenNumber(MatPartitioning,PetscInt*); 124771306c60Sjroman 124871306c60Sjroman #define MP_PARTY_OPT "opt" 124971306c60Sjroman #define MP_PARTY_LIN "lin" 125071306c60Sjroman #define MP_PARTY_SCA "sca" 125171306c60Sjroman #define MP_PARTY_RAN "ran" 125271306c60Sjroman #define MP_PARTY_GBF "gbf" 125371306c60Sjroman #define MP_PARTY_GCF "gcf" 125471306c60Sjroman #define MP_PARTY_BUB "bub" 125571306c60Sjroman #define MP_PARTY_DEF "def" 1256014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetGlobal(MatPartitioning,const char*); 125771306c60Sjroman #define MP_PARTY_HELPFUL_SETS "hs" 125871306c60Sjroman #define MP_PARTY_KERNIGHAN_LIN "kl" 125971306c60Sjroman #define MP_PARTY_NONE "no" 1260014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetLocal(MatPartitioning,const char*); 1261014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal); 1262014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetBipart(MatPartitioning,PetscBool); 1263014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscBool); 126471306c60Sjroman 126550d91057SBarry Smith typedef enum { MP_PTSCOTCH_QUALITY,MP_PTSCOTCH_SPEED,MP_PTSCOTCH_BALANCE,MP_PTSCOTCH_SAFETY,MP_PTSCOTCH_SCALABILITY } MPPTScotchStrategyType; 12666a6fc655SJed Brown PETSC_EXTERN const char *const MPPTScotchStrategyTypes[]; 1267e0f1cffaSJose Roman 1268014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetImbalance(MatPartitioning,PetscReal); 1269014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetImbalance(MatPartitioning,PetscReal*); 1270014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetStrategy(MatPartitioning,MPPTScotchStrategyType); 1271014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetStrategy(MatPartitioning,MPPTScotchStrategyType*); 127271306c60Sjroman 1273b43b03e9SMark F. Adams /* 1274b43b03e9SMark F. Adams These routines are for coarsening matrices: 1275b43b03e9SMark F. Adams */ 1276b43b03e9SMark F. Adams 1277b43b03e9SMark F. Adams /*S 1278b43b03e9SMark F. Adams MatCoarsen - Object for managing the coarsening of a graph (symmetric matrix) 1279b43b03e9SMark F. Adams 1280b43b03e9SMark F. Adams Level: beginner 1281b43b03e9SMark F. Adams 1282b43b03e9SMark F. Adams Concepts: coarsen 1283b43b03e9SMark F. Adams 1284b43b03e9SMark F. Adams .seealso: MatCoarsenCreate), MatCoarsenType 1285b43b03e9SMark F. Adams S*/ 1286b43b03e9SMark F. Adams typedef struct _p_MatCoarsen* MatCoarsen; 1287b43b03e9SMark F. Adams 1288b43b03e9SMark F. Adams /*J 12898f6c3df8SBarry Smith MatCoarsenType - String with the name of a PETSc matrix coarsen 1290b43b03e9SMark F. Adams 1291b43b03e9SMark F. Adams Level: beginner 12928f6c3df8SBarry Smith 1293b43b03e9SMark F. Adams .seealso: MatCoarsenCreate(), MatCoarsen 1294b43b03e9SMark F. Adams J*/ 129519fd82e9SBarry Smith typedef const char* MatCoarsenType; 1296b43b03e9SMark F. Adams #define MATCOARSENMIS "mis" 12979057884aSMark F. Adams #define MATCOARSENHEM "hem" 1298b43b03e9SMark F. Adams 12990cbbd2e1SMark F. Adams /* linked list for aggregates */ 130041b27cdeSMark F. Adams typedef struct _PetscCDIntNd{ 130141b27cdeSMark F. Adams struct _PetscCDIntNd *next; 13020cbbd2e1SMark F. Adams PetscInt gid; 130341b27cdeSMark F. Adams }PetscCDIntNd; 13040cbbd2e1SMark F. Adams 13050cbbd2e1SMark F. Adams /* only used by node pool */ 130641b27cdeSMark F. Adams typedef struct _PetscCDArrNd{ 130741b27cdeSMark F. Adams struct _PetscCDArrNd *next; 130841b27cdeSMark F. Adams struct _PetscCDIntNd *array; 130941b27cdeSMark F. Adams }PetscCDArrNd; 13100cbbd2e1SMark F. Adams 13110cbbd2e1SMark F. Adams typedef struct _PetscCoarsenData{ 131282c86c8fSBarry Smith PetscCDArrNd pool_list; /* node pool */ 131341b27cdeSMark F. Adams PetscCDIntNd *new_node; 13140cbbd2e1SMark F. Adams PetscInt new_left; 13150cbbd2e1SMark F. Adams PetscInt chk_sz; 131641b27cdeSMark F. Adams PetscCDIntNd *extra_nodes; 131782c86c8fSBarry Smith PetscCDIntNd **array; /* Array of lists */ 13180cbbd2e1SMark F. Adams PetscInt size; 131982c86c8fSBarry Smith Mat mat; /* cache a Mat for communication data */ 13200cbbd2e1SMark F. Adams }PetscCoarsenData; 13210cbbd2e1SMark F. Adams 1322014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenCreate(MPI_Comm,MatCoarsen*); 132319fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenSetType(MatCoarsen,MatCoarsenType); 1324014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetAdjacency(MatCoarsen,Mat); 1325014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetGreedyOrdering(MatCoarsen,const IS); 1326014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetStrictAggs(MatCoarsen,PetscBool); 1327014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetVerbose(MatCoarsen,PetscInt); 1328014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenGetData( MatCoarsen, PetscCoarsenData ** ); 1329014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenApply(MatCoarsen); 1330014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenDestroy(MatCoarsen*); 1331b43b03e9SMark F. Adams 1332bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenRegister(const char[],PetscErrorCode (*)(MatCoarsen)); 1333b43b03e9SMark F. Adams 1334014dd563SJed Brown PETSC_EXTERN PetscBool MatCoarsenRegisterAllCalled; 1335b43b03e9SMark F. Adams 1336607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenRegisterAll(void); 1337b43b03e9SMark F. Adams 1338014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenView(MatCoarsen,PetscViewer); 1339014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetFromOptions(MatCoarsen); 134019fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenGetType(MatCoarsen,MatCoarsenType*); 1341ce1779c8SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatCoarsenViewFromOptions(MatCoarsen A,const char prefix[],const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,prefix,name);} 1342b43b03e9SMark F. Adams 1343014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMeshToVertexGraph(Mat,PetscInt,Mat*); 1344014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMeshToCellGraph(Mat,PetscInt,Mat*); 1345591294e4SBarry Smith 13460752156aSBarry Smith /* 13470a835dfdSSatish Balay If you add entries here you must also add them to finclude/petscmat.h 1348d4fbbf0eSBarry Smith */ 13491c1c02c0SLois Curfman McInnes typedef enum { MATOP_SET_VALUES=0, 13501c1c02c0SLois Curfman McInnes MATOP_GET_ROW=1, 13511c1c02c0SLois Curfman McInnes MATOP_RESTORE_ROW=2, 13521c1c02c0SLois Curfman McInnes MATOP_MULT=3, 13531c1c02c0SLois Curfman McInnes MATOP_MULT_ADD=4, 13547c922b88SBarry Smith MATOP_MULT_TRANSPOSE=5, 13557c922b88SBarry Smith MATOP_MULT_TRANSPOSE_ADD=6, 13561c1c02c0SLois Curfman McInnes MATOP_SOLVE=7, 13571c1c02c0SLois Curfman McInnes MATOP_SOLVE_ADD=8, 13587c922b88SBarry Smith MATOP_SOLVE_TRANSPOSE=9, 13597c922b88SBarry Smith MATOP_SOLVE_TRANSPOSE_ADD=10, 13601c1c02c0SLois Curfman McInnes MATOP_LUFACTOR=11, 13611c1c02c0SLois Curfman McInnes MATOP_CHOLESKYFACTOR=12, 1362a714c06dSBarry Smith MATOP_SOR=13, 13631c1c02c0SLois Curfman McInnes MATOP_TRANSPOSE=14, 13641c1c02c0SLois Curfman McInnes MATOP_GETINFO=15, 13651c1c02c0SLois Curfman McInnes MATOP_EQUAL=16, 13661c1c02c0SLois Curfman McInnes MATOP_GET_DIAGONAL=17, 13671c1c02c0SLois Curfman McInnes MATOP_DIAGONAL_SCALE=18, 13681c1c02c0SLois Curfman McInnes MATOP_NORM=19, 13691c1c02c0SLois Curfman McInnes MATOP_ASSEMBLY_BEGIN=20, 13701c1c02c0SLois Curfman McInnes MATOP_ASSEMBLY_END=21, 1371d519adbfSMatthew Knepley MATOP_SET_OPTION=22, 1372d519adbfSMatthew Knepley MATOP_ZERO_ENTRIES=23, 1373d519adbfSMatthew Knepley MATOP_ZERO_ROWS=24, 1374d519adbfSMatthew Knepley MATOP_LUFACTOR_SYMBOLIC=25, 1375d519adbfSMatthew Knepley MATOP_LUFACTOR_NUMERIC=26, 1376d519adbfSMatthew Knepley MATOP_CHOLESKY_FACTOR_SYMBOLIC=27, 1377d519adbfSMatthew Knepley MATOP_CHOLESKY_FACTOR_NUMERIC=28, 1378d519adbfSMatthew Knepley MATOP_SETUP_PREALLOCATION=29, 1379d519adbfSMatthew Knepley MATOP_ILUFACTOR_SYMBOLIC=30, 1380d519adbfSMatthew Knepley MATOP_ICCFACTOR_SYMBOLIC=31, 13819d39f69dSJed Brown /* MATOP_PLACEHOLDER_32=32, */ 13829d39f69dSJed Brown /* MATOP_PLACEHOLDER_33=33, */ 1383d519adbfSMatthew Knepley MATOP_DUPLICATE=34, 1384d519adbfSMatthew Knepley MATOP_FORWARD_SOLVE=35, 1385d519adbfSMatthew Knepley MATOP_BACKWARD_SOLVE=36, 1386d519adbfSMatthew Knepley MATOP_ILUFACTOR=37, 1387d519adbfSMatthew Knepley MATOP_ICCFACTOR=38, 1388d519adbfSMatthew Knepley MATOP_AXPY=39, 1389d519adbfSMatthew Knepley MATOP_GET_SUBMATRICES=40, 1390d519adbfSMatthew Knepley MATOP_INCREASE_OVERLAP=41, 1391d519adbfSMatthew Knepley MATOP_GET_VALUES=42, 1392d519adbfSMatthew Knepley MATOP_COPY=43, 1393d519adbfSMatthew Knepley MATOP_GET_ROW_MAX=44, 1394d519adbfSMatthew Knepley MATOP_SCALE=45, 1395d519adbfSMatthew Knepley MATOP_SHIFT=46, 139635153367SBarry Smith MATOP_DIAGONAL_SET=47, 13979d39f69dSJed Brown MATOP_ZERO_ROWS_COLUMNS=48, 13989d39f69dSJed Brown MATOP_SET_RANDOM=49, 1399d519adbfSMatthew Knepley MATOP_GET_ROW_IJ=50, 1400d519adbfSMatthew Knepley MATOP_RESTORE_ROW_IJ=51, 1401d519adbfSMatthew Knepley MATOP_GET_COLUMN_IJ=52, 1402d519adbfSMatthew Knepley MATOP_RESTORE_COLUMN_IJ=53, 1403d519adbfSMatthew Knepley MATOP_FDCOLORING_CREATE=54, 1404d519adbfSMatthew Knepley MATOP_COLORING_PATCH=55, 1405d519adbfSMatthew Knepley MATOP_SET_UNFACTORED=56, 1406d519adbfSMatthew Knepley MATOP_PERMUTE=57, 1407d519adbfSMatthew Knepley MATOP_SET_VALUES_BLOCKED=58, 1408d519adbfSMatthew Knepley MATOP_GET_SUBMATRIX=59, 1409d519adbfSMatthew Knepley MATOP_DESTROY=60, 1410d519adbfSMatthew Knepley MATOP_VIEW=61, 1411d519adbfSMatthew Knepley MATOP_CONVERT_FROM=62, 14127bab7c10SHong Zhang MATOP_MATMAT_MULT=63, 14137bab7c10SHong Zhang MATOP_MATMAT_MULT_SYMBOLIC=64, 14147bab7c10SHong Zhang MATOP_MATMAT_MULT_NUMERIC=65, 1415d519adbfSMatthew Knepley MATOP_SET_LOCAL_TO_GLOBAL_MAP=66, 1416d519adbfSMatthew Knepley MATOP_SET_VALUES_LOCAL=67, 1417d519adbfSMatthew Knepley MATOP_ZERO_ROWS_LOCAL=68, 1418d519adbfSMatthew Knepley MATOP_GET_ROW_MAX_ABS=69, 1419d519adbfSMatthew Knepley MATOP_GET_ROW_MIN_ABS=70, 1420d519adbfSMatthew Knepley MATOP_CONVERT=71, 1421d519adbfSMatthew Knepley MATOP_SET_COLORING=72, 14229d39f69dSJed Brown /* MATOP_PLACEHOLDER_73=73, */ 1423d519adbfSMatthew Knepley MATOP_SET_VALUES_ADIFOR=74, 1424d519adbfSMatthew Knepley MATOP_FD_COLORING_APPLY=75, 1425d519adbfSMatthew Knepley MATOP_SET_FROM_OPTIONS=76, 1426bda6c4cbSBarry Smith MATOP_MULT_CONSTRAINED=77, 1427bda6c4cbSBarry Smith MATOP_MULT_TRANSPOSE_CONSTRAIN=78, 14289d39f69dSJed Brown MATOP_FIND_ZERO_DIAGONALS=79, 1429d519adbfSMatthew Knepley MATOP_MULT_MULTIPLE=80, 1430d519adbfSMatthew Knepley MATOP_SOLVE_MULTIPLE=81, 1431d519adbfSMatthew Knepley MATOP_GET_INERTIA=82, 1432d519adbfSMatthew Knepley MATOP_LOAD=83, 1433d519adbfSMatthew Knepley MATOP_IS_SYMMETRIC=84, 1434d519adbfSMatthew Knepley MATOP_IS_HERMITIAN=85, 1435d519adbfSMatthew Knepley MATOP_IS_STRUCTURALLY_SYMMETRIC=86, 1436820469bcSHong Zhang MATOP_SET_VALUES_BLOCKEDLOCAL=87, 1437d519adbfSMatthew Knepley MATOP_GET_VECS=88, 1438d519adbfSMatthew Knepley MATOP_MAT_MULT=89, 1439d519adbfSMatthew Knepley MATOP_MAT_MULT_SYMBOLIC=90, 1440d519adbfSMatthew Knepley MATOP_MAT_MULT_NUMERIC=91, 1441d519adbfSMatthew Knepley MATOP_PTAP=92, 1442d519adbfSMatthew Knepley MATOP_PTAP_SYMBOLIC=93, 1443d519adbfSMatthew Knepley MATOP_PTAP_NUMERIC=94, 1444bda6c4cbSBarry Smith MATOP_MAT_TRANSPOSE_MULT=95, 1445bda6c4cbSBarry Smith MATOP_MAT_TRANSPOSE_MULT_SYMBO=96, 1446bda6c4cbSBarry Smith MATOP_MAT_TRANSPOSE_MULT_NUMER=97, 14479d39f69dSJed Brown /* MATOP_PLACEHOLDER_98=98, */ 14489d39f69dSJed Brown /* MATOP_PLACEHOLDER_99=99, */ 14499d39f69dSJed Brown /* MATOP_PLACEHOLDER_100=100, */ 14509d39f69dSJed Brown /* MATOP_PLACEHOLDER_101=101, */ 1451d519adbfSMatthew Knepley MATOP_CONJUGATE=102, 14529d39f69dSJed Brown /* MATOP_PLACEHOLDER_103=103, */ 1453d519adbfSMatthew Knepley MATOP_SET_VALUES_ROW=104, 1454d519adbfSMatthew Knepley MATOP_REAL_PART=105, 1455bda6c4cbSBarry Smith MATOP_IMAGINARY_PART=106, 1456bda6c4cbSBarry Smith MATOP_GET_ROW_UPPER_TRIANGULAR=107, 1457bda6c4cbSBarry Smith MATOP_RESTORE_ROW_UPPER_TRIANG=108, 1458bda6c4cbSBarry Smith MATOP_MAT_SOLVE=109, 14590e8d04f7SBarry Smith MATOP_GET_REDUNDANT_MATRIX=110, 1460d519adbfSMatthew Knepley MATOP_GET_ROW_MIN=111, 14610e8d04f7SBarry Smith MATOP_GET_COLUMN_VECTOR=112, 1462d519adbfSMatthew Knepley MATOP_MISSING_DIAGONAL=113, 14630e8d04f7SBarry Smith MATOP_GET_SEQ_NONZERO_STRUCTUR=114, 146489c6957cSBarry Smith MATOP_CREATE=115, 146589c6957cSBarry Smith MATOP_GET_GHOSTS=116, 14660e8d04f7SBarry Smith MATOP_GET_LOCAL_SUB_MATRIX=117, 14670e8d04f7SBarry Smith MATOP_RESTORE_LOCALSUB_MATRIX=118, 1468eeffb40dSHong Zhang MATOP_MULT_DIAGONAL_BLOCK=119, 14690e8d04f7SBarry Smith MATOP_HERMITIAN_TRANSPOSE=120, 14700e8d04f7SBarry Smith MATOP_MULT_HERMITIAN_TRANSPOSE=121, 14710e8d04f7SBarry Smith MATOP_MULT_HERMITIAN_TRANS_ADD=122, 14720e8d04f7SBarry Smith MATOP_GET_MULTI_PROC_BLOCK=123, 14739d39f69dSJed Brown MATOP_FIND_NONZERO_ROWS=124, 14740e8d04f7SBarry Smith MATOP_GET_COLUMN_NORMS=125, 14759d39f69dSJed Brown MATOP_INVERT_BLOCK_DIAGONAL=126, 14769d39f69dSJed Brown /* MATOP_PLACEHOLDER_127=127, */ 14770e8d04f7SBarry Smith MATOP_GET_SUB_MATRICES_PARALLE=128, 14782b8ad9a3SHong Zhang MATOP_SET_VALUES_BATCH=129, 14790e8d04f7SBarry Smith MATOP_TRANSPOSE_MAT_MULT=130, 14800e8d04f7SBarry Smith MATOP_TRANSPOSE_MAT_MULT_SYMBO=131, 14810e8d04f7SBarry Smith MATOP_TRANSPOSE_MAT_MULT_NUMER=132, 14820e8d04f7SBarry Smith MATOP_TRANSPOSE_COLORING_CREAT=133, 14830e8d04f7SBarry Smith MATOP_TRANS_COLORING_APPLY_SPT=134, 14840e8d04f7SBarry Smith MATOP_TRANS_COLORING_APPLY_DEN=135, 14850e8d04f7SBarry Smith MATOP_RART=136, 14860e8d04f7SBarry Smith MATOP_RART_SYMBOLIC=137, 14870e8d04f7SBarry Smith MATOP_RART_NUMERIC=138, 1488e09a3074SHong Zhang MATOP_SET_BLOCK_SIZES=139, 1489f9426fe0SMark Adams MATOP_AYPX=140, 14901919a2e2SJed Brown MATOP_RESIDUAL=141, 14919c8f2541SHong Zhang MATOP_FDCOLORING_SETUP=142, 14929c8f2541SHong Zhang MATOP_MPICONCATENATESEQ=144 1493fae171e0SBarry Smith } MatOperation; 1494014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatHasOperation(Mat,MatOperation,PetscBool *); 1495014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellSetOperation(Mat,MatOperation,void(*)(void)); 1496014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellGetOperation(Mat,MatOperation,void(**)(void)); 1497014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellSetContext(Mat,void*); 1498112a2221SBarry Smith 149990ace30eSBarry Smith /* 150090ace30eSBarry Smith Codes for matrices stored on disk. By default they are 150190ace30eSBarry Smith stored in a universal format. By changing the format with 15027973ad04SBarry Smith PetscViewerSetFormat(viewer,PETSC_VIEWER_NATIVE); the matrices will 150390ace30eSBarry Smith be stored in a way natural for the matrix, for example dense matrices 150490ace30eSBarry Smith would be stored as dense. Matrices stored this way may only be 1505f4403165SShri Abhyankar read into matrices of the same type. 150690ace30eSBarry Smith */ 150790ace30eSBarry Smith #define MATRIX_BINARY_FORMAT_DENSE -1 150890ace30eSBarry Smith 1509014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetHashTableFactor(Mat,PetscReal); 1510014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatISGetLocalMat(Mat,Mat*); 1511014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatISSetLocalMat(Mat,Mat); 1512b7ce53b6SStefano Zampini PETSC_EXTERN PetscErrorCode MatISGetMPIXAIJ(Mat,MatReuse,Mat*); 15131f4e1ec7SBarry Smith 1514d9274352SBarry Smith /*S 1515d9274352SBarry Smith MatNullSpace - Object that removes a null space from a vector, i.e. 1516d9274352SBarry Smith orthogonalizes the vector to a subsapce 1517d9274352SBarry Smith 1518f7a9e4ceSBarry Smith Level: advanced 1519d9274352SBarry Smith 1520d9274352SBarry Smith Concepts: matrix; linear operator, null space 1521d9274352SBarry Smith 15226e1639daSBarry Smith Users manual sections: 15236e1639daSBarry Smith . sec_singular 15246e1639daSBarry Smith 1525d9274352SBarry Smith .seealso: MatNullSpaceCreate() 1526d9274352SBarry Smith S*/ 152774637425SBarry Smith typedef struct _p_MatNullSpace* MatNullSpace; 1528d9274352SBarry Smith 1529014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceCreate(MPI_Comm,PetscBool ,PetscInt,const Vec[],MatNullSpace*); 1530014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(MatNullSpace,Vec,void*),void*); 1531014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceDestroy(MatNullSpace*); 1532d0195637SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceRemove(MatNullSpace,Vec); 1533014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNullSpace(Mat, MatNullSpace *); 1534014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNullSpace(Mat,MatNullSpace); 1535014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNearNullSpace(Mat,MatNullSpace); 1536014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNearNullSpace(Mat,MatNullSpace*); 1537014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceTest(MatNullSpace,Mat,PetscBool *); 1538014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceView(MatNullSpace,PetscViewer); 1539014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceGetVecs(MatNullSpace,PetscBool*,PetscInt*,const Vec**); 1540014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceCreateRigidBody(Vec,MatNullSpace*); 154174637425SBarry Smith 1542014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatReorderingSeqSBAIJ(Mat,IS); 1543014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetHashTableFactor(Mat,PetscReal); 1544014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetColumnIndices(Mat,PetscInt *); 1545014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJInvertBlockDiagonal(Mat); 15463f1d51d7SBarry Smith 1547014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMAIJ(Mat,PetscInt,Mat*); 1548014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMAIJRedimension(Mat,PetscInt,Mat*); 1549014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMAIJGetAIJ(Mat,Mat*); 1550c4f061fbSSatish Balay 1551014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatComputeExplicitOperator(Mat,Mat*); 1552b0a32e0cSBarry Smith 1553014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalScaleLocal(Mat,Vec); 155404f1ad80SBarry Smith 1555014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*); 1556014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetBase(Mat,Vec,Vec); 1557014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*); 1558014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*)); 1559014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec)); 1560014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDAddNullSpace(Mat,MatNullSpace); 1561014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt); 1562014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDResetHHistory(Mat); 1563014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctionError(Mat,PetscReal); 1564014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetPeriod(Mat,PetscInt); 1565014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDGetH(Mat,PetscScalar *); 1566014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetOptionsPrefix(Mat,const char[]); 1567014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*); 1568014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*); 1569e884886fSBarry Smith 15706370053bSBarry Smith /*S 15716370053bSBarry Smith MatMFFD - A data structured used to manage the computation of the h differencing parameter for matrix-free 15726370053bSBarry Smith Jacobian vector products 1573e884886fSBarry Smith 15746370053bSBarry Smith Notes: MATMFFD is a specific MatType which uses the MatMFFD data structure 15756370053bSBarry Smith 15766370053bSBarry Smith MatMFFD*() methods actually take the Mat as their first argument. Not a MatMFFD data structure 15776370053bSBarry Smith 15786370053bSBarry Smith Level: developer 15796370053bSBarry Smith 15806370053bSBarry Smith .seealso: MATMFFD, MatCreateMFFD(), MatMFFDSetFuction(), MatMFFDSetType(), MatMFFDRegister() 15816370053bSBarry Smith S*/ 1582e884886fSBarry Smith typedef struct _p_MatMFFD* MatMFFD; 1583e884886fSBarry Smith 158476bdecfbSBarry Smith /*J 1585e884886fSBarry Smith MatMFFDType - algorithm used to compute the h used in computing matrix-vector products via differencing of the function 1586e884886fSBarry Smith 1587e884886fSBarry Smith Level: beginner 1588e884886fSBarry Smith 1589e884886fSBarry Smith .seealso: MatMFFDSetType(), MatMFFDRegister() 159076bdecfbSBarry Smith J*/ 159119fd82e9SBarry Smith typedef const char* MatMFFDType; 1592e884886fSBarry Smith #define MATMFFD_DS "ds" 1593e884886fSBarry Smith #define MATMFFD_WP "wp" 1594e884886fSBarry Smith 159519fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDSetType(Mat,MatMFFDType); 1596bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDRegister(const char[],PetscErrorCode (*)(MatMFFD)); 1597e884886fSBarry Smith 1598607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDRegisterAll(void); 1599014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDDSSetUmin(Mat,PetscReal); 1600014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDWPSetComputeNormU(Mat,PetscBool ); 1601e884886fSBarry Smith 1602014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *); 1603014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *); 16047dbadf16SMatthew Knepley 160597969023SHong Zhang /* 160697969023SHong Zhang PETSc interface to MUMPS 160797969023SHong Zhang */ 160897969023SHong Zhang #ifdef PETSC_HAVE_MUMPS 1609014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMumpsSetIcntl(Mat,PetscInt,PetscInt); 1610bc6112feSHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetIcntl(Mat,PetscInt,PetscInt*); 1611b9e9ea26SSatish Balay PETSC_EXTERN PetscErrorCode MatMumpsSetCntl(Mat,PetscInt,PetscReal); 1612bc6112feSHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetCntl(Mat,PetscInt,PetscReal*); 1613bc6112feSHong Zhang 1614ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetInfo(Mat,PetscInt,PetscInt*); 1615ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetInfog(Mat,PetscInt,PetscInt*); 1616ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetRinfo(Mat,PetscInt,PetscReal*); 1617ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetRinfog(Mat,PetscInt,PetscReal*); 161897969023SHong Zhang #endif 161923a5497aSJed Brown 1620d954db57SHong Zhang /* 1621d615f992SSatish Balay PETSc interface to Mkl_Pardiso 1622b500a6b7SJose David Bermeol */ 1623b500a6b7SJose David Bermeol #ifdef PETSC_HAVE_MKL_PARDISO 1624d615f992SSatish Balay PETSC_EXTERN PetscErrorCode MatMkl_PardisoSetCntl(Mat,PetscInt,PetscInt); 1625b500a6b7SJose David Bermeol #endif 1626b500a6b7SJose David Bermeol 1627b500a6b7SJose David Bermeol /* 1628d954db57SHong Zhang PETSc interface to SUPERLU 1629d954db57SHong Zhang */ 1630d954db57SHong Zhang #ifdef PETSC_HAVE_SUPERLU 1631014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSuperluSetILUDropTol(Mat,PetscReal); 1632d954db57SHong Zhang #endif 1633d954db57SHong Zhang 1634aa372e3fSPaul Mullowney #ifdef PETSC_HAVE_CUDA 16353f9c0db1SPaul Mullowney /*E 1636e057df02SPaul Mullowney MatCUSPARSEStorageFormat - indicates the storage format for CUSPARSE (GPU) 16372692e278SPaul Mullowney matrices. 1638e057df02SPaul Mullowney 1639e057df02SPaul Mullowney Not Collective 1640e057df02SPaul Mullowney 16418468deeeSKarl Rupp + MAT_CUSPARSE_CSR - Compressed Sparse Row 16422692e278SPaul Mullowney . MAT_CUSPARSE_ELL - Ellpack (requires CUDA 4.2 or later). 16432692e278SPaul Mullowney - MAT_CUSPARSE_HYB - Hybrid, a combination of Ellpack and Coordinate format (requires CUDA 4.2 or later). 1644e057df02SPaul Mullowney 1645e057df02SPaul Mullowney Level: intermediate 1646e057df02SPaul Mullowney 1647e057df02SPaul Mullowney Any additions/changes here MUST also be made in include/finclude/petscmat.h 1648e057df02SPaul Mullowney 1649e057df02SPaul Mullowney .seealso: MatCUSPARSESetFormat(), MatCUSPARSEFormatOperation 1650e057df02SPaul Mullowney E*/ 1651e057df02SPaul Mullowney 16523f9c0db1SPaul Mullowney typedef enum {MAT_CUSPARSE_CSR, MAT_CUSPARSE_ELL, MAT_CUSPARSE_HYB} MatCUSPARSEStorageFormat; 1653e057df02SPaul Mullowney 1654e057df02SPaul Mullowney /* these will be strings associated with enumerated type defined above */ 1655e057df02SPaul Mullowney PETSC_EXTERN const char *const MatCUSPARSEStorageFormats[]; 1656e057df02SPaul Mullowney 16573f9c0db1SPaul Mullowney /*E 1658e057df02SPaul Mullowney MatCUSPARSEFormatOperation - indicates the operation of CUSPARSE (GPU) 16592692e278SPaul Mullowney matrices whose operation should use a particular storage format. 1660e057df02SPaul Mullowney 1661e057df02SPaul Mullowney Not Collective 1662e057df02SPaul Mullowney 16638468deeeSKarl Rupp + MAT_CUSPARSE_MULT_DIAG - sets the storage format for the diagonal matrix in the parallel MatMult 16648468deeeSKarl Rupp . MAT_CUSPARSE_MULT_OFFDIAG - sets the storage format for the offdiagonal matrix in the parallel MatMult 16658468deeeSKarl Rupp . MAT_CUSPARSE_MULT - sets the storage format for the entire matrix in the serial (single GPU) MatMult 16668468deeeSKarl Rupp - MAT_CUSPARSE_ALL - sets the storage format for all CUSPARSE (GPU) matrices 1667e057df02SPaul Mullowney 1668e057df02SPaul Mullowney Level: intermediate 1669e057df02SPaul Mullowney 1670e057df02SPaul Mullowney .seealso: MatCUSPARSESetFormat(), MatCUSPARSEStorageFormat 1671e057df02SPaul Mullowney E*/ 167236d62e41SPaul Mullowney typedef enum {MAT_CUSPARSE_MULT_DIAG, MAT_CUSPARSE_MULT_OFFDIAG, MAT_CUSPARSE_MULT, MAT_CUSPARSE_ALL} MatCUSPARSEFormatOperation; 1673e057df02SPaul Mullowney 167410e9db80SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 167510e9db80SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCreateAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 1676e057df02SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCUSPARSESetFormat(Mat,MatCUSPARSEFormatOperation,MatCUSPARSEStorageFormat); 16779ae82921SPaul Mullowney #endif 16789ae82921SPaul Mullowney 167990c53902SBarry Smith #if defined(PETSC_HAVE_CUSP) 1680014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCUSP(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 1681014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateAIJCUSP(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 1682e057df02SPaul Mullowney 16833f9c0db1SPaul Mullowney /*E 1684e057df02SPaul Mullowney MatCUSPStorageFormat - indicates the storage format for CUSP (GPU) 168536d62e41SPaul Mullowney matrices. 1686e057df02SPaul Mullowney 1687e057df02SPaul Mullowney Not Collective 1688e057df02SPaul Mullowney 16898468deeeSKarl Rupp + MAT_CUSP_CSR - Compressed Sparse Row 16908468deeeSKarl Rupp . MAT_CUSP_DIA - Diagonal 16918468deeeSKarl Rupp - MAT_CUSP_ELL - Ellpack 1692e057df02SPaul Mullowney 1693e057df02SPaul Mullowney Level: intermediate 1694e057df02SPaul Mullowney 1695e057df02SPaul Mullowney Any additions/changes here MUST also be made in include/finclude/petscmat.h 1696e057df02SPaul Mullowney 1697e057df02SPaul Mullowney .seealso: MatCUSPSetFormat(), MatCUSPFormatOperation 1698e057df02SPaul Mullowney E*/ 16993f9c0db1SPaul Mullowney typedef enum {MAT_CUSP_CSR, MAT_CUSP_DIA, MAT_CUSP_ELL} MatCUSPStorageFormat; 1700e057df02SPaul Mullowney 1701e057df02SPaul Mullowney /* these will be strings associated with enumerated type defined above */ 1702e057df02SPaul Mullowney PETSC_EXTERN const char *const MatCUSPStorageFormats[]; 1703e057df02SPaul Mullowney 17043f9c0db1SPaul Mullowney /*E 1705e057df02SPaul Mullowney MatCUSPFormatOperation - indicates the operation of CUSP (GPU) 170636d62e41SPaul Mullowney matrices whose operation should use a particular storage format. 1707e057df02SPaul Mullowney 1708e057df02SPaul Mullowney Not Collective 1709e057df02SPaul Mullowney 17108468deeeSKarl Rupp + MAT_CUSP_MULT_DIAG - sets the storage format for the diagonal matrix in the parallel MatMult 17118468deeeSKarl Rupp . MAT_CUSP_MULT_OFFDIAG - sets the storage format for the offdiagonal matrix in the parallel MatMult 17128468deeeSKarl Rupp . MAT_CUSP_MULT - sets the storage format for the entire matrix in the serial (single GPU) MatMult 17138468deeeSKarl Rupp - MAT_CUSP_ALL - sets the storage format for all CUSP (GPU) matrices 1714e057df02SPaul Mullowney 1715e057df02SPaul Mullowney Level: intermediate 1716e057df02SPaul Mullowney 1717e057df02SPaul Mullowney Any additions/changes here MUST also be made in include/finclude/petscmat.h 1718e057df02SPaul Mullowney 1719e057df02SPaul Mullowney .seealso: MatCUSPSetFormat(), MatCUSPStorageFormat 1720e057df02SPaul Mullowney E*/ 172136d62e41SPaul Mullowney typedef enum {MAT_CUSP_MULT_DIAG, MAT_CUSP_MULT_OFFDIAG, MAT_CUSP_MULT, MAT_CUSP_ALL} MatCUSPFormatOperation; 1722e057df02SPaul Mullowney 1723e057df02SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCUSPSetFormat(Mat,MatCUSPFormatOperation,MatCUSPStorageFormat); 1724e057df02SPaul Mullowney #endif 172590c53902SBarry Smith 1726d67ff14aSKarl Rupp #if defined(PETSC_HAVE_VIENNACL) 1727d67ff14aSKarl Rupp PETSC_EXTERN PetscErrorCode MatCreateSeqAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 1728d67ff14aSKarl Rupp PETSC_EXTERN PetscErrorCode MatCreateAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 1729d67ff14aSKarl Rupp #endif 1730d67ff14aSKarl Rupp 173154efbe56SHong Zhang /* 173254efbe56SHong Zhang PETSc interface to FFTW 173354efbe56SHong Zhang */ 173454efbe56SHong Zhang #if defined(PETSC_HAVE_FFTW) 1735014dd563SJed Brown PETSC_EXTERN PetscErrorCode VecScatterPetscToFFTW(Mat,Vec,Vec); 1736014dd563SJed Brown PETSC_EXTERN PetscErrorCode VecScatterFFTWToPetsc(Mat,Vec,Vec); 17372a7a6963SBarry Smith PETSC_EXTERN PetscErrorCode MatCreateVecsFFTW(Mat,Vec*,Vec*,Vec*); 173854efbe56SHong Zhang #endif 173954efbe56SHong Zhang 1740db31f6deSJed Brown #if defined(PETSC_HAVE_ELEMENTAL) 1741607a6623SBarry Smith PETSC_EXTERN PetscErrorCode PetscElementalInitializePackage(void); 1742db31f6deSJed Brown PETSC_EXTERN PetscErrorCode PetscElementalFinalizePackage(void); 1743db31f6deSJed Brown #endif 1744db31f6deSJed Brown 1745014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateNest(MPI_Comm,PetscInt,const IS[],PetscInt,const IS[],const Mat[],Mat*); 1746014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSize(Mat,PetscInt*,PetscInt*); 1747014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetISs(Mat,IS[],IS[]); 1748014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetLocalISs(Mat,IS[],IS[]); 1749014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSubMats(Mat,PetscInt*,PetscInt*,Mat***); 1750014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSubMat(Mat,PetscInt,PetscInt,Mat*); 175119fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatNestSetVecType(Mat,VecType); 1752014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestSetSubMats(Mat,PetscInt,const IS[],PetscInt,const IS[],const Mat[]); 1753014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestSetSubMat(Mat,PetscInt,PetscInt,Mat); 1754d8588912SDave May 17554325cce7SMatthew G Knepley PETSC_EXTERN PetscErrorCode MatChop(Mat,PetscReal); 1756e0a58c10SMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatComputeBandwidth(Mat,PetscReal,PetscInt*); 17574325cce7SMatthew G Knepley 175823a5497aSJed Brown #endif 1759