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" 35273d9f13SBarry Smith #define MATMPIAIJ "mpiaij" 365a11e1b2SBarry Smith #define MATAIJCRL "aijcrl" 375a11e1b2SBarry Smith #define MATSEQAIJCRL "seqaijcrl" 385a11e1b2SBarry Smith #define MATMPIAIJCRL "mpiaijcrl" 398154be41SBarry Smith #define MATAIJCUSP "aijcusp" 408154be41SBarry Smith #define MATSEQAIJCUSP "seqaijcusp" 418154be41SBarry Smith #define MATMPIAIJCUSP "mpiaijcusp" 429ae82921SPaul Mullowney #define MATAIJCUSPARSE "aijcusparse" 439ae82921SPaul Mullowney #define MATSEQAIJCUSPARSE "seqaijcusparse" 449ae82921SPaul Mullowney #define MATMPIAIJCUSPARSE "mpiaijcusparse" 45d67ff14aSKarl Rupp #define MATAIJVIENNACL "aijviennacl" 46d67ff14aSKarl Rupp #define MATSEQAIJVIENNACL "seqaijviennacl" 47d67ff14aSKarl Rupp #define MATMPIAIJVIENNACL "mpiaijviennacl" 485a11e1b2SBarry Smith #define MATAIJPERM "aijperm" 495a11e1b2SBarry Smith #define MATSEQAIJPERM "seqaijperm" 505a11e1b2SBarry Smith #define MATMPIAIJPERM "mpiaijperm" 51273d9f13SBarry Smith #define MATSHELL "shell" 525a11e1b2SBarry Smith #define MATDENSE "dense" 53209238afSKris Buschelman #define MATSEQDENSE "seqdense" 54273d9f13SBarry Smith #define MATMPIDENSE "mpidense" 55db31f6deSJed Brown #define MATELEMENTAL "elemental" 565a11e1b2SBarry Smith #define MATBAIJ "baij" 57273d9f13SBarry Smith #define MATSEQBAIJ "seqbaij" 58273d9f13SBarry Smith #define MATMPIBAIJ "mpibaij" 59273d9f13SBarry Smith #define MATMPIADJ "mpiadj" 605a11e1b2SBarry Smith #define MATSBAIJ "sbaij" 61273d9f13SBarry Smith #define MATSEQSBAIJ "seqsbaij" 62273d9f13SBarry Smith #define MATMPISBAIJ "mpisbaij" 63cebc7f6cSBarry Smith #define MATDAAD "daad" 64cebc7f6cSBarry Smith #define MATMFFD "mffd" 65c8a8475eSBarry Smith #define MATNORMAL "normal" 66c5e4d11fSDmitry Karpeev #define MATNORMALHERMITIAN "normalh" 67ab92ecdeSBarry Smith #define MATLRC "lrc" 682a6744ebSBarry Smith #define MATSCATTER "scatter" 69421e10b8SBarry Smith #define MATBLOCKMAT "blockmat" 70793850ffSBarry Smith #define MATCOMPOSITE "composite" 711f8c7532SHong Zhang #define MATFFT "fft" 721f8c7532SHong Zhang #define MATFFTW "fftw" 73e133240eSMatthew G Knepley #define MATSEQCUFFT "seqcufft" 74557cca28SSatish Balay #define MATTRANSPOSEMAT "transpose" 7572ca8751SBarry Smith #define MATSCHURCOMPLEMENT "schurcomplement" 761d6018f0SLisandro Dalcin #define MATPYTHON "python" 7763c07aadSStefano Zampini #define MATHYPRE "hypre" 78f91d8e95SBarry Smith #define MATHYPRESTRUCT "hyprestruct" 79a9e6138eSGlenn Hammond #define MATHYPRESSTRUCT "hypresstruct" 80ee1cef2cSJed Brown #define MATSUBMATRIX "submatrix" 812c0dbf93SJed Brown #define MATLOCALREF "localref" 82d8588912SDave May #define MATNEST "nest" 83c094ef40SMatthew G. Knepley #define MATPREALLOCATOR "preallocator" 84a3b2e22bSHong Zhang #define MATDUMMY "dummy" 85773941d6SBarry Smith 8676bdecfbSBarry Smith /*J 87c7393fdbSBarry Smith MatSolverPackage - String with the name of a PETSc matrix solver type. 889989ab13SBarry Smith 89c2b89b5dSBarry Smith For example: "petsc" indicates what PETSc provides, "superlu_dist" the parallel SuperLU_DIST package etc 909989ab13SBarry Smith 919989ab13SBarry Smith Level: beginner 929989ab13SBarry Smith 935c9eb25fSBarry Smith .seealso: MatGetFactor(), Mat, MatSetType(), MatType 9476bdecfbSBarry Smith J*/ 95c7393fdbSBarry Smith #define MatSolverPackage char* 962692d6eeSBarry Smith #define MATSOLVERSUPERLU "superlu" 972692d6eeSBarry Smith #define MATSOLVERSUPERLU_DIST "superlu_dist" 9808f5efcfSPieter Ghysels #define MATSOLVERSTRUMPACK "strumpack" 992692d6eeSBarry Smith #define MATSOLVERUMFPACK "umfpack" 10020db9a53SJed Brown #define MATSOLVERCHOLMOD "cholmod" 101d89f5e7aSBarry Smith #define MATSOLVERKLU "klu" 102418810c4SBarry Smith #define MATSOLVERSPARSEELEMENTAL "sparseelemental" 103d89f5e7aSBarry Smith #define MATSOLVERELEMENTAL "elemental" 1042692d6eeSBarry Smith #define MATSOLVERESSL "essl" 1052692d6eeSBarry Smith #define MATSOLVERLUSOL "lusol" 1062692d6eeSBarry Smith #define MATSOLVERMUMPS "mumps" 107d615f992SSatish Balay #define MATSOLVERMKL_PARDISO "mkl_pardiso" 108d305a81bSVasiliy Kozyrev #define MATSOLVERMKL_CPARDISO "mkl_cpardiso" 1092692d6eeSBarry Smith #define MATSOLVERPASTIX "pastix" 1102692d6eeSBarry Smith #define MATSOLVERMATLAB "matlab" 1112692d6eeSBarry Smith #define MATSOLVERPETSC "petsc" 1122692d6eeSBarry Smith #define MATSOLVERBAS "bas" 1139ae82921SPaul Mullowney #define MATSOLVERCUSPARSE "cusparse" 114c0cdd4a1SDahai Guo 115b24902e0SBarry Smith /*E 116b24902e0SBarry Smith MatFactorType - indicates what type of factorization is requested 117b24902e0SBarry Smith 118b24902e0SBarry Smith Level: beginner 119b24902e0SBarry Smith 120af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 121b24902e0SBarry Smith 122c7393fdbSBarry Smith .seealso: MatSolverPackage, MatGetFactor() 123b24902e0SBarry Smith E*/ 124599ef60dSHong Zhang typedef enum {MAT_FACTOR_NONE, MAT_FACTOR_LU, MAT_FACTOR_CHOLESKY, MAT_FACTOR_ILU, MAT_FACTOR_ICC,MAT_FACTOR_ILUDT} MatFactorType; 125014dd563SJed Brown PETSC_EXTERN const char *const MatFactorTypes[]; 126e92e720dSBarry Smith 127014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactor(Mat,const MatSolverPackage,MatFactorType,Mat*); 128014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactorAvailable(Mat,const MatSolverPackage,MatFactorType,PetscBool *); 129014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFactorGetSolverPackage(Mat,const MatSolverPackage*); 130014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactorType(Mat,MatFactorType*); 13118713533SBarry Smith PETSC_EXTERN PetscErrorCode MatSolverPackageRegister(const MatSolverPackage,const MatType,MatFactorType,PetscErrorCode(*)(Mat,MatFactorType,Mat*)); 13242c9c57cSBarry Smith PETSC_EXTERN PetscErrorCode MatSolverPackageGet(const MatSolverPackage,const MatType,MatFactorType,PetscBool*,PetscBool*,PetscErrorCode (**)(Mat,MatFactorType,Mat*)); 1339989ab13SBarry Smith 134c06d978dSMatthew Knepley /* Logging support */ 1350700a824SBarry Smith #define MAT_FILE_CLASSID 1211216 /* used to indicate matrices in binary files */ 136014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_CLASSID; 137335efc43SPeter Brune PETSC_EXTERN PetscClassId MAT_COLORING_CLASSID; 138014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_FDCOLORING_CLASSID; 139014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_TRANSPOSECOLORING_CLASSID; 140014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_PARTITIONING_CLASSID; 141014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_COARSEN_CLASSID; 142014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_NULLSPACE_CLASSID; 143014dd563SJed Brown PETSC_EXTERN PetscClassId MATMFFD_CLASSID; 144c06d978dSMatthew Knepley 145ceb03754SKris Buschelman /*E 1467dae84e0SHong Zhang MatReuse - Indicates if matrices obtained from a previous call to MatCreateSubMatrices(), MatCreateSubMatrix(), MatConvert() or several other functions 147cf37664fSBarry Smith are to be reused to store the new matrix values. 148cf37664fSBarry Smith 149cf37664fSBarry Smith $ MAT_INITIAL_MATRIX - create a new matrix 150cf37664fSBarry Smith $ MAT_REUSE_MATRIX - reuse the matrix created with a previous call that used MAT_INITIAL_MATRIX 151cf37664fSBarry Smith $ MAT_INPLACE_MATRIX - replace the first input matrix with the new matrix (not applicable to all functions) 152cf37664fSBarry Smith $ MAT_IGNORE_MATRIX - do not create a new matrix or reuse a give matrix, just ignore that matrix argument (not applicable to all functions) 153ceb03754SKris Buschelman 154ceb03754SKris Buschelman Level: beginner 155ceb03754SKris Buschelman 156af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 157ceb03754SKris Buschelman 1587dae84e0SHong Zhang .seealso: MatCreateSubMatrices(), MatCreateSubMatrix(), MatDestroyMatrices(), MatConvert() 159ceb03754SKris Buschelman E*/ 160511c6705SHong Zhang typedef enum {MAT_INITIAL_MATRIX,MAT_REUSE_MATRIX,MAT_IGNORE_MATRIX,MAT_INPLACE_MATRIX} MatReuse; 161ceb03754SKris Buschelman 1625494a064SHong Zhang /*E 1637dae84e0SHong Zhang MatCreateSubMatrixOption - Indicates if matrices obtained from a call to MatCreateSubMatrices() 164829201f2SHong Zhang include the matrix values. Currently it is only used by MatGetSeqNonzerostructure(). 1655494a064SHong Zhang 1665494a064SHong Zhang Level: beginner 1675494a064SHong Zhang 168829201f2SHong Zhang .seealso: MatGetSeqNonzerostructure() 1695494a064SHong Zhang E*/ 1707dae84e0SHong Zhang typedef enum {MAT_DO_NOT_GET_VALUES,MAT_GET_VALUES} MatCreateSubMatrixOption; 1715494a064SHong Zhang 172607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatInitializePackage(void); 173c06d978dSMatthew Knepley 174014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreate(MPI_Comm,Mat*); 175014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetSizes(Mat,PetscInt,PetscInt,PetscInt,PetscInt); 17619fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatSetType(Mat,MatType); 177014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetFromOptions(Mat); 178685405a1SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatViewFromOptions(Mat A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);} 179bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatRegister(const char[],PetscErrorCode(*)(Mat)); 180014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRegisterBaseName(const char[],const char[],const char[]); 181014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetOptionsPrefix(Mat,const char[]); 182014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAppendOptionsPrefix(Mat,const char[]); 183014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOptionsPrefix(Mat,const char*[]); 18484d44b13SHong Zhang PETSC_EXTERN PetscErrorCode MatSetErrorIfFailure(Mat,PetscBool); 185f69a0ea3SMatthew Knepley 186140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatList; 187140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatColoringList; 188140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatPartitioningList; 18928988994SBarry Smith 1903b224e63SBarry Smith /*E 191aa6c7ce3SBarry Smith MatStructure - Indicates if two matrices have the same nonzero structure 1923b224e63SBarry Smith 1933b224e63SBarry Smith Level: beginner 1943b224e63SBarry Smith 195af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 1963b224e63SBarry Smith 197aa6c7ce3SBarry Smith .seealso: MatCopy(), MatAXPY() 1983b224e63SBarry Smith E*/ 199aa6c7ce3SBarry Smith typedef enum {DIFFERENT_NONZERO_PATTERN,SUBSET_NONZERO_PATTERN,SAME_NONZERO_PATTERN} MatStructure; 2003b224e63SBarry Smith 201014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqDense(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*); 202014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateDense(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*); 203014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 204014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 205014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *); 206014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithSplitArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],PetscInt[],PetscInt[],PetscScalar[],Mat*); 2073b00a383SHong Zhang PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithSeqAIJ(MPI_Comm,Mat,Mat,const PetscInt[],Mat*); 2083b00a383SHong Zhang 209014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 210014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 211014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat*); 212d21a29f3SJed Brown 213014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAdj(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscInt[],Mat*); 214014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 215d21a29f3SJed Brown 216014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 217014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPISBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *); 21838f409ebSLisandro Dalcin PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 219014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 2204cce697cSJed Brown PETSC_EXTERN PetscErrorCode MatXAIJSetPreallocation(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscInt[],const PetscInt[]); 221dfb205c3SBarry Smith 222014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateShell(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,void *,Mat*); 223014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateNormal(Mat,Mat*); 224c5e4d11fSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatCreateNormalHermitian(Mat,Mat*); 225986c4d72SJose E. Roman PETSC_EXTERN PetscErrorCode MatCreateLRC(Mat,Mat,Vec,Mat,Mat*); 226267ca84cSJose E. Roman PETSC_EXTERN PetscErrorCode MatLRCGetMats(Mat,Mat*,Mat*,Vec*,Mat*); 227c5e4d11fSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatCreateIS(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,ISLocalToGlobalMapping,ISLocalToGlobalMapping,Mat*); 228014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 229014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 230c0cdd4a1SDahai Guo 231014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateScatter(MPI_Comm,VecScatter,Mat*); 232014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScatterSetVecScatter(Mat,VecScatter); 233014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScatterGetVecScatter(Mat,VecScatter*); 234014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateBlockMat(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,Mat*); 235014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeAddMat(Mat,Mat); 236014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeMerge(Mat); 237014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateComposite(MPI_Comm,PetscInt,const Mat*,Mat*); 2386d7c1e57SBarry Smith typedef enum {MAT_COMPOSITE_ADDITIVE,MAT_COMPOSITE_MULTIPLICATIVE} MatCompositeType; 239014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeSetType(Mat,MatCompositeType); 2406d7c1e57SBarry Smith 24119fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatCreateFFT(MPI_Comm,PetscInt,const PetscInt[],MatType,Mat*); 242014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqCUFFT(MPI_Comm,PetscInt,const PetscInt[],Mat*); 243dedccee8SHong Zhang 244014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateTranspose(Mat,Mat*); 2458060fb66Sstefano_zampini PETSC_EXTERN PetscErrorCode MatTransposeGetMat(Mat,Mat*); 246d0de2241SAndrew Spott PETSC_EXTERN PetscErrorCode MatCreateHermitianTranspose(Mat,Mat*); 24754e05e6cSHong Zhang PETSC_EXTERN PetscErrorCode MatCreateSubMatrixVirtual(Mat,IS,IS,Mat*); 24854e05e6cSHong Zhang PETSC_EXTERN PetscErrorCode MatSubMatrixVirtualUpdate(Mat,Mat,IS,IS); 249014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateLocalRef(Mat,IS,IS,Mat*); 2501472f72bSBarry Smith 251978814f1SStefano Zampini #if defined(PETSC_HAVE_HYPRE) 252d975228cSstefano_zampini PETSC_EXTERN PetscErrorCode MatHYPRESetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 253978814f1SStefano Zampini #endif 254978814f1SStefano Zampini 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 27862ef0227SBarry Smith column of a matrix as indexed on an associated grid. These are arguments to MatSetStencil() and MatSetBlockStencil() 27962ef0227SBarry Smith 28062ef0227SBarry Smith The i,j, and k represent the logical coordinates over the entire grid (for 2 and 1 dimensional problems the k and j entries are ignored). 28162ef0227SBarry Smith The c represents the the degrees of freedom at each grid point (the dof argument to DMDASetDOF()). If dof is 1 then this entry is ignored. 28262ef0227SBarry Smith 28362ef0227SBarry Smith For stencil access to vectors see DMDAVecGetArray(), DMDAVecGetArrayF90(). 284be479b30SJed Brown 285be479b30SJed Brown Fortran usage is different, see MatSetValuesStencil() for details. 2862ef4de8bSBarry Smith 2872ef4de8bSBarry Smith Level: beginner 2882ef4de8bSBarry Smith 2892ef4de8bSBarry Smith Concepts: matrix; linear operator 2902ef4de8bSBarry Smith 29162ef0227SBarry Smith .seealso: MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockedStencil(), DMDAVecGetArray(), DMDAVecGetArrayF90() 2922ef4de8bSBarry Smith S*/ 293435da068SBarry Smith typedef struct { 294c1ac3661SBarry Smith PetscInt k,j,i,c; 295435da068SBarry Smith } MatStencil; 2962ef4de8bSBarry Smith 297014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode); 298014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode); 299014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt); 300435da068SBarry Smith 301d91e6319SBarry Smith /*E 302d91e6319SBarry Smith MatAssemblyType - Indicates if the matrix is now to be used, or if you plan 30362ef0227SBarry Smith to continue to add or insert values to it 304d91e6319SBarry Smith 305d91e6319SBarry Smith Level: beginner 306d91e6319SBarry Smith 307d91e6319SBarry Smith .seealso: MatAssemblyBegin(), MatAssemblyEnd() 308d91e6319SBarry Smith E*/ 3096d4a8577SBarry Smith typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType; 310014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssemblyBegin(Mat,MatAssemblyType); 311014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssemblyEnd(Mat,MatAssemblyType); 312014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssembled(Mat,PetscBool *); 3134f9c727eSBarry Smith 3141d73ed98SBarry Smith 31530de9b25SBarry Smith 316d91e6319SBarry Smith /*E 317d91e6319SBarry Smith MatOption - Options that may be set for a matrix and its behavior or storage 318d91e6319SBarry Smith 319d91e6319SBarry Smith Level: beginner 320d91e6319SBarry Smith 321af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 322d91e6319SBarry Smith 323382164b0SBarry Smith Developer Notes: Entries that are negative need not be called collectively by all processes. 324382164b0SBarry Smith 325d91e6319SBarry Smith .seealso: MatSetOption() 326d91e6319SBarry Smith E*/ 327c5e4d11fSDmitry Karpeev typedef enum {MAT_OPTION_MIN = -3, 328c5e4d11fSDmitry Karpeev MAT_UNUSED_NONZERO_LOCATION_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, 340c5e4d11fSDmitry Karpeev MAT_NEW_NONZERO_LOCATION_ERR = 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, 348c5e4d11fSDmitry Karpeev MAT_NEW_NONZERO_ALLOCATION_ERR = 19, 349c5e4d11fSDmitry Karpeev MAT_SUBSET_OFF_PROC_ENTRIES = 20, 350c10200c1SHong Zhang MAT_SUBMAT_SINGLEIS = 21, 351957cac9fSHong Zhang MAT_STRUCTURE_ONLY = 22, 352957cac9fSHong Zhang MAT_OPTION_MAX = 23} MatOption; 353e057df02SPaul Mullowney 354014dd563SJed Brown PETSC_EXTERN const char *MatOptions[]; 355014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetOption(Mat,MatOption,PetscBool); 3567d68702bSBarry Smith PETSC_EXTERN PetscErrorCode MatGetOption(Mat,MatOption,PetscBool*); 35719fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatGetType(Mat,MatType*); 35884cb2905SBarry Smith 359014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]); 360014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 361014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 362014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowUpperTriangular(Mat); 363014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRowUpperTriangular(Mat); 364014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnVector(Mat,Vec,PetscInt); 3658c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJGetArray(Mat,PetscScalar *[]); 3668c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJRestoreArray(Mat,PetscScalar *[]); 36721e72a00SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJGetMaxRowNonzeros(Mat,PetscInt*); 368bd04181cSBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJSetValuesLocalFast(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 3698397e458SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqSBAIJGetArray(Mat,PetscScalar *[]); 3708397e458SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqSBAIJRestoreArray(Mat,PetscScalar *[]); 3718c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatDenseGetArray(Mat,PetscScalar *[]); 3728c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatDenseRestoreArray(Mat,PetscScalar *[]); 373*d3042a70SBarry Smith PETSC_EXTERN PetscErrorCode MatDensePlaceArray(Mat,const PetscScalar[]); 374*d3042a70SBarry Smith PETSC_EXTERN PetscErrorCode MatDenseResetArray(Mat); 375014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBlockSize(Mat,PetscInt *); 376014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSize(Mat,PetscInt); 377014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBlockSizes(Mat,PetscInt *,PetscInt *); 378014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSizes(Mat,PetscInt,PetscInt); 37933d57670SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSizesFromMats(Mat,Mat,Mat); 3801620fd73SBarry Smith 381014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMult(Mat,Vec,Vec); 382014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultDiagonalBlock(Mat,Vec,Vec); 383014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultAdd(Mat,Vec,Vec,Vec); 384014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTranspose(Mat,Vec,Vec); 385014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultHermitianTranspose(Mat,Vec,Vec); 386014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsTranspose(Mat,Mat,PetscReal,PetscBool *); 387014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitianTranspose(Mat,Mat,PetscReal,PetscBool *); 388014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeAdd(Mat,Vec,Vec,Vec); 389014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultHermitianTransposeAdd(Mat,Vec,Vec,Vec); 390014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultConstrained(Mat,Vec,Vec); 391014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeConstrained(Mat,Vec,Vec); 392014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatSolve(Mat,Mat,Mat); 393f9426fe0SMark Adams PETSC_EXTERN PetscErrorCode MatResidual(Mat,Vec,Vec,Vec); 394f5edf698SHong Zhang 395d91e6319SBarry Smith /*E 396d91e6319SBarry Smith MatDuplicateOption - Indicates if a duplicated sparse matrix should have 397d91e6319SBarry Smith its numerical values copied over or just its nonzero structure. 398d91e6319SBarry Smith 399d91e6319SBarry Smith Level: beginner 400d91e6319SBarry Smith 401af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 402d91e6319SBarry Smith 40370dcbbb9SBarry Smith $ MAT_SHARE_NONZERO_PATTERN - the i and j arrays in the new matrix will be shared with the original matrix 40470dcbbb9SBarry Smith $ this also triggers the MAT_DO_NOT_COPY_VALUES option. This is used when you 40570dcbbb9SBarry Smith $ have several matrices with the same nonzero pattern. 40670dcbbb9SBarry Smith 407d91e6319SBarry Smith .seealso: MatDuplicate() 408d91e6319SBarry Smith E*/ 40970dcbbb9SBarry Smith typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES,MAT_SHARE_NONZERO_PATTERN} MatDuplicateOption; 4102e8a6d31SBarry Smith 41119fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatConvert(Mat,MatType,MatReuse,Mat*); 412014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDuplicate(Mat,MatDuplicateOption,Mat*); 41394a9d846SBarry Smith 414cb5b572fSBarry Smith 415014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCopy(Mat,Mat,MatStructure); 416014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatView(Mat,PetscViewer); 417014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsSymmetric(Mat,PetscReal,PetscBool *); 418014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsStructurallySymmetric(Mat,PetscBool *); 419014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitian(Mat,PetscReal,PetscBool *); 420014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsSymmetricKnown(Mat,PetscBool *,PetscBool *); 421014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitianKnown(Mat,PetscBool *,PetscBool *); 422014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMissingDiagonal(Mat,PetscBool *,PetscInt *); 423014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLoad(Mat, PetscViewer); 4247b80b807SBarry Smith 4251a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool *); 4261a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool *); 4271a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool *); 4281a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool *); 429d4fbbf0eSBarry Smith 430d91e6319SBarry Smith /*S 431d91e6319SBarry Smith MatInfo - Context of matrix information, used with MatGetInfo() 432d91e6319SBarry Smith 433d91e6319SBarry Smith In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE 434d91e6319SBarry Smith 435d91e6319SBarry Smith Level: intermediate 436d91e6319SBarry Smith 437d91e6319SBarry Smith Concepts: matrix^nonzero information 438d91e6319SBarry Smith 439d9274352SBarry Smith .seealso: MatGetInfo(), MatInfoType 440d91e6319SBarry Smith S*/ 4414e220ebcSLois Curfman McInnes typedef struct { 442b0a32e0cSBarry Smith PetscLogDouble block_size; /* block size */ 443b0a32e0cSBarry Smith PetscLogDouble nz_allocated,nz_used,nz_unneeded; /* number of nonzeros */ 444b0a32e0cSBarry Smith PetscLogDouble memory; /* memory allocated */ 445b0a32e0cSBarry Smith PetscLogDouble assemblies; /* number of matrix assemblies called */ 446b0a32e0cSBarry Smith PetscLogDouble mallocs; /* number of mallocs during MatSetValues() */ 447b0a32e0cSBarry Smith PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */ 448b0a32e0cSBarry Smith PetscLogDouble factor_mallocs; /* number of mallocs during factorization */ 4494e220ebcSLois Curfman McInnes } MatInfo; 4504e220ebcSLois Curfman McInnes 451d9274352SBarry Smith /*E 452d9274352SBarry Smith MatInfoType - Indicates if you want information about the local part of the matrix, 453d9274352SBarry Smith the entire parallel matrix or the maximum over all the local parts. 454d9274352SBarry Smith 455d9274352SBarry Smith Level: beginner 456d9274352SBarry Smith 457af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 458d9274352SBarry Smith 459d9274352SBarry Smith .seealso: MatGetInfo(), MatInfo 460d9274352SBarry Smith E*/ 4617b80b807SBarry Smith typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType; 462014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetInfo(Mat,MatInfoType,MatInfo*); 463014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetDiagonal(Mat,Vec); 464014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMax(Mat,Vec,PetscInt[]); 465014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMin(Mat,Vec,PetscInt[]); 466014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMaxAbs(Mat,Vec,PetscInt[]); 467014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMinAbs(Mat,Vec,PetscInt[]); 468014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowSum(Mat,Vec); 469014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTranspose(Mat,MatReuse,Mat*); 470014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatHermitianTranspose(Mat,MatReuse,Mat*); 471014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPermute(Mat,IS,IS,Mat*); 472014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalScale(Mat,Vec,Vec); 473014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalSet(Mat,Vec,InsertMode); 474a52676ecSHong Zhang 475014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatEqual(Mat,Mat,PetscBool*); 476014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultEqual(Mat,Mat,PetscInt,PetscBool*); 477014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultAddEqual(Mat,Mat,PetscInt,PetscBool*); 478014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeEqual(Mat,Mat,PetscInt,PetscBool*); 479014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscBool*); 480a52676ecSHong Zhang PETSC_EXTERN PetscErrorCode MatMatMultEqual(Mat,Mat,Mat,PetscInt,PetscBool*); 481a52676ecSHong Zhang PETSC_EXTERN PetscErrorCode MatTransposeMatMultEqual(Mat,Mat,Mat,PetscInt,PetscBool*); 4827b80b807SBarry Smith 483014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNorm(Mat,NormType,PetscReal*); 484014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnNorms(Mat,NormType,PetscReal*); 485014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroEntries(Mat); 486014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 487014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsIS(Mat,IS,PetscScalar,Vec,Vec); 488014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsStencil(Mat,PetscInt,const MatStencil [],PetscScalar,Vec,Vec); 489014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsStencil(Mat,PetscInt,const MatStencil[],PetscScalar,Vec,Vec); 490014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumns(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 491014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsIS(Mat,IS,PetscScalar,Vec,Vec); 4927b80b807SBarry Smith 493014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSize(Mat,PetscInt*,PetscInt*); 494014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalSize(Mat,PetscInt*,PetscInt*); 495014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRange(Mat,PetscInt*,PetscInt*); 496014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRanges(Mat,const PetscInt**); 497014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRangeColumn(Mat,PetscInt*,PetscInt*); 498014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRangesColumn(Mat,const PetscInt**); 499566876eaSJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipIS(Mat,IS*,IS*); 5005ef9f2a5SBarry Smith 5017dae84e0SHong Zhang PETSC_EXTERN PetscErrorCode MatCreateSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]); 502cd2534ddSHong Zhang PETSC_DEPRECATED("Use MatCreateSubMatrices()") PETSC_STATIC_INLINE PetscErrorCode MatGetSubMatrices(Mat mat,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *submat[]) {return MatCreateSubMatrices(mat,n,irow,icol,scall,submat);} 5037dae84e0SHong Zhang PETSC_EXTERN PetscErrorCode MatCreateSubMatricesMPI(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]); 504cd2534ddSHong Zhang PETSC_DEPRECATED("Use MatCreateSubMatricesMPI()") PETSC_STATIC_INLINE PetscErrorCode MatGetSubMatricesMPI(Mat mat,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *submat[]) {return MatCreateSubMatricesMPI(mat,n,irow,icol,scall,submat);} 505014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroyMatrices(PetscInt,Mat *[]); 506df750dc8SHong Zhang PETSC_EXTERN PetscErrorCode MatDestroySubMatrices(PetscInt,Mat *[]); 5077dae84e0SHong Zhang PETSC_EXTERN PetscErrorCode MatCreateSubMatrix(Mat,IS,IS,MatReuse,Mat *); 508cd2534ddSHong Zhang PETSC_DEPRECATED("Use MatCreateSubMatrix()") PETSC_STATIC_INLINE PetscErrorCode MatGetSubMatrix(Mat mat,IS isrow,IS iscol,MatReuse cll,Mat *newmat) {return MatCreateSubMatrix(mat,isrow,iscol,cll,newmat);} 509014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalSubMatrix(Mat,IS,IS,Mat*); 510014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreLocalSubMatrix(Mat,IS,IS,Mat*); 511014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSeqNonzeroStructure(Mat,Mat*); 512014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroySeqNonzeroStructure(Mat*); 5137b80b807SBarry Smith 514014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJ(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*); 515014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJSymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*); 516014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJNumeric(Mat,Mat); 517014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMat(Mat,MatReuse,Mat*); 518014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*); 519014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,Mat*); 520014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetGhosts(Mat, PetscInt *,const PetscInt *[]); 5218efafbd8SBarry Smith 522014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt); 523aa1e27eaSFande Kong PETSC_EXTERN PetscErrorCode MatIncreaseOverlapSplit(Mat mat,PetscInt n,IS is[],PetscInt ov); 524d2b2a242SBarry Smith PETSC_EXTERN PetscErrorCode MatMPIAIJSetUseScalableIncreaseOverlap(Mat,PetscBool); 5257b80b807SBarry Smith 526014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*); 527014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMultSymbolic(Mat,Mat,PetscReal,Mat*); 528014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMultNumeric(Mat,Mat,Mat); 52922440eb1SKris Buschelman 5307bab7c10SHong Zhang PETSC_EXTERN PetscErrorCode MatMatMatMult(Mat,Mat,Mat,MatReuse,PetscReal,Mat*); 5312df6c5c3SPatrick Farrell PETSC_EXTERN PetscErrorCode MatGalerkin(Mat,Mat,Mat,MatReuse,PetscReal,Mat*); 5327bab7c10SHong Zhang 533014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*); 534014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAPSymbolic(Mat,Mat,PetscReal,Mat*); 535014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAPNumeric(Mat,Mat,Mat); 536014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARt(Mat,Mat,MatReuse,PetscReal,Mat*); 537014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARtSymbolic(Mat,Mat,PetscReal,Mat*); 538014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARtNumeric(Mat,Mat,Mat); 53922440eb1SKris Buschelman 540014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeMatMult(Mat,Mat,MatReuse,PetscReal,Mat*); 541014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposetMatMultSymbolic(Mat,Mat,PetscReal,Mat*); 542014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposetMatMultNumeric(Mat,Mat,Mat); 543014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMult(Mat,Mat,MatReuse,PetscReal,Mat*); 544014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMultSymbolic(Mat,Mat,PetscReal,Mat*); 545014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMultNumeric(Mat,Mat,Mat); 546bc011b1eSHong Zhang 547014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAXPY(Mat,PetscScalar,Mat,MatStructure); 548014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAYPX(Mat,PetscScalar,Mat,MatStructure); 5491c741599SHong Zhang 550014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScale(Mat,PetscScalar); 551014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShift(Mat,PetscScalar); 5527b80b807SBarry Smith 553014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping); 554014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*); 555a93ff8c4SPeter Brune PETSC_EXTERN PetscErrorCode MatGetLayouts(Mat,PetscLayout*,PetscLayout*); 556014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 557014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsLocalIS(Mat,IS,PetscScalar,Vec,Vec); 558014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 559014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocalIS(Mat,IS,PetscScalar,Vec,Vec); 560014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 561014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 562052efed2SBarry Smith 563014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStashSetInitialSize(Mat,PetscInt,PetscInt); 564014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*); 56590f02eecSBarry Smith 566014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInterpolate(Mat,Vec,Vec); 567014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInterpolateAdd(Mat,Vec,Vec,Vec); 568014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestrict(Mat,Vec,Vec); 5692a7a6963SBarry Smith PETSC_EXTERN PetscErrorCode MatCreateVecs(Mat,Vec*,Vec*); 57084cc2a5aSBarry Smith PETSC_DEPRECATED("Use MatCreateVecs()") PETSC_STATIC_INLINE PetscErrorCode MatGetVecs(Mat mat,Vec *x,Vec *y) {return MatCreateVecs(mat,x,y);} 57153cd1579SHong Zhang PETSC_EXTERN PetscErrorCode MatCreateRedundantMatrix(Mat,PetscInt,MPI_Comm,MatReuse,Mat*); 572014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetMultiProcBlock(Mat,MPI_Comm,MatReuse,Mat*); 573014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFindZeroDiagonals(Mat,IS*); 5743a062f41SBarry Smith PETSC_EXTERN PetscErrorCode MatFindOffBlockDiagonalEntries(Mat,IS*); 5759c8f2541SHong Zhang PETSC_EXTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat(MPI_Comm,Mat,PetscInt,MatReuse,Mat*); 576bd481603SBarry Smith 577bd481603SBarry Smith /*MC 5781620fd73SBarry Smith MatSetValue - Set a single entry into a matrix. 5791620fd73SBarry Smith 5801620fd73SBarry Smith Not collective 5811620fd73SBarry Smith 582a9834a7dSSatish Balay Synopsis: 583a9834a7dSSatish Balay #include <petscmat.h> 584a9834a7dSSatish Balay PetscErrorCode MatSetValue(Mat m,PetscInt row,PetscInt col,PetscScalar value,InsertMode mode) 585a9834a7dSSatish Balay 5861620fd73SBarry Smith Input Parameters: 5871620fd73SBarry Smith + m - the matrix 5881620fd73SBarry Smith . row - the row location of the entry 5891620fd73SBarry Smith . col - the column location of the entry 5901620fd73SBarry Smith . value - the value to insert 5911620fd73SBarry Smith - mode - either INSERT_VALUES or ADD_VALUES 5921620fd73SBarry Smith 5931620fd73SBarry Smith Notes: 5941620fd73SBarry Smith For efficiency one should use MatSetValues() and set several or many 5951620fd73SBarry Smith values simultaneously if possible. 5961620fd73SBarry Smith 5971620fd73SBarry Smith Level: beginner 5981620fd73SBarry Smith 5991620fd73SBarry Smith .seealso: MatSetValues(), MatSetValueLocal() 6001620fd73SBarry Smith M*/ 6011620fd73SBarry 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);} 6021620fd73SBarry Smith 6031620fd73SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatGetValue(Mat v,PetscInt i,PetscInt j,PetscScalar *va) {return MatGetValues(v,1,&i,1,&j,va);} 6041620fd73SBarry Smith 6051620fd73SBarry 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);} 6061620fd73SBarry Smith 6071620fd73SBarry Smith /*MC 6080d2a7ff2SSatish Balay MatPreallocateInitialize - Begins the block of code that will count the number of nonzeros per 609bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 610bd481603SBarry Smith 611bd481603SBarry Smith Synopsis: 612aaa7dc30SBarry Smith #include <petscmat.h> 613c1ac3661SBarry Smith PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz) 614bd481603SBarry Smith 615bd481603SBarry Smith Collective on MPI_Comm 616bd481603SBarry Smith 617bd481603SBarry Smith Input Parameters: 618bd481603SBarry Smith + comm - the communicator that will share the eventually allocated matrix 619859fcb39SBarry Smith . nrows - the number of LOCAL rows in the matrix 620859fcb39SBarry Smith - ncols - the number of LOCAL columns in the matrix 621bd481603SBarry Smith 622bd481603SBarry Smith Output Parameters: 623bd481603SBarry Smith + dnz - the array that will be passed to the matrix preallocation routines 624bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 625bd481603SBarry Smith 626bd481603SBarry Smith Level: intermediate 627bd481603SBarry Smith 628bd481603SBarry Smith Notes: 629a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 630bd481603SBarry Smith 6311d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 632bd481603SBarry Smith 633bd481603SBarry Smith Use MatPreallocateInitializeSymmetric() for symmetric matrices (MPISBAIJ matrices) 634bd481603SBarry Smith 6351620fd73SBarry Smith This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize(). 6361620fd73SBarry Smith 637bd481603SBarry Smith Concepts: preallocation^Matrix 638bd481603SBarry Smith 639d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(), 640d6e23781SBarry Smith MatPreallocateInitializeSymmetric(), MatPreallocateSymmetricSetLocalBlock() 641bd481603SBarry Smith M*/ 642c4f061fbSSatish Balay #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \ 6437c922b88SBarry Smith { \ 644a89bc211SBarry Smith PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__start,__end; \ 6451795a4d1SJed Brown _4_ierr = PetscCalloc2(__nrows,&dnz,__nrows,&onz);CHKERRQ(_4_ierr); \ 6461795a4d1SJed Brown __start = 0; __end = __start; \ 647c1ac3661SBarry Smith _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ctmp;\ 648a89bc211SBarry Smith _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows; 6497c922b88SBarry Smith 650bd481603SBarry Smith /*MC 651bd481603SBarry Smith MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be 652bd481603SBarry Smith inserted using a local number of the rows and columns 653bd481603SBarry Smith 654bd481603SBarry Smith Synopsis: 655aaa7dc30SBarry Smith #include <petscmat.h> 656c1ac3661SBarry Smith PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 657bd481603SBarry Smith 658bd481603SBarry Smith Not Collective 659bd481603SBarry Smith 660bd481603SBarry Smith Input Parameters: 661784ac674SJed Brown + map - the row mapping from local numbering to global numbering 662bd481603SBarry Smith . nrows - the number of rows indicated 6631d73ed98SBarry Smith . rows - the indices of the rows 664784ac674SJed Brown . cmap - the column mapping from local to global numbering 665bd481603SBarry Smith . ncols - the number of columns in the matrix 666bd481603SBarry Smith . cols - the columns indicated 667bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 668bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 669bd481603SBarry Smith 670bd481603SBarry Smith Level: intermediate 671bd481603SBarry Smith 672bd481603SBarry Smith Notes: 673a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 674bd481603SBarry Smith 6751d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 676bd481603SBarry Smith 677bd481603SBarry Smith Concepts: preallocation^Matrix 678bd481603SBarry Smith 679d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(), 680c1154cd5SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocalRemoveDups() 681bd481603SBarry Smith M*/ 682784ac674SJed Brown #define MatPreallocateSetLocal(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \ 683c4f061fbSSatish Balay {\ 684c1ac3661SBarry Smith PetscInt __l;\ 685784ac674SJed Brown _4_ierr = ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\ 686784ac674SJed Brown _4_ierr = ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\ 687c4f061fbSSatish Balay for (__l=0;__l<nrows;__l++) {\ 688ef66eb69SBarry Smith _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 689c4f061fbSSatish Balay }\ 690c4f061fbSSatish Balay } 691c4f061fbSSatish Balay 692bd481603SBarry Smith /*MC 693c1154cd5SBarry Smith MatPreallocateSetLocalRemoveDups - Indicates the locations (rows and columns) in the matrix where nonzeros will be 694c1154cd5SBarry Smith inserted using a local number of the rows and columns. This version removes any duplicate columns in cols 695c1154cd5SBarry Smith 696c1154cd5SBarry Smith Synopsis: 697c1154cd5SBarry Smith #include <petscmat.h> 698c1154cd5SBarry Smith PetscErrorCode MatPreallocateSetLocalRemoveDups(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 699c1154cd5SBarry Smith 700c1154cd5SBarry Smith Not Collective 701c1154cd5SBarry Smith 702c1154cd5SBarry Smith Input Parameters: 703c1154cd5SBarry Smith + map - the row mapping from local numbering to global numbering 704c1154cd5SBarry Smith . nrows - the number of rows indicated 705c1154cd5SBarry Smith . rows - the indices of the rows (these values are mapped to the global values) 706c1154cd5SBarry Smith . cmap - the column mapping from local to global numbering 707c1154cd5SBarry Smith . ncols - the number of columns in the matrix (this value will be changed if duplicate columns are found) 708c1154cd5SBarry Smith . cols - the columns indicated (these values are mapped to the global values, they are then sorted and duplicates removed) 709c1154cd5SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 710c1154cd5SBarry Smith - ozn - the other array passed to the matrix preallocation routines 711c1154cd5SBarry Smith 712c1154cd5SBarry Smith Level: intermediate 713c1154cd5SBarry Smith 714c1154cd5SBarry Smith Notes: 715c1154cd5SBarry Smith See Users-Manual: ch_performance for more details. 716c1154cd5SBarry Smith 717c1154cd5SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 718c1154cd5SBarry Smith 719c1154cd5SBarry Smith Concepts: preallocation^Matrix 720c1154cd5SBarry Smith 721c1154cd5SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(), 722c1154cd5SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocal() 723c1154cd5SBarry Smith M*/ 724c1154cd5SBarry Smith #define MatPreallocateSetLocalRemoveDups(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \ 725c1154cd5SBarry Smith {\ 726c1154cd5SBarry Smith PetscInt __l;\ 727c1154cd5SBarry Smith _4_ierr = ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\ 728c1154cd5SBarry Smith _4_ierr = ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\ 729c1154cd5SBarry Smith _4_ierr = PetscSortRemoveDupsInt(&ncols,cols);CHKERRQ(_4_ierr);\ 730c1154cd5SBarry Smith for (__l=0;__l<nrows;__l++) {\ 731c1154cd5SBarry Smith _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 732c1154cd5SBarry Smith }\ 733c1154cd5SBarry Smith } 734c1154cd5SBarry Smith 735c1154cd5SBarry Smith /*MC 736d6e23781SBarry Smith MatPreallocateSetLocalBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be 737bd481603SBarry Smith inserted using a local number of the rows and columns 738bd481603SBarry Smith 739bd481603SBarry Smith Synopsis: 740aaa7dc30SBarry Smith #include <petscmat.h> 741d6e23781SBarry Smith PetscErrorCode MatPreallocateSetLocalBlock(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 742d6e23781SBarry Smith 743d6e23781SBarry Smith Not Collective 744d6e23781SBarry Smith 745d6e23781SBarry Smith Input Parameters: 746d6e23781SBarry Smith + map - the row mapping from local numbering to global numbering 747d6e23781SBarry Smith . nrows - the number of rows indicated 748d6e23781SBarry Smith . rows - the indices of the rows 749d6e23781SBarry Smith . cmap - the column mapping from local to global numbering 750d6e23781SBarry Smith . ncols - the number of columns in the matrix 751d6e23781SBarry Smith . cols - the columns indicated 752d6e23781SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 753d6e23781SBarry Smith - ozn - the other array passed to the matrix preallocation routines 754d6e23781SBarry Smith 755d6e23781SBarry Smith Level: intermediate 756d6e23781SBarry Smith 757d6e23781SBarry Smith Notes: 758d6e23781SBarry Smith See Users-Manual: ch_performance for more details. 759d6e23781SBarry Smith 760d6e23781SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 761d6e23781SBarry Smith 762d6e23781SBarry Smith Concepts: preallocation^Matrix 763d6e23781SBarry Smith 764d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(), 765d6e23781SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock() 766d6e23781SBarry Smith M*/ 767d6e23781SBarry Smith #define MatPreallocateSetLocalBlock(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \ 768d6e23781SBarry Smith {\ 769d6e23781SBarry Smith PetscInt __l;\ 770d6e23781SBarry Smith _4_ierr = ISLocalToGlobalMappingApplyBlock(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\ 771d6e23781SBarry Smith _4_ierr = ISLocalToGlobalMappingApplyBlock(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\ 772d6e23781SBarry Smith for (__l=0;__l<nrows;__l++) {\ 773d6e23781SBarry Smith _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 774d6e23781SBarry Smith }\ 775d6e23781SBarry Smith } 776d6e23781SBarry Smith 777d6e23781SBarry Smith /*MC 778d6e23781SBarry Smith MatPreallocateSymmetricSetLocalBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be 779d6e23781SBarry Smith inserted using a local number of the rows and columns 780d6e23781SBarry Smith 781d6e23781SBarry Smith Synopsis: 782d6e23781SBarry Smith #include <petscmat.h> 783d6e23781SBarry Smith PetscErrorCode MatPreallocateSymmetricSetLocalBlock(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 784bd481603SBarry Smith 785bd481603SBarry Smith Not Collective 786bd481603SBarry Smith 787bd481603SBarry Smith Input Parameters: 788bd481603SBarry Smith + map - the mapping between local numbering and global numbering 789bd481603SBarry Smith . nrows - the number of rows indicated 7901d73ed98SBarry Smith . rows - the indices of the rows 791bd481603SBarry Smith . ncols - the number of columns in the matrix 792bd481603SBarry Smith . cols - the columns indicated 793bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 794bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 795bd481603SBarry Smith 796bd481603SBarry Smith Level: intermediate 797bd481603SBarry Smith 798bd481603SBarry Smith Notes: 799a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 800bd481603SBarry Smith 801bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 802bd481603SBarry Smith 803bd481603SBarry Smith Concepts: preallocation^Matrix 804bd481603SBarry Smith 805d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateInitialize(), 806d6e23781SBarry Smith MatPreallocateInitialize(), MatPreallocateSetLocal() 807bd481603SBarry Smith M*/ 808d6e23781SBarry Smith #define MatPreallocateSymmetricSetLocalBlock(map,nrows,rows,ncols,cols,dnz,onz) 0;\ 809d3d32019SBarry Smith {\ 810c1ac3661SBarry Smith PetscInt __l;\ 811d6e23781SBarry Smith _4_ierr = ISLocalToGlobalMappingApplyBlock(map,nrows,rows,rows);CHKERRQ(_4_ierr);\ 812d6e23781SBarry Smith _4_ierr = ISLocalToGlobalMappingApplyBlock(map,ncols,cols,cols);CHKERRQ(_4_ierr);\ 813d3d32019SBarry Smith for (__l=0;__l<nrows;__l++) {\ 814d6e23781SBarry Smith _4_ierr = MatPreallocateSymmetricSetBlock((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 815d3d32019SBarry Smith }\ 816d3d32019SBarry Smith } 817bd481603SBarry Smith /*MC 818bd481603SBarry Smith MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be 819bd481603SBarry Smith inserted using a local number of the rows and columns 820bd481603SBarry Smith 821bd481603SBarry Smith Synopsis: 822aaa7dc30SBarry Smith #include <petscmat.h> 823c1ac3661SBarry Smith PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 824bd481603SBarry Smith 825bd481603SBarry Smith Not Collective 826bd481603SBarry Smith 827bd481603SBarry Smith Input Parameters: 82864075487SBarry Smith + row - the row 829bd481603SBarry Smith . ncols - the number of columns in the matrix 830a50f8bf6SHong Zhang - cols - the columns indicated 831a50f8bf6SHong Zhang 832a50f8bf6SHong Zhang Output Parameters: 833a50f8bf6SHong Zhang + dnz - the array that will be passed to the matrix preallocation routines 834bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 835bd481603SBarry Smith 836bd481603SBarry Smith Level: intermediate 837bd481603SBarry Smith 838bd481603SBarry Smith Notes: 839a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 840bd481603SBarry Smith 841bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 842bd481603SBarry Smith 8431620fd73SBarry Smith This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize(). 8441620fd73SBarry Smith 845bd481603SBarry Smith Concepts: preallocation^Matrix 846bd481603SBarry Smith 847d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(), 848d6e23781SBarry Smith MatPreallocateInitialize(), MatPreallocateSetLocal() 849bd481603SBarry Smith M*/ 850c4f061fbSSatish Balay #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\ 851c1ac3661SBarry Smith { PetscInt __i; \ 852e32f2f54SBarry 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);\ 853e32f2f54SBarry 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);\ 8547c922b88SBarry Smith for (__i=0; __i<nc; __i++) {\ 85564075487SBarry Smith if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \ 8567cc688d7SBarry Smith else dnz[row - __rstart]++;\ 8577c922b88SBarry Smith }\ 8587c922b88SBarry Smith } 8597c922b88SBarry Smith 860bd481603SBarry Smith /*MC 861d6e23781SBarry Smith MatPreallocateSymmetricSetBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be 862bd481603SBarry Smith inserted using a local number of the rows and columns 863bd481603SBarry Smith 864bd481603SBarry Smith Synopsis: 865aaa7dc30SBarry Smith #include <petscmat.h> 866d6e23781SBarry Smith PetscErrorCode MatPreallocateSymmetricSetBlock(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 867bd481603SBarry Smith 868bd481603SBarry Smith Not Collective 869bd481603SBarry Smith 870bd481603SBarry Smith Input Parameters: 871bd481603SBarry Smith + nrows - the number of rows indicated 8721d73ed98SBarry Smith . rows - the indices of the rows 873bd481603SBarry Smith . ncols - the number of columns in the matrix 874bd481603SBarry Smith . cols - the columns indicated 875bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 876bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 877bd481603SBarry Smith 878bd481603SBarry Smith Level: intermediate 879bd481603SBarry Smith 880bd481603SBarry Smith Notes: 881a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 882bd481603SBarry Smith 883bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 884bd481603SBarry Smith 8851620fd73SBarry Smith This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize(). 8861620fd73SBarry Smith 887bd481603SBarry Smith Concepts: preallocation^Matrix 888bd481603SBarry Smith 889d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateInitialize(), 890d6e23781SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocal() 891bd481603SBarry Smith M*/ 892d6e23781SBarry Smith #define MatPreallocateSymmetricSetBlock(row,nc,cols,dnz,onz) 0;\ 893c1ac3661SBarry Smith { PetscInt __i; \ 894d3d32019SBarry Smith for (__i=0; __i<nc; __i++) {\ 895d3d32019SBarry Smith if (cols[__i] >= __end) onz[row - __rstart]++; \ 896d3d32019SBarry Smith else if (cols[__i] >= row) dnz[row - __rstart]++;\ 897d3d32019SBarry Smith }\ 898d3d32019SBarry Smith } 899d3d32019SBarry Smith 900bd481603SBarry Smith /*MC 90116371a99SBarry Smith MatPreallocateLocation - An alternative to MatPreallocationSet() that puts the nonzero locations into the matrix if it exists 90216371a99SBarry Smith 90316371a99SBarry Smith Synopsis: 904aaa7dc30SBarry Smith #include <petscmat.h> 90516371a99SBarry Smith PetscErrorCode MatPreallocateLocations(Mat A,PetscInt row,PetscInt ncols,PetscInt *cols,PetscInt *dnz,PetscInt *onz) 90616371a99SBarry Smith 90716371a99SBarry Smith Not Collective 90816371a99SBarry Smith 90916371a99SBarry Smith Input Parameters: 91016371a99SBarry Smith . A - matrix 91116371a99SBarry Smith . row - row where values exist (must be local to this process) 91216371a99SBarry Smith . ncols - number of columns 91316371a99SBarry Smith . cols - columns with nonzeros 91416371a99SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 91516371a99SBarry Smith - ozn - the other array passed to the matrix preallocation routines 91616371a99SBarry Smith 91716371a99SBarry Smith Level: intermediate 91816371a99SBarry Smith 91916371a99SBarry Smith Notes: 920a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 92116371a99SBarry Smith 92216371a99SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 92316371a99SBarry Smith 92416371a99SBarry Smith This is a MACRO not a function because it uses a bunch of variables private to the MatPreallocation.... routines. 92516371a99SBarry Smith 92616371a99SBarry Smith Concepts: preallocation^Matrix 92716371a99SBarry Smith 928d6e23781SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(), 929d6e23781SBarry Smith MatPreallocateSymmetricSetLocalBlock() 93016371a99SBarry Smith M*/ 9310298fd71SBarry 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);} 93216371a99SBarry Smith 93316371a99SBarry Smith 93416371a99SBarry Smith /*MC 9350d2a7ff2SSatish Balay MatPreallocateFinalize - Ends the block of code that will count the number of nonzeros per 936bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 937bd481603SBarry Smith 938bd481603SBarry Smith Synopsis: 939aaa7dc30SBarry Smith #include <petscmat.h> 940c1ac3661SBarry Smith PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz) 941bd481603SBarry Smith 942bd481603SBarry Smith Collective on MPI_Comm 943bd481603SBarry Smith 944bd481603SBarry Smith Input Parameters: 94516371a99SBarry Smith + dnz - the array that was be passed to the matrix preallocation routines 946bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 947bd481603SBarry Smith 948bd481603SBarry Smith Level: intermediate 949bd481603SBarry Smith 950bd481603SBarry Smith Notes: 951a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 952bd481603SBarry Smith 953bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 954bd481603SBarry Smith 9551620fd73SBarry Smith This is a MACRO not a function because it closes the { started in MatPreallocateInitialize(). 9561620fd73SBarry Smith 957bd481603SBarry Smith Concepts: preallocation^Matrix 958bd481603SBarry Smith 959d6e23781SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(), 960d6e23781SBarry Smith MatPreallocateSymmetricSetLocalBlock() 961bd481603SBarry Smith M*/ 962a89bc211SBarry Smith #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree2(dnz,onz);CHKERRQ(_4_ierr);} 9637c922b88SBarry Smith 9647b80b807SBarry Smith /* Routines unique to particular data structures */ 965014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellGetContext(Mat,void *); 966698d4c6aSKris Buschelman 967014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInodeAdjustForInodes(Mat,IS*,IS*); 968014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *); 969698d4c6aSKris Buschelman 970014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetColumnIndices(Mat,PetscInt[]); 971014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetColumnIndices(Mat,PetscInt[]); 972014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 973014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 974014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 975014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJFromTriple(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*,PetscInt,PetscBool); 9767b80b807SBarry Smith 977a96a251dSBarry Smith #define MAT_SKIP_ALLOCATION -4 978a96a251dSBarry Smith 979014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]); 980014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]); 981014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]); 982273d9f13SBarry Smith 983014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 984014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 985014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 986014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []); 987014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 988014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]); 989014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 990014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]); 991014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIDenseSetPreallocation(Mat,PetscScalar[]); 992014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqDenseSetPreallocation(Mat,PetscScalar[]); 9939230625dSJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,const PetscInt*[]); 9949230625dSJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,const PetscInt*[]); 995014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAdjCreateNonemptySubcommMat(Mat,Mat*); 996273d9f13SBarry Smith 9972e1947a5SStefano Zampini PETSC_EXTERN PetscErrorCode MatISSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 998cf0a3239SStefano Zampini PETSC_EXTERN PetscErrorCode MatISSetUpSF(Mat); 999014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqDenseSetLDA(Mat,PetscInt); 1000014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDenseGetLocalMatrix(Mat,Mat*); 10011b807ce4Svictorle 1002014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStoreValues(Mat); 1003014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRetrieveValues(Mat); 10042e8a6d31SBarry Smith 1005014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDAADSetCtx(Mat,void*); 10063a7fca6bSBarry Smith 1007014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFindNonzeroRows(Mat,IS*); 10087cfe41e4SPatrick Farrell PETSC_EXTERN PetscErrorCode MatFindZeroRows(Mat,IS*); 10097b80b807SBarry Smith /* 10107b80b807SBarry Smith These routines are not usually accessed directly, rather solving is 101194b7f48cSBarry Smith done through the KSP and PC interfaces. 10127b80b807SBarry Smith */ 10137b80b807SBarry Smith 101476bdecfbSBarry Smith /*J 10158f6c3df8SBarry Smith MatOrderingType - String with the name of a PETSc matrix ordering 1016d9274352SBarry Smith 1017d9274352SBarry Smith Level: beginner 1018d9274352SBarry Smith 1019d9274352SBarry Smith .seealso: MatGetOrdering() 102076bdecfbSBarry Smith J*/ 102119fd82e9SBarry Smith typedef const char* MatOrderingType; 10222692d6eeSBarry Smith #define MATORDERINGNATURAL "natural" 10232692d6eeSBarry Smith #define MATORDERINGND "nd" 10242692d6eeSBarry Smith #define MATORDERING1WD "1wd" 10252692d6eeSBarry Smith #define MATORDERINGRCM "rcm" 10262692d6eeSBarry Smith #define MATORDERINGQMD "qmd" 10272692d6eeSBarry Smith #define MATORDERINGROWLENGTH "rowlength" 1028510184a7SJed Brown #define MATORDERINGWBM "wbm" 1029fc1bef75SJed Brown #define MATORDERINGSPECTRAL "spectral" 1030312e372aSBarry Smith #define MATORDERINGAMD "amd" /* only works if UMFPACK is installed with PETSc */ 1031b12f92e5SBarry Smith 103219fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatGetOrdering(Mat,MatOrderingType,IS*,IS*); 1033140e18c1SBarry Smith PETSC_EXTERN PetscErrorCode MatGetOrderingList(PetscFunctionList*); 1034bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatOrderingRegister(const char[],PetscErrorCode(*)(Mat,MatOrderingType,IS*,IS*)); 1035140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatOrderingList; 1036d4fbbf0eSBarry Smith 1037014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS); 1038fc1bef75SJed Brown PETSC_EXTERN PetscErrorCode MatCreateLaplacian(Mat,PetscReal,PetscBool,Mat*); 1039a2ce50c7SBarry Smith 1040d91e6319SBarry Smith /*S 1041d90ac83dSHong Zhang MatFactorShiftType - Numeric Shift. 1042d90ac83dSHong Zhang 1043d90ac83dSHong Zhang Level: beginner 1044d90ac83dSHong Zhang 1045d90ac83dSHong Zhang S*/ 1046d90ac83dSHong Zhang typedef enum {MAT_SHIFT_NONE,MAT_SHIFT_NONZERO,MAT_SHIFT_POSITIVE_DEFINITE,MAT_SHIFT_INBLOCKS} MatFactorShiftType; 10476a6fc655SJed Brown PETSC_EXTERN const char *const MatFactorShiftTypes[]; 10485e9742b9SJed Brown PETSC_EXTERN const char *const MatFactorShiftTypesDetail[]; 1049d90ac83dSHong Zhang 1050d90ac83dSHong Zhang /*S 10519aa7eafdSHong Zhang MatFactorError - indicates what type of error in matrix factor 10529aa7eafdSHong Zhang 10539aa7eafdSHong Zhang Level: beginner 10549aa7eafdSHong Zhang 105500e125f8SBarry Smith Developer Notes: Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 105600e125f8SBarry Smith 105700e125f8SBarry Smith .seealso: MatGetFactor() 10589aa7eafdSHong Zhang S*/ 10595cd7cf9dSHong Zhang typedef enum {MAT_FACTOR_NOERROR,MAT_FACTOR_STRUCT_ZEROPIVOT,MAT_FACTOR_NUMERIC_ZEROPIVOT,MAT_FACTOR_OUTMEMORY,MAT_FACTOR_OTHER} MatFactorError; 10609aa7eafdSHong Zhang 106100e125f8SBarry Smith PETSC_EXTERN PetscErrorCode MatFactorGetError(Mat,MatFactorError*); 1062b8b68cfdSBarry Smith PETSC_EXTERN PetscErrorCode MatFactorClearError(Mat); 10637b6c816cSBarry Smith PETSC_EXTERN PetscErrorCode MatFactorGetErrorZeroPivot(Mat,PetscReal*,PetscInt*); 106400e125f8SBarry Smith 10659aa7eafdSHong Zhang /*S 1066422a814eSBarry Smith MatFactorInfo - Data passed into the matrix factorization routines, and information about the resulting factorization 1067b00f7748SHong Zhang 106861cad860SBarry Smith In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE, that is use 106961cad860SBarry Smith $ MatFactorInfo info(MAT_FACTORINFO_SIZE) 1070b00f7748SHong Zhang 107115e8a5b3SHong Zhang Notes: These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC. 1072b00f7748SHong Zhang 107361cad860SBarry Smith You can use MatFactorInfoInitialize() to set default values. 107461cad860SBarry Smith 1075b00f7748SHong Zhang Level: developer 1076b00f7748SHong Zhang 1077d7d82daaSBarry Smith .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(), 1078d7d82daaSBarry Smith MatFactorInfoInitialize() 1079b00f7748SHong Zhang 1080b00f7748SHong Zhang S*/ 1081b00f7748SHong Zhang typedef struct { 108215e8a5b3SHong Zhang PetscReal diagonal_fill; /* force diagonal to fill in if initially not filled */ 108385317021SBarry Smith PetscReal usedt; 108415e8a5b3SHong Zhang PetscReal dt; /* drop tolerance */ 1085b00f7748SHong Zhang PetscReal dtcol; /* tolerance for pivoting */ 108615e8a5b3SHong Zhang PetscReal dtcount; /* maximum nonzeros to be allowed per row */ 108767e28bfeSBarry Smith PetscReal fill; /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */ 1088348344bbSBarry Smith PetscReal levels; /* ICC/ILU(levels) */ 1089bcd9e38bSBarry Smith PetscReal pivotinblocks; /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0 1090bcd9e38bSBarry Smith factorization may be faster if do not pivot */ 109115e8a5b3SHong Zhang PetscReal zeropivot; /* pivot is called zero if less than this */ 1092f4db908eSBarry Smith PetscReal shifttype; /* type of shift added to matrix factor to prevent zero pivots */ 1093f4db908eSBarry Smith PetscReal shiftamount; /* how large the shift is */ 109415e8a5b3SHong Zhang } MatFactorInfo; 1095ffa6d0a5SLois Curfman McInnes 1096014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFactorInfoInitialize(MatFactorInfo*); 10979a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatCholeskyFactor(Mat,IS,const MatFactorInfo*); 10989a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatCholeskyFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*); 10999a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatCholeskyFactorNumeric(Mat,Mat,const MatFactorInfo*); 11009a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatLUFactor(Mat,IS,IS,const MatFactorInfo*); 11019a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatILUFactor(Mat,IS,IS,const MatFactorInfo*); 11029a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatLUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*); 11039a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatILUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*); 11049a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatICCFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*); 11059a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatICCFactor(Mat,IS,const MatFactorInfo*); 11069a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatLUFactorNumeric(Mat,Mat,const MatFactorInfo*); 1107014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*); 1108014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolve(Mat,Vec,Vec); 1109014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatForwardSolve(Mat,Vec,Vec); 1110014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatBackwardSolve(Mat,Vec,Vec); 1111014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveAdd(Mat,Vec,Vec,Vec); 1112014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveTranspose(Mat,Vec,Vec); 1113014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveTransposeAdd(Mat,Vec,Vec,Vec); 1114014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolves(Mat,Vecs,Vecs); 11158e7ba810SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorSetSchurIS(Mat,IS); 11165a05ddb0SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorGetSchurComplement(Mat,Mat*); 11175a05ddb0SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorRestoreSchurComplement(Mat,Mat*); 11185a05ddb0SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorInvertSchurComplement(Mat); 11195a05ddb0SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorCreateSchurComplement(Mat,Mat*); 11205a05ddb0SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorSolveSchurComplement(Mat,Vec,Vec); 11215a05ddb0SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorSolveSchurComplementTranspose(Mat,Vec,Vec); 11226dba178dSStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorFactorizeSchurComplement(Mat); 1123e8ade678SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorSetSchurComplementSolverType(Mat,PetscInt); 1124014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetUnfactored(Mat); 1125bb5a7306SBarry Smith 1126d91e6319SBarry Smith /*E 1127d91e6319SBarry Smith MatSORType - What type of (S)SOR to perform 1128bb1eb677SSatish Balay 1129d91e6319SBarry Smith Level: beginner 1130d91e6319SBarry Smith 1131d9274352SBarry Smith May be bitwise ORd together 1132d9274352SBarry Smith 1133af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 1134d91e6319SBarry Smith 11354e7234bfSBarry Smith MatSORType may be bitwise ORd together, so do not change the numbers 11364e7234bfSBarry Smith 113741f059aeSBarry Smith .seealso: MatSOR() 1138d91e6319SBarry Smith E*/ 1139ee50ffe9SBarry Smith typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3, 1140ee50ffe9SBarry Smith SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8, 1141ee50ffe9SBarry Smith SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16, 114284cb2905SBarry Smith SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType; 1143014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSOR(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec); 11448ed539a5SBarry Smith 1145d4fbbf0eSBarry Smith /* 1146639f9d9dSBarry Smith These routines are for efficiently computing Jacobians via finite differences. 1147639f9d9dSBarry Smith */ 1148b12f92e5SBarry Smith 11497086a01eSPeter Brune /*S 11507086a01eSPeter Brune MatColoring - Object for managing the coloring of matrices. 11517086a01eSPeter Brune 11527086a01eSPeter Brune Level: beginner 11537086a01eSPeter Brune 11547086a01eSPeter Brune Concepts: matrix, coloring 11557086a01eSPeter Brune 11567086a01eSPeter Brune .seealso: MatFDColoringCreate() ISColoring MatFDColoring 11577086a01eSPeter Brune S*/ 11587086a01eSPeter Brune typedef struct _p_MatColoring* MatColoring; 115976bdecfbSBarry Smith /*J 11608f6c3df8SBarry Smith MatColoringType - String with the name of a PETSc matrix coloring 1161d9274352SBarry Smith 1162d9274352SBarry Smith Level: beginner 1163d9274352SBarry Smith 11647086a01eSPeter Brune .seealso: MatColoringSetType(), MatColoring 116576bdecfbSBarry Smith J*/ 11667086a01eSPeter Brune 116719fd82e9SBarry Smith typedef const char* MatColoringType; 11685567d71aSPeter Brune #define MATCOLORINGJP "jp" 1169b9acb617SPeter Brune #define MATCOLORINGPOWER "power" 11702692d6eeSBarry Smith #define MATCOLORINGNATURAL "natural" 11712692d6eeSBarry Smith #define MATCOLORINGSL "sl" 11722692d6eeSBarry Smith #define MATCOLORINGLF "lf" 11732692d6eeSBarry Smith #define MATCOLORINGID "id" 11748121bdceSPeter Brune #define MATCOLORINGGREEDY "greedy" 1175b12f92e5SBarry Smith 11768ac6da0aSPeter Brune /*E 11778ac6da0aSPeter Brune MatColoringWeightType - Type of weight scheme 11788ac6da0aSPeter Brune 11798ac6da0aSPeter Brune Not Collective 11808ac6da0aSPeter Brune 11818ac6da0aSPeter Brune + MAT_COLORING_RANDOM - Random weights 11828ac6da0aSPeter Brune . MAT_COLORING_LEXICAL - Lexical weighting based upon global numbering. 11838ac6da0aSPeter Brune - MAT_COLORING_LF - Last-first weighting. 11848ac6da0aSPeter Brune 11858ac6da0aSPeter Brune Level: intermediate 11868ac6da0aSPeter Brune 1187af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 11888ac6da0aSPeter Brune 11898ac6da0aSPeter Brune .seealso: MatCUSPSetFormat(), MatCUSPFormatOperation 11908ac6da0aSPeter Brune E*/ 1191301c140bSPeter Brune typedef enum {MAT_COLORING_WEIGHT_RANDOM,MAT_COLORING_WEIGHT_LEXICAL,MAT_COLORING_WEIGHT_LF,MAT_COLORING_WEIGHT_SL} MatColoringWeightType; 11928ac6da0aSPeter Brune 1193335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringCreate(Mat,MatColoring*); 1194d7f2a307SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringGetDegrees(Mat,PetscInt,PetscInt*); 1195335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringDestroy(MatColoring*); 1196335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringView(MatColoring,PetscViewer); 1197335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetType(MatColoring,MatColoringType); 1198335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetFromOptions(MatColoring); 1199335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetDistance(MatColoring,PetscInt); 1200335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringGetDistance(MatColoring,PetscInt*); 1201335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetMaxColors(MatColoring,PetscInt); 1202335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringGetMaxColors(MatColoring,PetscInt*); 1203335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringApply(MatColoring,ISColoring*); 1204335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringRegister(const char[],PetscErrorCode(*)(MatColoring)); 1205014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*); 12068ac6da0aSPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetWeightType(MatColoring,MatColoringWeightType); 1207c29971fcSPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetWeights(MatColoring,PetscReal*,PetscInt*); 12088ac6da0aSPeter Brune PETSC_EXTERN PetscErrorCode MatColoringCreateWeights(MatColoring,PetscReal **,PetscInt **lperm); 1209cb394dc5SBarry Smith PETSC_EXTERN PetscErrorCode MatColoringTestValid(MatColoring,ISColoring); 1210639f9d9dSBarry Smith 1211d9274352SBarry Smith /*S 1212d9274352SBarry Smith MatFDColoring - Object for computing a sparse Jacobian via finite differences 1213d9274352SBarry Smith and coloring 1214639f9d9dSBarry Smith 1215d9274352SBarry Smith Level: beginner 1216d9274352SBarry Smith 1217d9274352SBarry Smith Concepts: coloring, sparse Jacobian, finite differences 1218d9274352SBarry Smith 1219d9274352SBarry Smith .seealso: MatFDColoringCreate() 1220d9274352SBarry Smith S*/ 1221e2a1c21fSSatish Balay typedef struct _p_MatFDColoring* MatFDColoring; 1222639f9d9dSBarry Smith 1223014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringCreate(Mat,ISColoring,MatFDColoring *); 1224014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringDestroy(MatFDColoring*); 1225014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringView(MatFDColoring,PetscViewer); 1226014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*); 1227014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**); 1228014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal); 1229014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetFromOptions(MatFDColoring); 1230d1e9a80fSBarry Smith PETSC_EXTERN PetscErrorCode MatFDColoringApply(Mat,MatFDColoring,Vec,void *); 1231014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetF(MatFDColoring,Vec); 1232014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,PetscInt*[]); 1233f86b9fbaSHong Zhang PETSC_EXTERN PetscErrorCode MatFDColoringSetUp(Mat,ISColoring,MatFDColoring); 1234f86b9fbaSHong Zhang PETSC_EXTERN PetscErrorCode MatFDColoringSetBlockSize(MatFDColoring,PetscInt,PetscInt); 1235f86b9fbaSHong Zhang 1236b1683b59SHong Zhang 1237b1683b59SHong Zhang /*S 1238b9af6bddSHong Zhang MatTransposeColoring - Object for computing a sparse matrix product C=A*B^T via coloring 1239b1683b59SHong Zhang 1240b1683b59SHong Zhang Level: beginner 1241b1683b59SHong Zhang 1242b1683b59SHong Zhang Concepts: coloring, sparse matrix product 1243b1683b59SHong Zhang 1244b9af6bddSHong Zhang .seealso: MatTransposeColoringCreate() 1245b1683b59SHong Zhang S*/ 1246b9af6bddSHong Zhang typedef struct _p_MatTransposeColoring* MatTransposeColoring; 1247b1683b59SHong Zhang 1248014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeColoringCreate(Mat,ISColoring,MatTransposeColoring *); 1249014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransColoringApplySpToDen(MatTransposeColoring,Mat,Mat); 1250014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransColoringApplyDenToSp(MatTransposeColoring,Mat,Mat); 1251014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeColoringDestroy(MatTransposeColoring*); 1252b1683b59SHong Zhang 1253639f9d9dSBarry Smith /* 12540752156aSBarry Smith These routines are for partitioning matrices: currently used only 12553eda8832SBarry Smith for adjacency matrix, MatCreateMPIAdj(). 12560752156aSBarry Smith */ 1257ca161407SBarry Smith 1258d9274352SBarry Smith /*S 1259d9274352SBarry Smith MatPartitioning - Object for managing the partitioning of a matrix or graph 1260d9274352SBarry Smith 1261d9274352SBarry Smith Level: beginner 1262d9274352SBarry Smith 1263d9274352SBarry Smith Concepts: partitioning 1264d9274352SBarry Smith 1265743a1026Svictorle .seealso: MatPartitioningCreate(), MatPartitioningType 1266d9274352SBarry Smith S*/ 126791e9ee9fSBarry Smith typedef struct _p_MatPartitioning* MatPartitioning; 1268d9274352SBarry Smith 126976bdecfbSBarry Smith /*J 12708f6c3df8SBarry Smith MatPartitioningType - String with the name of a PETSc matrix partitioning 1271d9274352SBarry Smith 1272d9274352SBarry Smith Level: beginner 12730d04baf8SBarry Smith dm 1274b547a13bSvictorle .seealso: MatPartitioningCreate(), MatPartitioning 127576bdecfbSBarry Smith J*/ 127619fd82e9SBarry Smith typedef const char* MatPartitioningType; 12772692d6eeSBarry Smith #define MATPARTITIONINGCURRENT "current" 1278aa1e27eaSFande Kong #define MATPARTITIONINGAVERAGE "average" 12792692d6eeSBarry Smith #define MATPARTITIONINGSQUARE "square" 12802692d6eeSBarry Smith #define MATPARTITIONINGPARMETIS "parmetis" 12812692d6eeSBarry Smith #define MATPARTITIONINGCHACO "chaco" 12822692d6eeSBarry Smith #define MATPARTITIONINGPARTY "party" 128350d91057SBarry Smith #define MATPARTITIONINGPTSCOTCH "ptscotch" 128488d2ac2bSFande Kong #define MATPARTITIONINGHIERARCH "hierarch" 128571306c60Sjroman 1286ca161407SBarry Smith 1287014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningCreate(MPI_Comm,MatPartitioning*); 128819fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningSetType(MatPartitioning,MatPartitioningType); 1289014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetNParts(MatPartitioning,PetscInt); 1290014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetAdjacency(MatPartitioning,Mat); 1291014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]); 1292014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []); 1293014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningApply(MatPartitioning,IS*); 1294014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningDestroy(MatPartitioning*); 12952aabb6bbSBarry Smith 1296bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningRegister(const char[],PetscErrorCode (*)(MatPartitioning)); 129730de9b25SBarry Smith 1298f1144a30SSatish Balay 12992bad1931SBarry Smith 1300014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningView(MatPartitioning,PetscViewer); 1301014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetFromOptions(MatPartitioning); 130219fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningGetType(MatPartitioning,MatPartitioningType*); 1303ca161407SBarry Smith 130427538973SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetRepartition(MatPartitioning part); 1305014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetCoarseSequential(MatPartitioning); 1306014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningParmetisGetEdgeCut(MatPartitioning, PetscInt *); 13070752156aSBarry Smith 1308b6956c12SJose Roman typedef enum { MP_CHACO_MULTILEVEL=1,MP_CHACO_SPECTRAL=2,MP_CHACO_LINEAR=4,MP_CHACO_RANDOM=5,MP_CHACO_SCATTERED=6 } MPChacoGlobalType; 13096a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoGlobalTypes[]; 1310b6956c12SJose Roman typedef enum { MP_CHACO_KERNIGHAN=1,MP_CHACO_NONE=2 } MPChacoLocalType; 13116a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoLocalTypes[]; 1312b6956c12SJose Roman typedef enum { MP_CHACO_LANCZOS=0,MP_CHACO_RQI=1 } MPChacoEigenType; 13136a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoEigenTypes[]; 1314b6956c12SJose Roman 1315014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetGlobal(MatPartitioning,MPChacoGlobalType); 1316014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetGlobal(MatPartitioning,MPChacoGlobalType*); 1317014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetLocal(MatPartitioning,MPChacoLocalType); 1318014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetLocal(MatPartitioning,MPChacoLocalType*); 1319014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal); 1320014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType); 1321014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenSolver(MatPartitioning,MPChacoEigenType*); 1322014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenTol(MatPartitioning,PetscReal); 1323014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenTol(MatPartitioning,PetscReal*); 1324014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenNumber(MatPartitioning,PetscInt); 1325014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenNumber(MatPartitioning,PetscInt*); 132671306c60Sjroman 132771306c60Sjroman #define MP_PARTY_OPT "opt" 132871306c60Sjroman #define MP_PARTY_LIN "lin" 132971306c60Sjroman #define MP_PARTY_SCA "sca" 133071306c60Sjroman #define MP_PARTY_RAN "ran" 133171306c60Sjroman #define MP_PARTY_GBF "gbf" 133271306c60Sjroman #define MP_PARTY_GCF "gcf" 133371306c60Sjroman #define MP_PARTY_BUB "bub" 133471306c60Sjroman #define MP_PARTY_DEF "def" 1335014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetGlobal(MatPartitioning,const char*); 133671306c60Sjroman #define MP_PARTY_HELPFUL_SETS "hs" 133771306c60Sjroman #define MP_PARTY_KERNIGHAN_LIN "kl" 133871306c60Sjroman #define MP_PARTY_NONE "no" 1339014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetLocal(MatPartitioning,const char*); 1340014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal); 1341014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetBipart(MatPartitioning,PetscBool); 1342014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscBool); 134371306c60Sjroman 134450d91057SBarry Smith typedef enum { MP_PTSCOTCH_QUALITY,MP_PTSCOTCH_SPEED,MP_PTSCOTCH_BALANCE,MP_PTSCOTCH_SAFETY,MP_PTSCOTCH_SCALABILITY } MPPTScotchStrategyType; 13456a6fc655SJed Brown PETSC_EXTERN const char *const MPPTScotchStrategyTypes[]; 1346e0f1cffaSJose Roman 1347014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetImbalance(MatPartitioning,PetscReal); 1348014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetImbalance(MatPartitioning,PetscReal*); 1349014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetStrategy(MatPartitioning,MPPTScotchStrategyType); 1350014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetStrategy(MatPartitioning,MPPTScotchStrategyType*); 135171306c60Sjroman 1352b43b03e9SMark F. Adams /* 13536294fa2bSFande Kong * hierarchical partitioning 13546294fa2bSFande Kong */ 13554e713f55SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalGetFineparts(MatPartitioning,IS*); 13564e713f55SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalGetCoarseparts(MatPartitioning,IS*); 13574e713f55SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalSetNcoarseparts(MatPartitioning,PetscInt); 13584e713f55SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalSetNfineparts(MatPartitioning, PetscInt); 13596294fa2bSFande Kong 1360014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMeshToVertexGraph(Mat,PetscInt,Mat*); 1361014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMeshToCellGraph(Mat,PetscInt,Mat*); 1362591294e4SBarry Smith 13630752156aSBarry Smith /* 1364af0996ceSBarry Smith If you add entries here you must also add them to petsc/finclude/petscmat.h 1365d4fbbf0eSBarry Smith */ 13661c1c02c0SLois Curfman McInnes typedef enum { MATOP_SET_VALUES=0, 13671c1c02c0SLois Curfman McInnes MATOP_GET_ROW=1, 13681c1c02c0SLois Curfman McInnes MATOP_RESTORE_ROW=2, 13691c1c02c0SLois Curfman McInnes MATOP_MULT=3, 13701c1c02c0SLois Curfman McInnes MATOP_MULT_ADD=4, 13717c922b88SBarry Smith MATOP_MULT_TRANSPOSE=5, 13727c922b88SBarry Smith MATOP_MULT_TRANSPOSE_ADD=6, 13731c1c02c0SLois Curfman McInnes MATOP_SOLVE=7, 13741c1c02c0SLois Curfman McInnes MATOP_SOLVE_ADD=8, 13757c922b88SBarry Smith MATOP_SOLVE_TRANSPOSE=9, 13767c922b88SBarry Smith MATOP_SOLVE_TRANSPOSE_ADD=10, 13771c1c02c0SLois Curfman McInnes MATOP_LUFACTOR=11, 13781c1c02c0SLois Curfman McInnes MATOP_CHOLESKYFACTOR=12, 1379a714c06dSBarry Smith MATOP_SOR=13, 13801c1c02c0SLois Curfman McInnes MATOP_TRANSPOSE=14, 13811c1c02c0SLois Curfman McInnes MATOP_GETINFO=15, 13821c1c02c0SLois Curfman McInnes MATOP_EQUAL=16, 13831c1c02c0SLois Curfman McInnes MATOP_GET_DIAGONAL=17, 13841c1c02c0SLois Curfman McInnes MATOP_DIAGONAL_SCALE=18, 13851c1c02c0SLois Curfman McInnes MATOP_NORM=19, 13861c1c02c0SLois Curfman McInnes MATOP_ASSEMBLY_BEGIN=20, 13871c1c02c0SLois Curfman McInnes MATOP_ASSEMBLY_END=21, 1388d519adbfSMatthew Knepley MATOP_SET_OPTION=22, 1389d519adbfSMatthew Knepley MATOP_ZERO_ENTRIES=23, 1390d519adbfSMatthew Knepley MATOP_ZERO_ROWS=24, 1391d519adbfSMatthew Knepley MATOP_LUFACTOR_SYMBOLIC=25, 1392d519adbfSMatthew Knepley MATOP_LUFACTOR_NUMERIC=26, 1393d519adbfSMatthew Knepley MATOP_CHOLESKY_FACTOR_SYMBOLIC=27, 1394d519adbfSMatthew Knepley MATOP_CHOLESKY_FACTOR_NUMERIC=28, 1395d519adbfSMatthew Knepley MATOP_SETUP_PREALLOCATION=29, 1396d519adbfSMatthew Knepley MATOP_ILUFACTOR_SYMBOLIC=30, 1397d519adbfSMatthew Knepley MATOP_ICCFACTOR_SYMBOLIC=31, 1398a5b7ff6bSBarry Smith MATOP_GET_DIAGONAL_BLOCK=32, 13999d39f69dSJed Brown /* MATOP_PLACEHOLDER_33=33, */ 1400d519adbfSMatthew Knepley MATOP_DUPLICATE=34, 1401d519adbfSMatthew Knepley MATOP_FORWARD_SOLVE=35, 1402d519adbfSMatthew Knepley MATOP_BACKWARD_SOLVE=36, 1403d519adbfSMatthew Knepley MATOP_ILUFACTOR=37, 1404d519adbfSMatthew Knepley MATOP_ICCFACTOR=38, 1405d519adbfSMatthew Knepley MATOP_AXPY=39, 14067dae84e0SHong Zhang MATOP_CREATE_SUBMATRICES=40, 1407d519adbfSMatthew Knepley MATOP_INCREASE_OVERLAP=41, 1408d519adbfSMatthew Knepley MATOP_GET_VALUES=42, 1409d519adbfSMatthew Knepley MATOP_COPY=43, 1410d519adbfSMatthew Knepley MATOP_GET_ROW_MAX=44, 1411d519adbfSMatthew Knepley MATOP_SCALE=45, 1412d519adbfSMatthew Knepley MATOP_SHIFT=46, 141335153367SBarry Smith MATOP_DIAGONAL_SET=47, 14149d39f69dSJed Brown MATOP_ZERO_ROWS_COLUMNS=48, 14159d39f69dSJed Brown MATOP_SET_RANDOM=49, 1416d519adbfSMatthew Knepley MATOP_GET_ROW_IJ=50, 1417d519adbfSMatthew Knepley MATOP_RESTORE_ROW_IJ=51, 1418d519adbfSMatthew Knepley MATOP_GET_COLUMN_IJ=52, 1419d519adbfSMatthew Knepley MATOP_RESTORE_COLUMN_IJ=53, 1420d519adbfSMatthew Knepley MATOP_FDCOLORING_CREATE=54, 1421d519adbfSMatthew Knepley MATOP_COLORING_PATCH=55, 1422d519adbfSMatthew Knepley MATOP_SET_UNFACTORED=56, 1423d519adbfSMatthew Knepley MATOP_PERMUTE=57, 1424d519adbfSMatthew Knepley MATOP_SET_VALUES_BLOCKED=58, 14257dae84e0SHong Zhang MATOP_CREATE_SUBMATRIX=59, 1426d519adbfSMatthew Knepley MATOP_DESTROY=60, 1427d519adbfSMatthew Knepley MATOP_VIEW=61, 1428d519adbfSMatthew Knepley MATOP_CONVERT_FROM=62, 14297bab7c10SHong Zhang MATOP_MATMAT_MULT=63, 14307bab7c10SHong Zhang MATOP_MATMAT_MULT_SYMBOLIC=64, 14317bab7c10SHong Zhang MATOP_MATMAT_MULT_NUMERIC=65, 1432d519adbfSMatthew Knepley MATOP_SET_LOCAL_TO_GLOBAL_MAP=66, 1433d519adbfSMatthew Knepley MATOP_SET_VALUES_LOCAL=67, 1434d519adbfSMatthew Knepley MATOP_ZERO_ROWS_LOCAL=68, 1435d519adbfSMatthew Knepley MATOP_GET_ROW_MAX_ABS=69, 1436d519adbfSMatthew Knepley MATOP_GET_ROW_MIN_ABS=70, 1437d519adbfSMatthew Knepley MATOP_CONVERT=71, 1438d519adbfSMatthew Knepley MATOP_SET_COLORING=72, 14399d39f69dSJed Brown /* MATOP_PLACEHOLDER_73=73, */ 1440d519adbfSMatthew Knepley MATOP_SET_VALUES_ADIFOR=74, 1441d519adbfSMatthew Knepley MATOP_FD_COLORING_APPLY=75, 1442d519adbfSMatthew Knepley MATOP_SET_FROM_OPTIONS=76, 1443bda6c4cbSBarry Smith MATOP_MULT_CONSTRAINED=77, 1444bda6c4cbSBarry Smith MATOP_MULT_TRANSPOSE_CONSTRAIN=78, 14459d39f69dSJed Brown MATOP_FIND_ZERO_DIAGONALS=79, 1446d519adbfSMatthew Knepley MATOP_MULT_MULTIPLE=80, 1447d519adbfSMatthew Knepley MATOP_SOLVE_MULTIPLE=81, 1448d519adbfSMatthew Knepley MATOP_GET_INERTIA=82, 1449d519adbfSMatthew Knepley MATOP_LOAD=83, 1450d519adbfSMatthew Knepley MATOP_IS_SYMMETRIC=84, 1451d519adbfSMatthew Knepley MATOP_IS_HERMITIAN=85, 1452d519adbfSMatthew Knepley MATOP_IS_STRUCTURALLY_SYMMETRIC=86, 1453820469bcSHong Zhang MATOP_SET_VALUES_BLOCKEDLOCAL=87, 1454d519adbfSMatthew Knepley MATOP_GET_VECS=88, 1455d519adbfSMatthew Knepley MATOP_MAT_MULT=89, 1456d519adbfSMatthew Knepley MATOP_MAT_MULT_SYMBOLIC=90, 1457d519adbfSMatthew Knepley MATOP_MAT_MULT_NUMERIC=91, 1458d519adbfSMatthew Knepley MATOP_PTAP=92, 1459d519adbfSMatthew Knepley MATOP_PTAP_SYMBOLIC=93, 1460d519adbfSMatthew Knepley MATOP_PTAP_NUMERIC=94, 1461bda6c4cbSBarry Smith MATOP_MAT_TRANSPOSE_MULT=95, 1462bda6c4cbSBarry Smith MATOP_MAT_TRANSPOSE_MULT_SYMBO=96, 1463bda6c4cbSBarry Smith MATOP_MAT_TRANSPOSE_MULT_NUMER=97, 14649d39f69dSJed Brown /* MATOP_PLACEHOLDER_98=98, */ 14659d39f69dSJed Brown /* MATOP_PLACEHOLDER_99=99, */ 14669d39f69dSJed Brown /* MATOP_PLACEHOLDER_100=100, */ 14679d39f69dSJed Brown /* MATOP_PLACEHOLDER_101=101, */ 1468d519adbfSMatthew Knepley MATOP_CONJUGATE=102, 14699d39f69dSJed Brown /* MATOP_PLACEHOLDER_103=103, */ 1470d519adbfSMatthew Knepley MATOP_SET_VALUES_ROW=104, 1471d519adbfSMatthew Knepley MATOP_REAL_PART=105, 1472bda6c4cbSBarry Smith MATOP_IMAGINARY_PART=106, 1473bda6c4cbSBarry Smith MATOP_GET_ROW_UPPER_TRIANGULAR=107, 1474bda6c4cbSBarry Smith MATOP_RESTORE_ROW_UPPER_TRIANG=108, 1475bda6c4cbSBarry Smith MATOP_MAT_SOLVE=109, 14760e8d04f7SBarry Smith MATOP_GET_REDUNDANT_MATRIX=110, 1477d519adbfSMatthew Knepley MATOP_GET_ROW_MIN=111, 14780e8d04f7SBarry Smith MATOP_GET_COLUMN_VECTOR=112, 1479d519adbfSMatthew Knepley MATOP_MISSING_DIAGONAL=113, 14800e8d04f7SBarry Smith MATOP_GET_SEQ_NONZERO_STRUCTUR=114, 148189c6957cSBarry Smith MATOP_CREATE=115, 148289c6957cSBarry Smith MATOP_GET_GHOSTS=116, 14830e8d04f7SBarry Smith MATOP_GET_LOCAL_SUB_MATRIX=117, 14840e8d04f7SBarry Smith MATOP_RESTORE_LOCALSUB_MATRIX=118, 1485eeffb40dSHong Zhang MATOP_MULT_DIAGONAL_BLOCK=119, 14860e8d04f7SBarry Smith MATOP_HERMITIAN_TRANSPOSE=120, 14870e8d04f7SBarry Smith MATOP_MULT_HERMITIAN_TRANSPOSE=121, 14880e8d04f7SBarry Smith MATOP_MULT_HERMITIAN_TRANS_ADD=122, 14890e8d04f7SBarry Smith MATOP_GET_MULTI_PROC_BLOCK=123, 14909d39f69dSJed Brown MATOP_FIND_NONZERO_ROWS=124, 14910e8d04f7SBarry Smith MATOP_GET_COLUMN_NORMS=125, 14929d39f69dSJed Brown MATOP_INVERT_BLOCK_DIAGONAL=126, 14939d39f69dSJed Brown /* MATOP_PLACEHOLDER_127=127, */ 14940e8d04f7SBarry Smith MATOP_GET_SUB_MATRICES_PARALLE=128, 14952b8ad9a3SHong Zhang MATOP_SET_VALUES_BATCH=129, 14960e8d04f7SBarry Smith MATOP_TRANSPOSE_MAT_MULT=130, 14970e8d04f7SBarry Smith MATOP_TRANSPOSE_MAT_MULT_SYMBO=131, 1498e9b602ebSSatish Balay MATOP_TRANSPOSE_MAT_MULT_NUMER=132, 14990e8d04f7SBarry Smith MATOP_TRANSPOSE_COLORING_CREAT=133, 15000e8d04f7SBarry Smith MATOP_TRANS_COLORING_APPLY_SPT=134, 15010e8d04f7SBarry Smith MATOP_TRANS_COLORING_APPLY_DEN=135, 15020e8d04f7SBarry Smith MATOP_RART=136, 15030e8d04f7SBarry Smith MATOP_RART_SYMBOLIC=137, 15040e8d04f7SBarry Smith MATOP_RART_NUMERIC=138, 1505e09a3074SHong Zhang MATOP_SET_BLOCK_SIZES=139, 1506f9426fe0SMark Adams MATOP_AYPX=140, 15071919a2e2SJed Brown MATOP_RESIDUAL=141, 15089c8f2541SHong Zhang MATOP_FDCOLORING_SETUP=142, 15099c8f2541SHong Zhang MATOP_MPICONCATENATESEQ=144 1510fae171e0SBarry Smith } MatOperation; 1511014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatHasOperation(Mat,MatOperation,PetscBool *); 1512014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellSetOperation(Mat,MatOperation,void(*)(void)); 1513014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellGetOperation(Mat,MatOperation,void(**)(void)); 1514014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellSetContext(Mat,void*); 1515112a2221SBarry Smith 151690ace30eSBarry Smith /* 151790ace30eSBarry Smith Codes for matrices stored on disk. By default they are 151890ace30eSBarry Smith stored in a universal format. By changing the format with 15196a9046bcSBarry Smith PetscViewerPushFormat(viewer,PETSC_VIEWER_NATIVE); the matrices will 152090ace30eSBarry Smith be stored in a way natural for the matrix, for example dense matrices 152190ace30eSBarry Smith would be stored as dense. Matrices stored this way may only be 1522f4403165SShri Abhyankar read into matrices of the same type. 152390ace30eSBarry Smith */ 152490ace30eSBarry Smith #define MATRIX_BINARY_FORMAT_DENSE -1 152590ace30eSBarry Smith 1526014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetHashTableFactor(Mat,PetscReal); 1527014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatISGetLocalMat(Mat,Mat*); 1528014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatISSetLocalMat(Mat,Mat); 1529b7ce53b6SStefano Zampini PETSC_EXTERN PetscErrorCode MatISGetMPIXAIJ(Mat,MatReuse,Mat*); 15301f4e1ec7SBarry Smith 1531d9274352SBarry Smith /*S 1532d9274352SBarry Smith MatNullSpace - Object that removes a null space from a vector, i.e. 1533d9274352SBarry Smith orthogonalizes the vector to a subsapce 1534d9274352SBarry Smith 1535f7a9e4ceSBarry Smith Level: advanced 1536d9274352SBarry Smith 1537d9274352SBarry Smith Concepts: matrix; linear operator, null space 1538d9274352SBarry Smith 15396e1639daSBarry Smith Users manual sections: 15406e1639daSBarry Smith . sec_singular 15416e1639daSBarry Smith 1542d9274352SBarry Smith .seealso: MatNullSpaceCreate() 1543d9274352SBarry Smith S*/ 154474637425SBarry Smith typedef struct _p_MatNullSpace* MatNullSpace; 1545d9274352SBarry Smith 1546014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceCreate(MPI_Comm,PetscBool ,PetscInt,const Vec[],MatNullSpace*); 1547014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(MatNullSpace,Vec,void*),void*); 1548014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceDestroy(MatNullSpace*); 1549d0195637SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceRemove(MatNullSpace,Vec); 1550014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNullSpace(Mat, MatNullSpace *); 15515fa7ec2dSBarry Smith PETSC_EXTERN PetscErrorCode MatGetTransposeNullSpace(Mat, MatNullSpace *); 15525fa7ec2dSBarry Smith PETSC_EXTERN PetscErrorCode MatSetTransposeNullSpace(Mat,MatNullSpace); 1553014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNullSpace(Mat,MatNullSpace); 1554014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNearNullSpace(Mat,MatNullSpace); 1555014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNearNullSpace(Mat,MatNullSpace*); 1556014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceTest(MatNullSpace,Mat,PetscBool *); 1557014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceView(MatNullSpace,PetscViewer); 1558014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceGetVecs(MatNullSpace,PetscBool*,PetscInt*,const Vec**); 1559014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceCreateRigidBody(Vec,MatNullSpace*); 156074637425SBarry Smith 1561014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatReorderingSeqSBAIJ(Mat,IS); 1562014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetHashTableFactor(Mat,PetscReal); 1563014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetColumnIndices(Mat,PetscInt *); 1564014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJInvertBlockDiagonal(Mat); 15653f1d51d7SBarry Smith 1566014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMAIJ(Mat,PetscInt,Mat*); 1567014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMAIJRedimension(Mat,PetscInt,Mat*); 1568014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMAIJGetAIJ(Mat,Mat*); 1569c4f061fbSSatish Balay 1570014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatComputeExplicitOperator(Mat,Mat*); 1571b0a32e0cSBarry Smith 1572014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalScaleLocal(Mat,Vec); 157304f1ad80SBarry Smith 1574014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*); 1575014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetBase(Mat,Vec,Vec); 1576014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*); 1577014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*)); 1578014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec)); 1579014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt); 1580014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDResetHHistory(Mat); 1581014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctionError(Mat,PetscReal); 1582014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetPeriod(Mat,PetscInt); 1583014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDGetH(Mat,PetscScalar *); 1584014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetOptionsPrefix(Mat,const char[]); 1585014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*); 1586014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*); 1587e884886fSBarry Smith 15886370053bSBarry Smith /*S 15896370053bSBarry Smith MatMFFD - A data structured used to manage the computation of the h differencing parameter for matrix-free 15906370053bSBarry Smith Jacobian vector products 1591e884886fSBarry Smith 15926370053bSBarry Smith Notes: MATMFFD is a specific MatType which uses the MatMFFD data structure 15936370053bSBarry Smith 15946370053bSBarry Smith MatMFFD*() methods actually take the Mat as their first argument. Not a MatMFFD data structure 15956370053bSBarry Smith 15966370053bSBarry Smith Level: developer 15976370053bSBarry Smith 15986370053bSBarry Smith .seealso: MATMFFD, MatCreateMFFD(), MatMFFDSetFuction(), MatMFFDSetType(), MatMFFDRegister() 15996370053bSBarry Smith S*/ 1600e884886fSBarry Smith typedef struct _p_MatMFFD* MatMFFD; 1601e884886fSBarry Smith 160276bdecfbSBarry Smith /*J 1603e884886fSBarry Smith MatMFFDType - algorithm used to compute the h used in computing matrix-vector products via differencing of the function 1604e884886fSBarry Smith 1605e884886fSBarry Smith Level: beginner 1606e884886fSBarry Smith 1607e884886fSBarry Smith .seealso: MatMFFDSetType(), MatMFFDRegister() 160876bdecfbSBarry Smith J*/ 160919fd82e9SBarry Smith typedef const char* MatMFFDType; 1610e884886fSBarry Smith #define MATMFFD_DS "ds" 1611e884886fSBarry Smith #define MATMFFD_WP "wp" 1612e884886fSBarry Smith 161319fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDSetType(Mat,MatMFFDType); 1614bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDRegister(const char[],PetscErrorCode (*)(MatMFFD)); 1615e884886fSBarry Smith 1616014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDDSSetUmin(Mat,PetscReal); 1617014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDWPSetComputeNormU(Mat,PetscBool ); 1618e884886fSBarry Smith 161961ab5df0SBarry Smith PETSC_EXTERN PetscErrorCode MatFDColoringSetType(MatFDColoring,MatMFFDType); 162061ab5df0SBarry Smith 1621014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *); 1622014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *); 16237dbadf16SMatthew Knepley 162497969023SHong Zhang /* 162597969023SHong Zhang PETSc interface to MUMPS 162697969023SHong Zhang */ 162797969023SHong Zhang #ifdef PETSC_HAVE_MUMPS 1628014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMumpsSetIcntl(Mat,PetscInt,PetscInt); 1629bc6112feSHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetIcntl(Mat,PetscInt,PetscInt*); 1630b9e9ea26SSatish Balay PETSC_EXTERN PetscErrorCode MatMumpsSetCntl(Mat,PetscInt,PetscReal); 1631bc6112feSHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetCntl(Mat,PetscInt,PetscReal*); 1632bc6112feSHong Zhang 1633ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetInfo(Mat,PetscInt,PetscInt*); 1634ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetInfog(Mat,PetscInt,PetscInt*); 1635ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetRinfo(Mat,PetscInt,PetscReal*); 1636ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetRinfog(Mat,PetscInt,PetscReal*); 163797969023SHong Zhang #endif 163823a5497aSJed Brown 1639d954db57SHong Zhang /* 1640d615f992SSatish Balay PETSc interface to Mkl_Pardiso 1641b500a6b7SJose David Bermeol */ 1642b500a6b7SJose David Bermeol #ifdef PETSC_HAVE_MKL_PARDISO 1643d615f992SSatish Balay PETSC_EXTERN PetscErrorCode MatMkl_PardisoSetCntl(Mat,PetscInt,PetscInt); 1644b500a6b7SJose David Bermeol #endif 1645b500a6b7SJose David Bermeol 1646b500a6b7SJose David Bermeol /* 1647d305a81bSVasiliy Kozyrev PETSc interface to Mkl_CPardiso 1648d305a81bSVasiliy Kozyrev */ 1649d305a81bSVasiliy Kozyrev #ifdef PETSC_HAVE_MKL_CPARDISO 1650d305a81bSVasiliy Kozyrev PETSC_EXTERN PetscErrorCode MatMkl_CPardisoSetCntl(Mat,PetscInt,PetscInt); 1651d305a81bSVasiliy Kozyrev #endif 1652d305a81bSVasiliy Kozyrev 1653d305a81bSVasiliy Kozyrev /* 1654d954db57SHong Zhang PETSc interface to SUPERLU 1655d954db57SHong Zhang */ 1656d954db57SHong Zhang #ifdef PETSC_HAVE_SUPERLU 1657014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSuperluSetILUDropTol(Mat,PetscReal); 1658d954db57SHong Zhang #endif 1659d954db57SHong Zhang 1660fb785984SHong Zhang /* 1661fb785984SHong Zhang PETSc interface to SUPERLU_DIST 1662fb785984SHong Zhang */ 1663fb785984SHong Zhang #ifdef PETSC_HAVE_SUPERLU_DIST 1664fb785984SHong Zhang PETSC_EXTERN PetscErrorCode MatSuperluDistGetDiagU(Mat,PetscScalar*); 1665fb785984SHong Zhang #endif 1666fb785984SHong Zhang 1667575704cbSPieter Ghysels /* 1668575704cbSPieter Ghysels PETSc interface to STRUMPACK 1669575704cbSPieter Ghysels */ 1670575704cbSPieter Ghysels #ifdef PETSC_HAVE_STRUMPACK 1671575704cbSPieter Ghysels PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSRelCompTol(Mat,PetscReal); 1672575704cbSPieter Ghysels PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSMinSize(Mat,PetscInt); 1673575704cbSPieter Ghysels PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetColPerm(Mat,PetscBool); 1674575704cbSPieter Ghysels #endif 1675575704cbSPieter Ghysels 1676575704cbSPieter Ghysels 1677aa372e3fSPaul Mullowney #ifdef PETSC_HAVE_CUDA 16783f9c0db1SPaul Mullowney /*E 1679e057df02SPaul Mullowney MatCUSPARSEStorageFormat - indicates the storage format for CUSPARSE (GPU) 16802692e278SPaul Mullowney matrices. 1681e057df02SPaul Mullowney 1682e057df02SPaul Mullowney Not Collective 1683e057df02SPaul Mullowney 16848468deeeSKarl Rupp + MAT_CUSPARSE_CSR - Compressed Sparse Row 16852692e278SPaul Mullowney . MAT_CUSPARSE_ELL - Ellpack (requires CUDA 4.2 or later). 16862692e278SPaul Mullowney - MAT_CUSPARSE_HYB - Hybrid, a combination of Ellpack and Coordinate format (requires CUDA 4.2 or later). 1687e057df02SPaul Mullowney 1688e057df02SPaul Mullowney Level: intermediate 1689e057df02SPaul Mullowney 1690af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 1691e057df02SPaul Mullowney 1692e057df02SPaul Mullowney .seealso: MatCUSPARSESetFormat(), MatCUSPARSEFormatOperation 1693e057df02SPaul Mullowney E*/ 1694e057df02SPaul Mullowney 16953f9c0db1SPaul Mullowney typedef enum {MAT_CUSPARSE_CSR, MAT_CUSPARSE_ELL, MAT_CUSPARSE_HYB} MatCUSPARSEStorageFormat; 1696e057df02SPaul Mullowney 1697e057df02SPaul Mullowney /* these will be strings associated with enumerated type defined above */ 1698e057df02SPaul Mullowney PETSC_EXTERN const char *const MatCUSPARSEStorageFormats[]; 1699e057df02SPaul Mullowney 17003f9c0db1SPaul Mullowney /*E 1701e057df02SPaul Mullowney MatCUSPARSEFormatOperation - indicates the operation of CUSPARSE (GPU) 17022692e278SPaul Mullowney matrices whose operation should use a particular storage format. 1703e057df02SPaul Mullowney 1704e057df02SPaul Mullowney Not Collective 1705e057df02SPaul Mullowney 17068468deeeSKarl Rupp + MAT_CUSPARSE_MULT_DIAG - sets the storage format for the diagonal matrix in the parallel MatMult 17078468deeeSKarl Rupp . MAT_CUSPARSE_MULT_OFFDIAG - sets the storage format for the offdiagonal matrix in the parallel MatMult 17088468deeeSKarl Rupp . MAT_CUSPARSE_MULT - sets the storage format for the entire matrix in the serial (single GPU) MatMult 17098468deeeSKarl Rupp - MAT_CUSPARSE_ALL - sets the storage format for all CUSPARSE (GPU) matrices 1710e057df02SPaul Mullowney 1711e057df02SPaul Mullowney Level: intermediate 1712e057df02SPaul Mullowney 1713e057df02SPaul Mullowney .seealso: MatCUSPARSESetFormat(), MatCUSPARSEStorageFormat 1714e057df02SPaul Mullowney E*/ 171536d62e41SPaul Mullowney typedef enum {MAT_CUSPARSE_MULT_DIAG, MAT_CUSPARSE_MULT_OFFDIAG, MAT_CUSPARSE_MULT, MAT_CUSPARSE_ALL} MatCUSPARSEFormatOperation; 1716e057df02SPaul Mullowney 171710e9db80SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 171810e9db80SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCreateAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 1719e057df02SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCUSPARSESetFormat(Mat,MatCUSPARSEFormatOperation,MatCUSPARSEStorageFormat); 17209ae82921SPaul Mullowney #endif 17219ae82921SPaul Mullowney 172290c53902SBarry Smith #if defined(PETSC_HAVE_CUSP) 1723014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCUSP(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 1724014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateAIJCUSP(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 1725e057df02SPaul Mullowney 17263f9c0db1SPaul Mullowney /*E 1727e057df02SPaul Mullowney MatCUSPStorageFormat - indicates the storage format for CUSP (GPU) 172836d62e41SPaul Mullowney matrices. 1729e057df02SPaul Mullowney 1730e057df02SPaul Mullowney Not Collective 1731e057df02SPaul Mullowney 17328468deeeSKarl Rupp + MAT_CUSP_CSR - Compressed Sparse Row 17338468deeeSKarl Rupp . MAT_CUSP_DIA - Diagonal 17348468deeeSKarl Rupp - MAT_CUSP_ELL - Ellpack 1735e057df02SPaul Mullowney 1736e057df02SPaul Mullowney Level: intermediate 1737e057df02SPaul Mullowney 1738af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 1739e057df02SPaul Mullowney 1740e057df02SPaul Mullowney .seealso: MatCUSPSetFormat(), MatCUSPFormatOperation 1741e057df02SPaul Mullowney E*/ 17423f9c0db1SPaul Mullowney typedef enum {MAT_CUSP_CSR, MAT_CUSP_DIA, MAT_CUSP_ELL} MatCUSPStorageFormat; 1743e057df02SPaul Mullowney 1744e057df02SPaul Mullowney /* these will be strings associated with enumerated type defined above */ 1745e057df02SPaul Mullowney PETSC_EXTERN const char *const MatCUSPStorageFormats[]; 1746e057df02SPaul Mullowney 17473f9c0db1SPaul Mullowney /*E 1748e057df02SPaul Mullowney MatCUSPFormatOperation - indicates the operation of CUSP (GPU) 174936d62e41SPaul Mullowney matrices whose operation should use a particular storage format. 1750e057df02SPaul Mullowney 1751e057df02SPaul Mullowney Not Collective 1752e057df02SPaul Mullowney 17538468deeeSKarl Rupp + MAT_CUSP_MULT_DIAG - sets the storage format for the diagonal matrix in the parallel MatMult 17548468deeeSKarl Rupp . MAT_CUSP_MULT_OFFDIAG - sets the storage format for the offdiagonal matrix in the parallel MatMult 17558468deeeSKarl Rupp . MAT_CUSP_MULT - sets the storage format for the entire matrix in the serial (single GPU) MatMult 17568468deeeSKarl Rupp - MAT_CUSP_ALL - sets the storage format for all CUSP (GPU) matrices 1757e057df02SPaul Mullowney 1758e057df02SPaul Mullowney Level: intermediate 1759e057df02SPaul Mullowney 1760af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 1761e057df02SPaul Mullowney 1762e057df02SPaul Mullowney .seealso: MatCUSPSetFormat(), MatCUSPStorageFormat 1763e057df02SPaul Mullowney E*/ 176436d62e41SPaul Mullowney typedef enum {MAT_CUSP_MULT_DIAG, MAT_CUSP_MULT_OFFDIAG, MAT_CUSP_MULT, MAT_CUSP_ALL} MatCUSPFormatOperation; 1765e057df02SPaul Mullowney 1766e057df02SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCUSPSetFormat(Mat,MatCUSPFormatOperation,MatCUSPStorageFormat); 1767e057df02SPaul Mullowney #endif 176890c53902SBarry Smith 1769d67ff14aSKarl Rupp #if defined(PETSC_HAVE_VIENNACL) 1770d67ff14aSKarl Rupp PETSC_EXTERN PetscErrorCode MatCreateSeqAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 1771d67ff14aSKarl Rupp PETSC_EXTERN PetscErrorCode MatCreateAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 1772d67ff14aSKarl Rupp #endif 1773d67ff14aSKarl Rupp 177454efbe56SHong Zhang /* 177554efbe56SHong Zhang PETSc interface to FFTW 177654efbe56SHong Zhang */ 177754efbe56SHong Zhang #if defined(PETSC_HAVE_FFTW) 1778014dd563SJed Brown PETSC_EXTERN PetscErrorCode VecScatterPetscToFFTW(Mat,Vec,Vec); 1779014dd563SJed Brown PETSC_EXTERN PetscErrorCode VecScatterFFTWToPetsc(Mat,Vec,Vec); 17802a7a6963SBarry Smith PETSC_EXTERN PetscErrorCode MatCreateVecsFFTW(Mat,Vec*,Vec*,Vec*); 178154efbe56SHong Zhang #endif 178254efbe56SHong Zhang 1783014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateNest(MPI_Comm,PetscInt,const IS[],PetscInt,const IS[],const Mat[],Mat*); 1784014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSize(Mat,PetscInt*,PetscInt*); 1785014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetISs(Mat,IS[],IS[]); 1786014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetLocalISs(Mat,IS[],IS[]); 1787014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSubMats(Mat,PetscInt*,PetscInt*,Mat***); 1788014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSubMat(Mat,PetscInt,PetscInt,Mat*); 178919fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatNestSetVecType(Mat,VecType); 1790014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestSetSubMats(Mat,PetscInt,const IS[],PetscInt,const IS[],const Mat[]); 1791014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestSetSubMat(Mat,PetscInt,PetscInt,Mat); 1792d8588912SDave May 17934325cce7SMatthew G Knepley PETSC_EXTERN PetscErrorCode MatChop(Mat,PetscReal); 1794e0a58c10SMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatComputeBandwidth(Mat,PetscReal,PetscInt*); 17954325cce7SMatthew G Knepley 1796b541e6a4SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSubdomainsCreateCoalesce(Mat,PetscInt,PetscInt*,IS**); 179753dd7562SDmitry Karpeev 1798c094ef40SMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatPreallocatorPreallocate(Mat,PetscBool,Mat); 1799c094ef40SMatthew G. Knepley 1800539c167fSBarry Smith PETSC_INTERN PetscErrorCode MatHeaderMerge(Mat,Mat*); 1801539c167fSBarry Smith PETSC_EXTERN PetscErrorCode MatHeaderReplace(Mat,Mat*); 1802539c167fSBarry Smith 180323a5497aSJed Brown #endif 1804