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*); 2078d7a6e47SBarry Smith 208014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 209014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 210014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat*); 211d21a29f3SJed Brown 212014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAdj(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscInt[],Mat*); 213014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 214d21a29f3SJed Brown 215014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 216014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPISBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *); 21738f409ebSLisandro Dalcin PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 218014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 2194cce697cSJed Brown PETSC_EXTERN PetscErrorCode MatXAIJSetPreallocation(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscInt[],const PetscInt[]); 220dfb205c3SBarry Smith 221014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateShell(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,void *,Mat*); 222014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateNormal(Mat,Mat*); 223c5e4d11fSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatCreateNormalHermitian(Mat,Mat*); 224986c4d72SJose E. Roman PETSC_EXTERN PetscErrorCode MatCreateLRC(Mat,Mat,Vec,Mat,Mat*); 225267ca84cSJose E. Roman PETSC_EXTERN PetscErrorCode MatLRCGetMats(Mat,Mat*,Mat*,Vec*,Mat*); 226c5e4d11fSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatCreateIS(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,ISLocalToGlobalMapping,ISLocalToGlobalMapping,Mat*); 227014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 228014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 229c0cdd4a1SDahai Guo 230014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateScatter(MPI_Comm,VecScatter,Mat*); 231014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScatterSetVecScatter(Mat,VecScatter); 232014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScatterGetVecScatter(Mat,VecScatter*); 233014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateBlockMat(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,Mat*); 234014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeAddMat(Mat,Mat); 235014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeMerge(Mat); 236014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateComposite(MPI_Comm,PetscInt,const Mat*,Mat*); 2376d7c1e57SBarry Smith typedef enum {MAT_COMPOSITE_ADDITIVE,MAT_COMPOSITE_MULTIPLICATIVE} MatCompositeType; 238014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeSetType(Mat,MatCompositeType); 2396d7c1e57SBarry Smith 24019fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatCreateFFT(MPI_Comm,PetscInt,const PetscInt[],MatType,Mat*); 241014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqCUFFT(MPI_Comm,PetscInt,const PetscInt[],Mat*); 242dedccee8SHong Zhang 243014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateTranspose(Mat,Mat*); 2448060fb66Sstefano_zampini PETSC_EXTERN PetscErrorCode MatTransposeGetMat(Mat,Mat*); 245d0de2241SAndrew Spott PETSC_EXTERN PetscErrorCode MatCreateHermitianTranspose(Mat,Mat*); 24654e05e6cSHong Zhang PETSC_EXTERN PetscErrorCode MatCreateSubMatrixVirtual(Mat,IS,IS,Mat*); 24754e05e6cSHong Zhang PETSC_EXTERN PetscErrorCode MatSubMatrixVirtualUpdate(Mat,Mat,IS,IS); 248014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateLocalRef(Mat,IS,IS,Mat*); 2491472f72bSBarry Smith 250978814f1SStefano Zampini #if defined(PETSC_HAVE_HYPRE) 251d975228cSstefano_zampini PETSC_EXTERN PetscErrorCode MatHYPRESetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 252978814f1SStefano Zampini #endif 253978814f1SStefano Zampini 254014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPythonSetType(Mat,const char[]); 2551d6018f0SLisandro Dalcin 256014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetUp(Mat); 257014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroy(Mat*); 258e56f5c9eSBarry Smith PETSC_EXTERN PetscErrorCode MatGetNonzeroState(Mat,PetscObjectState*); 25921c89e3eSBarry Smith 260014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatConjugate(Mat); 261014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRealPart(Mat); 262014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatImaginaryPart(Mat); 263014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetDiagonalBlock(Mat,Mat*); 264014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetTrace(Mat,PetscScalar*); 265713ccfa9SJed Brown PETSC_EXTERN PetscErrorCode MatInvertBlockDiagonal(Mat,const PetscScalar **); 26699cafbc1SBarry Smith 2678ed539a5SBarry Smith /* ------------------------------------------------------------*/ 268014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 269014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 270014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesRow(Mat,PetscInt,const PetscScalar[]); 271014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesRowLocal(Mat,PetscInt,const PetscScalar[]); 272014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBatch(Mat,PetscInt,PetscInt,PetscInt[],const PetscScalar[]); 27373a71a0fSBarry Smith PETSC_EXTERN PetscErrorCode MatSetRandom(Mat,PetscRandom); 27484cb2905SBarry Smith 2752ef4de8bSBarry Smith /*S 2762ef4de8bSBarry Smith MatStencil - Data structure (C struct) for storing information about a single row or 27762ef0227SBarry Smith column of a matrix as indexed on an associated grid. These are arguments to MatSetStencil() and MatSetBlockStencil() 27862ef0227SBarry Smith 27962ef0227SBarry 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). 28062ef0227SBarry 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. 28162ef0227SBarry Smith 28262ef0227SBarry Smith For stencil access to vectors see DMDAVecGetArray(), DMDAVecGetArrayF90(). 283be479b30SJed Brown 284be479b30SJed Brown Fortran usage is different, see MatSetValuesStencil() for details. 2852ef4de8bSBarry Smith 2862ef4de8bSBarry Smith Level: beginner 2872ef4de8bSBarry Smith 2882ef4de8bSBarry Smith Concepts: matrix; linear operator 2892ef4de8bSBarry Smith 29062ef0227SBarry Smith .seealso: MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockedStencil(), DMDAVecGetArray(), DMDAVecGetArrayF90() 2912ef4de8bSBarry Smith S*/ 292435da068SBarry Smith typedef struct { 293c1ac3661SBarry Smith PetscInt k,j,i,c; 294435da068SBarry Smith } MatStencil; 2952ef4de8bSBarry Smith 296014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode); 297014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode); 298014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt); 299435da068SBarry Smith 300d91e6319SBarry Smith /*E 301d91e6319SBarry Smith MatAssemblyType - Indicates if the matrix is now to be used, or if you plan 30262ef0227SBarry Smith to continue to add or insert values to it 303d91e6319SBarry Smith 304d91e6319SBarry Smith Level: beginner 305d91e6319SBarry Smith 306d91e6319SBarry Smith .seealso: MatAssemblyBegin(), MatAssemblyEnd() 307d91e6319SBarry Smith E*/ 3086d4a8577SBarry Smith typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType; 309014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssemblyBegin(Mat,MatAssemblyType); 310014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssemblyEnd(Mat,MatAssemblyType); 311014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssembled(Mat,PetscBool *); 3124f9c727eSBarry Smith 3131d73ed98SBarry Smith 31430de9b25SBarry Smith 315d91e6319SBarry Smith /*E 316d91e6319SBarry Smith MatOption - Options that may be set for a matrix and its behavior or storage 317d91e6319SBarry Smith 318d91e6319SBarry Smith Level: beginner 319d91e6319SBarry Smith 320af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 321d91e6319SBarry Smith 322382164b0SBarry Smith Developer Notes: Entries that are negative need not be called collectively by all processes. 323382164b0SBarry Smith 324d91e6319SBarry Smith .seealso: MatSetOption() 325d91e6319SBarry Smith E*/ 326c5e4d11fSDmitry Karpeev typedef enum {MAT_OPTION_MIN = -3, 327c5e4d11fSDmitry Karpeev MAT_UNUSED_NONZERO_LOCATION_ERR = -2, 32892d4d306SBarry Smith MAT_ROW_ORIENTED = -1, 329382164b0SBarry Smith MAT_SYMMETRIC = 1, 330382164b0SBarry Smith MAT_STRUCTURALLY_SYMMETRIC = 2, 331382164b0SBarry Smith MAT_NEW_DIAGONALS = 3, 332382164b0SBarry Smith MAT_IGNORE_OFF_PROC_ENTRIES = 4, 333382164b0SBarry Smith MAT_USE_HASH_TABLE = 5, 334382164b0SBarry Smith MAT_KEEP_NONZERO_PATTERN = 6, 335382164b0SBarry Smith MAT_IGNORE_ZERO_ENTRIES = 7, 336382164b0SBarry Smith MAT_USE_INODES = 8, 337382164b0SBarry Smith MAT_HERMITIAN = 9, 338382164b0SBarry Smith MAT_SYMMETRY_ETERNAL = 10, 339c5e4d11fSDmitry Karpeev MAT_NEW_NONZERO_LOCATION_ERR = 11, 340382164b0SBarry Smith MAT_IGNORE_LOWER_TRIANGULAR = 12, 341382164b0SBarry Smith MAT_ERROR_LOWER_TRIANGULAR = 13, 342382164b0SBarry Smith MAT_GETROW_UPPERTRIANGULAR = 14, 343382164b0SBarry Smith MAT_SPD = 15, 34492d4d306SBarry Smith MAT_NO_OFF_PROC_ZERO_ROWS = 16, 34592d4d306SBarry Smith MAT_NO_OFF_PROC_ENTRIES = 17, 34692d4d306SBarry Smith MAT_NEW_NONZERO_LOCATIONS = 18, 347c5e4d11fSDmitry Karpeev MAT_NEW_NONZERO_ALLOCATION_ERR = 19, 348c5e4d11fSDmitry Karpeev MAT_SUBSET_OFF_PROC_ENTRIES = 20, 349c10200c1SHong Zhang MAT_SUBMAT_SINGLEIS = 21, 350c10200c1SHong Zhang MAT_OPTION_MAX = 22} MatOption; 351e057df02SPaul Mullowney 352014dd563SJed Brown PETSC_EXTERN const char *MatOptions[]; 353014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetOption(Mat,MatOption,PetscBool); 3547d68702bSBarry Smith PETSC_EXTERN PetscErrorCode MatGetOption(Mat,MatOption,PetscBool*); 35519fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatGetType(Mat,MatType*); 35684cb2905SBarry Smith 357014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]); 358014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 359014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 360014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowUpperTriangular(Mat); 361014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRowUpperTriangular(Mat); 362014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnVector(Mat,Vec,PetscInt); 3638c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJGetArray(Mat,PetscScalar *[]); 3648c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJRestoreArray(Mat,PetscScalar *[]); 36521e72a00SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJGetMaxRowNonzeros(Mat,PetscInt*); 366bd04181cSBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJSetValuesLocalFast(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 3678397e458SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqSBAIJGetArray(Mat,PetscScalar *[]); 3688397e458SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqSBAIJRestoreArray(Mat,PetscScalar *[]); 3698c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatDenseGetArray(Mat,PetscScalar *[]); 3708c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatDenseRestoreArray(Mat,PetscScalar *[]); 371014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBlockSize(Mat,PetscInt *); 372014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSize(Mat,PetscInt); 373014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBlockSizes(Mat,PetscInt *,PetscInt *); 374014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSizes(Mat,PetscInt,PetscInt); 37533d57670SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSizesFromMats(Mat,Mat,Mat); 3761620fd73SBarry Smith 377014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMult(Mat,Vec,Vec); 378014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultDiagonalBlock(Mat,Vec,Vec); 379014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultAdd(Mat,Vec,Vec,Vec); 380014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTranspose(Mat,Vec,Vec); 381014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultHermitianTranspose(Mat,Vec,Vec); 382014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsTranspose(Mat,Mat,PetscReal,PetscBool *); 383014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitianTranspose(Mat,Mat,PetscReal,PetscBool *); 384014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeAdd(Mat,Vec,Vec,Vec); 385014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultHermitianTransposeAdd(Mat,Vec,Vec,Vec); 386014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultConstrained(Mat,Vec,Vec); 387014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeConstrained(Mat,Vec,Vec); 388014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatSolve(Mat,Mat,Mat); 389f9426fe0SMark Adams PETSC_EXTERN PetscErrorCode MatResidual(Mat,Vec,Vec,Vec); 390f5edf698SHong Zhang 391d91e6319SBarry Smith /*E 392d91e6319SBarry Smith MatDuplicateOption - Indicates if a duplicated sparse matrix should have 393d91e6319SBarry Smith its numerical values copied over or just its nonzero structure. 394d91e6319SBarry Smith 395d91e6319SBarry Smith Level: beginner 396d91e6319SBarry Smith 397af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 398d91e6319SBarry Smith 39970dcbbb9SBarry Smith $ MAT_SHARE_NONZERO_PATTERN - the i and j arrays in the new matrix will be shared with the original matrix 40070dcbbb9SBarry Smith $ this also triggers the MAT_DO_NOT_COPY_VALUES option. This is used when you 40170dcbbb9SBarry Smith $ have several matrices with the same nonzero pattern. 40270dcbbb9SBarry Smith 403d91e6319SBarry Smith .seealso: MatDuplicate() 404d91e6319SBarry Smith E*/ 40570dcbbb9SBarry Smith typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES,MAT_SHARE_NONZERO_PATTERN} MatDuplicateOption; 4062e8a6d31SBarry Smith 40719fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatConvert(Mat,MatType,MatReuse,Mat*); 408014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDuplicate(Mat,MatDuplicateOption,Mat*); 40994a9d846SBarry Smith 410cb5b572fSBarry Smith 411014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCopy(Mat,Mat,MatStructure); 412014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatView(Mat,PetscViewer); 413014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsSymmetric(Mat,PetscReal,PetscBool *); 414014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsStructurallySymmetric(Mat,PetscBool *); 415014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitian(Mat,PetscReal,PetscBool *); 416014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsSymmetricKnown(Mat,PetscBool *,PetscBool *); 417014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitianKnown(Mat,PetscBool *,PetscBool *); 418014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMissingDiagonal(Mat,PetscBool *,PetscInt *); 419014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLoad(Mat, PetscViewer); 4207b80b807SBarry Smith 4211a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool *); 4221a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool *); 4231a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool *); 4241a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool *); 425d4fbbf0eSBarry Smith 426d91e6319SBarry Smith /*S 427d91e6319SBarry Smith MatInfo - Context of matrix information, used with MatGetInfo() 428d91e6319SBarry Smith 429d91e6319SBarry Smith In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE 430d91e6319SBarry Smith 431d91e6319SBarry Smith Level: intermediate 432d91e6319SBarry Smith 433d91e6319SBarry Smith Concepts: matrix^nonzero information 434d91e6319SBarry Smith 435d9274352SBarry Smith .seealso: MatGetInfo(), MatInfoType 436d91e6319SBarry Smith S*/ 4374e220ebcSLois Curfman McInnes typedef struct { 438b0a32e0cSBarry Smith PetscLogDouble block_size; /* block size */ 439b0a32e0cSBarry Smith PetscLogDouble nz_allocated,nz_used,nz_unneeded; /* number of nonzeros */ 440b0a32e0cSBarry Smith PetscLogDouble memory; /* memory allocated */ 441b0a32e0cSBarry Smith PetscLogDouble assemblies; /* number of matrix assemblies called */ 442b0a32e0cSBarry Smith PetscLogDouble mallocs; /* number of mallocs during MatSetValues() */ 443b0a32e0cSBarry Smith PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */ 444b0a32e0cSBarry Smith PetscLogDouble factor_mallocs; /* number of mallocs during factorization */ 4454e220ebcSLois Curfman McInnes } MatInfo; 4464e220ebcSLois Curfman McInnes 447d9274352SBarry Smith /*E 448d9274352SBarry Smith MatInfoType - Indicates if you want information about the local part of the matrix, 449d9274352SBarry Smith the entire parallel matrix or the maximum over all the local parts. 450d9274352SBarry Smith 451d9274352SBarry Smith Level: beginner 452d9274352SBarry Smith 453af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 454d9274352SBarry Smith 455d9274352SBarry Smith .seealso: MatGetInfo(), MatInfo 456d9274352SBarry Smith E*/ 4577b80b807SBarry Smith typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType; 458014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetInfo(Mat,MatInfoType,MatInfo*); 459014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetDiagonal(Mat,Vec); 460014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMax(Mat,Vec,PetscInt[]); 461014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMin(Mat,Vec,PetscInt[]); 462014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMaxAbs(Mat,Vec,PetscInt[]); 463014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMinAbs(Mat,Vec,PetscInt[]); 464014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowSum(Mat,Vec); 465014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTranspose(Mat,MatReuse,Mat*); 466014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatHermitianTranspose(Mat,MatReuse,Mat*); 467014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPermute(Mat,IS,IS,Mat*); 468014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalScale(Mat,Vec,Vec); 469014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalSet(Mat,Vec,InsertMode); 470a52676ecSHong Zhang 471014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatEqual(Mat,Mat,PetscBool*); 472014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultEqual(Mat,Mat,PetscInt,PetscBool*); 473014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultAddEqual(Mat,Mat,PetscInt,PetscBool*); 474014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeEqual(Mat,Mat,PetscInt,PetscBool*); 475014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscBool*); 476a52676ecSHong Zhang PETSC_EXTERN PetscErrorCode MatMatMultEqual(Mat,Mat,Mat,PetscInt,PetscBool*); 477a52676ecSHong Zhang PETSC_EXTERN PetscErrorCode MatTransposeMatMultEqual(Mat,Mat,Mat,PetscInt,PetscBool*); 4787b80b807SBarry Smith 479014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNorm(Mat,NormType,PetscReal*); 480014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnNorms(Mat,NormType,PetscReal*); 481014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroEntries(Mat); 482014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 483014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsIS(Mat,IS,PetscScalar,Vec,Vec); 484014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsStencil(Mat,PetscInt,const MatStencil [],PetscScalar,Vec,Vec); 485014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsStencil(Mat,PetscInt,const MatStencil[],PetscScalar,Vec,Vec); 486014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumns(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 487014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsIS(Mat,IS,PetscScalar,Vec,Vec); 4887b80b807SBarry Smith 489014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSize(Mat,PetscInt*,PetscInt*); 490014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalSize(Mat,PetscInt*,PetscInt*); 491014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRange(Mat,PetscInt*,PetscInt*); 492014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRanges(Mat,const PetscInt**); 493014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRangeColumn(Mat,PetscInt*,PetscInt*); 494014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRangesColumn(Mat,const PetscInt**); 495566876eaSJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipIS(Mat,IS*,IS*); 4965ef9f2a5SBarry Smith 4977dae84e0SHong Zhang PETSC_EXTERN PetscErrorCode MatCreateSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]); 498cd2534ddSHong 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);} 4997dae84e0SHong Zhang PETSC_EXTERN PetscErrorCode MatCreateSubMatricesMPI(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]); 500cd2534ddSHong 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);} 501014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroyMatrices(PetscInt,Mat *[]); 502df750dc8SHong Zhang PETSC_EXTERN PetscErrorCode MatDestroySubMatrices(PetscInt,Mat *[]); 5037dae84e0SHong Zhang PETSC_EXTERN PetscErrorCode MatCreateSubMatrix(Mat,IS,IS,MatReuse,Mat *); 504cd2534ddSHong 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);} 505014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalSubMatrix(Mat,IS,IS,Mat*); 506014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreLocalSubMatrix(Mat,IS,IS,Mat*); 507014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSeqNonzeroStructure(Mat,Mat*); 508014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroySeqNonzeroStructure(Mat*); 5097b80b807SBarry Smith 510014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJ(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*); 511014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJSymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*); 512014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJNumeric(Mat,Mat); 513014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMat(Mat,MatReuse,Mat*); 514014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*); 515014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,Mat*); 516014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetGhosts(Mat, PetscInt *,const PetscInt *[]); 5178efafbd8SBarry Smith 518014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt); 519aa1e27eaSFande Kong PETSC_EXTERN PetscErrorCode MatIncreaseOverlapSplit(Mat mat,PetscInt n,IS is[],PetscInt ov); 520d2b2a242SBarry Smith PETSC_EXTERN PetscErrorCode MatMPIAIJSetUseScalableIncreaseOverlap(Mat,PetscBool); 5217b80b807SBarry Smith 522014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*); 523014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMultSymbolic(Mat,Mat,PetscReal,Mat*); 524014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMultNumeric(Mat,Mat,Mat); 52522440eb1SKris Buschelman 5267bab7c10SHong Zhang PETSC_EXTERN PetscErrorCode MatMatMatMult(Mat,Mat,Mat,MatReuse,PetscReal,Mat*); 5272df6c5c3SPatrick Farrell PETSC_EXTERN PetscErrorCode MatGalerkin(Mat,Mat,Mat,MatReuse,PetscReal,Mat*); 5287bab7c10SHong Zhang 529014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*); 530014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAPSymbolic(Mat,Mat,PetscReal,Mat*); 531014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAPNumeric(Mat,Mat,Mat); 532014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARt(Mat,Mat,MatReuse,PetscReal,Mat*); 533014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARtSymbolic(Mat,Mat,PetscReal,Mat*); 534014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARtNumeric(Mat,Mat,Mat); 53522440eb1SKris Buschelman 536014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeMatMult(Mat,Mat,MatReuse,PetscReal,Mat*); 537014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposetMatMultSymbolic(Mat,Mat,PetscReal,Mat*); 538014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposetMatMultNumeric(Mat,Mat,Mat); 539014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMult(Mat,Mat,MatReuse,PetscReal,Mat*); 540014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMultSymbolic(Mat,Mat,PetscReal,Mat*); 541014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMultNumeric(Mat,Mat,Mat); 542bc011b1eSHong Zhang 543014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAXPY(Mat,PetscScalar,Mat,MatStructure); 544014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAYPX(Mat,PetscScalar,Mat,MatStructure); 5451c741599SHong Zhang 546014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScale(Mat,PetscScalar); 547014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShift(Mat,PetscScalar); 5487b80b807SBarry Smith 549014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping); 550014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*); 551a93ff8c4SPeter Brune PETSC_EXTERN PetscErrorCode MatGetLayouts(Mat,PetscLayout*,PetscLayout*); 552014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 553014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsLocalIS(Mat,IS,PetscScalar,Vec,Vec); 554014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 555014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocalIS(Mat,IS,PetscScalar,Vec,Vec); 556014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 557014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 558052efed2SBarry Smith 559014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStashSetInitialSize(Mat,PetscInt,PetscInt); 560014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*); 56190f02eecSBarry Smith 562014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInterpolate(Mat,Vec,Vec); 563014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInterpolateAdd(Mat,Vec,Vec,Vec); 564014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestrict(Mat,Vec,Vec); 5652a7a6963SBarry Smith PETSC_EXTERN PetscErrorCode MatCreateVecs(Mat,Vec*,Vec*); 56684cc2a5aSBarry Smith PETSC_DEPRECATED("Use MatCreateVecs()") PETSC_STATIC_INLINE PetscErrorCode MatGetVecs(Mat mat,Vec *x,Vec *y) {return MatCreateVecs(mat,x,y);} 56753cd1579SHong Zhang PETSC_EXTERN PetscErrorCode MatCreateRedundantMatrix(Mat,PetscInt,MPI_Comm,MatReuse,Mat*); 568014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetMultiProcBlock(Mat,MPI_Comm,MatReuse,Mat*); 569014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFindZeroDiagonals(Mat,IS*); 5703a062f41SBarry Smith PETSC_EXTERN PetscErrorCode MatFindOffBlockDiagonalEntries(Mat,IS*); 5719c8f2541SHong Zhang PETSC_EXTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat(MPI_Comm,Mat,PetscInt,MatReuse,Mat*); 572bd481603SBarry Smith 573bd481603SBarry Smith /*MC 5741620fd73SBarry Smith MatSetValue - Set a single entry into a matrix. 5751620fd73SBarry Smith 5761620fd73SBarry Smith Not collective 5771620fd73SBarry Smith 578a9834a7dSSatish Balay Synopsis: 579a9834a7dSSatish Balay #include <petscmat.h> 580a9834a7dSSatish Balay PetscErrorCode MatSetValue(Mat m,PetscInt row,PetscInt col,PetscScalar value,InsertMode mode) 581a9834a7dSSatish Balay 5821620fd73SBarry Smith Input Parameters: 5831620fd73SBarry Smith + m - the matrix 5841620fd73SBarry Smith . row - the row location of the entry 5851620fd73SBarry Smith . col - the column location of the entry 5861620fd73SBarry Smith . value - the value to insert 5871620fd73SBarry Smith - mode - either INSERT_VALUES or ADD_VALUES 5881620fd73SBarry Smith 5891620fd73SBarry Smith Notes: 5901620fd73SBarry Smith For efficiency one should use MatSetValues() and set several or many 5911620fd73SBarry Smith values simultaneously if possible. 5921620fd73SBarry Smith 5931620fd73SBarry Smith Level: beginner 5941620fd73SBarry Smith 5951620fd73SBarry Smith .seealso: MatSetValues(), MatSetValueLocal() 5961620fd73SBarry Smith M*/ 5971620fd73SBarry 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);} 5981620fd73SBarry Smith 5991620fd73SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatGetValue(Mat v,PetscInt i,PetscInt j,PetscScalar *va) {return MatGetValues(v,1,&i,1,&j,va);} 6001620fd73SBarry Smith 6011620fd73SBarry 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);} 6021620fd73SBarry Smith 6031620fd73SBarry Smith /*MC 6040d2a7ff2SSatish Balay MatPreallocateInitialize - Begins the block of code that will count the number of nonzeros per 605bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 606bd481603SBarry Smith 607bd481603SBarry Smith Synopsis: 608aaa7dc30SBarry Smith #include <petscmat.h> 609c1ac3661SBarry Smith PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz) 610bd481603SBarry Smith 611bd481603SBarry Smith Collective on MPI_Comm 612bd481603SBarry Smith 613bd481603SBarry Smith Input Parameters: 614bd481603SBarry Smith + comm - the communicator that will share the eventually allocated matrix 615859fcb39SBarry Smith . nrows - the number of LOCAL rows in the matrix 616859fcb39SBarry Smith - ncols - the number of LOCAL columns in the matrix 617bd481603SBarry Smith 618bd481603SBarry Smith Output Parameters: 619bd481603SBarry Smith + dnz - the array that will be passed to the matrix preallocation routines 620bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 621bd481603SBarry Smith 622bd481603SBarry Smith Level: intermediate 623bd481603SBarry Smith 624bd481603SBarry Smith Notes: 625a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 626bd481603SBarry Smith 6271d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 628bd481603SBarry Smith 629bd481603SBarry Smith Use MatPreallocateInitializeSymmetric() for symmetric matrices (MPISBAIJ matrices) 630bd481603SBarry Smith 6311620fd73SBarry Smith This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize(). 6321620fd73SBarry Smith 633bd481603SBarry Smith Concepts: preallocation^Matrix 634bd481603SBarry Smith 635d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(), 636d6e23781SBarry Smith MatPreallocateInitializeSymmetric(), MatPreallocateSymmetricSetLocalBlock() 637bd481603SBarry Smith M*/ 638c4f061fbSSatish Balay #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \ 6397c922b88SBarry Smith { \ 640a89bc211SBarry Smith PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__start,__end; \ 6411795a4d1SJed Brown _4_ierr = PetscCalloc2(__nrows,&dnz,__nrows,&onz);CHKERRQ(_4_ierr); \ 6421795a4d1SJed Brown __start = 0; __end = __start; \ 643c1ac3661SBarry Smith _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ctmp;\ 644a89bc211SBarry Smith _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows; 6457c922b88SBarry Smith 646bd481603SBarry Smith /*MC 647bd481603SBarry Smith MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be 648bd481603SBarry Smith inserted using a local number of the rows and columns 649bd481603SBarry Smith 650bd481603SBarry Smith Synopsis: 651aaa7dc30SBarry Smith #include <petscmat.h> 652c1ac3661SBarry Smith PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 653bd481603SBarry Smith 654bd481603SBarry Smith Not Collective 655bd481603SBarry Smith 656bd481603SBarry Smith Input Parameters: 657784ac674SJed Brown + map - the row mapping from local numbering to global numbering 658bd481603SBarry Smith . nrows - the number of rows indicated 6591d73ed98SBarry Smith . rows - the indices of the rows 660784ac674SJed Brown . cmap - the column mapping from local to global numbering 661bd481603SBarry Smith . ncols - the number of columns in the matrix 662bd481603SBarry Smith . cols - the columns indicated 663bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 664bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 665bd481603SBarry Smith 666bd481603SBarry Smith Level: intermediate 667bd481603SBarry Smith 668bd481603SBarry Smith Notes: 669a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 670bd481603SBarry Smith 6711d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 672bd481603SBarry Smith 673bd481603SBarry Smith Concepts: preallocation^Matrix 674bd481603SBarry Smith 675d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(), 676c1154cd5SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocalRemoveDups() 677bd481603SBarry Smith M*/ 678784ac674SJed Brown #define MatPreallocateSetLocal(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \ 679c4f061fbSSatish Balay {\ 680c1ac3661SBarry Smith PetscInt __l;\ 681784ac674SJed Brown _4_ierr = ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\ 682784ac674SJed Brown _4_ierr = ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\ 683c4f061fbSSatish Balay for (__l=0;__l<nrows;__l++) {\ 684ef66eb69SBarry Smith _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 685c4f061fbSSatish Balay }\ 686c4f061fbSSatish Balay } 687c4f061fbSSatish Balay 688bd481603SBarry Smith /*MC 689c1154cd5SBarry Smith MatPreallocateSetLocalRemoveDups - Indicates the locations (rows and columns) in the matrix where nonzeros will be 690c1154cd5SBarry Smith inserted using a local number of the rows and columns. This version removes any duplicate columns in cols 691c1154cd5SBarry Smith 692c1154cd5SBarry Smith Synopsis: 693c1154cd5SBarry Smith #include <petscmat.h> 694c1154cd5SBarry Smith PetscErrorCode MatPreallocateSetLocalRemoveDups(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 695c1154cd5SBarry Smith 696c1154cd5SBarry Smith Not Collective 697c1154cd5SBarry Smith 698c1154cd5SBarry Smith Input Parameters: 699c1154cd5SBarry Smith + map - the row mapping from local numbering to global numbering 700c1154cd5SBarry Smith . nrows - the number of rows indicated 701c1154cd5SBarry Smith . rows - the indices of the rows (these values are mapped to the global values) 702c1154cd5SBarry Smith . cmap - the column mapping from local to global numbering 703c1154cd5SBarry Smith . ncols - the number of columns in the matrix (this value will be changed if duplicate columns are found) 704c1154cd5SBarry Smith . cols - the columns indicated (these values are mapped to the global values, they are then sorted and duplicates removed) 705c1154cd5SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 706c1154cd5SBarry Smith - ozn - the other array passed to the matrix preallocation routines 707c1154cd5SBarry Smith 708c1154cd5SBarry Smith Level: intermediate 709c1154cd5SBarry Smith 710c1154cd5SBarry Smith Notes: 711c1154cd5SBarry Smith See Users-Manual: ch_performance for more details. 712c1154cd5SBarry Smith 713c1154cd5SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 714c1154cd5SBarry Smith 715c1154cd5SBarry Smith Concepts: preallocation^Matrix 716c1154cd5SBarry Smith 717c1154cd5SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(), 718c1154cd5SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocal() 719c1154cd5SBarry Smith M*/ 720c1154cd5SBarry Smith #define MatPreallocateSetLocalRemoveDups(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \ 721c1154cd5SBarry Smith {\ 722c1154cd5SBarry Smith PetscInt __l;\ 723c1154cd5SBarry Smith _4_ierr = ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\ 724c1154cd5SBarry Smith _4_ierr = ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\ 725c1154cd5SBarry Smith _4_ierr = PetscSortRemoveDupsInt(&ncols,cols);CHKERRQ(_4_ierr);\ 726c1154cd5SBarry Smith for (__l=0;__l<nrows;__l++) {\ 727c1154cd5SBarry Smith _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 728c1154cd5SBarry Smith }\ 729c1154cd5SBarry Smith } 730c1154cd5SBarry Smith 731c1154cd5SBarry Smith /*MC 732d6e23781SBarry Smith MatPreallocateSetLocalBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be 733bd481603SBarry Smith inserted using a local number of the rows and columns 734bd481603SBarry Smith 735bd481603SBarry Smith Synopsis: 736aaa7dc30SBarry Smith #include <petscmat.h> 737d6e23781SBarry Smith PetscErrorCode MatPreallocateSetLocalBlock(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 738d6e23781SBarry Smith 739d6e23781SBarry Smith Not Collective 740d6e23781SBarry Smith 741d6e23781SBarry Smith Input Parameters: 742d6e23781SBarry Smith + map - the row mapping from local numbering to global numbering 743d6e23781SBarry Smith . nrows - the number of rows indicated 744d6e23781SBarry Smith . rows - the indices of the rows 745d6e23781SBarry Smith . cmap - the column mapping from local to global numbering 746d6e23781SBarry Smith . ncols - the number of columns in the matrix 747d6e23781SBarry Smith . cols - the columns indicated 748d6e23781SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 749d6e23781SBarry Smith - ozn - the other array passed to the matrix preallocation routines 750d6e23781SBarry Smith 751d6e23781SBarry Smith Level: intermediate 752d6e23781SBarry Smith 753d6e23781SBarry Smith Notes: 754d6e23781SBarry Smith See Users-Manual: ch_performance for more details. 755d6e23781SBarry Smith 756d6e23781SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 757d6e23781SBarry Smith 758d6e23781SBarry Smith Concepts: preallocation^Matrix 759d6e23781SBarry Smith 760d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(), 761d6e23781SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock() 762d6e23781SBarry Smith M*/ 763d6e23781SBarry Smith #define MatPreallocateSetLocalBlock(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \ 764d6e23781SBarry Smith {\ 765d6e23781SBarry Smith PetscInt __l;\ 766d6e23781SBarry Smith _4_ierr = ISLocalToGlobalMappingApplyBlock(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\ 767d6e23781SBarry Smith _4_ierr = ISLocalToGlobalMappingApplyBlock(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\ 768d6e23781SBarry Smith for (__l=0;__l<nrows;__l++) {\ 769d6e23781SBarry Smith _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 770d6e23781SBarry Smith }\ 771d6e23781SBarry Smith } 772d6e23781SBarry Smith 773d6e23781SBarry Smith /*MC 774d6e23781SBarry Smith MatPreallocateSymmetricSetLocalBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be 775d6e23781SBarry Smith inserted using a local number of the rows and columns 776d6e23781SBarry Smith 777d6e23781SBarry Smith Synopsis: 778d6e23781SBarry Smith #include <petscmat.h> 779d6e23781SBarry Smith PetscErrorCode MatPreallocateSymmetricSetLocalBlock(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 780bd481603SBarry Smith 781bd481603SBarry Smith Not Collective 782bd481603SBarry Smith 783bd481603SBarry Smith Input Parameters: 784bd481603SBarry Smith + map - the mapping between local numbering and global numbering 785bd481603SBarry Smith . nrows - the number of rows indicated 7861d73ed98SBarry Smith . rows - the indices of the rows 787bd481603SBarry Smith . ncols - the number of columns in the matrix 788bd481603SBarry Smith . cols - the columns indicated 789bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 790bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 791bd481603SBarry Smith 792bd481603SBarry Smith Level: intermediate 793bd481603SBarry Smith 794bd481603SBarry Smith Notes: 795a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 796bd481603SBarry Smith 797bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 798bd481603SBarry Smith 799bd481603SBarry Smith Concepts: preallocation^Matrix 800bd481603SBarry Smith 801d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateInitialize(), 802d6e23781SBarry Smith MatPreallocateInitialize(), MatPreallocateSetLocal() 803bd481603SBarry Smith M*/ 804d6e23781SBarry Smith #define MatPreallocateSymmetricSetLocalBlock(map,nrows,rows,ncols,cols,dnz,onz) 0;\ 805d3d32019SBarry Smith {\ 806c1ac3661SBarry Smith PetscInt __l;\ 807d6e23781SBarry Smith _4_ierr = ISLocalToGlobalMappingApplyBlock(map,nrows,rows,rows);CHKERRQ(_4_ierr);\ 808d6e23781SBarry Smith _4_ierr = ISLocalToGlobalMappingApplyBlock(map,ncols,cols,cols);CHKERRQ(_4_ierr);\ 809d3d32019SBarry Smith for (__l=0;__l<nrows;__l++) {\ 810d6e23781SBarry Smith _4_ierr = MatPreallocateSymmetricSetBlock((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 811d3d32019SBarry Smith }\ 812d3d32019SBarry Smith } 813bd481603SBarry Smith /*MC 814bd481603SBarry Smith MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be 815bd481603SBarry Smith inserted using a local number of the rows and columns 816bd481603SBarry Smith 817bd481603SBarry Smith Synopsis: 818aaa7dc30SBarry Smith #include <petscmat.h> 819c1ac3661SBarry Smith PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 820bd481603SBarry Smith 821bd481603SBarry Smith Not Collective 822bd481603SBarry Smith 823bd481603SBarry Smith Input Parameters: 82464075487SBarry Smith + row - the row 825bd481603SBarry Smith . ncols - the number of columns in the matrix 826a50f8bf6SHong Zhang - cols - the columns indicated 827a50f8bf6SHong Zhang 828a50f8bf6SHong Zhang Output Parameters: 829a50f8bf6SHong Zhang + dnz - the array that will be passed to the matrix preallocation routines 830bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 831bd481603SBarry Smith 832bd481603SBarry Smith Level: intermediate 833bd481603SBarry Smith 834bd481603SBarry Smith Notes: 835a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 836bd481603SBarry Smith 837bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 838bd481603SBarry Smith 8391620fd73SBarry Smith This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize(). 8401620fd73SBarry Smith 841bd481603SBarry Smith Concepts: preallocation^Matrix 842bd481603SBarry Smith 843d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(), 844d6e23781SBarry Smith MatPreallocateInitialize(), MatPreallocateSetLocal() 845bd481603SBarry Smith M*/ 846c4f061fbSSatish Balay #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\ 847c1ac3661SBarry Smith { PetscInt __i; \ 848e32f2f54SBarry 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);\ 849e32f2f54SBarry 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);\ 8507c922b88SBarry Smith for (__i=0; __i<nc; __i++) {\ 85164075487SBarry Smith if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \ 8527cc688d7SBarry Smith else dnz[row - __rstart]++;\ 8537c922b88SBarry Smith }\ 8547c922b88SBarry Smith } 8557c922b88SBarry Smith 856bd481603SBarry Smith /*MC 857d6e23781SBarry Smith MatPreallocateSymmetricSetBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be 858bd481603SBarry Smith inserted using a local number of the rows and columns 859bd481603SBarry Smith 860bd481603SBarry Smith Synopsis: 861aaa7dc30SBarry Smith #include <petscmat.h> 862d6e23781SBarry Smith PetscErrorCode MatPreallocateSymmetricSetBlock(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 863bd481603SBarry Smith 864bd481603SBarry Smith Not Collective 865bd481603SBarry Smith 866bd481603SBarry Smith Input Parameters: 867bd481603SBarry Smith + nrows - the number of rows indicated 8681d73ed98SBarry Smith . rows - the indices of the rows 869bd481603SBarry Smith . ncols - the number of columns in the matrix 870bd481603SBarry Smith . cols - the columns indicated 871bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 872bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 873bd481603SBarry Smith 874bd481603SBarry Smith Level: intermediate 875bd481603SBarry Smith 876bd481603SBarry Smith Notes: 877a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 878bd481603SBarry Smith 879bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 880bd481603SBarry Smith 8811620fd73SBarry Smith This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize(). 8821620fd73SBarry Smith 883bd481603SBarry Smith Concepts: preallocation^Matrix 884bd481603SBarry Smith 885d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateInitialize(), 886d6e23781SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocal() 887bd481603SBarry Smith M*/ 888d6e23781SBarry Smith #define MatPreallocateSymmetricSetBlock(row,nc,cols,dnz,onz) 0;\ 889c1ac3661SBarry Smith { PetscInt __i; \ 890d3d32019SBarry Smith for (__i=0; __i<nc; __i++) {\ 891d3d32019SBarry Smith if (cols[__i] >= __end) onz[row - __rstart]++; \ 892d3d32019SBarry Smith else if (cols[__i] >= row) dnz[row - __rstart]++;\ 893d3d32019SBarry Smith }\ 894d3d32019SBarry Smith } 895d3d32019SBarry Smith 896bd481603SBarry Smith /*MC 89716371a99SBarry Smith MatPreallocateLocation - An alternative to MatPreallocationSet() that puts the nonzero locations into the matrix if it exists 89816371a99SBarry Smith 89916371a99SBarry Smith Synopsis: 900aaa7dc30SBarry Smith #include <petscmat.h> 90116371a99SBarry Smith PetscErrorCode MatPreallocateLocations(Mat A,PetscInt row,PetscInt ncols,PetscInt *cols,PetscInt *dnz,PetscInt *onz) 90216371a99SBarry Smith 90316371a99SBarry Smith Not Collective 90416371a99SBarry Smith 90516371a99SBarry Smith Input Parameters: 90616371a99SBarry Smith . A - matrix 90716371a99SBarry Smith . row - row where values exist (must be local to this process) 90816371a99SBarry Smith . ncols - number of columns 90916371a99SBarry Smith . cols - columns with nonzeros 91016371a99SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 91116371a99SBarry Smith - ozn - the other array passed to the matrix preallocation routines 91216371a99SBarry Smith 91316371a99SBarry Smith Level: intermediate 91416371a99SBarry Smith 91516371a99SBarry Smith Notes: 916a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 91716371a99SBarry Smith 91816371a99SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 91916371a99SBarry Smith 92016371a99SBarry Smith This is a MACRO not a function because it uses a bunch of variables private to the MatPreallocation.... routines. 92116371a99SBarry Smith 92216371a99SBarry Smith Concepts: preallocation^Matrix 92316371a99SBarry Smith 924d6e23781SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(), 925d6e23781SBarry Smith MatPreallocateSymmetricSetLocalBlock() 92616371a99SBarry Smith M*/ 9270298fd71SBarry 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);} 92816371a99SBarry Smith 92916371a99SBarry Smith 93016371a99SBarry Smith /*MC 9310d2a7ff2SSatish Balay MatPreallocateFinalize - Ends the block of code that will count the number of nonzeros per 932bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 933bd481603SBarry Smith 934bd481603SBarry Smith Synopsis: 935aaa7dc30SBarry Smith #include <petscmat.h> 936c1ac3661SBarry Smith PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz) 937bd481603SBarry Smith 938bd481603SBarry Smith Collective on MPI_Comm 939bd481603SBarry Smith 940bd481603SBarry Smith Input Parameters: 94116371a99SBarry Smith + dnz - the array that was be passed to the matrix preallocation routines 942bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 943bd481603SBarry Smith 944bd481603SBarry Smith Level: intermediate 945bd481603SBarry Smith 946bd481603SBarry Smith Notes: 947a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 948bd481603SBarry Smith 949bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 950bd481603SBarry Smith 9511620fd73SBarry Smith This is a MACRO not a function because it closes the { started in MatPreallocateInitialize(). 9521620fd73SBarry Smith 953bd481603SBarry Smith Concepts: preallocation^Matrix 954bd481603SBarry Smith 955d6e23781SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(), 956d6e23781SBarry Smith MatPreallocateSymmetricSetLocalBlock() 957bd481603SBarry Smith M*/ 958a89bc211SBarry Smith #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree2(dnz,onz);CHKERRQ(_4_ierr);} 9597c922b88SBarry Smith 9607b80b807SBarry Smith /* Routines unique to particular data structures */ 961014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellGetContext(Mat,void *); 962698d4c6aSKris Buschelman 963014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInodeAdjustForInodes(Mat,IS*,IS*); 964014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *); 965698d4c6aSKris Buschelman 966014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetColumnIndices(Mat,PetscInt[]); 967014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetColumnIndices(Mat,PetscInt[]); 968014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 969014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 970014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 971014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJFromTriple(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*,PetscInt,PetscBool); 9727b80b807SBarry Smith 973a96a251dSBarry Smith #define MAT_SKIP_ALLOCATION -4 974a96a251dSBarry Smith 975014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]); 976014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]); 977014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]); 978273d9f13SBarry Smith 979014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 980014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 981014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 982014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []); 983014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 984014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]); 985014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 986014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]); 987014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIDenseSetPreallocation(Mat,PetscScalar[]); 988014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqDenseSetPreallocation(Mat,PetscScalar[]); 9899230625dSJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,const PetscInt*[]); 9909230625dSJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,const PetscInt*[]); 991014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAdjCreateNonemptySubcommMat(Mat,Mat*); 992273d9f13SBarry Smith 9932e1947a5SStefano Zampini PETSC_EXTERN PetscErrorCode MatISSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 994cf0a3239SStefano Zampini PETSC_EXTERN PetscErrorCode MatISSetUpSF(Mat); 995014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqDenseSetLDA(Mat,PetscInt); 996014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDenseGetLocalMatrix(Mat,Mat*); 9971b807ce4Svictorle 998014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStoreValues(Mat); 999014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRetrieveValues(Mat); 10002e8a6d31SBarry Smith 1001014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDAADSetCtx(Mat,void*); 10023a7fca6bSBarry Smith 1003014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFindNonzeroRows(Mat,IS*); 10047cfe41e4SPatrick Farrell PETSC_EXTERN PetscErrorCode MatFindZeroRows(Mat,IS*); 10057b80b807SBarry Smith /* 10067b80b807SBarry Smith These routines are not usually accessed directly, rather solving is 100794b7f48cSBarry Smith done through the KSP and PC interfaces. 10087b80b807SBarry Smith */ 10097b80b807SBarry Smith 101076bdecfbSBarry Smith /*J 10118f6c3df8SBarry Smith MatOrderingType - String with the name of a PETSc matrix ordering 1012d9274352SBarry Smith 1013d9274352SBarry Smith Level: beginner 1014d9274352SBarry Smith 1015d9274352SBarry Smith .seealso: MatGetOrdering() 101676bdecfbSBarry Smith J*/ 101719fd82e9SBarry Smith typedef const char* MatOrderingType; 10182692d6eeSBarry Smith #define MATORDERINGNATURAL "natural" 10192692d6eeSBarry Smith #define MATORDERINGND "nd" 10202692d6eeSBarry Smith #define MATORDERING1WD "1wd" 10212692d6eeSBarry Smith #define MATORDERINGRCM "rcm" 10222692d6eeSBarry Smith #define MATORDERINGQMD "qmd" 10232692d6eeSBarry Smith #define MATORDERINGROWLENGTH "rowlength" 1024510184a7SJed Brown #define MATORDERINGWBM "wbm" 1025fc1bef75SJed Brown #define MATORDERINGSPECTRAL "spectral" 1026312e372aSBarry Smith #define MATORDERINGAMD "amd" /* only works if UMFPACK is installed with PETSc */ 1027b12f92e5SBarry Smith 102819fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatGetOrdering(Mat,MatOrderingType,IS*,IS*); 1029140e18c1SBarry Smith PETSC_EXTERN PetscErrorCode MatGetOrderingList(PetscFunctionList*); 1030bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatOrderingRegister(const char[],PetscErrorCode(*)(Mat,MatOrderingType,IS*,IS*)); 1031140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatOrderingList; 1032d4fbbf0eSBarry Smith 1033014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS); 1034fc1bef75SJed Brown PETSC_EXTERN PetscErrorCode MatCreateLaplacian(Mat,PetscReal,PetscBool,Mat*); 1035a2ce50c7SBarry Smith 1036d91e6319SBarry Smith /*S 1037d90ac83dSHong Zhang MatFactorShiftType - Numeric Shift. 1038d90ac83dSHong Zhang 1039d90ac83dSHong Zhang Level: beginner 1040d90ac83dSHong Zhang 1041d90ac83dSHong Zhang S*/ 1042d90ac83dSHong Zhang typedef enum {MAT_SHIFT_NONE,MAT_SHIFT_NONZERO,MAT_SHIFT_POSITIVE_DEFINITE,MAT_SHIFT_INBLOCKS} MatFactorShiftType; 10436a6fc655SJed Brown PETSC_EXTERN const char *const MatFactorShiftTypes[]; 10445e9742b9SJed Brown PETSC_EXTERN const char *const MatFactorShiftTypesDetail[]; 1045d90ac83dSHong Zhang 1046d90ac83dSHong Zhang /*S 10479aa7eafdSHong Zhang MatFactorError - indicates what type of error in matrix factor 10489aa7eafdSHong Zhang 10499aa7eafdSHong Zhang Level: beginner 10509aa7eafdSHong Zhang 105100e125f8SBarry Smith Developer Notes: Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 105200e125f8SBarry Smith 105300e125f8SBarry Smith .seealso: MatGetFactor() 10549aa7eafdSHong Zhang S*/ 10555cd7cf9dSHong Zhang typedef enum {MAT_FACTOR_NOERROR,MAT_FACTOR_STRUCT_ZEROPIVOT,MAT_FACTOR_NUMERIC_ZEROPIVOT,MAT_FACTOR_OUTMEMORY,MAT_FACTOR_OTHER} MatFactorError; 10569aa7eafdSHong Zhang 105700e125f8SBarry Smith PETSC_EXTERN PetscErrorCode MatFactorGetError(Mat,MatFactorError*); 1058b8b68cfdSBarry Smith PETSC_EXTERN PetscErrorCode MatFactorClearError(Mat); 10597b6c816cSBarry Smith PETSC_EXTERN PetscErrorCode MatFactorGetErrorZeroPivot(Mat,PetscReal*,PetscInt*); 106000e125f8SBarry Smith 10619aa7eafdSHong Zhang /*S 1062422a814eSBarry Smith MatFactorInfo - Data passed into the matrix factorization routines, and information about the resulting factorization 1063b00f7748SHong Zhang 106461cad860SBarry Smith In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE, that is use 106561cad860SBarry Smith $ MatFactorInfo info(MAT_FACTORINFO_SIZE) 1066b00f7748SHong Zhang 106715e8a5b3SHong Zhang Notes: These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC. 1068b00f7748SHong Zhang 106961cad860SBarry Smith You can use MatFactorInfoInitialize() to set default values. 107061cad860SBarry Smith 1071b00f7748SHong Zhang Level: developer 1072b00f7748SHong Zhang 1073d7d82daaSBarry Smith .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(), 1074d7d82daaSBarry Smith MatFactorInfoInitialize() 1075b00f7748SHong Zhang 1076b00f7748SHong Zhang S*/ 1077b00f7748SHong Zhang typedef struct { 107815e8a5b3SHong Zhang PetscReal diagonal_fill; /* force diagonal to fill in if initially not filled */ 107985317021SBarry Smith PetscReal usedt; 108015e8a5b3SHong Zhang PetscReal dt; /* drop tolerance */ 1081b00f7748SHong Zhang PetscReal dtcol; /* tolerance for pivoting */ 108215e8a5b3SHong Zhang PetscReal dtcount; /* maximum nonzeros to be allowed per row */ 108367e28bfeSBarry Smith PetscReal fill; /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */ 1084348344bbSBarry Smith PetscReal levels; /* ICC/ILU(levels) */ 1085bcd9e38bSBarry Smith PetscReal pivotinblocks; /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0 1086bcd9e38bSBarry Smith factorization may be faster if do not pivot */ 108715e8a5b3SHong Zhang PetscReal zeropivot; /* pivot is called zero if less than this */ 1088f4db908eSBarry Smith PetscReal shifttype; /* type of shift added to matrix factor to prevent zero pivots */ 1089f4db908eSBarry Smith PetscReal shiftamount; /* how large the shift is */ 109015e8a5b3SHong Zhang } MatFactorInfo; 1091ffa6d0a5SLois Curfman McInnes 1092014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFactorInfoInitialize(MatFactorInfo*); 10939a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatCholeskyFactor(Mat,IS,const MatFactorInfo*); 10949a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatCholeskyFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*); 10959a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatCholeskyFactorNumeric(Mat,Mat,const MatFactorInfo*); 10969a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatLUFactor(Mat,IS,IS,const MatFactorInfo*); 10979a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatILUFactor(Mat,IS,IS,const MatFactorInfo*); 10989a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatLUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*); 10999a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatILUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*); 11009a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatICCFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*); 11019a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatICCFactor(Mat,IS,const MatFactorInfo*); 11029a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatLUFactorNumeric(Mat,Mat,const MatFactorInfo*); 1103014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*); 1104014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolve(Mat,Vec,Vec); 1105014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatForwardSolve(Mat,Vec,Vec); 1106014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatBackwardSolve(Mat,Vec,Vec); 1107014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveAdd(Mat,Vec,Vec,Vec); 1108014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveTranspose(Mat,Vec,Vec); 1109014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveTransposeAdd(Mat,Vec,Vec,Vec); 1110014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolves(Mat,Vecs,Vecs); 1111*7c2f51b8SStefano Zampini PETSC_EXTERN PetscErrorCode MatSetUnfactored(Mat); 1112*7c2f51b8SStefano Zampini 1113*7c2f51b8SStefano Zampini typedef enum {MAT_FACTOR_SCHUR_UNFACTORED, MAT_FACTOR_SCHUR_FACTORED, MAT_FACTOR_SCHUR_INVERTED} MatFactorSchurStatus; 11148e7ba810SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorSetSchurIS(Mat,IS); 1115*7c2f51b8SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorGetSchurComplement(Mat,Mat*,MatFactorSchurStatus*); 1116*7c2f51b8SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorRestoreSchurComplement(Mat,Mat*,MatFactorSchurStatus); 11175a05ddb0SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorInvertSchurComplement(Mat); 1118*7c2f51b8SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorCreateSchurComplement(Mat,Mat*,MatFactorSchurStatus*); 11195a05ddb0SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorSolveSchurComplement(Mat,Vec,Vec); 11205a05ddb0SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorSolveSchurComplementTranspose(Mat,Vec,Vec); 11216dba178dSStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorFactorizeSchurComplement(Mat); 1122e8ade678SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorSetSchurComplementSolverType(Mat,PetscInt); 1123bb5a7306SBarry Smith 1124d91e6319SBarry Smith /*E 1125d91e6319SBarry Smith MatSORType - What type of (S)SOR to perform 1126bb1eb677SSatish Balay 1127d91e6319SBarry Smith Level: beginner 1128d91e6319SBarry Smith 1129d9274352SBarry Smith May be bitwise ORd together 1130d9274352SBarry Smith 1131af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 1132d91e6319SBarry Smith 11334e7234bfSBarry Smith MatSORType may be bitwise ORd together, so do not change the numbers 11344e7234bfSBarry Smith 113541f059aeSBarry Smith .seealso: MatSOR() 1136d91e6319SBarry Smith E*/ 1137ee50ffe9SBarry Smith typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3, 1138ee50ffe9SBarry Smith SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8, 1139ee50ffe9SBarry Smith SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16, 114084cb2905SBarry Smith SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType; 1141014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSOR(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec); 11428ed539a5SBarry Smith 1143d4fbbf0eSBarry Smith /* 1144639f9d9dSBarry Smith These routines are for efficiently computing Jacobians via finite differences. 1145639f9d9dSBarry Smith */ 1146b12f92e5SBarry Smith 11477086a01eSPeter Brune /*S 11487086a01eSPeter Brune MatColoring - Object for managing the coloring of matrices. 11497086a01eSPeter Brune 11507086a01eSPeter Brune Level: beginner 11517086a01eSPeter Brune 11527086a01eSPeter Brune Concepts: matrix, coloring 11537086a01eSPeter Brune 11547086a01eSPeter Brune .seealso: MatFDColoringCreate() ISColoring MatFDColoring 11557086a01eSPeter Brune S*/ 11567086a01eSPeter Brune typedef struct _p_MatColoring* MatColoring; 115776bdecfbSBarry Smith /*J 11588f6c3df8SBarry Smith MatColoringType - String with the name of a PETSc matrix coloring 1159d9274352SBarry Smith 1160d9274352SBarry Smith Level: beginner 1161d9274352SBarry Smith 11627086a01eSPeter Brune .seealso: MatColoringSetType(), MatColoring 116376bdecfbSBarry Smith J*/ 11647086a01eSPeter Brune 116519fd82e9SBarry Smith typedef const char* MatColoringType; 11665567d71aSPeter Brune #define MATCOLORINGJP "jp" 1167b9acb617SPeter Brune #define MATCOLORINGPOWER "power" 11682692d6eeSBarry Smith #define MATCOLORINGNATURAL "natural" 11692692d6eeSBarry Smith #define MATCOLORINGSL "sl" 11702692d6eeSBarry Smith #define MATCOLORINGLF "lf" 11712692d6eeSBarry Smith #define MATCOLORINGID "id" 11728121bdceSPeter Brune #define MATCOLORINGGREEDY "greedy" 1173b12f92e5SBarry Smith 11748ac6da0aSPeter Brune /*E 11758ac6da0aSPeter Brune MatColoringWeightType - Type of weight scheme 11768ac6da0aSPeter Brune 11778ac6da0aSPeter Brune Not Collective 11788ac6da0aSPeter Brune 11798ac6da0aSPeter Brune + MAT_COLORING_RANDOM - Random weights 11808ac6da0aSPeter Brune . MAT_COLORING_LEXICAL - Lexical weighting based upon global numbering. 11818ac6da0aSPeter Brune - MAT_COLORING_LF - Last-first weighting. 11828ac6da0aSPeter Brune 11838ac6da0aSPeter Brune Level: intermediate 11848ac6da0aSPeter Brune 1185af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 11868ac6da0aSPeter Brune 11878ac6da0aSPeter Brune .seealso: MatCUSPSetFormat(), MatCUSPFormatOperation 11888ac6da0aSPeter Brune E*/ 1189301c140bSPeter Brune typedef enum {MAT_COLORING_WEIGHT_RANDOM,MAT_COLORING_WEIGHT_LEXICAL,MAT_COLORING_WEIGHT_LF,MAT_COLORING_WEIGHT_SL} MatColoringWeightType; 11908ac6da0aSPeter Brune 1191335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringCreate(Mat,MatColoring*); 1192d7f2a307SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringGetDegrees(Mat,PetscInt,PetscInt*); 1193335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringDestroy(MatColoring*); 1194335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringView(MatColoring,PetscViewer); 1195335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetType(MatColoring,MatColoringType); 1196335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetFromOptions(MatColoring); 1197335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetDistance(MatColoring,PetscInt); 1198335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringGetDistance(MatColoring,PetscInt*); 1199335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetMaxColors(MatColoring,PetscInt); 1200335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringGetMaxColors(MatColoring,PetscInt*); 1201335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringApply(MatColoring,ISColoring*); 1202335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringRegister(const char[],PetscErrorCode(*)(MatColoring)); 1203014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*); 12048ac6da0aSPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetWeightType(MatColoring,MatColoringWeightType); 1205c29971fcSPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetWeights(MatColoring,PetscReal*,PetscInt*); 12068ac6da0aSPeter Brune PETSC_EXTERN PetscErrorCode MatColoringCreateWeights(MatColoring,PetscReal **,PetscInt **lperm); 1207cb394dc5SBarry Smith PETSC_EXTERN PetscErrorCode MatColoringTestValid(MatColoring,ISColoring); 1208639f9d9dSBarry Smith 1209d9274352SBarry Smith /*S 1210d9274352SBarry Smith MatFDColoring - Object for computing a sparse Jacobian via finite differences 1211d9274352SBarry Smith and coloring 1212639f9d9dSBarry Smith 1213d9274352SBarry Smith Level: beginner 1214d9274352SBarry Smith 1215d9274352SBarry Smith Concepts: coloring, sparse Jacobian, finite differences 1216d9274352SBarry Smith 1217d9274352SBarry Smith .seealso: MatFDColoringCreate() 1218d9274352SBarry Smith S*/ 1219e2a1c21fSSatish Balay typedef struct _p_MatFDColoring* MatFDColoring; 1220639f9d9dSBarry Smith 1221014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringCreate(Mat,ISColoring,MatFDColoring *); 1222014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringDestroy(MatFDColoring*); 1223014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringView(MatFDColoring,PetscViewer); 1224014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*); 1225014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**); 1226014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal); 1227014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetFromOptions(MatFDColoring); 1228d1e9a80fSBarry Smith PETSC_EXTERN PetscErrorCode MatFDColoringApply(Mat,MatFDColoring,Vec,void *); 1229014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetF(MatFDColoring,Vec); 1230014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,PetscInt*[]); 1231f86b9fbaSHong Zhang PETSC_EXTERN PetscErrorCode MatFDColoringSetUp(Mat,ISColoring,MatFDColoring); 1232f86b9fbaSHong Zhang PETSC_EXTERN PetscErrorCode MatFDColoringSetBlockSize(MatFDColoring,PetscInt,PetscInt); 1233f86b9fbaSHong Zhang 1234b1683b59SHong Zhang 1235b1683b59SHong Zhang /*S 1236b9af6bddSHong Zhang MatTransposeColoring - Object for computing a sparse matrix product C=A*B^T via coloring 1237b1683b59SHong Zhang 1238b1683b59SHong Zhang Level: beginner 1239b1683b59SHong Zhang 1240b1683b59SHong Zhang Concepts: coloring, sparse matrix product 1241b1683b59SHong Zhang 1242b9af6bddSHong Zhang .seealso: MatTransposeColoringCreate() 1243b1683b59SHong Zhang S*/ 1244b9af6bddSHong Zhang typedef struct _p_MatTransposeColoring* MatTransposeColoring; 1245b1683b59SHong Zhang 1246014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeColoringCreate(Mat,ISColoring,MatTransposeColoring *); 1247014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransColoringApplySpToDen(MatTransposeColoring,Mat,Mat); 1248014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransColoringApplyDenToSp(MatTransposeColoring,Mat,Mat); 1249014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeColoringDestroy(MatTransposeColoring*); 1250b1683b59SHong Zhang 1251639f9d9dSBarry Smith /* 12520752156aSBarry Smith These routines are for partitioning matrices: currently used only 12533eda8832SBarry Smith for adjacency matrix, MatCreateMPIAdj(). 12540752156aSBarry Smith */ 1255ca161407SBarry Smith 1256d9274352SBarry Smith /*S 1257d9274352SBarry Smith MatPartitioning - Object for managing the partitioning of a matrix or graph 1258d9274352SBarry Smith 1259d9274352SBarry Smith Level: beginner 1260d9274352SBarry Smith 1261d9274352SBarry Smith Concepts: partitioning 1262d9274352SBarry Smith 1263743a1026Svictorle .seealso: MatPartitioningCreate(), MatPartitioningType 1264d9274352SBarry Smith S*/ 126591e9ee9fSBarry Smith typedef struct _p_MatPartitioning* MatPartitioning; 1266d9274352SBarry Smith 126776bdecfbSBarry Smith /*J 12688f6c3df8SBarry Smith MatPartitioningType - String with the name of a PETSc matrix partitioning 1269d9274352SBarry Smith 1270d9274352SBarry Smith Level: beginner 12710d04baf8SBarry Smith dm 1272b547a13bSvictorle .seealso: MatPartitioningCreate(), MatPartitioning 127376bdecfbSBarry Smith J*/ 127419fd82e9SBarry Smith typedef const char* MatPartitioningType; 12752692d6eeSBarry Smith #define MATPARTITIONINGCURRENT "current" 1276aa1e27eaSFande Kong #define MATPARTITIONINGAVERAGE "average" 12772692d6eeSBarry Smith #define MATPARTITIONINGSQUARE "square" 12782692d6eeSBarry Smith #define MATPARTITIONINGPARMETIS "parmetis" 12792692d6eeSBarry Smith #define MATPARTITIONINGCHACO "chaco" 12802692d6eeSBarry Smith #define MATPARTITIONINGPARTY "party" 128150d91057SBarry Smith #define MATPARTITIONINGPTSCOTCH "ptscotch" 128288d2ac2bSFande Kong #define MATPARTITIONINGHIERARCH "hierarch" 128371306c60Sjroman 1284ca161407SBarry Smith 1285014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningCreate(MPI_Comm,MatPartitioning*); 128619fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningSetType(MatPartitioning,MatPartitioningType); 1287014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetNParts(MatPartitioning,PetscInt); 1288014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetAdjacency(MatPartitioning,Mat); 1289014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]); 1290014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []); 1291014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningApply(MatPartitioning,IS*); 1292014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningDestroy(MatPartitioning*); 12932aabb6bbSBarry Smith 1294bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningRegister(const char[],PetscErrorCode (*)(MatPartitioning)); 129530de9b25SBarry Smith 1296f1144a30SSatish Balay 12972bad1931SBarry Smith 1298014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningView(MatPartitioning,PetscViewer); 1299014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetFromOptions(MatPartitioning); 130019fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningGetType(MatPartitioning,MatPartitioningType*); 1301ca161407SBarry Smith 130227538973SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetRepartition(MatPartitioning part); 1303014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetCoarseSequential(MatPartitioning); 1304014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningParmetisGetEdgeCut(MatPartitioning, PetscInt *); 13050752156aSBarry Smith 1306b6956c12SJose Roman typedef enum { MP_CHACO_MULTILEVEL=1,MP_CHACO_SPECTRAL=2,MP_CHACO_LINEAR=4,MP_CHACO_RANDOM=5,MP_CHACO_SCATTERED=6 } MPChacoGlobalType; 13076a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoGlobalTypes[]; 1308b6956c12SJose Roman typedef enum { MP_CHACO_KERNIGHAN=1,MP_CHACO_NONE=2 } MPChacoLocalType; 13096a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoLocalTypes[]; 1310b6956c12SJose Roman typedef enum { MP_CHACO_LANCZOS=0,MP_CHACO_RQI=1 } MPChacoEigenType; 13116a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoEigenTypes[]; 1312b6956c12SJose Roman 1313014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetGlobal(MatPartitioning,MPChacoGlobalType); 1314014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetGlobal(MatPartitioning,MPChacoGlobalType*); 1315014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetLocal(MatPartitioning,MPChacoLocalType); 1316014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetLocal(MatPartitioning,MPChacoLocalType*); 1317014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal); 1318014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType); 1319014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenSolver(MatPartitioning,MPChacoEigenType*); 1320014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenTol(MatPartitioning,PetscReal); 1321014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenTol(MatPartitioning,PetscReal*); 1322014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenNumber(MatPartitioning,PetscInt); 1323014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenNumber(MatPartitioning,PetscInt*); 132471306c60Sjroman 132571306c60Sjroman #define MP_PARTY_OPT "opt" 132671306c60Sjroman #define MP_PARTY_LIN "lin" 132771306c60Sjroman #define MP_PARTY_SCA "sca" 132871306c60Sjroman #define MP_PARTY_RAN "ran" 132971306c60Sjroman #define MP_PARTY_GBF "gbf" 133071306c60Sjroman #define MP_PARTY_GCF "gcf" 133171306c60Sjroman #define MP_PARTY_BUB "bub" 133271306c60Sjroman #define MP_PARTY_DEF "def" 1333014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetGlobal(MatPartitioning,const char*); 133471306c60Sjroman #define MP_PARTY_HELPFUL_SETS "hs" 133571306c60Sjroman #define MP_PARTY_KERNIGHAN_LIN "kl" 133671306c60Sjroman #define MP_PARTY_NONE "no" 1337014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetLocal(MatPartitioning,const char*); 1338014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal); 1339014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetBipart(MatPartitioning,PetscBool); 1340014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscBool); 134171306c60Sjroman 134250d91057SBarry Smith typedef enum { MP_PTSCOTCH_QUALITY,MP_PTSCOTCH_SPEED,MP_PTSCOTCH_BALANCE,MP_PTSCOTCH_SAFETY,MP_PTSCOTCH_SCALABILITY } MPPTScotchStrategyType; 13436a6fc655SJed Brown PETSC_EXTERN const char *const MPPTScotchStrategyTypes[]; 1344e0f1cffaSJose Roman 1345014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetImbalance(MatPartitioning,PetscReal); 1346014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetImbalance(MatPartitioning,PetscReal*); 1347014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetStrategy(MatPartitioning,MPPTScotchStrategyType); 1348014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetStrategy(MatPartitioning,MPPTScotchStrategyType*); 134971306c60Sjroman 1350b43b03e9SMark F. Adams /* 13516294fa2bSFande Kong * hierarchical partitioning 13526294fa2bSFande Kong */ 13534e713f55SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalGetFineparts(MatPartitioning,IS*); 13544e713f55SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalGetCoarseparts(MatPartitioning,IS*); 13554e713f55SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalSetNcoarseparts(MatPartitioning,PetscInt); 13564e713f55SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalSetNfineparts(MatPartitioning, PetscInt); 13576294fa2bSFande Kong 1358014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMeshToVertexGraph(Mat,PetscInt,Mat*); 1359014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMeshToCellGraph(Mat,PetscInt,Mat*); 1360591294e4SBarry Smith 13610752156aSBarry Smith /* 1362af0996ceSBarry Smith If you add entries here you must also add them to petsc/finclude/petscmat.h 1363d4fbbf0eSBarry Smith */ 13641c1c02c0SLois Curfman McInnes typedef enum { MATOP_SET_VALUES=0, 13651c1c02c0SLois Curfman McInnes MATOP_GET_ROW=1, 13661c1c02c0SLois Curfman McInnes MATOP_RESTORE_ROW=2, 13671c1c02c0SLois Curfman McInnes MATOP_MULT=3, 13681c1c02c0SLois Curfman McInnes MATOP_MULT_ADD=4, 13697c922b88SBarry Smith MATOP_MULT_TRANSPOSE=5, 13707c922b88SBarry Smith MATOP_MULT_TRANSPOSE_ADD=6, 13711c1c02c0SLois Curfman McInnes MATOP_SOLVE=7, 13721c1c02c0SLois Curfman McInnes MATOP_SOLVE_ADD=8, 13737c922b88SBarry Smith MATOP_SOLVE_TRANSPOSE=9, 13747c922b88SBarry Smith MATOP_SOLVE_TRANSPOSE_ADD=10, 13751c1c02c0SLois Curfman McInnes MATOP_LUFACTOR=11, 13761c1c02c0SLois Curfman McInnes MATOP_CHOLESKYFACTOR=12, 1377a714c06dSBarry Smith MATOP_SOR=13, 13781c1c02c0SLois Curfman McInnes MATOP_TRANSPOSE=14, 13791c1c02c0SLois Curfman McInnes MATOP_GETINFO=15, 13801c1c02c0SLois Curfman McInnes MATOP_EQUAL=16, 13811c1c02c0SLois Curfman McInnes MATOP_GET_DIAGONAL=17, 13821c1c02c0SLois Curfman McInnes MATOP_DIAGONAL_SCALE=18, 13831c1c02c0SLois Curfman McInnes MATOP_NORM=19, 13841c1c02c0SLois Curfman McInnes MATOP_ASSEMBLY_BEGIN=20, 13851c1c02c0SLois Curfman McInnes MATOP_ASSEMBLY_END=21, 1386d519adbfSMatthew Knepley MATOP_SET_OPTION=22, 1387d519adbfSMatthew Knepley MATOP_ZERO_ENTRIES=23, 1388d519adbfSMatthew Knepley MATOP_ZERO_ROWS=24, 1389d519adbfSMatthew Knepley MATOP_LUFACTOR_SYMBOLIC=25, 1390d519adbfSMatthew Knepley MATOP_LUFACTOR_NUMERIC=26, 1391d519adbfSMatthew Knepley MATOP_CHOLESKY_FACTOR_SYMBOLIC=27, 1392d519adbfSMatthew Knepley MATOP_CHOLESKY_FACTOR_NUMERIC=28, 1393d519adbfSMatthew Knepley MATOP_SETUP_PREALLOCATION=29, 1394d519adbfSMatthew Knepley MATOP_ILUFACTOR_SYMBOLIC=30, 1395d519adbfSMatthew Knepley MATOP_ICCFACTOR_SYMBOLIC=31, 1396a5b7ff6bSBarry Smith MATOP_GET_DIAGONAL_BLOCK=32, 13979d39f69dSJed Brown /* MATOP_PLACEHOLDER_33=33, */ 1398d519adbfSMatthew Knepley MATOP_DUPLICATE=34, 1399d519adbfSMatthew Knepley MATOP_FORWARD_SOLVE=35, 1400d519adbfSMatthew Knepley MATOP_BACKWARD_SOLVE=36, 1401d519adbfSMatthew Knepley MATOP_ILUFACTOR=37, 1402d519adbfSMatthew Knepley MATOP_ICCFACTOR=38, 1403d519adbfSMatthew Knepley MATOP_AXPY=39, 14047dae84e0SHong Zhang MATOP_CREATE_SUBMATRICES=40, 1405d519adbfSMatthew Knepley MATOP_INCREASE_OVERLAP=41, 1406d519adbfSMatthew Knepley MATOP_GET_VALUES=42, 1407d519adbfSMatthew Knepley MATOP_COPY=43, 1408d519adbfSMatthew Knepley MATOP_GET_ROW_MAX=44, 1409d519adbfSMatthew Knepley MATOP_SCALE=45, 1410d519adbfSMatthew Knepley MATOP_SHIFT=46, 141135153367SBarry Smith MATOP_DIAGONAL_SET=47, 14129d39f69dSJed Brown MATOP_ZERO_ROWS_COLUMNS=48, 14139d39f69dSJed Brown MATOP_SET_RANDOM=49, 1414d519adbfSMatthew Knepley MATOP_GET_ROW_IJ=50, 1415d519adbfSMatthew Knepley MATOP_RESTORE_ROW_IJ=51, 1416d519adbfSMatthew Knepley MATOP_GET_COLUMN_IJ=52, 1417d519adbfSMatthew Knepley MATOP_RESTORE_COLUMN_IJ=53, 1418d519adbfSMatthew Knepley MATOP_FDCOLORING_CREATE=54, 1419d519adbfSMatthew Knepley MATOP_COLORING_PATCH=55, 1420d519adbfSMatthew Knepley MATOP_SET_UNFACTORED=56, 1421d519adbfSMatthew Knepley MATOP_PERMUTE=57, 1422d519adbfSMatthew Knepley MATOP_SET_VALUES_BLOCKED=58, 14237dae84e0SHong Zhang MATOP_CREATE_SUBMATRIX=59, 1424d519adbfSMatthew Knepley MATOP_DESTROY=60, 1425d519adbfSMatthew Knepley MATOP_VIEW=61, 1426d519adbfSMatthew Knepley MATOP_CONVERT_FROM=62, 14277bab7c10SHong Zhang MATOP_MATMAT_MULT=63, 14287bab7c10SHong Zhang MATOP_MATMAT_MULT_SYMBOLIC=64, 14297bab7c10SHong Zhang MATOP_MATMAT_MULT_NUMERIC=65, 1430d519adbfSMatthew Knepley MATOP_SET_LOCAL_TO_GLOBAL_MAP=66, 1431d519adbfSMatthew Knepley MATOP_SET_VALUES_LOCAL=67, 1432d519adbfSMatthew Knepley MATOP_ZERO_ROWS_LOCAL=68, 1433d519adbfSMatthew Knepley MATOP_GET_ROW_MAX_ABS=69, 1434d519adbfSMatthew Knepley MATOP_GET_ROW_MIN_ABS=70, 1435d519adbfSMatthew Knepley MATOP_CONVERT=71, 1436d519adbfSMatthew Knepley MATOP_SET_COLORING=72, 14379d39f69dSJed Brown /* MATOP_PLACEHOLDER_73=73, */ 1438d519adbfSMatthew Knepley MATOP_SET_VALUES_ADIFOR=74, 1439d519adbfSMatthew Knepley MATOP_FD_COLORING_APPLY=75, 1440d519adbfSMatthew Knepley MATOP_SET_FROM_OPTIONS=76, 1441bda6c4cbSBarry Smith MATOP_MULT_CONSTRAINED=77, 1442bda6c4cbSBarry Smith MATOP_MULT_TRANSPOSE_CONSTRAIN=78, 14439d39f69dSJed Brown MATOP_FIND_ZERO_DIAGONALS=79, 1444d519adbfSMatthew Knepley MATOP_MULT_MULTIPLE=80, 1445d519adbfSMatthew Knepley MATOP_SOLVE_MULTIPLE=81, 1446d519adbfSMatthew Knepley MATOP_GET_INERTIA=82, 1447d519adbfSMatthew Knepley MATOP_LOAD=83, 1448d519adbfSMatthew Knepley MATOP_IS_SYMMETRIC=84, 1449d519adbfSMatthew Knepley MATOP_IS_HERMITIAN=85, 1450d519adbfSMatthew Knepley MATOP_IS_STRUCTURALLY_SYMMETRIC=86, 1451820469bcSHong Zhang MATOP_SET_VALUES_BLOCKEDLOCAL=87, 1452d519adbfSMatthew Knepley MATOP_GET_VECS=88, 1453d519adbfSMatthew Knepley MATOP_MAT_MULT=89, 1454d519adbfSMatthew Knepley MATOP_MAT_MULT_SYMBOLIC=90, 1455d519adbfSMatthew Knepley MATOP_MAT_MULT_NUMERIC=91, 1456d519adbfSMatthew Knepley MATOP_PTAP=92, 1457d519adbfSMatthew Knepley MATOP_PTAP_SYMBOLIC=93, 1458d519adbfSMatthew Knepley MATOP_PTAP_NUMERIC=94, 1459bda6c4cbSBarry Smith MATOP_MAT_TRANSPOSE_MULT=95, 1460bda6c4cbSBarry Smith MATOP_MAT_TRANSPOSE_MULT_SYMBO=96, 1461bda6c4cbSBarry Smith MATOP_MAT_TRANSPOSE_MULT_NUMER=97, 14629d39f69dSJed Brown /* MATOP_PLACEHOLDER_98=98, */ 14639d39f69dSJed Brown /* MATOP_PLACEHOLDER_99=99, */ 14649d39f69dSJed Brown /* MATOP_PLACEHOLDER_100=100, */ 14659d39f69dSJed Brown /* MATOP_PLACEHOLDER_101=101, */ 1466d519adbfSMatthew Knepley MATOP_CONJUGATE=102, 14679d39f69dSJed Brown /* MATOP_PLACEHOLDER_103=103, */ 1468d519adbfSMatthew Knepley MATOP_SET_VALUES_ROW=104, 1469d519adbfSMatthew Knepley MATOP_REAL_PART=105, 1470bda6c4cbSBarry Smith MATOP_IMAGINARY_PART=106, 1471bda6c4cbSBarry Smith MATOP_GET_ROW_UPPER_TRIANGULAR=107, 1472bda6c4cbSBarry Smith MATOP_RESTORE_ROW_UPPER_TRIANG=108, 1473bda6c4cbSBarry Smith MATOP_MAT_SOLVE=109, 14740e8d04f7SBarry Smith MATOP_GET_REDUNDANT_MATRIX=110, 1475d519adbfSMatthew Knepley MATOP_GET_ROW_MIN=111, 14760e8d04f7SBarry Smith MATOP_GET_COLUMN_VECTOR=112, 1477d519adbfSMatthew Knepley MATOP_MISSING_DIAGONAL=113, 14780e8d04f7SBarry Smith MATOP_GET_SEQ_NONZERO_STRUCTUR=114, 147989c6957cSBarry Smith MATOP_CREATE=115, 148089c6957cSBarry Smith MATOP_GET_GHOSTS=116, 14810e8d04f7SBarry Smith MATOP_GET_LOCAL_SUB_MATRIX=117, 14820e8d04f7SBarry Smith MATOP_RESTORE_LOCALSUB_MATRIX=118, 1483eeffb40dSHong Zhang MATOP_MULT_DIAGONAL_BLOCK=119, 14840e8d04f7SBarry Smith MATOP_HERMITIAN_TRANSPOSE=120, 14850e8d04f7SBarry Smith MATOP_MULT_HERMITIAN_TRANSPOSE=121, 14860e8d04f7SBarry Smith MATOP_MULT_HERMITIAN_TRANS_ADD=122, 14870e8d04f7SBarry Smith MATOP_GET_MULTI_PROC_BLOCK=123, 14889d39f69dSJed Brown MATOP_FIND_NONZERO_ROWS=124, 14890e8d04f7SBarry Smith MATOP_GET_COLUMN_NORMS=125, 14909d39f69dSJed Brown MATOP_INVERT_BLOCK_DIAGONAL=126, 14919d39f69dSJed Brown /* MATOP_PLACEHOLDER_127=127, */ 14920e8d04f7SBarry Smith MATOP_GET_SUB_MATRICES_PARALLE=128, 14932b8ad9a3SHong Zhang MATOP_SET_VALUES_BATCH=129, 14940e8d04f7SBarry Smith MATOP_TRANSPOSE_MAT_MULT=130, 14950e8d04f7SBarry Smith MATOP_TRANSPOSE_MAT_MULT_SYMBO=131, 1496e9b602ebSSatish Balay MATOP_TRANSPOSE_MAT_MULT_NUMER=132, 14970e8d04f7SBarry Smith MATOP_TRANSPOSE_COLORING_CREAT=133, 14980e8d04f7SBarry Smith MATOP_TRANS_COLORING_APPLY_SPT=134, 14990e8d04f7SBarry Smith MATOP_TRANS_COLORING_APPLY_DEN=135, 15000e8d04f7SBarry Smith MATOP_RART=136, 15010e8d04f7SBarry Smith MATOP_RART_SYMBOLIC=137, 15020e8d04f7SBarry Smith MATOP_RART_NUMERIC=138, 1503e09a3074SHong Zhang MATOP_SET_BLOCK_SIZES=139, 1504f9426fe0SMark Adams MATOP_AYPX=140, 15051919a2e2SJed Brown MATOP_RESIDUAL=141, 15069c8f2541SHong Zhang MATOP_FDCOLORING_SETUP=142, 15079c8f2541SHong Zhang MATOP_MPICONCATENATESEQ=144 1508fae171e0SBarry Smith } MatOperation; 1509014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatHasOperation(Mat,MatOperation,PetscBool *); 1510014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellSetOperation(Mat,MatOperation,void(*)(void)); 1511014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellGetOperation(Mat,MatOperation,void(**)(void)); 1512014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellSetContext(Mat,void*); 1513112a2221SBarry Smith 151490ace30eSBarry Smith /* 151590ace30eSBarry Smith Codes for matrices stored on disk. By default they are 151690ace30eSBarry Smith stored in a universal format. By changing the format with 15176a9046bcSBarry Smith PetscViewerPushFormat(viewer,PETSC_VIEWER_NATIVE); the matrices will 151890ace30eSBarry Smith be stored in a way natural for the matrix, for example dense matrices 151990ace30eSBarry Smith would be stored as dense. Matrices stored this way may only be 1520f4403165SShri Abhyankar read into matrices of the same type. 152190ace30eSBarry Smith */ 152290ace30eSBarry Smith #define MATRIX_BINARY_FORMAT_DENSE -1 152390ace30eSBarry Smith 1524014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetHashTableFactor(Mat,PetscReal); 1525014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatISGetLocalMat(Mat,Mat*); 1526014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatISSetLocalMat(Mat,Mat); 1527b7ce53b6SStefano Zampini PETSC_EXTERN PetscErrorCode MatISGetMPIXAIJ(Mat,MatReuse,Mat*); 15281f4e1ec7SBarry Smith 1529d9274352SBarry Smith /*S 1530d9274352SBarry Smith MatNullSpace - Object that removes a null space from a vector, i.e. 1531d9274352SBarry Smith orthogonalizes the vector to a subsapce 1532d9274352SBarry Smith 1533f7a9e4ceSBarry Smith Level: advanced 1534d9274352SBarry Smith 1535d9274352SBarry Smith Concepts: matrix; linear operator, null space 1536d9274352SBarry Smith 15376e1639daSBarry Smith Users manual sections: 15386e1639daSBarry Smith . sec_singular 15396e1639daSBarry Smith 1540d9274352SBarry Smith .seealso: MatNullSpaceCreate() 1541d9274352SBarry Smith S*/ 154274637425SBarry Smith typedef struct _p_MatNullSpace* MatNullSpace; 1543d9274352SBarry Smith 1544014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceCreate(MPI_Comm,PetscBool ,PetscInt,const Vec[],MatNullSpace*); 1545014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(MatNullSpace,Vec,void*),void*); 1546014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceDestroy(MatNullSpace*); 1547d0195637SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceRemove(MatNullSpace,Vec); 1548014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNullSpace(Mat, MatNullSpace *); 15495fa7ec2dSBarry Smith PETSC_EXTERN PetscErrorCode MatGetTransposeNullSpace(Mat, MatNullSpace *); 15505fa7ec2dSBarry Smith PETSC_EXTERN PetscErrorCode MatSetTransposeNullSpace(Mat,MatNullSpace); 1551014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNullSpace(Mat,MatNullSpace); 1552014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNearNullSpace(Mat,MatNullSpace); 1553014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNearNullSpace(Mat,MatNullSpace*); 1554014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceTest(MatNullSpace,Mat,PetscBool *); 1555014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceView(MatNullSpace,PetscViewer); 1556014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceGetVecs(MatNullSpace,PetscBool*,PetscInt*,const Vec**); 1557014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceCreateRigidBody(Vec,MatNullSpace*); 155874637425SBarry Smith 1559014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatReorderingSeqSBAIJ(Mat,IS); 1560014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetHashTableFactor(Mat,PetscReal); 1561014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetColumnIndices(Mat,PetscInt *); 1562014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJInvertBlockDiagonal(Mat); 15633f1d51d7SBarry Smith 1564014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMAIJ(Mat,PetscInt,Mat*); 1565014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMAIJRedimension(Mat,PetscInt,Mat*); 1566014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMAIJGetAIJ(Mat,Mat*); 1567c4f061fbSSatish Balay 1568014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatComputeExplicitOperator(Mat,Mat*); 1569b0a32e0cSBarry Smith 1570014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalScaleLocal(Mat,Vec); 157104f1ad80SBarry Smith 1572014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*); 1573014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetBase(Mat,Vec,Vec); 1574014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*); 1575014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*)); 1576014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec)); 1577014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt); 1578014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDResetHHistory(Mat); 1579014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctionError(Mat,PetscReal); 1580014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetPeriod(Mat,PetscInt); 1581014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDGetH(Mat,PetscScalar *); 1582014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetOptionsPrefix(Mat,const char[]); 1583014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*); 1584014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*); 1585e884886fSBarry Smith 15866370053bSBarry Smith /*S 15876370053bSBarry Smith MatMFFD - A data structured used to manage the computation of the h differencing parameter for matrix-free 15886370053bSBarry Smith Jacobian vector products 1589e884886fSBarry Smith 15906370053bSBarry Smith Notes: MATMFFD is a specific MatType which uses the MatMFFD data structure 15916370053bSBarry Smith 15926370053bSBarry Smith MatMFFD*() methods actually take the Mat as their first argument. Not a MatMFFD data structure 15936370053bSBarry Smith 15946370053bSBarry Smith Level: developer 15956370053bSBarry Smith 15966370053bSBarry Smith .seealso: MATMFFD, MatCreateMFFD(), MatMFFDSetFuction(), MatMFFDSetType(), MatMFFDRegister() 15976370053bSBarry Smith S*/ 1598e884886fSBarry Smith typedef struct _p_MatMFFD* MatMFFD; 1599e884886fSBarry Smith 160076bdecfbSBarry Smith /*J 1601e884886fSBarry Smith MatMFFDType - algorithm used to compute the h used in computing matrix-vector products via differencing of the function 1602e884886fSBarry Smith 1603e884886fSBarry Smith Level: beginner 1604e884886fSBarry Smith 1605e884886fSBarry Smith .seealso: MatMFFDSetType(), MatMFFDRegister() 160676bdecfbSBarry Smith J*/ 160719fd82e9SBarry Smith typedef const char* MatMFFDType; 1608e884886fSBarry Smith #define MATMFFD_DS "ds" 1609e884886fSBarry Smith #define MATMFFD_WP "wp" 1610e884886fSBarry Smith 161119fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDSetType(Mat,MatMFFDType); 1612bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDRegister(const char[],PetscErrorCode (*)(MatMFFD)); 1613e884886fSBarry Smith 1614014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDDSSetUmin(Mat,PetscReal); 1615014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDWPSetComputeNormU(Mat,PetscBool ); 1616e884886fSBarry Smith 161761ab5df0SBarry Smith PETSC_EXTERN PetscErrorCode MatFDColoringSetType(MatFDColoring,MatMFFDType); 161861ab5df0SBarry Smith 1619014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *); 1620014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *); 16217dbadf16SMatthew Knepley 162297969023SHong Zhang /* 162397969023SHong Zhang PETSc interface to MUMPS 162497969023SHong Zhang */ 162597969023SHong Zhang #ifdef PETSC_HAVE_MUMPS 1626014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMumpsSetIcntl(Mat,PetscInt,PetscInt); 1627bc6112feSHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetIcntl(Mat,PetscInt,PetscInt*); 1628b9e9ea26SSatish Balay PETSC_EXTERN PetscErrorCode MatMumpsSetCntl(Mat,PetscInt,PetscReal); 1629bc6112feSHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetCntl(Mat,PetscInt,PetscReal*); 1630bc6112feSHong Zhang 1631ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetInfo(Mat,PetscInt,PetscInt*); 1632ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetInfog(Mat,PetscInt,PetscInt*); 1633ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetRinfo(Mat,PetscInt,PetscReal*); 1634ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetRinfog(Mat,PetscInt,PetscReal*); 163597969023SHong Zhang #endif 163623a5497aSJed Brown 1637d954db57SHong Zhang /* 1638d615f992SSatish Balay PETSc interface to Mkl_Pardiso 1639b500a6b7SJose David Bermeol */ 1640b500a6b7SJose David Bermeol #ifdef PETSC_HAVE_MKL_PARDISO 1641d615f992SSatish Balay PETSC_EXTERN PetscErrorCode MatMkl_PardisoSetCntl(Mat,PetscInt,PetscInt); 1642b500a6b7SJose David Bermeol #endif 1643b500a6b7SJose David Bermeol 1644b500a6b7SJose David Bermeol /* 1645d305a81bSVasiliy Kozyrev PETSc interface to Mkl_CPardiso 1646d305a81bSVasiliy Kozyrev */ 1647d305a81bSVasiliy Kozyrev #ifdef PETSC_HAVE_MKL_CPARDISO 1648d305a81bSVasiliy Kozyrev PETSC_EXTERN PetscErrorCode MatMkl_CPardisoSetCntl(Mat,PetscInt,PetscInt); 1649d305a81bSVasiliy Kozyrev #endif 1650d305a81bSVasiliy Kozyrev 1651d305a81bSVasiliy Kozyrev /* 1652d954db57SHong Zhang PETSc interface to SUPERLU 1653d954db57SHong Zhang */ 1654d954db57SHong Zhang #ifdef PETSC_HAVE_SUPERLU 1655014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSuperluSetILUDropTol(Mat,PetscReal); 1656d954db57SHong Zhang #endif 1657d954db57SHong Zhang 1658fb785984SHong Zhang /* 1659fb785984SHong Zhang PETSc interface to SUPERLU_DIST 1660fb785984SHong Zhang */ 1661fb785984SHong Zhang #ifdef PETSC_HAVE_SUPERLU_DIST 1662fb785984SHong Zhang PETSC_EXTERN PetscErrorCode MatSuperluDistGetDiagU(Mat,PetscScalar*); 1663fb785984SHong Zhang #endif 1664fb785984SHong Zhang 1665575704cbSPieter Ghysels /* 1666575704cbSPieter Ghysels PETSc interface to STRUMPACK 1667575704cbSPieter Ghysels */ 1668575704cbSPieter Ghysels #ifdef PETSC_HAVE_STRUMPACK 1669575704cbSPieter Ghysels PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSRelCompTol(Mat,PetscReal); 1670575704cbSPieter Ghysels PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSMinSize(Mat,PetscInt); 1671575704cbSPieter Ghysels PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetColPerm(Mat,PetscBool); 1672575704cbSPieter Ghysels #endif 1673575704cbSPieter Ghysels 1674575704cbSPieter Ghysels 1675aa372e3fSPaul Mullowney #ifdef PETSC_HAVE_CUDA 16763f9c0db1SPaul Mullowney /*E 1677e057df02SPaul Mullowney MatCUSPARSEStorageFormat - indicates the storage format for CUSPARSE (GPU) 16782692e278SPaul Mullowney matrices. 1679e057df02SPaul Mullowney 1680e057df02SPaul Mullowney Not Collective 1681e057df02SPaul Mullowney 16828468deeeSKarl Rupp + MAT_CUSPARSE_CSR - Compressed Sparse Row 16832692e278SPaul Mullowney . MAT_CUSPARSE_ELL - Ellpack (requires CUDA 4.2 or later). 16842692e278SPaul Mullowney - MAT_CUSPARSE_HYB - Hybrid, a combination of Ellpack and Coordinate format (requires CUDA 4.2 or later). 1685e057df02SPaul Mullowney 1686e057df02SPaul Mullowney Level: intermediate 1687e057df02SPaul Mullowney 1688af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 1689e057df02SPaul Mullowney 1690e057df02SPaul Mullowney .seealso: MatCUSPARSESetFormat(), MatCUSPARSEFormatOperation 1691e057df02SPaul Mullowney E*/ 1692e057df02SPaul Mullowney 16933f9c0db1SPaul Mullowney typedef enum {MAT_CUSPARSE_CSR, MAT_CUSPARSE_ELL, MAT_CUSPARSE_HYB} MatCUSPARSEStorageFormat; 1694e057df02SPaul Mullowney 1695e057df02SPaul Mullowney /* these will be strings associated with enumerated type defined above */ 1696e057df02SPaul Mullowney PETSC_EXTERN const char *const MatCUSPARSEStorageFormats[]; 1697e057df02SPaul Mullowney 16983f9c0db1SPaul Mullowney /*E 1699e057df02SPaul Mullowney MatCUSPARSEFormatOperation - indicates the operation of CUSPARSE (GPU) 17002692e278SPaul Mullowney matrices whose operation should use a particular storage format. 1701e057df02SPaul Mullowney 1702e057df02SPaul Mullowney Not Collective 1703e057df02SPaul Mullowney 17048468deeeSKarl Rupp + MAT_CUSPARSE_MULT_DIAG - sets the storage format for the diagonal matrix in the parallel MatMult 17058468deeeSKarl Rupp . MAT_CUSPARSE_MULT_OFFDIAG - sets the storage format for the offdiagonal matrix in the parallel MatMult 17068468deeeSKarl Rupp . MAT_CUSPARSE_MULT - sets the storage format for the entire matrix in the serial (single GPU) MatMult 17078468deeeSKarl Rupp - MAT_CUSPARSE_ALL - sets the storage format for all CUSPARSE (GPU) matrices 1708e057df02SPaul Mullowney 1709e057df02SPaul Mullowney Level: intermediate 1710e057df02SPaul Mullowney 1711e057df02SPaul Mullowney .seealso: MatCUSPARSESetFormat(), MatCUSPARSEStorageFormat 1712e057df02SPaul Mullowney E*/ 171336d62e41SPaul Mullowney typedef enum {MAT_CUSPARSE_MULT_DIAG, MAT_CUSPARSE_MULT_OFFDIAG, MAT_CUSPARSE_MULT, MAT_CUSPARSE_ALL} MatCUSPARSEFormatOperation; 1714e057df02SPaul Mullowney 171510e9db80SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 171610e9db80SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCreateAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 1717e057df02SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCUSPARSESetFormat(Mat,MatCUSPARSEFormatOperation,MatCUSPARSEStorageFormat); 17189ae82921SPaul Mullowney #endif 17199ae82921SPaul Mullowney 172090c53902SBarry Smith #if defined(PETSC_HAVE_CUSP) 1721014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCUSP(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 1722014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateAIJCUSP(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 1723e057df02SPaul Mullowney 17243f9c0db1SPaul Mullowney /*E 1725e057df02SPaul Mullowney MatCUSPStorageFormat - indicates the storage format for CUSP (GPU) 172636d62e41SPaul Mullowney matrices. 1727e057df02SPaul Mullowney 1728e057df02SPaul Mullowney Not Collective 1729e057df02SPaul Mullowney 17308468deeeSKarl Rupp + MAT_CUSP_CSR - Compressed Sparse Row 17318468deeeSKarl Rupp . MAT_CUSP_DIA - Diagonal 17328468deeeSKarl Rupp - MAT_CUSP_ELL - Ellpack 1733e057df02SPaul Mullowney 1734e057df02SPaul Mullowney Level: intermediate 1735e057df02SPaul Mullowney 1736af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 1737e057df02SPaul Mullowney 1738e057df02SPaul Mullowney .seealso: MatCUSPSetFormat(), MatCUSPFormatOperation 1739e057df02SPaul Mullowney E*/ 17403f9c0db1SPaul Mullowney typedef enum {MAT_CUSP_CSR, MAT_CUSP_DIA, MAT_CUSP_ELL} MatCUSPStorageFormat; 1741e057df02SPaul Mullowney 1742e057df02SPaul Mullowney /* these will be strings associated with enumerated type defined above */ 1743e057df02SPaul Mullowney PETSC_EXTERN const char *const MatCUSPStorageFormats[]; 1744e057df02SPaul Mullowney 17453f9c0db1SPaul Mullowney /*E 1746e057df02SPaul Mullowney MatCUSPFormatOperation - indicates the operation of CUSP (GPU) 174736d62e41SPaul Mullowney matrices whose operation should use a particular storage format. 1748e057df02SPaul Mullowney 1749e057df02SPaul Mullowney Not Collective 1750e057df02SPaul Mullowney 17518468deeeSKarl Rupp + MAT_CUSP_MULT_DIAG - sets the storage format for the diagonal matrix in the parallel MatMult 17528468deeeSKarl Rupp . MAT_CUSP_MULT_OFFDIAG - sets the storage format for the offdiagonal matrix in the parallel MatMult 17538468deeeSKarl Rupp . MAT_CUSP_MULT - sets the storage format for the entire matrix in the serial (single GPU) MatMult 17548468deeeSKarl Rupp - MAT_CUSP_ALL - sets the storage format for all CUSP (GPU) matrices 1755e057df02SPaul Mullowney 1756e057df02SPaul Mullowney Level: intermediate 1757e057df02SPaul Mullowney 1758af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 1759e057df02SPaul Mullowney 1760e057df02SPaul Mullowney .seealso: MatCUSPSetFormat(), MatCUSPStorageFormat 1761e057df02SPaul Mullowney E*/ 176236d62e41SPaul Mullowney typedef enum {MAT_CUSP_MULT_DIAG, MAT_CUSP_MULT_OFFDIAG, MAT_CUSP_MULT, MAT_CUSP_ALL} MatCUSPFormatOperation; 1763e057df02SPaul Mullowney 1764e057df02SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCUSPSetFormat(Mat,MatCUSPFormatOperation,MatCUSPStorageFormat); 1765e057df02SPaul Mullowney #endif 176690c53902SBarry Smith 1767d67ff14aSKarl Rupp #if defined(PETSC_HAVE_VIENNACL) 1768d67ff14aSKarl Rupp PETSC_EXTERN PetscErrorCode MatCreateSeqAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 1769d67ff14aSKarl Rupp PETSC_EXTERN PetscErrorCode MatCreateAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 1770d67ff14aSKarl Rupp #endif 1771d67ff14aSKarl Rupp 177254efbe56SHong Zhang /* 177354efbe56SHong Zhang PETSc interface to FFTW 177454efbe56SHong Zhang */ 177554efbe56SHong Zhang #if defined(PETSC_HAVE_FFTW) 1776014dd563SJed Brown PETSC_EXTERN PetscErrorCode VecScatterPetscToFFTW(Mat,Vec,Vec); 1777014dd563SJed Brown PETSC_EXTERN PetscErrorCode VecScatterFFTWToPetsc(Mat,Vec,Vec); 17782a7a6963SBarry Smith PETSC_EXTERN PetscErrorCode MatCreateVecsFFTW(Mat,Vec*,Vec*,Vec*); 177954efbe56SHong Zhang #endif 178054efbe56SHong Zhang 1781014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateNest(MPI_Comm,PetscInt,const IS[],PetscInt,const IS[],const Mat[],Mat*); 1782014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSize(Mat,PetscInt*,PetscInt*); 1783014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetISs(Mat,IS[],IS[]); 1784014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetLocalISs(Mat,IS[],IS[]); 1785014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSubMats(Mat,PetscInt*,PetscInt*,Mat***); 1786014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSubMat(Mat,PetscInt,PetscInt,Mat*); 178719fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatNestSetVecType(Mat,VecType); 1788014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestSetSubMats(Mat,PetscInt,const IS[],PetscInt,const IS[],const Mat[]); 1789014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestSetSubMat(Mat,PetscInt,PetscInt,Mat); 1790d8588912SDave May 17914325cce7SMatthew G Knepley PETSC_EXTERN PetscErrorCode MatChop(Mat,PetscReal); 1792e0a58c10SMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatComputeBandwidth(Mat,PetscReal,PetscInt*); 17934325cce7SMatthew G Knepley 1794b541e6a4SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSubdomainsCreateCoalesce(Mat,PetscInt,PetscInt*,IS**); 179553dd7562SDmitry Karpeev 1796c094ef40SMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatPreallocatorPreallocate(Mat,PetscBool,Mat); 1797c094ef40SMatthew G. Knepley 1798539c167fSBarry Smith PETSC_INTERN PetscErrorCode MatHeaderMerge(Mat,Mat*); 1799539c167fSBarry Smith PETSC_EXTERN PetscErrorCode MatHeaderReplace(Mat,Mat*); 1800539c167fSBarry Smith 180123a5497aSJed Brown #endif 1802