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 3220f8fb01aSBarry Smith Any additions/changes here must also be made in src/mat/interface/dlregismat.c in MatOptions[] 323d91e6319SBarry Smith 324382164b0SBarry Smith Developer Notes: Entries that are negative need not be called collectively by all processes. 325382164b0SBarry Smith 326d91e6319SBarry Smith .seealso: MatSetOption() 327d91e6319SBarry Smith E*/ 328c5e4d11fSDmitry Karpeev typedef enum {MAT_OPTION_MIN = -3, 329c5e4d11fSDmitry Karpeev MAT_UNUSED_NONZERO_LOCATION_ERR = -2, 33092d4d306SBarry Smith MAT_ROW_ORIENTED = -1, 331382164b0SBarry Smith MAT_SYMMETRIC = 1, 332382164b0SBarry Smith MAT_STRUCTURALLY_SYMMETRIC = 2, 333382164b0SBarry Smith MAT_NEW_DIAGONALS = 3, 334382164b0SBarry Smith MAT_IGNORE_OFF_PROC_ENTRIES = 4, 335382164b0SBarry Smith MAT_USE_HASH_TABLE = 5, 336382164b0SBarry Smith MAT_KEEP_NONZERO_PATTERN = 6, 337382164b0SBarry Smith MAT_IGNORE_ZERO_ENTRIES = 7, 338382164b0SBarry Smith MAT_USE_INODES = 8, 339382164b0SBarry Smith MAT_HERMITIAN = 9, 340382164b0SBarry Smith MAT_SYMMETRY_ETERNAL = 10, 341c5e4d11fSDmitry Karpeev MAT_NEW_NONZERO_LOCATION_ERR = 11, 342382164b0SBarry Smith MAT_IGNORE_LOWER_TRIANGULAR = 12, 343382164b0SBarry Smith MAT_ERROR_LOWER_TRIANGULAR = 13, 344382164b0SBarry Smith MAT_GETROW_UPPERTRIANGULAR = 14, 345382164b0SBarry Smith MAT_SPD = 15, 34692d4d306SBarry Smith MAT_NO_OFF_PROC_ZERO_ROWS = 16, 34792d4d306SBarry Smith MAT_NO_OFF_PROC_ENTRIES = 17, 34892d4d306SBarry Smith MAT_NEW_NONZERO_LOCATIONS = 18, 349c5e4d11fSDmitry Karpeev MAT_NEW_NONZERO_ALLOCATION_ERR = 19, 350c5e4d11fSDmitry Karpeev MAT_SUBSET_OFF_PROC_ENTRIES = 20, 351c10200c1SHong Zhang MAT_SUBMAT_SINGLEIS = 21, 352957cac9fSHong Zhang MAT_STRUCTURE_ONLY = 22, 353957cac9fSHong Zhang MAT_OPTION_MAX = 23} MatOption; 354e057df02SPaul Mullowney 3550f8fb01aSBarry Smith PETSC_EXTERN const char *const *MatOptions; 356014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetOption(Mat,MatOption,PetscBool); 3577d68702bSBarry Smith PETSC_EXTERN PetscErrorCode MatGetOption(Mat,MatOption,PetscBool*); 35819fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatGetType(Mat,MatType*); 35984cb2905SBarry Smith 360014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]); 361014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 362014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 363014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowUpperTriangular(Mat); 364014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRowUpperTriangular(Mat); 365014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnVector(Mat,Vec,PetscInt); 3668c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJGetArray(Mat,PetscScalar *[]); 3678c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJRestoreArray(Mat,PetscScalar *[]); 36821e72a00SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJGetMaxRowNonzeros(Mat,PetscInt*); 369bd04181cSBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJSetValuesLocalFast(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 3704099cc6bSBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJSetType(Mat,MatType); 3714099cc6bSBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJRegister(const char[],PetscErrorCode (*)(Mat,const MatType,MatReuse,Mat *)); 3724099cc6bSBarry Smith PETSC_EXTERN PetscFunctionList MatSeqAIJList; 3738397e458SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqSBAIJGetArray(Mat,PetscScalar *[]); 3748397e458SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqSBAIJRestoreArray(Mat,PetscScalar *[]); 3758c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatDenseGetArray(Mat,PetscScalar *[]); 3768c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatDenseRestoreArray(Mat,PetscScalar *[]); 377d3042a70SBarry Smith PETSC_EXTERN PetscErrorCode MatDensePlaceArray(Mat,const PetscScalar[]); 378d3042a70SBarry Smith PETSC_EXTERN PetscErrorCode MatDenseResetArray(Mat); 379014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBlockSize(Mat,PetscInt *); 380014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSize(Mat,PetscInt); 381014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBlockSizes(Mat,PetscInt *,PetscInt *); 382014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSizes(Mat,PetscInt,PetscInt); 38333d57670SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSizesFromMats(Mat,Mat,Mat); 3841620fd73SBarry Smith 385014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMult(Mat,Vec,Vec); 386014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultDiagonalBlock(Mat,Vec,Vec); 387014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultAdd(Mat,Vec,Vec,Vec); 388014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTranspose(Mat,Vec,Vec); 389014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultHermitianTranspose(Mat,Vec,Vec); 390014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsTranspose(Mat,Mat,PetscReal,PetscBool *); 391014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitianTranspose(Mat,Mat,PetscReal,PetscBool *); 392014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeAdd(Mat,Vec,Vec,Vec); 393014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultHermitianTransposeAdd(Mat,Vec,Vec,Vec); 394014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultConstrained(Mat,Vec,Vec); 395014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeConstrained(Mat,Vec,Vec); 396014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatSolve(Mat,Mat,Mat); 397bdc285e1SStefano Zampini PETSC_EXTERN PetscErrorCode MatMatSolveTranspose(Mat,Mat,Mat); 398f9426fe0SMark Adams PETSC_EXTERN PetscErrorCode MatResidual(Mat,Vec,Vec,Vec); 399f5edf698SHong Zhang 400d91e6319SBarry Smith /*E 401d91e6319SBarry Smith MatDuplicateOption - Indicates if a duplicated sparse matrix should have 402d91e6319SBarry Smith its numerical values copied over or just its nonzero structure. 403d91e6319SBarry Smith 404d91e6319SBarry Smith Level: beginner 405d91e6319SBarry Smith 406af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 407d91e6319SBarry Smith 4080018da39SRichard Tran Mills $ MAT_DO_NOT_COPY_VALUES - Create a matrix using the same nonzero pattern as the original matrix, 4090018da39SRichard Tran Mills $ with zeros for the numerical values. 4100018da39SRichard Tran Mills $ MAT_COPY_VALUES - Create a matrix with the same nonzero pattern as the original matrix 4110018da39SRichard Tran Mills $ and with the same numerical values. 4120018da39SRichard Tran Mills $ MAT_SHARE_NONZERO_PATTERN - Create a matrix that shares the nonzero structure with the previous matrix 4130018da39SRichard Tran Mills $ and does not copy it, using zeros for the numerical values. The parent and 4140018da39SRichard Tran Mills $ child matrices will share their index (i and j) arrays, and you cannot 4150018da39SRichard Tran Mills $ insert new nonzero entries into either matrix. 4160018da39SRichard Tran Mills 4170018da39SRichard Tran Mills Notes: Many matrix types (including SeqAIJ) do not support the MAT_SHARE_NONZERO_PATTERN optimization; in 4180018da39SRichard Tran Mills this case the behavior is as if MAT_DO_NOT_COPY_VALUES has been specified. 41970dcbbb9SBarry Smith 420d91e6319SBarry Smith .seealso: MatDuplicate() 421d91e6319SBarry Smith E*/ 42270dcbbb9SBarry Smith typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES,MAT_SHARE_NONZERO_PATTERN} MatDuplicateOption; 4232e8a6d31SBarry Smith 42419fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatConvert(Mat,MatType,MatReuse,Mat*); 425014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDuplicate(Mat,MatDuplicateOption,Mat*); 42694a9d846SBarry Smith 427cb5b572fSBarry Smith 428014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCopy(Mat,Mat,MatStructure); 429014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatView(Mat,PetscViewer); 430014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsSymmetric(Mat,PetscReal,PetscBool *); 431014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsStructurallySymmetric(Mat,PetscBool *); 432014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitian(Mat,PetscReal,PetscBool *); 433014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsSymmetricKnown(Mat,PetscBool *,PetscBool *); 434014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitianKnown(Mat,PetscBool *,PetscBool *); 435014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMissingDiagonal(Mat,PetscBool *,PetscInt *); 436014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLoad(Mat, PetscViewer); 4377b80b807SBarry Smith 4381a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool *); 4391a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool *); 4401a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool *); 4411a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool *); 442d4fbbf0eSBarry Smith 443d91e6319SBarry Smith /*S 444d91e6319SBarry Smith MatInfo - Context of matrix information, used with MatGetInfo() 445d91e6319SBarry Smith 446d91e6319SBarry Smith In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE 447d91e6319SBarry Smith 448d91e6319SBarry Smith Level: intermediate 449d91e6319SBarry Smith 450d91e6319SBarry Smith Concepts: matrix^nonzero information 451d91e6319SBarry Smith 452d9274352SBarry Smith .seealso: MatGetInfo(), MatInfoType 453d91e6319SBarry Smith S*/ 4544e220ebcSLois Curfman McInnes typedef struct { 455b0a32e0cSBarry Smith PetscLogDouble block_size; /* block size */ 456b0a32e0cSBarry Smith PetscLogDouble nz_allocated,nz_used,nz_unneeded; /* number of nonzeros */ 457b0a32e0cSBarry Smith PetscLogDouble memory; /* memory allocated */ 458b0a32e0cSBarry Smith PetscLogDouble assemblies; /* number of matrix assemblies called */ 459b0a32e0cSBarry Smith PetscLogDouble mallocs; /* number of mallocs during MatSetValues() */ 460b0a32e0cSBarry Smith PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */ 461b0a32e0cSBarry Smith PetscLogDouble factor_mallocs; /* number of mallocs during factorization */ 4624e220ebcSLois Curfman McInnes } MatInfo; 4634e220ebcSLois Curfman McInnes 464d9274352SBarry Smith /*E 465d9274352SBarry Smith MatInfoType - Indicates if you want information about the local part of the matrix, 466d9274352SBarry Smith the entire parallel matrix or the maximum over all the local parts. 467d9274352SBarry Smith 468d9274352SBarry Smith Level: beginner 469d9274352SBarry Smith 470af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 471d9274352SBarry Smith 472d9274352SBarry Smith .seealso: MatGetInfo(), MatInfo 473d9274352SBarry Smith E*/ 4747b80b807SBarry Smith typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType; 475014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetInfo(Mat,MatInfoType,MatInfo*); 476014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetDiagonal(Mat,Vec); 477014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMax(Mat,Vec,PetscInt[]); 478014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMin(Mat,Vec,PetscInt[]); 479014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMaxAbs(Mat,Vec,PetscInt[]); 480014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMinAbs(Mat,Vec,PetscInt[]); 481014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowSum(Mat,Vec); 482014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTranspose(Mat,MatReuse,Mat*); 483014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatHermitianTranspose(Mat,MatReuse,Mat*); 484014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPermute(Mat,IS,IS,Mat*); 485014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalScale(Mat,Vec,Vec); 486014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalSet(Mat,Vec,InsertMode); 487a52676ecSHong Zhang 488014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatEqual(Mat,Mat,PetscBool*); 489014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultEqual(Mat,Mat,PetscInt,PetscBool*); 490014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultAddEqual(Mat,Mat,PetscInt,PetscBool*); 491014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeEqual(Mat,Mat,PetscInt,PetscBool*); 492014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscBool*); 493a52676ecSHong Zhang PETSC_EXTERN PetscErrorCode MatMatMultEqual(Mat,Mat,Mat,PetscInt,PetscBool*); 494a52676ecSHong Zhang PETSC_EXTERN PetscErrorCode MatTransposeMatMultEqual(Mat,Mat,Mat,PetscInt,PetscBool*); 4957b80b807SBarry Smith 496014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNorm(Mat,NormType,PetscReal*); 497014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnNorms(Mat,NormType,PetscReal*); 498014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroEntries(Mat); 499014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 500014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsIS(Mat,IS,PetscScalar,Vec,Vec); 501014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsStencil(Mat,PetscInt,const MatStencil [],PetscScalar,Vec,Vec); 502014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsStencil(Mat,PetscInt,const MatStencil[],PetscScalar,Vec,Vec); 503014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumns(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 504014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsIS(Mat,IS,PetscScalar,Vec,Vec); 5057b80b807SBarry Smith 506014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSize(Mat,PetscInt*,PetscInt*); 507014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalSize(Mat,PetscInt*,PetscInt*); 508014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRange(Mat,PetscInt*,PetscInt*); 509014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRanges(Mat,const PetscInt**); 510014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRangeColumn(Mat,PetscInt*,PetscInt*); 511014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRangesColumn(Mat,const PetscInt**); 512566876eaSJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipIS(Mat,IS*,IS*); 5135ef9f2a5SBarry Smith 5147dae84e0SHong Zhang PETSC_EXTERN PetscErrorCode MatCreateSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]); 515cd2534ddSHong 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);} 5167dae84e0SHong Zhang PETSC_EXTERN PetscErrorCode MatCreateSubMatricesMPI(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]); 517cd2534ddSHong 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);} 518014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroyMatrices(PetscInt,Mat *[]); 519df750dc8SHong Zhang PETSC_EXTERN PetscErrorCode MatDestroySubMatrices(PetscInt,Mat *[]); 5207dae84e0SHong Zhang PETSC_EXTERN PetscErrorCode MatCreateSubMatrix(Mat,IS,IS,MatReuse,Mat *); 521cd2534ddSHong 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);} 522014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalSubMatrix(Mat,IS,IS,Mat*); 523014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreLocalSubMatrix(Mat,IS,IS,Mat*); 524014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSeqNonzeroStructure(Mat,Mat*); 525014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroySeqNonzeroStructure(Mat*); 5267b80b807SBarry Smith 527014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJ(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*); 528014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJSymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*); 529014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJNumeric(Mat,Mat); 530014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMat(Mat,MatReuse,Mat*); 531014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*); 532014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,Mat*); 533014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetGhosts(Mat, PetscInt *,const PetscInt *[]); 5348efafbd8SBarry Smith 535014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt); 536aa1e27eaSFande Kong PETSC_EXTERN PetscErrorCode MatIncreaseOverlapSplit(Mat mat,PetscInt n,IS is[],PetscInt ov); 537d2b2a242SBarry Smith PETSC_EXTERN PetscErrorCode MatMPIAIJSetUseScalableIncreaseOverlap(Mat,PetscBool); 5387b80b807SBarry Smith 539014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*); 540014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMultSymbolic(Mat,Mat,PetscReal,Mat*); 541014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMultNumeric(Mat,Mat,Mat); 54222440eb1SKris Buschelman 5437bab7c10SHong Zhang PETSC_EXTERN PetscErrorCode MatMatMatMult(Mat,Mat,Mat,MatReuse,PetscReal,Mat*); 5442df6c5c3SPatrick Farrell PETSC_EXTERN PetscErrorCode MatGalerkin(Mat,Mat,Mat,MatReuse,PetscReal,Mat*); 5457bab7c10SHong Zhang 546014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*); 547014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAPSymbolic(Mat,Mat,PetscReal,Mat*); 548014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAPNumeric(Mat,Mat,Mat); 549014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARt(Mat,Mat,MatReuse,PetscReal,Mat*); 550014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARtSymbolic(Mat,Mat,PetscReal,Mat*); 551014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARtNumeric(Mat,Mat,Mat); 55222440eb1SKris Buschelman 553014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeMatMult(Mat,Mat,MatReuse,PetscReal,Mat*); 554014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposetMatMultSymbolic(Mat,Mat,PetscReal,Mat*); 555014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposetMatMultNumeric(Mat,Mat,Mat); 556014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMult(Mat,Mat,MatReuse,PetscReal,Mat*); 557014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMultSymbolic(Mat,Mat,PetscReal,Mat*); 558014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMultNumeric(Mat,Mat,Mat); 559bc011b1eSHong Zhang 560014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAXPY(Mat,PetscScalar,Mat,MatStructure); 561014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAYPX(Mat,PetscScalar,Mat,MatStructure); 5621c741599SHong Zhang 563014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScale(Mat,PetscScalar); 564014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShift(Mat,PetscScalar); 5657b80b807SBarry Smith 566014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping); 567014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*); 568a93ff8c4SPeter Brune PETSC_EXTERN PetscErrorCode MatGetLayouts(Mat,PetscLayout*,PetscLayout*); 569014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 570014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsLocalIS(Mat,IS,PetscScalar,Vec,Vec); 571014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 572014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocalIS(Mat,IS,PetscScalar,Vec,Vec); 573014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 574014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 575052efed2SBarry Smith 576014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStashSetInitialSize(Mat,PetscInt,PetscInt); 577014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*); 57890f02eecSBarry Smith 579014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInterpolate(Mat,Vec,Vec); 580014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInterpolateAdd(Mat,Vec,Vec,Vec); 581014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestrict(Mat,Vec,Vec); 5822a7a6963SBarry Smith PETSC_EXTERN PetscErrorCode MatCreateVecs(Mat,Vec*,Vec*); 58384cc2a5aSBarry Smith PETSC_DEPRECATED("Use MatCreateVecs()") PETSC_STATIC_INLINE PetscErrorCode MatGetVecs(Mat mat,Vec *x,Vec *y) {return MatCreateVecs(mat,x,y);} 58453cd1579SHong Zhang PETSC_EXTERN PetscErrorCode MatCreateRedundantMatrix(Mat,PetscInt,MPI_Comm,MatReuse,Mat*); 585014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetMultiProcBlock(Mat,MPI_Comm,MatReuse,Mat*); 586014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFindZeroDiagonals(Mat,IS*); 5873a062f41SBarry Smith PETSC_EXTERN PetscErrorCode MatFindOffBlockDiagonalEntries(Mat,IS*); 5889c8f2541SHong Zhang PETSC_EXTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat(MPI_Comm,Mat,PetscInt,MatReuse,Mat*); 589bd481603SBarry Smith 590bd481603SBarry Smith /*MC 5911620fd73SBarry Smith MatSetValue - Set a single entry into a matrix. 5921620fd73SBarry Smith 5931620fd73SBarry Smith Not collective 5941620fd73SBarry Smith 595a9834a7dSSatish Balay Synopsis: 596a9834a7dSSatish Balay #include <petscmat.h> 597a9834a7dSSatish Balay PetscErrorCode MatSetValue(Mat m,PetscInt row,PetscInt col,PetscScalar value,InsertMode mode) 598a9834a7dSSatish Balay 5991620fd73SBarry Smith Input Parameters: 6001620fd73SBarry Smith + m - the matrix 6011620fd73SBarry Smith . row - the row location of the entry 6021620fd73SBarry Smith . col - the column location of the entry 6031620fd73SBarry Smith . value - the value to insert 6041620fd73SBarry Smith - mode - either INSERT_VALUES or ADD_VALUES 6051620fd73SBarry Smith 6061620fd73SBarry Smith Notes: 6071620fd73SBarry Smith For efficiency one should use MatSetValues() and set several or many 6081620fd73SBarry Smith values simultaneously if possible. 6091620fd73SBarry Smith 6101620fd73SBarry Smith Level: beginner 6111620fd73SBarry Smith 6121620fd73SBarry Smith .seealso: MatSetValues(), MatSetValueLocal() 6131620fd73SBarry Smith M*/ 6141620fd73SBarry 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);} 6151620fd73SBarry Smith 6161620fd73SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatGetValue(Mat v,PetscInt i,PetscInt j,PetscScalar *va) {return MatGetValues(v,1,&i,1,&j,va);} 6171620fd73SBarry Smith 6181620fd73SBarry 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);} 6191620fd73SBarry Smith 6201620fd73SBarry Smith /*MC 6210d2a7ff2SSatish Balay MatPreallocateInitialize - Begins the block of code that will count the number of nonzeros per 622bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 623bd481603SBarry Smith 624bd481603SBarry Smith Synopsis: 625aaa7dc30SBarry Smith #include <petscmat.h> 626c1ac3661SBarry Smith PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz) 627bd481603SBarry Smith 628bd481603SBarry Smith Collective on MPI_Comm 629bd481603SBarry Smith 630bd481603SBarry Smith Input Parameters: 631bd481603SBarry Smith + comm - the communicator that will share the eventually allocated matrix 632859fcb39SBarry Smith . nrows - the number of LOCAL rows in the matrix 633859fcb39SBarry Smith - ncols - the number of LOCAL columns in the matrix 634bd481603SBarry Smith 635bd481603SBarry Smith Output Parameters: 636bd481603SBarry Smith + dnz - the array that will be passed to the matrix preallocation routines 637bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 638bd481603SBarry Smith 639bd481603SBarry Smith Level: intermediate 640bd481603SBarry Smith 641bd481603SBarry Smith Notes: 642a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 643bd481603SBarry Smith 6441d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 645bd481603SBarry Smith 646bd481603SBarry Smith Use MatPreallocateInitializeSymmetric() for symmetric matrices (MPISBAIJ matrices) 647bd481603SBarry Smith 6481620fd73SBarry Smith This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize(). 6491620fd73SBarry Smith 650bd481603SBarry Smith Concepts: preallocation^Matrix 651bd481603SBarry Smith 652d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(), 653d6e23781SBarry Smith MatPreallocateInitializeSymmetric(), MatPreallocateSymmetricSetLocalBlock() 654bd481603SBarry Smith M*/ 655c4f061fbSSatish Balay #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \ 6567c922b88SBarry Smith { \ 657a89bc211SBarry Smith PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__start,__end; \ 6587282cfc1SLisandro Dalcin _4_ierr = PetscCalloc2((size_t)__nrows,&dnz,(size_t)__nrows,&onz);CHKERRQ(_4_ierr); \ 6591795a4d1SJed Brown __start = 0; __end = __start; \ 660c1ac3661SBarry Smith _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ctmp;\ 661a89bc211SBarry Smith _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows; 6627c922b88SBarry Smith 663bd481603SBarry Smith /*MC 664bd481603SBarry Smith MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be 665bd481603SBarry Smith inserted using a local number of the rows and columns 666bd481603SBarry Smith 667bd481603SBarry Smith Synopsis: 668aaa7dc30SBarry Smith #include <petscmat.h> 669c1ac3661SBarry Smith PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 670bd481603SBarry Smith 671bd481603SBarry Smith Not Collective 672bd481603SBarry Smith 673bd481603SBarry Smith Input Parameters: 674784ac674SJed Brown + map - the row mapping from local numbering to global numbering 675bd481603SBarry Smith . nrows - the number of rows indicated 6761d73ed98SBarry Smith . rows - the indices of the rows 677784ac674SJed Brown . cmap - the column mapping from local to global numbering 678bd481603SBarry Smith . ncols - the number of columns in the matrix 679bd481603SBarry Smith . cols - the columns indicated 680bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 681bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 682bd481603SBarry Smith 683bd481603SBarry Smith Level: intermediate 684bd481603SBarry Smith 685bd481603SBarry Smith Notes: 686a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 687bd481603SBarry Smith 6881d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 689bd481603SBarry Smith 690bd481603SBarry Smith Concepts: preallocation^Matrix 691bd481603SBarry Smith 692d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(), 693c1154cd5SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocalRemoveDups() 694bd481603SBarry Smith M*/ 695784ac674SJed Brown #define MatPreallocateSetLocal(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \ 696c4f061fbSSatish Balay {\ 697c1ac3661SBarry Smith PetscInt __l;\ 698784ac674SJed Brown _4_ierr = ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\ 699784ac674SJed Brown _4_ierr = ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\ 700c4f061fbSSatish Balay for (__l=0;__l<nrows;__l++) {\ 701ef66eb69SBarry Smith _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 702c4f061fbSSatish Balay }\ 703c4f061fbSSatish Balay } 704c4f061fbSSatish Balay 705bd481603SBarry Smith /*MC 706c1154cd5SBarry Smith MatPreallocateSetLocalRemoveDups - Indicates the locations (rows and columns) in the matrix where nonzeros will be 707c1154cd5SBarry Smith inserted using a local number of the rows and columns. This version removes any duplicate columns in cols 708c1154cd5SBarry Smith 709c1154cd5SBarry Smith Synopsis: 710c1154cd5SBarry Smith #include <petscmat.h> 711c1154cd5SBarry Smith PetscErrorCode MatPreallocateSetLocalRemoveDups(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 712c1154cd5SBarry Smith 713c1154cd5SBarry Smith Not Collective 714c1154cd5SBarry Smith 715c1154cd5SBarry Smith Input Parameters: 716c1154cd5SBarry Smith + map - the row mapping from local numbering to global numbering 717c1154cd5SBarry Smith . nrows - the number of rows indicated 718c1154cd5SBarry Smith . rows - the indices of the rows (these values are mapped to the global values) 719c1154cd5SBarry Smith . cmap - the column mapping from local to global numbering 720c1154cd5SBarry Smith . ncols - the number of columns in the matrix (this value will be changed if duplicate columns are found) 721c1154cd5SBarry Smith . cols - the columns indicated (these values are mapped to the global values, they are then sorted and duplicates removed) 722c1154cd5SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 723c1154cd5SBarry Smith - ozn - the other array passed to the matrix preallocation routines 724c1154cd5SBarry Smith 725c1154cd5SBarry Smith Level: intermediate 726c1154cd5SBarry Smith 727c1154cd5SBarry Smith Notes: 728c1154cd5SBarry Smith See Users-Manual: ch_performance for more details. 729c1154cd5SBarry Smith 730c1154cd5SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 731c1154cd5SBarry Smith 732c1154cd5SBarry Smith Concepts: preallocation^Matrix 733c1154cd5SBarry Smith 734c1154cd5SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(), 735c1154cd5SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocal() 736c1154cd5SBarry Smith M*/ 737c1154cd5SBarry Smith #define MatPreallocateSetLocalRemoveDups(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \ 738c1154cd5SBarry Smith {\ 739c1154cd5SBarry Smith PetscInt __l;\ 740c1154cd5SBarry Smith _4_ierr = ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\ 741c1154cd5SBarry Smith _4_ierr = ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\ 742c1154cd5SBarry Smith _4_ierr = PetscSortRemoveDupsInt(&ncols,cols);CHKERRQ(_4_ierr);\ 743c1154cd5SBarry Smith for (__l=0;__l<nrows;__l++) {\ 744c1154cd5SBarry Smith _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 745c1154cd5SBarry Smith }\ 746c1154cd5SBarry Smith } 747c1154cd5SBarry Smith 748c1154cd5SBarry Smith /*MC 749d6e23781SBarry Smith MatPreallocateSetLocalBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be 750bd481603SBarry Smith inserted using a local number of the rows and columns 751bd481603SBarry Smith 752bd481603SBarry Smith Synopsis: 753aaa7dc30SBarry Smith #include <petscmat.h> 754d6e23781SBarry Smith PetscErrorCode MatPreallocateSetLocalBlock(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 755d6e23781SBarry Smith 756d6e23781SBarry Smith Not Collective 757d6e23781SBarry Smith 758d6e23781SBarry Smith Input Parameters: 759d6e23781SBarry Smith + map - the row mapping from local numbering to global numbering 760d6e23781SBarry Smith . nrows - the number of rows indicated 761d6e23781SBarry Smith . rows - the indices of the rows 762d6e23781SBarry Smith . cmap - the column mapping from local to global numbering 763d6e23781SBarry Smith . ncols - the number of columns in the matrix 764d6e23781SBarry Smith . cols - the columns indicated 765d6e23781SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 766d6e23781SBarry Smith - ozn - the other array passed to the matrix preallocation routines 767d6e23781SBarry Smith 768d6e23781SBarry Smith Level: intermediate 769d6e23781SBarry Smith 770d6e23781SBarry Smith Notes: 771d6e23781SBarry Smith See Users-Manual: ch_performance for more details. 772d6e23781SBarry Smith 773d6e23781SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 774d6e23781SBarry Smith 775d6e23781SBarry Smith Concepts: preallocation^Matrix 776d6e23781SBarry Smith 777d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(), 778d6e23781SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock() 779d6e23781SBarry Smith M*/ 780d6e23781SBarry Smith #define MatPreallocateSetLocalBlock(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \ 781d6e23781SBarry Smith {\ 782d6e23781SBarry Smith PetscInt __l;\ 783d6e23781SBarry Smith _4_ierr = ISLocalToGlobalMappingApplyBlock(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\ 784d6e23781SBarry Smith _4_ierr = ISLocalToGlobalMappingApplyBlock(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\ 785d6e23781SBarry Smith for (__l=0;__l<nrows;__l++) {\ 786d6e23781SBarry Smith _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 787d6e23781SBarry Smith }\ 788d6e23781SBarry Smith } 789d6e23781SBarry Smith 790d6e23781SBarry Smith /*MC 791d6e23781SBarry Smith MatPreallocateSymmetricSetLocalBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be 792d6e23781SBarry Smith inserted using a local number of the rows and columns 793d6e23781SBarry Smith 794d6e23781SBarry Smith Synopsis: 795d6e23781SBarry Smith #include <petscmat.h> 796d6e23781SBarry Smith PetscErrorCode MatPreallocateSymmetricSetLocalBlock(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 797bd481603SBarry Smith 798bd481603SBarry Smith Not Collective 799bd481603SBarry Smith 800bd481603SBarry Smith Input Parameters: 801bd481603SBarry Smith + map - the mapping between local numbering and global numbering 802bd481603SBarry Smith . nrows - the number of rows indicated 8031d73ed98SBarry Smith . rows - the indices of the rows 804bd481603SBarry Smith . ncols - the number of columns in the matrix 805bd481603SBarry Smith . cols - the columns indicated 806bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 807bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 808bd481603SBarry Smith 809bd481603SBarry Smith Level: intermediate 810bd481603SBarry Smith 811bd481603SBarry Smith Notes: 812a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 813bd481603SBarry Smith 814bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 815bd481603SBarry Smith 816bd481603SBarry Smith Concepts: preallocation^Matrix 817bd481603SBarry Smith 818d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateInitialize(), 819d6e23781SBarry Smith MatPreallocateInitialize(), MatPreallocateSetLocal() 820bd481603SBarry Smith M*/ 821d6e23781SBarry Smith #define MatPreallocateSymmetricSetLocalBlock(map,nrows,rows,ncols,cols,dnz,onz) 0;\ 822d3d32019SBarry Smith {\ 823c1ac3661SBarry Smith PetscInt __l;\ 824d6e23781SBarry Smith _4_ierr = ISLocalToGlobalMappingApplyBlock(map,nrows,rows,rows);CHKERRQ(_4_ierr);\ 825d6e23781SBarry Smith _4_ierr = ISLocalToGlobalMappingApplyBlock(map,ncols,cols,cols);CHKERRQ(_4_ierr);\ 826d3d32019SBarry Smith for (__l=0;__l<nrows;__l++) {\ 827d6e23781SBarry Smith _4_ierr = MatPreallocateSymmetricSetBlock((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 828d3d32019SBarry Smith }\ 829d3d32019SBarry Smith } 830bd481603SBarry Smith /*MC 831bd481603SBarry Smith MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be 832bd481603SBarry Smith inserted using a local number of the rows and columns 833bd481603SBarry Smith 834bd481603SBarry Smith Synopsis: 835aaa7dc30SBarry Smith #include <petscmat.h> 836c1ac3661SBarry Smith PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 837bd481603SBarry Smith 838bd481603SBarry Smith Not Collective 839bd481603SBarry Smith 840bd481603SBarry Smith Input Parameters: 84164075487SBarry Smith + row - the row 842bd481603SBarry Smith . ncols - the number of columns in the matrix 843a50f8bf6SHong Zhang - cols - the columns indicated 844a50f8bf6SHong Zhang 845a50f8bf6SHong Zhang Output Parameters: 846a50f8bf6SHong Zhang + dnz - the array that will be passed to the matrix preallocation routines 847bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 848bd481603SBarry Smith 849bd481603SBarry Smith Level: intermediate 850bd481603SBarry Smith 851bd481603SBarry Smith Notes: 852a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 853bd481603SBarry Smith 854bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 855bd481603SBarry Smith 8561620fd73SBarry Smith This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize(). 8571620fd73SBarry Smith 858bd481603SBarry Smith Concepts: preallocation^Matrix 859bd481603SBarry Smith 860d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(), 861d6e23781SBarry Smith MatPreallocateInitialize(), MatPreallocateSetLocal() 862bd481603SBarry Smith M*/ 863c4f061fbSSatish Balay #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\ 864c1ac3661SBarry Smith { PetscInt __i; \ 865e32f2f54SBarry 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);\ 866e32f2f54SBarry 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);\ 8677c922b88SBarry Smith for (__i=0; __i<nc; __i++) {\ 86864075487SBarry Smith if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \ 8697cc688d7SBarry Smith else dnz[row - __rstart]++;\ 8707c922b88SBarry Smith }\ 8717c922b88SBarry Smith } 8727c922b88SBarry Smith 873bd481603SBarry Smith /*MC 874d6e23781SBarry Smith MatPreallocateSymmetricSetBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be 875bd481603SBarry Smith inserted using a local number of the rows and columns 876bd481603SBarry Smith 877bd481603SBarry Smith Synopsis: 878aaa7dc30SBarry Smith #include <petscmat.h> 879d6e23781SBarry Smith PetscErrorCode MatPreallocateSymmetricSetBlock(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 880bd481603SBarry Smith 881bd481603SBarry Smith Not Collective 882bd481603SBarry Smith 883bd481603SBarry Smith Input Parameters: 884bd481603SBarry Smith + nrows - the number of rows indicated 8851d73ed98SBarry Smith . rows - the indices of the rows 886bd481603SBarry Smith . ncols - the number of columns in the matrix 887bd481603SBarry Smith . cols - the columns indicated 888bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 889bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 890bd481603SBarry Smith 891bd481603SBarry Smith Level: intermediate 892bd481603SBarry Smith 893bd481603SBarry Smith Notes: 894a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 895bd481603SBarry Smith 896bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 897bd481603SBarry Smith 8981620fd73SBarry Smith This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize(). 8991620fd73SBarry Smith 900bd481603SBarry Smith Concepts: preallocation^Matrix 901bd481603SBarry Smith 902d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateInitialize(), 903d6e23781SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocal() 904bd481603SBarry Smith M*/ 905d6e23781SBarry Smith #define MatPreallocateSymmetricSetBlock(row,nc,cols,dnz,onz) 0;\ 906c1ac3661SBarry Smith { PetscInt __i; \ 907d3d32019SBarry Smith for (__i=0; __i<nc; __i++) {\ 908d3d32019SBarry Smith if (cols[__i] >= __end) onz[row - __rstart]++; \ 909d3d32019SBarry Smith else if (cols[__i] >= row) dnz[row - __rstart]++;\ 910d3d32019SBarry Smith }\ 911d3d32019SBarry Smith } 912d3d32019SBarry Smith 913bd481603SBarry Smith /*MC 91416371a99SBarry Smith MatPreallocateLocation - An alternative to MatPreallocationSet() that puts the nonzero locations into the matrix if it exists 91516371a99SBarry Smith 91616371a99SBarry Smith Synopsis: 917aaa7dc30SBarry Smith #include <petscmat.h> 91816371a99SBarry Smith PetscErrorCode MatPreallocateLocations(Mat A,PetscInt row,PetscInt ncols,PetscInt *cols,PetscInt *dnz,PetscInt *onz) 91916371a99SBarry Smith 92016371a99SBarry Smith Not Collective 92116371a99SBarry Smith 92216371a99SBarry Smith Input Parameters: 92316371a99SBarry Smith . A - matrix 92416371a99SBarry Smith . row - row where values exist (must be local to this process) 92516371a99SBarry Smith . ncols - number of columns 92616371a99SBarry Smith . cols - columns with nonzeros 92716371a99SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 92816371a99SBarry Smith - ozn - the other array passed to the matrix preallocation routines 92916371a99SBarry Smith 93016371a99SBarry Smith Level: intermediate 93116371a99SBarry Smith 93216371a99SBarry Smith Notes: 933a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 93416371a99SBarry Smith 93516371a99SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 93616371a99SBarry Smith 93716371a99SBarry Smith This is a MACRO not a function because it uses a bunch of variables private to the MatPreallocation.... routines. 93816371a99SBarry Smith 93916371a99SBarry Smith Concepts: preallocation^Matrix 94016371a99SBarry Smith 941d6e23781SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(), 942d6e23781SBarry Smith MatPreallocateSymmetricSetLocalBlock() 94316371a99SBarry Smith M*/ 9440298fd71SBarry 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);} 94516371a99SBarry Smith 94616371a99SBarry Smith 94716371a99SBarry Smith /*MC 9480d2a7ff2SSatish Balay MatPreallocateFinalize - Ends the block of code that will count the number of nonzeros per 949bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 950bd481603SBarry Smith 951bd481603SBarry Smith Synopsis: 952aaa7dc30SBarry Smith #include <petscmat.h> 953c1ac3661SBarry Smith PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz) 954bd481603SBarry Smith 955bd481603SBarry Smith Collective on MPI_Comm 956bd481603SBarry Smith 957bd481603SBarry Smith Input Parameters: 95816371a99SBarry Smith + dnz - the array that was be passed to the matrix preallocation routines 959bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 960bd481603SBarry Smith 961bd481603SBarry Smith Level: intermediate 962bd481603SBarry Smith 963bd481603SBarry Smith Notes: 964a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 965bd481603SBarry Smith 966bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 967bd481603SBarry Smith 9681620fd73SBarry Smith This is a MACRO not a function because it closes the { started in MatPreallocateInitialize(). 9691620fd73SBarry Smith 970bd481603SBarry Smith Concepts: preallocation^Matrix 971bd481603SBarry Smith 972d6e23781SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(), 973d6e23781SBarry Smith MatPreallocateSymmetricSetLocalBlock() 974bd481603SBarry Smith M*/ 975a89bc211SBarry Smith #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree2(dnz,onz);CHKERRQ(_4_ierr);} 9767c922b88SBarry Smith 9777b80b807SBarry Smith /* Routines unique to particular data structures */ 978014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellGetContext(Mat,void *); 979698d4c6aSKris Buschelman 980014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInodeAdjustForInodes(Mat,IS*,IS*); 981014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *); 982698d4c6aSKris Buschelman 983014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetColumnIndices(Mat,PetscInt[]); 984014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetColumnIndices(Mat,PetscInt[]); 985014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 986014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 987014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 988014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJFromTriple(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*,PetscInt,PetscBool); 9897b80b807SBarry Smith 990a96a251dSBarry Smith #define MAT_SKIP_ALLOCATION -4 991a96a251dSBarry Smith 992014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]); 993014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]); 994014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]); 995273d9f13SBarry Smith 996014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 997014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 998014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 999014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []); 1000014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 1001014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]); 1002014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 1003014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]); 1004014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIDenseSetPreallocation(Mat,PetscScalar[]); 1005014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqDenseSetPreallocation(Mat,PetscScalar[]); 10069230625dSJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,const PetscInt*[]); 10079230625dSJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,const PetscInt*[]); 1008014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAdjCreateNonemptySubcommMat(Mat,Mat*); 1009273d9f13SBarry Smith 10102e1947a5SStefano Zampini PETSC_EXTERN PetscErrorCode MatISSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 1011cf0a3239SStefano Zampini PETSC_EXTERN PetscErrorCode MatISSetUpSF(Mat); 1012014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqDenseSetLDA(Mat,PetscInt); 1013014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDenseGetLocalMatrix(Mat,Mat*); 10141b807ce4Svictorle 1015014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStoreValues(Mat); 1016014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRetrieveValues(Mat); 10172e8a6d31SBarry Smith 1018014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDAADSetCtx(Mat,void*); 10193a7fca6bSBarry Smith 1020014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFindNonzeroRows(Mat,IS*); 10217cfe41e4SPatrick Farrell PETSC_EXTERN PetscErrorCode MatFindZeroRows(Mat,IS*); 10227b80b807SBarry Smith /* 10237b80b807SBarry Smith These routines are not usually accessed directly, rather solving is 102494b7f48cSBarry Smith done through the KSP and PC interfaces. 10257b80b807SBarry Smith */ 10267b80b807SBarry Smith 102776bdecfbSBarry Smith /*J 10288f6c3df8SBarry Smith MatOrderingType - String with the name of a PETSc matrix ordering 1029d9274352SBarry Smith 1030d9274352SBarry Smith Level: beginner 1031d9274352SBarry Smith 1032d9274352SBarry Smith .seealso: MatGetOrdering() 103376bdecfbSBarry Smith J*/ 103419fd82e9SBarry Smith typedef const char* MatOrderingType; 10352692d6eeSBarry Smith #define MATORDERINGNATURAL "natural" 10362692d6eeSBarry Smith #define MATORDERINGND "nd" 10372692d6eeSBarry Smith #define MATORDERING1WD "1wd" 10382692d6eeSBarry Smith #define MATORDERINGRCM "rcm" 10392692d6eeSBarry Smith #define MATORDERINGQMD "qmd" 10402692d6eeSBarry Smith #define MATORDERINGROWLENGTH "rowlength" 1041510184a7SJed Brown #define MATORDERINGWBM "wbm" 1042fc1bef75SJed Brown #define MATORDERINGSPECTRAL "spectral" 1043312e372aSBarry Smith #define MATORDERINGAMD "amd" /* only works if UMFPACK is installed with PETSc */ 1044b12f92e5SBarry Smith 104519fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatGetOrdering(Mat,MatOrderingType,IS*,IS*); 1046140e18c1SBarry Smith PETSC_EXTERN PetscErrorCode MatGetOrderingList(PetscFunctionList*); 1047bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatOrderingRegister(const char[],PetscErrorCode(*)(Mat,MatOrderingType,IS*,IS*)); 1048140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatOrderingList; 1049d4fbbf0eSBarry Smith 1050014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS); 1051fc1bef75SJed Brown PETSC_EXTERN PetscErrorCode MatCreateLaplacian(Mat,PetscReal,PetscBool,Mat*); 1052a2ce50c7SBarry Smith 1053d91e6319SBarry Smith /*S 1054d90ac83dSHong Zhang MatFactorShiftType - Numeric Shift. 1055d90ac83dSHong Zhang 1056d90ac83dSHong Zhang Level: beginner 1057d90ac83dSHong Zhang 1058d90ac83dSHong Zhang S*/ 1059d90ac83dSHong Zhang typedef enum {MAT_SHIFT_NONE,MAT_SHIFT_NONZERO,MAT_SHIFT_POSITIVE_DEFINITE,MAT_SHIFT_INBLOCKS} MatFactorShiftType; 10606a6fc655SJed Brown PETSC_EXTERN const char *const MatFactorShiftTypes[]; 10615e9742b9SJed Brown PETSC_EXTERN const char *const MatFactorShiftTypesDetail[]; 1062d90ac83dSHong Zhang 1063d90ac83dSHong Zhang /*S 10649aa7eafdSHong Zhang MatFactorError - indicates what type of error in matrix factor 10659aa7eafdSHong Zhang 10669aa7eafdSHong Zhang Level: beginner 10679aa7eafdSHong Zhang 106800e125f8SBarry Smith Developer Notes: Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 106900e125f8SBarry Smith 107000e125f8SBarry Smith .seealso: MatGetFactor() 10719aa7eafdSHong Zhang S*/ 10725cd7cf9dSHong Zhang typedef enum {MAT_FACTOR_NOERROR,MAT_FACTOR_STRUCT_ZEROPIVOT,MAT_FACTOR_NUMERIC_ZEROPIVOT,MAT_FACTOR_OUTMEMORY,MAT_FACTOR_OTHER} MatFactorError; 10739aa7eafdSHong Zhang 107400e125f8SBarry Smith PETSC_EXTERN PetscErrorCode MatFactorGetError(Mat,MatFactorError*); 1075b8b68cfdSBarry Smith PETSC_EXTERN PetscErrorCode MatFactorClearError(Mat); 10767b6c816cSBarry Smith PETSC_EXTERN PetscErrorCode MatFactorGetErrorZeroPivot(Mat,PetscReal*,PetscInt*); 107700e125f8SBarry Smith 10789aa7eafdSHong Zhang /*S 1079422a814eSBarry Smith MatFactorInfo - Data passed into the matrix factorization routines, and information about the resulting factorization 1080b00f7748SHong Zhang 108161cad860SBarry Smith In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE, that is use 108261cad860SBarry Smith $ MatFactorInfo info(MAT_FACTORINFO_SIZE) 1083b00f7748SHong Zhang 108415e8a5b3SHong Zhang Notes: These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC. 1085b00f7748SHong Zhang 108661cad860SBarry Smith You can use MatFactorInfoInitialize() to set default values. 108761cad860SBarry Smith 1088b00f7748SHong Zhang Level: developer 1089b00f7748SHong Zhang 1090d7d82daaSBarry Smith .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(), 1091d7d82daaSBarry Smith MatFactorInfoInitialize() 1092b00f7748SHong Zhang 1093b00f7748SHong Zhang S*/ 1094b00f7748SHong Zhang typedef struct { 109515e8a5b3SHong Zhang PetscReal diagonal_fill; /* force diagonal to fill in if initially not filled */ 109685317021SBarry Smith PetscReal usedt; 109715e8a5b3SHong Zhang PetscReal dt; /* drop tolerance */ 1098b00f7748SHong Zhang PetscReal dtcol; /* tolerance for pivoting */ 109915e8a5b3SHong Zhang PetscReal dtcount; /* maximum nonzeros to be allowed per row */ 110067e28bfeSBarry Smith PetscReal fill; /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */ 1101348344bbSBarry Smith PetscReal levels; /* ICC/ILU(levels) */ 1102bcd9e38bSBarry Smith PetscReal pivotinblocks; /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0 1103bcd9e38bSBarry Smith factorization may be faster if do not pivot */ 110415e8a5b3SHong Zhang PetscReal zeropivot; /* pivot is called zero if less than this */ 1105f4db908eSBarry Smith PetscReal shifttype; /* type of shift added to matrix factor to prevent zero pivots */ 1106f4db908eSBarry Smith PetscReal shiftamount; /* how large the shift is */ 110715e8a5b3SHong Zhang } MatFactorInfo; 1108ffa6d0a5SLois Curfman McInnes 1109014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFactorInfoInitialize(MatFactorInfo*); 11109a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatCholeskyFactor(Mat,IS,const MatFactorInfo*); 11119a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatCholeskyFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*); 11129a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatCholeskyFactorNumeric(Mat,Mat,const MatFactorInfo*); 11139a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatLUFactor(Mat,IS,IS,const MatFactorInfo*); 11149a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatILUFactor(Mat,IS,IS,const MatFactorInfo*); 11159a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatLUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*); 11169a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatILUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*); 11179a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatICCFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*); 11189a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatICCFactor(Mat,IS,const MatFactorInfo*); 11199a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatLUFactorNumeric(Mat,Mat,const MatFactorInfo*); 1120014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*); 1121014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolve(Mat,Vec,Vec); 1122014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatForwardSolve(Mat,Vec,Vec); 1123014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatBackwardSolve(Mat,Vec,Vec); 1124014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveAdd(Mat,Vec,Vec,Vec); 1125014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveTranspose(Mat,Vec,Vec); 1126014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveTransposeAdd(Mat,Vec,Vec,Vec); 1127014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolves(Mat,Vecs,Vecs); 11287c2f51b8SStefano Zampini PETSC_EXTERN PetscErrorCode MatSetUnfactored(Mat); 11297c2f51b8SStefano Zampini 11307c2f51b8SStefano Zampini typedef enum {MAT_FACTOR_SCHUR_UNFACTORED, MAT_FACTOR_SCHUR_FACTORED, MAT_FACTOR_SCHUR_INVERTED} MatFactorSchurStatus; 11318e7ba810SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorSetSchurIS(Mat,IS); 11327c2f51b8SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorGetSchurComplement(Mat,Mat*,MatFactorSchurStatus*); 11337c2f51b8SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorRestoreSchurComplement(Mat,Mat*,MatFactorSchurStatus); 11345a05ddb0SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorInvertSchurComplement(Mat); 11357c2f51b8SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorCreateSchurComplement(Mat,Mat*,MatFactorSchurStatus*); 11365a05ddb0SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorSolveSchurComplement(Mat,Vec,Vec); 11375a05ddb0SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorSolveSchurComplementTranspose(Mat,Vec,Vec); 11386dba178dSStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorFactorizeSchurComplement(Mat); 1139bb5a7306SBarry Smith 1140d91e6319SBarry Smith /*E 1141d91e6319SBarry Smith MatSORType - What type of (S)SOR to perform 1142bb1eb677SSatish Balay 1143d91e6319SBarry Smith Level: beginner 1144d91e6319SBarry Smith 1145d9274352SBarry Smith May be bitwise ORd together 1146d9274352SBarry Smith 1147af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 1148d91e6319SBarry Smith 11494e7234bfSBarry Smith MatSORType may be bitwise ORd together, so do not change the numbers 11504e7234bfSBarry Smith 115141f059aeSBarry Smith .seealso: MatSOR() 1152d91e6319SBarry Smith E*/ 1153ee50ffe9SBarry Smith typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3, 1154ee50ffe9SBarry Smith SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8, 1155ee50ffe9SBarry Smith SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16, 115684cb2905SBarry Smith SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType; 1157014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSOR(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec); 11588ed539a5SBarry Smith 1159d4fbbf0eSBarry Smith /* 1160639f9d9dSBarry Smith These routines are for efficiently computing Jacobians via finite differences. 1161639f9d9dSBarry Smith */ 1162b12f92e5SBarry Smith 11637086a01eSPeter Brune /*S 11647086a01eSPeter Brune MatColoring - Object for managing the coloring of matrices. 11657086a01eSPeter Brune 11667086a01eSPeter Brune Level: beginner 11677086a01eSPeter Brune 11687086a01eSPeter Brune Concepts: matrix, coloring 11697086a01eSPeter Brune 11707086a01eSPeter Brune .seealso: MatFDColoringCreate() ISColoring MatFDColoring 11717086a01eSPeter Brune S*/ 11727086a01eSPeter Brune typedef struct _p_MatColoring* MatColoring; 117376bdecfbSBarry Smith /*J 11748f6c3df8SBarry Smith MatColoringType - String with the name of a PETSc matrix coloring 1175d9274352SBarry Smith 1176d9274352SBarry Smith Level: beginner 1177d9274352SBarry Smith 11787086a01eSPeter Brune .seealso: MatColoringSetType(), MatColoring 117976bdecfbSBarry Smith J*/ 11807086a01eSPeter Brune 118119fd82e9SBarry Smith typedef const char* MatColoringType; 11825567d71aSPeter Brune #define MATCOLORINGJP "jp" 1183b9acb617SPeter Brune #define MATCOLORINGPOWER "power" 11842692d6eeSBarry Smith #define MATCOLORINGNATURAL "natural" 11852692d6eeSBarry Smith #define MATCOLORINGSL "sl" 11862692d6eeSBarry Smith #define MATCOLORINGLF "lf" 11872692d6eeSBarry Smith #define MATCOLORINGID "id" 11888121bdceSPeter Brune #define MATCOLORINGGREEDY "greedy" 1189b12f92e5SBarry Smith 11908ac6da0aSPeter Brune /*E 11918ac6da0aSPeter Brune MatColoringWeightType - Type of weight scheme 11928ac6da0aSPeter Brune 11938ac6da0aSPeter Brune Not Collective 11948ac6da0aSPeter Brune 11958ac6da0aSPeter Brune + MAT_COLORING_RANDOM - Random weights 11968ac6da0aSPeter Brune . MAT_COLORING_LEXICAL - Lexical weighting based upon global numbering. 11978ac6da0aSPeter Brune - MAT_COLORING_LF - Last-first weighting. 11988ac6da0aSPeter Brune 11998ac6da0aSPeter Brune Level: intermediate 12008ac6da0aSPeter Brune 1201af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 12028ac6da0aSPeter Brune 12038ac6da0aSPeter Brune .seealso: MatCUSPSetFormat(), MatCUSPFormatOperation 12048ac6da0aSPeter Brune E*/ 1205301c140bSPeter Brune typedef enum {MAT_COLORING_WEIGHT_RANDOM,MAT_COLORING_WEIGHT_LEXICAL,MAT_COLORING_WEIGHT_LF,MAT_COLORING_WEIGHT_SL} MatColoringWeightType; 12068ac6da0aSPeter Brune 1207335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringCreate(Mat,MatColoring*); 1208d7f2a307SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringGetDegrees(Mat,PetscInt,PetscInt*); 1209335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringDestroy(MatColoring*); 1210335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringView(MatColoring,PetscViewer); 1211335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetType(MatColoring,MatColoringType); 1212335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetFromOptions(MatColoring); 1213335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetDistance(MatColoring,PetscInt); 1214335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringGetDistance(MatColoring,PetscInt*); 1215335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetMaxColors(MatColoring,PetscInt); 1216335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringGetMaxColors(MatColoring,PetscInt*); 1217335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringApply(MatColoring,ISColoring*); 1218335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringRegister(const char[],PetscErrorCode(*)(MatColoring)); 1219014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*); 12208ac6da0aSPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetWeightType(MatColoring,MatColoringWeightType); 1221c29971fcSPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetWeights(MatColoring,PetscReal*,PetscInt*); 12228ac6da0aSPeter Brune PETSC_EXTERN PetscErrorCode MatColoringCreateWeights(MatColoring,PetscReal **,PetscInt **lperm); 1223cb394dc5SBarry Smith PETSC_EXTERN PetscErrorCode MatColoringTestValid(MatColoring,ISColoring); 1224639f9d9dSBarry Smith 1225d9274352SBarry Smith /*S 1226d9274352SBarry Smith MatFDColoring - Object for computing a sparse Jacobian via finite differences 1227d9274352SBarry Smith and coloring 1228639f9d9dSBarry Smith 1229d9274352SBarry Smith Level: beginner 1230d9274352SBarry Smith 1231d9274352SBarry Smith Concepts: coloring, sparse Jacobian, finite differences 1232d9274352SBarry Smith 1233d9274352SBarry Smith .seealso: MatFDColoringCreate() 1234d9274352SBarry Smith S*/ 1235e2a1c21fSSatish Balay typedef struct _p_MatFDColoring* MatFDColoring; 1236639f9d9dSBarry Smith 1237014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringCreate(Mat,ISColoring,MatFDColoring *); 1238014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringDestroy(MatFDColoring*); 1239014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringView(MatFDColoring,PetscViewer); 1240014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*); 1241014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**); 1242014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal); 1243014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetFromOptions(MatFDColoring); 1244d1e9a80fSBarry Smith PETSC_EXTERN PetscErrorCode MatFDColoringApply(Mat,MatFDColoring,Vec,void *); 1245014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetF(MatFDColoring,Vec); 1246edaa7c33SBarry Smith PETSC_EXTERN PetscErrorCode MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,const PetscInt*[]); 1247f86b9fbaSHong Zhang PETSC_EXTERN PetscErrorCode MatFDColoringSetUp(Mat,ISColoring,MatFDColoring); 1248f86b9fbaSHong Zhang PETSC_EXTERN PetscErrorCode MatFDColoringSetBlockSize(MatFDColoring,PetscInt,PetscInt); 1249f86b9fbaSHong Zhang 1250b1683b59SHong Zhang 1251b1683b59SHong Zhang /*S 1252b9af6bddSHong Zhang MatTransposeColoring - Object for computing a sparse matrix product C=A*B^T via coloring 1253b1683b59SHong Zhang 1254b1683b59SHong Zhang Level: beginner 1255b1683b59SHong Zhang 1256b1683b59SHong Zhang Concepts: coloring, sparse matrix product 1257b1683b59SHong Zhang 1258b9af6bddSHong Zhang .seealso: MatTransposeColoringCreate() 1259b1683b59SHong Zhang S*/ 1260b9af6bddSHong Zhang typedef struct _p_MatTransposeColoring* MatTransposeColoring; 1261b1683b59SHong Zhang 1262014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeColoringCreate(Mat,ISColoring,MatTransposeColoring *); 1263014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransColoringApplySpToDen(MatTransposeColoring,Mat,Mat); 1264014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransColoringApplyDenToSp(MatTransposeColoring,Mat,Mat); 1265014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeColoringDestroy(MatTransposeColoring*); 1266b1683b59SHong Zhang 1267639f9d9dSBarry Smith /* 12680752156aSBarry Smith These routines are for partitioning matrices: currently used only 12693eda8832SBarry Smith for adjacency matrix, MatCreateMPIAdj(). 12700752156aSBarry Smith */ 1271ca161407SBarry Smith 1272d9274352SBarry Smith /*S 1273d9274352SBarry Smith MatPartitioning - Object for managing the partitioning of a matrix or graph 1274d9274352SBarry Smith 1275d9274352SBarry Smith Level: beginner 1276d9274352SBarry Smith 1277d9274352SBarry Smith Concepts: partitioning 1278d9274352SBarry Smith 1279743a1026Svictorle .seealso: MatPartitioningCreate(), MatPartitioningType 1280d9274352SBarry Smith S*/ 128191e9ee9fSBarry Smith typedef struct _p_MatPartitioning* MatPartitioning; 1282d9274352SBarry Smith 128376bdecfbSBarry Smith /*J 12848f6c3df8SBarry Smith MatPartitioningType - String with the name of a PETSc matrix partitioning 1285d9274352SBarry Smith 1286d9274352SBarry Smith Level: beginner 12870d04baf8SBarry Smith dm 1288b547a13bSvictorle .seealso: MatPartitioningCreate(), MatPartitioning 128976bdecfbSBarry Smith J*/ 129019fd82e9SBarry Smith typedef const char* MatPartitioningType; 12912692d6eeSBarry Smith #define MATPARTITIONINGCURRENT "current" 1292aa1e27eaSFande Kong #define MATPARTITIONINGAVERAGE "average" 12932692d6eeSBarry Smith #define MATPARTITIONINGSQUARE "square" 12942692d6eeSBarry Smith #define MATPARTITIONINGPARMETIS "parmetis" 12952692d6eeSBarry Smith #define MATPARTITIONINGCHACO "chaco" 12962692d6eeSBarry Smith #define MATPARTITIONINGPARTY "party" 129750d91057SBarry Smith #define MATPARTITIONINGPTSCOTCH "ptscotch" 129888d2ac2bSFande Kong #define MATPARTITIONINGHIERARCH "hierarch" 129971306c60Sjroman 1300ca161407SBarry Smith 1301014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningCreate(MPI_Comm,MatPartitioning*); 130219fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningSetType(MatPartitioning,MatPartitioningType); 1303014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetNParts(MatPartitioning,PetscInt); 1304014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetAdjacency(MatPartitioning,Mat); 1305014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]); 1306014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []); 1307014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningApply(MatPartitioning,IS*); 1308014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningDestroy(MatPartitioning*); 13092aabb6bbSBarry Smith 1310bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningRegister(const char[],PetscErrorCode (*)(MatPartitioning)); 131130de9b25SBarry Smith 1312f1144a30SSatish Balay 13132bad1931SBarry Smith 1314014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningView(MatPartitioning,PetscViewer); 1315014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetFromOptions(MatPartitioning); 131619fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningGetType(MatPartitioning,MatPartitioningType*); 1317ca161407SBarry Smith 131827538973SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetRepartition(MatPartitioning part); 1319014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetCoarseSequential(MatPartitioning); 1320014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningParmetisGetEdgeCut(MatPartitioning, PetscInt *); 13210752156aSBarry Smith 1322b6956c12SJose Roman typedef enum { MP_CHACO_MULTILEVEL=1,MP_CHACO_SPECTRAL=2,MP_CHACO_LINEAR=4,MP_CHACO_RANDOM=5,MP_CHACO_SCATTERED=6 } MPChacoGlobalType; 13236a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoGlobalTypes[]; 1324b6956c12SJose Roman typedef enum { MP_CHACO_KERNIGHAN=1,MP_CHACO_NONE=2 } MPChacoLocalType; 13256a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoLocalTypes[]; 1326b6956c12SJose Roman typedef enum { MP_CHACO_LANCZOS=0,MP_CHACO_RQI=1 } MPChacoEigenType; 13276a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoEigenTypes[]; 1328b6956c12SJose Roman 1329014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetGlobal(MatPartitioning,MPChacoGlobalType); 1330014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetGlobal(MatPartitioning,MPChacoGlobalType*); 1331014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetLocal(MatPartitioning,MPChacoLocalType); 1332014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetLocal(MatPartitioning,MPChacoLocalType*); 1333014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal); 1334014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType); 1335014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenSolver(MatPartitioning,MPChacoEigenType*); 1336014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenTol(MatPartitioning,PetscReal); 1337014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenTol(MatPartitioning,PetscReal*); 1338014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenNumber(MatPartitioning,PetscInt); 1339014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenNumber(MatPartitioning,PetscInt*); 134071306c60Sjroman 134171306c60Sjroman #define MP_PARTY_OPT "opt" 134271306c60Sjroman #define MP_PARTY_LIN "lin" 134371306c60Sjroman #define MP_PARTY_SCA "sca" 134471306c60Sjroman #define MP_PARTY_RAN "ran" 134571306c60Sjroman #define MP_PARTY_GBF "gbf" 134671306c60Sjroman #define MP_PARTY_GCF "gcf" 134771306c60Sjroman #define MP_PARTY_BUB "bub" 134871306c60Sjroman #define MP_PARTY_DEF "def" 1349014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetGlobal(MatPartitioning,const char*); 135071306c60Sjroman #define MP_PARTY_HELPFUL_SETS "hs" 135171306c60Sjroman #define MP_PARTY_KERNIGHAN_LIN "kl" 135271306c60Sjroman #define MP_PARTY_NONE "no" 1353014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetLocal(MatPartitioning,const char*); 1354014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal); 1355014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetBipart(MatPartitioning,PetscBool); 1356014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscBool); 135771306c60Sjroman 135850d91057SBarry Smith typedef enum { MP_PTSCOTCH_QUALITY,MP_PTSCOTCH_SPEED,MP_PTSCOTCH_BALANCE,MP_PTSCOTCH_SAFETY,MP_PTSCOTCH_SCALABILITY } MPPTScotchStrategyType; 13596a6fc655SJed Brown PETSC_EXTERN const char *const MPPTScotchStrategyTypes[]; 1360e0f1cffaSJose Roman 1361014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetImbalance(MatPartitioning,PetscReal); 1362014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetImbalance(MatPartitioning,PetscReal*); 1363014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetStrategy(MatPartitioning,MPPTScotchStrategyType); 1364014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetStrategy(MatPartitioning,MPPTScotchStrategyType*); 136571306c60Sjroman 1366b43b03e9SMark F. Adams /* 13676294fa2bSFande Kong * hierarchical partitioning 13686294fa2bSFande Kong */ 13694e713f55SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalGetFineparts(MatPartitioning,IS*); 13704e713f55SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalGetCoarseparts(MatPartitioning,IS*); 13714e713f55SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalSetNcoarseparts(MatPartitioning,PetscInt); 13724e713f55SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalSetNfineparts(MatPartitioning, PetscInt); 13736294fa2bSFande Kong 1374014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMeshToVertexGraph(Mat,PetscInt,Mat*); 1375014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMeshToCellGraph(Mat,PetscInt,Mat*); 1376591294e4SBarry Smith 13770752156aSBarry Smith /* 1378af0996ceSBarry Smith If you add entries here you must also add them to petsc/finclude/petscmat.h 1379d4fbbf0eSBarry Smith */ 13801c1c02c0SLois Curfman McInnes typedef enum { MATOP_SET_VALUES=0, 13811c1c02c0SLois Curfman McInnes MATOP_GET_ROW=1, 13821c1c02c0SLois Curfman McInnes MATOP_RESTORE_ROW=2, 13831c1c02c0SLois Curfman McInnes MATOP_MULT=3, 13841c1c02c0SLois Curfman McInnes MATOP_MULT_ADD=4, 13857c922b88SBarry Smith MATOP_MULT_TRANSPOSE=5, 13867c922b88SBarry Smith MATOP_MULT_TRANSPOSE_ADD=6, 13871c1c02c0SLois Curfman McInnes MATOP_SOLVE=7, 13881c1c02c0SLois Curfman McInnes MATOP_SOLVE_ADD=8, 13897c922b88SBarry Smith MATOP_SOLVE_TRANSPOSE=9, 13907c922b88SBarry Smith MATOP_SOLVE_TRANSPOSE_ADD=10, 13911c1c02c0SLois Curfman McInnes MATOP_LUFACTOR=11, 13921c1c02c0SLois Curfman McInnes MATOP_CHOLESKYFACTOR=12, 1393a714c06dSBarry Smith MATOP_SOR=13, 13941c1c02c0SLois Curfman McInnes MATOP_TRANSPOSE=14, 13951c1c02c0SLois Curfman McInnes MATOP_GETINFO=15, 13961c1c02c0SLois Curfman McInnes MATOP_EQUAL=16, 13971c1c02c0SLois Curfman McInnes MATOP_GET_DIAGONAL=17, 13981c1c02c0SLois Curfman McInnes MATOP_DIAGONAL_SCALE=18, 13991c1c02c0SLois Curfman McInnes MATOP_NORM=19, 14001c1c02c0SLois Curfman McInnes MATOP_ASSEMBLY_BEGIN=20, 14011c1c02c0SLois Curfman McInnes MATOP_ASSEMBLY_END=21, 1402d519adbfSMatthew Knepley MATOP_SET_OPTION=22, 1403d519adbfSMatthew Knepley MATOP_ZERO_ENTRIES=23, 1404d519adbfSMatthew Knepley MATOP_ZERO_ROWS=24, 1405d519adbfSMatthew Knepley MATOP_LUFACTOR_SYMBOLIC=25, 1406d519adbfSMatthew Knepley MATOP_LUFACTOR_NUMERIC=26, 1407d519adbfSMatthew Knepley MATOP_CHOLESKY_FACTOR_SYMBOLIC=27, 1408d519adbfSMatthew Knepley MATOP_CHOLESKY_FACTOR_NUMERIC=28, 1409d519adbfSMatthew Knepley MATOP_SETUP_PREALLOCATION=29, 1410d519adbfSMatthew Knepley MATOP_ILUFACTOR_SYMBOLIC=30, 1411d519adbfSMatthew Knepley MATOP_ICCFACTOR_SYMBOLIC=31, 1412a5b7ff6bSBarry Smith MATOP_GET_DIAGONAL_BLOCK=32, 14139d39f69dSJed Brown /* MATOP_PLACEHOLDER_33=33, */ 1414d519adbfSMatthew Knepley MATOP_DUPLICATE=34, 1415d519adbfSMatthew Knepley MATOP_FORWARD_SOLVE=35, 1416d519adbfSMatthew Knepley MATOP_BACKWARD_SOLVE=36, 1417d519adbfSMatthew Knepley MATOP_ILUFACTOR=37, 1418d519adbfSMatthew Knepley MATOP_ICCFACTOR=38, 1419d519adbfSMatthew Knepley MATOP_AXPY=39, 14207dae84e0SHong Zhang MATOP_CREATE_SUBMATRICES=40, 1421d519adbfSMatthew Knepley MATOP_INCREASE_OVERLAP=41, 1422d519adbfSMatthew Knepley MATOP_GET_VALUES=42, 1423d519adbfSMatthew Knepley MATOP_COPY=43, 1424d519adbfSMatthew Knepley MATOP_GET_ROW_MAX=44, 1425d519adbfSMatthew Knepley MATOP_SCALE=45, 1426d519adbfSMatthew Knepley MATOP_SHIFT=46, 142735153367SBarry Smith MATOP_DIAGONAL_SET=47, 14289d39f69dSJed Brown MATOP_ZERO_ROWS_COLUMNS=48, 14299d39f69dSJed Brown MATOP_SET_RANDOM=49, 1430d519adbfSMatthew Knepley MATOP_GET_ROW_IJ=50, 1431d519adbfSMatthew Knepley MATOP_RESTORE_ROW_IJ=51, 1432d519adbfSMatthew Knepley MATOP_GET_COLUMN_IJ=52, 1433d519adbfSMatthew Knepley MATOP_RESTORE_COLUMN_IJ=53, 1434d519adbfSMatthew Knepley MATOP_FDCOLORING_CREATE=54, 1435d519adbfSMatthew Knepley MATOP_COLORING_PATCH=55, 1436d519adbfSMatthew Knepley MATOP_SET_UNFACTORED=56, 1437d519adbfSMatthew Knepley MATOP_PERMUTE=57, 1438d519adbfSMatthew Knepley MATOP_SET_VALUES_BLOCKED=58, 14397dae84e0SHong Zhang MATOP_CREATE_SUBMATRIX=59, 1440d519adbfSMatthew Knepley MATOP_DESTROY=60, 1441d519adbfSMatthew Knepley MATOP_VIEW=61, 1442d519adbfSMatthew Knepley MATOP_CONVERT_FROM=62, 14437bab7c10SHong Zhang MATOP_MATMAT_MULT=63, 14447bab7c10SHong Zhang MATOP_MATMAT_MULT_SYMBOLIC=64, 14457bab7c10SHong Zhang MATOP_MATMAT_MULT_NUMERIC=65, 1446d519adbfSMatthew Knepley MATOP_SET_LOCAL_TO_GLOBAL_MAP=66, 1447d519adbfSMatthew Knepley MATOP_SET_VALUES_LOCAL=67, 1448d519adbfSMatthew Knepley MATOP_ZERO_ROWS_LOCAL=68, 1449d519adbfSMatthew Knepley MATOP_GET_ROW_MAX_ABS=69, 1450d519adbfSMatthew Knepley MATOP_GET_ROW_MIN_ABS=70, 1451d519adbfSMatthew Knepley MATOP_CONVERT=71, 1452d519adbfSMatthew Knepley MATOP_SET_COLORING=72, 14539d39f69dSJed Brown /* MATOP_PLACEHOLDER_73=73, */ 1454d519adbfSMatthew Knepley MATOP_SET_VALUES_ADIFOR=74, 1455d519adbfSMatthew Knepley MATOP_FD_COLORING_APPLY=75, 1456d519adbfSMatthew Knepley MATOP_SET_FROM_OPTIONS=76, 1457bda6c4cbSBarry Smith MATOP_MULT_CONSTRAINED=77, 1458bda6c4cbSBarry Smith MATOP_MULT_TRANSPOSE_CONSTRAIN=78, 14599d39f69dSJed Brown MATOP_FIND_ZERO_DIAGONALS=79, 1460d519adbfSMatthew Knepley MATOP_MULT_MULTIPLE=80, 1461d519adbfSMatthew Knepley MATOP_SOLVE_MULTIPLE=81, 1462d519adbfSMatthew Knepley MATOP_GET_INERTIA=82, 1463d519adbfSMatthew Knepley MATOP_LOAD=83, 1464d519adbfSMatthew Knepley MATOP_IS_SYMMETRIC=84, 1465d519adbfSMatthew Knepley MATOP_IS_HERMITIAN=85, 1466d519adbfSMatthew Knepley MATOP_IS_STRUCTURALLY_SYMMETRIC=86, 1467820469bcSHong Zhang MATOP_SET_VALUES_BLOCKEDLOCAL=87, 1468*b41ce5d5SBarry Smith MATOP_CREATE_VECS=88, 1469d519adbfSMatthew Knepley MATOP_MAT_MULT=89, 1470d519adbfSMatthew Knepley MATOP_MAT_MULT_SYMBOLIC=90, 1471d519adbfSMatthew Knepley MATOP_MAT_MULT_NUMERIC=91, 1472d519adbfSMatthew Knepley MATOP_PTAP=92, 1473d519adbfSMatthew Knepley MATOP_PTAP_SYMBOLIC=93, 1474d519adbfSMatthew Knepley MATOP_PTAP_NUMERIC=94, 1475bda6c4cbSBarry Smith MATOP_MAT_TRANSPOSE_MULT=95, 1476bda6c4cbSBarry Smith MATOP_MAT_TRANSPOSE_MULT_SYMBO=96, 1477bda6c4cbSBarry Smith MATOP_MAT_TRANSPOSE_MULT_NUMER=97, 14789d39f69dSJed Brown /* MATOP_PLACEHOLDER_98=98, */ 14799d39f69dSJed Brown /* MATOP_PLACEHOLDER_99=99, */ 14809d39f69dSJed Brown /* MATOP_PLACEHOLDER_100=100, */ 14819d39f69dSJed Brown /* MATOP_PLACEHOLDER_101=101, */ 1482d519adbfSMatthew Knepley MATOP_CONJUGATE=102, 14839d39f69dSJed Brown /* MATOP_PLACEHOLDER_103=103, */ 1484d519adbfSMatthew Knepley MATOP_SET_VALUES_ROW=104, 1485d519adbfSMatthew Knepley MATOP_REAL_PART=105, 1486bda6c4cbSBarry Smith MATOP_IMAGINARY_PART=106, 1487bda6c4cbSBarry Smith MATOP_GET_ROW_UPPER_TRIANGULAR=107, 1488bda6c4cbSBarry Smith MATOP_RESTORE_ROW_UPPER_TRIANG=108, 1489bda6c4cbSBarry Smith MATOP_MAT_SOLVE=109, 149029eadf9eSStefano Zampini MATOP_MAT_SOLVE_TRANSPOSE=110, 1491d519adbfSMatthew Knepley MATOP_GET_ROW_MIN=111, 14920e8d04f7SBarry Smith MATOP_GET_COLUMN_VECTOR=112, 1493d519adbfSMatthew Knepley MATOP_MISSING_DIAGONAL=113, 14940e8d04f7SBarry Smith MATOP_GET_SEQ_NONZERO_STRUCTUR=114, 149589c6957cSBarry Smith MATOP_CREATE=115, 149689c6957cSBarry Smith MATOP_GET_GHOSTS=116, 14970e8d04f7SBarry Smith MATOP_GET_LOCAL_SUB_MATRIX=117, 14980e8d04f7SBarry Smith MATOP_RESTORE_LOCALSUB_MATRIX=118, 1499eeffb40dSHong Zhang MATOP_MULT_DIAGONAL_BLOCK=119, 15000e8d04f7SBarry Smith MATOP_HERMITIAN_TRANSPOSE=120, 15010e8d04f7SBarry Smith MATOP_MULT_HERMITIAN_TRANSPOSE=121, 15020e8d04f7SBarry Smith MATOP_MULT_HERMITIAN_TRANS_ADD=122, 15030e8d04f7SBarry Smith MATOP_GET_MULTI_PROC_BLOCK=123, 15049d39f69dSJed Brown MATOP_FIND_NONZERO_ROWS=124, 15050e8d04f7SBarry Smith MATOP_GET_COLUMN_NORMS=125, 15069d39f69dSJed Brown MATOP_INVERT_BLOCK_DIAGONAL=126, 15079d39f69dSJed Brown /* MATOP_PLACEHOLDER_127=127, */ 1508*b41ce5d5SBarry Smith MATOP_CREATE_SUB_MATRICES_MPI=128, 15092b8ad9a3SHong Zhang MATOP_SET_VALUES_BATCH=129, 15100e8d04f7SBarry Smith MATOP_TRANSPOSE_MAT_MULT=130, 15110e8d04f7SBarry Smith MATOP_TRANSPOSE_MAT_MULT_SYMBO=131, 1512e9b602ebSSatish Balay MATOP_TRANSPOSE_MAT_MULT_NUMER=132, 15130e8d04f7SBarry Smith MATOP_TRANSPOSE_COLORING_CREAT=133, 15140e8d04f7SBarry Smith MATOP_TRANS_COLORING_APPLY_SPT=134, 15150e8d04f7SBarry Smith MATOP_TRANS_COLORING_APPLY_DEN=135, 15160e8d04f7SBarry Smith MATOP_RART=136, 15170e8d04f7SBarry Smith MATOP_RART_SYMBOLIC=137, 15180e8d04f7SBarry Smith MATOP_RART_NUMERIC=138, 1519e09a3074SHong Zhang MATOP_SET_BLOCK_SIZES=139, 1520f9426fe0SMark Adams MATOP_AYPX=140, 15211919a2e2SJed Brown MATOP_RESIDUAL=141, 15229c8f2541SHong Zhang MATOP_FDCOLORING_SETUP=142, 15232d033e1fSHong Zhang MATOP_MPICONCATENATESEQ=144, 15242d033e1fSHong Zhang MATOP_DESTROYSUBMATRICES=145 1525fae171e0SBarry Smith } MatOperation; 1526014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatHasOperation(Mat,MatOperation,PetscBool *); 1527014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellSetOperation(Mat,MatOperation,void(*)(void)); 1528014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellGetOperation(Mat,MatOperation,void(**)(void)); 1529014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellSetContext(Mat,void*); 1530112a2221SBarry Smith 153190ace30eSBarry Smith /* 153290ace30eSBarry Smith Codes for matrices stored on disk. By default they are 153390ace30eSBarry Smith stored in a universal format. By changing the format with 15346a9046bcSBarry Smith PetscViewerPushFormat(viewer,PETSC_VIEWER_NATIVE); the matrices will 153590ace30eSBarry Smith be stored in a way natural for the matrix, for example dense matrices 153690ace30eSBarry Smith would be stored as dense. Matrices stored this way may only be 1537f4403165SShri Abhyankar read into matrices of the same type. 153890ace30eSBarry Smith */ 153990ace30eSBarry Smith #define MATRIX_BINARY_FORMAT_DENSE -1 154090ace30eSBarry Smith 1541014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetHashTableFactor(Mat,PetscReal); 1542014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatISGetLocalMat(Mat,Mat*); 15433b3b1effSJed Brown PETSC_EXTERN PetscErrorCode MatISRestoreLocalMat(Mat,Mat*); 1544014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatISSetLocalMat(Mat,Mat); 1545b7ce53b6SStefano Zampini PETSC_EXTERN PetscErrorCode MatISGetMPIXAIJ(Mat,MatReuse,Mat*); 15461f4e1ec7SBarry Smith 1547d9274352SBarry Smith /*S 1548d9274352SBarry Smith MatNullSpace - Object that removes a null space from a vector, i.e. 1549d9274352SBarry Smith orthogonalizes the vector to a subsapce 1550d9274352SBarry Smith 1551f7a9e4ceSBarry Smith Level: advanced 1552d9274352SBarry Smith 1553d9274352SBarry Smith Concepts: matrix; linear operator, null space 1554d9274352SBarry Smith 15556e1639daSBarry Smith Users manual sections: 15566e1639daSBarry Smith . sec_singular 15576e1639daSBarry Smith 1558d9274352SBarry Smith .seealso: MatNullSpaceCreate() 1559d9274352SBarry Smith S*/ 156074637425SBarry Smith typedef struct _p_MatNullSpace* MatNullSpace; 1561d9274352SBarry Smith 1562014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceCreate(MPI_Comm,PetscBool ,PetscInt,const Vec[],MatNullSpace*); 1563014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(MatNullSpace,Vec,void*),void*); 1564014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceDestroy(MatNullSpace*); 1565d0195637SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceRemove(MatNullSpace,Vec); 1566014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNullSpace(Mat, MatNullSpace *); 15675fa7ec2dSBarry Smith PETSC_EXTERN PetscErrorCode MatGetTransposeNullSpace(Mat, MatNullSpace *); 15685fa7ec2dSBarry Smith PETSC_EXTERN PetscErrorCode MatSetTransposeNullSpace(Mat,MatNullSpace); 1569014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNullSpace(Mat,MatNullSpace); 1570014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNearNullSpace(Mat,MatNullSpace); 1571014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNearNullSpace(Mat,MatNullSpace*); 1572014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceTest(MatNullSpace,Mat,PetscBool *); 1573014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceView(MatNullSpace,PetscViewer); 1574014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceGetVecs(MatNullSpace,PetscBool*,PetscInt*,const Vec**); 1575014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceCreateRigidBody(Vec,MatNullSpace*); 157674637425SBarry Smith 1577014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatReorderingSeqSBAIJ(Mat,IS); 1578014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetHashTableFactor(Mat,PetscReal); 1579014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetColumnIndices(Mat,PetscInt *); 1580014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJInvertBlockDiagonal(Mat); 15813f1d51d7SBarry Smith 1582014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMAIJ(Mat,PetscInt,Mat*); 1583014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMAIJRedimension(Mat,PetscInt,Mat*); 1584014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMAIJGetAIJ(Mat,Mat*); 1585c4f061fbSSatish Balay 1586014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatComputeExplicitOperator(Mat,Mat*); 1587b0a32e0cSBarry Smith 1588014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalScaleLocal(Mat,Vec); 158904f1ad80SBarry Smith 1590014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*); 1591014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetBase(Mat,Vec,Vec); 1592014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*); 1593014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*)); 1594014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec)); 1595014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt); 1596014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDResetHHistory(Mat); 1597014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctionError(Mat,PetscReal); 1598014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetPeriod(Mat,PetscInt); 1599014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDGetH(Mat,PetscScalar *); 1600014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetOptionsPrefix(Mat,const char[]); 1601014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*); 1602014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*); 1603e884886fSBarry Smith 16046370053bSBarry Smith /*S 16056370053bSBarry Smith MatMFFD - A data structured used to manage the computation of the h differencing parameter for matrix-free 16066370053bSBarry Smith Jacobian vector products 1607e884886fSBarry Smith 16086370053bSBarry Smith Notes: MATMFFD is a specific MatType which uses the MatMFFD data structure 16096370053bSBarry Smith 16106370053bSBarry Smith MatMFFD*() methods actually take the Mat as their first argument. Not a MatMFFD data structure 16116370053bSBarry Smith 16126370053bSBarry Smith Level: developer 16136370053bSBarry Smith 16146370053bSBarry Smith .seealso: MATMFFD, MatCreateMFFD(), MatMFFDSetFuction(), MatMFFDSetType(), MatMFFDRegister() 16156370053bSBarry Smith S*/ 1616e884886fSBarry Smith typedef struct _p_MatMFFD* MatMFFD; 1617e884886fSBarry Smith 161876bdecfbSBarry Smith /*J 1619e884886fSBarry Smith MatMFFDType - algorithm used to compute the h used in computing matrix-vector products via differencing of the function 1620e884886fSBarry Smith 1621e884886fSBarry Smith Level: beginner 1622e884886fSBarry Smith 1623e884886fSBarry Smith .seealso: MatMFFDSetType(), MatMFFDRegister() 162476bdecfbSBarry Smith J*/ 162519fd82e9SBarry Smith typedef const char* MatMFFDType; 1626e884886fSBarry Smith #define MATMFFD_DS "ds" 1627e884886fSBarry Smith #define MATMFFD_WP "wp" 1628e884886fSBarry Smith 162919fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDSetType(Mat,MatMFFDType); 1630bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDRegister(const char[],PetscErrorCode (*)(MatMFFD)); 1631e884886fSBarry Smith 1632014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDDSSetUmin(Mat,PetscReal); 1633014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDWPSetComputeNormU(Mat,PetscBool ); 1634e884886fSBarry Smith 163561ab5df0SBarry Smith PETSC_EXTERN PetscErrorCode MatFDColoringSetType(MatFDColoring,MatMFFDType); 163661ab5df0SBarry Smith 1637014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *); 1638014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *); 16397dbadf16SMatthew Knepley 164097969023SHong Zhang /* 164197969023SHong Zhang PETSc interface to MUMPS 164297969023SHong Zhang */ 164397969023SHong Zhang #ifdef PETSC_HAVE_MUMPS 1644014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMumpsSetIcntl(Mat,PetscInt,PetscInt); 1645bc6112feSHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetIcntl(Mat,PetscInt,PetscInt*); 1646b9e9ea26SSatish Balay PETSC_EXTERN PetscErrorCode MatMumpsSetCntl(Mat,PetscInt,PetscReal); 1647bc6112feSHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetCntl(Mat,PetscInt,PetscReal*); 1648bc6112feSHong Zhang 1649ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetInfo(Mat,PetscInt,PetscInt*); 1650ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetInfog(Mat,PetscInt,PetscInt*); 1651ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetRinfo(Mat,PetscInt,PetscReal*); 1652ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetRinfog(Mat,PetscInt,PetscReal*); 165397969023SHong Zhang #endif 165423a5497aSJed Brown 1655d954db57SHong Zhang /* 1656d615f992SSatish Balay PETSc interface to Mkl_Pardiso 1657b500a6b7SJose David Bermeol */ 1658b500a6b7SJose David Bermeol #ifdef PETSC_HAVE_MKL_PARDISO 1659d615f992SSatish Balay PETSC_EXTERN PetscErrorCode MatMkl_PardisoSetCntl(Mat,PetscInt,PetscInt); 1660b500a6b7SJose David Bermeol #endif 1661b500a6b7SJose David Bermeol 1662b500a6b7SJose David Bermeol /* 1663d305a81bSVasiliy Kozyrev PETSc interface to Mkl_CPardiso 1664d305a81bSVasiliy Kozyrev */ 1665d305a81bSVasiliy Kozyrev #ifdef PETSC_HAVE_MKL_CPARDISO 1666d305a81bSVasiliy Kozyrev PETSC_EXTERN PetscErrorCode MatMkl_CPardisoSetCntl(Mat,PetscInt,PetscInt); 1667d305a81bSVasiliy Kozyrev #endif 1668d305a81bSVasiliy Kozyrev 1669d305a81bSVasiliy Kozyrev /* 1670d954db57SHong Zhang PETSc interface to SUPERLU 1671d954db57SHong Zhang */ 1672d954db57SHong Zhang #ifdef PETSC_HAVE_SUPERLU 1673014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSuperluSetILUDropTol(Mat,PetscReal); 1674d954db57SHong Zhang #endif 1675d954db57SHong Zhang 1676fb785984SHong Zhang /* 1677fb785984SHong Zhang PETSc interface to SUPERLU_DIST 1678fb785984SHong Zhang */ 1679fb785984SHong Zhang #ifdef PETSC_HAVE_SUPERLU_DIST 1680fb785984SHong Zhang PETSC_EXTERN PetscErrorCode MatSuperluDistGetDiagU(Mat,PetscScalar*); 1681fb785984SHong Zhang #endif 1682fb785984SHong Zhang 1683575704cbSPieter Ghysels /* 1684575704cbSPieter Ghysels PETSc interface to STRUMPACK 1685575704cbSPieter Ghysels */ 1686575704cbSPieter Ghysels #ifdef PETSC_HAVE_STRUMPACK 1687575704cbSPieter Ghysels PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSRelCompTol(Mat,PetscReal); 1688575704cbSPieter Ghysels PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSMinSize(Mat,PetscInt); 1689575704cbSPieter Ghysels PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetColPerm(Mat,PetscBool); 1690575704cbSPieter Ghysels #endif 1691575704cbSPieter Ghysels 1692575704cbSPieter Ghysels 1693aa372e3fSPaul Mullowney #ifdef PETSC_HAVE_CUDA 16943f9c0db1SPaul Mullowney /*E 1695e057df02SPaul Mullowney MatCUSPARSEStorageFormat - indicates the storage format for CUSPARSE (GPU) 16962692e278SPaul Mullowney matrices. 1697e057df02SPaul Mullowney 1698e057df02SPaul Mullowney Not Collective 1699e057df02SPaul Mullowney 17008468deeeSKarl Rupp + MAT_CUSPARSE_CSR - Compressed Sparse Row 17012692e278SPaul Mullowney . MAT_CUSPARSE_ELL - Ellpack (requires CUDA 4.2 or later). 17022692e278SPaul Mullowney - MAT_CUSPARSE_HYB - Hybrid, a combination of Ellpack and Coordinate format (requires CUDA 4.2 or later). 1703e057df02SPaul Mullowney 1704e057df02SPaul Mullowney Level: intermediate 1705e057df02SPaul Mullowney 1706af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 1707e057df02SPaul Mullowney 1708e057df02SPaul Mullowney .seealso: MatCUSPARSESetFormat(), MatCUSPARSEFormatOperation 1709e057df02SPaul Mullowney E*/ 1710e057df02SPaul Mullowney 17113f9c0db1SPaul Mullowney typedef enum {MAT_CUSPARSE_CSR, MAT_CUSPARSE_ELL, MAT_CUSPARSE_HYB} MatCUSPARSEStorageFormat; 1712e057df02SPaul Mullowney 1713e057df02SPaul Mullowney /* these will be strings associated with enumerated type defined above */ 1714e057df02SPaul Mullowney PETSC_EXTERN const char *const MatCUSPARSEStorageFormats[]; 1715e057df02SPaul Mullowney 17163f9c0db1SPaul Mullowney /*E 1717e057df02SPaul Mullowney MatCUSPARSEFormatOperation - indicates the operation of CUSPARSE (GPU) 17182692e278SPaul Mullowney matrices whose operation should use a particular storage format. 1719e057df02SPaul Mullowney 1720e057df02SPaul Mullowney Not Collective 1721e057df02SPaul Mullowney 17228468deeeSKarl Rupp + MAT_CUSPARSE_MULT_DIAG - sets the storage format for the diagonal matrix in the parallel MatMult 17238468deeeSKarl Rupp . MAT_CUSPARSE_MULT_OFFDIAG - sets the storage format for the offdiagonal matrix in the parallel MatMult 17248468deeeSKarl Rupp . MAT_CUSPARSE_MULT - sets the storage format for the entire matrix in the serial (single GPU) MatMult 17258468deeeSKarl Rupp - MAT_CUSPARSE_ALL - sets the storage format for all CUSPARSE (GPU) matrices 1726e057df02SPaul Mullowney 1727e057df02SPaul Mullowney Level: intermediate 1728e057df02SPaul Mullowney 1729e057df02SPaul Mullowney .seealso: MatCUSPARSESetFormat(), MatCUSPARSEStorageFormat 1730e057df02SPaul Mullowney E*/ 173136d62e41SPaul Mullowney typedef enum {MAT_CUSPARSE_MULT_DIAG, MAT_CUSPARSE_MULT_OFFDIAG, MAT_CUSPARSE_MULT, MAT_CUSPARSE_ALL} MatCUSPARSEFormatOperation; 1732e057df02SPaul Mullowney 173310e9db80SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 173410e9db80SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCreateAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 1735e057df02SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCUSPARSESetFormat(Mat,MatCUSPARSEFormatOperation,MatCUSPARSEStorageFormat); 17369ae82921SPaul Mullowney #endif 17379ae82921SPaul Mullowney 173890c53902SBarry Smith #if defined(PETSC_HAVE_CUSP) 1739014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCUSP(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 1740014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateAIJCUSP(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 1741e057df02SPaul Mullowney 17423f9c0db1SPaul Mullowney /*E 1743e057df02SPaul Mullowney MatCUSPStorageFormat - indicates the storage format for CUSP (GPU) 174436d62e41SPaul Mullowney matrices. 1745e057df02SPaul Mullowney 1746e057df02SPaul Mullowney Not Collective 1747e057df02SPaul Mullowney 17488468deeeSKarl Rupp + MAT_CUSP_CSR - Compressed Sparse Row 17498468deeeSKarl Rupp . MAT_CUSP_DIA - Diagonal 17508468deeeSKarl Rupp - MAT_CUSP_ELL - Ellpack 1751e057df02SPaul Mullowney 1752e057df02SPaul Mullowney Level: intermediate 1753e057df02SPaul Mullowney 1754af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 1755e057df02SPaul Mullowney 1756e057df02SPaul Mullowney .seealso: MatCUSPSetFormat(), MatCUSPFormatOperation 1757e057df02SPaul Mullowney E*/ 17583f9c0db1SPaul Mullowney typedef enum {MAT_CUSP_CSR, MAT_CUSP_DIA, MAT_CUSP_ELL} MatCUSPStorageFormat; 1759e057df02SPaul Mullowney 1760e057df02SPaul Mullowney /* these will be strings associated with enumerated type defined above */ 1761e057df02SPaul Mullowney PETSC_EXTERN const char *const MatCUSPStorageFormats[]; 1762e057df02SPaul Mullowney 17633f9c0db1SPaul Mullowney /*E 1764e057df02SPaul Mullowney MatCUSPFormatOperation - indicates the operation of CUSP (GPU) 176536d62e41SPaul Mullowney matrices whose operation should use a particular storage format. 1766e057df02SPaul Mullowney 1767e057df02SPaul Mullowney Not Collective 1768e057df02SPaul Mullowney 17698468deeeSKarl Rupp + MAT_CUSP_MULT_DIAG - sets the storage format for the diagonal matrix in the parallel MatMult 17708468deeeSKarl Rupp . MAT_CUSP_MULT_OFFDIAG - sets the storage format for the offdiagonal matrix in the parallel MatMult 17718468deeeSKarl Rupp . MAT_CUSP_MULT - sets the storage format for the entire matrix in the serial (single GPU) MatMult 17728468deeeSKarl Rupp - MAT_CUSP_ALL - sets the storage format for all CUSP (GPU) matrices 1773e057df02SPaul Mullowney 1774e057df02SPaul Mullowney Level: intermediate 1775e057df02SPaul Mullowney 1776af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 1777e057df02SPaul Mullowney 1778e057df02SPaul Mullowney .seealso: MatCUSPSetFormat(), MatCUSPStorageFormat 1779e057df02SPaul Mullowney E*/ 178036d62e41SPaul Mullowney typedef enum {MAT_CUSP_MULT_DIAG, MAT_CUSP_MULT_OFFDIAG, MAT_CUSP_MULT, MAT_CUSP_ALL} MatCUSPFormatOperation; 1781e057df02SPaul Mullowney 1782e057df02SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCUSPSetFormat(Mat,MatCUSPFormatOperation,MatCUSPStorageFormat); 1783e057df02SPaul Mullowney #endif 178490c53902SBarry Smith 1785d67ff14aSKarl Rupp #if defined(PETSC_HAVE_VIENNACL) 1786d67ff14aSKarl Rupp PETSC_EXTERN PetscErrorCode MatCreateSeqAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 1787d67ff14aSKarl Rupp PETSC_EXTERN PetscErrorCode MatCreateAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 1788d67ff14aSKarl Rupp #endif 1789d67ff14aSKarl Rupp 179054efbe56SHong Zhang /* 179154efbe56SHong Zhang PETSc interface to FFTW 179254efbe56SHong Zhang */ 179354efbe56SHong Zhang #if defined(PETSC_HAVE_FFTW) 1794014dd563SJed Brown PETSC_EXTERN PetscErrorCode VecScatterPetscToFFTW(Mat,Vec,Vec); 1795014dd563SJed Brown PETSC_EXTERN PetscErrorCode VecScatterFFTWToPetsc(Mat,Vec,Vec); 17962a7a6963SBarry Smith PETSC_EXTERN PetscErrorCode MatCreateVecsFFTW(Mat,Vec*,Vec*,Vec*); 179754efbe56SHong Zhang #endif 179854efbe56SHong Zhang 1799014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateNest(MPI_Comm,PetscInt,const IS[],PetscInt,const IS[],const Mat[],Mat*); 1800014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSize(Mat,PetscInt*,PetscInt*); 1801014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetISs(Mat,IS[],IS[]); 1802014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetLocalISs(Mat,IS[],IS[]); 1803014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSubMats(Mat,PetscInt*,PetscInt*,Mat***); 1804014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSubMat(Mat,PetscInt,PetscInt,Mat*); 180519fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatNestSetVecType(Mat,VecType); 1806014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestSetSubMats(Mat,PetscInt,const IS[],PetscInt,const IS[],const Mat[]); 1807014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestSetSubMat(Mat,PetscInt,PetscInt,Mat); 1808d8588912SDave May 18094325cce7SMatthew G Knepley PETSC_EXTERN PetscErrorCode MatChop(Mat,PetscReal); 1810e0a58c10SMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatComputeBandwidth(Mat,PetscReal,PetscInt*); 18114325cce7SMatthew G Knepley 1812b541e6a4SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSubdomainsCreateCoalesce(Mat,PetscInt,PetscInt*,IS**); 181353dd7562SDmitry Karpeev 1814c094ef40SMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatPreallocatorPreallocate(Mat,PetscBool,Mat); 1815c094ef40SMatthew G. Knepley 1816539c167fSBarry Smith PETSC_INTERN PetscErrorCode MatHeaderMerge(Mat,Mat*); 1817539c167fSBarry Smith PETSC_EXTERN PetscErrorCode MatHeaderReplace(Mat,Mat*); 1818539c167fSBarry Smith 181923a5497aSJed Brown #endif 1820