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" 357d6a0e61SBarry Smith #define MATSEQAIJPTHREAD "seqaijpthread" 36bf2c1783SBarry Smith #define MATAIJPTHREAD "aijpthread" 37273d9f13SBarry Smith #define MATMPIAIJ "mpiaij" 385a11e1b2SBarry Smith #define MATAIJCRL "aijcrl" 395a11e1b2SBarry Smith #define MATSEQAIJCRL "seqaijcrl" 405a11e1b2SBarry Smith #define MATMPIAIJCRL "mpiaijcrl" 418154be41SBarry Smith #define MATAIJCUSP "aijcusp" 428154be41SBarry Smith #define MATSEQAIJCUSP "seqaijcusp" 438154be41SBarry Smith #define MATMPIAIJCUSP "mpiaijcusp" 449ae82921SPaul Mullowney #define MATAIJCUSPARSE "aijcusparse" 459ae82921SPaul Mullowney #define MATSEQAIJCUSPARSE "seqaijcusparse" 469ae82921SPaul Mullowney #define MATMPIAIJCUSPARSE "mpiaijcusparse" 47d67ff14aSKarl Rupp #define MATAIJVIENNACL "aijviennacl" 48d67ff14aSKarl Rupp #define MATSEQAIJVIENNACL "seqaijviennacl" 49d67ff14aSKarl Rupp #define MATMPIAIJVIENNACL "mpiaijviennacl" 505a11e1b2SBarry Smith #define MATAIJPERM "aijperm" 515a11e1b2SBarry Smith #define MATSEQAIJPERM "seqaijperm" 525a11e1b2SBarry Smith #define MATMPIAIJPERM "mpiaijperm" 53273d9f13SBarry Smith #define MATSHELL "shell" 545a11e1b2SBarry Smith #define MATDENSE "dense" 55209238afSKris Buschelman #define MATSEQDENSE "seqdense" 56273d9f13SBarry Smith #define MATMPIDENSE "mpidense" 57db31f6deSJed Brown #define MATELEMENTAL "elemental" 585a11e1b2SBarry Smith #define MATBAIJ "baij" 59273d9f13SBarry Smith #define MATSEQBAIJ "seqbaij" 60273d9f13SBarry Smith #define MATMPIBAIJ "mpibaij" 61273d9f13SBarry Smith #define MATMPIADJ "mpiadj" 625a11e1b2SBarry Smith #define MATSBAIJ "sbaij" 63273d9f13SBarry Smith #define MATSEQSBAIJ "seqsbaij" 64273d9f13SBarry Smith #define MATMPISBAIJ "mpisbaij" 65c0cdd4a1SDahai Guo #define MATSEQBSTRM "seqbstrm" 66c0cdd4a1SDahai Guo #define MATMPIBSTRM "mpibstrm" 67c0cdd4a1SDahai Guo #define MATBSTRM "bstrm" 68aa5a9175SDahai Guo #define MATSEQSBSTRM "seqsbstrm" 69aa5a9175SDahai Guo #define MATMPISBSTRM "mpisbstrm" 70aa5a9175SDahai Guo #define MATSBSTRM "sbstrm" 71cebc7f6cSBarry Smith #define MATDAAD "daad" 72cebc7f6cSBarry Smith #define MATMFFD "mffd" 73c8a8475eSBarry Smith #define MATNORMAL "normal" 74ab92ecdeSBarry Smith #define MATLRC "lrc" 752a6744ebSBarry Smith #define MATSCATTER "scatter" 76421e10b8SBarry Smith #define MATBLOCKMAT "blockmat" 77793850ffSBarry Smith #define MATCOMPOSITE "composite" 781f8c7532SHong Zhang #define MATFFT "fft" 791f8c7532SHong Zhang #define MATFFTW "fftw" 80e133240eSMatthew G Knepley #define MATSEQCUFFT "seqcufft" 81557cca28SSatish Balay #define MATTRANSPOSEMAT "transpose" 8272ca8751SBarry Smith #define MATSCHURCOMPLEMENT "schurcomplement" 831d6018f0SLisandro Dalcin #define MATPYTHON "python" 84f91d8e95SBarry Smith #define MATHYPRESTRUCT "hyprestruct" 85a9e6138eSGlenn Hammond #define MATHYPRESSTRUCT "hypresstruct" 86ee1cef2cSJed Brown #define MATSUBMATRIX "submatrix" 872c0dbf93SJed Brown #define MATLOCALREF "localref" 88d8588912SDave May #define MATNEST "nest" 89773941d6SBarry Smith 9076bdecfbSBarry Smith /*J 91c7393fdbSBarry Smith MatSolverPackage - String with the name of a PETSc matrix solver type. 929989ab13SBarry Smith 939989ab13SBarry Smith For example: "petsc" indicates what PETSc provides, "superlu" indicates either 949989ab13SBarry Smith SuperLU or SuperLU_Dist etc. 959989ab13SBarry Smith 969989ab13SBarry Smith 979989ab13SBarry Smith Level: beginner 989989ab13SBarry Smith 995c9eb25fSBarry Smith .seealso: MatGetFactor(), Mat, MatSetType(), MatType 10076bdecfbSBarry Smith J*/ 101c7393fdbSBarry Smith #define MatSolverPackage char* 1022692d6eeSBarry Smith #define MATSOLVERSUPERLU "superlu" 1032692d6eeSBarry Smith #define MATSOLVERSUPERLU_DIST "superlu_dist" 1042692d6eeSBarry Smith #define MATSOLVERUMFPACK "umfpack" 10520db9a53SJed Brown #define MATSOLVERCHOLMOD "cholmod" 1062692d6eeSBarry Smith #define MATSOLVERESSL "essl" 1072692d6eeSBarry Smith #define MATSOLVERLUSOL "lusol" 1082692d6eeSBarry Smith #define MATSOLVERMUMPS "mumps" 1092692d6eeSBarry Smith #define MATSOLVERPASTIX "pastix" 1102692d6eeSBarry Smith #define MATSOLVERMATLAB "matlab" 1112692d6eeSBarry Smith #define MATSOLVERPETSC "petsc" 1122692d6eeSBarry Smith #define MATSOLVERBAS "bas" 1139ae82921SPaul Mullowney #define MATSOLVERCUSPARSE "cusparse" 114aa5a9175SDahai Guo #define MATSOLVERBSTRM "bstrm" 115aa5a9175SDahai Guo #define MATSOLVERSBSTRM "sbstrm" 11615767789SHong Zhang #define MATSOLVERELEMENTAL "elemental" 11717f1a0eaSHong Zhang #define MATSOLVERCLIQUE "clique" 11869e15a41SShri Abhyankar #define MATSOLVERKLU "klu" 119c0cdd4a1SDahai Guo 120b24902e0SBarry Smith /*E 121b24902e0SBarry Smith MatFactorType - indicates what type of factorization is requested 122b24902e0SBarry Smith 123b24902e0SBarry Smith Level: beginner 124b24902e0SBarry Smith 125b24902e0SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 126b24902e0SBarry Smith 127c7393fdbSBarry Smith .seealso: MatSolverPackage, MatGetFactor() 128b24902e0SBarry Smith E*/ 129599ef60dSHong Zhang typedef enum {MAT_FACTOR_NONE, MAT_FACTOR_LU, MAT_FACTOR_CHOLESKY, MAT_FACTOR_ILU, MAT_FACTOR_ICC,MAT_FACTOR_ILUDT} MatFactorType; 130014dd563SJed Brown PETSC_EXTERN const char *const MatFactorTypes[]; 131e92e720dSBarry Smith 132014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactor(Mat,const MatSolverPackage,MatFactorType,Mat*); 133014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactorAvailable(Mat,const MatSolverPackage,MatFactorType,PetscBool *); 134014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFactorGetSolverPackage(Mat,const MatSolverPackage*); 135014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactorType(Mat,MatFactorType*); 1369989ab13SBarry Smith 137c06d978dSMatthew Knepley /* Logging support */ 1380700a824SBarry Smith #define MAT_FILE_CLASSID 1211216 /* used to indicate matrices in binary files */ 139014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_CLASSID; 140335efc43SPeter Brune PETSC_EXTERN PetscClassId MAT_COLORING_CLASSID; 141014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_FDCOLORING_CLASSID; 142014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_TRANSPOSECOLORING_CLASSID; 143014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_PARTITIONING_CLASSID; 144014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_COARSEN_CLASSID; 145014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_NULLSPACE_CLASSID; 146014dd563SJed Brown PETSC_EXTERN PetscClassId MATMFFD_CLASSID; 147c06d978dSMatthew Knepley 148ceb03754SKris Buschelman /*E 149ceb03754SKris Buschelman MatReuse - Indicates if matrices obtained from a previous call to MatGetSubMatrices() 150d6eff37eSBarry Smith or MatGetSubMatrix() are to be reused to store the new matrix values. For MatConvert() is used to indicate 151d6eff37eSBarry Smith that the input matrix is to be replaced with the converted matrix. 152ceb03754SKris Buschelman 153ceb03754SKris Buschelman Level: beginner 154ceb03754SKris Buschelman 155ceb03754SKris Buschelman Any additions/changes here MUST also be made in include/finclude/petscmat.h 156ceb03754SKris Buschelman 1570c451bc4SBarry Smith .seealso: MatGetSubMatrices(), MatGetSubMatrix(), MatDestroyMatrices(), MatConvert() 158ceb03754SKris Buschelman E*/ 159dfe085dbSJed Brown typedef enum {MAT_INITIAL_MATRIX,MAT_REUSE_MATRIX,MAT_IGNORE_MATRIX} MatReuse; 160ceb03754SKris Buschelman 1615494a064SHong Zhang /*E 1625494a064SHong Zhang MatGetSubMatrixOption - Indicates if matrices obtained from a call to MatGetSubMatrices() 163829201f2SHong Zhang include the matrix values. Currently it is only used by MatGetSeqNonzerostructure(). 1645494a064SHong Zhang 1655494a064SHong Zhang Level: beginner 1665494a064SHong Zhang 167829201f2SHong Zhang .seealso: MatGetSeqNonzerostructure() 1685494a064SHong Zhang E*/ 1695494a064SHong Zhang typedef enum {MAT_DO_NOT_GET_VALUES,MAT_GET_VALUES} MatGetSubMatrixOption; 1705494a064SHong Zhang 171607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatInitializePackage(void); 172c06d978dSMatthew Knepley 173014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreate(MPI_Comm,Mat*); 174014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetSizes(Mat,PetscInt,PetscInt,PetscInt,PetscInt); 17519fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatSetType(Mat,MatType); 176014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetFromOptions(Mat); 177ce1779c8SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatViewFromOptions(Mat A,const char prefix[],const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,prefix,name);} 178607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatRegisterAll(void); 179bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatRegister(const char[],PetscErrorCode(*)(Mat)); 180014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRegisterBaseName(const char[],const char[],const char[]); 181014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetOptionsPrefix(Mat,const char[]); 182014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAppendOptionsPrefix(Mat,const char[]); 183014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOptionsPrefix(Mat,const char*[]); 184f69a0ea3SMatthew Knepley 185014dd563SJed Brown PETSC_EXTERN PetscBool MatRegisterAllCalled; 186140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatList; 187140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatColoringList; 188140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatPartitioningList; 189140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatCoarsenList; 19028988994SBarry Smith 1913b224e63SBarry Smith /*E 192aa6c7ce3SBarry Smith MatStructure - Indicates if two matrices have the same nonzero structure 1933b224e63SBarry Smith 1943b224e63SBarry Smith Level: beginner 1953b224e63SBarry Smith 1963b224e63SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 1973b224e63SBarry Smith 198aa6c7ce3SBarry Smith .seealso: MatCopy(), MatAXPY() 1993b224e63SBarry Smith E*/ 200aa6c7ce3SBarry Smith typedef enum {DIFFERENT_NONZERO_PATTERN,SUBSET_NONZERO_PATTERN,SAME_NONZERO_PATTERN} MatStructure; 2013b224e63SBarry Smith 202014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqDense(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*); 203014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateDense(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*); 204014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 205014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 206014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *); 207014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithSplitArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],PetscInt[],PetscInt[],PetscScalar[],Mat*); 2088d7a6e47SBarry Smith 209014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 210014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 211014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat*); 212d21a29f3SJed Brown 213014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAdj(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscInt[],Mat*); 214014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 215d21a29f3SJed Brown 216014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 217014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPISBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *); 21838f409ebSLisandro Dalcin PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 219014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 2204cce697cSJed Brown PETSC_EXTERN PetscErrorCode MatXAIJSetPreallocation(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscInt[],const PetscInt[]); 221dfb205c3SBarry Smith 222014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateShell(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,void *,Mat*); 223014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateNormal(Mat,Mat*); 224014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateLRC(Mat,Mat,Mat,Mat*); 225014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateIS(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,ISLocalToGlobalMapping,Mat*); 226014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 227014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 228c0cdd4a1SDahai Guo 229014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 230014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 231014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqSBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 232014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPISBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 233c0cdd4a1SDahai Guo 234014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateScatter(MPI_Comm,VecScatter,Mat*); 235014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScatterSetVecScatter(Mat,VecScatter); 236014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScatterGetVecScatter(Mat,VecScatter*); 237014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateBlockMat(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,Mat*); 238014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeAddMat(Mat,Mat); 239014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeMerge(Mat); 240014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateComposite(MPI_Comm,PetscInt,const Mat*,Mat*); 2416d7c1e57SBarry Smith typedef enum {MAT_COMPOSITE_ADDITIVE,MAT_COMPOSITE_MULTIPLICATIVE} MatCompositeType; 242014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeSetType(Mat,MatCompositeType); 2436d7c1e57SBarry Smith 24419fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatCreateFFT(MPI_Comm,PetscInt,const PetscInt[],MatType,Mat*); 245014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqCUFFT(MPI_Comm,PetscInt,const PetscInt[],Mat*); 246dedccee8SHong Zhang 247014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateTranspose(Mat,Mat*); 248014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSubMatrix(Mat,IS,IS,Mat*); 249014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSubMatrixUpdate(Mat,Mat,IS,IS); 250014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateLocalRef(Mat,IS,IS,Mat*); 2511472f72bSBarry Smith 252014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPythonSetType(Mat,const char[]); 2531d6018f0SLisandro Dalcin 254014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetUp(Mat); 255014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroy(Mat*); 256e56f5c9eSBarry Smith PETSC_EXTERN PetscErrorCode MatGetNonzeroState(Mat,PetscObjectState*); 25721c89e3eSBarry Smith 258014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatConjugate(Mat); 259014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRealPart(Mat); 260014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatImaginaryPart(Mat); 261014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetDiagonalBlock(Mat,Mat*); 262014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetTrace(Mat,PetscScalar*); 263713ccfa9SJed Brown PETSC_EXTERN PetscErrorCode MatInvertBlockDiagonal(Mat,const PetscScalar **); 26499cafbc1SBarry Smith 2658ed539a5SBarry Smith /* ------------------------------------------------------------*/ 266014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 267014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 268014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesRow(Mat,PetscInt,const PetscScalar[]); 269014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesRowLocal(Mat,PetscInt,const PetscScalar[]); 270014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBatch(Mat,PetscInt,PetscInt,PetscInt[],const PetscScalar[]); 27173a71a0fSBarry Smith PETSC_EXTERN PetscErrorCode MatSetRandom(Mat,PetscRandom); 27284cb2905SBarry Smith 2732ef4de8bSBarry Smith /*S 2742ef4de8bSBarry Smith MatStencil - Data structure (C struct) for storing information about a single row or 275be479b30SJed Brown column of a matrix as indexed on an associated grid. 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 283be479b30SJed Brown .seealso: MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockedStencil() 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 293014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetColoring(Mat,ISColoring); 294014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesAdifor(Mat,PetscInt,void*); 2953a7fca6bSBarry Smith 296d91e6319SBarry Smith /*E 297d91e6319SBarry Smith MatAssemblyType - Indicates if the matrix is now to be used, or if you plan 298d91e6319SBarry Smith to continue to add values to it 299d91e6319SBarry Smith 300d91e6319SBarry Smith Level: beginner 301d91e6319SBarry Smith 302d91e6319SBarry Smith .seealso: MatAssemblyBegin(), MatAssemblyEnd() 303d91e6319SBarry Smith E*/ 3046d4a8577SBarry Smith typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType; 305014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssemblyBegin(Mat,MatAssemblyType); 306014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssemblyEnd(Mat,MatAssemblyType); 307014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssembled(Mat,PetscBool *); 3084f9c727eSBarry Smith 3091d73ed98SBarry Smith 31030de9b25SBarry Smith 311d91e6319SBarry Smith /*E 312d91e6319SBarry Smith MatOption - Options that may be set for a matrix and its behavior or storage 313d91e6319SBarry Smith 314d91e6319SBarry Smith Level: beginner 315d91e6319SBarry Smith 3160a835dfdSSatish Balay Any additions/changes here MUST also be made in include/finclude/petscmat.h 317d91e6319SBarry Smith 318382164b0SBarry Smith Developer Notes: Entries that are negative need not be called collectively by all processes. 319382164b0SBarry Smith 320d91e6319SBarry Smith .seealso: MatSetOption() 321d91e6319SBarry Smith E*/ 322382164b0SBarry Smith typedef enum {MAT_OPTION_MIN = -8, 323382164b0SBarry Smith MAT_NEW_NONZERO_LOCATION_ERR = -7, 324382164b0SBarry Smith MAT_NO_OFF_PROC_ZERO_ROWS = -6, 325382164b0SBarry Smith MAT_NO_OFF_PROC_ENTRIES = -5, 326382164b0SBarry Smith MAT_UNUSED_NONZERO_LOCATION_ERR = -4, 327382164b0SBarry Smith MAT_NEW_NONZERO_ALLOCATION_ERR = -3, 328382164b0SBarry Smith MAT_ROW_ORIENTED = -2, 329382164b0SBarry Smith MAT_NEW_NONZERO_LOCATIONS = -1, 330382164b0SBarry Smith MAT_SYMMETRIC = 1, 331382164b0SBarry Smith MAT_STRUCTURALLY_SYMMETRIC = 2, 332382164b0SBarry Smith MAT_NEW_DIAGONALS = 3, 333382164b0SBarry Smith MAT_IGNORE_OFF_PROC_ENTRIES = 4, 334382164b0SBarry Smith MAT_USE_HASH_TABLE = 5, 335382164b0SBarry Smith MAT_KEEP_NONZERO_PATTERN = 6, 336382164b0SBarry Smith MAT_IGNORE_ZERO_ENTRIES = 7, 337382164b0SBarry Smith MAT_USE_INODES = 8, 338382164b0SBarry Smith MAT_HERMITIAN = 9, 339382164b0SBarry Smith MAT_SYMMETRY_ETERNAL = 10, 34011e456e1SBarry Smith MAT_DUMMY = 11, 341382164b0SBarry Smith MAT_IGNORE_LOWER_TRIANGULAR = 12, 342382164b0SBarry Smith MAT_ERROR_LOWER_TRIANGULAR = 13, 343382164b0SBarry Smith MAT_GETROW_UPPERTRIANGULAR = 14, 344382164b0SBarry Smith MAT_SPD = 15, 345382164b0SBarry Smith MAT_OPTION_MAX = 16} MatOption; 346e057df02SPaul Mullowney 347014dd563SJed Brown PETSC_EXTERN const char *MatOptions[]; 348014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetOption(Mat,MatOption,PetscBool ); 34919fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatGetType(Mat,MatType*); 35084cb2905SBarry Smith 351014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]); 352014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 353014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 354014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowUpperTriangular(Mat); 355014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRowUpperTriangular(Mat); 356014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 357014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 358014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnVector(Mat,Vec,PetscInt); 3598c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJGetArray(Mat,PetscScalar *[]); 3608c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJRestoreArray(Mat,PetscScalar *[]); 361bd04181cSBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJSetValuesLocalFast(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 3628c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatDenseGetArray(Mat,PetscScalar *[]); 3638c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatDenseRestoreArray(Mat,PetscScalar *[]); 364014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBlockSize(Mat,PetscInt *); 365014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSize(Mat,PetscInt); 366014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBlockSizes(Mat,PetscInt *,PetscInt *); 367014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSizes(Mat,PetscInt,PetscInt); 36833d57670SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSizesFromMats(Mat,Mat,Mat); 369014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNThreads(Mat,PetscInt); 370014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNThreads(Mat,PetscInt*); 3711620fd73SBarry Smith 372014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMult(Mat,Vec,Vec); 373014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultDiagonalBlock(Mat,Vec,Vec); 374014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultAdd(Mat,Vec,Vec,Vec); 375014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTranspose(Mat,Vec,Vec); 376014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultHermitianTranspose(Mat,Vec,Vec); 377014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsTranspose(Mat,Mat,PetscReal,PetscBool *); 378014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitianTranspose(Mat,Mat,PetscReal,PetscBool *); 379014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeAdd(Mat,Vec,Vec,Vec); 380014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultHermitianTransposeAdd(Mat,Vec,Vec,Vec); 381014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultConstrained(Mat,Vec,Vec); 382014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeConstrained(Mat,Vec,Vec); 383014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatSolve(Mat,Mat,Mat); 384f9426fe0SMark Adams PETSC_EXTERN PetscErrorCode MatResidual(Mat,Vec,Vec,Vec); 385f5edf698SHong Zhang 386d91e6319SBarry Smith /*E 387d91e6319SBarry Smith MatDuplicateOption - Indicates if a duplicated sparse matrix should have 388d91e6319SBarry Smith its numerical values copied over or just its nonzero structure. 389d91e6319SBarry Smith 390d91e6319SBarry Smith Level: beginner 391d91e6319SBarry Smith 392d91e6319SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 393d91e6319SBarry Smith 39470dcbbb9SBarry Smith $ MAT_SHARE_NONZERO_PATTERN - the i and j arrays in the new matrix will be shared with the original matrix 39570dcbbb9SBarry Smith $ this also triggers the MAT_DO_NOT_COPY_VALUES option. This is used when you 39670dcbbb9SBarry Smith $ have several matrices with the same nonzero pattern. 39770dcbbb9SBarry Smith 398d91e6319SBarry Smith .seealso: MatDuplicate() 399d91e6319SBarry Smith E*/ 40070dcbbb9SBarry Smith typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES,MAT_SHARE_NONZERO_PATTERN} MatDuplicateOption; 4012e8a6d31SBarry Smith 40219fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatConvert(Mat,MatType,MatReuse,Mat*); 403014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDuplicate(Mat,MatDuplicateOption,Mat*); 40494a9d846SBarry Smith 405cb5b572fSBarry Smith 406014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCopy(Mat,Mat,MatStructure); 407014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatView(Mat,PetscViewer); 408014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsSymmetric(Mat,PetscReal,PetscBool *); 409014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsStructurallySymmetric(Mat,PetscBool *); 410014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitian(Mat,PetscReal,PetscBool *); 411014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsSymmetricKnown(Mat,PetscBool *,PetscBool *); 412014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitianKnown(Mat,PetscBool *,PetscBool *); 413014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMissingDiagonal(Mat,PetscBool *,PetscInt *); 414014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLoad(Mat, PetscViewer); 4157b80b807SBarry Smith 4161a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool *); 4171a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool *); 4181a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool *); 4191a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool *); 420d4fbbf0eSBarry Smith 421d91e6319SBarry Smith /*S 422d91e6319SBarry Smith MatInfo - Context of matrix information, used with MatGetInfo() 423d91e6319SBarry Smith 424d91e6319SBarry Smith In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE 425d91e6319SBarry Smith 426d91e6319SBarry Smith Level: intermediate 427d91e6319SBarry Smith 428d91e6319SBarry Smith Concepts: matrix^nonzero information 429d91e6319SBarry Smith 430d9274352SBarry Smith .seealso: MatGetInfo(), MatInfoType 431d91e6319SBarry Smith S*/ 4324e220ebcSLois Curfman McInnes typedef struct { 433b0a32e0cSBarry Smith PetscLogDouble block_size; /* block size */ 434b0a32e0cSBarry Smith PetscLogDouble nz_allocated,nz_used,nz_unneeded; /* number of nonzeros */ 435b0a32e0cSBarry Smith PetscLogDouble memory; /* memory allocated */ 436b0a32e0cSBarry Smith PetscLogDouble assemblies; /* number of matrix assemblies called */ 437b0a32e0cSBarry Smith PetscLogDouble mallocs; /* number of mallocs during MatSetValues() */ 438b0a32e0cSBarry Smith PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */ 439b0a32e0cSBarry Smith PetscLogDouble factor_mallocs; /* number of mallocs during factorization */ 4404e220ebcSLois Curfman McInnes } MatInfo; 4414e220ebcSLois Curfman McInnes 442d9274352SBarry Smith /*E 443d9274352SBarry Smith MatInfoType - Indicates if you want information about the local part of the matrix, 444d9274352SBarry Smith the entire parallel matrix or the maximum over all the local parts. 445d9274352SBarry Smith 446d9274352SBarry Smith Level: beginner 447d9274352SBarry Smith 448d9274352SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 449d9274352SBarry Smith 450d9274352SBarry Smith .seealso: MatGetInfo(), MatInfo 451d9274352SBarry Smith E*/ 4527b80b807SBarry Smith typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType; 453014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetInfo(Mat,MatInfoType,MatInfo*); 454014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetDiagonal(Mat,Vec); 455014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMax(Mat,Vec,PetscInt[]); 456014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMin(Mat,Vec,PetscInt[]); 457014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMaxAbs(Mat,Vec,PetscInt[]); 458014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMinAbs(Mat,Vec,PetscInt[]); 459014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowSum(Mat,Vec); 460014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTranspose(Mat,MatReuse,Mat*); 461014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatHermitianTranspose(Mat,MatReuse,Mat*); 462014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPermute(Mat,IS,IS,Mat *); 463014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalScale(Mat,Vec,Vec); 464014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalSet(Mat,Vec,InsertMode); 465014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatEqual(Mat,Mat,PetscBool *); 466014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultEqual(Mat,Mat,PetscInt,PetscBool *); 467014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultAddEqual(Mat,Mat,PetscInt,PetscBool *); 468014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeEqual(Mat,Mat,PetscInt,PetscBool *); 469014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscBool *); 4707b80b807SBarry Smith 471014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNorm(Mat,NormType,PetscReal *); 472014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnNorms(Mat,NormType,PetscReal *); 473014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroEntries(Mat); 474014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 475014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsIS(Mat,IS,PetscScalar,Vec,Vec); 476014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsStencil(Mat,PetscInt,const MatStencil [],PetscScalar,Vec,Vec); 477014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsStencil(Mat,PetscInt,const MatStencil[],PetscScalar,Vec,Vec); 478014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumns(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 479014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsIS(Mat,IS,PetscScalar,Vec,Vec); 4807b80b807SBarry Smith 481014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSize(Mat,PetscInt*,PetscInt*); 482014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalSize(Mat,PetscInt*,PetscInt*); 483014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRange(Mat,PetscInt*,PetscInt*); 484014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRanges(Mat,const PetscInt**); 485014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRangeColumn(Mat,PetscInt*,PetscInt*); 486014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRangesColumn(Mat,const PetscInt**); 487566876eaSJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipIS(Mat,IS*,IS*); 4885ef9f2a5SBarry Smith 489014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]); 490014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSubMatricesParallel(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]); 491014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroyMatrices(PetscInt,Mat *[]); 492014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSubMatrix(Mat,IS,IS,MatReuse,Mat *); 493014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalSubMatrix(Mat,IS,IS,Mat*); 494014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreLocalSubMatrix(Mat,IS,IS,Mat*); 495014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSeqNonzeroStructure(Mat,Mat*); 496014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroySeqNonzeroStructure(Mat*); 4977b80b807SBarry Smith 498014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJConcatenateSeqAIJ(MPI_Comm,Mat,PetscInt,MatReuse,Mat*); 499014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJConcatenateSeqAIJSymbolic(MPI_Comm,Mat,PetscInt,Mat*); 500014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJConcatenateSeqAIJNumeric(MPI_Comm,Mat,PetscInt,Mat); 501014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJ(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*); 502014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJSymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*); 503014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJNumeric(Mat,Mat); 504014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMat(Mat,MatReuse,Mat*); 505014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*); 506014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,Mat*); 507014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetGhosts(Mat, PetscInt *,const PetscInt *[]); 5088efafbd8SBarry Smith 509014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt); 5107b80b807SBarry Smith 511014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*); 512014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMultSymbolic(Mat,Mat,PetscReal,Mat*); 513014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMultNumeric(Mat,Mat,Mat); 51422440eb1SKris Buschelman 5157bab7c10SHong Zhang PETSC_EXTERN PetscErrorCode MatMatMatMult(Mat,Mat,Mat,MatReuse,PetscReal,Mat*); 5167bab7c10SHong Zhang PETSC_EXTERN PetscErrorCode MatMatMatMultSymbolic(Mat,Mat,Mat,PetscReal,Mat*); 5177bab7c10SHong Zhang PETSC_EXTERN PetscErrorCode MatMatMatMultNumeric(Mat,Mat,Mat,Mat); 5187bab7c10SHong Zhang 519014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*); 520014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAPSymbolic(Mat,Mat,PetscReal,Mat*); 521014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAPNumeric(Mat,Mat,Mat); 522014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARt(Mat,Mat,MatReuse,PetscReal,Mat*); 523014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARtSymbolic(Mat,Mat,PetscReal,Mat*); 524014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARtNumeric(Mat,Mat,Mat); 52522440eb1SKris Buschelman 526014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeMatMult(Mat,Mat,MatReuse,PetscReal,Mat*); 527014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposetMatMultSymbolic(Mat,Mat,PetscReal,Mat*); 528014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposetMatMultNumeric(Mat,Mat,Mat); 529014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMult(Mat,Mat,MatReuse,PetscReal,Mat*); 530014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMultSymbolic(Mat,Mat,PetscReal,Mat*); 531014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMultNumeric(Mat,Mat,Mat); 532bc011b1eSHong Zhang 533014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAXPY(Mat,PetscScalar,Mat,MatStructure); 534014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAYPX(Mat,PetscScalar,Mat,MatStructure); 5351c741599SHong Zhang 536014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScale(Mat,PetscScalar); 537014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShift(Mat,PetscScalar); 5387b80b807SBarry Smith 539014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping); 540014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetLocalToGlobalMappingBlock(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping); 541014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*); 542014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalToGlobalMappingBlock(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*); 543a93ff8c4SPeter Brune PETSC_EXTERN PetscErrorCode MatGetLayouts(Mat,PetscLayout*,PetscLayout*); 544014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 545014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsLocalIS(Mat,IS,PetscScalar,Vec,Vec); 546014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 547014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocalIS(Mat,IS,PetscScalar,Vec,Vec); 548014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 549014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 550052efed2SBarry Smith 551014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStashSetInitialSize(Mat,PetscInt,PetscInt); 552014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*); 55390f02eecSBarry Smith 554014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInterpolate(Mat,Vec,Vec); 555014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInterpolateAdd(Mat,Vec,Vec,Vec); 556014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestrict(Mat,Vec,Vec); 557014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetVecs(Mat,Vec*,Vec*); 558b2bf6370SHong Zhang PETSC_EXTERN PetscErrorCode MatGetRedundantMatrix(Mat,PetscInt,MPI_Comm,MatReuse,Mat*); 559014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetMultiProcBlock(Mat,MPI_Comm,MatReuse,Mat*); 560014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFindZeroDiagonals(Mat,IS*); 5613a062f41SBarry Smith PETSC_EXTERN PetscErrorCode MatFindOffBlockDiagonalEntries(Mat,IS*); 562bd481603SBarry Smith 563bd481603SBarry Smith /*MC 5641620fd73SBarry Smith MatSetValue - Set a single entry into a matrix. 5651620fd73SBarry Smith 5661620fd73SBarry Smith Not collective 5671620fd73SBarry Smith 5681620fd73SBarry Smith Input Parameters: 5691620fd73SBarry Smith + m - the matrix 5701620fd73SBarry Smith . row - the row location of the entry 5711620fd73SBarry Smith . col - the column location of the entry 5721620fd73SBarry Smith . value - the value to insert 5731620fd73SBarry Smith - mode - either INSERT_VALUES or ADD_VALUES 5741620fd73SBarry Smith 5751620fd73SBarry Smith Notes: 5761620fd73SBarry Smith For efficiency one should use MatSetValues() and set several or many 5771620fd73SBarry Smith values simultaneously if possible. 5781620fd73SBarry Smith 5791620fd73SBarry Smith Level: beginner 5801620fd73SBarry Smith 5811620fd73SBarry Smith .seealso: MatSetValues(), MatSetValueLocal() 5821620fd73SBarry Smith M*/ 5831620fd73SBarry 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);} 5841620fd73SBarry Smith 5851620fd73SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatGetValue(Mat v,PetscInt i,PetscInt j,PetscScalar *va) {return MatGetValues(v,1,&i,1,&j,va);} 5861620fd73SBarry Smith 5871620fd73SBarry 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);} 5881620fd73SBarry Smith 5891620fd73SBarry Smith /*MC 5900d2a7ff2SSatish Balay MatPreallocateInitialize - Begins the block of code that will count the number of nonzeros per 591bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 592bd481603SBarry Smith 593bd481603SBarry Smith Synopsis: 594aaa7dc30SBarry Smith #include <petscmat.h> 595c1ac3661SBarry Smith PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz) 596bd481603SBarry Smith 597bd481603SBarry Smith Collective on MPI_Comm 598bd481603SBarry Smith 599bd481603SBarry Smith Input Parameters: 600bd481603SBarry Smith + comm - the communicator that will share the eventually allocated matrix 601859fcb39SBarry Smith . nrows - the number of LOCAL rows in the matrix 602859fcb39SBarry Smith - ncols - the number of LOCAL columns in the matrix 603bd481603SBarry Smith 604bd481603SBarry Smith Output Parameters: 605bd481603SBarry Smith + dnz - the array that will be passed to the matrix preallocation routines 606bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 607bd481603SBarry Smith 608bd481603SBarry Smith Level: intermediate 609bd481603SBarry Smith 610bd481603SBarry Smith Notes: 6110598bfebSBarry Smith See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details. 612bd481603SBarry Smith 6131d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 614bd481603SBarry Smith 615bd481603SBarry Smith Use MatPreallocateInitializeSymmetric() for symmetric matrices (MPISBAIJ matrices) 616bd481603SBarry Smith 6171620fd73SBarry Smith This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize(). 6181620fd73SBarry Smith 619bd481603SBarry Smith Concepts: preallocation^Matrix 620bd481603SBarry Smith 621bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(), 622bd481603SBarry Smith MatPreallocateInitializeSymmetric(), MatPreallocateSymmetricSetLocal() 623bd481603SBarry Smith M*/ 624c4f061fbSSatish Balay #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \ 6257c922b88SBarry Smith { \ 626a89bc211SBarry Smith PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__start,__end; \ 6271795a4d1SJed Brown _4_ierr = PetscCalloc2(__nrows,&dnz,__nrows,&onz);CHKERRQ(_4_ierr); \ 6281795a4d1SJed Brown __start = 0; __end = __start; \ 629c1ac3661SBarry Smith _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ctmp;\ 630a89bc211SBarry Smith _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows; 6317c922b88SBarry Smith 632bd481603SBarry Smith /*MC 633bd481603SBarry Smith MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be 634bd481603SBarry Smith inserted using a local number of the rows and columns 635bd481603SBarry Smith 636bd481603SBarry Smith Synopsis: 637aaa7dc30SBarry Smith #include <petscmat.h> 638c1ac3661SBarry Smith PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 639bd481603SBarry Smith 640bd481603SBarry Smith Not Collective 641bd481603SBarry Smith 642bd481603SBarry Smith Input Parameters: 643784ac674SJed Brown + map - the row mapping from local numbering to global numbering 644bd481603SBarry Smith . nrows - the number of rows indicated 6451d73ed98SBarry Smith . rows - the indices of the rows 646784ac674SJed Brown . cmap - the column mapping from local to global numbering 647bd481603SBarry Smith . ncols - the number of columns in the matrix 648bd481603SBarry Smith . cols - the columns indicated 649bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 650bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 651bd481603SBarry Smith 652bd481603SBarry Smith Level: intermediate 653bd481603SBarry Smith 654bd481603SBarry Smith Notes: 6550598bfebSBarry Smith See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details. 656bd481603SBarry Smith 6571d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 658bd481603SBarry Smith 659bd481603SBarry Smith Concepts: preallocation^Matrix 660bd481603SBarry Smith 661bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(), 662bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal() 663bd481603SBarry Smith M*/ 664784ac674SJed Brown #define MatPreallocateSetLocal(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \ 665c4f061fbSSatish Balay {\ 666c1ac3661SBarry Smith PetscInt __l;\ 667784ac674SJed Brown _4_ierr = ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\ 668784ac674SJed Brown _4_ierr = ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\ 669c4f061fbSSatish Balay for (__l=0;__l<nrows;__l++) {\ 670ef66eb69SBarry Smith _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 671c4f061fbSSatish Balay }\ 672c4f061fbSSatish Balay } 673c4f061fbSSatish Balay 674bd481603SBarry Smith /*MC 675bd481603SBarry Smith MatPreallocateSymmetricSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be 676bd481603SBarry Smith inserted using a local number of the rows and columns 677bd481603SBarry Smith 678bd481603SBarry Smith Synopsis: 679aaa7dc30SBarry Smith #include <petscmat.h> 680c1ac3661SBarry Smith PetscErrorCode MatPreallocateSymmetricSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 681bd481603SBarry Smith 682bd481603SBarry Smith Not Collective 683bd481603SBarry Smith 684bd481603SBarry Smith Input Parameters: 685bd481603SBarry Smith + map - the mapping between local numbering and global numbering 686bd481603SBarry Smith . nrows - the number of rows indicated 6871d73ed98SBarry Smith . rows - the indices of the rows 688bd481603SBarry Smith . ncols - the number of columns in the matrix 689bd481603SBarry Smith . cols - the columns indicated 690bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 691bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 692bd481603SBarry Smith 693bd481603SBarry Smith Level: intermediate 694bd481603SBarry Smith 695bd481603SBarry Smith Notes: 6960598bfebSBarry Smith See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details. 697bd481603SBarry Smith 698bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 699bd481603SBarry Smith 700bd481603SBarry Smith Concepts: preallocation^Matrix 701bd481603SBarry Smith 702bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(), 703bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal() 704bd481603SBarry Smith M*/ 705d3d32019SBarry Smith #define MatPreallocateSymmetricSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\ 706d3d32019SBarry Smith {\ 707c1ac3661SBarry Smith PetscInt __l;\ 708d3d32019SBarry Smith _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\ 709d3d32019SBarry Smith _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\ 710d3d32019SBarry Smith for (__l=0;__l<nrows;__l++) {\ 711d3d32019SBarry Smith _4_ierr = MatPreallocateSymmetricSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 712d3d32019SBarry Smith }\ 713d3d32019SBarry Smith } 714d3d32019SBarry Smith 715bd481603SBarry Smith /*MC 716bd481603SBarry Smith MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be 717bd481603SBarry Smith inserted using a local number of the rows and columns 718bd481603SBarry Smith 719bd481603SBarry Smith Synopsis: 720aaa7dc30SBarry Smith #include <petscmat.h> 721c1ac3661SBarry Smith PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 722bd481603SBarry Smith 723bd481603SBarry Smith Not Collective 724bd481603SBarry Smith 725bd481603SBarry Smith Input Parameters: 72664075487SBarry Smith + row - the row 727bd481603SBarry Smith . ncols - the number of columns in the matrix 728a50f8bf6SHong Zhang - cols - the columns indicated 729a50f8bf6SHong Zhang 730a50f8bf6SHong Zhang Output Parameters: 731a50f8bf6SHong Zhang + 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: 7370598bfebSBarry Smith See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details. 738bd481603SBarry Smith 739bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 740bd481603SBarry Smith 7411620fd73SBarry Smith This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize(). 7421620fd73SBarry Smith 743bd481603SBarry Smith Concepts: preallocation^Matrix 744bd481603SBarry Smith 745bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(), 746bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal() 747bd481603SBarry Smith M*/ 748c4f061fbSSatish Balay #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\ 749c1ac3661SBarry Smith { PetscInt __i; \ 750e32f2f54SBarry 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);\ 751e32f2f54SBarry 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);\ 7527c922b88SBarry Smith for (__i=0; __i<nc; __i++) {\ 75364075487SBarry Smith if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \ 7547cc688d7SBarry Smith else dnz[row - __rstart]++;\ 7557c922b88SBarry Smith }\ 7567c922b88SBarry Smith } 7577c922b88SBarry Smith 758bd481603SBarry Smith /*MC 759bd481603SBarry Smith MatPreallocateSymmetricSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be 760bd481603SBarry Smith inserted using a local number of the rows and columns 761bd481603SBarry Smith 762bd481603SBarry Smith Synopsis: 763aaa7dc30SBarry Smith #include <petscmat.h> 764c1ac3661SBarry Smith PetscErrorCode MatPreallocateSymmetricSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 765bd481603SBarry Smith 766bd481603SBarry Smith Not Collective 767bd481603SBarry Smith 768bd481603SBarry Smith Input Parameters: 769bd481603SBarry Smith + nrows - the number of rows indicated 7701d73ed98SBarry Smith . rows - the indices of the rows 771bd481603SBarry Smith . ncols - the number of columns in the matrix 772bd481603SBarry Smith . cols - the columns indicated 773bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 774bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 775bd481603SBarry Smith 776bd481603SBarry Smith Level: intermediate 777bd481603SBarry Smith 778bd481603SBarry Smith Notes: 7790598bfebSBarry Smith See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details. 780bd481603SBarry Smith 781bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 782bd481603SBarry Smith 7831620fd73SBarry Smith This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize(). 7841620fd73SBarry Smith 785bd481603SBarry Smith Concepts: preallocation^Matrix 786bd481603SBarry Smith 787bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(), 788bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal() 789bd481603SBarry Smith M*/ 790d3d32019SBarry Smith #define MatPreallocateSymmetricSet(row,nc,cols,dnz,onz) 0;\ 791c1ac3661SBarry Smith { PetscInt __i; \ 792d3d32019SBarry Smith for (__i=0; __i<nc; __i++) {\ 793d3d32019SBarry Smith if (cols[__i] >= __end) onz[row - __rstart]++; \ 794d3d32019SBarry Smith else if (cols[__i] >= row) dnz[row - __rstart]++;\ 795d3d32019SBarry Smith }\ 796d3d32019SBarry Smith } 797d3d32019SBarry Smith 798bd481603SBarry Smith /*MC 79916371a99SBarry Smith MatPreallocateLocation - An alternative to MatPreallocationSet() that puts the nonzero locations into the matrix if it exists 80016371a99SBarry Smith 80116371a99SBarry Smith Synopsis: 802aaa7dc30SBarry Smith #include <petscmat.h> 80316371a99SBarry Smith PetscErrorCode MatPreallocateLocations(Mat A,PetscInt row,PetscInt ncols,PetscInt *cols,PetscInt *dnz,PetscInt *onz) 80416371a99SBarry Smith 80516371a99SBarry Smith Not Collective 80616371a99SBarry Smith 80716371a99SBarry Smith Input Parameters: 80816371a99SBarry Smith . A - matrix 80916371a99SBarry Smith . row - row where values exist (must be local to this process) 81016371a99SBarry Smith . ncols - number of columns 81116371a99SBarry Smith . cols - columns with nonzeros 81216371a99SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 81316371a99SBarry Smith - ozn - the other array passed to the matrix preallocation routines 81416371a99SBarry Smith 81516371a99SBarry Smith Level: intermediate 81616371a99SBarry Smith 81716371a99SBarry Smith Notes: 8180598bfebSBarry Smith See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details. 81916371a99SBarry Smith 82016371a99SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 82116371a99SBarry Smith 82216371a99SBarry Smith This is a MACRO not a function because it uses a bunch of variables private to the MatPreallocation.... routines. 82316371a99SBarry Smith 82416371a99SBarry Smith Concepts: preallocation^Matrix 82516371a99SBarry Smith 82616371a99SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(), 827eabe889fSLisandro Dalcin MatPreallocateSymmetricSetLocal() 82816371a99SBarry Smith M*/ 8290298fd71SBarry 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);} 83016371a99SBarry Smith 83116371a99SBarry Smith 83216371a99SBarry Smith /*MC 8330d2a7ff2SSatish Balay MatPreallocateFinalize - Ends the block of code that will count the number of nonzeros per 834bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 835bd481603SBarry Smith 836bd481603SBarry Smith Synopsis: 837aaa7dc30SBarry Smith #include <petscmat.h> 838c1ac3661SBarry Smith PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz) 839bd481603SBarry Smith 840bd481603SBarry Smith Collective on MPI_Comm 841bd481603SBarry Smith 842bd481603SBarry Smith Input Parameters: 84316371a99SBarry Smith + dnz - the array that was be passed to the matrix preallocation routines 844bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 845bd481603SBarry Smith 846bd481603SBarry Smith Level: intermediate 847bd481603SBarry Smith 848bd481603SBarry Smith Notes: 8490598bfebSBarry Smith See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details. 850bd481603SBarry Smith 851bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 852bd481603SBarry Smith 8531620fd73SBarry Smith This is a MACRO not a function because it closes the { started in MatPreallocateInitialize(). 8541620fd73SBarry Smith 855bd481603SBarry Smith Concepts: preallocation^Matrix 856bd481603SBarry Smith 857bd481603SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(), 858eabe889fSLisandro Dalcin MatPreallocateSymmetricSetLocal() 859bd481603SBarry Smith M*/ 860a89bc211SBarry Smith #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree2(dnz,onz);CHKERRQ(_4_ierr);} 8617c922b88SBarry Smith 8627b80b807SBarry Smith /* Routines unique to particular data structures */ 863014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellGetContext(Mat,void *); 864698d4c6aSKris Buschelman 865014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInodeAdjustForInodes(Mat,IS*,IS*); 866014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *); 867698d4c6aSKris Buschelman 868014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetColumnIndices(Mat,PetscInt[]); 869014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetColumnIndices(Mat,PetscInt[]); 870014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 871014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 872014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 873014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJFromTriple(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*,PetscInt,PetscBool); 8747b80b807SBarry Smith 875a96a251dSBarry Smith #define MAT_SKIP_ALLOCATION -4 876a96a251dSBarry Smith 877014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]); 878014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]); 879014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]); 880273d9f13SBarry Smith 881014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 882014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 883014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 884014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []); 885014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 886014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]); 887014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 888014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]); 889014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIDenseSetPreallocation(Mat,PetscScalar[]); 890014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqDenseSetPreallocation(Mat,PetscScalar[]); 8919230625dSJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,const PetscInt*[]); 8929230625dSJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,const PetscInt*[]); 893014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAdjCreateNonemptySubcommMat(Mat,Mat*); 894273d9f13SBarry Smith 895014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqDenseSetLDA(Mat,PetscInt); 896014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDenseGetLocalMatrix(Mat,Mat*); 8971b807ce4Svictorle 898014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStoreValues(Mat); 899014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRetrieveValues(Mat); 9002e8a6d31SBarry Smith 901014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDAADSetCtx(Mat,void*); 9023a7fca6bSBarry Smith 903014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFindNonzeroRows(Mat,IS*); 9047b80b807SBarry Smith /* 9057b80b807SBarry Smith These routines are not usually accessed directly, rather solving is 90694b7f48cSBarry Smith done through the KSP and PC interfaces. 9077b80b807SBarry Smith */ 9087b80b807SBarry Smith 90976bdecfbSBarry Smith /*J 9108f6c3df8SBarry Smith MatOrderingType - String with the name of a PETSc matrix ordering 911d9274352SBarry Smith 912d9274352SBarry Smith Level: beginner 913d9274352SBarry Smith 914d9274352SBarry Smith .seealso: MatGetOrdering() 91576bdecfbSBarry Smith J*/ 91619fd82e9SBarry Smith typedef const char* MatOrderingType; 9172692d6eeSBarry Smith #define MATORDERINGNATURAL "natural" 9182692d6eeSBarry Smith #define MATORDERINGND "nd" 9192692d6eeSBarry Smith #define MATORDERING1WD "1wd" 9202692d6eeSBarry Smith #define MATORDERINGRCM "rcm" 9212692d6eeSBarry Smith #define MATORDERINGQMD "qmd" 9222692d6eeSBarry Smith #define MATORDERINGROWLENGTH "rowlength" 923510184a7SJed Brown #define MATORDERINGWBM "wbm" 924fc1bef75SJed Brown #define MATORDERINGSPECTRAL "spectral" 925312e372aSBarry Smith #define MATORDERINGAMD "amd" /* only works if UMFPACK is installed with PETSc */ 926b12f92e5SBarry Smith 92719fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatGetOrdering(Mat,MatOrderingType,IS*,IS*); 928140e18c1SBarry Smith PETSC_EXTERN PetscErrorCode MatGetOrderingList(PetscFunctionList*); 929bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatOrderingRegister(const char[],PetscErrorCode(*)(Mat,MatOrderingType,IS*,IS*)); 930607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatOrderingRegisterAll(void); 931014dd563SJed Brown PETSC_EXTERN PetscBool MatOrderingRegisterAllCalled; 932140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatOrderingList; 933d4fbbf0eSBarry Smith 934014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS); 935fc1bef75SJed Brown PETSC_EXTERN PetscErrorCode MatCreateLaplacian(Mat,PetscReal,PetscBool,Mat*); 936a2ce50c7SBarry Smith 937d91e6319SBarry Smith /*S 938d90ac83dSHong Zhang MatFactorShiftType - Numeric Shift. 939d90ac83dSHong Zhang 940d90ac83dSHong Zhang Level: beginner 941d90ac83dSHong Zhang 942d90ac83dSHong Zhang S*/ 943d90ac83dSHong Zhang typedef enum {MAT_SHIFT_NONE,MAT_SHIFT_NONZERO,MAT_SHIFT_POSITIVE_DEFINITE,MAT_SHIFT_INBLOCKS} MatFactorShiftType; 9446a6fc655SJed Brown PETSC_EXTERN const char *const MatFactorShiftTypes[]; 9455e9742b9SJed Brown PETSC_EXTERN const char *const MatFactorShiftTypesDetail[]; 946d90ac83dSHong Zhang 947d90ac83dSHong Zhang /*S 9482401956bSBarry Smith MatFactorInfo - Data passed into the matrix factorization routines 949b00f7748SHong Zhang 95061cad860SBarry Smith In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE, that is use 95161cad860SBarry Smith $ MatFactorInfo info(MAT_FACTORINFO_SIZE) 952b00f7748SHong Zhang 95315e8a5b3SHong Zhang Notes: These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC. 954b00f7748SHong Zhang 95561cad860SBarry Smith You can use MatFactorInfoInitialize() to set default values. 95661cad860SBarry Smith 957b00f7748SHong Zhang Level: developer 958b00f7748SHong Zhang 959d7d82daaSBarry Smith .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(), 960d7d82daaSBarry Smith MatFactorInfoInitialize() 961b00f7748SHong Zhang 962b00f7748SHong Zhang S*/ 963b00f7748SHong Zhang typedef struct { 96415e8a5b3SHong Zhang PetscReal diagonal_fill; /* force diagonal to fill in if initially not filled */ 96585317021SBarry Smith PetscReal usedt; 96615e8a5b3SHong Zhang PetscReal dt; /* drop tolerance */ 967b00f7748SHong Zhang PetscReal dtcol; /* tolerance for pivoting */ 96815e8a5b3SHong Zhang PetscReal dtcount; /* maximum nonzeros to be allowed per row */ 96967e28bfeSBarry Smith PetscReal fill; /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */ 970348344bbSBarry Smith PetscReal levels; /* ICC/ILU(levels) */ 971bcd9e38bSBarry Smith PetscReal pivotinblocks; /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0 972bcd9e38bSBarry Smith factorization may be faster if do not pivot */ 97315e8a5b3SHong Zhang PetscReal zeropivot; /* pivot is called zero if less than this */ 974f4db908eSBarry Smith PetscReal shifttype; /* type of shift added to matrix factor to prevent zero pivots */ 975f4db908eSBarry Smith PetscReal shiftamount; /* how large the shift is */ 97615e8a5b3SHong Zhang } MatFactorInfo; 977ffa6d0a5SLois Curfman McInnes 978014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFactorInfoInitialize(MatFactorInfo*); 979014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCholeskyFactor(Mat,IS,const MatFactorInfo*); 980014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCholeskyFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*); 981014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCholeskyFactorNumeric(Mat,Mat,const MatFactorInfo*); 982014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLUFactor(Mat,IS,IS,const MatFactorInfo*); 983014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatILUFactor(Mat,IS,IS,const MatFactorInfo*); 984014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*); 985014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatILUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*); 986014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatICCFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*); 987014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatICCFactor(Mat,IS,const MatFactorInfo*); 988014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLUFactorNumeric(Mat,Mat,const MatFactorInfo*); 989014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*); 990014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolve(Mat,Vec,Vec); 991014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatForwardSolve(Mat,Vec,Vec); 992014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatBackwardSolve(Mat,Vec,Vec); 993014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveAdd(Mat,Vec,Vec,Vec); 994014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveTranspose(Mat,Vec,Vec); 995014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveTransposeAdd(Mat,Vec,Vec,Vec); 996014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolves(Mat,Vecs,Vecs); 9978ed539a5SBarry Smith 998014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetUnfactored(Mat); 999bb5a7306SBarry Smith 1000d91e6319SBarry Smith /*E 1001d91e6319SBarry Smith MatSORType - What type of (S)SOR to perform 1002bb1eb677SSatish Balay 1003d91e6319SBarry Smith Level: beginner 1004d91e6319SBarry Smith 1005d9274352SBarry Smith May be bitwise ORd together 1006d9274352SBarry Smith 1007d91e6319SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 1008d91e6319SBarry Smith 10094e7234bfSBarry Smith MatSORType may be bitwise ORd together, so do not change the numbers 10104e7234bfSBarry Smith 101141f059aeSBarry Smith .seealso: MatSOR() 1012d91e6319SBarry Smith E*/ 1013ee50ffe9SBarry Smith typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3, 1014ee50ffe9SBarry Smith SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8, 1015ee50ffe9SBarry Smith SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16, 101684cb2905SBarry Smith SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType; 1017014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSOR(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec); 10188ed539a5SBarry Smith 1019d4fbbf0eSBarry Smith /* 1020639f9d9dSBarry Smith These routines are for efficiently computing Jacobians via finite differences. 1021639f9d9dSBarry Smith */ 1022b12f92e5SBarry Smith 10237086a01eSPeter Brune /*S 10247086a01eSPeter Brune MatColoring - Object for managing the coloring of matrices. 10257086a01eSPeter Brune 10267086a01eSPeter Brune Level: beginner 10277086a01eSPeter Brune 10287086a01eSPeter Brune Concepts: matrix, coloring 10297086a01eSPeter Brune 10307086a01eSPeter Brune .seealso: MatFDColoringCreate() ISColoring MatFDColoring 10317086a01eSPeter Brune S*/ 10327086a01eSPeter Brune typedef struct _p_MatColoring* MatColoring; 103376bdecfbSBarry Smith /*J 10348f6c3df8SBarry Smith MatColoringType - String with the name of a PETSc matrix coloring 1035d9274352SBarry Smith 1036d9274352SBarry Smith Level: beginner 1037d9274352SBarry Smith 10387086a01eSPeter Brune .seealso: MatColoringSetType(), MatColoring 103976bdecfbSBarry Smith J*/ 10407086a01eSPeter Brune 104119fd82e9SBarry Smith typedef const char* MatColoringType; 10425567d71aSPeter Brune #define MATCOLORINGJP "jp" 10432aa85f04SPeter Brune #define MATCOLORINGMIS "mis" 10442692d6eeSBarry Smith #define MATCOLORINGNATURAL "natural" 10452692d6eeSBarry Smith #define MATCOLORINGSL "sl" 10462692d6eeSBarry Smith #define MATCOLORINGLF "lf" 10472692d6eeSBarry Smith #define MATCOLORINGID "id" 10488121bdceSPeter Brune #define MATCOLORINGGREEDY "greedy" 1049b12f92e5SBarry Smith 10508ac6da0aSPeter Brune /*E 10518ac6da0aSPeter Brune MatColoringWeightType - Type of weight scheme 10528ac6da0aSPeter Brune 10538ac6da0aSPeter Brune Not Collective 10548ac6da0aSPeter Brune 10558ac6da0aSPeter Brune + MAT_COLORING_RANDOM - Random weights 10568ac6da0aSPeter Brune . MAT_COLORING_LEXICAL - Lexical weighting based upon global numbering. 10578ac6da0aSPeter Brune - MAT_COLORING_LF - Last-first weighting. 10588ac6da0aSPeter Brune 10598ac6da0aSPeter Brune Level: intermediate 10608ac6da0aSPeter Brune 10618ac6da0aSPeter Brune Any additions/changes here MUST also be made in include/finclude/petscmat.h 10628ac6da0aSPeter Brune 10638ac6da0aSPeter Brune .seealso: MatCUSPSetFormat(), MatCUSPFormatOperation 10648ac6da0aSPeter Brune E*/ 1065301c140bSPeter Brune typedef enum {MAT_COLORING_WEIGHT_RANDOM,MAT_COLORING_WEIGHT_LEXICAL,MAT_COLORING_WEIGHT_LF,MAT_COLORING_WEIGHT_SL} MatColoringWeightType; 10668ac6da0aSPeter Brune 1067335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringCreate(Mat,MatColoring*); 1068d7f2a307SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringGetDegrees(Mat,PetscInt,PetscInt*); 1069335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringDestroy(MatColoring*); 1070335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringView(MatColoring,PetscViewer); 1071335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetType(MatColoring,MatColoringType); 1072335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetFromOptions(MatColoring); 1073335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetDistance(MatColoring,PetscInt); 1074335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringGetDistance(MatColoring,PetscInt*); 1075335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetMaxColors(MatColoring,PetscInt); 1076335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringGetMaxColors(MatColoring,PetscInt*); 1077335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringApply(MatColoring,ISColoring*); 1078607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatColoringRegisterAll(void); 1079335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringRegister(const char[],PetscErrorCode(*)(MatColoring)); 1080335efc43SPeter Brune PETSC_EXTERN PetscBool MatColoringRegisterAllCalled; 1081014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*); 10828ac6da0aSPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetWeightType(MatColoring,MatColoringWeightType); 1083*c29971fcSPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetWeights(MatColoring,PetscReal*,PetscInt*); 10848ac6da0aSPeter Brune PETSC_EXTERN PetscErrorCode MatColoringCreateWeights(MatColoring,PetscReal **,PetscInt **lperm); 1085639f9d9dSBarry Smith 1086d9274352SBarry Smith /*S 1087d9274352SBarry Smith MatFDColoring - Object for computing a sparse Jacobian via finite differences 1088d9274352SBarry Smith and coloring 1089639f9d9dSBarry Smith 1090d9274352SBarry Smith Level: beginner 1091d9274352SBarry Smith 1092d9274352SBarry Smith Concepts: coloring, sparse Jacobian, finite differences 1093d9274352SBarry Smith 1094d9274352SBarry Smith .seealso: MatFDColoringCreate() 1095d9274352SBarry Smith S*/ 1096e2a1c21fSSatish Balay typedef struct _p_MatFDColoring* MatFDColoring; 1097639f9d9dSBarry Smith 1098014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringCreate(Mat,ISColoring,MatFDColoring *); 1099014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringDestroy(MatFDColoring*); 1100014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringView(MatFDColoring,PetscViewer); 1101014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*); 1102014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**); 1103014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal); 1104014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetFromOptions(MatFDColoring); 1105d1e9a80fSBarry Smith PETSC_EXTERN PetscErrorCode MatFDColoringApply(Mat,MatFDColoring,Vec,void *); 1106014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetF(MatFDColoring,Vec); 1107014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,PetscInt*[]); 1108f86b9fbaSHong Zhang PETSC_EXTERN PetscErrorCode MatFDColoringSetUp(Mat,ISColoring,MatFDColoring); 1109f86b9fbaSHong Zhang PETSC_EXTERN PetscErrorCode MatFDColoringSetBlockSize(MatFDColoring,PetscInt,PetscInt); 1110f86b9fbaSHong Zhang 1111b1683b59SHong Zhang 1112b1683b59SHong Zhang /*S 1113b9af6bddSHong Zhang MatTransposeColoring - Object for computing a sparse matrix product C=A*B^T via coloring 1114b1683b59SHong Zhang 1115b1683b59SHong Zhang Level: beginner 1116b1683b59SHong Zhang 1117b1683b59SHong Zhang Concepts: coloring, sparse matrix product 1118b1683b59SHong Zhang 1119b9af6bddSHong Zhang .seealso: MatTransposeColoringCreate() 1120b1683b59SHong Zhang S*/ 1121b9af6bddSHong Zhang typedef struct _p_MatTransposeColoring* MatTransposeColoring; 1122b1683b59SHong Zhang 1123014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeColoringCreate(Mat,ISColoring,MatTransposeColoring *); 1124014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransColoringApplySpToDen(MatTransposeColoring,Mat,Mat); 1125014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransColoringApplyDenToSp(MatTransposeColoring,Mat,Mat); 1126014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeColoringDestroy(MatTransposeColoring*); 1127b1683b59SHong Zhang 1128639f9d9dSBarry Smith /* 11290752156aSBarry Smith These routines are for partitioning matrices: currently used only 11303eda8832SBarry Smith for adjacency matrix, MatCreateMPIAdj(). 11310752156aSBarry Smith */ 1132ca161407SBarry Smith 1133d9274352SBarry Smith /*S 1134d9274352SBarry Smith MatPartitioning - Object for managing the partitioning of a matrix or graph 1135d9274352SBarry Smith 1136d9274352SBarry Smith Level: beginner 1137d9274352SBarry Smith 1138d9274352SBarry Smith Concepts: partitioning 1139d9274352SBarry Smith 1140743a1026Svictorle .seealso: MatPartitioningCreate(), MatPartitioningType 1141d9274352SBarry Smith S*/ 114291e9ee9fSBarry Smith typedef struct _p_MatPartitioning* MatPartitioning; 1143d9274352SBarry Smith 114476bdecfbSBarry Smith /*J 11458f6c3df8SBarry Smith MatPartitioningType - String with the name of a PETSc matrix partitioning 1146d9274352SBarry Smith 1147d9274352SBarry Smith Level: beginner 11480d04baf8SBarry Smith dm 1149b547a13bSvictorle .seealso: MatPartitioningCreate(), MatPartitioning 115076bdecfbSBarry Smith J*/ 115119fd82e9SBarry Smith typedef const char* MatPartitioningType; 11522692d6eeSBarry Smith #define MATPARTITIONINGCURRENT "current" 11532692d6eeSBarry Smith #define MATPARTITIONINGSQUARE "square" 11542692d6eeSBarry Smith #define MATPARTITIONINGPARMETIS "parmetis" 11552692d6eeSBarry Smith #define MATPARTITIONINGCHACO "chaco" 11562692d6eeSBarry Smith #define MATPARTITIONINGPARTY "party" 115750d91057SBarry Smith #define MATPARTITIONINGPTSCOTCH "ptscotch" 115871306c60Sjroman 1159ca161407SBarry Smith 1160014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningCreate(MPI_Comm,MatPartitioning*); 116119fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningSetType(MatPartitioning,MatPartitioningType); 1162014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetNParts(MatPartitioning,PetscInt); 1163014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetAdjacency(MatPartitioning,Mat); 1164014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]); 1165014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []); 1166014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningApply(MatPartitioning,IS*); 1167014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningDestroy(MatPartitioning*); 11682aabb6bbSBarry Smith 1169bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningRegister(const char[],PetscErrorCode (*)(MatPartitioning)); 117030de9b25SBarry Smith 1171014dd563SJed Brown PETSC_EXTERN PetscBool MatPartitioningRegisterAllCalled; 1172f1144a30SSatish Balay 1173607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningRegisterAll(void); 11742bad1931SBarry Smith 1175014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningView(MatPartitioning,PetscViewer); 1176014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetFromOptions(MatPartitioning); 117719fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningGetType(MatPartitioning,MatPartitioningType*); 1178ca161407SBarry Smith 1179014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetCoarseSequential(MatPartitioning); 1180014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningParmetisGetEdgeCut(MatPartitioning, PetscInt *); 11810752156aSBarry Smith 1182b6956c12SJose Roman typedef enum { MP_CHACO_MULTILEVEL=1,MP_CHACO_SPECTRAL=2,MP_CHACO_LINEAR=4,MP_CHACO_RANDOM=5,MP_CHACO_SCATTERED=6 } MPChacoGlobalType; 11836a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoGlobalTypes[]; 1184b6956c12SJose Roman typedef enum { MP_CHACO_KERNIGHAN=1,MP_CHACO_NONE=2 } MPChacoLocalType; 11856a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoLocalTypes[]; 1186b6956c12SJose Roman typedef enum { MP_CHACO_LANCZOS=0,MP_CHACO_RQI=1 } MPChacoEigenType; 11876a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoEigenTypes[]; 1188b6956c12SJose Roman 1189014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetGlobal(MatPartitioning,MPChacoGlobalType); 1190014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetGlobal(MatPartitioning,MPChacoGlobalType*); 1191014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetLocal(MatPartitioning,MPChacoLocalType); 1192014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetLocal(MatPartitioning,MPChacoLocalType*); 1193014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal); 1194014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType); 1195014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenSolver(MatPartitioning,MPChacoEigenType*); 1196014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenTol(MatPartitioning,PetscReal); 1197014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenTol(MatPartitioning,PetscReal*); 1198014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenNumber(MatPartitioning,PetscInt); 1199014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenNumber(MatPartitioning,PetscInt*); 120071306c60Sjroman 120171306c60Sjroman #define MP_PARTY_OPT "opt" 120271306c60Sjroman #define MP_PARTY_LIN "lin" 120371306c60Sjroman #define MP_PARTY_SCA "sca" 120471306c60Sjroman #define MP_PARTY_RAN "ran" 120571306c60Sjroman #define MP_PARTY_GBF "gbf" 120671306c60Sjroman #define MP_PARTY_GCF "gcf" 120771306c60Sjroman #define MP_PARTY_BUB "bub" 120871306c60Sjroman #define MP_PARTY_DEF "def" 1209014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetGlobal(MatPartitioning,const char*); 121071306c60Sjroman #define MP_PARTY_HELPFUL_SETS "hs" 121171306c60Sjroman #define MP_PARTY_KERNIGHAN_LIN "kl" 121271306c60Sjroman #define MP_PARTY_NONE "no" 1213014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetLocal(MatPartitioning,const char*); 1214014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal); 1215014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetBipart(MatPartitioning,PetscBool); 1216014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscBool); 121771306c60Sjroman 121850d91057SBarry Smith typedef enum { MP_PTSCOTCH_QUALITY,MP_PTSCOTCH_SPEED,MP_PTSCOTCH_BALANCE,MP_PTSCOTCH_SAFETY,MP_PTSCOTCH_SCALABILITY } MPPTScotchStrategyType; 12196a6fc655SJed Brown PETSC_EXTERN const char *const MPPTScotchStrategyTypes[]; 1220e0f1cffaSJose Roman 1221014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetImbalance(MatPartitioning,PetscReal); 1222014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetImbalance(MatPartitioning,PetscReal*); 1223014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetStrategy(MatPartitioning,MPPTScotchStrategyType); 1224014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetStrategy(MatPartitioning,MPPTScotchStrategyType*); 122571306c60Sjroman 1226b43b03e9SMark F. Adams /* 1227b43b03e9SMark F. Adams These routines are for coarsening matrices: 1228b43b03e9SMark F. Adams */ 1229b43b03e9SMark F. Adams 1230b43b03e9SMark F. Adams /*S 1231b43b03e9SMark F. Adams MatCoarsen - Object for managing the coarsening of a graph (symmetric matrix) 1232b43b03e9SMark F. Adams 1233b43b03e9SMark F. Adams Level: beginner 1234b43b03e9SMark F. Adams 1235b43b03e9SMark F. Adams Concepts: coarsen 1236b43b03e9SMark F. Adams 1237b43b03e9SMark F. Adams .seealso: MatCoarsenCreate), MatCoarsenType 1238b43b03e9SMark F. Adams S*/ 1239b43b03e9SMark F. Adams typedef struct _p_MatCoarsen* MatCoarsen; 1240b43b03e9SMark F. Adams 1241b43b03e9SMark F. Adams /*J 12428f6c3df8SBarry Smith MatCoarsenType - String with the name of a PETSc matrix coarsen 1243b43b03e9SMark F. Adams 1244b43b03e9SMark F. Adams Level: beginner 12458f6c3df8SBarry Smith 1246b43b03e9SMark F. Adams .seealso: MatCoarsenCreate(), MatCoarsen 1247b43b03e9SMark F. Adams J*/ 124819fd82e9SBarry Smith typedef const char* MatCoarsenType; 1249b43b03e9SMark F. Adams #define MATCOARSENMIS "mis" 12509057884aSMark F. Adams #define MATCOARSENHEM "hem" 1251b43b03e9SMark F. Adams 12520cbbd2e1SMark F. Adams /* linked list for aggregates */ 125341b27cdeSMark F. Adams typedef struct _PetscCDIntNd{ 125441b27cdeSMark F. Adams struct _PetscCDIntNd *next; 12550cbbd2e1SMark F. Adams PetscInt gid; 125641b27cdeSMark F. Adams }PetscCDIntNd; 12570cbbd2e1SMark F. Adams 12580cbbd2e1SMark F. Adams /* only used by node pool */ 125941b27cdeSMark F. Adams typedef struct _PetscCDArrNd{ 126041b27cdeSMark F. Adams struct _PetscCDArrNd *next; 126141b27cdeSMark F. Adams struct _PetscCDIntNd *array; 126241b27cdeSMark F. Adams }PetscCDArrNd; 12630cbbd2e1SMark F. Adams 12640cbbd2e1SMark F. Adams typedef struct _PetscCoarsenData{ 126582c86c8fSBarry Smith PetscCDArrNd pool_list; /* node pool */ 126641b27cdeSMark F. Adams PetscCDIntNd *new_node; 12670cbbd2e1SMark F. Adams PetscInt new_left; 12680cbbd2e1SMark F. Adams PetscInt chk_sz; 126941b27cdeSMark F. Adams PetscCDIntNd *extra_nodes; 127082c86c8fSBarry Smith PetscCDIntNd **array; /* Array of lists */ 12710cbbd2e1SMark F. Adams PetscInt size; 127282c86c8fSBarry Smith Mat mat; /* cache a Mat for communication data */ 12730cbbd2e1SMark F. Adams }PetscCoarsenData; 12740cbbd2e1SMark F. Adams 1275014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenCreate(MPI_Comm,MatCoarsen*); 127619fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenSetType(MatCoarsen,MatCoarsenType); 1277014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetAdjacency(MatCoarsen,Mat); 1278014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetGreedyOrdering(MatCoarsen,const IS); 1279014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetStrictAggs(MatCoarsen,PetscBool); 1280014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetVerbose(MatCoarsen,PetscInt); 1281014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenGetData( MatCoarsen, PetscCoarsenData ** ); 1282014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenApply(MatCoarsen); 1283014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenDestroy(MatCoarsen*); 1284b43b03e9SMark F. Adams 1285bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenRegister(const char[],PetscErrorCode (*)(MatCoarsen)); 1286b43b03e9SMark F. Adams 1287014dd563SJed Brown PETSC_EXTERN PetscBool MatCoarsenRegisterAllCalled; 1288b43b03e9SMark F. Adams 1289607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenRegisterAll(void); 1290b43b03e9SMark F. Adams 1291014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenView(MatCoarsen,PetscViewer); 1292014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetFromOptions(MatCoarsen); 129319fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenGetType(MatCoarsen,MatCoarsenType*); 1294ce1779c8SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatCoarsenViewFromOptions(MatCoarsen A,const char prefix[],const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,prefix,name);} 1295b43b03e9SMark F. Adams 1296014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMeshToVertexGraph(Mat,PetscInt,Mat*); 1297014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMeshToCellGraph(Mat,PetscInt,Mat*); 1298591294e4SBarry Smith 12990752156aSBarry Smith /* 13000a835dfdSSatish Balay If you add entries here you must also add them to finclude/petscmat.h 1301d4fbbf0eSBarry Smith */ 13021c1c02c0SLois Curfman McInnes typedef enum { MATOP_SET_VALUES=0, 13031c1c02c0SLois Curfman McInnes MATOP_GET_ROW=1, 13041c1c02c0SLois Curfman McInnes MATOP_RESTORE_ROW=2, 13051c1c02c0SLois Curfman McInnes MATOP_MULT=3, 13061c1c02c0SLois Curfman McInnes MATOP_MULT_ADD=4, 13077c922b88SBarry Smith MATOP_MULT_TRANSPOSE=5, 13087c922b88SBarry Smith MATOP_MULT_TRANSPOSE_ADD=6, 13091c1c02c0SLois Curfman McInnes MATOP_SOLVE=7, 13101c1c02c0SLois Curfman McInnes MATOP_SOLVE_ADD=8, 13117c922b88SBarry Smith MATOP_SOLVE_TRANSPOSE=9, 13127c922b88SBarry Smith MATOP_SOLVE_TRANSPOSE_ADD=10, 13131c1c02c0SLois Curfman McInnes MATOP_LUFACTOR=11, 13141c1c02c0SLois Curfman McInnes MATOP_CHOLESKYFACTOR=12, 1315a714c06dSBarry Smith MATOP_SOR=13, 13161c1c02c0SLois Curfman McInnes MATOP_TRANSPOSE=14, 13171c1c02c0SLois Curfman McInnes MATOP_GETINFO=15, 13181c1c02c0SLois Curfman McInnes MATOP_EQUAL=16, 13191c1c02c0SLois Curfman McInnes MATOP_GET_DIAGONAL=17, 13201c1c02c0SLois Curfman McInnes MATOP_DIAGONAL_SCALE=18, 13211c1c02c0SLois Curfman McInnes MATOP_NORM=19, 13221c1c02c0SLois Curfman McInnes MATOP_ASSEMBLY_BEGIN=20, 13231c1c02c0SLois Curfman McInnes MATOP_ASSEMBLY_END=21, 1324d519adbfSMatthew Knepley MATOP_SET_OPTION=22, 1325d519adbfSMatthew Knepley MATOP_ZERO_ENTRIES=23, 1326d519adbfSMatthew Knepley MATOP_ZERO_ROWS=24, 1327d519adbfSMatthew Knepley MATOP_LUFACTOR_SYMBOLIC=25, 1328d519adbfSMatthew Knepley MATOP_LUFACTOR_NUMERIC=26, 1329d519adbfSMatthew Knepley MATOP_CHOLESKY_FACTOR_SYMBOLIC=27, 1330d519adbfSMatthew Knepley MATOP_CHOLESKY_FACTOR_NUMERIC=28, 1331d519adbfSMatthew Knepley MATOP_SETUP_PREALLOCATION=29, 1332d519adbfSMatthew Knepley MATOP_ILUFACTOR_SYMBOLIC=30, 1333d519adbfSMatthew Knepley MATOP_ICCFACTOR_SYMBOLIC=31, 13349d39f69dSJed Brown /* MATOP_PLACEHOLDER_32=32, */ 13359d39f69dSJed Brown /* MATOP_PLACEHOLDER_33=33, */ 1336d519adbfSMatthew Knepley MATOP_DUPLICATE=34, 1337d519adbfSMatthew Knepley MATOP_FORWARD_SOLVE=35, 1338d519adbfSMatthew Knepley MATOP_BACKWARD_SOLVE=36, 1339d519adbfSMatthew Knepley MATOP_ILUFACTOR=37, 1340d519adbfSMatthew Knepley MATOP_ICCFACTOR=38, 1341d519adbfSMatthew Knepley MATOP_AXPY=39, 1342d519adbfSMatthew Knepley MATOP_GET_SUBMATRICES=40, 1343d519adbfSMatthew Knepley MATOP_INCREASE_OVERLAP=41, 1344d519adbfSMatthew Knepley MATOP_GET_VALUES=42, 1345d519adbfSMatthew Knepley MATOP_COPY=43, 1346d519adbfSMatthew Knepley MATOP_GET_ROW_MAX=44, 1347d519adbfSMatthew Knepley MATOP_SCALE=45, 1348d519adbfSMatthew Knepley MATOP_SHIFT=46, 134935153367SBarry Smith MATOP_DIAGONAL_SET=47, 13509d39f69dSJed Brown MATOP_ZERO_ROWS_COLUMNS=48, 13519d39f69dSJed Brown MATOP_SET_RANDOM=49, 1352d519adbfSMatthew Knepley MATOP_GET_ROW_IJ=50, 1353d519adbfSMatthew Knepley MATOP_RESTORE_ROW_IJ=51, 1354d519adbfSMatthew Knepley MATOP_GET_COLUMN_IJ=52, 1355d519adbfSMatthew Knepley MATOP_RESTORE_COLUMN_IJ=53, 1356d519adbfSMatthew Knepley MATOP_FDCOLORING_CREATE=54, 1357d519adbfSMatthew Knepley MATOP_COLORING_PATCH=55, 1358d519adbfSMatthew Knepley MATOP_SET_UNFACTORED=56, 1359d519adbfSMatthew Knepley MATOP_PERMUTE=57, 1360d519adbfSMatthew Knepley MATOP_SET_VALUES_BLOCKED=58, 1361d519adbfSMatthew Knepley MATOP_GET_SUBMATRIX=59, 1362d519adbfSMatthew Knepley MATOP_DESTROY=60, 1363d519adbfSMatthew Knepley MATOP_VIEW=61, 1364d519adbfSMatthew Knepley MATOP_CONVERT_FROM=62, 13657bab7c10SHong Zhang MATOP_MATMAT_MULT=63, 13667bab7c10SHong Zhang MATOP_MATMAT_MULT_SYMBOLIC=64, 13677bab7c10SHong Zhang MATOP_MATMAT_MULT_NUMERIC=65, 1368d519adbfSMatthew Knepley MATOP_SET_LOCAL_TO_GLOBAL_MAP=66, 1369d519adbfSMatthew Knepley MATOP_SET_VALUES_LOCAL=67, 1370d519adbfSMatthew Knepley MATOP_ZERO_ROWS_LOCAL=68, 1371d519adbfSMatthew Knepley MATOP_GET_ROW_MAX_ABS=69, 1372d519adbfSMatthew Knepley MATOP_GET_ROW_MIN_ABS=70, 1373d519adbfSMatthew Knepley MATOP_CONVERT=71, 1374d519adbfSMatthew Knepley MATOP_SET_COLORING=72, 13759d39f69dSJed Brown /* MATOP_PLACEHOLDER_73=73, */ 1376d519adbfSMatthew Knepley MATOP_SET_VALUES_ADIFOR=74, 1377d519adbfSMatthew Knepley MATOP_FD_COLORING_APPLY=75, 1378d519adbfSMatthew Knepley MATOP_SET_FROM_OPTIONS=76, 1379bda6c4cbSBarry Smith MATOP_MULT_CONSTRAINED=77, 1380bda6c4cbSBarry Smith MATOP_MULT_TRANSPOSE_CONSTRAIN=78, 13819d39f69dSJed Brown MATOP_FIND_ZERO_DIAGONALS=79, 1382d519adbfSMatthew Knepley MATOP_MULT_MULTIPLE=80, 1383d519adbfSMatthew Knepley MATOP_SOLVE_MULTIPLE=81, 1384d519adbfSMatthew Knepley MATOP_GET_INERTIA=82, 1385d519adbfSMatthew Knepley MATOP_LOAD=83, 1386d519adbfSMatthew Knepley MATOP_IS_SYMMETRIC=84, 1387d519adbfSMatthew Knepley MATOP_IS_HERMITIAN=85, 1388d519adbfSMatthew Knepley MATOP_IS_STRUCTURALLY_SYMMETRIC=86, 1389820469bcSHong Zhang MATOP_SET_VALUES_BLOCKEDLOCAL=87, 1390d519adbfSMatthew Knepley MATOP_GET_VECS=88, 1391d519adbfSMatthew Knepley MATOP_MAT_MULT=89, 1392d519adbfSMatthew Knepley MATOP_MAT_MULT_SYMBOLIC=90, 1393d519adbfSMatthew Knepley MATOP_MAT_MULT_NUMERIC=91, 1394d519adbfSMatthew Knepley MATOP_PTAP=92, 1395d519adbfSMatthew Knepley MATOP_PTAP_SYMBOLIC=93, 1396d519adbfSMatthew Knepley MATOP_PTAP_NUMERIC=94, 1397bda6c4cbSBarry Smith MATOP_MAT_TRANSPOSE_MULT=95, 1398bda6c4cbSBarry Smith MATOP_MAT_TRANSPOSE_MULT_SYMBO=96, 1399bda6c4cbSBarry Smith MATOP_MAT_TRANSPOSE_MULT_NUMER=97, 14009d39f69dSJed Brown /* MATOP_PLACEHOLDER_98=98, */ 14019d39f69dSJed Brown /* MATOP_PLACEHOLDER_99=99, */ 14029d39f69dSJed Brown /* MATOP_PLACEHOLDER_100=100, */ 14039d39f69dSJed Brown /* MATOP_PLACEHOLDER_101=101, */ 1404d519adbfSMatthew Knepley MATOP_CONJUGATE=102, 14059d39f69dSJed Brown /* MATOP_PLACEHOLDER_103=103, */ 1406d519adbfSMatthew Knepley MATOP_SET_VALUES_ROW=104, 1407d519adbfSMatthew Knepley MATOP_REAL_PART=105, 1408bda6c4cbSBarry Smith MATOP_IMAGINARY_PART=106, 1409bda6c4cbSBarry Smith MATOP_GET_ROW_UPPER_TRIANGULAR=107, 1410bda6c4cbSBarry Smith MATOP_RESTORE_ROW_UPPER_TRIANG=108, 1411bda6c4cbSBarry Smith MATOP_MAT_SOLVE=109, 14120e8d04f7SBarry Smith MATOP_GET_REDUNDANT_MATRIX=110, 1413d519adbfSMatthew Knepley MATOP_GET_ROW_MIN=111, 14140e8d04f7SBarry Smith MATOP_GET_COLUMN_VECTOR=112, 1415d519adbfSMatthew Knepley MATOP_MISSING_DIAGONAL=113, 14160e8d04f7SBarry Smith MATOP_GET_SEQ_NONZERO_STRUCTUR=114, 141789c6957cSBarry Smith MATOP_CREATE=115, 141889c6957cSBarry Smith MATOP_GET_GHOSTS=116, 14190e8d04f7SBarry Smith MATOP_GET_LOCAL_SUB_MATRIX=117, 14200e8d04f7SBarry Smith MATOP_RESTORE_LOCALSUB_MATRIX=118, 1421eeffb40dSHong Zhang MATOP_MULT_DIAGONAL_BLOCK=119, 14220e8d04f7SBarry Smith MATOP_HERMITIAN_TRANSPOSE=120, 14230e8d04f7SBarry Smith MATOP_MULT_HERMITIAN_TRANSPOSE=121, 14240e8d04f7SBarry Smith MATOP_MULT_HERMITIAN_TRANS_ADD=122, 14250e8d04f7SBarry Smith MATOP_GET_MULTI_PROC_BLOCK=123, 14269d39f69dSJed Brown MATOP_FIND_NONZERO_ROWS=124, 14270e8d04f7SBarry Smith MATOP_GET_COLUMN_NORMS=125, 14289d39f69dSJed Brown MATOP_INVERT_BLOCK_DIAGONAL=126, 14299d39f69dSJed Brown /* MATOP_PLACEHOLDER_127=127, */ 14300e8d04f7SBarry Smith MATOP_GET_SUB_MATRICES_PARALLE=128, 14312b8ad9a3SHong Zhang MATOP_SET_VALUES_BATCH=129, 14320e8d04f7SBarry Smith MATOP_TRANSPOSE_MAT_MULT=130, 14330e8d04f7SBarry Smith MATOP_TRANSPOSE_MAT_MULT_SYMBO=131, 14340e8d04f7SBarry Smith MATOP_TRANSPOSE_MAT_MULT_NUMER=132, 14350e8d04f7SBarry Smith MATOP_TRANSPOSE_COLORING_CREAT=133, 14360e8d04f7SBarry Smith MATOP_TRANS_COLORING_APPLY_SPT=134, 14370e8d04f7SBarry Smith MATOP_TRANS_COLORING_APPLY_DEN=135, 14380e8d04f7SBarry Smith MATOP_RART=136, 14390e8d04f7SBarry Smith MATOP_RART_SYMBOLIC=137, 14400e8d04f7SBarry Smith MATOP_RART_NUMERIC=138, 1441e09a3074SHong Zhang MATOP_SET_BLOCK_SIZES=139, 1442f9426fe0SMark Adams MATOP_AYPX=140, 14431919a2e2SJed Brown MATOP_RESIDUAL=141, 14441919a2e2SJed Brown MATOP_FDCOLORING_SETUP= 142 1445fae171e0SBarry Smith } MatOperation; 1446014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatHasOperation(Mat,MatOperation,PetscBool *); 1447014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellSetOperation(Mat,MatOperation,void(*)(void)); 1448014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellGetOperation(Mat,MatOperation,void(**)(void)); 1449014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellSetContext(Mat,void*); 1450112a2221SBarry Smith 145190ace30eSBarry Smith /* 145290ace30eSBarry Smith Codes for matrices stored on disk. By default they are 145390ace30eSBarry Smith stored in a universal format. By changing the format with 14547973ad04SBarry Smith PetscViewerSetFormat(viewer,PETSC_VIEWER_NATIVE); the matrices will 145590ace30eSBarry Smith be stored in a way natural for the matrix, for example dense matrices 145690ace30eSBarry Smith would be stored as dense. Matrices stored this way may only be 1457f4403165SShri Abhyankar read into matrices of the same type. 145890ace30eSBarry Smith */ 145990ace30eSBarry Smith #define MATRIX_BINARY_FORMAT_DENSE -1 146090ace30eSBarry Smith 1461014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetHashTableFactor(Mat,PetscReal); 1462014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatISGetLocalMat(Mat,Mat*); 1463014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatISSetLocalMat(Mat,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 *); 1485014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNullSpace(Mat,MatNullSpace); 1486014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNearNullSpace(Mat,MatNullSpace); 1487014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNearNullSpace(Mat,MatNullSpace*); 1488014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceTest(MatNullSpace,Mat,PetscBool *); 1489014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceView(MatNullSpace,PetscViewer); 1490014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceGetVecs(MatNullSpace,PetscBool*,PetscInt*,const Vec**); 1491014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceCreateRigidBody(Vec,MatNullSpace*); 149274637425SBarry Smith 1493014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatReorderingSeqSBAIJ(Mat,IS); 1494014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetHashTableFactor(Mat,PetscReal); 1495014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetColumnIndices(Mat,PetscInt *); 1496014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJInvertBlockDiagonal(Mat); 14973f1d51d7SBarry Smith 1498014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMAIJ(Mat,PetscInt,Mat*); 1499014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMAIJRedimension(Mat,PetscInt,Mat*); 1500014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMAIJGetAIJ(Mat,Mat*); 1501c4f061fbSSatish Balay 1502014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatComputeExplicitOperator(Mat,Mat*); 1503b0a32e0cSBarry Smith 1504014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalScaleLocal(Mat,Vec); 150504f1ad80SBarry Smith 1506014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*); 1507014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetBase(Mat,Vec,Vec); 1508014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*); 1509014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*)); 1510014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec)); 1511014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDAddNullSpace(Mat,MatNullSpace); 1512014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt); 1513014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDResetHHistory(Mat); 1514014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctionError(Mat,PetscReal); 1515014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetPeriod(Mat,PetscInt); 1516014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDGetH(Mat,PetscScalar *); 1517014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetOptionsPrefix(Mat,const char[]); 1518014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*); 1519014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*); 1520e884886fSBarry Smith 15216370053bSBarry Smith /*S 15226370053bSBarry Smith MatMFFD - A data structured used to manage the computation of the h differencing parameter for matrix-free 15236370053bSBarry Smith Jacobian vector products 1524e884886fSBarry Smith 15256370053bSBarry Smith Notes: MATMFFD is a specific MatType which uses the MatMFFD data structure 15266370053bSBarry Smith 15276370053bSBarry Smith MatMFFD*() methods actually take the Mat as their first argument. Not a MatMFFD data structure 15286370053bSBarry Smith 15296370053bSBarry Smith Level: developer 15306370053bSBarry Smith 15316370053bSBarry Smith .seealso: MATMFFD, MatCreateMFFD(), MatMFFDSetFuction(), MatMFFDSetType(), MatMFFDRegister() 15326370053bSBarry Smith S*/ 1533e884886fSBarry Smith typedef struct _p_MatMFFD* MatMFFD; 1534e884886fSBarry Smith 153576bdecfbSBarry Smith /*J 1536e884886fSBarry Smith MatMFFDType - algorithm used to compute the h used in computing matrix-vector products via differencing of the function 1537e884886fSBarry Smith 1538e884886fSBarry Smith Level: beginner 1539e884886fSBarry Smith 1540e884886fSBarry Smith .seealso: MatMFFDSetType(), MatMFFDRegister() 154176bdecfbSBarry Smith J*/ 154219fd82e9SBarry Smith typedef const char* MatMFFDType; 1543e884886fSBarry Smith #define MATMFFD_DS "ds" 1544e884886fSBarry Smith #define MATMFFD_WP "wp" 1545e884886fSBarry Smith 154619fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDSetType(Mat,MatMFFDType); 1547bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDRegister(const char[],PetscErrorCode (*)(MatMFFD)); 1548e884886fSBarry Smith 1549607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDRegisterAll(void); 1550014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDDSSetUmin(Mat,PetscReal); 1551014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDWPSetComputeNormU(Mat,PetscBool ); 1552e884886fSBarry Smith 1553014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *); 1554014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *); 15557dbadf16SMatthew Knepley 155697969023SHong Zhang /* 155797969023SHong Zhang PETSc interface to MUMPS 155897969023SHong Zhang */ 155997969023SHong Zhang #ifdef PETSC_HAVE_MUMPS 1560014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMumpsSetIcntl(Mat,PetscInt,PetscInt); 1561b9e9ea26SSatish Balay PETSC_EXTERN PetscErrorCode MatMumpsSetCntl(Mat,PetscInt,PetscReal); 156297969023SHong Zhang #endif 156323a5497aSJed Brown 1564d954db57SHong Zhang /* 1565d954db57SHong Zhang PETSc interface to SUPERLU 1566d954db57SHong Zhang */ 1567d954db57SHong Zhang #ifdef PETSC_HAVE_SUPERLU 1568014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSuperluSetILUDropTol(Mat,PetscReal); 1569d954db57SHong Zhang #endif 1570d954db57SHong Zhang 1571aa372e3fSPaul Mullowney #ifdef PETSC_HAVE_CUDA 15723f9c0db1SPaul Mullowney /*E 1573e057df02SPaul Mullowney MatCUSPARSEStorageFormat - indicates the storage format for CUSPARSE (GPU) 15742692e278SPaul Mullowney matrices. 1575e057df02SPaul Mullowney 1576e057df02SPaul Mullowney Not Collective 1577e057df02SPaul Mullowney 15788468deeeSKarl Rupp + MAT_CUSPARSE_CSR - Compressed Sparse Row 15792692e278SPaul Mullowney . MAT_CUSPARSE_ELL - Ellpack (requires CUDA 4.2 or later). 15802692e278SPaul Mullowney - MAT_CUSPARSE_HYB - Hybrid, a combination of Ellpack and Coordinate format (requires CUDA 4.2 or later). 1581e057df02SPaul Mullowney 1582e057df02SPaul Mullowney Level: intermediate 1583e057df02SPaul Mullowney 1584e057df02SPaul Mullowney Any additions/changes here MUST also be made in include/finclude/petscmat.h 1585e057df02SPaul Mullowney 1586e057df02SPaul Mullowney .seealso: MatCUSPARSESetFormat(), MatCUSPARSEFormatOperation 1587e057df02SPaul Mullowney E*/ 1588e057df02SPaul Mullowney 15893f9c0db1SPaul Mullowney typedef enum {MAT_CUSPARSE_CSR, MAT_CUSPARSE_ELL, MAT_CUSPARSE_HYB} MatCUSPARSEStorageFormat; 1590e057df02SPaul Mullowney 1591e057df02SPaul Mullowney /* these will be strings associated with enumerated type defined above */ 1592e057df02SPaul Mullowney PETSC_EXTERN const char *const MatCUSPARSEStorageFormats[]; 1593e057df02SPaul Mullowney 15943f9c0db1SPaul Mullowney /*E 1595e057df02SPaul Mullowney MatCUSPARSEFormatOperation - indicates the operation of CUSPARSE (GPU) 15962692e278SPaul Mullowney matrices whose operation should use a particular storage format. 1597e057df02SPaul Mullowney 1598e057df02SPaul Mullowney Not Collective 1599e057df02SPaul Mullowney 16008468deeeSKarl Rupp + MAT_CUSPARSE_MULT_DIAG - sets the storage format for the diagonal matrix in the parallel MatMult 16018468deeeSKarl Rupp . MAT_CUSPARSE_MULT_OFFDIAG - sets the storage format for the offdiagonal matrix in the parallel MatMult 16028468deeeSKarl Rupp . MAT_CUSPARSE_MULT - sets the storage format for the entire matrix in the serial (single GPU) MatMult 16038468deeeSKarl Rupp - MAT_CUSPARSE_ALL - sets the storage format for all CUSPARSE (GPU) matrices 1604e057df02SPaul Mullowney 1605e057df02SPaul Mullowney Level: intermediate 1606e057df02SPaul Mullowney 1607e057df02SPaul Mullowney .seealso: MatCUSPARSESetFormat(), MatCUSPARSEStorageFormat 1608e057df02SPaul Mullowney E*/ 160936d62e41SPaul Mullowney typedef enum {MAT_CUSPARSE_MULT_DIAG, MAT_CUSPARSE_MULT_OFFDIAG, MAT_CUSPARSE_MULT, MAT_CUSPARSE_ALL} MatCUSPARSEFormatOperation; 1610e057df02SPaul Mullowney 161110e9db80SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 161210e9db80SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCreateAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 1613e057df02SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCUSPARSESetFormat(Mat,MatCUSPARSEFormatOperation,MatCUSPARSEStorageFormat); 16149ae82921SPaul Mullowney #endif 16159ae82921SPaul Mullowney 161690c53902SBarry Smith #if defined(PETSC_HAVE_CUSP) 1617014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCUSP(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 1618014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateAIJCUSP(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 1619e057df02SPaul Mullowney 16203f9c0db1SPaul Mullowney /*E 1621e057df02SPaul Mullowney MatCUSPStorageFormat - indicates the storage format for CUSP (GPU) 162236d62e41SPaul Mullowney matrices. 1623e057df02SPaul Mullowney 1624e057df02SPaul Mullowney Not Collective 1625e057df02SPaul Mullowney 16268468deeeSKarl Rupp + MAT_CUSP_CSR - Compressed Sparse Row 16278468deeeSKarl Rupp . MAT_CUSP_DIA - Diagonal 16288468deeeSKarl Rupp - MAT_CUSP_ELL - Ellpack 1629e057df02SPaul Mullowney 1630e057df02SPaul Mullowney Level: intermediate 1631e057df02SPaul Mullowney 1632e057df02SPaul Mullowney Any additions/changes here MUST also be made in include/finclude/petscmat.h 1633e057df02SPaul Mullowney 1634e057df02SPaul Mullowney .seealso: MatCUSPSetFormat(), MatCUSPFormatOperation 1635e057df02SPaul Mullowney E*/ 16363f9c0db1SPaul Mullowney typedef enum {MAT_CUSP_CSR, MAT_CUSP_DIA, MAT_CUSP_ELL} MatCUSPStorageFormat; 1637e057df02SPaul Mullowney 1638e057df02SPaul Mullowney /* these will be strings associated with enumerated type defined above */ 1639e057df02SPaul Mullowney PETSC_EXTERN const char *const MatCUSPStorageFormats[]; 1640e057df02SPaul Mullowney 16413f9c0db1SPaul Mullowney /*E 1642e057df02SPaul Mullowney MatCUSPFormatOperation - indicates the operation of CUSP (GPU) 164336d62e41SPaul Mullowney matrices whose operation should use a particular storage format. 1644e057df02SPaul Mullowney 1645e057df02SPaul Mullowney Not Collective 1646e057df02SPaul Mullowney 16478468deeeSKarl Rupp + MAT_CUSP_MULT_DIAG - sets the storage format for the diagonal matrix in the parallel MatMult 16488468deeeSKarl Rupp . MAT_CUSP_MULT_OFFDIAG - sets the storage format for the offdiagonal matrix in the parallel MatMult 16498468deeeSKarl Rupp . MAT_CUSP_MULT - sets the storage format for the entire matrix in the serial (single GPU) MatMult 16508468deeeSKarl Rupp - MAT_CUSP_ALL - sets the storage format for all CUSP (GPU) matrices 1651e057df02SPaul Mullowney 1652e057df02SPaul Mullowney Level: intermediate 1653e057df02SPaul Mullowney 1654e057df02SPaul Mullowney Any additions/changes here MUST also be made in include/finclude/petscmat.h 1655e057df02SPaul Mullowney 1656e057df02SPaul Mullowney .seealso: MatCUSPSetFormat(), MatCUSPStorageFormat 1657e057df02SPaul Mullowney E*/ 165836d62e41SPaul Mullowney typedef enum {MAT_CUSP_MULT_DIAG, MAT_CUSP_MULT_OFFDIAG, MAT_CUSP_MULT, MAT_CUSP_ALL} MatCUSPFormatOperation; 1659e057df02SPaul Mullowney 1660e057df02SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCUSPSetFormat(Mat,MatCUSPFormatOperation,MatCUSPStorageFormat); 1661e057df02SPaul Mullowney #endif 166290c53902SBarry Smith 1663d67ff14aSKarl Rupp #if defined(PETSC_HAVE_VIENNACL) 1664d67ff14aSKarl Rupp PETSC_EXTERN PetscErrorCode MatCreateSeqAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 1665d67ff14aSKarl Rupp PETSC_EXTERN PetscErrorCode MatCreateAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 1666d67ff14aSKarl Rupp #endif 1667d67ff14aSKarl Rupp 166854efbe56SHong Zhang /* 166954efbe56SHong Zhang PETSc interface to FFTW 167054efbe56SHong Zhang */ 167154efbe56SHong Zhang #if defined(PETSC_HAVE_FFTW) 1672014dd563SJed Brown PETSC_EXTERN PetscErrorCode VecScatterPetscToFFTW(Mat,Vec,Vec); 1673014dd563SJed Brown PETSC_EXTERN PetscErrorCode VecScatterFFTWToPetsc(Mat,Vec,Vec); 1674014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetVecsFFTW(Mat,Vec*,Vec*,Vec*); 167554efbe56SHong Zhang #endif 167654efbe56SHong Zhang 1677db31f6deSJed Brown #if defined(PETSC_HAVE_ELEMENTAL) 1678607a6623SBarry Smith PETSC_EXTERN PetscErrorCode PetscElementalInitializePackage(void); 1679db31f6deSJed Brown PETSC_EXTERN PetscErrorCode PetscElementalFinalizePackage(void); 1680db31f6deSJed Brown #endif 1681db31f6deSJed Brown 1682014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateNest(MPI_Comm,PetscInt,const IS[],PetscInt,const IS[],const Mat[],Mat*); 1683014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSize(Mat,PetscInt*,PetscInt*); 1684014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetISs(Mat,IS[],IS[]); 1685014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetLocalISs(Mat,IS[],IS[]); 1686014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSubMats(Mat,PetscInt*,PetscInt*,Mat***); 1687014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSubMat(Mat,PetscInt,PetscInt,Mat*); 168819fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatNestSetVecType(Mat,VecType); 1689014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestSetSubMats(Mat,PetscInt,const IS[],PetscInt,const IS[],const Mat[]); 1690014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestSetSubMat(Mat,PetscInt,PetscInt,Mat); 1691d8588912SDave May 16924325cce7SMatthew G Knepley PETSC_EXTERN PetscErrorCode MatChop(Mat,PetscReal); 1693e0a58c10SMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatComputeBandwidth(Mat,PetscReal,PetscInt*); 16944325cce7SMatthew G Knepley 169523a5497aSJed Brown #endif 1696