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" 84773941d6SBarry Smith 8576bdecfbSBarry Smith /*J 86c7393fdbSBarry Smith MatSolverPackage - String with the name of a PETSc matrix solver type. 879989ab13SBarry Smith 88c2b89b5dSBarry Smith For example: "petsc" indicates what PETSc provides, "superlu_dist" the parallel SuperLU_DIST package etc 899989ab13SBarry Smith 909989ab13SBarry Smith Level: beginner 919989ab13SBarry Smith 925c9eb25fSBarry Smith .seealso: MatGetFactor(), Mat, MatSetType(), MatType 9376bdecfbSBarry Smith J*/ 94c7393fdbSBarry Smith #define MatSolverPackage char* 952692d6eeSBarry Smith #define MATSOLVERSUPERLU "superlu" 962692d6eeSBarry Smith #define MATSOLVERSUPERLU_DIST "superlu_dist" 9708f5efcfSPieter Ghysels #define MATSOLVERSTRUMPACK "strumpack" 982692d6eeSBarry Smith #define MATSOLVERUMFPACK "umfpack" 9920db9a53SJed Brown #define MATSOLVERCHOLMOD "cholmod" 100d89f5e7aSBarry Smith #define MATSOLVERKLU "klu" 101ecddaf3cSBarry Smith #define MATSOLVERCLIQUE "clique" 102d89f5e7aSBarry Smith #define MATSOLVERELEMENTAL "elemental" 1032692d6eeSBarry Smith #define MATSOLVERESSL "essl" 1042692d6eeSBarry Smith #define MATSOLVERLUSOL "lusol" 1052692d6eeSBarry Smith #define MATSOLVERMUMPS "mumps" 106d615f992SSatish Balay #define MATSOLVERMKL_PARDISO "mkl_pardiso" 107d305a81bSVasiliy Kozyrev #define MATSOLVERMKL_CPARDISO "mkl_cpardiso" 1082692d6eeSBarry Smith #define MATSOLVERPASTIX "pastix" 1092692d6eeSBarry Smith #define MATSOLVERMATLAB "matlab" 1102692d6eeSBarry Smith #define MATSOLVERPETSC "petsc" 1112692d6eeSBarry Smith #define MATSOLVERBAS "bas" 1129ae82921SPaul Mullowney #define MATSOLVERCUSPARSE "cusparse" 113c0cdd4a1SDahai Guo 114b24902e0SBarry Smith /*E 115b24902e0SBarry Smith MatFactorType - indicates what type of factorization is requested 116b24902e0SBarry Smith 117b24902e0SBarry Smith Level: beginner 118b24902e0SBarry Smith 119af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 120b24902e0SBarry Smith 121c7393fdbSBarry Smith .seealso: MatSolverPackage, MatGetFactor() 122b24902e0SBarry Smith E*/ 123599ef60dSHong Zhang typedef enum {MAT_FACTOR_NONE, MAT_FACTOR_LU, MAT_FACTOR_CHOLESKY, MAT_FACTOR_ILU, MAT_FACTOR_ICC,MAT_FACTOR_ILUDT} MatFactorType; 124014dd563SJed Brown PETSC_EXTERN const char *const MatFactorTypes[]; 125e92e720dSBarry Smith 126014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactor(Mat,const MatSolverPackage,MatFactorType,Mat*); 127014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactorAvailable(Mat,const MatSolverPackage,MatFactorType,PetscBool *); 128014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFactorGetSolverPackage(Mat,const MatSolverPackage*); 129014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactorType(Mat,MatFactorType*); 13018713533SBarry Smith PETSC_EXTERN PetscErrorCode MatSolverPackageRegister(const MatSolverPackage,const MatType,MatFactorType,PetscErrorCode(*)(Mat,MatFactorType,Mat*)); 13142c9c57cSBarry Smith PETSC_EXTERN PetscErrorCode MatSolverPackageGet(const MatSolverPackage,const MatType,MatFactorType,PetscBool*,PetscBool*,PetscErrorCode (**)(Mat,MatFactorType,Mat*)); 1329989ab13SBarry Smith 133c06d978dSMatthew Knepley /* Logging support */ 1340700a824SBarry Smith #define MAT_FILE_CLASSID 1211216 /* used to indicate matrices in binary files */ 135014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_CLASSID; 136335efc43SPeter Brune PETSC_EXTERN PetscClassId MAT_COLORING_CLASSID; 137014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_FDCOLORING_CLASSID; 138014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_TRANSPOSECOLORING_CLASSID; 139014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_PARTITIONING_CLASSID; 140014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_COARSEN_CLASSID; 141014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_NULLSPACE_CLASSID; 142014dd563SJed Brown PETSC_EXTERN PetscClassId MATMFFD_CLASSID; 143c06d978dSMatthew Knepley 144ceb03754SKris Buschelman /*E 145*cf37664fSBarry Smith MatReuse - Indicates if matrices obtained from a previous call to MatGetSubMatrices(), MatGetSubMatrix(), MatConvert() or several other functions 146*cf37664fSBarry Smith are to be reused to store the new matrix values. 147*cf37664fSBarry Smith 148*cf37664fSBarry Smith $ MAT_INITIAL_MATRIX - create a new matrix 149*cf37664fSBarry Smith $ MAT_REUSE_MATRIX - reuse the matrix created with a previous call that used MAT_INITIAL_MATRIX 150*cf37664fSBarry Smith $ MAT_INPLACE_MATRIX - replace the first input matrix with the new matrix (not applicable to all functions) 151*cf37664fSBarry 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) 152ceb03754SKris Buschelman 153ceb03754SKris Buschelman Level: beginner 154ceb03754SKris Buschelman 155af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 156ceb03754SKris Buschelman 1570c451bc4SBarry Smith .seealso: MatGetSubMatrices(), MatGetSubMatrix(), MatDestroyMatrices(), MatConvert() 158ceb03754SKris Buschelman E*/ 159511c6705SHong Zhang typedef enum {MAT_INITIAL_MATRIX,MAT_REUSE_MATRIX,MAT_IGNORE_MATRIX,MAT_INPLACE_MATRIX} MatReuse; 160ceb03754SKris Buschelman 1615494a064SHong Zhang /*E 1625494a064SHong Zhang MatGetSubMatrixOption - Indicates if matrices obtained from a call to MatGetSubMatrices() 163829201f2SHong Zhang include the matrix values. Currently it is only used by MatGetSeqNonzerostructure(). 1645494a064SHong Zhang 1655494a064SHong Zhang Level: beginner 1665494a064SHong Zhang 167829201f2SHong Zhang .seealso: MatGetSeqNonzerostructure() 1685494a064SHong Zhang E*/ 1695494a064SHong Zhang typedef enum {MAT_DO_NOT_GET_VALUES,MAT_GET_VALUES} MatGetSubMatrixOption; 1705494a064SHong Zhang 171607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatInitializePackage(void); 172c06d978dSMatthew Knepley 173014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreate(MPI_Comm,Mat*); 174014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetSizes(Mat,PetscInt,PetscInt,PetscInt,PetscInt); 17519fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatSetType(Mat,MatType); 176014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetFromOptions(Mat); 177685405a1SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatViewFromOptions(Mat A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);} 178bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatRegister(const char[],PetscErrorCode(*)(Mat)); 179014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRegisterBaseName(const char[],const char[],const char[]); 180014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetOptionsPrefix(Mat,const char[]); 181014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAppendOptionsPrefix(Mat,const char[]); 182014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOptionsPrefix(Mat,const char*[]); 18384d44b13SHong Zhang PETSC_EXTERN PetscErrorCode MatSetErrorIfFailure(Mat,PetscBool); 184f69a0ea3SMatthew Knepley 185140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatList; 186140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatColoringList; 187140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatPartitioningList; 18828988994SBarry Smith 1893b224e63SBarry Smith /*E 190aa6c7ce3SBarry Smith MatStructure - Indicates if two matrices have the same nonzero structure 1913b224e63SBarry Smith 1923b224e63SBarry Smith Level: beginner 1933b224e63SBarry Smith 194af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 1953b224e63SBarry Smith 196aa6c7ce3SBarry Smith .seealso: MatCopy(), MatAXPY() 1973b224e63SBarry Smith E*/ 198aa6c7ce3SBarry Smith typedef enum {DIFFERENT_NONZERO_PATTERN,SUBSET_NONZERO_PATTERN,SAME_NONZERO_PATTERN} MatStructure; 1993b224e63SBarry Smith 200014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqDense(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*); 201014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateDense(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*); 202014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 203014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 204014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *); 205014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithSplitArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],PetscInt[],PetscInt[],PetscScalar[],Mat*); 2068d7a6e47SBarry Smith 207014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 208014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 209014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat*); 210d21a29f3SJed Brown 211014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAdj(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscInt[],Mat*); 212014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 213d21a29f3SJed Brown 214014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 215014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPISBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *); 21638f409ebSLisandro Dalcin PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 217014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 2184cce697cSJed Brown PETSC_EXTERN PetscErrorCode MatXAIJSetPreallocation(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscInt[],const PetscInt[]); 219dfb205c3SBarry Smith 220014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateShell(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,void *,Mat*); 221014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateNormal(Mat,Mat*); 222c5e4d11fSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatCreateNormalHermitian(Mat,Mat*); 223986c4d72SJose E. Roman PETSC_EXTERN PetscErrorCode MatCreateLRC(Mat,Mat,Vec,Mat,Mat*); 224267ca84cSJose E. Roman PETSC_EXTERN PetscErrorCode MatLRCGetMats(Mat,Mat*,Mat*,Vec*,Mat*); 225c5e4d11fSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatCreateIS(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,ISLocalToGlobalMapping,ISLocalToGlobalMapping,Mat*); 226014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 227014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 228c0cdd4a1SDahai Guo 229014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateScatter(MPI_Comm,VecScatter,Mat*); 230014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScatterSetVecScatter(Mat,VecScatter); 231014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScatterGetVecScatter(Mat,VecScatter*); 232014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateBlockMat(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,Mat*); 233014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeAddMat(Mat,Mat); 234014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeMerge(Mat); 235014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateComposite(MPI_Comm,PetscInt,const Mat*,Mat*); 2366d7c1e57SBarry Smith typedef enum {MAT_COMPOSITE_ADDITIVE,MAT_COMPOSITE_MULTIPLICATIVE} MatCompositeType; 237014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeSetType(Mat,MatCompositeType); 2386d7c1e57SBarry Smith 23919fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatCreateFFT(MPI_Comm,PetscInt,const PetscInt[],MatType,Mat*); 240014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqCUFFT(MPI_Comm,PetscInt,const PetscInt[],Mat*); 241dedccee8SHong Zhang 242014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateTranspose(Mat,Mat*); 2438060fb66Sstefano_zampini PETSC_EXTERN PetscErrorCode MatTransposeGetMat(Mat,Mat*); 244d0de2241SAndrew Spott PETSC_EXTERN PetscErrorCode MatCreateHermitianTranspose(Mat,Mat*); 245014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSubMatrix(Mat,IS,IS,Mat*); 246014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSubMatrixUpdate(Mat,Mat,IS,IS); 247014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateLocalRef(Mat,IS,IS,Mat*); 2481472f72bSBarry Smith 249978814f1SStefano Zampini #if defined(PETSC_HAVE_HYPRE) 250d975228cSstefano_zampini PETSC_EXTERN PetscErrorCode MatHYPRESetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 251978814f1SStefano Zampini #endif 252978814f1SStefano Zampini 253014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPythonSetType(Mat,const char[]); 2541d6018f0SLisandro Dalcin 255014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetUp(Mat); 256014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroy(Mat*); 257e56f5c9eSBarry Smith PETSC_EXTERN PetscErrorCode MatGetNonzeroState(Mat,PetscObjectState*); 25821c89e3eSBarry Smith 259014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatConjugate(Mat); 260014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRealPart(Mat); 261014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatImaginaryPart(Mat); 262014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetDiagonalBlock(Mat,Mat*); 263014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetTrace(Mat,PetscScalar*); 264713ccfa9SJed Brown PETSC_EXTERN PetscErrorCode MatInvertBlockDiagonal(Mat,const PetscScalar **); 26599cafbc1SBarry Smith 2668ed539a5SBarry Smith /* ------------------------------------------------------------*/ 267014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 268014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 269014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesRow(Mat,PetscInt,const PetscScalar[]); 270014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesRowLocal(Mat,PetscInt,const PetscScalar[]); 271014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBatch(Mat,PetscInt,PetscInt,PetscInt[],const PetscScalar[]); 27273a71a0fSBarry Smith PETSC_EXTERN PetscErrorCode MatSetRandom(Mat,PetscRandom); 27384cb2905SBarry Smith 2742ef4de8bSBarry Smith /*S 2752ef4de8bSBarry Smith MatStencil - Data structure (C struct) for storing information about a single row or 27662ef0227SBarry Smith column of a matrix as indexed on an associated grid. These are arguments to MatSetStencil() and MatSetBlockStencil() 27762ef0227SBarry Smith 27862ef0227SBarry 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). 27962ef0227SBarry 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. 28062ef0227SBarry Smith 28162ef0227SBarry Smith For stencil access to vectors see DMDAVecGetArray(), DMDAVecGetArrayF90(). 282be479b30SJed Brown 283be479b30SJed Brown Fortran usage is different, see MatSetValuesStencil() for details. 2842ef4de8bSBarry Smith 2852ef4de8bSBarry Smith Level: beginner 2862ef4de8bSBarry Smith 2872ef4de8bSBarry Smith Concepts: matrix; linear operator 2882ef4de8bSBarry Smith 28962ef0227SBarry Smith .seealso: MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockedStencil(), DMDAVecGetArray(), DMDAVecGetArrayF90() 2902ef4de8bSBarry Smith S*/ 291435da068SBarry Smith typedef struct { 292c1ac3661SBarry Smith PetscInt k,j,i,c; 293435da068SBarry Smith } MatStencil; 2942ef4de8bSBarry Smith 295014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode); 296014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode); 297014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt); 298435da068SBarry Smith 299d91e6319SBarry Smith /*E 300d91e6319SBarry Smith MatAssemblyType - Indicates if the matrix is now to be used, or if you plan 30162ef0227SBarry Smith to continue to add or insert values to it 302d91e6319SBarry Smith 303d91e6319SBarry Smith Level: beginner 304d91e6319SBarry Smith 305d91e6319SBarry Smith .seealso: MatAssemblyBegin(), MatAssemblyEnd() 306d91e6319SBarry Smith E*/ 3076d4a8577SBarry Smith typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType; 308014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssemblyBegin(Mat,MatAssemblyType); 309014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssemblyEnd(Mat,MatAssemblyType); 310014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssembled(Mat,PetscBool *); 3114f9c727eSBarry Smith 3121d73ed98SBarry Smith 31330de9b25SBarry Smith 314d91e6319SBarry Smith /*E 315d91e6319SBarry Smith MatOption - Options that may be set for a matrix and its behavior or storage 316d91e6319SBarry Smith 317d91e6319SBarry Smith Level: beginner 318d91e6319SBarry Smith 319af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 320d91e6319SBarry Smith 321382164b0SBarry Smith Developer Notes: Entries that are negative need not be called collectively by all processes. 322382164b0SBarry Smith 323d91e6319SBarry Smith .seealso: MatSetOption() 324d91e6319SBarry Smith E*/ 325c5e4d11fSDmitry Karpeev typedef enum {MAT_OPTION_MIN = -3, 326c5e4d11fSDmitry Karpeev MAT_UNUSED_NONZERO_LOCATION_ERR = -2, 32792d4d306SBarry Smith MAT_ROW_ORIENTED = -1, 328382164b0SBarry Smith MAT_SYMMETRIC = 1, 329382164b0SBarry Smith MAT_STRUCTURALLY_SYMMETRIC = 2, 330382164b0SBarry Smith MAT_NEW_DIAGONALS = 3, 331382164b0SBarry Smith MAT_IGNORE_OFF_PROC_ENTRIES = 4, 332382164b0SBarry Smith MAT_USE_HASH_TABLE = 5, 333382164b0SBarry Smith MAT_KEEP_NONZERO_PATTERN = 6, 334382164b0SBarry Smith MAT_IGNORE_ZERO_ENTRIES = 7, 335382164b0SBarry Smith MAT_USE_INODES = 8, 336382164b0SBarry Smith MAT_HERMITIAN = 9, 337382164b0SBarry Smith MAT_SYMMETRY_ETERNAL = 10, 338c5e4d11fSDmitry Karpeev MAT_NEW_NONZERO_LOCATION_ERR = 11, 339382164b0SBarry Smith MAT_IGNORE_LOWER_TRIANGULAR = 12, 340382164b0SBarry Smith MAT_ERROR_LOWER_TRIANGULAR = 13, 341382164b0SBarry Smith MAT_GETROW_UPPERTRIANGULAR = 14, 342382164b0SBarry Smith MAT_SPD = 15, 34392d4d306SBarry Smith MAT_NO_OFF_PROC_ZERO_ROWS = 16, 34492d4d306SBarry Smith MAT_NO_OFF_PROC_ENTRIES = 17, 34592d4d306SBarry Smith MAT_NEW_NONZERO_LOCATIONS = 18, 346c5e4d11fSDmitry Karpeev MAT_NEW_NONZERO_ALLOCATION_ERR = 19, 347c5e4d11fSDmitry Karpeev MAT_SUBSET_OFF_PROC_ENTRIES = 20, 348c10200c1SHong Zhang MAT_SUBMAT_SINGLEIS = 21, 349c10200c1SHong Zhang MAT_OPTION_MAX = 22} MatOption; 350e057df02SPaul Mullowney 351014dd563SJed Brown PETSC_EXTERN const char *MatOptions[]; 352014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetOption(Mat,MatOption,PetscBool); 3537d68702bSBarry Smith PETSC_EXTERN PetscErrorCode MatGetOption(Mat,MatOption,PetscBool*); 35419fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatGetType(Mat,MatType*); 35584cb2905SBarry Smith 356014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]); 357014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 358014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 359014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowUpperTriangular(Mat); 360014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRowUpperTriangular(Mat); 361014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnVector(Mat,Vec,PetscInt); 3628c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJGetArray(Mat,PetscScalar *[]); 3638c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJRestoreArray(Mat,PetscScalar *[]); 36421e72a00SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJGetMaxRowNonzeros(Mat,PetscInt*); 365bd04181cSBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJSetValuesLocalFast(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 3668c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatDenseGetArray(Mat,PetscScalar *[]); 3678c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatDenseRestoreArray(Mat,PetscScalar *[]); 368014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBlockSize(Mat,PetscInt *); 369014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSize(Mat,PetscInt); 370014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBlockSizes(Mat,PetscInt *,PetscInt *); 371014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSizes(Mat,PetscInt,PetscInt); 37233d57670SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSizesFromMats(Mat,Mat,Mat); 3731620fd73SBarry Smith 374014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMult(Mat,Vec,Vec); 375014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultDiagonalBlock(Mat,Vec,Vec); 376014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultAdd(Mat,Vec,Vec,Vec); 377014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTranspose(Mat,Vec,Vec); 378014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultHermitianTranspose(Mat,Vec,Vec); 379014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsTranspose(Mat,Mat,PetscReal,PetscBool *); 380014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitianTranspose(Mat,Mat,PetscReal,PetscBool *); 381014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeAdd(Mat,Vec,Vec,Vec); 382014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultHermitianTransposeAdd(Mat,Vec,Vec,Vec); 383014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultConstrained(Mat,Vec,Vec); 384014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeConstrained(Mat,Vec,Vec); 385014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatSolve(Mat,Mat,Mat); 386f9426fe0SMark Adams PETSC_EXTERN PetscErrorCode MatResidual(Mat,Vec,Vec,Vec); 387f5edf698SHong Zhang 388d91e6319SBarry Smith /*E 389d91e6319SBarry Smith MatDuplicateOption - Indicates if a duplicated sparse matrix should have 390d91e6319SBarry Smith its numerical values copied over or just its nonzero structure. 391d91e6319SBarry Smith 392d91e6319SBarry Smith Level: beginner 393d91e6319SBarry Smith 394af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 395d91e6319SBarry Smith 39670dcbbb9SBarry Smith $ MAT_SHARE_NONZERO_PATTERN - the i and j arrays in the new matrix will be shared with the original matrix 39770dcbbb9SBarry Smith $ this also triggers the MAT_DO_NOT_COPY_VALUES option. This is used when you 39870dcbbb9SBarry Smith $ have several matrices with the same nonzero pattern. 39970dcbbb9SBarry Smith 400d91e6319SBarry Smith .seealso: MatDuplicate() 401d91e6319SBarry Smith E*/ 40270dcbbb9SBarry Smith typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES,MAT_SHARE_NONZERO_PATTERN} MatDuplicateOption; 4032e8a6d31SBarry Smith 40419fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatConvert(Mat,MatType,MatReuse,Mat*); 405014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDuplicate(Mat,MatDuplicateOption,Mat*); 40694a9d846SBarry Smith 407cb5b572fSBarry Smith 408014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCopy(Mat,Mat,MatStructure); 409014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatView(Mat,PetscViewer); 410014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsSymmetric(Mat,PetscReal,PetscBool *); 411014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsStructurallySymmetric(Mat,PetscBool *); 412014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitian(Mat,PetscReal,PetscBool *); 413014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsSymmetricKnown(Mat,PetscBool *,PetscBool *); 414014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitianKnown(Mat,PetscBool *,PetscBool *); 415014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMissingDiagonal(Mat,PetscBool *,PetscInt *); 416014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLoad(Mat, PetscViewer); 4177b80b807SBarry Smith 4181a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool *); 4191a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool *); 4201a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool *); 4211a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool *); 422d4fbbf0eSBarry Smith 423d91e6319SBarry Smith /*S 424d91e6319SBarry Smith MatInfo - Context of matrix information, used with MatGetInfo() 425d91e6319SBarry Smith 426d91e6319SBarry Smith In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE 427d91e6319SBarry Smith 428d91e6319SBarry Smith Level: intermediate 429d91e6319SBarry Smith 430d91e6319SBarry Smith Concepts: matrix^nonzero information 431d91e6319SBarry Smith 432d9274352SBarry Smith .seealso: MatGetInfo(), MatInfoType 433d91e6319SBarry Smith S*/ 4344e220ebcSLois Curfman McInnes typedef struct { 435b0a32e0cSBarry Smith PetscLogDouble block_size; /* block size */ 436b0a32e0cSBarry Smith PetscLogDouble nz_allocated,nz_used,nz_unneeded; /* number of nonzeros */ 437b0a32e0cSBarry Smith PetscLogDouble memory; /* memory allocated */ 438b0a32e0cSBarry Smith PetscLogDouble assemblies; /* number of matrix assemblies called */ 439b0a32e0cSBarry Smith PetscLogDouble mallocs; /* number of mallocs during MatSetValues() */ 440b0a32e0cSBarry Smith PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */ 441b0a32e0cSBarry Smith PetscLogDouble factor_mallocs; /* number of mallocs during factorization */ 4424e220ebcSLois Curfman McInnes } MatInfo; 4434e220ebcSLois Curfman McInnes 444d9274352SBarry Smith /*E 445d9274352SBarry Smith MatInfoType - Indicates if you want information about the local part of the matrix, 446d9274352SBarry Smith the entire parallel matrix or the maximum over all the local parts. 447d9274352SBarry Smith 448d9274352SBarry Smith Level: beginner 449d9274352SBarry Smith 450af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 451d9274352SBarry Smith 452d9274352SBarry Smith .seealso: MatGetInfo(), MatInfo 453d9274352SBarry Smith E*/ 4547b80b807SBarry Smith typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType; 455014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetInfo(Mat,MatInfoType,MatInfo*); 456014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetDiagonal(Mat,Vec); 457014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMax(Mat,Vec,PetscInt[]); 458014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMin(Mat,Vec,PetscInt[]); 459014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMaxAbs(Mat,Vec,PetscInt[]); 460014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMinAbs(Mat,Vec,PetscInt[]); 461014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowSum(Mat,Vec); 462014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTranspose(Mat,MatReuse,Mat*); 463014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatHermitianTranspose(Mat,MatReuse,Mat*); 464014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPermute(Mat,IS,IS,Mat*); 465014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalScale(Mat,Vec,Vec); 466014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalSet(Mat,Vec,InsertMode); 467a52676ecSHong Zhang 468014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatEqual(Mat,Mat,PetscBool*); 469014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultEqual(Mat,Mat,PetscInt,PetscBool*); 470014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultAddEqual(Mat,Mat,PetscInt,PetscBool*); 471014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeEqual(Mat,Mat,PetscInt,PetscBool*); 472014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscBool*); 473a52676ecSHong Zhang PETSC_EXTERN PetscErrorCode MatMatMultEqual(Mat,Mat,Mat,PetscInt,PetscBool*); 474a52676ecSHong Zhang PETSC_EXTERN PetscErrorCode MatTransposeMatMultEqual(Mat,Mat,Mat,PetscInt,PetscBool*); 4757b80b807SBarry Smith 476014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNorm(Mat,NormType,PetscReal*); 477014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnNorms(Mat,NormType,PetscReal*); 478014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroEntries(Mat); 479014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 480014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsIS(Mat,IS,PetscScalar,Vec,Vec); 481014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsStencil(Mat,PetscInt,const MatStencil [],PetscScalar,Vec,Vec); 482014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsStencil(Mat,PetscInt,const MatStencil[],PetscScalar,Vec,Vec); 483014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumns(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 484014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsIS(Mat,IS,PetscScalar,Vec,Vec); 4857b80b807SBarry Smith 486014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSize(Mat,PetscInt*,PetscInt*); 487014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalSize(Mat,PetscInt*,PetscInt*); 488014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRange(Mat,PetscInt*,PetscInt*); 489014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRanges(Mat,const PetscInt**); 490014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRangeColumn(Mat,PetscInt*,PetscInt*); 491014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRangesColumn(Mat,const PetscInt**); 492566876eaSJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipIS(Mat,IS*,IS*); 4935ef9f2a5SBarry Smith 494014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]); 49553dd7562SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatGetSubMatricesMPI(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]); 496014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroyMatrices(PetscInt,Mat *[]); 497014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSubMatrix(Mat,IS,IS,MatReuse,Mat *); 498014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalSubMatrix(Mat,IS,IS,Mat*); 499014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreLocalSubMatrix(Mat,IS,IS,Mat*); 500014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSeqNonzeroStructure(Mat,Mat*); 501014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroySeqNonzeroStructure(Mat*); 5027b80b807SBarry Smith 503014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJ(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*); 504014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJSymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*); 505014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJNumeric(Mat,Mat); 506014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMat(Mat,MatReuse,Mat*); 507014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*); 508014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,Mat*); 509014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetGhosts(Mat, PetscInt *,const PetscInt *[]); 5108efafbd8SBarry Smith 511014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt); 512aa1e27eaSFande Kong PETSC_EXTERN PetscErrorCode MatIncreaseOverlapSplit(Mat mat,PetscInt n,IS is[],PetscInt ov); 513d2b2a242SBarry Smith PETSC_EXTERN PetscErrorCode MatMPIAIJSetUseScalableIncreaseOverlap(Mat,PetscBool); 5147b80b807SBarry Smith 515014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*); 516014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMultSymbolic(Mat,Mat,PetscReal,Mat*); 517014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMultNumeric(Mat,Mat,Mat); 51822440eb1SKris Buschelman 5197bab7c10SHong Zhang PETSC_EXTERN PetscErrorCode MatMatMatMult(Mat,Mat,Mat,MatReuse,PetscReal,Mat*); 5202df6c5c3SPatrick Farrell PETSC_EXTERN PetscErrorCode MatGalerkin(Mat,Mat,Mat,MatReuse,PetscReal,Mat*); 5217bab7c10SHong Zhang 522014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*); 523014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAPSymbolic(Mat,Mat,PetscReal,Mat*); 524014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAPNumeric(Mat,Mat,Mat); 525014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARt(Mat,Mat,MatReuse,PetscReal,Mat*); 526014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARtSymbolic(Mat,Mat,PetscReal,Mat*); 527014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARtNumeric(Mat,Mat,Mat); 52822440eb1SKris Buschelman 529014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeMatMult(Mat,Mat,MatReuse,PetscReal,Mat*); 530014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposetMatMultSymbolic(Mat,Mat,PetscReal,Mat*); 531014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposetMatMultNumeric(Mat,Mat,Mat); 532014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMult(Mat,Mat,MatReuse,PetscReal,Mat*); 533014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMultSymbolic(Mat,Mat,PetscReal,Mat*); 534014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMultNumeric(Mat,Mat,Mat); 535bc011b1eSHong Zhang 536014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAXPY(Mat,PetscScalar,Mat,MatStructure); 537014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAYPX(Mat,PetscScalar,Mat,MatStructure); 5381c741599SHong Zhang 539014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScale(Mat,PetscScalar); 540014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShift(Mat,PetscScalar); 5417b80b807SBarry Smith 542014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping); 543014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*); 544a93ff8c4SPeter Brune PETSC_EXTERN PetscErrorCode MatGetLayouts(Mat,PetscLayout*,PetscLayout*); 545014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 546014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsLocalIS(Mat,IS,PetscScalar,Vec,Vec); 547014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 548014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocalIS(Mat,IS,PetscScalar,Vec,Vec); 549014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 550014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 551052efed2SBarry Smith 552014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStashSetInitialSize(Mat,PetscInt,PetscInt); 553014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*); 55490f02eecSBarry Smith 555014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInterpolate(Mat,Vec,Vec); 556014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInterpolateAdd(Mat,Vec,Vec,Vec); 557014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestrict(Mat,Vec,Vec); 5582a7a6963SBarry Smith PETSC_EXTERN PetscErrorCode MatCreateVecs(Mat,Vec*,Vec*); 55984cc2a5aSBarry Smith PETSC_DEPRECATED("Use MatCreateVecs()") PETSC_STATIC_INLINE PetscErrorCode MatGetVecs(Mat mat,Vec *x,Vec *y) {return MatCreateVecs(mat,x,y);} 56053cd1579SHong Zhang PETSC_EXTERN PetscErrorCode MatCreateRedundantMatrix(Mat,PetscInt,MPI_Comm,MatReuse,Mat*); 561014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetMultiProcBlock(Mat,MPI_Comm,MatReuse,Mat*); 562014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFindZeroDiagonals(Mat,IS*); 5633a062f41SBarry Smith PETSC_EXTERN PetscErrorCode MatFindOffBlockDiagonalEntries(Mat,IS*); 5649c8f2541SHong Zhang PETSC_EXTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat(MPI_Comm,Mat,PetscInt,MatReuse,Mat*); 565bd481603SBarry Smith 566bd481603SBarry Smith /*MC 5671620fd73SBarry Smith MatSetValue - Set a single entry into a matrix. 5681620fd73SBarry Smith 5691620fd73SBarry Smith Not collective 5701620fd73SBarry Smith 571a9834a7dSSatish Balay Synopsis: 572a9834a7dSSatish Balay #include <petscmat.h> 573a9834a7dSSatish Balay PetscErrorCode MatSetValue(Mat m,PetscInt row,PetscInt col,PetscScalar value,InsertMode mode) 574a9834a7dSSatish Balay 5751620fd73SBarry Smith Input Parameters: 5761620fd73SBarry Smith + m - the matrix 5771620fd73SBarry Smith . row - the row location of the entry 5781620fd73SBarry Smith . col - the column location of the entry 5791620fd73SBarry Smith . value - the value to insert 5801620fd73SBarry Smith - mode - either INSERT_VALUES or ADD_VALUES 5811620fd73SBarry Smith 5821620fd73SBarry Smith Notes: 5831620fd73SBarry Smith For efficiency one should use MatSetValues() and set several or many 5841620fd73SBarry Smith values simultaneously if possible. 5851620fd73SBarry Smith 5861620fd73SBarry Smith Level: beginner 5871620fd73SBarry Smith 5881620fd73SBarry Smith .seealso: MatSetValues(), MatSetValueLocal() 5891620fd73SBarry Smith M*/ 5901620fd73SBarry 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);} 5911620fd73SBarry Smith 5921620fd73SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatGetValue(Mat v,PetscInt i,PetscInt j,PetscScalar *va) {return MatGetValues(v,1,&i,1,&j,va);} 5931620fd73SBarry Smith 5941620fd73SBarry 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);} 5951620fd73SBarry Smith 5961620fd73SBarry Smith /*MC 5970d2a7ff2SSatish Balay MatPreallocateInitialize - Begins the block of code that will count the number of nonzeros per 598bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 599bd481603SBarry Smith 600bd481603SBarry Smith Synopsis: 601aaa7dc30SBarry Smith #include <petscmat.h> 602c1ac3661SBarry Smith PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz) 603bd481603SBarry Smith 604bd481603SBarry Smith Collective on MPI_Comm 605bd481603SBarry Smith 606bd481603SBarry Smith Input Parameters: 607bd481603SBarry Smith + comm - the communicator that will share the eventually allocated matrix 608859fcb39SBarry Smith . nrows - the number of LOCAL rows in the matrix 609859fcb39SBarry Smith - ncols - the number of LOCAL columns in the matrix 610bd481603SBarry Smith 611bd481603SBarry Smith Output Parameters: 612bd481603SBarry Smith + dnz - the array that will be passed to the matrix preallocation routines 613bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 614bd481603SBarry Smith 615bd481603SBarry Smith Level: intermediate 616bd481603SBarry Smith 617bd481603SBarry Smith Notes: 618a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 619bd481603SBarry Smith 6201d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 621bd481603SBarry Smith 622bd481603SBarry Smith Use MatPreallocateInitializeSymmetric() for symmetric matrices (MPISBAIJ matrices) 623bd481603SBarry Smith 6241620fd73SBarry Smith This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize(). 6251620fd73SBarry Smith 626bd481603SBarry Smith Concepts: preallocation^Matrix 627bd481603SBarry Smith 628d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(), 629d6e23781SBarry Smith MatPreallocateInitializeSymmetric(), MatPreallocateSymmetricSetLocalBlock() 630bd481603SBarry Smith M*/ 631c4f061fbSSatish Balay #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \ 6327c922b88SBarry Smith { \ 633a89bc211SBarry Smith PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__start,__end; \ 6341795a4d1SJed Brown _4_ierr = PetscCalloc2(__nrows,&dnz,__nrows,&onz);CHKERRQ(_4_ierr); \ 6351795a4d1SJed Brown __start = 0; __end = __start; \ 636c1ac3661SBarry Smith _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ctmp;\ 637a89bc211SBarry Smith _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows; 6387c922b88SBarry Smith 639bd481603SBarry Smith /*MC 640bd481603SBarry Smith MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be 641bd481603SBarry Smith inserted using a local number of the rows and columns 642bd481603SBarry Smith 643bd481603SBarry Smith Synopsis: 644aaa7dc30SBarry Smith #include <petscmat.h> 645c1ac3661SBarry Smith PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 646bd481603SBarry Smith 647bd481603SBarry Smith Not Collective 648bd481603SBarry Smith 649bd481603SBarry Smith Input Parameters: 650784ac674SJed Brown + map - the row mapping from local numbering to global numbering 651bd481603SBarry Smith . nrows - the number of rows indicated 6521d73ed98SBarry Smith . rows - the indices of the rows 653784ac674SJed Brown . cmap - the column mapping from local to global numbering 654bd481603SBarry Smith . ncols - the number of columns in the matrix 655bd481603SBarry Smith . cols - the columns indicated 656bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 657bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 658bd481603SBarry Smith 659bd481603SBarry Smith Level: intermediate 660bd481603SBarry Smith 661bd481603SBarry Smith Notes: 662a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 663bd481603SBarry Smith 6641d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 665bd481603SBarry Smith 666bd481603SBarry Smith Concepts: preallocation^Matrix 667bd481603SBarry Smith 668d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(), 669c1154cd5SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocalRemoveDups() 670bd481603SBarry Smith M*/ 671784ac674SJed Brown #define MatPreallocateSetLocal(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \ 672c4f061fbSSatish Balay {\ 673c1ac3661SBarry Smith PetscInt __l;\ 674784ac674SJed Brown _4_ierr = ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\ 675784ac674SJed Brown _4_ierr = ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\ 676c4f061fbSSatish Balay for (__l=0;__l<nrows;__l++) {\ 677ef66eb69SBarry Smith _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 678c4f061fbSSatish Balay }\ 679c4f061fbSSatish Balay } 680c4f061fbSSatish Balay 681bd481603SBarry Smith /*MC 682c1154cd5SBarry Smith MatPreallocateSetLocalRemoveDups - Indicates the locations (rows and columns) in the matrix where nonzeros will be 683c1154cd5SBarry Smith inserted using a local number of the rows and columns. This version removes any duplicate columns in cols 684c1154cd5SBarry Smith 685c1154cd5SBarry Smith Synopsis: 686c1154cd5SBarry Smith #include <petscmat.h> 687c1154cd5SBarry Smith PetscErrorCode MatPreallocateSetLocalRemoveDups(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 688c1154cd5SBarry Smith 689c1154cd5SBarry Smith Not Collective 690c1154cd5SBarry Smith 691c1154cd5SBarry Smith Input Parameters: 692c1154cd5SBarry Smith + map - the row mapping from local numbering to global numbering 693c1154cd5SBarry Smith . nrows - the number of rows indicated 694c1154cd5SBarry Smith . rows - the indices of the rows (these values are mapped to the global values) 695c1154cd5SBarry Smith . cmap - the column mapping from local to global numbering 696c1154cd5SBarry Smith . ncols - the number of columns in the matrix (this value will be changed if duplicate columns are found) 697c1154cd5SBarry Smith . cols - the columns indicated (these values are mapped to the global values, they are then sorted and duplicates removed) 698c1154cd5SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 699c1154cd5SBarry Smith - ozn - the other array passed to the matrix preallocation routines 700c1154cd5SBarry Smith 701c1154cd5SBarry Smith Level: intermediate 702c1154cd5SBarry Smith 703c1154cd5SBarry Smith Notes: 704c1154cd5SBarry Smith See Users-Manual: ch_performance for more details. 705c1154cd5SBarry Smith 706c1154cd5SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 707c1154cd5SBarry Smith 708c1154cd5SBarry Smith Concepts: preallocation^Matrix 709c1154cd5SBarry Smith 710c1154cd5SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(), 711c1154cd5SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocal() 712c1154cd5SBarry Smith M*/ 713c1154cd5SBarry Smith #define MatPreallocateSetLocalRemoveDups(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \ 714c1154cd5SBarry Smith {\ 715c1154cd5SBarry Smith PetscInt __l;\ 716c1154cd5SBarry Smith _4_ierr = ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\ 717c1154cd5SBarry Smith _4_ierr = ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\ 718c1154cd5SBarry Smith _4_ierr = PetscSortRemoveDupsInt(&ncols,cols);CHKERRQ(_4_ierr);\ 719c1154cd5SBarry Smith for (__l=0;__l<nrows;__l++) {\ 720c1154cd5SBarry Smith _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 721c1154cd5SBarry Smith }\ 722c1154cd5SBarry Smith } 723c1154cd5SBarry Smith 724c1154cd5SBarry Smith /*MC 725d6e23781SBarry Smith MatPreallocateSetLocalBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be 726bd481603SBarry Smith inserted using a local number of the rows and columns 727bd481603SBarry Smith 728bd481603SBarry Smith Synopsis: 729aaa7dc30SBarry Smith #include <petscmat.h> 730d6e23781SBarry Smith PetscErrorCode MatPreallocateSetLocalBlock(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 731d6e23781SBarry Smith 732d6e23781SBarry Smith Not Collective 733d6e23781SBarry Smith 734d6e23781SBarry Smith Input Parameters: 735d6e23781SBarry Smith + map - the row mapping from local numbering to global numbering 736d6e23781SBarry Smith . nrows - the number of rows indicated 737d6e23781SBarry Smith . rows - the indices of the rows 738d6e23781SBarry Smith . cmap - the column mapping from local to global numbering 739d6e23781SBarry Smith . ncols - the number of columns in the matrix 740d6e23781SBarry Smith . cols - the columns indicated 741d6e23781SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 742d6e23781SBarry Smith - ozn - the other array passed to the matrix preallocation routines 743d6e23781SBarry Smith 744d6e23781SBarry Smith Level: intermediate 745d6e23781SBarry Smith 746d6e23781SBarry Smith Notes: 747d6e23781SBarry Smith See Users-Manual: ch_performance for more details. 748d6e23781SBarry Smith 749d6e23781SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 750d6e23781SBarry Smith 751d6e23781SBarry Smith Concepts: preallocation^Matrix 752d6e23781SBarry Smith 753d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(), 754d6e23781SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock() 755d6e23781SBarry Smith M*/ 756d6e23781SBarry Smith #define MatPreallocateSetLocalBlock(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \ 757d6e23781SBarry Smith {\ 758d6e23781SBarry Smith PetscInt __l;\ 759d6e23781SBarry Smith _4_ierr = ISLocalToGlobalMappingApplyBlock(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\ 760d6e23781SBarry Smith _4_ierr = ISLocalToGlobalMappingApplyBlock(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\ 761d6e23781SBarry Smith for (__l=0;__l<nrows;__l++) {\ 762d6e23781SBarry Smith _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 763d6e23781SBarry Smith }\ 764d6e23781SBarry Smith } 765d6e23781SBarry Smith 766d6e23781SBarry Smith /*MC 767d6e23781SBarry Smith MatPreallocateSymmetricSetLocalBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be 768d6e23781SBarry Smith inserted using a local number of the rows and columns 769d6e23781SBarry Smith 770d6e23781SBarry Smith Synopsis: 771d6e23781SBarry Smith #include <petscmat.h> 772d6e23781SBarry Smith PetscErrorCode MatPreallocateSymmetricSetLocalBlock(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 773bd481603SBarry Smith 774bd481603SBarry Smith Not Collective 775bd481603SBarry Smith 776bd481603SBarry Smith Input Parameters: 777bd481603SBarry Smith + map - the mapping between local numbering and global numbering 778bd481603SBarry Smith . nrows - the number of rows indicated 7791d73ed98SBarry Smith . rows - the indices of the rows 780bd481603SBarry Smith . ncols - the number of columns in the matrix 781bd481603SBarry Smith . cols - the columns indicated 782bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 783bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 784bd481603SBarry Smith 785bd481603SBarry Smith Level: intermediate 786bd481603SBarry Smith 787bd481603SBarry Smith Notes: 788a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 789bd481603SBarry Smith 790bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 791bd481603SBarry Smith 792bd481603SBarry Smith Concepts: preallocation^Matrix 793bd481603SBarry Smith 794d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateInitialize(), 795d6e23781SBarry Smith MatPreallocateInitialize(), MatPreallocateSetLocal() 796bd481603SBarry Smith M*/ 797d6e23781SBarry Smith #define MatPreallocateSymmetricSetLocalBlock(map,nrows,rows,ncols,cols,dnz,onz) 0;\ 798d3d32019SBarry Smith {\ 799c1ac3661SBarry Smith PetscInt __l;\ 800d6e23781SBarry Smith _4_ierr = ISLocalToGlobalMappingApplyBlock(map,nrows,rows,rows);CHKERRQ(_4_ierr);\ 801d6e23781SBarry Smith _4_ierr = ISLocalToGlobalMappingApplyBlock(map,ncols,cols,cols);CHKERRQ(_4_ierr);\ 802d3d32019SBarry Smith for (__l=0;__l<nrows;__l++) {\ 803d6e23781SBarry Smith _4_ierr = MatPreallocateSymmetricSetBlock((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 804d3d32019SBarry Smith }\ 805d3d32019SBarry Smith } 806bd481603SBarry Smith /*MC 807bd481603SBarry Smith MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be 808bd481603SBarry Smith inserted using a local number of the rows and columns 809bd481603SBarry Smith 810bd481603SBarry Smith Synopsis: 811aaa7dc30SBarry Smith #include <petscmat.h> 812c1ac3661SBarry Smith PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 813bd481603SBarry Smith 814bd481603SBarry Smith Not Collective 815bd481603SBarry Smith 816bd481603SBarry Smith Input Parameters: 81764075487SBarry Smith + row - the row 818bd481603SBarry Smith . ncols - the number of columns in the matrix 819a50f8bf6SHong Zhang - cols - the columns indicated 820a50f8bf6SHong Zhang 821a50f8bf6SHong Zhang Output Parameters: 822a50f8bf6SHong Zhang + dnz - the array that will be passed to the matrix preallocation routines 823bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 824bd481603SBarry Smith 825bd481603SBarry Smith Level: intermediate 826bd481603SBarry Smith 827bd481603SBarry Smith Notes: 828a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 829bd481603SBarry Smith 830bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 831bd481603SBarry Smith 8321620fd73SBarry Smith This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize(). 8331620fd73SBarry Smith 834bd481603SBarry Smith Concepts: preallocation^Matrix 835bd481603SBarry Smith 836d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(), 837d6e23781SBarry Smith MatPreallocateInitialize(), MatPreallocateSetLocal() 838bd481603SBarry Smith M*/ 839c4f061fbSSatish Balay #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\ 840c1ac3661SBarry Smith { PetscInt __i; \ 841e32f2f54SBarry 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);\ 842e32f2f54SBarry 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);\ 8437c922b88SBarry Smith for (__i=0; __i<nc; __i++) {\ 84464075487SBarry Smith if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \ 8457cc688d7SBarry Smith else dnz[row - __rstart]++;\ 8467c922b88SBarry Smith }\ 8477c922b88SBarry Smith } 8487c922b88SBarry Smith 849bd481603SBarry Smith /*MC 850d6e23781SBarry Smith MatPreallocateSymmetricSetBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be 851bd481603SBarry Smith inserted using a local number of the rows and columns 852bd481603SBarry Smith 853bd481603SBarry Smith Synopsis: 854aaa7dc30SBarry Smith #include <petscmat.h> 855d6e23781SBarry Smith PetscErrorCode MatPreallocateSymmetricSetBlock(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 856bd481603SBarry Smith 857bd481603SBarry Smith Not Collective 858bd481603SBarry Smith 859bd481603SBarry Smith Input Parameters: 860bd481603SBarry Smith + nrows - the number of rows indicated 8611d73ed98SBarry Smith . rows - the indices of the rows 862bd481603SBarry Smith . ncols - the number of columns in the matrix 863bd481603SBarry Smith . cols - the columns indicated 864bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 865bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 866bd481603SBarry Smith 867bd481603SBarry Smith Level: intermediate 868bd481603SBarry Smith 869bd481603SBarry Smith Notes: 870a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 871bd481603SBarry Smith 872bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 873bd481603SBarry Smith 8741620fd73SBarry Smith This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize(). 8751620fd73SBarry Smith 876bd481603SBarry Smith Concepts: preallocation^Matrix 877bd481603SBarry Smith 878d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateInitialize(), 879d6e23781SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocal() 880bd481603SBarry Smith M*/ 881d6e23781SBarry Smith #define MatPreallocateSymmetricSetBlock(row,nc,cols,dnz,onz) 0;\ 882c1ac3661SBarry Smith { PetscInt __i; \ 883d3d32019SBarry Smith for (__i=0; __i<nc; __i++) {\ 884d3d32019SBarry Smith if (cols[__i] >= __end) onz[row - __rstart]++; \ 885d3d32019SBarry Smith else if (cols[__i] >= row) dnz[row - __rstart]++;\ 886d3d32019SBarry Smith }\ 887d3d32019SBarry Smith } 888d3d32019SBarry Smith 889bd481603SBarry Smith /*MC 89016371a99SBarry Smith MatPreallocateLocation - An alternative to MatPreallocationSet() that puts the nonzero locations into the matrix if it exists 89116371a99SBarry Smith 89216371a99SBarry Smith Synopsis: 893aaa7dc30SBarry Smith #include <petscmat.h> 89416371a99SBarry Smith PetscErrorCode MatPreallocateLocations(Mat A,PetscInt row,PetscInt ncols,PetscInt *cols,PetscInt *dnz,PetscInt *onz) 89516371a99SBarry Smith 89616371a99SBarry Smith Not Collective 89716371a99SBarry Smith 89816371a99SBarry Smith Input Parameters: 89916371a99SBarry Smith . A - matrix 90016371a99SBarry Smith . row - row where values exist (must be local to this process) 90116371a99SBarry Smith . ncols - number of columns 90216371a99SBarry Smith . cols - columns with nonzeros 90316371a99SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 90416371a99SBarry Smith - ozn - the other array passed to the matrix preallocation routines 90516371a99SBarry Smith 90616371a99SBarry Smith Level: intermediate 90716371a99SBarry Smith 90816371a99SBarry Smith Notes: 909a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 91016371a99SBarry Smith 91116371a99SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 91216371a99SBarry Smith 91316371a99SBarry Smith This is a MACRO not a function because it uses a bunch of variables private to the MatPreallocation.... routines. 91416371a99SBarry Smith 91516371a99SBarry Smith Concepts: preallocation^Matrix 91616371a99SBarry Smith 917d6e23781SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(), 918d6e23781SBarry Smith MatPreallocateSymmetricSetLocalBlock() 91916371a99SBarry Smith M*/ 9200298fd71SBarry 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);} 92116371a99SBarry Smith 92216371a99SBarry Smith 92316371a99SBarry Smith /*MC 9240d2a7ff2SSatish Balay MatPreallocateFinalize - Ends the block of code that will count the number of nonzeros per 925bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 926bd481603SBarry Smith 927bd481603SBarry Smith Synopsis: 928aaa7dc30SBarry Smith #include <petscmat.h> 929c1ac3661SBarry Smith PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz) 930bd481603SBarry Smith 931bd481603SBarry Smith Collective on MPI_Comm 932bd481603SBarry Smith 933bd481603SBarry Smith Input Parameters: 93416371a99SBarry Smith + dnz - the array that was be passed to the matrix preallocation routines 935bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 936bd481603SBarry Smith 937bd481603SBarry Smith Level: intermediate 938bd481603SBarry Smith 939bd481603SBarry Smith Notes: 940a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 941bd481603SBarry Smith 942bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 943bd481603SBarry Smith 9441620fd73SBarry Smith This is a MACRO not a function because it closes the { started in MatPreallocateInitialize(). 9451620fd73SBarry Smith 946bd481603SBarry Smith Concepts: preallocation^Matrix 947bd481603SBarry Smith 948d6e23781SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(), 949d6e23781SBarry Smith MatPreallocateSymmetricSetLocalBlock() 950bd481603SBarry Smith M*/ 951a89bc211SBarry Smith #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree2(dnz,onz);CHKERRQ(_4_ierr);} 9527c922b88SBarry Smith 9537b80b807SBarry Smith /* Routines unique to particular data structures */ 954014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellGetContext(Mat,void *); 955698d4c6aSKris Buschelman 956014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInodeAdjustForInodes(Mat,IS*,IS*); 957014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *); 958698d4c6aSKris Buschelman 959014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetColumnIndices(Mat,PetscInt[]); 960014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetColumnIndices(Mat,PetscInt[]); 961014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 962014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 963014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 964014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJFromTriple(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*,PetscInt,PetscBool); 9657b80b807SBarry Smith 966a96a251dSBarry Smith #define MAT_SKIP_ALLOCATION -4 967a96a251dSBarry Smith 968014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]); 969014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]); 970014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]); 971273d9f13SBarry Smith 972014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 973014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 974014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 975014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []); 976014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 977014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]); 978014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 979014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]); 980014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIDenseSetPreallocation(Mat,PetscScalar[]); 981014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqDenseSetPreallocation(Mat,PetscScalar[]); 9829230625dSJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,const PetscInt*[]); 9839230625dSJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,const PetscInt*[]); 984014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAdjCreateNonemptySubcommMat(Mat,Mat*); 985273d9f13SBarry Smith 9862e1947a5SStefano Zampini PETSC_EXTERN PetscErrorCode MatISSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 987cf0a3239SStefano Zampini PETSC_EXTERN PetscErrorCode MatISSetUpSF(Mat); 988014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqDenseSetLDA(Mat,PetscInt); 989014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDenseGetLocalMatrix(Mat,Mat*); 9901b807ce4Svictorle 991014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStoreValues(Mat); 992014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRetrieveValues(Mat); 9932e8a6d31SBarry Smith 994014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDAADSetCtx(Mat,void*); 9953a7fca6bSBarry Smith 996014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFindNonzeroRows(Mat,IS*); 9977cfe41e4SPatrick Farrell PETSC_EXTERN PetscErrorCode MatFindZeroRows(Mat,IS*); 9987b80b807SBarry Smith /* 9997b80b807SBarry Smith These routines are not usually accessed directly, rather solving is 100094b7f48cSBarry Smith done through the KSP and PC interfaces. 10017b80b807SBarry Smith */ 10027b80b807SBarry Smith 100376bdecfbSBarry Smith /*J 10048f6c3df8SBarry Smith MatOrderingType - String with the name of a PETSc matrix ordering 1005d9274352SBarry Smith 1006d9274352SBarry Smith Level: beginner 1007d9274352SBarry Smith 1008d9274352SBarry Smith .seealso: MatGetOrdering() 100976bdecfbSBarry Smith J*/ 101019fd82e9SBarry Smith typedef const char* MatOrderingType; 10112692d6eeSBarry Smith #define MATORDERINGNATURAL "natural" 10122692d6eeSBarry Smith #define MATORDERINGND "nd" 10132692d6eeSBarry Smith #define MATORDERING1WD "1wd" 10142692d6eeSBarry Smith #define MATORDERINGRCM "rcm" 10152692d6eeSBarry Smith #define MATORDERINGQMD "qmd" 10162692d6eeSBarry Smith #define MATORDERINGROWLENGTH "rowlength" 1017510184a7SJed Brown #define MATORDERINGWBM "wbm" 1018fc1bef75SJed Brown #define MATORDERINGSPECTRAL "spectral" 1019312e372aSBarry Smith #define MATORDERINGAMD "amd" /* only works if UMFPACK is installed with PETSc */ 1020b12f92e5SBarry Smith 102119fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatGetOrdering(Mat,MatOrderingType,IS*,IS*); 1022140e18c1SBarry Smith PETSC_EXTERN PetscErrorCode MatGetOrderingList(PetscFunctionList*); 1023bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatOrderingRegister(const char[],PetscErrorCode(*)(Mat,MatOrderingType,IS*,IS*)); 1024140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatOrderingList; 1025d4fbbf0eSBarry Smith 1026014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS); 1027fc1bef75SJed Brown PETSC_EXTERN PetscErrorCode MatCreateLaplacian(Mat,PetscReal,PetscBool,Mat*); 1028a2ce50c7SBarry Smith 1029d91e6319SBarry Smith /*S 1030d90ac83dSHong Zhang MatFactorShiftType - Numeric Shift. 1031d90ac83dSHong Zhang 1032d90ac83dSHong Zhang Level: beginner 1033d90ac83dSHong Zhang 1034d90ac83dSHong Zhang S*/ 1035d90ac83dSHong Zhang typedef enum {MAT_SHIFT_NONE,MAT_SHIFT_NONZERO,MAT_SHIFT_POSITIVE_DEFINITE,MAT_SHIFT_INBLOCKS} MatFactorShiftType; 10366a6fc655SJed Brown PETSC_EXTERN const char *const MatFactorShiftTypes[]; 10375e9742b9SJed Brown PETSC_EXTERN const char *const MatFactorShiftTypesDetail[]; 1038d90ac83dSHong Zhang 1039d90ac83dSHong Zhang /*S 10409aa7eafdSHong Zhang MatFactorError - indicates what type of error in matrix factor 10419aa7eafdSHong Zhang 10429aa7eafdSHong Zhang Level: beginner 10439aa7eafdSHong Zhang 104400e125f8SBarry Smith Developer Notes: Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 104500e125f8SBarry Smith 104600e125f8SBarry Smith .seealso: MatGetFactor() 10479aa7eafdSHong Zhang S*/ 10485cd7cf9dSHong Zhang typedef enum {MAT_FACTOR_NOERROR,MAT_FACTOR_STRUCT_ZEROPIVOT,MAT_FACTOR_NUMERIC_ZEROPIVOT,MAT_FACTOR_OUTMEMORY,MAT_FACTOR_OTHER} MatFactorError; 10499aa7eafdSHong Zhang 105000e125f8SBarry Smith PETSC_EXTERN PetscErrorCode MatFactorGetError(Mat,MatFactorError*); 1051b8b68cfdSBarry Smith PETSC_EXTERN PetscErrorCode MatFactorClearError(Mat); 10527b6c816cSBarry Smith PETSC_EXTERN PetscErrorCode MatFactorGetErrorZeroPivot(Mat,PetscReal*,PetscInt*); 105300e125f8SBarry Smith 10549aa7eafdSHong Zhang /*S 1055422a814eSBarry Smith MatFactorInfo - Data passed into the matrix factorization routines, and information about the resulting factorization 1056b00f7748SHong Zhang 105761cad860SBarry Smith In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE, that is use 105861cad860SBarry Smith $ MatFactorInfo info(MAT_FACTORINFO_SIZE) 1059b00f7748SHong Zhang 106015e8a5b3SHong Zhang Notes: These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC. 1061b00f7748SHong Zhang 106261cad860SBarry Smith You can use MatFactorInfoInitialize() to set default values. 106361cad860SBarry Smith 1064b00f7748SHong Zhang Level: developer 1065b00f7748SHong Zhang 1066d7d82daaSBarry Smith .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(), 1067d7d82daaSBarry Smith MatFactorInfoInitialize() 1068b00f7748SHong Zhang 1069b00f7748SHong Zhang S*/ 1070b00f7748SHong Zhang typedef struct { 107115e8a5b3SHong Zhang PetscReal diagonal_fill; /* force diagonal to fill in if initially not filled */ 107285317021SBarry Smith PetscReal usedt; 107315e8a5b3SHong Zhang PetscReal dt; /* drop tolerance */ 1074b00f7748SHong Zhang PetscReal dtcol; /* tolerance for pivoting */ 107515e8a5b3SHong Zhang PetscReal dtcount; /* maximum nonzeros to be allowed per row */ 107667e28bfeSBarry Smith PetscReal fill; /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */ 1077348344bbSBarry Smith PetscReal levels; /* ICC/ILU(levels) */ 1078bcd9e38bSBarry Smith PetscReal pivotinblocks; /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0 1079bcd9e38bSBarry Smith factorization may be faster if do not pivot */ 108015e8a5b3SHong Zhang PetscReal zeropivot; /* pivot is called zero if less than this */ 1081f4db908eSBarry Smith PetscReal shifttype; /* type of shift added to matrix factor to prevent zero pivots */ 1082f4db908eSBarry Smith PetscReal shiftamount; /* how large the shift is */ 108315e8a5b3SHong Zhang } MatFactorInfo; 1084ffa6d0a5SLois Curfman McInnes 1085014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFactorInfoInitialize(MatFactorInfo*); 10869a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatCholeskyFactor(Mat,IS,const MatFactorInfo*); 10879a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatCholeskyFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*); 10889a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatCholeskyFactorNumeric(Mat,Mat,const MatFactorInfo*); 10899a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatLUFactor(Mat,IS,IS,const MatFactorInfo*); 10909a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatILUFactor(Mat,IS,IS,const MatFactorInfo*); 10919a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatLUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*); 10929a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatILUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*); 10939a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatICCFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*); 10949a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatICCFactor(Mat,IS,const MatFactorInfo*); 10959a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatLUFactorNumeric(Mat,Mat,const MatFactorInfo*); 1096014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*); 1097014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolve(Mat,Vec,Vec); 1098014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatForwardSolve(Mat,Vec,Vec); 1099014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatBackwardSolve(Mat,Vec,Vec); 1100014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveAdd(Mat,Vec,Vec,Vec); 1101014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveTranspose(Mat,Vec,Vec); 1102014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveTransposeAdd(Mat,Vec,Vec,Vec); 1103014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolves(Mat,Vecs,Vecs); 11048e7ba810SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorSetSchurIS(Mat,IS); 11055a05ddb0SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorGetSchurComplement(Mat,Mat*); 11065a05ddb0SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorRestoreSchurComplement(Mat,Mat*); 11075a05ddb0SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorInvertSchurComplement(Mat); 11085a05ddb0SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorCreateSchurComplement(Mat,Mat*); 11095a05ddb0SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorSolveSchurComplement(Mat,Vec,Vec); 11105a05ddb0SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorSolveSchurComplementTranspose(Mat,Vec,Vec); 11116dba178dSStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorFactorizeSchurComplement(Mat); 1112e8ade678SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorSetSchurComplementSolverType(Mat,PetscInt); 1113014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetUnfactored(Mat); 1114bb5a7306SBarry Smith 1115d91e6319SBarry Smith /*E 1116d91e6319SBarry Smith MatSORType - What type of (S)SOR to perform 1117bb1eb677SSatish Balay 1118d91e6319SBarry Smith Level: beginner 1119d91e6319SBarry Smith 1120d9274352SBarry Smith May be bitwise ORd together 1121d9274352SBarry Smith 1122af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 1123d91e6319SBarry Smith 11244e7234bfSBarry Smith MatSORType may be bitwise ORd together, so do not change the numbers 11254e7234bfSBarry Smith 112641f059aeSBarry Smith .seealso: MatSOR() 1127d91e6319SBarry Smith E*/ 1128ee50ffe9SBarry Smith typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3, 1129ee50ffe9SBarry Smith SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8, 1130ee50ffe9SBarry Smith SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16, 113184cb2905SBarry Smith SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType; 1132014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSOR(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec); 11338ed539a5SBarry Smith 1134d4fbbf0eSBarry Smith /* 1135639f9d9dSBarry Smith These routines are for efficiently computing Jacobians via finite differences. 1136639f9d9dSBarry Smith */ 1137b12f92e5SBarry Smith 11387086a01eSPeter Brune /*S 11397086a01eSPeter Brune MatColoring - Object for managing the coloring of matrices. 11407086a01eSPeter Brune 11417086a01eSPeter Brune Level: beginner 11427086a01eSPeter Brune 11437086a01eSPeter Brune Concepts: matrix, coloring 11447086a01eSPeter Brune 11457086a01eSPeter Brune .seealso: MatFDColoringCreate() ISColoring MatFDColoring 11467086a01eSPeter Brune S*/ 11477086a01eSPeter Brune typedef struct _p_MatColoring* MatColoring; 114876bdecfbSBarry Smith /*J 11498f6c3df8SBarry Smith MatColoringType - String with the name of a PETSc matrix coloring 1150d9274352SBarry Smith 1151d9274352SBarry Smith Level: beginner 1152d9274352SBarry Smith 11537086a01eSPeter Brune .seealso: MatColoringSetType(), MatColoring 115476bdecfbSBarry Smith J*/ 11557086a01eSPeter Brune 115619fd82e9SBarry Smith typedef const char* MatColoringType; 11575567d71aSPeter Brune #define MATCOLORINGJP "jp" 1158b9acb617SPeter Brune #define MATCOLORINGPOWER "power" 11592692d6eeSBarry Smith #define MATCOLORINGNATURAL "natural" 11602692d6eeSBarry Smith #define MATCOLORINGSL "sl" 11612692d6eeSBarry Smith #define MATCOLORINGLF "lf" 11622692d6eeSBarry Smith #define MATCOLORINGID "id" 11638121bdceSPeter Brune #define MATCOLORINGGREEDY "greedy" 1164b12f92e5SBarry Smith 11658ac6da0aSPeter Brune /*E 11668ac6da0aSPeter Brune MatColoringWeightType - Type of weight scheme 11678ac6da0aSPeter Brune 11688ac6da0aSPeter Brune Not Collective 11698ac6da0aSPeter Brune 11708ac6da0aSPeter Brune + MAT_COLORING_RANDOM - Random weights 11718ac6da0aSPeter Brune . MAT_COLORING_LEXICAL - Lexical weighting based upon global numbering. 11728ac6da0aSPeter Brune - MAT_COLORING_LF - Last-first weighting. 11738ac6da0aSPeter Brune 11748ac6da0aSPeter Brune Level: intermediate 11758ac6da0aSPeter Brune 1176af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 11778ac6da0aSPeter Brune 11788ac6da0aSPeter Brune .seealso: MatCUSPSetFormat(), MatCUSPFormatOperation 11798ac6da0aSPeter Brune E*/ 1180301c140bSPeter Brune typedef enum {MAT_COLORING_WEIGHT_RANDOM,MAT_COLORING_WEIGHT_LEXICAL,MAT_COLORING_WEIGHT_LF,MAT_COLORING_WEIGHT_SL} MatColoringWeightType; 11818ac6da0aSPeter Brune 1182335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringCreate(Mat,MatColoring*); 1183d7f2a307SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringGetDegrees(Mat,PetscInt,PetscInt*); 1184335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringDestroy(MatColoring*); 1185335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringView(MatColoring,PetscViewer); 1186335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetType(MatColoring,MatColoringType); 1187335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetFromOptions(MatColoring); 1188335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetDistance(MatColoring,PetscInt); 1189335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringGetDistance(MatColoring,PetscInt*); 1190335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetMaxColors(MatColoring,PetscInt); 1191335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringGetMaxColors(MatColoring,PetscInt*); 1192335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringApply(MatColoring,ISColoring*); 1193335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringRegister(const char[],PetscErrorCode(*)(MatColoring)); 1194014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*); 11958ac6da0aSPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetWeightType(MatColoring,MatColoringWeightType); 1196c29971fcSPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetWeights(MatColoring,PetscReal*,PetscInt*); 11978ac6da0aSPeter Brune PETSC_EXTERN PetscErrorCode MatColoringCreateWeights(MatColoring,PetscReal **,PetscInt **lperm); 1198cb394dc5SBarry Smith PETSC_EXTERN PetscErrorCode MatColoringTestValid(MatColoring,ISColoring); 1199639f9d9dSBarry Smith 1200d9274352SBarry Smith /*S 1201d9274352SBarry Smith MatFDColoring - Object for computing a sparse Jacobian via finite differences 1202d9274352SBarry Smith and coloring 1203639f9d9dSBarry Smith 1204d9274352SBarry Smith Level: beginner 1205d9274352SBarry Smith 1206d9274352SBarry Smith Concepts: coloring, sparse Jacobian, finite differences 1207d9274352SBarry Smith 1208d9274352SBarry Smith .seealso: MatFDColoringCreate() 1209d9274352SBarry Smith S*/ 1210e2a1c21fSSatish Balay typedef struct _p_MatFDColoring* MatFDColoring; 1211639f9d9dSBarry Smith 1212014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringCreate(Mat,ISColoring,MatFDColoring *); 1213014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringDestroy(MatFDColoring*); 1214014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringView(MatFDColoring,PetscViewer); 1215014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*); 1216014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**); 1217014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal); 1218014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetFromOptions(MatFDColoring); 1219d1e9a80fSBarry Smith PETSC_EXTERN PetscErrorCode MatFDColoringApply(Mat,MatFDColoring,Vec,void *); 1220014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetF(MatFDColoring,Vec); 1221014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,PetscInt*[]); 1222f86b9fbaSHong Zhang PETSC_EXTERN PetscErrorCode MatFDColoringSetUp(Mat,ISColoring,MatFDColoring); 1223f86b9fbaSHong Zhang PETSC_EXTERN PetscErrorCode MatFDColoringSetBlockSize(MatFDColoring,PetscInt,PetscInt); 1224f86b9fbaSHong Zhang 1225b1683b59SHong Zhang 1226b1683b59SHong Zhang /*S 1227b9af6bddSHong Zhang MatTransposeColoring - Object for computing a sparse matrix product C=A*B^T via coloring 1228b1683b59SHong Zhang 1229b1683b59SHong Zhang Level: beginner 1230b1683b59SHong Zhang 1231b1683b59SHong Zhang Concepts: coloring, sparse matrix product 1232b1683b59SHong Zhang 1233b9af6bddSHong Zhang .seealso: MatTransposeColoringCreate() 1234b1683b59SHong Zhang S*/ 1235b9af6bddSHong Zhang typedef struct _p_MatTransposeColoring* MatTransposeColoring; 1236b1683b59SHong Zhang 1237014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeColoringCreate(Mat,ISColoring,MatTransposeColoring *); 1238014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransColoringApplySpToDen(MatTransposeColoring,Mat,Mat); 1239014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransColoringApplyDenToSp(MatTransposeColoring,Mat,Mat); 1240014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeColoringDestroy(MatTransposeColoring*); 1241b1683b59SHong Zhang 1242639f9d9dSBarry Smith /* 12430752156aSBarry Smith These routines are for partitioning matrices: currently used only 12443eda8832SBarry Smith for adjacency matrix, MatCreateMPIAdj(). 12450752156aSBarry Smith */ 1246ca161407SBarry Smith 1247d9274352SBarry Smith /*S 1248d9274352SBarry Smith MatPartitioning - Object for managing the partitioning of a matrix or graph 1249d9274352SBarry Smith 1250d9274352SBarry Smith Level: beginner 1251d9274352SBarry Smith 1252d9274352SBarry Smith Concepts: partitioning 1253d9274352SBarry Smith 1254743a1026Svictorle .seealso: MatPartitioningCreate(), MatPartitioningType 1255d9274352SBarry Smith S*/ 125691e9ee9fSBarry Smith typedef struct _p_MatPartitioning* MatPartitioning; 1257d9274352SBarry Smith 125876bdecfbSBarry Smith /*J 12598f6c3df8SBarry Smith MatPartitioningType - String with the name of a PETSc matrix partitioning 1260d9274352SBarry Smith 1261d9274352SBarry Smith Level: beginner 12620d04baf8SBarry Smith dm 1263b547a13bSvictorle .seealso: MatPartitioningCreate(), MatPartitioning 126476bdecfbSBarry Smith J*/ 126519fd82e9SBarry Smith typedef const char* MatPartitioningType; 12662692d6eeSBarry Smith #define MATPARTITIONINGCURRENT "current" 1267aa1e27eaSFande Kong #define MATPARTITIONINGAVERAGE "average" 12682692d6eeSBarry Smith #define MATPARTITIONINGSQUARE "square" 12692692d6eeSBarry Smith #define MATPARTITIONINGPARMETIS "parmetis" 12702692d6eeSBarry Smith #define MATPARTITIONINGCHACO "chaco" 12712692d6eeSBarry Smith #define MATPARTITIONINGPARTY "party" 127250d91057SBarry Smith #define MATPARTITIONINGPTSCOTCH "ptscotch" 127388d2ac2bSFande Kong #define MATPARTITIONINGHIERARCH "hierarch" 127471306c60Sjroman 1275ca161407SBarry Smith 1276014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningCreate(MPI_Comm,MatPartitioning*); 127719fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningSetType(MatPartitioning,MatPartitioningType); 1278014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetNParts(MatPartitioning,PetscInt); 1279014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetAdjacency(MatPartitioning,Mat); 1280014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]); 1281014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []); 1282014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningApply(MatPartitioning,IS*); 1283014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningDestroy(MatPartitioning*); 12842aabb6bbSBarry Smith 1285bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningRegister(const char[],PetscErrorCode (*)(MatPartitioning)); 128630de9b25SBarry Smith 1287f1144a30SSatish Balay 12882bad1931SBarry Smith 1289014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningView(MatPartitioning,PetscViewer); 1290014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetFromOptions(MatPartitioning); 129119fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningGetType(MatPartitioning,MatPartitioningType*); 1292ca161407SBarry Smith 129327538973SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetRepartition(MatPartitioning part); 1294014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetCoarseSequential(MatPartitioning); 1295014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningParmetisGetEdgeCut(MatPartitioning, PetscInt *); 12960752156aSBarry Smith 1297b6956c12SJose Roman typedef enum { MP_CHACO_MULTILEVEL=1,MP_CHACO_SPECTRAL=2,MP_CHACO_LINEAR=4,MP_CHACO_RANDOM=5,MP_CHACO_SCATTERED=6 } MPChacoGlobalType; 12986a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoGlobalTypes[]; 1299b6956c12SJose Roman typedef enum { MP_CHACO_KERNIGHAN=1,MP_CHACO_NONE=2 } MPChacoLocalType; 13006a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoLocalTypes[]; 1301b6956c12SJose Roman typedef enum { MP_CHACO_LANCZOS=0,MP_CHACO_RQI=1 } MPChacoEigenType; 13026a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoEigenTypes[]; 1303b6956c12SJose Roman 1304014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetGlobal(MatPartitioning,MPChacoGlobalType); 1305014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetGlobal(MatPartitioning,MPChacoGlobalType*); 1306014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetLocal(MatPartitioning,MPChacoLocalType); 1307014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetLocal(MatPartitioning,MPChacoLocalType*); 1308014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal); 1309014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType); 1310014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenSolver(MatPartitioning,MPChacoEigenType*); 1311014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenTol(MatPartitioning,PetscReal); 1312014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenTol(MatPartitioning,PetscReal*); 1313014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenNumber(MatPartitioning,PetscInt); 1314014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenNumber(MatPartitioning,PetscInt*); 131571306c60Sjroman 131671306c60Sjroman #define MP_PARTY_OPT "opt" 131771306c60Sjroman #define MP_PARTY_LIN "lin" 131871306c60Sjroman #define MP_PARTY_SCA "sca" 131971306c60Sjroman #define MP_PARTY_RAN "ran" 132071306c60Sjroman #define MP_PARTY_GBF "gbf" 132171306c60Sjroman #define MP_PARTY_GCF "gcf" 132271306c60Sjroman #define MP_PARTY_BUB "bub" 132371306c60Sjroman #define MP_PARTY_DEF "def" 1324014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetGlobal(MatPartitioning,const char*); 132571306c60Sjroman #define MP_PARTY_HELPFUL_SETS "hs" 132671306c60Sjroman #define MP_PARTY_KERNIGHAN_LIN "kl" 132771306c60Sjroman #define MP_PARTY_NONE "no" 1328014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetLocal(MatPartitioning,const char*); 1329014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal); 1330014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetBipart(MatPartitioning,PetscBool); 1331014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscBool); 133271306c60Sjroman 133350d91057SBarry Smith typedef enum { MP_PTSCOTCH_QUALITY,MP_PTSCOTCH_SPEED,MP_PTSCOTCH_BALANCE,MP_PTSCOTCH_SAFETY,MP_PTSCOTCH_SCALABILITY } MPPTScotchStrategyType; 13346a6fc655SJed Brown PETSC_EXTERN const char *const MPPTScotchStrategyTypes[]; 1335e0f1cffaSJose Roman 1336014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetImbalance(MatPartitioning,PetscReal); 1337014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetImbalance(MatPartitioning,PetscReal*); 1338014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetStrategy(MatPartitioning,MPPTScotchStrategyType); 1339014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetStrategy(MatPartitioning,MPPTScotchStrategyType*); 134071306c60Sjroman 1341b43b03e9SMark F. Adams /* 13426294fa2bSFande Kong * hierarchical partitioning 13436294fa2bSFande Kong */ 13444e713f55SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalGetFineparts(MatPartitioning,IS*); 13454e713f55SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalGetCoarseparts(MatPartitioning,IS*); 13464e713f55SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalSetNcoarseparts(MatPartitioning,PetscInt); 13474e713f55SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalSetNfineparts(MatPartitioning, PetscInt); 13486294fa2bSFande Kong 1349014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMeshToVertexGraph(Mat,PetscInt,Mat*); 1350014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMeshToCellGraph(Mat,PetscInt,Mat*); 1351591294e4SBarry Smith 13520752156aSBarry Smith /* 1353af0996ceSBarry Smith If you add entries here you must also add them to petsc/finclude/petscmat.h 1354d4fbbf0eSBarry Smith */ 13551c1c02c0SLois Curfman McInnes typedef enum { MATOP_SET_VALUES=0, 13561c1c02c0SLois Curfman McInnes MATOP_GET_ROW=1, 13571c1c02c0SLois Curfman McInnes MATOP_RESTORE_ROW=2, 13581c1c02c0SLois Curfman McInnes MATOP_MULT=3, 13591c1c02c0SLois Curfman McInnes MATOP_MULT_ADD=4, 13607c922b88SBarry Smith MATOP_MULT_TRANSPOSE=5, 13617c922b88SBarry Smith MATOP_MULT_TRANSPOSE_ADD=6, 13621c1c02c0SLois Curfman McInnes MATOP_SOLVE=7, 13631c1c02c0SLois Curfman McInnes MATOP_SOLVE_ADD=8, 13647c922b88SBarry Smith MATOP_SOLVE_TRANSPOSE=9, 13657c922b88SBarry Smith MATOP_SOLVE_TRANSPOSE_ADD=10, 13661c1c02c0SLois Curfman McInnes MATOP_LUFACTOR=11, 13671c1c02c0SLois Curfman McInnes MATOP_CHOLESKYFACTOR=12, 1368a714c06dSBarry Smith MATOP_SOR=13, 13691c1c02c0SLois Curfman McInnes MATOP_TRANSPOSE=14, 13701c1c02c0SLois Curfman McInnes MATOP_GETINFO=15, 13711c1c02c0SLois Curfman McInnes MATOP_EQUAL=16, 13721c1c02c0SLois Curfman McInnes MATOP_GET_DIAGONAL=17, 13731c1c02c0SLois Curfman McInnes MATOP_DIAGONAL_SCALE=18, 13741c1c02c0SLois Curfman McInnes MATOP_NORM=19, 13751c1c02c0SLois Curfman McInnes MATOP_ASSEMBLY_BEGIN=20, 13761c1c02c0SLois Curfman McInnes MATOP_ASSEMBLY_END=21, 1377d519adbfSMatthew Knepley MATOP_SET_OPTION=22, 1378d519adbfSMatthew Knepley MATOP_ZERO_ENTRIES=23, 1379d519adbfSMatthew Knepley MATOP_ZERO_ROWS=24, 1380d519adbfSMatthew Knepley MATOP_LUFACTOR_SYMBOLIC=25, 1381d519adbfSMatthew Knepley MATOP_LUFACTOR_NUMERIC=26, 1382d519adbfSMatthew Knepley MATOP_CHOLESKY_FACTOR_SYMBOLIC=27, 1383d519adbfSMatthew Knepley MATOP_CHOLESKY_FACTOR_NUMERIC=28, 1384d519adbfSMatthew Knepley MATOP_SETUP_PREALLOCATION=29, 1385d519adbfSMatthew Knepley MATOP_ILUFACTOR_SYMBOLIC=30, 1386d519adbfSMatthew Knepley MATOP_ICCFACTOR_SYMBOLIC=31, 1387a5b7ff6bSBarry Smith MATOP_GET_DIAGONAL_BLOCK=32, 13889d39f69dSJed Brown /* MATOP_PLACEHOLDER_33=33, */ 1389d519adbfSMatthew Knepley MATOP_DUPLICATE=34, 1390d519adbfSMatthew Knepley MATOP_FORWARD_SOLVE=35, 1391d519adbfSMatthew Knepley MATOP_BACKWARD_SOLVE=36, 1392d519adbfSMatthew Knepley MATOP_ILUFACTOR=37, 1393d519adbfSMatthew Knepley MATOP_ICCFACTOR=38, 1394d519adbfSMatthew Knepley MATOP_AXPY=39, 1395d519adbfSMatthew Knepley MATOP_GET_SUBMATRICES=40, 1396d519adbfSMatthew Knepley MATOP_INCREASE_OVERLAP=41, 1397d519adbfSMatthew Knepley MATOP_GET_VALUES=42, 1398d519adbfSMatthew Knepley MATOP_COPY=43, 1399d519adbfSMatthew Knepley MATOP_GET_ROW_MAX=44, 1400d519adbfSMatthew Knepley MATOP_SCALE=45, 1401d519adbfSMatthew Knepley MATOP_SHIFT=46, 140235153367SBarry Smith MATOP_DIAGONAL_SET=47, 14039d39f69dSJed Brown MATOP_ZERO_ROWS_COLUMNS=48, 14049d39f69dSJed Brown MATOP_SET_RANDOM=49, 1405d519adbfSMatthew Knepley MATOP_GET_ROW_IJ=50, 1406d519adbfSMatthew Knepley MATOP_RESTORE_ROW_IJ=51, 1407d519adbfSMatthew Knepley MATOP_GET_COLUMN_IJ=52, 1408d519adbfSMatthew Knepley MATOP_RESTORE_COLUMN_IJ=53, 1409d519adbfSMatthew Knepley MATOP_FDCOLORING_CREATE=54, 1410d519adbfSMatthew Knepley MATOP_COLORING_PATCH=55, 1411d519adbfSMatthew Knepley MATOP_SET_UNFACTORED=56, 1412d519adbfSMatthew Knepley MATOP_PERMUTE=57, 1413d519adbfSMatthew Knepley MATOP_SET_VALUES_BLOCKED=58, 1414d519adbfSMatthew Knepley MATOP_GET_SUBMATRIX=59, 1415d519adbfSMatthew Knepley MATOP_DESTROY=60, 1416d519adbfSMatthew Knepley MATOP_VIEW=61, 1417d519adbfSMatthew Knepley MATOP_CONVERT_FROM=62, 14187bab7c10SHong Zhang MATOP_MATMAT_MULT=63, 14197bab7c10SHong Zhang MATOP_MATMAT_MULT_SYMBOLIC=64, 14207bab7c10SHong Zhang MATOP_MATMAT_MULT_NUMERIC=65, 1421d519adbfSMatthew Knepley MATOP_SET_LOCAL_TO_GLOBAL_MAP=66, 1422d519adbfSMatthew Knepley MATOP_SET_VALUES_LOCAL=67, 1423d519adbfSMatthew Knepley MATOP_ZERO_ROWS_LOCAL=68, 1424d519adbfSMatthew Knepley MATOP_GET_ROW_MAX_ABS=69, 1425d519adbfSMatthew Knepley MATOP_GET_ROW_MIN_ABS=70, 1426d519adbfSMatthew Knepley MATOP_CONVERT=71, 1427d519adbfSMatthew Knepley MATOP_SET_COLORING=72, 14289d39f69dSJed Brown /* MATOP_PLACEHOLDER_73=73, */ 1429d519adbfSMatthew Knepley MATOP_SET_VALUES_ADIFOR=74, 1430d519adbfSMatthew Knepley MATOP_FD_COLORING_APPLY=75, 1431d519adbfSMatthew Knepley MATOP_SET_FROM_OPTIONS=76, 1432bda6c4cbSBarry Smith MATOP_MULT_CONSTRAINED=77, 1433bda6c4cbSBarry Smith MATOP_MULT_TRANSPOSE_CONSTRAIN=78, 14349d39f69dSJed Brown MATOP_FIND_ZERO_DIAGONALS=79, 1435d519adbfSMatthew Knepley MATOP_MULT_MULTIPLE=80, 1436d519adbfSMatthew Knepley MATOP_SOLVE_MULTIPLE=81, 1437d519adbfSMatthew Knepley MATOP_GET_INERTIA=82, 1438d519adbfSMatthew Knepley MATOP_LOAD=83, 1439d519adbfSMatthew Knepley MATOP_IS_SYMMETRIC=84, 1440d519adbfSMatthew Knepley MATOP_IS_HERMITIAN=85, 1441d519adbfSMatthew Knepley MATOP_IS_STRUCTURALLY_SYMMETRIC=86, 1442820469bcSHong Zhang MATOP_SET_VALUES_BLOCKEDLOCAL=87, 1443d519adbfSMatthew Knepley MATOP_GET_VECS=88, 1444d519adbfSMatthew Knepley MATOP_MAT_MULT=89, 1445d519adbfSMatthew Knepley MATOP_MAT_MULT_SYMBOLIC=90, 1446d519adbfSMatthew Knepley MATOP_MAT_MULT_NUMERIC=91, 1447d519adbfSMatthew Knepley MATOP_PTAP=92, 1448d519adbfSMatthew Knepley MATOP_PTAP_SYMBOLIC=93, 1449d519adbfSMatthew Knepley MATOP_PTAP_NUMERIC=94, 1450bda6c4cbSBarry Smith MATOP_MAT_TRANSPOSE_MULT=95, 1451bda6c4cbSBarry Smith MATOP_MAT_TRANSPOSE_MULT_SYMBO=96, 1452bda6c4cbSBarry Smith MATOP_MAT_TRANSPOSE_MULT_NUMER=97, 14539d39f69dSJed Brown /* MATOP_PLACEHOLDER_98=98, */ 14549d39f69dSJed Brown /* MATOP_PLACEHOLDER_99=99, */ 14559d39f69dSJed Brown /* MATOP_PLACEHOLDER_100=100, */ 14569d39f69dSJed Brown /* MATOP_PLACEHOLDER_101=101, */ 1457d519adbfSMatthew Knepley MATOP_CONJUGATE=102, 14589d39f69dSJed Brown /* MATOP_PLACEHOLDER_103=103, */ 1459d519adbfSMatthew Knepley MATOP_SET_VALUES_ROW=104, 1460d519adbfSMatthew Knepley MATOP_REAL_PART=105, 1461bda6c4cbSBarry Smith MATOP_IMAGINARY_PART=106, 1462bda6c4cbSBarry Smith MATOP_GET_ROW_UPPER_TRIANGULAR=107, 1463bda6c4cbSBarry Smith MATOP_RESTORE_ROW_UPPER_TRIANG=108, 1464bda6c4cbSBarry Smith MATOP_MAT_SOLVE=109, 14650e8d04f7SBarry Smith MATOP_GET_REDUNDANT_MATRIX=110, 1466d519adbfSMatthew Knepley MATOP_GET_ROW_MIN=111, 14670e8d04f7SBarry Smith MATOP_GET_COLUMN_VECTOR=112, 1468d519adbfSMatthew Knepley MATOP_MISSING_DIAGONAL=113, 14690e8d04f7SBarry Smith MATOP_GET_SEQ_NONZERO_STRUCTUR=114, 147089c6957cSBarry Smith MATOP_CREATE=115, 147189c6957cSBarry Smith MATOP_GET_GHOSTS=116, 14720e8d04f7SBarry Smith MATOP_GET_LOCAL_SUB_MATRIX=117, 14730e8d04f7SBarry Smith MATOP_RESTORE_LOCALSUB_MATRIX=118, 1474eeffb40dSHong Zhang MATOP_MULT_DIAGONAL_BLOCK=119, 14750e8d04f7SBarry Smith MATOP_HERMITIAN_TRANSPOSE=120, 14760e8d04f7SBarry Smith MATOP_MULT_HERMITIAN_TRANSPOSE=121, 14770e8d04f7SBarry Smith MATOP_MULT_HERMITIAN_TRANS_ADD=122, 14780e8d04f7SBarry Smith MATOP_GET_MULTI_PROC_BLOCK=123, 14799d39f69dSJed Brown MATOP_FIND_NONZERO_ROWS=124, 14800e8d04f7SBarry Smith MATOP_GET_COLUMN_NORMS=125, 14819d39f69dSJed Brown MATOP_INVERT_BLOCK_DIAGONAL=126, 14829d39f69dSJed Brown /* MATOP_PLACEHOLDER_127=127, */ 14830e8d04f7SBarry Smith MATOP_GET_SUB_MATRICES_PARALLE=128, 14842b8ad9a3SHong Zhang MATOP_SET_VALUES_BATCH=129, 14850e8d04f7SBarry Smith MATOP_TRANSPOSE_MAT_MULT=130, 14860e8d04f7SBarry Smith MATOP_TRANSPOSE_MAT_MULT_SYMBO=131, 1487e9b602ebSSatish Balay MATOP_TRANSPOSE_MAT_MULT_NUMER=132, 14880e8d04f7SBarry Smith MATOP_TRANSPOSE_COLORING_CREAT=133, 14890e8d04f7SBarry Smith MATOP_TRANS_COLORING_APPLY_SPT=134, 14900e8d04f7SBarry Smith MATOP_TRANS_COLORING_APPLY_DEN=135, 14910e8d04f7SBarry Smith MATOP_RART=136, 14920e8d04f7SBarry Smith MATOP_RART_SYMBOLIC=137, 14930e8d04f7SBarry Smith MATOP_RART_NUMERIC=138, 1494e09a3074SHong Zhang MATOP_SET_BLOCK_SIZES=139, 1495f9426fe0SMark Adams MATOP_AYPX=140, 14961919a2e2SJed Brown MATOP_RESIDUAL=141, 14979c8f2541SHong Zhang MATOP_FDCOLORING_SETUP=142, 14989c8f2541SHong Zhang MATOP_MPICONCATENATESEQ=144 1499fae171e0SBarry Smith } MatOperation; 1500014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatHasOperation(Mat,MatOperation,PetscBool *); 1501014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellSetOperation(Mat,MatOperation,void(*)(void)); 1502014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellGetOperation(Mat,MatOperation,void(**)(void)); 1503014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellSetContext(Mat,void*); 1504112a2221SBarry Smith 150590ace30eSBarry Smith /* 150690ace30eSBarry Smith Codes for matrices stored on disk. By default they are 150790ace30eSBarry Smith stored in a universal format. By changing the format with 15086a9046bcSBarry Smith PetscViewerPushFormat(viewer,PETSC_VIEWER_NATIVE); the matrices will 150990ace30eSBarry Smith be stored in a way natural for the matrix, for example dense matrices 151090ace30eSBarry Smith would be stored as dense. Matrices stored this way may only be 1511f4403165SShri Abhyankar read into matrices of the same type. 151290ace30eSBarry Smith */ 151390ace30eSBarry Smith #define MATRIX_BINARY_FORMAT_DENSE -1 151490ace30eSBarry Smith 1515014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetHashTableFactor(Mat,PetscReal); 1516014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatISGetLocalMat(Mat,Mat*); 1517014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatISSetLocalMat(Mat,Mat); 1518b7ce53b6SStefano Zampini PETSC_EXTERN PetscErrorCode MatISGetMPIXAIJ(Mat,MatReuse,Mat*); 15191f4e1ec7SBarry Smith 1520d9274352SBarry Smith /*S 1521d9274352SBarry Smith MatNullSpace - Object that removes a null space from a vector, i.e. 1522d9274352SBarry Smith orthogonalizes the vector to a subsapce 1523d9274352SBarry Smith 1524f7a9e4ceSBarry Smith Level: advanced 1525d9274352SBarry Smith 1526d9274352SBarry Smith Concepts: matrix; linear operator, null space 1527d9274352SBarry Smith 15286e1639daSBarry Smith Users manual sections: 15296e1639daSBarry Smith . sec_singular 15306e1639daSBarry Smith 1531d9274352SBarry Smith .seealso: MatNullSpaceCreate() 1532d9274352SBarry Smith S*/ 153374637425SBarry Smith typedef struct _p_MatNullSpace* MatNullSpace; 1534d9274352SBarry Smith 1535014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceCreate(MPI_Comm,PetscBool ,PetscInt,const Vec[],MatNullSpace*); 1536014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(MatNullSpace,Vec,void*),void*); 1537014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceDestroy(MatNullSpace*); 1538d0195637SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceRemove(MatNullSpace,Vec); 1539014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNullSpace(Mat, MatNullSpace *); 15405fa7ec2dSBarry Smith PETSC_EXTERN PetscErrorCode MatGetTransposeNullSpace(Mat, MatNullSpace *); 15415fa7ec2dSBarry Smith PETSC_EXTERN PetscErrorCode MatSetTransposeNullSpace(Mat,MatNullSpace); 1542014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNullSpace(Mat,MatNullSpace); 1543014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNearNullSpace(Mat,MatNullSpace); 1544014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNearNullSpace(Mat,MatNullSpace*); 1545014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceTest(MatNullSpace,Mat,PetscBool *); 1546014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceView(MatNullSpace,PetscViewer); 1547014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceGetVecs(MatNullSpace,PetscBool*,PetscInt*,const Vec**); 1548014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceCreateRigidBody(Vec,MatNullSpace*); 154974637425SBarry Smith 1550014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatReorderingSeqSBAIJ(Mat,IS); 1551014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetHashTableFactor(Mat,PetscReal); 1552014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetColumnIndices(Mat,PetscInt *); 1553014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJInvertBlockDiagonal(Mat); 15543f1d51d7SBarry Smith 1555014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMAIJ(Mat,PetscInt,Mat*); 1556014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMAIJRedimension(Mat,PetscInt,Mat*); 1557014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMAIJGetAIJ(Mat,Mat*); 1558c4f061fbSSatish Balay 1559014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatComputeExplicitOperator(Mat,Mat*); 1560b0a32e0cSBarry Smith 1561014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalScaleLocal(Mat,Vec); 156204f1ad80SBarry Smith 1563014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*); 1564014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetBase(Mat,Vec,Vec); 1565014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*); 1566014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*)); 1567014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec)); 1568014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt); 1569014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDResetHHistory(Mat); 1570014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctionError(Mat,PetscReal); 1571014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetPeriod(Mat,PetscInt); 1572014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDGetH(Mat,PetscScalar *); 1573014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetOptionsPrefix(Mat,const char[]); 1574014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*); 1575014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*); 1576e884886fSBarry Smith 15776370053bSBarry Smith /*S 15786370053bSBarry Smith MatMFFD - A data structured used to manage the computation of the h differencing parameter for matrix-free 15796370053bSBarry Smith Jacobian vector products 1580e884886fSBarry Smith 15816370053bSBarry Smith Notes: MATMFFD is a specific MatType which uses the MatMFFD data structure 15826370053bSBarry Smith 15836370053bSBarry Smith MatMFFD*() methods actually take the Mat as their first argument. Not a MatMFFD data structure 15846370053bSBarry Smith 15856370053bSBarry Smith Level: developer 15866370053bSBarry Smith 15876370053bSBarry Smith .seealso: MATMFFD, MatCreateMFFD(), MatMFFDSetFuction(), MatMFFDSetType(), MatMFFDRegister() 15886370053bSBarry Smith S*/ 1589e884886fSBarry Smith typedef struct _p_MatMFFD* MatMFFD; 1590e884886fSBarry Smith 159176bdecfbSBarry Smith /*J 1592e884886fSBarry Smith MatMFFDType - algorithm used to compute the h used in computing matrix-vector products via differencing of the function 1593e884886fSBarry Smith 1594e884886fSBarry Smith Level: beginner 1595e884886fSBarry Smith 1596e884886fSBarry Smith .seealso: MatMFFDSetType(), MatMFFDRegister() 159776bdecfbSBarry Smith J*/ 159819fd82e9SBarry Smith typedef const char* MatMFFDType; 1599e884886fSBarry Smith #define MATMFFD_DS "ds" 1600e884886fSBarry Smith #define MATMFFD_WP "wp" 1601e884886fSBarry Smith 160219fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDSetType(Mat,MatMFFDType); 1603bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDRegister(const char[],PetscErrorCode (*)(MatMFFD)); 1604e884886fSBarry Smith 1605014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDDSSetUmin(Mat,PetscReal); 1606014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDWPSetComputeNormU(Mat,PetscBool ); 1607e884886fSBarry Smith 160861ab5df0SBarry Smith PETSC_EXTERN PetscErrorCode MatFDColoringSetType(MatFDColoring,MatMFFDType); 160961ab5df0SBarry Smith 1610014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *); 1611014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *); 16127dbadf16SMatthew Knepley 161397969023SHong Zhang /* 161497969023SHong Zhang PETSc interface to MUMPS 161597969023SHong Zhang */ 161697969023SHong Zhang #ifdef PETSC_HAVE_MUMPS 1617014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMumpsSetIcntl(Mat,PetscInt,PetscInt); 1618bc6112feSHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetIcntl(Mat,PetscInt,PetscInt*); 1619b9e9ea26SSatish Balay PETSC_EXTERN PetscErrorCode MatMumpsSetCntl(Mat,PetscInt,PetscReal); 1620bc6112feSHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetCntl(Mat,PetscInt,PetscReal*); 1621bc6112feSHong Zhang 1622ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetInfo(Mat,PetscInt,PetscInt*); 1623ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetInfog(Mat,PetscInt,PetscInt*); 1624ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetRinfo(Mat,PetscInt,PetscReal*); 1625ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetRinfog(Mat,PetscInt,PetscReal*); 162697969023SHong Zhang #endif 162723a5497aSJed Brown 1628d954db57SHong Zhang /* 1629d615f992SSatish Balay PETSc interface to Mkl_Pardiso 1630b500a6b7SJose David Bermeol */ 1631b500a6b7SJose David Bermeol #ifdef PETSC_HAVE_MKL_PARDISO 1632d615f992SSatish Balay PETSC_EXTERN PetscErrorCode MatMkl_PardisoSetCntl(Mat,PetscInt,PetscInt); 1633b500a6b7SJose David Bermeol #endif 1634b500a6b7SJose David Bermeol 1635b500a6b7SJose David Bermeol /* 1636d305a81bSVasiliy Kozyrev PETSc interface to Mkl_CPardiso 1637d305a81bSVasiliy Kozyrev */ 1638d305a81bSVasiliy Kozyrev #ifdef PETSC_HAVE_MKL_CPARDISO 1639d305a81bSVasiliy Kozyrev PETSC_EXTERN PetscErrorCode MatMkl_CPardisoSetCntl(Mat,PetscInt,PetscInt); 1640d305a81bSVasiliy Kozyrev #endif 1641d305a81bSVasiliy Kozyrev 1642d305a81bSVasiliy Kozyrev /* 1643d954db57SHong Zhang PETSc interface to SUPERLU 1644d954db57SHong Zhang */ 1645d954db57SHong Zhang #ifdef PETSC_HAVE_SUPERLU 1646014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSuperluSetILUDropTol(Mat,PetscReal); 1647d954db57SHong Zhang #endif 1648d954db57SHong Zhang 1649fb785984SHong Zhang /* 1650fb785984SHong Zhang PETSc interface to SUPERLU_DIST 1651fb785984SHong Zhang */ 1652fb785984SHong Zhang #ifdef PETSC_HAVE_SUPERLU_DIST 1653fb785984SHong Zhang PETSC_EXTERN PetscErrorCode MatSuperluDistGetDiagU(Mat,PetscScalar*); 1654fb785984SHong Zhang #endif 1655fb785984SHong Zhang 1656575704cbSPieter Ghysels /* 1657575704cbSPieter Ghysels PETSc interface to STRUMPACK 1658575704cbSPieter Ghysels */ 1659575704cbSPieter Ghysels #ifdef PETSC_HAVE_STRUMPACK 1660575704cbSPieter Ghysels PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSRelCompTol(Mat,PetscReal); 1661575704cbSPieter Ghysels PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSMinSize(Mat,PetscInt); 1662575704cbSPieter Ghysels PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetColPerm(Mat,PetscBool); 1663575704cbSPieter Ghysels #endif 1664575704cbSPieter Ghysels 1665575704cbSPieter Ghysels 1666aa372e3fSPaul Mullowney #ifdef PETSC_HAVE_CUDA 16673f9c0db1SPaul Mullowney /*E 1668e057df02SPaul Mullowney MatCUSPARSEStorageFormat - indicates the storage format for CUSPARSE (GPU) 16692692e278SPaul Mullowney matrices. 1670e057df02SPaul Mullowney 1671e057df02SPaul Mullowney Not Collective 1672e057df02SPaul Mullowney 16738468deeeSKarl Rupp + MAT_CUSPARSE_CSR - Compressed Sparse Row 16742692e278SPaul Mullowney . MAT_CUSPARSE_ELL - Ellpack (requires CUDA 4.2 or later). 16752692e278SPaul Mullowney - MAT_CUSPARSE_HYB - Hybrid, a combination of Ellpack and Coordinate format (requires CUDA 4.2 or later). 1676e057df02SPaul Mullowney 1677e057df02SPaul Mullowney Level: intermediate 1678e057df02SPaul Mullowney 1679af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 1680e057df02SPaul Mullowney 1681e057df02SPaul Mullowney .seealso: MatCUSPARSESetFormat(), MatCUSPARSEFormatOperation 1682e057df02SPaul Mullowney E*/ 1683e057df02SPaul Mullowney 16843f9c0db1SPaul Mullowney typedef enum {MAT_CUSPARSE_CSR, MAT_CUSPARSE_ELL, MAT_CUSPARSE_HYB} MatCUSPARSEStorageFormat; 1685e057df02SPaul Mullowney 1686e057df02SPaul Mullowney /* these will be strings associated with enumerated type defined above */ 1687e057df02SPaul Mullowney PETSC_EXTERN const char *const MatCUSPARSEStorageFormats[]; 1688e057df02SPaul Mullowney 16893f9c0db1SPaul Mullowney /*E 1690e057df02SPaul Mullowney MatCUSPARSEFormatOperation - indicates the operation of CUSPARSE (GPU) 16912692e278SPaul Mullowney matrices whose operation should use a particular storage format. 1692e057df02SPaul Mullowney 1693e057df02SPaul Mullowney Not Collective 1694e057df02SPaul Mullowney 16958468deeeSKarl Rupp + MAT_CUSPARSE_MULT_DIAG - sets the storage format for the diagonal matrix in the parallel MatMult 16968468deeeSKarl Rupp . MAT_CUSPARSE_MULT_OFFDIAG - sets the storage format for the offdiagonal matrix in the parallel MatMult 16978468deeeSKarl Rupp . MAT_CUSPARSE_MULT - sets the storage format for the entire matrix in the serial (single GPU) MatMult 16988468deeeSKarl Rupp - MAT_CUSPARSE_ALL - sets the storage format for all CUSPARSE (GPU) matrices 1699e057df02SPaul Mullowney 1700e057df02SPaul Mullowney Level: intermediate 1701e057df02SPaul Mullowney 1702e057df02SPaul Mullowney .seealso: MatCUSPARSESetFormat(), MatCUSPARSEStorageFormat 1703e057df02SPaul Mullowney E*/ 170436d62e41SPaul Mullowney typedef enum {MAT_CUSPARSE_MULT_DIAG, MAT_CUSPARSE_MULT_OFFDIAG, MAT_CUSPARSE_MULT, MAT_CUSPARSE_ALL} MatCUSPARSEFormatOperation; 1705e057df02SPaul Mullowney 170610e9db80SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 170710e9db80SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCreateAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 1708e057df02SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCUSPARSESetFormat(Mat,MatCUSPARSEFormatOperation,MatCUSPARSEStorageFormat); 17099ae82921SPaul Mullowney #endif 17109ae82921SPaul Mullowney 171190c53902SBarry Smith #if defined(PETSC_HAVE_CUSP) 1712014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCUSP(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 1713014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateAIJCUSP(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 1714e057df02SPaul Mullowney 17153f9c0db1SPaul Mullowney /*E 1716e057df02SPaul Mullowney MatCUSPStorageFormat - indicates the storage format for CUSP (GPU) 171736d62e41SPaul Mullowney matrices. 1718e057df02SPaul Mullowney 1719e057df02SPaul Mullowney Not Collective 1720e057df02SPaul Mullowney 17218468deeeSKarl Rupp + MAT_CUSP_CSR - Compressed Sparse Row 17228468deeeSKarl Rupp . MAT_CUSP_DIA - Diagonal 17238468deeeSKarl Rupp - MAT_CUSP_ELL - Ellpack 1724e057df02SPaul Mullowney 1725e057df02SPaul Mullowney Level: intermediate 1726e057df02SPaul Mullowney 1727af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 1728e057df02SPaul Mullowney 1729e057df02SPaul Mullowney .seealso: MatCUSPSetFormat(), MatCUSPFormatOperation 1730e057df02SPaul Mullowney E*/ 17313f9c0db1SPaul Mullowney typedef enum {MAT_CUSP_CSR, MAT_CUSP_DIA, MAT_CUSP_ELL} MatCUSPStorageFormat; 1732e057df02SPaul Mullowney 1733e057df02SPaul Mullowney /* these will be strings associated with enumerated type defined above */ 1734e057df02SPaul Mullowney PETSC_EXTERN const char *const MatCUSPStorageFormats[]; 1735e057df02SPaul Mullowney 17363f9c0db1SPaul Mullowney /*E 1737e057df02SPaul Mullowney MatCUSPFormatOperation - indicates the operation of CUSP (GPU) 173836d62e41SPaul Mullowney matrices whose operation should use a particular storage format. 1739e057df02SPaul Mullowney 1740e057df02SPaul Mullowney Not Collective 1741e057df02SPaul Mullowney 17428468deeeSKarl Rupp + MAT_CUSP_MULT_DIAG - sets the storage format for the diagonal matrix in the parallel MatMult 17438468deeeSKarl Rupp . MAT_CUSP_MULT_OFFDIAG - sets the storage format for the offdiagonal matrix in the parallel MatMult 17448468deeeSKarl Rupp . MAT_CUSP_MULT - sets the storage format for the entire matrix in the serial (single GPU) MatMult 17458468deeeSKarl Rupp - MAT_CUSP_ALL - sets the storage format for all CUSP (GPU) matrices 1746e057df02SPaul Mullowney 1747e057df02SPaul Mullowney Level: intermediate 1748e057df02SPaul Mullowney 1749af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 1750e057df02SPaul Mullowney 1751e057df02SPaul Mullowney .seealso: MatCUSPSetFormat(), MatCUSPStorageFormat 1752e057df02SPaul Mullowney E*/ 175336d62e41SPaul Mullowney typedef enum {MAT_CUSP_MULT_DIAG, MAT_CUSP_MULT_OFFDIAG, MAT_CUSP_MULT, MAT_CUSP_ALL} MatCUSPFormatOperation; 1754e057df02SPaul Mullowney 1755e057df02SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCUSPSetFormat(Mat,MatCUSPFormatOperation,MatCUSPStorageFormat); 1756e057df02SPaul Mullowney #endif 175790c53902SBarry Smith 1758d67ff14aSKarl Rupp #if defined(PETSC_HAVE_VIENNACL) 1759d67ff14aSKarl Rupp PETSC_EXTERN PetscErrorCode MatCreateSeqAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 1760d67ff14aSKarl Rupp PETSC_EXTERN PetscErrorCode MatCreateAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 1761d67ff14aSKarl Rupp #endif 1762d67ff14aSKarl Rupp 176354efbe56SHong Zhang /* 176454efbe56SHong Zhang PETSc interface to FFTW 176554efbe56SHong Zhang */ 176654efbe56SHong Zhang #if defined(PETSC_HAVE_FFTW) 1767014dd563SJed Brown PETSC_EXTERN PetscErrorCode VecScatterPetscToFFTW(Mat,Vec,Vec); 1768014dd563SJed Brown PETSC_EXTERN PetscErrorCode VecScatterFFTWToPetsc(Mat,Vec,Vec); 17692a7a6963SBarry Smith PETSC_EXTERN PetscErrorCode MatCreateVecsFFTW(Mat,Vec*,Vec*,Vec*); 177054efbe56SHong Zhang #endif 177154efbe56SHong Zhang 1772014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateNest(MPI_Comm,PetscInt,const IS[],PetscInt,const IS[],const Mat[],Mat*); 1773014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSize(Mat,PetscInt*,PetscInt*); 1774014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetISs(Mat,IS[],IS[]); 1775014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetLocalISs(Mat,IS[],IS[]); 1776014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSubMats(Mat,PetscInt*,PetscInt*,Mat***); 1777014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSubMat(Mat,PetscInt,PetscInt,Mat*); 177819fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatNestSetVecType(Mat,VecType); 1779014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestSetSubMats(Mat,PetscInt,const IS[],PetscInt,const IS[],const Mat[]); 1780014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestSetSubMat(Mat,PetscInt,PetscInt,Mat); 1781d8588912SDave May 17824325cce7SMatthew G Knepley PETSC_EXTERN PetscErrorCode MatChop(Mat,PetscReal); 1783e0a58c10SMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatComputeBandwidth(Mat,PetscReal,PetscInt*); 17844325cce7SMatthew G Knepley 1785b541e6a4SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSubdomainsCreateCoalesce(Mat,PetscInt,PetscInt*,IS**); 178653dd7562SDmitry Karpeev 1787c094ef40SMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatPreallocatorPreallocate(Mat,PetscBool,Mat); 1788c094ef40SMatthew G. Knepley 1789539c167fSBarry Smith PETSC_INTERN PetscErrorCode MatHeaderMerge(Mat,Mat*); 1790539c167fSBarry Smith PETSC_EXTERN PetscErrorCode MatHeaderReplace(Mat,Mat*); 1791539c167fSBarry Smith 179223a5497aSJed Brown #endif 1793