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 145ceb03754SKris Buschelman MatReuse - Indicates if matrices obtained from a previous call to MatGetSubMatrices() 146d6eff37eSBarry Smith or MatGetSubMatrix() are to be reused to store the new matrix values. For MatConvert() is used to indicate 147d6eff37eSBarry Smith that the input matrix is to be replaced with the converted matrix. 148ceb03754SKris Buschelman 149ceb03754SKris Buschelman Level: beginner 150ceb03754SKris Buschelman 151af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 152ceb03754SKris Buschelman 1530c451bc4SBarry Smith .seealso: MatGetSubMatrices(), MatGetSubMatrix(), MatDestroyMatrices(), MatConvert() 154ceb03754SKris Buschelman E*/ 155511c6705SHong Zhang typedef enum {MAT_INITIAL_MATRIX,MAT_REUSE_MATRIX,MAT_IGNORE_MATRIX,MAT_INPLACE_MATRIX} MatReuse; 156ceb03754SKris Buschelman 1575494a064SHong Zhang /*E 1585494a064SHong Zhang MatGetSubMatrixOption - Indicates if matrices obtained from a call to MatGetSubMatrices() 159829201f2SHong Zhang include the matrix values. Currently it is only used by MatGetSeqNonzerostructure(). 1605494a064SHong Zhang 1615494a064SHong Zhang Level: beginner 1625494a064SHong Zhang 163829201f2SHong Zhang .seealso: MatGetSeqNonzerostructure() 1645494a064SHong Zhang E*/ 1655494a064SHong Zhang typedef enum {MAT_DO_NOT_GET_VALUES,MAT_GET_VALUES} MatGetSubMatrixOption; 1665494a064SHong Zhang 167607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatInitializePackage(void); 168c06d978dSMatthew Knepley 169014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreate(MPI_Comm,Mat*); 170014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetSizes(Mat,PetscInt,PetscInt,PetscInt,PetscInt); 17119fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatSetType(Mat,MatType); 172014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetFromOptions(Mat); 173685405a1SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatViewFromOptions(Mat A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);} 174bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatRegister(const char[],PetscErrorCode(*)(Mat)); 175014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRegisterBaseName(const char[],const char[],const char[]); 176014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetOptionsPrefix(Mat,const char[]); 177014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAppendOptionsPrefix(Mat,const char[]); 178014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOptionsPrefix(Mat,const char*[]); 17984d44b13SHong Zhang PETSC_EXTERN PetscErrorCode MatSetErrorIfFailure(Mat,PetscBool); 180f69a0ea3SMatthew Knepley 181140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatList; 182140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatColoringList; 183140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatPartitioningList; 18428988994SBarry Smith 1853b224e63SBarry Smith /*E 186aa6c7ce3SBarry Smith MatStructure - Indicates if two matrices have the same nonzero structure 1873b224e63SBarry Smith 1883b224e63SBarry Smith Level: beginner 1893b224e63SBarry Smith 190af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 1913b224e63SBarry Smith 192aa6c7ce3SBarry Smith .seealso: MatCopy(), MatAXPY() 1933b224e63SBarry Smith E*/ 194aa6c7ce3SBarry Smith typedef enum {DIFFERENT_NONZERO_PATTERN,SUBSET_NONZERO_PATTERN,SAME_NONZERO_PATTERN} MatStructure; 1953b224e63SBarry Smith 196014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqDense(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*); 197014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateDense(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*); 198014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 199014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 200014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *); 201014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithSplitArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],PetscInt[],PetscInt[],PetscScalar[],Mat*); 2028d7a6e47SBarry Smith 203014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 204014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 205014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat*); 206d21a29f3SJed Brown 207014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAdj(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscInt[],Mat*); 208014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 209d21a29f3SJed Brown 210014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 211014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPISBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *); 21238f409ebSLisandro Dalcin PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 213014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 2144cce697cSJed Brown PETSC_EXTERN PetscErrorCode MatXAIJSetPreallocation(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscInt[],const PetscInt[]); 215dfb205c3SBarry Smith 216014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateShell(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,void *,Mat*); 217014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateNormal(Mat,Mat*); 218c5e4d11fSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatCreateNormalHermitian(Mat,Mat*); 219014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateLRC(Mat,Mat,Mat,Mat*); 220c5e4d11fSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatCreateIS(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,ISLocalToGlobalMapping,ISLocalToGlobalMapping,Mat*); 221014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 222014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 223c0cdd4a1SDahai Guo 224014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateScatter(MPI_Comm,VecScatter,Mat*); 225014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScatterSetVecScatter(Mat,VecScatter); 226014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScatterGetVecScatter(Mat,VecScatter*); 227014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateBlockMat(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,Mat*); 228014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeAddMat(Mat,Mat); 229014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeMerge(Mat); 230014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateComposite(MPI_Comm,PetscInt,const Mat*,Mat*); 2316d7c1e57SBarry Smith typedef enum {MAT_COMPOSITE_ADDITIVE,MAT_COMPOSITE_MULTIPLICATIVE} MatCompositeType; 232014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeSetType(Mat,MatCompositeType); 2336d7c1e57SBarry Smith 23419fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatCreateFFT(MPI_Comm,PetscInt,const PetscInt[],MatType,Mat*); 235014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqCUFFT(MPI_Comm,PetscInt,const PetscInt[],Mat*); 236dedccee8SHong Zhang 237014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateTranspose(Mat,Mat*); 238d0de2241SAndrew Spott PETSC_EXTERN PetscErrorCode MatCreateHermitianTranspose(Mat,Mat*); 239014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSubMatrix(Mat,IS,IS,Mat*); 240014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSubMatrixUpdate(Mat,Mat,IS,IS); 241014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateLocalRef(Mat,IS,IS,Mat*); 2421472f72bSBarry Smith 243978814f1SStefano Zampini #if defined(PETSC_HAVE_HYPRE) 244*225daaf8SStefano Zampini PETSC_EXTERN PetscErrorCode MatCreateFromParCSR(void*,MatType,PetscCopyMode,Mat*); 245978814f1SStefano Zampini #endif 246978814f1SStefano Zampini 247014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPythonSetType(Mat,const char[]); 2481d6018f0SLisandro Dalcin 249014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetUp(Mat); 250014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroy(Mat*); 251e56f5c9eSBarry Smith PETSC_EXTERN PetscErrorCode MatGetNonzeroState(Mat,PetscObjectState*); 25221c89e3eSBarry Smith 253014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatConjugate(Mat); 254014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRealPart(Mat); 255014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatImaginaryPart(Mat); 256014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetDiagonalBlock(Mat,Mat*); 257014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetTrace(Mat,PetscScalar*); 258713ccfa9SJed Brown PETSC_EXTERN PetscErrorCode MatInvertBlockDiagonal(Mat,const PetscScalar **); 25999cafbc1SBarry Smith 2608ed539a5SBarry Smith /* ------------------------------------------------------------*/ 261014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 262014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 263014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesRow(Mat,PetscInt,const PetscScalar[]); 264014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesRowLocal(Mat,PetscInt,const PetscScalar[]); 265014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBatch(Mat,PetscInt,PetscInt,PetscInt[],const PetscScalar[]); 26673a71a0fSBarry Smith PETSC_EXTERN PetscErrorCode MatSetRandom(Mat,PetscRandom); 26784cb2905SBarry Smith 2682ef4de8bSBarry Smith /*S 2692ef4de8bSBarry Smith MatStencil - Data structure (C struct) for storing information about a single row or 27062ef0227SBarry Smith column of a matrix as indexed on an associated grid. These are arguments to MatSetStencil() and MatSetBlockStencil() 27162ef0227SBarry Smith 27262ef0227SBarry 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). 27362ef0227SBarry 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. 27462ef0227SBarry Smith 27562ef0227SBarry Smith For stencil access to vectors see DMDAVecGetArray(), DMDAVecGetArrayF90(). 276be479b30SJed Brown 277be479b30SJed Brown Fortran usage is different, see MatSetValuesStencil() for details. 2782ef4de8bSBarry Smith 2792ef4de8bSBarry Smith Level: beginner 2802ef4de8bSBarry Smith 2812ef4de8bSBarry Smith Concepts: matrix; linear operator 2822ef4de8bSBarry Smith 28362ef0227SBarry Smith .seealso: MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockedStencil(), DMDAVecGetArray(), DMDAVecGetArrayF90() 2842ef4de8bSBarry Smith S*/ 285435da068SBarry Smith typedef struct { 286c1ac3661SBarry Smith PetscInt k,j,i,c; 287435da068SBarry Smith } MatStencil; 2882ef4de8bSBarry Smith 289014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode); 290014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode); 291014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt); 292435da068SBarry Smith 293d91e6319SBarry Smith /*E 294d91e6319SBarry Smith MatAssemblyType - Indicates if the matrix is now to be used, or if you plan 29562ef0227SBarry Smith to continue to add or insert values to it 296d91e6319SBarry Smith 297d91e6319SBarry Smith Level: beginner 298d91e6319SBarry Smith 299d91e6319SBarry Smith .seealso: MatAssemblyBegin(), MatAssemblyEnd() 300d91e6319SBarry Smith E*/ 3016d4a8577SBarry Smith typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType; 302014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssemblyBegin(Mat,MatAssemblyType); 303014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssemblyEnd(Mat,MatAssemblyType); 304014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssembled(Mat,PetscBool *); 3054f9c727eSBarry Smith 3061d73ed98SBarry Smith 30730de9b25SBarry Smith 308d91e6319SBarry Smith /*E 309d91e6319SBarry Smith MatOption - Options that may be set for a matrix and its behavior or storage 310d91e6319SBarry Smith 311d91e6319SBarry Smith Level: beginner 312d91e6319SBarry Smith 313af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 314d91e6319SBarry Smith 315382164b0SBarry Smith Developer Notes: Entries that are negative need not be called collectively by all processes. 316382164b0SBarry Smith 317d91e6319SBarry Smith .seealso: MatSetOption() 318d91e6319SBarry Smith E*/ 319c5e4d11fSDmitry Karpeev typedef enum {MAT_OPTION_MIN = -3, 320c5e4d11fSDmitry Karpeev MAT_UNUSED_NONZERO_LOCATION_ERR = -2, 32192d4d306SBarry Smith MAT_ROW_ORIENTED = -1, 322382164b0SBarry Smith MAT_SYMMETRIC = 1, 323382164b0SBarry Smith MAT_STRUCTURALLY_SYMMETRIC = 2, 324382164b0SBarry Smith MAT_NEW_DIAGONALS = 3, 325382164b0SBarry Smith MAT_IGNORE_OFF_PROC_ENTRIES = 4, 326382164b0SBarry Smith MAT_USE_HASH_TABLE = 5, 327382164b0SBarry Smith MAT_KEEP_NONZERO_PATTERN = 6, 328382164b0SBarry Smith MAT_IGNORE_ZERO_ENTRIES = 7, 329382164b0SBarry Smith MAT_USE_INODES = 8, 330382164b0SBarry Smith MAT_HERMITIAN = 9, 331382164b0SBarry Smith MAT_SYMMETRY_ETERNAL = 10, 332c5e4d11fSDmitry Karpeev MAT_NEW_NONZERO_LOCATION_ERR = 11, 333382164b0SBarry Smith MAT_IGNORE_LOWER_TRIANGULAR = 12, 334382164b0SBarry Smith MAT_ERROR_LOWER_TRIANGULAR = 13, 335382164b0SBarry Smith MAT_GETROW_UPPERTRIANGULAR = 14, 336382164b0SBarry Smith MAT_SPD = 15, 33792d4d306SBarry Smith MAT_NO_OFF_PROC_ZERO_ROWS = 16, 33892d4d306SBarry Smith MAT_NO_OFF_PROC_ENTRIES = 17, 33992d4d306SBarry Smith MAT_NEW_NONZERO_LOCATIONS = 18, 340c5e4d11fSDmitry Karpeev MAT_NEW_NONZERO_ALLOCATION_ERR = 19, 341c5e4d11fSDmitry Karpeev MAT_SUBSET_OFF_PROC_ENTRIES = 20, 342c5e4d11fSDmitry Karpeev MAT_OPTION_MAX = 21} MatOption; 343e057df02SPaul Mullowney 344014dd563SJed Brown PETSC_EXTERN const char *MatOptions[]; 345014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetOption(Mat,MatOption,PetscBool); 3467d68702bSBarry Smith PETSC_EXTERN PetscErrorCode MatGetOption(Mat,MatOption,PetscBool*); 34719fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatGetType(Mat,MatType*); 34884cb2905SBarry Smith 349014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]); 350014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 351014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 352014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowUpperTriangular(Mat); 353014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRowUpperTriangular(Mat); 354014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnVector(Mat,Vec,PetscInt); 3558c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJGetArray(Mat,PetscScalar *[]); 3568c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJRestoreArray(Mat,PetscScalar *[]); 35721e72a00SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJGetMaxRowNonzeros(Mat,PetscInt*); 358bd04181cSBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJSetValuesLocalFast(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 3598c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatDenseGetArray(Mat,PetscScalar *[]); 3608c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatDenseRestoreArray(Mat,PetscScalar *[]); 361014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBlockSize(Mat,PetscInt *); 362014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSize(Mat,PetscInt); 363014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBlockSizes(Mat,PetscInt *,PetscInt *); 364014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSizes(Mat,PetscInt,PetscInt); 36533d57670SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSizesFromMats(Mat,Mat,Mat); 3661620fd73SBarry Smith 367014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMult(Mat,Vec,Vec); 368014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultDiagonalBlock(Mat,Vec,Vec); 369014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultAdd(Mat,Vec,Vec,Vec); 370014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTranspose(Mat,Vec,Vec); 371014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultHermitianTranspose(Mat,Vec,Vec); 372014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsTranspose(Mat,Mat,PetscReal,PetscBool *); 373014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitianTranspose(Mat,Mat,PetscReal,PetscBool *); 374014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeAdd(Mat,Vec,Vec,Vec); 375014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultHermitianTransposeAdd(Mat,Vec,Vec,Vec); 376014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultConstrained(Mat,Vec,Vec); 377014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeConstrained(Mat,Vec,Vec); 378014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatSolve(Mat,Mat,Mat); 379f9426fe0SMark Adams PETSC_EXTERN PetscErrorCode MatResidual(Mat,Vec,Vec,Vec); 380f5edf698SHong Zhang 381d91e6319SBarry Smith /*E 382d91e6319SBarry Smith MatDuplicateOption - Indicates if a duplicated sparse matrix should have 383d91e6319SBarry Smith its numerical values copied over or just its nonzero structure. 384d91e6319SBarry Smith 385d91e6319SBarry Smith Level: beginner 386d91e6319SBarry Smith 387af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 388d91e6319SBarry Smith 38970dcbbb9SBarry Smith $ MAT_SHARE_NONZERO_PATTERN - the i and j arrays in the new matrix will be shared with the original matrix 39070dcbbb9SBarry Smith $ this also triggers the MAT_DO_NOT_COPY_VALUES option. This is used when you 39170dcbbb9SBarry Smith $ have several matrices with the same nonzero pattern. 39270dcbbb9SBarry Smith 393d91e6319SBarry Smith .seealso: MatDuplicate() 394d91e6319SBarry Smith E*/ 39570dcbbb9SBarry Smith typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES,MAT_SHARE_NONZERO_PATTERN} MatDuplicateOption; 3962e8a6d31SBarry Smith 39719fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatConvert(Mat,MatType,MatReuse,Mat*); 398014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDuplicate(Mat,MatDuplicateOption,Mat*); 39994a9d846SBarry Smith 400cb5b572fSBarry Smith 401014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCopy(Mat,Mat,MatStructure); 402014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatView(Mat,PetscViewer); 403014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsSymmetric(Mat,PetscReal,PetscBool *); 404014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsStructurallySymmetric(Mat,PetscBool *); 405014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitian(Mat,PetscReal,PetscBool *); 406014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsSymmetricKnown(Mat,PetscBool *,PetscBool *); 407014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitianKnown(Mat,PetscBool *,PetscBool *); 408014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMissingDiagonal(Mat,PetscBool *,PetscInt *); 409014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLoad(Mat, PetscViewer); 4107b80b807SBarry Smith 4111a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool *); 4121a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool *); 4131a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool *); 4141a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool *); 415d4fbbf0eSBarry Smith 416d91e6319SBarry Smith /*S 417d91e6319SBarry Smith MatInfo - Context of matrix information, used with MatGetInfo() 418d91e6319SBarry Smith 419d91e6319SBarry Smith In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE 420d91e6319SBarry Smith 421d91e6319SBarry Smith Level: intermediate 422d91e6319SBarry Smith 423d91e6319SBarry Smith Concepts: matrix^nonzero information 424d91e6319SBarry Smith 425d9274352SBarry Smith .seealso: MatGetInfo(), MatInfoType 426d91e6319SBarry Smith S*/ 4274e220ebcSLois Curfman McInnes typedef struct { 428b0a32e0cSBarry Smith PetscLogDouble block_size; /* block size */ 429b0a32e0cSBarry Smith PetscLogDouble nz_allocated,nz_used,nz_unneeded; /* number of nonzeros */ 430b0a32e0cSBarry Smith PetscLogDouble memory; /* memory allocated */ 431b0a32e0cSBarry Smith PetscLogDouble assemblies; /* number of matrix assemblies called */ 432b0a32e0cSBarry Smith PetscLogDouble mallocs; /* number of mallocs during MatSetValues() */ 433b0a32e0cSBarry Smith PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */ 434b0a32e0cSBarry Smith PetscLogDouble factor_mallocs; /* number of mallocs during factorization */ 4354e220ebcSLois Curfman McInnes } MatInfo; 4364e220ebcSLois Curfman McInnes 437d9274352SBarry Smith /*E 438d9274352SBarry Smith MatInfoType - Indicates if you want information about the local part of the matrix, 439d9274352SBarry Smith the entire parallel matrix or the maximum over all the local parts. 440d9274352SBarry Smith 441d9274352SBarry Smith Level: beginner 442d9274352SBarry Smith 443af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 444d9274352SBarry Smith 445d9274352SBarry Smith .seealso: MatGetInfo(), MatInfo 446d9274352SBarry Smith E*/ 4477b80b807SBarry Smith typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType; 448014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetInfo(Mat,MatInfoType,MatInfo*); 449014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetDiagonal(Mat,Vec); 450014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMax(Mat,Vec,PetscInt[]); 451014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMin(Mat,Vec,PetscInt[]); 452014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMaxAbs(Mat,Vec,PetscInt[]); 453014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMinAbs(Mat,Vec,PetscInt[]); 454014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowSum(Mat,Vec); 455014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTranspose(Mat,MatReuse,Mat*); 456014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatHermitianTranspose(Mat,MatReuse,Mat*); 457014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPermute(Mat,IS,IS,Mat*); 458014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalScale(Mat,Vec,Vec); 459014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalSet(Mat,Vec,InsertMode); 460a52676ecSHong Zhang 461014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatEqual(Mat,Mat,PetscBool*); 462014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultEqual(Mat,Mat,PetscInt,PetscBool*); 463014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultAddEqual(Mat,Mat,PetscInt,PetscBool*); 464014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeEqual(Mat,Mat,PetscInt,PetscBool*); 465014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscBool*); 466a52676ecSHong Zhang PETSC_EXTERN PetscErrorCode MatMatMultEqual(Mat,Mat,Mat,PetscInt,PetscBool*); 467a52676ecSHong Zhang PETSC_EXTERN PetscErrorCode MatTransposeMatMultEqual(Mat,Mat,Mat,PetscInt,PetscBool*); 4687b80b807SBarry Smith 469014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNorm(Mat,NormType,PetscReal*); 470014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnNorms(Mat,NormType,PetscReal*); 471014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroEntries(Mat); 472014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 473014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsIS(Mat,IS,PetscScalar,Vec,Vec); 474014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsStencil(Mat,PetscInt,const MatStencil [],PetscScalar,Vec,Vec); 475014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsStencil(Mat,PetscInt,const MatStencil[],PetscScalar,Vec,Vec); 476014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumns(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 477014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsIS(Mat,IS,PetscScalar,Vec,Vec); 4787b80b807SBarry Smith 479014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSize(Mat,PetscInt*,PetscInt*); 480014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalSize(Mat,PetscInt*,PetscInt*); 481014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRange(Mat,PetscInt*,PetscInt*); 482014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRanges(Mat,const PetscInt**); 483014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRangeColumn(Mat,PetscInt*,PetscInt*); 484014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRangesColumn(Mat,const PetscInt**); 485566876eaSJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipIS(Mat,IS*,IS*); 4865ef9f2a5SBarry Smith 487014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]); 48853dd7562SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatGetSubMatricesMPI(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]); 489014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroyMatrices(PetscInt,Mat *[]); 490014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSubMatrix(Mat,IS,IS,MatReuse,Mat *); 491014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalSubMatrix(Mat,IS,IS,Mat*); 492014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreLocalSubMatrix(Mat,IS,IS,Mat*); 493014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSeqNonzeroStructure(Mat,Mat*); 494014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroySeqNonzeroStructure(Mat*); 4957b80b807SBarry Smith 496014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJ(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*); 497014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJSymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*); 498014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJNumeric(Mat,Mat); 499014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMat(Mat,MatReuse,Mat*); 500014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*); 501014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,Mat*); 502014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetGhosts(Mat, PetscInt *,const PetscInt *[]); 5038efafbd8SBarry Smith 504014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt); 505aa1e27eaSFande Kong PETSC_EXTERN PetscErrorCode MatIncreaseOverlapSplit(Mat mat,PetscInt n,IS is[],PetscInt ov); 506d2b2a242SBarry Smith PETSC_EXTERN PetscErrorCode MatMPIAIJSetUseScalableIncreaseOverlap(Mat,PetscBool); 5077b80b807SBarry Smith 508014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*); 509014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMultSymbolic(Mat,Mat,PetscReal,Mat*); 510014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMultNumeric(Mat,Mat,Mat); 51122440eb1SKris Buschelman 5127bab7c10SHong Zhang PETSC_EXTERN PetscErrorCode MatMatMatMult(Mat,Mat,Mat,MatReuse,PetscReal,Mat*); 5137bab7c10SHong Zhang 514014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*); 515014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAPSymbolic(Mat,Mat,PetscReal,Mat*); 516014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAPNumeric(Mat,Mat,Mat); 517014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARt(Mat,Mat,MatReuse,PetscReal,Mat*); 518014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARtSymbolic(Mat,Mat,PetscReal,Mat*); 519014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARtNumeric(Mat,Mat,Mat); 52022440eb1SKris Buschelman 521014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeMatMult(Mat,Mat,MatReuse,PetscReal,Mat*); 522014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposetMatMultSymbolic(Mat,Mat,PetscReal,Mat*); 523014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposetMatMultNumeric(Mat,Mat,Mat); 524014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMult(Mat,Mat,MatReuse,PetscReal,Mat*); 525014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMultSymbolic(Mat,Mat,PetscReal,Mat*); 526014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMultNumeric(Mat,Mat,Mat); 527bc011b1eSHong Zhang 528014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAXPY(Mat,PetscScalar,Mat,MatStructure); 529014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAYPX(Mat,PetscScalar,Mat,MatStructure); 5301c741599SHong Zhang 531014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScale(Mat,PetscScalar); 532014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShift(Mat,PetscScalar); 5337b80b807SBarry Smith 534014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping); 535014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*); 536a93ff8c4SPeter Brune PETSC_EXTERN PetscErrorCode MatGetLayouts(Mat,PetscLayout*,PetscLayout*); 537014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 538014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsLocalIS(Mat,IS,PetscScalar,Vec,Vec); 539014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 540014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocalIS(Mat,IS,PetscScalar,Vec,Vec); 541014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 542014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 543052efed2SBarry Smith 544014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStashSetInitialSize(Mat,PetscInt,PetscInt); 545014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*); 54690f02eecSBarry Smith 547014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInterpolate(Mat,Vec,Vec); 548014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInterpolateAdd(Mat,Vec,Vec,Vec); 549014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestrict(Mat,Vec,Vec); 5502a7a6963SBarry Smith PETSC_EXTERN PetscErrorCode MatCreateVecs(Mat,Vec*,Vec*); 55184cc2a5aSBarry Smith PETSC_DEPRECATED("Use MatCreateVecs()") PETSC_STATIC_INLINE PetscErrorCode MatGetVecs(Mat mat,Vec *x,Vec *y) {return MatCreateVecs(mat,x,y);} 55253cd1579SHong Zhang PETSC_EXTERN PetscErrorCode MatCreateRedundantMatrix(Mat,PetscInt,MPI_Comm,MatReuse,Mat*); 553014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetMultiProcBlock(Mat,MPI_Comm,MatReuse,Mat*); 554014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFindZeroDiagonals(Mat,IS*); 5553a062f41SBarry Smith PETSC_EXTERN PetscErrorCode MatFindOffBlockDiagonalEntries(Mat,IS*); 5569c8f2541SHong Zhang PETSC_EXTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat(MPI_Comm,Mat,PetscInt,MatReuse,Mat*); 557bd481603SBarry Smith 558bd481603SBarry Smith /*MC 5591620fd73SBarry Smith MatSetValue - Set a single entry into a matrix. 5601620fd73SBarry Smith 5611620fd73SBarry Smith Not collective 5621620fd73SBarry Smith 563a9834a7dSSatish Balay Synopsis: 564a9834a7dSSatish Balay #include <petscmat.h> 565a9834a7dSSatish Balay PetscErrorCode MatSetValue(Mat m,PetscInt row,PetscInt col,PetscScalar value,InsertMode mode) 566a9834a7dSSatish Balay 5671620fd73SBarry Smith Input Parameters: 5681620fd73SBarry Smith + m - the matrix 5691620fd73SBarry Smith . row - the row location of the entry 5701620fd73SBarry Smith . col - the column location of the entry 5711620fd73SBarry Smith . value - the value to insert 5721620fd73SBarry Smith - mode - either INSERT_VALUES or ADD_VALUES 5731620fd73SBarry Smith 5741620fd73SBarry Smith Notes: 5751620fd73SBarry Smith For efficiency one should use MatSetValues() and set several or many 5761620fd73SBarry Smith values simultaneously if possible. 5771620fd73SBarry Smith 5781620fd73SBarry Smith Level: beginner 5791620fd73SBarry Smith 5801620fd73SBarry Smith .seealso: MatSetValues(), MatSetValueLocal() 5811620fd73SBarry Smith M*/ 5821620fd73SBarry 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);} 5831620fd73SBarry Smith 5841620fd73SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatGetValue(Mat v,PetscInt i,PetscInt j,PetscScalar *va) {return MatGetValues(v,1,&i,1,&j,va);} 5851620fd73SBarry Smith 5861620fd73SBarry 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);} 5871620fd73SBarry Smith 5881620fd73SBarry Smith /*MC 5890d2a7ff2SSatish Balay MatPreallocateInitialize - Begins the block of code that will count the number of nonzeros per 590bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 591bd481603SBarry Smith 592bd481603SBarry Smith Synopsis: 593aaa7dc30SBarry Smith #include <petscmat.h> 594c1ac3661SBarry Smith PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz) 595bd481603SBarry Smith 596bd481603SBarry Smith Collective on MPI_Comm 597bd481603SBarry Smith 598bd481603SBarry Smith Input Parameters: 599bd481603SBarry Smith + comm - the communicator that will share the eventually allocated matrix 600859fcb39SBarry Smith . nrows - the number of LOCAL rows in the matrix 601859fcb39SBarry Smith - ncols - the number of LOCAL columns in the matrix 602bd481603SBarry Smith 603bd481603SBarry Smith Output Parameters: 604bd481603SBarry Smith + dnz - the array that will be passed to the matrix preallocation routines 605bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 606bd481603SBarry Smith 607bd481603SBarry Smith Level: intermediate 608bd481603SBarry Smith 609bd481603SBarry Smith Notes: 610a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 611bd481603SBarry Smith 6121d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 613bd481603SBarry Smith 614bd481603SBarry Smith Use MatPreallocateInitializeSymmetric() for symmetric matrices (MPISBAIJ matrices) 615bd481603SBarry Smith 6161620fd73SBarry Smith This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize(). 6171620fd73SBarry Smith 618bd481603SBarry Smith Concepts: preallocation^Matrix 619bd481603SBarry Smith 620d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(), 621d6e23781SBarry Smith MatPreallocateInitializeSymmetric(), MatPreallocateSymmetricSetLocalBlock() 622bd481603SBarry Smith M*/ 623c4f061fbSSatish Balay #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \ 6247c922b88SBarry Smith { \ 625a89bc211SBarry Smith PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__start,__end; \ 6261795a4d1SJed Brown _4_ierr = PetscCalloc2(__nrows,&dnz,__nrows,&onz);CHKERRQ(_4_ierr); \ 6271795a4d1SJed Brown __start = 0; __end = __start; \ 628c1ac3661SBarry Smith _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ctmp;\ 629a89bc211SBarry Smith _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows; 6307c922b88SBarry Smith 631bd481603SBarry Smith /*MC 632bd481603SBarry Smith MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be 633bd481603SBarry Smith inserted using a local number of the rows and columns 634bd481603SBarry Smith 635bd481603SBarry Smith Synopsis: 636aaa7dc30SBarry Smith #include <petscmat.h> 637c1ac3661SBarry Smith PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 638bd481603SBarry Smith 639bd481603SBarry Smith Not Collective 640bd481603SBarry Smith 641bd481603SBarry Smith Input Parameters: 642784ac674SJed Brown + map - the row mapping from local numbering to global numbering 643bd481603SBarry Smith . nrows - the number of rows indicated 6441d73ed98SBarry Smith . rows - the indices of the rows 645784ac674SJed Brown . cmap - the column mapping from local to global numbering 646bd481603SBarry Smith . ncols - the number of columns in the matrix 647bd481603SBarry Smith . cols - the columns indicated 648bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 649bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 650bd481603SBarry Smith 651bd481603SBarry Smith Level: intermediate 652bd481603SBarry Smith 653bd481603SBarry Smith Notes: 654a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 655bd481603SBarry Smith 6561d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 657bd481603SBarry Smith 658bd481603SBarry Smith Concepts: preallocation^Matrix 659bd481603SBarry Smith 660d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(), 661d6e23781SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock() 662bd481603SBarry Smith M*/ 663784ac674SJed Brown #define MatPreallocateSetLocal(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \ 664c4f061fbSSatish Balay {\ 665c1ac3661SBarry Smith PetscInt __l;\ 666784ac674SJed Brown _4_ierr = ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\ 667784ac674SJed Brown _4_ierr = ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\ 668c4f061fbSSatish Balay for (__l=0;__l<nrows;__l++) {\ 669ef66eb69SBarry Smith _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 670c4f061fbSSatish Balay }\ 671c4f061fbSSatish Balay } 672c4f061fbSSatish Balay 673bd481603SBarry Smith /*MC 674d6e23781SBarry Smith MatPreallocateSetLocalBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be 675bd481603SBarry Smith inserted using a local number of the rows and columns 676bd481603SBarry Smith 677bd481603SBarry Smith Synopsis: 678aaa7dc30SBarry Smith #include <petscmat.h> 679d6e23781SBarry Smith PetscErrorCode MatPreallocateSetLocalBlock(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 680d6e23781SBarry Smith 681d6e23781SBarry Smith Not Collective 682d6e23781SBarry Smith 683d6e23781SBarry Smith Input Parameters: 684d6e23781SBarry Smith + map - the row mapping from local numbering to global numbering 685d6e23781SBarry Smith . nrows - the number of rows indicated 686d6e23781SBarry Smith . rows - the indices of the rows 687d6e23781SBarry Smith . cmap - the column mapping from local to global numbering 688d6e23781SBarry Smith . ncols - the number of columns in the matrix 689d6e23781SBarry Smith . cols - the columns indicated 690d6e23781SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 691d6e23781SBarry Smith - ozn - the other array passed to the matrix preallocation routines 692d6e23781SBarry Smith 693d6e23781SBarry Smith Level: intermediate 694d6e23781SBarry Smith 695d6e23781SBarry Smith Notes: 696d6e23781SBarry Smith See Users-Manual: ch_performance for more details. 697d6e23781SBarry Smith 698d6e23781SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 699d6e23781SBarry Smith 700d6e23781SBarry Smith Concepts: preallocation^Matrix 701d6e23781SBarry Smith 702d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(), 703d6e23781SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock() 704d6e23781SBarry Smith M*/ 705d6e23781SBarry Smith #define MatPreallocateSetLocalBlock(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \ 706d6e23781SBarry Smith {\ 707d6e23781SBarry Smith PetscInt __l;\ 708d6e23781SBarry Smith _4_ierr = ISLocalToGlobalMappingApplyBlock(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\ 709d6e23781SBarry Smith _4_ierr = ISLocalToGlobalMappingApplyBlock(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\ 710d6e23781SBarry Smith for (__l=0;__l<nrows;__l++) {\ 711d6e23781SBarry Smith _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 712d6e23781SBarry Smith }\ 713d6e23781SBarry Smith } 714d6e23781SBarry Smith 715d6e23781SBarry Smith /*MC 716d6e23781SBarry Smith MatPreallocateSymmetricSetLocalBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be 717d6e23781SBarry Smith inserted using a local number of the rows and columns 718d6e23781SBarry Smith 719d6e23781SBarry Smith Synopsis: 720d6e23781SBarry Smith #include <petscmat.h> 721d6e23781SBarry Smith PetscErrorCode MatPreallocateSymmetricSetLocalBlock(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 722bd481603SBarry Smith 723bd481603SBarry Smith Not Collective 724bd481603SBarry Smith 725bd481603SBarry Smith Input Parameters: 726bd481603SBarry Smith + map - the mapping between local numbering and global numbering 727bd481603SBarry Smith . nrows - the number of rows indicated 7281d73ed98SBarry Smith . rows - the indices of the rows 729bd481603SBarry Smith . ncols - the number of columns in the matrix 730bd481603SBarry Smith . cols - the columns indicated 731bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 732bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 733bd481603SBarry Smith 734bd481603SBarry Smith Level: intermediate 735bd481603SBarry Smith 736bd481603SBarry Smith Notes: 737a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 738bd481603SBarry Smith 739bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 740bd481603SBarry Smith 741bd481603SBarry Smith Concepts: preallocation^Matrix 742bd481603SBarry Smith 743d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateInitialize(), 744d6e23781SBarry Smith MatPreallocateInitialize(), MatPreallocateSetLocal() 745bd481603SBarry Smith M*/ 746d6e23781SBarry Smith #define MatPreallocateSymmetricSetLocalBlock(map,nrows,rows,ncols,cols,dnz,onz) 0;\ 747d3d32019SBarry Smith {\ 748c1ac3661SBarry Smith PetscInt __l;\ 749d6e23781SBarry Smith _4_ierr = ISLocalToGlobalMappingApplyBlock(map,nrows,rows,rows);CHKERRQ(_4_ierr);\ 750d6e23781SBarry Smith _4_ierr = ISLocalToGlobalMappingApplyBlock(map,ncols,cols,cols);CHKERRQ(_4_ierr);\ 751d3d32019SBarry Smith for (__l=0;__l<nrows;__l++) {\ 752d6e23781SBarry Smith _4_ierr = MatPreallocateSymmetricSetBlock((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 753d3d32019SBarry Smith }\ 754d3d32019SBarry Smith } 755bd481603SBarry Smith /*MC 756bd481603SBarry Smith MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be 757bd481603SBarry Smith inserted using a local number of the rows and columns 758bd481603SBarry Smith 759bd481603SBarry Smith Synopsis: 760aaa7dc30SBarry Smith #include <petscmat.h> 761c1ac3661SBarry Smith PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 762bd481603SBarry Smith 763bd481603SBarry Smith Not Collective 764bd481603SBarry Smith 765bd481603SBarry Smith Input Parameters: 76664075487SBarry Smith + row - the row 767bd481603SBarry Smith . ncols - the number of columns in the matrix 768a50f8bf6SHong Zhang - cols - the columns indicated 769a50f8bf6SHong Zhang 770a50f8bf6SHong Zhang Output Parameters: 771a50f8bf6SHong Zhang + dnz - the array that will be passed to the matrix preallocation routines 772bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 773bd481603SBarry Smith 774bd481603SBarry Smith Level: intermediate 775bd481603SBarry Smith 776bd481603SBarry Smith Notes: 777a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 778bd481603SBarry Smith 779bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 780bd481603SBarry Smith 7811620fd73SBarry Smith This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize(). 7821620fd73SBarry Smith 783bd481603SBarry Smith Concepts: preallocation^Matrix 784bd481603SBarry Smith 785d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(), 786d6e23781SBarry Smith MatPreallocateInitialize(), MatPreallocateSetLocal() 787bd481603SBarry Smith M*/ 788c4f061fbSSatish Balay #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\ 789c1ac3661SBarry Smith { PetscInt __i; \ 790e32f2f54SBarry 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);\ 791e32f2f54SBarry 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);\ 7927c922b88SBarry Smith for (__i=0; __i<nc; __i++) {\ 79364075487SBarry Smith if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \ 7947cc688d7SBarry Smith else dnz[row - __rstart]++;\ 7957c922b88SBarry Smith }\ 7967c922b88SBarry Smith } 7977c922b88SBarry Smith 798bd481603SBarry Smith /*MC 799d6e23781SBarry Smith MatPreallocateSymmetricSetBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be 800bd481603SBarry Smith inserted using a local number of the rows and columns 801bd481603SBarry Smith 802bd481603SBarry Smith Synopsis: 803aaa7dc30SBarry Smith #include <petscmat.h> 804d6e23781SBarry Smith PetscErrorCode MatPreallocateSymmetricSetBlock(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 805bd481603SBarry Smith 806bd481603SBarry Smith Not Collective 807bd481603SBarry Smith 808bd481603SBarry Smith Input Parameters: 809bd481603SBarry Smith + nrows - the number of rows indicated 8101d73ed98SBarry Smith . rows - the indices of the rows 811bd481603SBarry Smith . ncols - the number of columns in the matrix 812bd481603SBarry Smith . cols - the columns indicated 813bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 814bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 815bd481603SBarry Smith 816bd481603SBarry Smith Level: intermediate 817bd481603SBarry Smith 818bd481603SBarry Smith Notes: 819a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 820bd481603SBarry Smith 821bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 822bd481603SBarry Smith 8231620fd73SBarry Smith This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize(). 8241620fd73SBarry Smith 825bd481603SBarry Smith Concepts: preallocation^Matrix 826bd481603SBarry Smith 827d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateInitialize(), 828d6e23781SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocal() 829bd481603SBarry Smith M*/ 830d6e23781SBarry Smith #define MatPreallocateSymmetricSetBlock(row,nc,cols,dnz,onz) 0;\ 831c1ac3661SBarry Smith { PetscInt __i; \ 832d3d32019SBarry Smith for (__i=0; __i<nc; __i++) {\ 833d3d32019SBarry Smith if (cols[__i] >= __end) onz[row - __rstart]++; \ 834d3d32019SBarry Smith else if (cols[__i] >= row) dnz[row - __rstart]++;\ 835d3d32019SBarry Smith }\ 836d3d32019SBarry Smith } 837d3d32019SBarry Smith 838bd481603SBarry Smith /*MC 83916371a99SBarry Smith MatPreallocateLocation - An alternative to MatPreallocationSet() that puts the nonzero locations into the matrix if it exists 84016371a99SBarry Smith 84116371a99SBarry Smith Synopsis: 842aaa7dc30SBarry Smith #include <petscmat.h> 84316371a99SBarry Smith PetscErrorCode MatPreallocateLocations(Mat A,PetscInt row,PetscInt ncols,PetscInt *cols,PetscInt *dnz,PetscInt *onz) 84416371a99SBarry Smith 84516371a99SBarry Smith Not Collective 84616371a99SBarry Smith 84716371a99SBarry Smith Input Parameters: 84816371a99SBarry Smith . A - matrix 84916371a99SBarry Smith . row - row where values exist (must be local to this process) 85016371a99SBarry Smith . ncols - number of columns 85116371a99SBarry Smith . cols - columns with nonzeros 85216371a99SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 85316371a99SBarry Smith - ozn - the other array passed to the matrix preallocation routines 85416371a99SBarry Smith 85516371a99SBarry Smith Level: intermediate 85616371a99SBarry Smith 85716371a99SBarry Smith Notes: 858a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 85916371a99SBarry Smith 86016371a99SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 86116371a99SBarry Smith 86216371a99SBarry Smith This is a MACRO not a function because it uses a bunch of variables private to the MatPreallocation.... routines. 86316371a99SBarry Smith 86416371a99SBarry Smith Concepts: preallocation^Matrix 86516371a99SBarry Smith 866d6e23781SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(), 867d6e23781SBarry Smith MatPreallocateSymmetricSetLocalBlock() 86816371a99SBarry Smith M*/ 8690298fd71SBarry 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);} 87016371a99SBarry Smith 87116371a99SBarry Smith 87216371a99SBarry Smith /*MC 8730d2a7ff2SSatish Balay MatPreallocateFinalize - Ends the block of code that will count the number of nonzeros per 874bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 875bd481603SBarry Smith 876bd481603SBarry Smith Synopsis: 877aaa7dc30SBarry Smith #include <petscmat.h> 878c1ac3661SBarry Smith PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz) 879bd481603SBarry Smith 880bd481603SBarry Smith Collective on MPI_Comm 881bd481603SBarry Smith 882bd481603SBarry Smith Input Parameters: 88316371a99SBarry Smith + dnz - the array that was be passed to the matrix preallocation routines 884bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 885bd481603SBarry Smith 886bd481603SBarry Smith Level: intermediate 887bd481603SBarry Smith 888bd481603SBarry Smith Notes: 889a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 890bd481603SBarry Smith 891bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 892bd481603SBarry Smith 8931620fd73SBarry Smith This is a MACRO not a function because it closes the { started in MatPreallocateInitialize(). 8941620fd73SBarry Smith 895bd481603SBarry Smith Concepts: preallocation^Matrix 896bd481603SBarry Smith 897d6e23781SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(), 898d6e23781SBarry Smith MatPreallocateSymmetricSetLocalBlock() 899bd481603SBarry Smith M*/ 900a89bc211SBarry Smith #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree2(dnz,onz);CHKERRQ(_4_ierr);} 9017c922b88SBarry Smith 9027b80b807SBarry Smith /* Routines unique to particular data structures */ 903014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellGetContext(Mat,void *); 904698d4c6aSKris Buschelman 905014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInodeAdjustForInodes(Mat,IS*,IS*); 906014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *); 907698d4c6aSKris Buschelman 908014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetColumnIndices(Mat,PetscInt[]); 909014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetColumnIndices(Mat,PetscInt[]); 910014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 911014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 912014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 913014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJFromTriple(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*,PetscInt,PetscBool); 9147b80b807SBarry Smith 915a96a251dSBarry Smith #define MAT_SKIP_ALLOCATION -4 916a96a251dSBarry Smith 917014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]); 918014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]); 919014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]); 920273d9f13SBarry Smith 921014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 922014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 923014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 924014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []); 925014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 926014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]); 927014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 928014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]); 929014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIDenseSetPreallocation(Mat,PetscScalar[]); 930014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqDenseSetPreallocation(Mat,PetscScalar[]); 9319230625dSJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,const PetscInt*[]); 9329230625dSJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,const PetscInt*[]); 933014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAdjCreateNonemptySubcommMat(Mat,Mat*); 934273d9f13SBarry Smith 9352e1947a5SStefano Zampini PETSC_EXTERN PetscErrorCode MatISSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 936014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqDenseSetLDA(Mat,PetscInt); 937014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDenseGetLocalMatrix(Mat,Mat*); 9381b807ce4Svictorle 939014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStoreValues(Mat); 940014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRetrieveValues(Mat); 9412e8a6d31SBarry Smith 942014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDAADSetCtx(Mat,void*); 9433a7fca6bSBarry Smith 944014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFindNonzeroRows(Mat,IS*); 9457b80b807SBarry Smith /* 9467b80b807SBarry Smith These routines are not usually accessed directly, rather solving is 94794b7f48cSBarry Smith done through the KSP and PC interfaces. 9487b80b807SBarry Smith */ 9497b80b807SBarry Smith 95076bdecfbSBarry Smith /*J 9518f6c3df8SBarry Smith MatOrderingType - String with the name of a PETSc matrix ordering 952d9274352SBarry Smith 953d9274352SBarry Smith Level: beginner 954d9274352SBarry Smith 955d9274352SBarry Smith .seealso: MatGetOrdering() 95676bdecfbSBarry Smith J*/ 95719fd82e9SBarry Smith typedef const char* MatOrderingType; 9582692d6eeSBarry Smith #define MATORDERINGNATURAL "natural" 9592692d6eeSBarry Smith #define MATORDERINGND "nd" 9602692d6eeSBarry Smith #define MATORDERING1WD "1wd" 9612692d6eeSBarry Smith #define MATORDERINGRCM "rcm" 9622692d6eeSBarry Smith #define MATORDERINGQMD "qmd" 9632692d6eeSBarry Smith #define MATORDERINGROWLENGTH "rowlength" 964510184a7SJed Brown #define MATORDERINGWBM "wbm" 965fc1bef75SJed Brown #define MATORDERINGSPECTRAL "spectral" 966312e372aSBarry Smith #define MATORDERINGAMD "amd" /* only works if UMFPACK is installed with PETSc */ 967b12f92e5SBarry Smith 96819fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatGetOrdering(Mat,MatOrderingType,IS*,IS*); 969140e18c1SBarry Smith PETSC_EXTERN PetscErrorCode MatGetOrderingList(PetscFunctionList*); 970bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatOrderingRegister(const char[],PetscErrorCode(*)(Mat,MatOrderingType,IS*,IS*)); 971140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatOrderingList; 972d4fbbf0eSBarry Smith 973014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS); 974fc1bef75SJed Brown PETSC_EXTERN PetscErrorCode MatCreateLaplacian(Mat,PetscReal,PetscBool,Mat*); 975a2ce50c7SBarry Smith 976d91e6319SBarry Smith /*S 977d90ac83dSHong Zhang MatFactorShiftType - Numeric Shift. 978d90ac83dSHong Zhang 979d90ac83dSHong Zhang Level: beginner 980d90ac83dSHong Zhang 981d90ac83dSHong Zhang S*/ 982d90ac83dSHong Zhang typedef enum {MAT_SHIFT_NONE,MAT_SHIFT_NONZERO,MAT_SHIFT_POSITIVE_DEFINITE,MAT_SHIFT_INBLOCKS} MatFactorShiftType; 9836a6fc655SJed Brown PETSC_EXTERN const char *const MatFactorShiftTypes[]; 9845e9742b9SJed Brown PETSC_EXTERN const char *const MatFactorShiftTypesDetail[]; 985d90ac83dSHong Zhang 986d90ac83dSHong Zhang /*S 9879aa7eafdSHong Zhang MatFactorError - indicates what type of error in matrix factor 9889aa7eafdSHong Zhang 9899aa7eafdSHong Zhang Level: beginner 9909aa7eafdSHong Zhang 99100e125f8SBarry Smith Developer Notes: Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 99200e125f8SBarry Smith 99300e125f8SBarry Smith .seealso: MatGetFactor() 9949aa7eafdSHong Zhang S*/ 9955cd7cf9dSHong Zhang typedef enum {MAT_FACTOR_NOERROR,MAT_FACTOR_STRUCT_ZEROPIVOT,MAT_FACTOR_NUMERIC_ZEROPIVOT,MAT_FACTOR_OUTMEMORY,MAT_FACTOR_OTHER} MatFactorError; 9969aa7eafdSHong Zhang 99700e125f8SBarry Smith PETSC_EXTERN PetscErrorCode MatFactorGetError(Mat,MatFactorError*); 998b8b68cfdSBarry Smith PETSC_EXTERN PetscErrorCode MatFactorClearError(Mat); 9997b6c816cSBarry Smith PETSC_EXTERN PetscErrorCode MatFactorGetErrorZeroPivot(Mat,PetscReal*,PetscInt*); 100000e125f8SBarry Smith 10019aa7eafdSHong Zhang /*S 1002422a814eSBarry Smith MatFactorInfo - Data passed into the matrix factorization routines, and information about the resulting factorization 1003b00f7748SHong Zhang 100461cad860SBarry Smith In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE, that is use 100561cad860SBarry Smith $ MatFactorInfo info(MAT_FACTORINFO_SIZE) 1006b00f7748SHong Zhang 100715e8a5b3SHong Zhang Notes: These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC. 1008b00f7748SHong Zhang 100961cad860SBarry Smith You can use MatFactorInfoInitialize() to set default values. 101061cad860SBarry Smith 1011b00f7748SHong Zhang Level: developer 1012b00f7748SHong Zhang 1013d7d82daaSBarry Smith .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(), 1014d7d82daaSBarry Smith MatFactorInfoInitialize() 1015b00f7748SHong Zhang 1016b00f7748SHong Zhang S*/ 1017b00f7748SHong Zhang typedef struct { 101815e8a5b3SHong Zhang PetscReal diagonal_fill; /* force diagonal to fill in if initially not filled */ 101985317021SBarry Smith PetscReal usedt; 102015e8a5b3SHong Zhang PetscReal dt; /* drop tolerance */ 1021b00f7748SHong Zhang PetscReal dtcol; /* tolerance for pivoting */ 102215e8a5b3SHong Zhang PetscReal dtcount; /* maximum nonzeros to be allowed per row */ 102367e28bfeSBarry Smith PetscReal fill; /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */ 1024348344bbSBarry Smith PetscReal levels; /* ICC/ILU(levels) */ 1025bcd9e38bSBarry Smith PetscReal pivotinblocks; /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0 1026bcd9e38bSBarry Smith factorization may be faster if do not pivot */ 102715e8a5b3SHong Zhang PetscReal zeropivot; /* pivot is called zero if less than this */ 1028f4db908eSBarry Smith PetscReal shifttype; /* type of shift added to matrix factor to prevent zero pivots */ 1029f4db908eSBarry Smith PetscReal shiftamount; /* how large the shift is */ 103015e8a5b3SHong Zhang } MatFactorInfo; 1031ffa6d0a5SLois Curfman McInnes 1032014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFactorInfoInitialize(MatFactorInfo*); 10339a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatCholeskyFactor(Mat,IS,const MatFactorInfo*); 10349a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatCholeskyFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*); 10359a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatCholeskyFactorNumeric(Mat,Mat,const MatFactorInfo*); 10369a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatLUFactor(Mat,IS,IS,const MatFactorInfo*); 10379a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatILUFactor(Mat,IS,IS,const MatFactorInfo*); 10389a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatLUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*); 10399a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatILUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*); 10409a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatICCFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*); 10419a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatICCFactor(Mat,IS,const MatFactorInfo*); 10429a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatLUFactorNumeric(Mat,Mat,const MatFactorInfo*); 1043014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*); 1044014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolve(Mat,Vec,Vec); 1045014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatForwardSolve(Mat,Vec,Vec); 1046014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatBackwardSolve(Mat,Vec,Vec); 1047014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveAdd(Mat,Vec,Vec,Vec); 1048014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveTranspose(Mat,Vec,Vec); 1049014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveTransposeAdd(Mat,Vec,Vec,Vec); 1050014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolves(Mat,Vecs,Vecs); 10518e7ba810SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorSetSchurIS(Mat,IS); 10525a05ddb0SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorGetSchurComplement(Mat,Mat*); 10535a05ddb0SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorRestoreSchurComplement(Mat,Mat*); 10545a05ddb0SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorInvertSchurComplement(Mat); 10555a05ddb0SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorCreateSchurComplement(Mat,Mat*); 10565a05ddb0SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorSolveSchurComplement(Mat,Vec,Vec); 10575a05ddb0SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorSolveSchurComplementTranspose(Mat,Vec,Vec); 1058014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetUnfactored(Mat); 1059bb5a7306SBarry Smith 1060d91e6319SBarry Smith /*E 1061d91e6319SBarry Smith MatSORType - What type of (S)SOR to perform 1062bb1eb677SSatish Balay 1063d91e6319SBarry Smith Level: beginner 1064d91e6319SBarry Smith 1065d9274352SBarry Smith May be bitwise ORd together 1066d9274352SBarry Smith 1067af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 1068d91e6319SBarry Smith 10694e7234bfSBarry Smith MatSORType may be bitwise ORd together, so do not change the numbers 10704e7234bfSBarry Smith 107141f059aeSBarry Smith .seealso: MatSOR() 1072d91e6319SBarry Smith E*/ 1073ee50ffe9SBarry Smith typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3, 1074ee50ffe9SBarry Smith SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8, 1075ee50ffe9SBarry Smith SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16, 107684cb2905SBarry Smith SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType; 1077014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSOR(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec); 10788ed539a5SBarry Smith 1079d4fbbf0eSBarry Smith /* 1080639f9d9dSBarry Smith These routines are for efficiently computing Jacobians via finite differences. 1081639f9d9dSBarry Smith */ 1082b12f92e5SBarry Smith 10837086a01eSPeter Brune /*S 10847086a01eSPeter Brune MatColoring - Object for managing the coloring of matrices. 10857086a01eSPeter Brune 10867086a01eSPeter Brune Level: beginner 10877086a01eSPeter Brune 10887086a01eSPeter Brune Concepts: matrix, coloring 10897086a01eSPeter Brune 10907086a01eSPeter Brune .seealso: MatFDColoringCreate() ISColoring MatFDColoring 10917086a01eSPeter Brune S*/ 10927086a01eSPeter Brune typedef struct _p_MatColoring* MatColoring; 109376bdecfbSBarry Smith /*J 10948f6c3df8SBarry Smith MatColoringType - String with the name of a PETSc matrix coloring 1095d9274352SBarry Smith 1096d9274352SBarry Smith Level: beginner 1097d9274352SBarry Smith 10987086a01eSPeter Brune .seealso: MatColoringSetType(), MatColoring 109976bdecfbSBarry Smith J*/ 11007086a01eSPeter Brune 110119fd82e9SBarry Smith typedef const char* MatColoringType; 11025567d71aSPeter Brune #define MATCOLORINGJP "jp" 1103b9acb617SPeter Brune #define MATCOLORINGPOWER "power" 11042692d6eeSBarry Smith #define MATCOLORINGNATURAL "natural" 11052692d6eeSBarry Smith #define MATCOLORINGSL "sl" 11062692d6eeSBarry Smith #define MATCOLORINGLF "lf" 11072692d6eeSBarry Smith #define MATCOLORINGID "id" 11088121bdceSPeter Brune #define MATCOLORINGGREEDY "greedy" 1109b12f92e5SBarry Smith 11108ac6da0aSPeter Brune /*E 11118ac6da0aSPeter Brune MatColoringWeightType - Type of weight scheme 11128ac6da0aSPeter Brune 11138ac6da0aSPeter Brune Not Collective 11148ac6da0aSPeter Brune 11158ac6da0aSPeter Brune + MAT_COLORING_RANDOM - Random weights 11168ac6da0aSPeter Brune . MAT_COLORING_LEXICAL - Lexical weighting based upon global numbering. 11178ac6da0aSPeter Brune - MAT_COLORING_LF - Last-first weighting. 11188ac6da0aSPeter Brune 11198ac6da0aSPeter Brune Level: intermediate 11208ac6da0aSPeter Brune 1121af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 11228ac6da0aSPeter Brune 11238ac6da0aSPeter Brune .seealso: MatCUSPSetFormat(), MatCUSPFormatOperation 11248ac6da0aSPeter Brune E*/ 1125301c140bSPeter Brune typedef enum {MAT_COLORING_WEIGHT_RANDOM,MAT_COLORING_WEIGHT_LEXICAL,MAT_COLORING_WEIGHT_LF,MAT_COLORING_WEIGHT_SL} MatColoringWeightType; 11268ac6da0aSPeter Brune 1127335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringCreate(Mat,MatColoring*); 1128d7f2a307SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringGetDegrees(Mat,PetscInt,PetscInt*); 1129335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringDestroy(MatColoring*); 1130335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringView(MatColoring,PetscViewer); 1131335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetType(MatColoring,MatColoringType); 1132335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetFromOptions(MatColoring); 1133335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetDistance(MatColoring,PetscInt); 1134335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringGetDistance(MatColoring,PetscInt*); 1135335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetMaxColors(MatColoring,PetscInt); 1136335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringGetMaxColors(MatColoring,PetscInt*); 1137335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringApply(MatColoring,ISColoring*); 1138335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringRegister(const char[],PetscErrorCode(*)(MatColoring)); 1139014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*); 11408ac6da0aSPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetWeightType(MatColoring,MatColoringWeightType); 1141c29971fcSPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetWeights(MatColoring,PetscReal*,PetscInt*); 11428ac6da0aSPeter Brune PETSC_EXTERN PetscErrorCode MatColoringCreateWeights(MatColoring,PetscReal **,PetscInt **lperm); 1143cb394dc5SBarry Smith PETSC_EXTERN PetscErrorCode MatColoringTestValid(MatColoring,ISColoring); 1144639f9d9dSBarry Smith 1145d9274352SBarry Smith /*S 1146d9274352SBarry Smith MatFDColoring - Object for computing a sparse Jacobian via finite differences 1147d9274352SBarry Smith and coloring 1148639f9d9dSBarry Smith 1149d9274352SBarry Smith Level: beginner 1150d9274352SBarry Smith 1151d9274352SBarry Smith Concepts: coloring, sparse Jacobian, finite differences 1152d9274352SBarry Smith 1153d9274352SBarry Smith .seealso: MatFDColoringCreate() 1154d9274352SBarry Smith S*/ 1155e2a1c21fSSatish Balay typedef struct _p_MatFDColoring* MatFDColoring; 1156639f9d9dSBarry Smith 1157014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringCreate(Mat,ISColoring,MatFDColoring *); 1158014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringDestroy(MatFDColoring*); 1159014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringView(MatFDColoring,PetscViewer); 1160014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*); 1161014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**); 1162014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal); 1163014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetFromOptions(MatFDColoring); 1164d1e9a80fSBarry Smith PETSC_EXTERN PetscErrorCode MatFDColoringApply(Mat,MatFDColoring,Vec,void *); 1165014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetF(MatFDColoring,Vec); 1166014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,PetscInt*[]); 1167f86b9fbaSHong Zhang PETSC_EXTERN PetscErrorCode MatFDColoringSetUp(Mat,ISColoring,MatFDColoring); 1168f86b9fbaSHong Zhang PETSC_EXTERN PetscErrorCode MatFDColoringSetBlockSize(MatFDColoring,PetscInt,PetscInt); 1169f86b9fbaSHong Zhang 1170b1683b59SHong Zhang 1171b1683b59SHong Zhang /*S 1172b9af6bddSHong Zhang MatTransposeColoring - Object for computing a sparse matrix product C=A*B^T via coloring 1173b1683b59SHong Zhang 1174b1683b59SHong Zhang Level: beginner 1175b1683b59SHong Zhang 1176b1683b59SHong Zhang Concepts: coloring, sparse matrix product 1177b1683b59SHong Zhang 1178b9af6bddSHong Zhang .seealso: MatTransposeColoringCreate() 1179b1683b59SHong Zhang S*/ 1180b9af6bddSHong Zhang typedef struct _p_MatTransposeColoring* MatTransposeColoring; 1181b1683b59SHong Zhang 1182014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeColoringCreate(Mat,ISColoring,MatTransposeColoring *); 1183014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransColoringApplySpToDen(MatTransposeColoring,Mat,Mat); 1184014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransColoringApplyDenToSp(MatTransposeColoring,Mat,Mat); 1185014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeColoringDestroy(MatTransposeColoring*); 1186b1683b59SHong Zhang 1187639f9d9dSBarry Smith /* 11880752156aSBarry Smith These routines are for partitioning matrices: currently used only 11893eda8832SBarry Smith for adjacency matrix, MatCreateMPIAdj(). 11900752156aSBarry Smith */ 1191ca161407SBarry Smith 1192d9274352SBarry Smith /*S 1193d9274352SBarry Smith MatPartitioning - Object for managing the partitioning of a matrix or graph 1194d9274352SBarry Smith 1195d9274352SBarry Smith Level: beginner 1196d9274352SBarry Smith 1197d9274352SBarry Smith Concepts: partitioning 1198d9274352SBarry Smith 1199743a1026Svictorle .seealso: MatPartitioningCreate(), MatPartitioningType 1200d9274352SBarry Smith S*/ 120191e9ee9fSBarry Smith typedef struct _p_MatPartitioning* MatPartitioning; 1202d9274352SBarry Smith 120376bdecfbSBarry Smith /*J 12048f6c3df8SBarry Smith MatPartitioningType - String with the name of a PETSc matrix partitioning 1205d9274352SBarry Smith 1206d9274352SBarry Smith Level: beginner 12070d04baf8SBarry Smith dm 1208b547a13bSvictorle .seealso: MatPartitioningCreate(), MatPartitioning 120976bdecfbSBarry Smith J*/ 121019fd82e9SBarry Smith typedef const char* MatPartitioningType; 12112692d6eeSBarry Smith #define MATPARTITIONINGCURRENT "current" 1212aa1e27eaSFande Kong #define MATPARTITIONINGAVERAGE "average" 12132692d6eeSBarry Smith #define MATPARTITIONINGSQUARE "square" 12142692d6eeSBarry Smith #define MATPARTITIONINGPARMETIS "parmetis" 12152692d6eeSBarry Smith #define MATPARTITIONINGCHACO "chaco" 12162692d6eeSBarry Smith #define MATPARTITIONINGPARTY "party" 121750d91057SBarry Smith #define MATPARTITIONINGPTSCOTCH "ptscotch" 121888d2ac2bSFande Kong #define MATPARTITIONINGHIERARCH "hierarch" 121971306c60Sjroman 1220ca161407SBarry Smith 1221014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningCreate(MPI_Comm,MatPartitioning*); 122219fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningSetType(MatPartitioning,MatPartitioningType); 1223014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetNParts(MatPartitioning,PetscInt); 1224014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetAdjacency(MatPartitioning,Mat); 1225014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]); 1226014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []); 1227014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningApply(MatPartitioning,IS*); 1228014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningDestroy(MatPartitioning*); 12292aabb6bbSBarry Smith 1230bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningRegister(const char[],PetscErrorCode (*)(MatPartitioning)); 123130de9b25SBarry Smith 1232f1144a30SSatish Balay 12332bad1931SBarry Smith 1234014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningView(MatPartitioning,PetscViewer); 1235014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetFromOptions(MatPartitioning); 123619fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningGetType(MatPartitioning,MatPartitioningType*); 1237ca161407SBarry Smith 123827538973SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetRepartition(MatPartitioning part); 1239014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetCoarseSequential(MatPartitioning); 1240014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningParmetisGetEdgeCut(MatPartitioning, PetscInt *); 12410752156aSBarry Smith 1242b6956c12SJose Roman typedef enum { MP_CHACO_MULTILEVEL=1,MP_CHACO_SPECTRAL=2,MP_CHACO_LINEAR=4,MP_CHACO_RANDOM=5,MP_CHACO_SCATTERED=6 } MPChacoGlobalType; 12436a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoGlobalTypes[]; 1244b6956c12SJose Roman typedef enum { MP_CHACO_KERNIGHAN=1,MP_CHACO_NONE=2 } MPChacoLocalType; 12456a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoLocalTypes[]; 1246b6956c12SJose Roman typedef enum { MP_CHACO_LANCZOS=0,MP_CHACO_RQI=1 } MPChacoEigenType; 12476a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoEigenTypes[]; 1248b6956c12SJose Roman 1249014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetGlobal(MatPartitioning,MPChacoGlobalType); 1250014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetGlobal(MatPartitioning,MPChacoGlobalType*); 1251014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetLocal(MatPartitioning,MPChacoLocalType); 1252014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetLocal(MatPartitioning,MPChacoLocalType*); 1253014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal); 1254014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType); 1255014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenSolver(MatPartitioning,MPChacoEigenType*); 1256014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenTol(MatPartitioning,PetscReal); 1257014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenTol(MatPartitioning,PetscReal*); 1258014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenNumber(MatPartitioning,PetscInt); 1259014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenNumber(MatPartitioning,PetscInt*); 126071306c60Sjroman 126171306c60Sjroman #define MP_PARTY_OPT "opt" 126271306c60Sjroman #define MP_PARTY_LIN "lin" 126371306c60Sjroman #define MP_PARTY_SCA "sca" 126471306c60Sjroman #define MP_PARTY_RAN "ran" 126571306c60Sjroman #define MP_PARTY_GBF "gbf" 126671306c60Sjroman #define MP_PARTY_GCF "gcf" 126771306c60Sjroman #define MP_PARTY_BUB "bub" 126871306c60Sjroman #define MP_PARTY_DEF "def" 1269014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetGlobal(MatPartitioning,const char*); 127071306c60Sjroman #define MP_PARTY_HELPFUL_SETS "hs" 127171306c60Sjroman #define MP_PARTY_KERNIGHAN_LIN "kl" 127271306c60Sjroman #define MP_PARTY_NONE "no" 1273014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetLocal(MatPartitioning,const char*); 1274014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal); 1275014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetBipart(MatPartitioning,PetscBool); 1276014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscBool); 127771306c60Sjroman 127850d91057SBarry Smith typedef enum { MP_PTSCOTCH_QUALITY,MP_PTSCOTCH_SPEED,MP_PTSCOTCH_BALANCE,MP_PTSCOTCH_SAFETY,MP_PTSCOTCH_SCALABILITY } MPPTScotchStrategyType; 12796a6fc655SJed Brown PETSC_EXTERN const char *const MPPTScotchStrategyTypes[]; 1280e0f1cffaSJose Roman 1281014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetImbalance(MatPartitioning,PetscReal); 1282014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetImbalance(MatPartitioning,PetscReal*); 1283014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetStrategy(MatPartitioning,MPPTScotchStrategyType); 1284014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetStrategy(MatPartitioning,MPPTScotchStrategyType*); 128571306c60Sjroman 1286b43b03e9SMark F. Adams /* 12876294fa2bSFande Kong * hierarchical partitioning 12886294fa2bSFande Kong */ 12894e713f55SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalGetFineparts(MatPartitioning,IS*); 12904e713f55SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalGetCoarseparts(MatPartitioning,IS*); 12914e713f55SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalSetNcoarseparts(MatPartitioning,PetscInt); 12924e713f55SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalSetNfineparts(MatPartitioning, PetscInt); 12936294fa2bSFande Kong 1294014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMeshToVertexGraph(Mat,PetscInt,Mat*); 1295014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMeshToCellGraph(Mat,PetscInt,Mat*); 1296591294e4SBarry Smith 12970752156aSBarry Smith /* 1298af0996ceSBarry Smith If you add entries here you must also add them to petsc/finclude/petscmat.h 1299d4fbbf0eSBarry Smith */ 13001c1c02c0SLois Curfman McInnes typedef enum { MATOP_SET_VALUES=0, 13011c1c02c0SLois Curfman McInnes MATOP_GET_ROW=1, 13021c1c02c0SLois Curfman McInnes MATOP_RESTORE_ROW=2, 13031c1c02c0SLois Curfman McInnes MATOP_MULT=3, 13041c1c02c0SLois Curfman McInnes MATOP_MULT_ADD=4, 13057c922b88SBarry Smith MATOP_MULT_TRANSPOSE=5, 13067c922b88SBarry Smith MATOP_MULT_TRANSPOSE_ADD=6, 13071c1c02c0SLois Curfman McInnes MATOP_SOLVE=7, 13081c1c02c0SLois Curfman McInnes MATOP_SOLVE_ADD=8, 13097c922b88SBarry Smith MATOP_SOLVE_TRANSPOSE=9, 13107c922b88SBarry Smith MATOP_SOLVE_TRANSPOSE_ADD=10, 13111c1c02c0SLois Curfman McInnes MATOP_LUFACTOR=11, 13121c1c02c0SLois Curfman McInnes MATOP_CHOLESKYFACTOR=12, 1313a714c06dSBarry Smith MATOP_SOR=13, 13141c1c02c0SLois Curfman McInnes MATOP_TRANSPOSE=14, 13151c1c02c0SLois Curfman McInnes MATOP_GETINFO=15, 13161c1c02c0SLois Curfman McInnes MATOP_EQUAL=16, 13171c1c02c0SLois Curfman McInnes MATOP_GET_DIAGONAL=17, 13181c1c02c0SLois Curfman McInnes MATOP_DIAGONAL_SCALE=18, 13191c1c02c0SLois Curfman McInnes MATOP_NORM=19, 13201c1c02c0SLois Curfman McInnes MATOP_ASSEMBLY_BEGIN=20, 13211c1c02c0SLois Curfman McInnes MATOP_ASSEMBLY_END=21, 1322d519adbfSMatthew Knepley MATOP_SET_OPTION=22, 1323d519adbfSMatthew Knepley MATOP_ZERO_ENTRIES=23, 1324d519adbfSMatthew Knepley MATOP_ZERO_ROWS=24, 1325d519adbfSMatthew Knepley MATOP_LUFACTOR_SYMBOLIC=25, 1326d519adbfSMatthew Knepley MATOP_LUFACTOR_NUMERIC=26, 1327d519adbfSMatthew Knepley MATOP_CHOLESKY_FACTOR_SYMBOLIC=27, 1328d519adbfSMatthew Knepley MATOP_CHOLESKY_FACTOR_NUMERIC=28, 1329d519adbfSMatthew Knepley MATOP_SETUP_PREALLOCATION=29, 1330d519adbfSMatthew Knepley MATOP_ILUFACTOR_SYMBOLIC=30, 1331d519adbfSMatthew Knepley MATOP_ICCFACTOR_SYMBOLIC=31, 1332a5b7ff6bSBarry Smith MATOP_GET_DIAGONAL_BLOCK=32, 13339d39f69dSJed Brown /* MATOP_PLACEHOLDER_33=33, */ 1334d519adbfSMatthew Knepley MATOP_DUPLICATE=34, 1335d519adbfSMatthew Knepley MATOP_FORWARD_SOLVE=35, 1336d519adbfSMatthew Knepley MATOP_BACKWARD_SOLVE=36, 1337d519adbfSMatthew Knepley MATOP_ILUFACTOR=37, 1338d519adbfSMatthew Knepley MATOP_ICCFACTOR=38, 1339d519adbfSMatthew Knepley MATOP_AXPY=39, 1340d519adbfSMatthew Knepley MATOP_GET_SUBMATRICES=40, 1341d519adbfSMatthew Knepley MATOP_INCREASE_OVERLAP=41, 1342d519adbfSMatthew Knepley MATOP_GET_VALUES=42, 1343d519adbfSMatthew Knepley MATOP_COPY=43, 1344d519adbfSMatthew Knepley MATOP_GET_ROW_MAX=44, 1345d519adbfSMatthew Knepley MATOP_SCALE=45, 1346d519adbfSMatthew Knepley MATOP_SHIFT=46, 134735153367SBarry Smith MATOP_DIAGONAL_SET=47, 13489d39f69dSJed Brown MATOP_ZERO_ROWS_COLUMNS=48, 13499d39f69dSJed Brown MATOP_SET_RANDOM=49, 1350d519adbfSMatthew Knepley MATOP_GET_ROW_IJ=50, 1351d519adbfSMatthew Knepley MATOP_RESTORE_ROW_IJ=51, 1352d519adbfSMatthew Knepley MATOP_GET_COLUMN_IJ=52, 1353d519adbfSMatthew Knepley MATOP_RESTORE_COLUMN_IJ=53, 1354d519adbfSMatthew Knepley MATOP_FDCOLORING_CREATE=54, 1355d519adbfSMatthew Knepley MATOP_COLORING_PATCH=55, 1356d519adbfSMatthew Knepley MATOP_SET_UNFACTORED=56, 1357d519adbfSMatthew Knepley MATOP_PERMUTE=57, 1358d519adbfSMatthew Knepley MATOP_SET_VALUES_BLOCKED=58, 1359d519adbfSMatthew Knepley MATOP_GET_SUBMATRIX=59, 1360d519adbfSMatthew Knepley MATOP_DESTROY=60, 1361d519adbfSMatthew Knepley MATOP_VIEW=61, 1362d519adbfSMatthew Knepley MATOP_CONVERT_FROM=62, 13637bab7c10SHong Zhang MATOP_MATMAT_MULT=63, 13647bab7c10SHong Zhang MATOP_MATMAT_MULT_SYMBOLIC=64, 13657bab7c10SHong Zhang MATOP_MATMAT_MULT_NUMERIC=65, 1366d519adbfSMatthew Knepley MATOP_SET_LOCAL_TO_GLOBAL_MAP=66, 1367d519adbfSMatthew Knepley MATOP_SET_VALUES_LOCAL=67, 1368d519adbfSMatthew Knepley MATOP_ZERO_ROWS_LOCAL=68, 1369d519adbfSMatthew Knepley MATOP_GET_ROW_MAX_ABS=69, 1370d519adbfSMatthew Knepley MATOP_GET_ROW_MIN_ABS=70, 1371d519adbfSMatthew Knepley MATOP_CONVERT=71, 1372d519adbfSMatthew Knepley MATOP_SET_COLORING=72, 13739d39f69dSJed Brown /* MATOP_PLACEHOLDER_73=73, */ 1374d519adbfSMatthew Knepley MATOP_SET_VALUES_ADIFOR=74, 1375d519adbfSMatthew Knepley MATOP_FD_COLORING_APPLY=75, 1376d519adbfSMatthew Knepley MATOP_SET_FROM_OPTIONS=76, 1377bda6c4cbSBarry Smith MATOP_MULT_CONSTRAINED=77, 1378bda6c4cbSBarry Smith MATOP_MULT_TRANSPOSE_CONSTRAIN=78, 13799d39f69dSJed Brown MATOP_FIND_ZERO_DIAGONALS=79, 1380d519adbfSMatthew Knepley MATOP_MULT_MULTIPLE=80, 1381d519adbfSMatthew Knepley MATOP_SOLVE_MULTIPLE=81, 1382d519adbfSMatthew Knepley MATOP_GET_INERTIA=82, 1383d519adbfSMatthew Knepley MATOP_LOAD=83, 1384d519adbfSMatthew Knepley MATOP_IS_SYMMETRIC=84, 1385d519adbfSMatthew Knepley MATOP_IS_HERMITIAN=85, 1386d519adbfSMatthew Knepley MATOP_IS_STRUCTURALLY_SYMMETRIC=86, 1387820469bcSHong Zhang MATOP_SET_VALUES_BLOCKEDLOCAL=87, 1388d519adbfSMatthew Knepley MATOP_GET_VECS=88, 1389d519adbfSMatthew Knepley MATOP_MAT_MULT=89, 1390d519adbfSMatthew Knepley MATOP_MAT_MULT_SYMBOLIC=90, 1391d519adbfSMatthew Knepley MATOP_MAT_MULT_NUMERIC=91, 1392d519adbfSMatthew Knepley MATOP_PTAP=92, 1393d519adbfSMatthew Knepley MATOP_PTAP_SYMBOLIC=93, 1394d519adbfSMatthew Knepley MATOP_PTAP_NUMERIC=94, 1395bda6c4cbSBarry Smith MATOP_MAT_TRANSPOSE_MULT=95, 1396bda6c4cbSBarry Smith MATOP_MAT_TRANSPOSE_MULT_SYMBO=96, 1397bda6c4cbSBarry Smith MATOP_MAT_TRANSPOSE_MULT_NUMER=97, 13989d39f69dSJed Brown /* MATOP_PLACEHOLDER_98=98, */ 13999d39f69dSJed Brown /* MATOP_PLACEHOLDER_99=99, */ 14009d39f69dSJed Brown /* MATOP_PLACEHOLDER_100=100, */ 14019d39f69dSJed Brown /* MATOP_PLACEHOLDER_101=101, */ 1402d519adbfSMatthew Knepley MATOP_CONJUGATE=102, 14039d39f69dSJed Brown /* MATOP_PLACEHOLDER_103=103, */ 1404d519adbfSMatthew Knepley MATOP_SET_VALUES_ROW=104, 1405d519adbfSMatthew Knepley MATOP_REAL_PART=105, 1406bda6c4cbSBarry Smith MATOP_IMAGINARY_PART=106, 1407bda6c4cbSBarry Smith MATOP_GET_ROW_UPPER_TRIANGULAR=107, 1408bda6c4cbSBarry Smith MATOP_RESTORE_ROW_UPPER_TRIANG=108, 1409bda6c4cbSBarry Smith MATOP_MAT_SOLVE=109, 14100e8d04f7SBarry Smith MATOP_GET_REDUNDANT_MATRIX=110, 1411d519adbfSMatthew Knepley MATOP_GET_ROW_MIN=111, 14120e8d04f7SBarry Smith MATOP_GET_COLUMN_VECTOR=112, 1413d519adbfSMatthew Knepley MATOP_MISSING_DIAGONAL=113, 14140e8d04f7SBarry Smith MATOP_GET_SEQ_NONZERO_STRUCTUR=114, 141589c6957cSBarry Smith MATOP_CREATE=115, 141689c6957cSBarry Smith MATOP_GET_GHOSTS=116, 14170e8d04f7SBarry Smith MATOP_GET_LOCAL_SUB_MATRIX=117, 14180e8d04f7SBarry Smith MATOP_RESTORE_LOCALSUB_MATRIX=118, 1419eeffb40dSHong Zhang MATOP_MULT_DIAGONAL_BLOCK=119, 14200e8d04f7SBarry Smith MATOP_HERMITIAN_TRANSPOSE=120, 14210e8d04f7SBarry Smith MATOP_MULT_HERMITIAN_TRANSPOSE=121, 14220e8d04f7SBarry Smith MATOP_MULT_HERMITIAN_TRANS_ADD=122, 14230e8d04f7SBarry Smith MATOP_GET_MULTI_PROC_BLOCK=123, 14249d39f69dSJed Brown MATOP_FIND_NONZERO_ROWS=124, 14250e8d04f7SBarry Smith MATOP_GET_COLUMN_NORMS=125, 14269d39f69dSJed Brown MATOP_INVERT_BLOCK_DIAGONAL=126, 14279d39f69dSJed Brown /* MATOP_PLACEHOLDER_127=127, */ 14280e8d04f7SBarry Smith MATOP_GET_SUB_MATRICES_PARALLE=128, 14292b8ad9a3SHong Zhang MATOP_SET_VALUES_BATCH=129, 14300e8d04f7SBarry Smith MATOP_TRANSPOSE_MAT_MULT=130, 14310e8d04f7SBarry Smith MATOP_TRANSPOSE_MAT_MULT_SYMBO=131, 1432e9b602ebSSatish Balay MATOP_TRANSPOSE_MAT_MULT_NUMER=132, 14330e8d04f7SBarry Smith MATOP_TRANSPOSE_COLORING_CREAT=133, 14340e8d04f7SBarry Smith MATOP_TRANS_COLORING_APPLY_SPT=134, 14350e8d04f7SBarry Smith MATOP_TRANS_COLORING_APPLY_DEN=135, 14360e8d04f7SBarry Smith MATOP_RART=136, 14370e8d04f7SBarry Smith MATOP_RART_SYMBOLIC=137, 14380e8d04f7SBarry Smith MATOP_RART_NUMERIC=138, 1439e09a3074SHong Zhang MATOP_SET_BLOCK_SIZES=139, 1440f9426fe0SMark Adams MATOP_AYPX=140, 14411919a2e2SJed Brown MATOP_RESIDUAL=141, 14429c8f2541SHong Zhang MATOP_FDCOLORING_SETUP=142, 14439c8f2541SHong Zhang MATOP_MPICONCATENATESEQ=144 1444fae171e0SBarry Smith } MatOperation; 1445014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatHasOperation(Mat,MatOperation,PetscBool *); 1446014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellSetOperation(Mat,MatOperation,void(*)(void)); 1447014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellGetOperation(Mat,MatOperation,void(**)(void)); 1448014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellSetContext(Mat,void*); 1449112a2221SBarry Smith 145090ace30eSBarry Smith /* 145190ace30eSBarry Smith Codes for matrices stored on disk. By default they are 145290ace30eSBarry Smith stored in a universal format. By changing the format with 14536a9046bcSBarry Smith PetscViewerPushFormat(viewer,PETSC_VIEWER_NATIVE); the matrices will 145490ace30eSBarry Smith be stored in a way natural for the matrix, for example dense matrices 145590ace30eSBarry Smith would be stored as dense. Matrices stored this way may only be 1456f4403165SShri Abhyankar read into matrices of the same type. 145790ace30eSBarry Smith */ 145890ace30eSBarry Smith #define MATRIX_BINARY_FORMAT_DENSE -1 145990ace30eSBarry Smith 1460014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetHashTableFactor(Mat,PetscReal); 1461014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatISGetLocalMat(Mat,Mat*); 1462014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatISSetLocalMat(Mat,Mat); 1463b7ce53b6SStefano Zampini PETSC_EXTERN PetscErrorCode MatISGetMPIXAIJ(Mat,MatReuse,Mat*); 14641f4e1ec7SBarry Smith 1465d9274352SBarry Smith /*S 1466d9274352SBarry Smith MatNullSpace - Object that removes a null space from a vector, i.e. 1467d9274352SBarry Smith orthogonalizes the vector to a subsapce 1468d9274352SBarry Smith 1469f7a9e4ceSBarry Smith Level: advanced 1470d9274352SBarry Smith 1471d9274352SBarry Smith Concepts: matrix; linear operator, null space 1472d9274352SBarry Smith 14736e1639daSBarry Smith Users manual sections: 14746e1639daSBarry Smith . sec_singular 14756e1639daSBarry Smith 1476d9274352SBarry Smith .seealso: MatNullSpaceCreate() 1477d9274352SBarry Smith S*/ 147874637425SBarry Smith typedef struct _p_MatNullSpace* MatNullSpace; 1479d9274352SBarry Smith 1480014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceCreate(MPI_Comm,PetscBool ,PetscInt,const Vec[],MatNullSpace*); 1481014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(MatNullSpace,Vec,void*),void*); 1482014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceDestroy(MatNullSpace*); 1483d0195637SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceRemove(MatNullSpace,Vec); 1484014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNullSpace(Mat, MatNullSpace *); 14855fa7ec2dSBarry Smith PETSC_EXTERN PetscErrorCode MatGetTransposeNullSpace(Mat, MatNullSpace *); 14865fa7ec2dSBarry Smith PETSC_EXTERN PetscErrorCode MatSetTransposeNullSpace(Mat,MatNullSpace); 1487014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNullSpace(Mat,MatNullSpace); 1488014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNearNullSpace(Mat,MatNullSpace); 1489014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNearNullSpace(Mat,MatNullSpace*); 1490014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceTest(MatNullSpace,Mat,PetscBool *); 1491014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceView(MatNullSpace,PetscViewer); 1492014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceGetVecs(MatNullSpace,PetscBool*,PetscInt*,const Vec**); 1493014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceCreateRigidBody(Vec,MatNullSpace*); 149474637425SBarry Smith 1495014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatReorderingSeqSBAIJ(Mat,IS); 1496014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetHashTableFactor(Mat,PetscReal); 1497014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetColumnIndices(Mat,PetscInt *); 1498014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJInvertBlockDiagonal(Mat); 14993f1d51d7SBarry Smith 1500014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMAIJ(Mat,PetscInt,Mat*); 1501014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMAIJRedimension(Mat,PetscInt,Mat*); 1502014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMAIJGetAIJ(Mat,Mat*); 1503c4f061fbSSatish Balay 1504014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatComputeExplicitOperator(Mat,Mat*); 1505b0a32e0cSBarry Smith 1506014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalScaleLocal(Mat,Vec); 150704f1ad80SBarry Smith 1508014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*); 1509014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetBase(Mat,Vec,Vec); 1510014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*); 1511014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*)); 1512014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec)); 1513014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt); 1514014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDResetHHistory(Mat); 1515014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctionError(Mat,PetscReal); 1516014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetPeriod(Mat,PetscInt); 1517014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDGetH(Mat,PetscScalar *); 1518014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetOptionsPrefix(Mat,const char[]); 1519014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*); 1520014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*); 1521e884886fSBarry Smith 15226370053bSBarry Smith /*S 15236370053bSBarry Smith MatMFFD - A data structured used to manage the computation of the h differencing parameter for matrix-free 15246370053bSBarry Smith Jacobian vector products 1525e884886fSBarry Smith 15266370053bSBarry Smith Notes: MATMFFD is a specific MatType which uses the MatMFFD data structure 15276370053bSBarry Smith 15286370053bSBarry Smith MatMFFD*() methods actually take the Mat as their first argument. Not a MatMFFD data structure 15296370053bSBarry Smith 15306370053bSBarry Smith Level: developer 15316370053bSBarry Smith 15326370053bSBarry Smith .seealso: MATMFFD, MatCreateMFFD(), MatMFFDSetFuction(), MatMFFDSetType(), MatMFFDRegister() 15336370053bSBarry Smith S*/ 1534e884886fSBarry Smith typedef struct _p_MatMFFD* MatMFFD; 1535e884886fSBarry Smith 153676bdecfbSBarry Smith /*J 1537e884886fSBarry Smith MatMFFDType - algorithm used to compute the h used in computing matrix-vector products via differencing of the function 1538e884886fSBarry Smith 1539e884886fSBarry Smith Level: beginner 1540e884886fSBarry Smith 1541e884886fSBarry Smith .seealso: MatMFFDSetType(), MatMFFDRegister() 154276bdecfbSBarry Smith J*/ 154319fd82e9SBarry Smith typedef const char* MatMFFDType; 1544e884886fSBarry Smith #define MATMFFD_DS "ds" 1545e884886fSBarry Smith #define MATMFFD_WP "wp" 1546e884886fSBarry Smith 154719fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDSetType(Mat,MatMFFDType); 1548bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDRegister(const char[],PetscErrorCode (*)(MatMFFD)); 1549e884886fSBarry Smith 1550014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDDSSetUmin(Mat,PetscReal); 1551014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDWPSetComputeNormU(Mat,PetscBool ); 1552e884886fSBarry Smith 155361ab5df0SBarry Smith PETSC_EXTERN PetscErrorCode MatFDColoringSetType(MatFDColoring,MatMFFDType); 155461ab5df0SBarry Smith 1555014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *); 1556014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *); 15577dbadf16SMatthew Knepley 155897969023SHong Zhang /* 155997969023SHong Zhang PETSc interface to MUMPS 156097969023SHong Zhang */ 156197969023SHong Zhang #ifdef PETSC_HAVE_MUMPS 1562014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMumpsSetIcntl(Mat,PetscInt,PetscInt); 1563bc6112feSHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetIcntl(Mat,PetscInt,PetscInt*); 1564b9e9ea26SSatish Balay PETSC_EXTERN PetscErrorCode MatMumpsSetCntl(Mat,PetscInt,PetscReal); 1565bc6112feSHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetCntl(Mat,PetscInt,PetscReal*); 1566bc6112feSHong Zhang 1567ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetInfo(Mat,PetscInt,PetscInt*); 1568ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetInfog(Mat,PetscInt,PetscInt*); 1569ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetRinfo(Mat,PetscInt,PetscReal*); 1570ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetRinfog(Mat,PetscInt,PetscReal*); 157197969023SHong Zhang #endif 157223a5497aSJed Brown 1573d954db57SHong Zhang /* 1574d615f992SSatish Balay PETSc interface to Mkl_Pardiso 1575b500a6b7SJose David Bermeol */ 1576b500a6b7SJose David Bermeol #ifdef PETSC_HAVE_MKL_PARDISO 1577d615f992SSatish Balay PETSC_EXTERN PetscErrorCode MatMkl_PardisoSetCntl(Mat,PetscInt,PetscInt); 1578b500a6b7SJose David Bermeol #endif 1579b500a6b7SJose David Bermeol 1580b500a6b7SJose David Bermeol /* 1581d305a81bSVasiliy Kozyrev PETSc interface to Mkl_CPardiso 1582d305a81bSVasiliy Kozyrev */ 1583d305a81bSVasiliy Kozyrev #ifdef PETSC_HAVE_MKL_CPARDISO 1584d305a81bSVasiliy Kozyrev PETSC_EXTERN PetscErrorCode MatMkl_CPardisoSetCntl(Mat,PetscInt,PetscInt); 1585d305a81bSVasiliy Kozyrev #endif 1586d305a81bSVasiliy Kozyrev 1587d305a81bSVasiliy Kozyrev /* 1588d954db57SHong Zhang PETSc interface to SUPERLU 1589d954db57SHong Zhang */ 1590d954db57SHong Zhang #ifdef PETSC_HAVE_SUPERLU 1591014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSuperluSetILUDropTol(Mat,PetscReal); 1592d954db57SHong Zhang #endif 1593d954db57SHong Zhang 1594fb785984SHong Zhang /* 1595fb785984SHong Zhang PETSc interface to SUPERLU_DIST 1596fb785984SHong Zhang */ 1597fb785984SHong Zhang #ifdef PETSC_HAVE_SUPERLU_DIST 1598fb785984SHong Zhang PETSC_EXTERN PetscErrorCode MatSuperluDistGetDiagU(Mat,PetscScalar*); 1599fb785984SHong Zhang #endif 1600fb785984SHong Zhang 1601575704cbSPieter Ghysels /* 1602575704cbSPieter Ghysels PETSc interface to STRUMPACK 1603575704cbSPieter Ghysels */ 1604575704cbSPieter Ghysels #ifdef PETSC_HAVE_STRUMPACK 1605575704cbSPieter Ghysels PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSRelCompTol(Mat,PetscReal); 1606575704cbSPieter Ghysels PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSMinSize(Mat,PetscInt); 1607575704cbSPieter Ghysels PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetColPerm(Mat,PetscBool); 1608575704cbSPieter Ghysels #endif 1609575704cbSPieter Ghysels 1610575704cbSPieter Ghysels 1611aa372e3fSPaul Mullowney #ifdef PETSC_HAVE_CUDA 16123f9c0db1SPaul Mullowney /*E 1613e057df02SPaul Mullowney MatCUSPARSEStorageFormat - indicates the storage format for CUSPARSE (GPU) 16142692e278SPaul Mullowney matrices. 1615e057df02SPaul Mullowney 1616e057df02SPaul Mullowney Not Collective 1617e057df02SPaul Mullowney 16188468deeeSKarl Rupp + MAT_CUSPARSE_CSR - Compressed Sparse Row 16192692e278SPaul Mullowney . MAT_CUSPARSE_ELL - Ellpack (requires CUDA 4.2 or later). 16202692e278SPaul Mullowney - MAT_CUSPARSE_HYB - Hybrid, a combination of Ellpack and Coordinate format (requires CUDA 4.2 or later). 1621e057df02SPaul Mullowney 1622e057df02SPaul Mullowney Level: intermediate 1623e057df02SPaul Mullowney 1624af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 1625e057df02SPaul Mullowney 1626e057df02SPaul Mullowney .seealso: MatCUSPARSESetFormat(), MatCUSPARSEFormatOperation 1627e057df02SPaul Mullowney E*/ 1628e057df02SPaul Mullowney 16293f9c0db1SPaul Mullowney typedef enum {MAT_CUSPARSE_CSR, MAT_CUSPARSE_ELL, MAT_CUSPARSE_HYB} MatCUSPARSEStorageFormat; 1630e057df02SPaul Mullowney 1631e057df02SPaul Mullowney /* these will be strings associated with enumerated type defined above */ 1632e057df02SPaul Mullowney PETSC_EXTERN const char *const MatCUSPARSEStorageFormats[]; 1633e057df02SPaul Mullowney 16343f9c0db1SPaul Mullowney /*E 1635e057df02SPaul Mullowney MatCUSPARSEFormatOperation - indicates the operation of CUSPARSE (GPU) 16362692e278SPaul Mullowney matrices whose operation should use a particular storage format. 1637e057df02SPaul Mullowney 1638e057df02SPaul Mullowney Not Collective 1639e057df02SPaul Mullowney 16408468deeeSKarl Rupp + MAT_CUSPARSE_MULT_DIAG - sets the storage format for the diagonal matrix in the parallel MatMult 16418468deeeSKarl Rupp . MAT_CUSPARSE_MULT_OFFDIAG - sets the storage format for the offdiagonal matrix in the parallel MatMult 16428468deeeSKarl Rupp . MAT_CUSPARSE_MULT - sets the storage format for the entire matrix in the serial (single GPU) MatMult 16438468deeeSKarl Rupp - MAT_CUSPARSE_ALL - sets the storage format for all CUSPARSE (GPU) matrices 1644e057df02SPaul Mullowney 1645e057df02SPaul Mullowney Level: intermediate 1646e057df02SPaul Mullowney 1647e057df02SPaul Mullowney .seealso: MatCUSPARSESetFormat(), MatCUSPARSEStorageFormat 1648e057df02SPaul Mullowney E*/ 164936d62e41SPaul Mullowney typedef enum {MAT_CUSPARSE_MULT_DIAG, MAT_CUSPARSE_MULT_OFFDIAG, MAT_CUSPARSE_MULT, MAT_CUSPARSE_ALL} MatCUSPARSEFormatOperation; 1650e057df02SPaul Mullowney 165110e9db80SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 165210e9db80SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCreateAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 1653e057df02SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCUSPARSESetFormat(Mat,MatCUSPARSEFormatOperation,MatCUSPARSEStorageFormat); 16549ae82921SPaul Mullowney #endif 16559ae82921SPaul Mullowney 165690c53902SBarry Smith #if defined(PETSC_HAVE_CUSP) 1657014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCUSP(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 1658014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateAIJCUSP(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 1659e057df02SPaul Mullowney 16603f9c0db1SPaul Mullowney /*E 1661e057df02SPaul Mullowney MatCUSPStorageFormat - indicates the storage format for CUSP (GPU) 166236d62e41SPaul Mullowney matrices. 1663e057df02SPaul Mullowney 1664e057df02SPaul Mullowney Not Collective 1665e057df02SPaul Mullowney 16668468deeeSKarl Rupp + MAT_CUSP_CSR - Compressed Sparse Row 16678468deeeSKarl Rupp . MAT_CUSP_DIA - Diagonal 16688468deeeSKarl Rupp - MAT_CUSP_ELL - Ellpack 1669e057df02SPaul Mullowney 1670e057df02SPaul Mullowney Level: intermediate 1671e057df02SPaul Mullowney 1672af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 1673e057df02SPaul Mullowney 1674e057df02SPaul Mullowney .seealso: MatCUSPSetFormat(), MatCUSPFormatOperation 1675e057df02SPaul Mullowney E*/ 16763f9c0db1SPaul Mullowney typedef enum {MAT_CUSP_CSR, MAT_CUSP_DIA, MAT_CUSP_ELL} MatCUSPStorageFormat; 1677e057df02SPaul Mullowney 1678e057df02SPaul Mullowney /* these will be strings associated with enumerated type defined above */ 1679e057df02SPaul Mullowney PETSC_EXTERN const char *const MatCUSPStorageFormats[]; 1680e057df02SPaul Mullowney 16813f9c0db1SPaul Mullowney /*E 1682e057df02SPaul Mullowney MatCUSPFormatOperation - indicates the operation of CUSP (GPU) 168336d62e41SPaul Mullowney matrices whose operation should use a particular storage format. 1684e057df02SPaul Mullowney 1685e057df02SPaul Mullowney Not Collective 1686e057df02SPaul Mullowney 16878468deeeSKarl Rupp + MAT_CUSP_MULT_DIAG - sets the storage format for the diagonal matrix in the parallel MatMult 16888468deeeSKarl Rupp . MAT_CUSP_MULT_OFFDIAG - sets the storage format for the offdiagonal matrix in the parallel MatMult 16898468deeeSKarl Rupp . MAT_CUSP_MULT - sets the storage format for the entire matrix in the serial (single GPU) MatMult 16908468deeeSKarl Rupp - MAT_CUSP_ALL - sets the storage format for all CUSP (GPU) matrices 1691e057df02SPaul Mullowney 1692e057df02SPaul Mullowney Level: intermediate 1693e057df02SPaul Mullowney 1694af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 1695e057df02SPaul Mullowney 1696e057df02SPaul Mullowney .seealso: MatCUSPSetFormat(), MatCUSPStorageFormat 1697e057df02SPaul Mullowney E*/ 169836d62e41SPaul Mullowney typedef enum {MAT_CUSP_MULT_DIAG, MAT_CUSP_MULT_OFFDIAG, MAT_CUSP_MULT, MAT_CUSP_ALL} MatCUSPFormatOperation; 1699e057df02SPaul Mullowney 1700e057df02SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCUSPSetFormat(Mat,MatCUSPFormatOperation,MatCUSPStorageFormat); 1701e057df02SPaul Mullowney #endif 170290c53902SBarry Smith 1703d67ff14aSKarl Rupp #if defined(PETSC_HAVE_VIENNACL) 1704d67ff14aSKarl Rupp PETSC_EXTERN PetscErrorCode MatCreateSeqAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 1705d67ff14aSKarl Rupp PETSC_EXTERN PetscErrorCode MatCreateAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 1706d67ff14aSKarl Rupp #endif 1707d67ff14aSKarl Rupp 170854efbe56SHong Zhang /* 170954efbe56SHong Zhang PETSc interface to FFTW 171054efbe56SHong Zhang */ 171154efbe56SHong Zhang #if defined(PETSC_HAVE_FFTW) 1712014dd563SJed Brown PETSC_EXTERN PetscErrorCode VecScatterPetscToFFTW(Mat,Vec,Vec); 1713014dd563SJed Brown PETSC_EXTERN PetscErrorCode VecScatterFFTWToPetsc(Mat,Vec,Vec); 17142a7a6963SBarry Smith PETSC_EXTERN PetscErrorCode MatCreateVecsFFTW(Mat,Vec*,Vec*,Vec*); 171554efbe56SHong Zhang #endif 171654efbe56SHong Zhang 1717382fd914SHong Zhang /* 1718382fd914SHong Zhang PETSc interface to ELEMENTAL 1719382fd914SHong Zhang */ 1720db31f6deSJed Brown #if defined(PETSC_HAVE_ELEMENTAL) 17212ef15aa2SHong Zhang #if defined(PETSC_USE_COMPLEX) 17222ef15aa2SHong Zhang typedef El::Complex<PetscReal> PetscElemScalar; 17232ef15aa2SHong Zhang #else 17242ef15aa2SHong Zhang typedef PetscScalar PetscElemScalar; 1725382fd914SHong Zhang #endif 1726607a6623SBarry Smith PETSC_EXTERN PetscErrorCode PetscElementalInitializePackage(void); 1727db31f6deSJed Brown PETSC_EXTERN PetscErrorCode PetscElementalFinalizePackage(void); 1728db31f6deSJed Brown #endif 1729db31f6deSJed Brown 1730014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateNest(MPI_Comm,PetscInt,const IS[],PetscInt,const IS[],const Mat[],Mat*); 1731014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSize(Mat,PetscInt*,PetscInt*); 1732014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetISs(Mat,IS[],IS[]); 1733014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetLocalISs(Mat,IS[],IS[]); 1734014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSubMats(Mat,PetscInt*,PetscInt*,Mat***); 1735014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSubMat(Mat,PetscInt,PetscInt,Mat*); 173619fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatNestSetVecType(Mat,VecType); 1737014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestSetSubMats(Mat,PetscInt,const IS[],PetscInt,const IS[],const Mat[]); 1738014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestSetSubMat(Mat,PetscInt,PetscInt,Mat); 1739d8588912SDave May 17404325cce7SMatthew G Knepley PETSC_EXTERN PetscErrorCode MatChop(Mat,PetscReal); 1741e0a58c10SMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatComputeBandwidth(Mat,PetscReal,PetscInt*); 17424325cce7SMatthew G Knepley 1743b541e6a4SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSubdomainsCreateCoalesce(Mat,PetscInt,PetscInt*,IS**); 174453dd7562SDmitry Karpeev 1745c094ef40SMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatPreallocatorPreallocate(Mat,PetscBool,Mat); 1746c094ef40SMatthew G. Knepley 1747539c167fSBarry Smith PETSC_INTERN PetscErrorCode MatHeaderMerge(Mat,Mat*); 1748539c167fSBarry Smith PETSC_EXTERN PetscErrorCode MatHeaderReplace(Mat,Mat*); 1749539c167fSBarry Smith 175023a5497aSJed Brown #endif 1751