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*/ 32292d4d306SBarry Smith typedef enum {MAT_OPTION_MIN = -5, 32392d4d306SBarry Smith MAT_NEW_NONZERO_LOCATION_ERR = -4, 32492d4d306SBarry Smith MAT_UNUSED_NONZERO_LOCATION_ERR = -3, 32592d4d306SBarry Smith MAT_NEW_NONZERO_ALLOCATION_ERR = -2, 32692d4d306SBarry Smith MAT_ROW_ORIENTED = -1, 327382164b0SBarry Smith MAT_SYMMETRIC = 1, 328382164b0SBarry Smith MAT_STRUCTURALLY_SYMMETRIC = 2, 329382164b0SBarry Smith MAT_NEW_DIAGONALS = 3, 330382164b0SBarry Smith MAT_IGNORE_OFF_PROC_ENTRIES = 4, 331382164b0SBarry Smith MAT_USE_HASH_TABLE = 5, 332382164b0SBarry Smith MAT_KEEP_NONZERO_PATTERN = 6, 333382164b0SBarry Smith MAT_IGNORE_ZERO_ENTRIES = 7, 334382164b0SBarry Smith MAT_USE_INODES = 8, 335382164b0SBarry Smith MAT_HERMITIAN = 9, 336382164b0SBarry Smith MAT_SYMMETRY_ETERNAL = 10, 33711e456e1SBarry Smith MAT_DUMMY = 11, 338382164b0SBarry Smith MAT_IGNORE_LOWER_TRIANGULAR = 12, 339382164b0SBarry Smith MAT_ERROR_LOWER_TRIANGULAR = 13, 340382164b0SBarry Smith MAT_GETROW_UPPERTRIANGULAR = 14, 341382164b0SBarry Smith MAT_SPD = 15, 34292d4d306SBarry Smith MAT_NO_OFF_PROC_ZERO_ROWS = 16, 34392d4d306SBarry Smith MAT_NO_OFF_PROC_ENTRIES = 17, 34492d4d306SBarry Smith MAT_NEW_NONZERO_LOCATIONS = 18, 34592d4d306SBarry Smith MAT_OPTION_MAX = 19} 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 MatGetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*); 541a93ff8c4SPeter Brune PETSC_EXTERN PetscErrorCode MatGetLayouts(Mat,PetscLayout*,PetscLayout*); 542014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 543014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsLocalIS(Mat,IS,PetscScalar,Vec,Vec); 544014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 545014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocalIS(Mat,IS,PetscScalar,Vec,Vec); 546014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 547014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 548052efed2SBarry Smith 549014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStashSetInitialSize(Mat,PetscInt,PetscInt); 550014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*); 55190f02eecSBarry Smith 552014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInterpolate(Mat,Vec,Vec); 553014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInterpolateAdd(Mat,Vec,Vec,Vec); 554014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestrict(Mat,Vec,Vec); 555014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetVecs(Mat,Vec*,Vec*); 556b2bf6370SHong Zhang PETSC_EXTERN PetscErrorCode MatGetRedundantMatrix(Mat,PetscInt,MPI_Comm,MatReuse,Mat*); 557014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetMultiProcBlock(Mat,MPI_Comm,MatReuse,Mat*); 558014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFindZeroDiagonals(Mat,IS*); 5593a062f41SBarry Smith PETSC_EXTERN PetscErrorCode MatFindOffBlockDiagonalEntries(Mat,IS*); 560bd481603SBarry Smith 561bd481603SBarry Smith /*MC 5621620fd73SBarry Smith MatSetValue - Set a single entry into a matrix. 5631620fd73SBarry Smith 5641620fd73SBarry Smith Not collective 5651620fd73SBarry Smith 5661620fd73SBarry Smith Input Parameters: 5671620fd73SBarry Smith + m - the matrix 5681620fd73SBarry Smith . row - the row location of the entry 5691620fd73SBarry Smith . col - the column location of the entry 5701620fd73SBarry Smith . value - the value to insert 5711620fd73SBarry Smith - mode - either INSERT_VALUES or ADD_VALUES 5721620fd73SBarry Smith 5731620fd73SBarry Smith Notes: 5741620fd73SBarry Smith For efficiency one should use MatSetValues() and set several or many 5751620fd73SBarry Smith values simultaneously if possible. 5761620fd73SBarry Smith 5771620fd73SBarry Smith Level: beginner 5781620fd73SBarry Smith 5791620fd73SBarry Smith .seealso: MatSetValues(), MatSetValueLocal() 5801620fd73SBarry Smith M*/ 5811620fd73SBarry 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);} 5821620fd73SBarry Smith 5831620fd73SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatGetValue(Mat v,PetscInt i,PetscInt j,PetscScalar *va) {return MatGetValues(v,1,&i,1,&j,va);} 5841620fd73SBarry Smith 5851620fd73SBarry 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);} 5861620fd73SBarry Smith 5871620fd73SBarry Smith /*MC 5880d2a7ff2SSatish Balay MatPreallocateInitialize - Begins the block of code that will count the number of nonzeros per 589bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 590bd481603SBarry Smith 591bd481603SBarry Smith Synopsis: 592aaa7dc30SBarry Smith #include <petscmat.h> 593c1ac3661SBarry Smith PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz) 594bd481603SBarry Smith 595bd481603SBarry Smith Collective on MPI_Comm 596bd481603SBarry Smith 597bd481603SBarry Smith Input Parameters: 598bd481603SBarry Smith + comm - the communicator that will share the eventually allocated matrix 599859fcb39SBarry Smith . nrows - the number of LOCAL rows in the matrix 600859fcb39SBarry Smith - ncols - the number of LOCAL columns in the matrix 601bd481603SBarry Smith 602bd481603SBarry Smith Output Parameters: 603bd481603SBarry Smith + dnz - the array that will be passed to the matrix preallocation routines 604bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 605bd481603SBarry Smith 606bd481603SBarry Smith Level: intermediate 607bd481603SBarry Smith 608bd481603SBarry Smith Notes: 609a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 610bd481603SBarry Smith 6111d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 612bd481603SBarry Smith 613bd481603SBarry Smith Use MatPreallocateInitializeSymmetric() for symmetric matrices (MPISBAIJ matrices) 614bd481603SBarry Smith 6151620fd73SBarry Smith This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize(). 6161620fd73SBarry Smith 617bd481603SBarry Smith Concepts: preallocation^Matrix 618bd481603SBarry Smith 619*d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(), 620*d6e23781SBarry Smith MatPreallocateInitializeSymmetric(), MatPreallocateSymmetricSetLocalBlock() 621bd481603SBarry Smith M*/ 622c4f061fbSSatish Balay #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \ 6237c922b88SBarry Smith { \ 624a89bc211SBarry Smith PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__start,__end; \ 6251795a4d1SJed Brown _4_ierr = PetscCalloc2(__nrows,&dnz,__nrows,&onz);CHKERRQ(_4_ierr); \ 6261795a4d1SJed Brown __start = 0; __end = __start; \ 627c1ac3661SBarry Smith _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ctmp;\ 628a89bc211SBarry Smith _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows; 6297c922b88SBarry Smith 630bd481603SBarry Smith /*MC 631bd481603SBarry Smith MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be 632bd481603SBarry Smith inserted using a local number of the rows and columns 633bd481603SBarry Smith 634bd481603SBarry Smith Synopsis: 635aaa7dc30SBarry Smith #include <petscmat.h> 636c1ac3661SBarry Smith PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 637bd481603SBarry Smith 638bd481603SBarry Smith Not Collective 639bd481603SBarry Smith 640bd481603SBarry Smith Input Parameters: 641784ac674SJed Brown + map - the row mapping from local numbering to global numbering 642bd481603SBarry Smith . nrows - the number of rows indicated 6431d73ed98SBarry Smith . rows - the indices of the rows 644784ac674SJed Brown . cmap - the column mapping from local to global numbering 645bd481603SBarry Smith . ncols - the number of columns in the matrix 646bd481603SBarry Smith . cols - the columns indicated 647bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 648bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 649bd481603SBarry Smith 650bd481603SBarry Smith Level: intermediate 651bd481603SBarry Smith 652bd481603SBarry Smith Notes: 653a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 654bd481603SBarry Smith 6551d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 656bd481603SBarry Smith 657bd481603SBarry Smith Concepts: preallocation^Matrix 658bd481603SBarry Smith 659*d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(), 660*d6e23781SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock() 661bd481603SBarry Smith M*/ 662784ac674SJed Brown #define MatPreallocateSetLocal(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \ 663c4f061fbSSatish Balay {\ 664c1ac3661SBarry Smith PetscInt __l;\ 665784ac674SJed Brown _4_ierr = ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\ 666784ac674SJed Brown _4_ierr = ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\ 667c4f061fbSSatish Balay for (__l=0;__l<nrows;__l++) {\ 668ef66eb69SBarry Smith _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 669c4f061fbSSatish Balay }\ 670c4f061fbSSatish Balay } 671c4f061fbSSatish Balay 672bd481603SBarry Smith /*MC 673*d6e23781SBarry Smith MatPreallocateSetLocalBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be 674bd481603SBarry Smith inserted using a local number of the rows and columns 675bd481603SBarry Smith 676bd481603SBarry Smith Synopsis: 677aaa7dc30SBarry Smith #include <petscmat.h> 678*d6e23781SBarry Smith PetscErrorCode MatPreallocateSetLocalBlock(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 679*d6e23781SBarry Smith 680*d6e23781SBarry Smith Not Collective 681*d6e23781SBarry Smith 682*d6e23781SBarry Smith Input Parameters: 683*d6e23781SBarry Smith + map - the row mapping from local numbering to global numbering 684*d6e23781SBarry Smith . nrows - the number of rows indicated 685*d6e23781SBarry Smith . rows - the indices of the rows 686*d6e23781SBarry Smith . cmap - the column mapping from local to global numbering 687*d6e23781SBarry Smith . ncols - the number of columns in the matrix 688*d6e23781SBarry Smith . cols - the columns indicated 689*d6e23781SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 690*d6e23781SBarry Smith - ozn - the other array passed to the matrix preallocation routines 691*d6e23781SBarry Smith 692*d6e23781SBarry Smith Level: intermediate 693*d6e23781SBarry Smith 694*d6e23781SBarry Smith Notes: 695*d6e23781SBarry Smith See Users-Manual: ch_performance for more details. 696*d6e23781SBarry Smith 697*d6e23781SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 698*d6e23781SBarry Smith 699*d6e23781SBarry Smith Concepts: preallocation^Matrix 700*d6e23781SBarry Smith 701*d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(), 702*d6e23781SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock() 703*d6e23781SBarry Smith M*/ 704*d6e23781SBarry Smith #define MatPreallocateSetLocalBlock(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \ 705*d6e23781SBarry Smith {\ 706*d6e23781SBarry Smith PetscInt __l;\ 707*d6e23781SBarry Smith _4_ierr = ISLocalToGlobalMappingApplyBlock(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\ 708*d6e23781SBarry Smith _4_ierr = ISLocalToGlobalMappingApplyBlock(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\ 709*d6e23781SBarry Smith for (__l=0;__l<nrows;__l++) {\ 710*d6e23781SBarry Smith _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 711*d6e23781SBarry Smith }\ 712*d6e23781SBarry Smith } 713*d6e23781SBarry Smith 714*d6e23781SBarry Smith /*MC 715*d6e23781SBarry Smith MatPreallocateSymmetricSetLocalBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be 716*d6e23781SBarry Smith inserted using a local number of the rows and columns 717*d6e23781SBarry Smith 718*d6e23781SBarry Smith Synopsis: 719*d6e23781SBarry Smith #include <petscmat.h> 720*d6e23781SBarry Smith PetscErrorCode MatPreallocateSymmetricSetLocalBlock(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 721bd481603SBarry Smith 722bd481603SBarry Smith Not Collective 723bd481603SBarry Smith 724bd481603SBarry Smith Input Parameters: 725bd481603SBarry Smith + map - the mapping between local numbering and global numbering 726bd481603SBarry Smith . nrows - the number of rows indicated 7271d73ed98SBarry Smith . rows - the indices of the rows 728bd481603SBarry Smith . ncols - the number of columns in the matrix 729bd481603SBarry Smith . cols - the columns indicated 730bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 731bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 732bd481603SBarry Smith 733bd481603SBarry Smith Level: intermediate 734bd481603SBarry Smith 735bd481603SBarry Smith Notes: 736a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 737bd481603SBarry Smith 738bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 739bd481603SBarry Smith 740bd481603SBarry Smith Concepts: preallocation^Matrix 741bd481603SBarry Smith 742*d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateInitialize(), 743*d6e23781SBarry Smith MatPreallocateInitialize(), MatPreallocateSetLocal() 744bd481603SBarry Smith M*/ 745*d6e23781SBarry Smith #define MatPreallocateSymmetricSetLocalBlock(map,nrows,rows,ncols,cols,dnz,onz) 0;\ 746d3d32019SBarry Smith {\ 747c1ac3661SBarry Smith PetscInt __l;\ 748*d6e23781SBarry Smith _4_ierr = ISLocalToGlobalMappingApplyBlock(map,nrows,rows,rows);CHKERRQ(_4_ierr);\ 749*d6e23781SBarry Smith _4_ierr = ISLocalToGlobalMappingApplyBlock(map,ncols,cols,cols);CHKERRQ(_4_ierr);\ 750d3d32019SBarry Smith for (__l=0;__l<nrows;__l++) {\ 751*d6e23781SBarry Smith _4_ierr = MatPreallocateSymmetricSetBlock((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 752d3d32019SBarry Smith }\ 753d3d32019SBarry Smith } 754bd481603SBarry Smith /*MC 755bd481603SBarry Smith MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be 756bd481603SBarry Smith inserted using a local number of the rows and columns 757bd481603SBarry Smith 758bd481603SBarry Smith Synopsis: 759aaa7dc30SBarry Smith #include <petscmat.h> 760c1ac3661SBarry Smith PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 761bd481603SBarry Smith 762bd481603SBarry Smith Not Collective 763bd481603SBarry Smith 764bd481603SBarry Smith Input Parameters: 76564075487SBarry Smith + row - the row 766bd481603SBarry Smith . ncols - the number of columns in the matrix 767a50f8bf6SHong Zhang - cols - the columns indicated 768a50f8bf6SHong Zhang 769a50f8bf6SHong Zhang Output Parameters: 770a50f8bf6SHong Zhang + dnz - the array that will be passed to the matrix preallocation routines 771bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 772bd481603SBarry Smith 773bd481603SBarry Smith Level: intermediate 774bd481603SBarry Smith 775bd481603SBarry Smith Notes: 776a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 777bd481603SBarry Smith 778bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 779bd481603SBarry Smith 7801620fd73SBarry Smith This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize(). 7811620fd73SBarry Smith 782bd481603SBarry Smith Concepts: preallocation^Matrix 783bd481603SBarry Smith 784*d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(), 785*d6e23781SBarry Smith MatPreallocateInitialize(), MatPreallocateSetLocal() 786bd481603SBarry Smith M*/ 787c4f061fbSSatish Balay #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\ 788c1ac3661SBarry Smith { PetscInt __i; \ 789e32f2f54SBarry 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);\ 790e32f2f54SBarry 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);\ 7917c922b88SBarry Smith for (__i=0; __i<nc; __i++) {\ 79264075487SBarry Smith if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \ 7937cc688d7SBarry Smith else dnz[row - __rstart]++;\ 7947c922b88SBarry Smith }\ 7957c922b88SBarry Smith } 7967c922b88SBarry Smith 797bd481603SBarry Smith /*MC 798*d6e23781SBarry Smith MatPreallocateSymmetricSetBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be 799bd481603SBarry Smith inserted using a local number of the rows and columns 800bd481603SBarry Smith 801bd481603SBarry Smith Synopsis: 802aaa7dc30SBarry Smith #include <petscmat.h> 803*d6e23781SBarry Smith PetscErrorCode MatPreallocateSymmetricSetBlock(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 804bd481603SBarry Smith 805bd481603SBarry Smith Not Collective 806bd481603SBarry Smith 807bd481603SBarry Smith Input Parameters: 808bd481603SBarry Smith + nrows - the number of rows indicated 8091d73ed98SBarry Smith . rows - the indices of the rows 810bd481603SBarry Smith . ncols - the number of columns in the matrix 811bd481603SBarry Smith . cols - the columns indicated 812bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 813bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 814bd481603SBarry Smith 815bd481603SBarry Smith Level: intermediate 816bd481603SBarry Smith 817bd481603SBarry Smith Notes: 818a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 819bd481603SBarry Smith 820bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 821bd481603SBarry Smith 8221620fd73SBarry Smith This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize(). 8231620fd73SBarry Smith 824bd481603SBarry Smith Concepts: preallocation^Matrix 825bd481603SBarry Smith 826*d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateInitialize(), 827*d6e23781SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocal() 828bd481603SBarry Smith M*/ 829*d6e23781SBarry Smith #define MatPreallocateSymmetricSetBlock(row,nc,cols,dnz,onz) 0;\ 830c1ac3661SBarry Smith { PetscInt __i; \ 831d3d32019SBarry Smith for (__i=0; __i<nc; __i++) {\ 832d3d32019SBarry Smith if (cols[__i] >= __end) onz[row - __rstart]++; \ 833d3d32019SBarry Smith else if (cols[__i] >= row) dnz[row - __rstart]++;\ 834d3d32019SBarry Smith }\ 835d3d32019SBarry Smith } 836d3d32019SBarry Smith 837bd481603SBarry Smith /*MC 83816371a99SBarry Smith MatPreallocateLocation - An alternative to MatPreallocationSet() that puts the nonzero locations into the matrix if it exists 83916371a99SBarry Smith 84016371a99SBarry Smith Synopsis: 841aaa7dc30SBarry Smith #include <petscmat.h> 84216371a99SBarry Smith PetscErrorCode MatPreallocateLocations(Mat A,PetscInt row,PetscInt ncols,PetscInt *cols,PetscInt *dnz,PetscInt *onz) 84316371a99SBarry Smith 84416371a99SBarry Smith Not Collective 84516371a99SBarry Smith 84616371a99SBarry Smith Input Parameters: 84716371a99SBarry Smith . A - matrix 84816371a99SBarry Smith . row - row where values exist (must be local to this process) 84916371a99SBarry Smith . ncols - number of columns 85016371a99SBarry Smith . cols - columns with nonzeros 85116371a99SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 85216371a99SBarry Smith - ozn - the other array passed to the matrix preallocation routines 85316371a99SBarry Smith 85416371a99SBarry Smith Level: intermediate 85516371a99SBarry Smith 85616371a99SBarry Smith Notes: 857a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 85816371a99SBarry Smith 85916371a99SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 86016371a99SBarry Smith 86116371a99SBarry Smith This is a MACRO not a function because it uses a bunch of variables private to the MatPreallocation.... routines. 86216371a99SBarry Smith 86316371a99SBarry Smith Concepts: preallocation^Matrix 86416371a99SBarry Smith 865*d6e23781SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(), 866*d6e23781SBarry Smith MatPreallocateSymmetricSetLocalBlock() 86716371a99SBarry Smith M*/ 8680298fd71SBarry 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);} 86916371a99SBarry Smith 87016371a99SBarry Smith 87116371a99SBarry Smith /*MC 8720d2a7ff2SSatish Balay MatPreallocateFinalize - Ends the block of code that will count the number of nonzeros per 873bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 874bd481603SBarry Smith 875bd481603SBarry Smith Synopsis: 876aaa7dc30SBarry Smith #include <petscmat.h> 877c1ac3661SBarry Smith PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz) 878bd481603SBarry Smith 879bd481603SBarry Smith Collective on MPI_Comm 880bd481603SBarry Smith 881bd481603SBarry Smith Input Parameters: 88216371a99SBarry Smith + dnz - the array that was be passed to the matrix preallocation routines 883bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 884bd481603SBarry Smith 885bd481603SBarry Smith Level: intermediate 886bd481603SBarry Smith 887bd481603SBarry Smith Notes: 888a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 889bd481603SBarry Smith 890bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 891bd481603SBarry Smith 8921620fd73SBarry Smith This is a MACRO not a function because it closes the { started in MatPreallocateInitialize(). 8931620fd73SBarry Smith 894bd481603SBarry Smith Concepts: preallocation^Matrix 895bd481603SBarry Smith 896*d6e23781SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(), 897*d6e23781SBarry Smith MatPreallocateSymmetricSetLocalBlock() 898bd481603SBarry Smith M*/ 899a89bc211SBarry Smith #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree2(dnz,onz);CHKERRQ(_4_ierr);} 9007c922b88SBarry Smith 9017b80b807SBarry Smith /* Routines unique to particular data structures */ 902014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellGetContext(Mat,void *); 903698d4c6aSKris Buschelman 904014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInodeAdjustForInodes(Mat,IS*,IS*); 905014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *); 906698d4c6aSKris Buschelman 907014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetColumnIndices(Mat,PetscInt[]); 908014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetColumnIndices(Mat,PetscInt[]); 909014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 910014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 911014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 912014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJFromTriple(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*,PetscInt,PetscBool); 9137b80b807SBarry Smith 914a96a251dSBarry Smith #define MAT_SKIP_ALLOCATION -4 915a96a251dSBarry Smith 916014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]); 917014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]); 918014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]); 919273d9f13SBarry Smith 920014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 921014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 922014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 923014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []); 924014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 925014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]); 926014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 927014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]); 928014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIDenseSetPreallocation(Mat,PetscScalar[]); 929014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqDenseSetPreallocation(Mat,PetscScalar[]); 9309230625dSJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,const PetscInt*[]); 9319230625dSJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,const PetscInt*[]); 932014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAdjCreateNonemptySubcommMat(Mat,Mat*); 933273d9f13SBarry Smith 934014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqDenseSetLDA(Mat,PetscInt); 935014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDenseGetLocalMatrix(Mat,Mat*); 9361b807ce4Svictorle 937014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStoreValues(Mat); 938014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRetrieveValues(Mat); 9392e8a6d31SBarry Smith 940014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDAADSetCtx(Mat,void*); 9413a7fca6bSBarry Smith 942014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFindNonzeroRows(Mat,IS*); 9437b80b807SBarry Smith /* 9447b80b807SBarry Smith These routines are not usually accessed directly, rather solving is 94594b7f48cSBarry Smith done through the KSP and PC interfaces. 9467b80b807SBarry Smith */ 9477b80b807SBarry Smith 94876bdecfbSBarry Smith /*J 9498f6c3df8SBarry Smith MatOrderingType - String with the name of a PETSc matrix ordering 950d9274352SBarry Smith 951d9274352SBarry Smith Level: beginner 952d9274352SBarry Smith 953d9274352SBarry Smith .seealso: MatGetOrdering() 95476bdecfbSBarry Smith J*/ 95519fd82e9SBarry Smith typedef const char* MatOrderingType; 9562692d6eeSBarry Smith #define MATORDERINGNATURAL "natural" 9572692d6eeSBarry Smith #define MATORDERINGND "nd" 9582692d6eeSBarry Smith #define MATORDERING1WD "1wd" 9592692d6eeSBarry Smith #define MATORDERINGRCM "rcm" 9602692d6eeSBarry Smith #define MATORDERINGQMD "qmd" 9612692d6eeSBarry Smith #define MATORDERINGROWLENGTH "rowlength" 962510184a7SJed Brown #define MATORDERINGWBM "wbm" 963fc1bef75SJed Brown #define MATORDERINGSPECTRAL "spectral" 964312e372aSBarry Smith #define MATORDERINGAMD "amd" /* only works if UMFPACK is installed with PETSc */ 965b12f92e5SBarry Smith 96619fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatGetOrdering(Mat,MatOrderingType,IS*,IS*); 967140e18c1SBarry Smith PETSC_EXTERN PetscErrorCode MatGetOrderingList(PetscFunctionList*); 968bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatOrderingRegister(const char[],PetscErrorCode(*)(Mat,MatOrderingType,IS*,IS*)); 969607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatOrderingRegisterAll(void); 970014dd563SJed Brown PETSC_EXTERN PetscBool MatOrderingRegisterAllCalled; 971140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatOrderingList; 972d4fbbf0eSBarry Smith 973014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS); 974fc1bef75SJed Brown PETSC_EXTERN PetscErrorCode MatCreateLaplacian(Mat,PetscReal,PetscBool,Mat*); 975a2ce50c7SBarry Smith 976d91e6319SBarry Smith /*S 977d90ac83dSHong Zhang MatFactorShiftType - Numeric Shift. 978d90ac83dSHong Zhang 979d90ac83dSHong Zhang Level: beginner 980d90ac83dSHong Zhang 981d90ac83dSHong Zhang S*/ 982d90ac83dSHong Zhang typedef enum {MAT_SHIFT_NONE,MAT_SHIFT_NONZERO,MAT_SHIFT_POSITIVE_DEFINITE,MAT_SHIFT_INBLOCKS} MatFactorShiftType; 9836a6fc655SJed Brown PETSC_EXTERN const char *const MatFactorShiftTypes[]; 9845e9742b9SJed Brown PETSC_EXTERN const char *const MatFactorShiftTypesDetail[]; 985d90ac83dSHong Zhang 986d90ac83dSHong Zhang /*S 9872401956bSBarry Smith MatFactorInfo - Data passed into the matrix factorization routines 988b00f7748SHong Zhang 98961cad860SBarry Smith In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE, that is use 99061cad860SBarry Smith $ MatFactorInfo info(MAT_FACTORINFO_SIZE) 991b00f7748SHong Zhang 99215e8a5b3SHong Zhang Notes: These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC. 993b00f7748SHong Zhang 99461cad860SBarry Smith You can use MatFactorInfoInitialize() to set default values. 99561cad860SBarry Smith 996b00f7748SHong Zhang Level: developer 997b00f7748SHong Zhang 998d7d82daaSBarry Smith .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(), 999d7d82daaSBarry Smith MatFactorInfoInitialize() 1000b00f7748SHong Zhang 1001b00f7748SHong Zhang S*/ 1002b00f7748SHong Zhang typedef struct { 100315e8a5b3SHong Zhang PetscReal diagonal_fill; /* force diagonal to fill in if initially not filled */ 100485317021SBarry Smith PetscReal usedt; 100515e8a5b3SHong Zhang PetscReal dt; /* drop tolerance */ 1006b00f7748SHong Zhang PetscReal dtcol; /* tolerance for pivoting */ 100715e8a5b3SHong Zhang PetscReal dtcount; /* maximum nonzeros to be allowed per row */ 100867e28bfeSBarry Smith PetscReal fill; /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */ 1009348344bbSBarry Smith PetscReal levels; /* ICC/ILU(levels) */ 1010bcd9e38bSBarry Smith PetscReal pivotinblocks; /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0 1011bcd9e38bSBarry Smith factorization may be faster if do not pivot */ 101215e8a5b3SHong Zhang PetscReal zeropivot; /* pivot is called zero if less than this */ 1013f4db908eSBarry Smith PetscReal shifttype; /* type of shift added to matrix factor to prevent zero pivots */ 1014f4db908eSBarry Smith PetscReal shiftamount; /* how large the shift is */ 101515e8a5b3SHong Zhang } MatFactorInfo; 1016ffa6d0a5SLois Curfman McInnes 1017014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFactorInfoInitialize(MatFactorInfo*); 1018014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCholeskyFactor(Mat,IS,const MatFactorInfo*); 1019014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCholeskyFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*); 1020014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCholeskyFactorNumeric(Mat,Mat,const MatFactorInfo*); 1021014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLUFactor(Mat,IS,IS,const MatFactorInfo*); 1022014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatILUFactor(Mat,IS,IS,const MatFactorInfo*); 1023014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*); 1024014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatILUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*); 1025014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatICCFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*); 1026014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatICCFactor(Mat,IS,const MatFactorInfo*); 1027014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLUFactorNumeric(Mat,Mat,const MatFactorInfo*); 1028014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*); 1029014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolve(Mat,Vec,Vec); 1030014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatForwardSolve(Mat,Vec,Vec); 1031014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatBackwardSolve(Mat,Vec,Vec); 1032014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveAdd(Mat,Vec,Vec,Vec); 1033014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveTranspose(Mat,Vec,Vec); 1034014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveTransposeAdd(Mat,Vec,Vec,Vec); 1035014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolves(Mat,Vecs,Vecs); 10368ed539a5SBarry Smith 1037014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetUnfactored(Mat); 1038bb5a7306SBarry Smith 1039d91e6319SBarry Smith /*E 1040d91e6319SBarry Smith MatSORType - What type of (S)SOR to perform 1041bb1eb677SSatish Balay 1042d91e6319SBarry Smith Level: beginner 1043d91e6319SBarry Smith 1044d9274352SBarry Smith May be bitwise ORd together 1045d9274352SBarry Smith 1046d91e6319SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 1047d91e6319SBarry Smith 10484e7234bfSBarry Smith MatSORType may be bitwise ORd together, so do not change the numbers 10494e7234bfSBarry Smith 105041f059aeSBarry Smith .seealso: MatSOR() 1051d91e6319SBarry Smith E*/ 1052ee50ffe9SBarry Smith typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3, 1053ee50ffe9SBarry Smith SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8, 1054ee50ffe9SBarry Smith SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16, 105584cb2905SBarry Smith SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType; 1056014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSOR(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec); 10578ed539a5SBarry Smith 1058d4fbbf0eSBarry Smith /* 1059639f9d9dSBarry Smith These routines are for efficiently computing Jacobians via finite differences. 1060639f9d9dSBarry Smith */ 1061b12f92e5SBarry Smith 10627086a01eSPeter Brune /*S 10637086a01eSPeter Brune MatColoring - Object for managing the coloring of matrices. 10647086a01eSPeter Brune 10657086a01eSPeter Brune Level: beginner 10667086a01eSPeter Brune 10677086a01eSPeter Brune Concepts: matrix, coloring 10687086a01eSPeter Brune 10697086a01eSPeter Brune .seealso: MatFDColoringCreate() ISColoring MatFDColoring 10707086a01eSPeter Brune S*/ 10717086a01eSPeter Brune typedef struct _p_MatColoring* MatColoring; 107276bdecfbSBarry Smith /*J 10738f6c3df8SBarry Smith MatColoringType - String with the name of a PETSc matrix coloring 1074d9274352SBarry Smith 1075d9274352SBarry Smith Level: beginner 1076d9274352SBarry Smith 10777086a01eSPeter Brune .seealso: MatColoringSetType(), MatColoring 107876bdecfbSBarry Smith J*/ 10797086a01eSPeter Brune 108019fd82e9SBarry Smith typedef const char* MatColoringType; 10815567d71aSPeter Brune #define MATCOLORINGJP "jp" 10822692d6eeSBarry Smith #define MATCOLORINGNATURAL "natural" 10832692d6eeSBarry Smith #define MATCOLORINGSL "sl" 10842692d6eeSBarry Smith #define MATCOLORINGLF "lf" 10852692d6eeSBarry Smith #define MATCOLORINGID "id" 10868121bdceSPeter Brune #define MATCOLORINGGREEDY "greedy" 1087b12f92e5SBarry Smith 10888ac6da0aSPeter Brune /*E 10898ac6da0aSPeter Brune MatColoringWeightType - Type of weight scheme 10908ac6da0aSPeter Brune 10918ac6da0aSPeter Brune Not Collective 10928ac6da0aSPeter Brune 10938ac6da0aSPeter Brune + MAT_COLORING_RANDOM - Random weights 10948ac6da0aSPeter Brune . MAT_COLORING_LEXICAL - Lexical weighting based upon global numbering. 10958ac6da0aSPeter Brune - MAT_COLORING_LF - Last-first weighting. 10968ac6da0aSPeter Brune 10978ac6da0aSPeter Brune Level: intermediate 10988ac6da0aSPeter Brune 10998ac6da0aSPeter Brune Any additions/changes here MUST also be made in include/finclude/petscmat.h 11008ac6da0aSPeter Brune 11018ac6da0aSPeter Brune .seealso: MatCUSPSetFormat(), MatCUSPFormatOperation 11028ac6da0aSPeter Brune E*/ 1103301c140bSPeter Brune typedef enum {MAT_COLORING_WEIGHT_RANDOM,MAT_COLORING_WEIGHT_LEXICAL,MAT_COLORING_WEIGHT_LF,MAT_COLORING_WEIGHT_SL} MatColoringWeightType; 11048ac6da0aSPeter Brune 1105335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringCreate(Mat,MatColoring*); 1106d7f2a307SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringGetDegrees(Mat,PetscInt,PetscInt*); 1107335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringDestroy(MatColoring*); 1108335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringView(MatColoring,PetscViewer); 1109335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetType(MatColoring,MatColoringType); 1110335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetFromOptions(MatColoring); 1111335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetDistance(MatColoring,PetscInt); 1112335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringGetDistance(MatColoring,PetscInt*); 1113335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetMaxColors(MatColoring,PetscInt); 1114335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringGetMaxColors(MatColoring,PetscInt*); 1115335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringApply(MatColoring,ISColoring*); 1116607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatColoringRegisterAll(void); 1117335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringRegister(const char[],PetscErrorCode(*)(MatColoring)); 1118335efc43SPeter Brune PETSC_EXTERN PetscBool MatColoringRegisterAllCalled; 1119014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*); 11208ac6da0aSPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetWeightType(MatColoring,MatColoringWeightType); 1121c29971fcSPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetWeights(MatColoring,PetscReal*,PetscInt*); 11228ac6da0aSPeter Brune PETSC_EXTERN PetscErrorCode MatColoringCreateWeights(MatColoring,PetscReal **,PetscInt **lperm); 1123639f9d9dSBarry Smith 1124d9274352SBarry Smith /*S 1125d9274352SBarry Smith MatFDColoring - Object for computing a sparse Jacobian via finite differences 1126d9274352SBarry Smith and coloring 1127639f9d9dSBarry Smith 1128d9274352SBarry Smith Level: beginner 1129d9274352SBarry Smith 1130d9274352SBarry Smith Concepts: coloring, sparse Jacobian, finite differences 1131d9274352SBarry Smith 1132d9274352SBarry Smith .seealso: MatFDColoringCreate() 1133d9274352SBarry Smith S*/ 1134e2a1c21fSSatish Balay typedef struct _p_MatFDColoring* MatFDColoring; 1135639f9d9dSBarry Smith 1136014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringCreate(Mat,ISColoring,MatFDColoring *); 1137014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringDestroy(MatFDColoring*); 1138014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringView(MatFDColoring,PetscViewer); 1139014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*); 1140014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**); 1141014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal); 1142014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetFromOptions(MatFDColoring); 1143d1e9a80fSBarry Smith PETSC_EXTERN PetscErrorCode MatFDColoringApply(Mat,MatFDColoring,Vec,void *); 1144014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetF(MatFDColoring,Vec); 1145014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,PetscInt*[]); 1146f86b9fbaSHong Zhang PETSC_EXTERN PetscErrorCode MatFDColoringSetUp(Mat,ISColoring,MatFDColoring); 1147f86b9fbaSHong Zhang PETSC_EXTERN PetscErrorCode MatFDColoringSetBlockSize(MatFDColoring,PetscInt,PetscInt); 1148f86b9fbaSHong Zhang 1149b1683b59SHong Zhang 1150b1683b59SHong Zhang /*S 1151b9af6bddSHong Zhang MatTransposeColoring - Object for computing a sparse matrix product C=A*B^T via coloring 1152b1683b59SHong Zhang 1153b1683b59SHong Zhang Level: beginner 1154b1683b59SHong Zhang 1155b1683b59SHong Zhang Concepts: coloring, sparse matrix product 1156b1683b59SHong Zhang 1157b9af6bddSHong Zhang .seealso: MatTransposeColoringCreate() 1158b1683b59SHong Zhang S*/ 1159b9af6bddSHong Zhang typedef struct _p_MatTransposeColoring* MatTransposeColoring; 1160b1683b59SHong Zhang 1161014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeColoringCreate(Mat,ISColoring,MatTransposeColoring *); 1162014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransColoringApplySpToDen(MatTransposeColoring,Mat,Mat); 1163014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransColoringApplyDenToSp(MatTransposeColoring,Mat,Mat); 1164014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeColoringDestroy(MatTransposeColoring*); 1165b1683b59SHong Zhang 1166639f9d9dSBarry Smith /* 11670752156aSBarry Smith These routines are for partitioning matrices: currently used only 11683eda8832SBarry Smith for adjacency matrix, MatCreateMPIAdj(). 11690752156aSBarry Smith */ 1170ca161407SBarry Smith 1171d9274352SBarry Smith /*S 1172d9274352SBarry Smith MatPartitioning - Object for managing the partitioning of a matrix or graph 1173d9274352SBarry Smith 1174d9274352SBarry Smith Level: beginner 1175d9274352SBarry Smith 1176d9274352SBarry Smith Concepts: partitioning 1177d9274352SBarry Smith 1178743a1026Svictorle .seealso: MatPartitioningCreate(), MatPartitioningType 1179d9274352SBarry Smith S*/ 118091e9ee9fSBarry Smith typedef struct _p_MatPartitioning* MatPartitioning; 1181d9274352SBarry Smith 118276bdecfbSBarry Smith /*J 11838f6c3df8SBarry Smith MatPartitioningType - String with the name of a PETSc matrix partitioning 1184d9274352SBarry Smith 1185d9274352SBarry Smith Level: beginner 11860d04baf8SBarry Smith dm 1187b547a13bSvictorle .seealso: MatPartitioningCreate(), MatPartitioning 118876bdecfbSBarry Smith J*/ 118919fd82e9SBarry Smith typedef const char* MatPartitioningType; 11902692d6eeSBarry Smith #define MATPARTITIONINGCURRENT "current" 11912692d6eeSBarry Smith #define MATPARTITIONINGSQUARE "square" 11922692d6eeSBarry Smith #define MATPARTITIONINGPARMETIS "parmetis" 11932692d6eeSBarry Smith #define MATPARTITIONINGCHACO "chaco" 11942692d6eeSBarry Smith #define MATPARTITIONINGPARTY "party" 119550d91057SBarry Smith #define MATPARTITIONINGPTSCOTCH "ptscotch" 119671306c60Sjroman 1197ca161407SBarry Smith 1198014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningCreate(MPI_Comm,MatPartitioning*); 119919fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningSetType(MatPartitioning,MatPartitioningType); 1200014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetNParts(MatPartitioning,PetscInt); 1201014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetAdjacency(MatPartitioning,Mat); 1202014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]); 1203014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []); 1204014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningApply(MatPartitioning,IS*); 1205014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningDestroy(MatPartitioning*); 12062aabb6bbSBarry Smith 1207bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningRegister(const char[],PetscErrorCode (*)(MatPartitioning)); 120830de9b25SBarry Smith 1209014dd563SJed Brown PETSC_EXTERN PetscBool MatPartitioningRegisterAllCalled; 1210f1144a30SSatish Balay 1211607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningRegisterAll(void); 12122bad1931SBarry Smith 1213014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningView(MatPartitioning,PetscViewer); 1214014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetFromOptions(MatPartitioning); 121519fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningGetType(MatPartitioning,MatPartitioningType*); 1216ca161407SBarry Smith 1217014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetCoarseSequential(MatPartitioning); 1218014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningParmetisGetEdgeCut(MatPartitioning, PetscInt *); 12190752156aSBarry Smith 1220b6956c12SJose Roman typedef enum { MP_CHACO_MULTILEVEL=1,MP_CHACO_SPECTRAL=2,MP_CHACO_LINEAR=4,MP_CHACO_RANDOM=5,MP_CHACO_SCATTERED=6 } MPChacoGlobalType; 12216a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoGlobalTypes[]; 1222b6956c12SJose Roman typedef enum { MP_CHACO_KERNIGHAN=1,MP_CHACO_NONE=2 } MPChacoLocalType; 12236a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoLocalTypes[]; 1224b6956c12SJose Roman typedef enum { MP_CHACO_LANCZOS=0,MP_CHACO_RQI=1 } MPChacoEigenType; 12256a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoEigenTypes[]; 1226b6956c12SJose Roman 1227014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetGlobal(MatPartitioning,MPChacoGlobalType); 1228014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetGlobal(MatPartitioning,MPChacoGlobalType*); 1229014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetLocal(MatPartitioning,MPChacoLocalType); 1230014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetLocal(MatPartitioning,MPChacoLocalType*); 1231014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal); 1232014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType); 1233014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenSolver(MatPartitioning,MPChacoEigenType*); 1234014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenTol(MatPartitioning,PetscReal); 1235014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenTol(MatPartitioning,PetscReal*); 1236014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenNumber(MatPartitioning,PetscInt); 1237014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenNumber(MatPartitioning,PetscInt*); 123871306c60Sjroman 123971306c60Sjroman #define MP_PARTY_OPT "opt" 124071306c60Sjroman #define MP_PARTY_LIN "lin" 124171306c60Sjroman #define MP_PARTY_SCA "sca" 124271306c60Sjroman #define MP_PARTY_RAN "ran" 124371306c60Sjroman #define MP_PARTY_GBF "gbf" 124471306c60Sjroman #define MP_PARTY_GCF "gcf" 124571306c60Sjroman #define MP_PARTY_BUB "bub" 124671306c60Sjroman #define MP_PARTY_DEF "def" 1247014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetGlobal(MatPartitioning,const char*); 124871306c60Sjroman #define MP_PARTY_HELPFUL_SETS "hs" 124971306c60Sjroman #define MP_PARTY_KERNIGHAN_LIN "kl" 125071306c60Sjroman #define MP_PARTY_NONE "no" 1251014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetLocal(MatPartitioning,const char*); 1252014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal); 1253014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetBipart(MatPartitioning,PetscBool); 1254014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscBool); 125571306c60Sjroman 125650d91057SBarry Smith typedef enum { MP_PTSCOTCH_QUALITY,MP_PTSCOTCH_SPEED,MP_PTSCOTCH_BALANCE,MP_PTSCOTCH_SAFETY,MP_PTSCOTCH_SCALABILITY } MPPTScotchStrategyType; 12576a6fc655SJed Brown PETSC_EXTERN const char *const MPPTScotchStrategyTypes[]; 1258e0f1cffaSJose Roman 1259014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetImbalance(MatPartitioning,PetscReal); 1260014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetImbalance(MatPartitioning,PetscReal*); 1261014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetStrategy(MatPartitioning,MPPTScotchStrategyType); 1262014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetStrategy(MatPartitioning,MPPTScotchStrategyType*); 126371306c60Sjroman 1264b43b03e9SMark F. Adams /* 1265b43b03e9SMark F. Adams These routines are for coarsening matrices: 1266b43b03e9SMark F. Adams */ 1267b43b03e9SMark F. Adams 1268b43b03e9SMark F. Adams /*S 1269b43b03e9SMark F. Adams MatCoarsen - Object for managing the coarsening of a graph (symmetric matrix) 1270b43b03e9SMark F. Adams 1271b43b03e9SMark F. Adams Level: beginner 1272b43b03e9SMark F. Adams 1273b43b03e9SMark F. Adams Concepts: coarsen 1274b43b03e9SMark F. Adams 1275b43b03e9SMark F. Adams .seealso: MatCoarsenCreate), MatCoarsenType 1276b43b03e9SMark F. Adams S*/ 1277b43b03e9SMark F. Adams typedef struct _p_MatCoarsen* MatCoarsen; 1278b43b03e9SMark F. Adams 1279b43b03e9SMark F. Adams /*J 12808f6c3df8SBarry Smith MatCoarsenType - String with the name of a PETSc matrix coarsen 1281b43b03e9SMark F. Adams 1282b43b03e9SMark F. Adams Level: beginner 12838f6c3df8SBarry Smith 1284b43b03e9SMark F. Adams .seealso: MatCoarsenCreate(), MatCoarsen 1285b43b03e9SMark F. Adams J*/ 128619fd82e9SBarry Smith typedef const char* MatCoarsenType; 1287b43b03e9SMark F. Adams #define MATCOARSENMIS "mis" 12889057884aSMark F. Adams #define MATCOARSENHEM "hem" 1289b43b03e9SMark F. Adams 12900cbbd2e1SMark F. Adams /* linked list for aggregates */ 129141b27cdeSMark F. Adams typedef struct _PetscCDIntNd{ 129241b27cdeSMark F. Adams struct _PetscCDIntNd *next; 12930cbbd2e1SMark F. Adams PetscInt gid; 129441b27cdeSMark F. Adams }PetscCDIntNd; 12950cbbd2e1SMark F. Adams 12960cbbd2e1SMark F. Adams /* only used by node pool */ 129741b27cdeSMark F. Adams typedef struct _PetscCDArrNd{ 129841b27cdeSMark F. Adams struct _PetscCDArrNd *next; 129941b27cdeSMark F. Adams struct _PetscCDIntNd *array; 130041b27cdeSMark F. Adams }PetscCDArrNd; 13010cbbd2e1SMark F. Adams 13020cbbd2e1SMark F. Adams typedef struct _PetscCoarsenData{ 130382c86c8fSBarry Smith PetscCDArrNd pool_list; /* node pool */ 130441b27cdeSMark F. Adams PetscCDIntNd *new_node; 13050cbbd2e1SMark F. Adams PetscInt new_left; 13060cbbd2e1SMark F. Adams PetscInt chk_sz; 130741b27cdeSMark F. Adams PetscCDIntNd *extra_nodes; 130882c86c8fSBarry Smith PetscCDIntNd **array; /* Array of lists */ 13090cbbd2e1SMark F. Adams PetscInt size; 131082c86c8fSBarry Smith Mat mat; /* cache a Mat for communication data */ 13110cbbd2e1SMark F. Adams }PetscCoarsenData; 13120cbbd2e1SMark F. Adams 1313014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenCreate(MPI_Comm,MatCoarsen*); 131419fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenSetType(MatCoarsen,MatCoarsenType); 1315014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetAdjacency(MatCoarsen,Mat); 1316014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetGreedyOrdering(MatCoarsen,const IS); 1317014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetStrictAggs(MatCoarsen,PetscBool); 1318014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetVerbose(MatCoarsen,PetscInt); 1319014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenGetData( MatCoarsen, PetscCoarsenData ** ); 1320014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenApply(MatCoarsen); 1321014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenDestroy(MatCoarsen*); 1322b43b03e9SMark F. Adams 1323bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenRegister(const char[],PetscErrorCode (*)(MatCoarsen)); 1324b43b03e9SMark F. Adams 1325014dd563SJed Brown PETSC_EXTERN PetscBool MatCoarsenRegisterAllCalled; 1326b43b03e9SMark F. Adams 1327607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenRegisterAll(void); 1328b43b03e9SMark F. Adams 1329014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenView(MatCoarsen,PetscViewer); 1330014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetFromOptions(MatCoarsen); 133119fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenGetType(MatCoarsen,MatCoarsenType*); 1332ce1779c8SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatCoarsenViewFromOptions(MatCoarsen A,const char prefix[],const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,prefix,name);} 1333b43b03e9SMark F. Adams 1334014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMeshToVertexGraph(Mat,PetscInt,Mat*); 1335014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMeshToCellGraph(Mat,PetscInt,Mat*); 1336591294e4SBarry Smith 13370752156aSBarry Smith /* 13380a835dfdSSatish Balay If you add entries here you must also add them to finclude/petscmat.h 1339d4fbbf0eSBarry Smith */ 13401c1c02c0SLois Curfman McInnes typedef enum { MATOP_SET_VALUES=0, 13411c1c02c0SLois Curfman McInnes MATOP_GET_ROW=1, 13421c1c02c0SLois Curfman McInnes MATOP_RESTORE_ROW=2, 13431c1c02c0SLois Curfman McInnes MATOP_MULT=3, 13441c1c02c0SLois Curfman McInnes MATOP_MULT_ADD=4, 13457c922b88SBarry Smith MATOP_MULT_TRANSPOSE=5, 13467c922b88SBarry Smith MATOP_MULT_TRANSPOSE_ADD=6, 13471c1c02c0SLois Curfman McInnes MATOP_SOLVE=7, 13481c1c02c0SLois Curfman McInnes MATOP_SOLVE_ADD=8, 13497c922b88SBarry Smith MATOP_SOLVE_TRANSPOSE=9, 13507c922b88SBarry Smith MATOP_SOLVE_TRANSPOSE_ADD=10, 13511c1c02c0SLois Curfman McInnes MATOP_LUFACTOR=11, 13521c1c02c0SLois Curfman McInnes MATOP_CHOLESKYFACTOR=12, 1353a714c06dSBarry Smith MATOP_SOR=13, 13541c1c02c0SLois Curfman McInnes MATOP_TRANSPOSE=14, 13551c1c02c0SLois Curfman McInnes MATOP_GETINFO=15, 13561c1c02c0SLois Curfman McInnes MATOP_EQUAL=16, 13571c1c02c0SLois Curfman McInnes MATOP_GET_DIAGONAL=17, 13581c1c02c0SLois Curfman McInnes MATOP_DIAGONAL_SCALE=18, 13591c1c02c0SLois Curfman McInnes MATOP_NORM=19, 13601c1c02c0SLois Curfman McInnes MATOP_ASSEMBLY_BEGIN=20, 13611c1c02c0SLois Curfman McInnes MATOP_ASSEMBLY_END=21, 1362d519adbfSMatthew Knepley MATOP_SET_OPTION=22, 1363d519adbfSMatthew Knepley MATOP_ZERO_ENTRIES=23, 1364d519adbfSMatthew Knepley MATOP_ZERO_ROWS=24, 1365d519adbfSMatthew Knepley MATOP_LUFACTOR_SYMBOLIC=25, 1366d519adbfSMatthew Knepley MATOP_LUFACTOR_NUMERIC=26, 1367d519adbfSMatthew Knepley MATOP_CHOLESKY_FACTOR_SYMBOLIC=27, 1368d519adbfSMatthew Knepley MATOP_CHOLESKY_FACTOR_NUMERIC=28, 1369d519adbfSMatthew Knepley MATOP_SETUP_PREALLOCATION=29, 1370d519adbfSMatthew Knepley MATOP_ILUFACTOR_SYMBOLIC=30, 1371d519adbfSMatthew Knepley MATOP_ICCFACTOR_SYMBOLIC=31, 13729d39f69dSJed Brown /* MATOP_PLACEHOLDER_32=32, */ 13739d39f69dSJed Brown /* MATOP_PLACEHOLDER_33=33, */ 1374d519adbfSMatthew Knepley MATOP_DUPLICATE=34, 1375d519adbfSMatthew Knepley MATOP_FORWARD_SOLVE=35, 1376d519adbfSMatthew Knepley MATOP_BACKWARD_SOLVE=36, 1377d519adbfSMatthew Knepley MATOP_ILUFACTOR=37, 1378d519adbfSMatthew Knepley MATOP_ICCFACTOR=38, 1379d519adbfSMatthew Knepley MATOP_AXPY=39, 1380d519adbfSMatthew Knepley MATOP_GET_SUBMATRICES=40, 1381d519adbfSMatthew Knepley MATOP_INCREASE_OVERLAP=41, 1382d519adbfSMatthew Knepley MATOP_GET_VALUES=42, 1383d519adbfSMatthew Knepley MATOP_COPY=43, 1384d519adbfSMatthew Knepley MATOP_GET_ROW_MAX=44, 1385d519adbfSMatthew Knepley MATOP_SCALE=45, 1386d519adbfSMatthew Knepley MATOP_SHIFT=46, 138735153367SBarry Smith MATOP_DIAGONAL_SET=47, 13889d39f69dSJed Brown MATOP_ZERO_ROWS_COLUMNS=48, 13899d39f69dSJed Brown MATOP_SET_RANDOM=49, 1390d519adbfSMatthew Knepley MATOP_GET_ROW_IJ=50, 1391d519adbfSMatthew Knepley MATOP_RESTORE_ROW_IJ=51, 1392d519adbfSMatthew Knepley MATOP_GET_COLUMN_IJ=52, 1393d519adbfSMatthew Knepley MATOP_RESTORE_COLUMN_IJ=53, 1394d519adbfSMatthew Knepley MATOP_FDCOLORING_CREATE=54, 1395d519adbfSMatthew Knepley MATOP_COLORING_PATCH=55, 1396d519adbfSMatthew Knepley MATOP_SET_UNFACTORED=56, 1397d519adbfSMatthew Knepley MATOP_PERMUTE=57, 1398d519adbfSMatthew Knepley MATOP_SET_VALUES_BLOCKED=58, 1399d519adbfSMatthew Knepley MATOP_GET_SUBMATRIX=59, 1400d519adbfSMatthew Knepley MATOP_DESTROY=60, 1401d519adbfSMatthew Knepley MATOP_VIEW=61, 1402d519adbfSMatthew Knepley MATOP_CONVERT_FROM=62, 14037bab7c10SHong Zhang MATOP_MATMAT_MULT=63, 14047bab7c10SHong Zhang MATOP_MATMAT_MULT_SYMBOLIC=64, 14057bab7c10SHong Zhang MATOP_MATMAT_MULT_NUMERIC=65, 1406d519adbfSMatthew Knepley MATOP_SET_LOCAL_TO_GLOBAL_MAP=66, 1407d519adbfSMatthew Knepley MATOP_SET_VALUES_LOCAL=67, 1408d519adbfSMatthew Knepley MATOP_ZERO_ROWS_LOCAL=68, 1409d519adbfSMatthew Knepley MATOP_GET_ROW_MAX_ABS=69, 1410d519adbfSMatthew Knepley MATOP_GET_ROW_MIN_ABS=70, 1411d519adbfSMatthew Knepley MATOP_CONVERT=71, 1412d519adbfSMatthew Knepley MATOP_SET_COLORING=72, 14139d39f69dSJed Brown /* MATOP_PLACEHOLDER_73=73, */ 1414d519adbfSMatthew Knepley MATOP_SET_VALUES_ADIFOR=74, 1415d519adbfSMatthew Knepley MATOP_FD_COLORING_APPLY=75, 1416d519adbfSMatthew Knepley MATOP_SET_FROM_OPTIONS=76, 1417bda6c4cbSBarry Smith MATOP_MULT_CONSTRAINED=77, 1418bda6c4cbSBarry Smith MATOP_MULT_TRANSPOSE_CONSTRAIN=78, 14199d39f69dSJed Brown MATOP_FIND_ZERO_DIAGONALS=79, 1420d519adbfSMatthew Knepley MATOP_MULT_MULTIPLE=80, 1421d519adbfSMatthew Knepley MATOP_SOLVE_MULTIPLE=81, 1422d519adbfSMatthew Knepley MATOP_GET_INERTIA=82, 1423d519adbfSMatthew Knepley MATOP_LOAD=83, 1424d519adbfSMatthew Knepley MATOP_IS_SYMMETRIC=84, 1425d519adbfSMatthew Knepley MATOP_IS_HERMITIAN=85, 1426d519adbfSMatthew Knepley MATOP_IS_STRUCTURALLY_SYMMETRIC=86, 1427820469bcSHong Zhang MATOP_SET_VALUES_BLOCKEDLOCAL=87, 1428d519adbfSMatthew Knepley MATOP_GET_VECS=88, 1429d519adbfSMatthew Knepley MATOP_MAT_MULT=89, 1430d519adbfSMatthew Knepley MATOP_MAT_MULT_SYMBOLIC=90, 1431d519adbfSMatthew Knepley MATOP_MAT_MULT_NUMERIC=91, 1432d519adbfSMatthew Knepley MATOP_PTAP=92, 1433d519adbfSMatthew Knepley MATOP_PTAP_SYMBOLIC=93, 1434d519adbfSMatthew Knepley MATOP_PTAP_NUMERIC=94, 1435bda6c4cbSBarry Smith MATOP_MAT_TRANSPOSE_MULT=95, 1436bda6c4cbSBarry Smith MATOP_MAT_TRANSPOSE_MULT_SYMBO=96, 1437bda6c4cbSBarry Smith MATOP_MAT_TRANSPOSE_MULT_NUMER=97, 14389d39f69dSJed Brown /* MATOP_PLACEHOLDER_98=98, */ 14399d39f69dSJed Brown /* MATOP_PLACEHOLDER_99=99, */ 14409d39f69dSJed Brown /* MATOP_PLACEHOLDER_100=100, */ 14419d39f69dSJed Brown /* MATOP_PLACEHOLDER_101=101, */ 1442d519adbfSMatthew Knepley MATOP_CONJUGATE=102, 14439d39f69dSJed Brown /* MATOP_PLACEHOLDER_103=103, */ 1444d519adbfSMatthew Knepley MATOP_SET_VALUES_ROW=104, 1445d519adbfSMatthew Knepley MATOP_REAL_PART=105, 1446bda6c4cbSBarry Smith MATOP_IMAGINARY_PART=106, 1447bda6c4cbSBarry Smith MATOP_GET_ROW_UPPER_TRIANGULAR=107, 1448bda6c4cbSBarry Smith MATOP_RESTORE_ROW_UPPER_TRIANG=108, 1449bda6c4cbSBarry Smith MATOP_MAT_SOLVE=109, 14500e8d04f7SBarry Smith MATOP_GET_REDUNDANT_MATRIX=110, 1451d519adbfSMatthew Knepley MATOP_GET_ROW_MIN=111, 14520e8d04f7SBarry Smith MATOP_GET_COLUMN_VECTOR=112, 1453d519adbfSMatthew Knepley MATOP_MISSING_DIAGONAL=113, 14540e8d04f7SBarry Smith MATOP_GET_SEQ_NONZERO_STRUCTUR=114, 145589c6957cSBarry Smith MATOP_CREATE=115, 145689c6957cSBarry Smith MATOP_GET_GHOSTS=116, 14570e8d04f7SBarry Smith MATOP_GET_LOCAL_SUB_MATRIX=117, 14580e8d04f7SBarry Smith MATOP_RESTORE_LOCALSUB_MATRIX=118, 1459eeffb40dSHong Zhang MATOP_MULT_DIAGONAL_BLOCK=119, 14600e8d04f7SBarry Smith MATOP_HERMITIAN_TRANSPOSE=120, 14610e8d04f7SBarry Smith MATOP_MULT_HERMITIAN_TRANSPOSE=121, 14620e8d04f7SBarry Smith MATOP_MULT_HERMITIAN_TRANS_ADD=122, 14630e8d04f7SBarry Smith MATOP_GET_MULTI_PROC_BLOCK=123, 14649d39f69dSJed Brown MATOP_FIND_NONZERO_ROWS=124, 14650e8d04f7SBarry Smith MATOP_GET_COLUMN_NORMS=125, 14669d39f69dSJed Brown MATOP_INVERT_BLOCK_DIAGONAL=126, 14679d39f69dSJed Brown /* MATOP_PLACEHOLDER_127=127, */ 14680e8d04f7SBarry Smith MATOP_GET_SUB_MATRICES_PARALLE=128, 14692b8ad9a3SHong Zhang MATOP_SET_VALUES_BATCH=129, 14700e8d04f7SBarry Smith MATOP_TRANSPOSE_MAT_MULT=130, 14710e8d04f7SBarry Smith MATOP_TRANSPOSE_MAT_MULT_SYMBO=131, 14720e8d04f7SBarry Smith MATOP_TRANSPOSE_MAT_MULT_NUMER=132, 14730e8d04f7SBarry Smith MATOP_TRANSPOSE_COLORING_CREAT=133, 14740e8d04f7SBarry Smith MATOP_TRANS_COLORING_APPLY_SPT=134, 14750e8d04f7SBarry Smith MATOP_TRANS_COLORING_APPLY_DEN=135, 14760e8d04f7SBarry Smith MATOP_RART=136, 14770e8d04f7SBarry Smith MATOP_RART_SYMBOLIC=137, 14780e8d04f7SBarry Smith MATOP_RART_NUMERIC=138, 1479e09a3074SHong Zhang MATOP_SET_BLOCK_SIZES=139, 1480f9426fe0SMark Adams MATOP_AYPX=140, 14811919a2e2SJed Brown MATOP_RESIDUAL=141, 14821919a2e2SJed Brown MATOP_FDCOLORING_SETUP= 142 1483fae171e0SBarry Smith } MatOperation; 1484014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatHasOperation(Mat,MatOperation,PetscBool *); 1485014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellSetOperation(Mat,MatOperation,void(*)(void)); 1486014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellGetOperation(Mat,MatOperation,void(**)(void)); 1487014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellSetContext(Mat,void*); 1488112a2221SBarry Smith 148990ace30eSBarry Smith /* 149090ace30eSBarry Smith Codes for matrices stored on disk. By default they are 149190ace30eSBarry Smith stored in a universal format. By changing the format with 14927973ad04SBarry Smith PetscViewerSetFormat(viewer,PETSC_VIEWER_NATIVE); the matrices will 149390ace30eSBarry Smith be stored in a way natural for the matrix, for example dense matrices 149490ace30eSBarry Smith would be stored as dense. Matrices stored this way may only be 1495f4403165SShri Abhyankar read into matrices of the same type. 149690ace30eSBarry Smith */ 149790ace30eSBarry Smith #define MATRIX_BINARY_FORMAT_DENSE -1 149890ace30eSBarry Smith 1499014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetHashTableFactor(Mat,PetscReal); 1500014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatISGetLocalMat(Mat,Mat*); 1501014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatISSetLocalMat(Mat,Mat); 15021f4e1ec7SBarry Smith 1503d9274352SBarry Smith /*S 1504d9274352SBarry Smith MatNullSpace - Object that removes a null space from a vector, i.e. 1505d9274352SBarry Smith orthogonalizes the vector to a subsapce 1506d9274352SBarry Smith 1507f7a9e4ceSBarry Smith Level: advanced 1508d9274352SBarry Smith 1509d9274352SBarry Smith Concepts: matrix; linear operator, null space 1510d9274352SBarry Smith 15116e1639daSBarry Smith Users manual sections: 15126e1639daSBarry Smith . sec_singular 15136e1639daSBarry Smith 1514d9274352SBarry Smith .seealso: MatNullSpaceCreate() 1515d9274352SBarry Smith S*/ 151674637425SBarry Smith typedef struct _p_MatNullSpace* MatNullSpace; 1517d9274352SBarry Smith 1518014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceCreate(MPI_Comm,PetscBool ,PetscInt,const Vec[],MatNullSpace*); 1519014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(MatNullSpace,Vec,void*),void*); 1520014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceDestroy(MatNullSpace*); 1521d0195637SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceRemove(MatNullSpace,Vec); 1522014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNullSpace(Mat, MatNullSpace *); 1523014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNullSpace(Mat,MatNullSpace); 1524014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNearNullSpace(Mat,MatNullSpace); 1525014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNearNullSpace(Mat,MatNullSpace*); 1526014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceTest(MatNullSpace,Mat,PetscBool *); 1527014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceView(MatNullSpace,PetscViewer); 1528014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceGetVecs(MatNullSpace,PetscBool*,PetscInt*,const Vec**); 1529014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceCreateRigidBody(Vec,MatNullSpace*); 153074637425SBarry Smith 1531014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatReorderingSeqSBAIJ(Mat,IS); 1532014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetHashTableFactor(Mat,PetscReal); 1533014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetColumnIndices(Mat,PetscInt *); 1534014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJInvertBlockDiagonal(Mat); 15353f1d51d7SBarry Smith 1536014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMAIJ(Mat,PetscInt,Mat*); 1537014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMAIJRedimension(Mat,PetscInt,Mat*); 1538014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMAIJGetAIJ(Mat,Mat*); 1539c4f061fbSSatish Balay 1540014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatComputeExplicitOperator(Mat,Mat*); 1541b0a32e0cSBarry Smith 1542014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalScaleLocal(Mat,Vec); 154304f1ad80SBarry Smith 1544014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*); 1545014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetBase(Mat,Vec,Vec); 1546014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*); 1547014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*)); 1548014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec)); 1549014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDAddNullSpace(Mat,MatNullSpace); 1550014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt); 1551014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDResetHHistory(Mat); 1552014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctionError(Mat,PetscReal); 1553014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetPeriod(Mat,PetscInt); 1554014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDGetH(Mat,PetscScalar *); 1555014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetOptionsPrefix(Mat,const char[]); 1556014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*); 1557014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*); 1558e884886fSBarry Smith 15596370053bSBarry Smith /*S 15606370053bSBarry Smith MatMFFD - A data structured used to manage the computation of the h differencing parameter for matrix-free 15616370053bSBarry Smith Jacobian vector products 1562e884886fSBarry Smith 15636370053bSBarry Smith Notes: MATMFFD is a specific MatType which uses the MatMFFD data structure 15646370053bSBarry Smith 15656370053bSBarry Smith MatMFFD*() methods actually take the Mat as their first argument. Not a MatMFFD data structure 15666370053bSBarry Smith 15676370053bSBarry Smith Level: developer 15686370053bSBarry Smith 15696370053bSBarry Smith .seealso: MATMFFD, MatCreateMFFD(), MatMFFDSetFuction(), MatMFFDSetType(), MatMFFDRegister() 15706370053bSBarry Smith S*/ 1571e884886fSBarry Smith typedef struct _p_MatMFFD* MatMFFD; 1572e884886fSBarry Smith 157376bdecfbSBarry Smith /*J 1574e884886fSBarry Smith MatMFFDType - algorithm used to compute the h used in computing matrix-vector products via differencing of the function 1575e884886fSBarry Smith 1576e884886fSBarry Smith Level: beginner 1577e884886fSBarry Smith 1578e884886fSBarry Smith .seealso: MatMFFDSetType(), MatMFFDRegister() 157976bdecfbSBarry Smith J*/ 158019fd82e9SBarry Smith typedef const char* MatMFFDType; 1581e884886fSBarry Smith #define MATMFFD_DS "ds" 1582e884886fSBarry Smith #define MATMFFD_WP "wp" 1583e884886fSBarry Smith 158419fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDSetType(Mat,MatMFFDType); 1585bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDRegister(const char[],PetscErrorCode (*)(MatMFFD)); 1586e884886fSBarry Smith 1587607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDRegisterAll(void); 1588014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDDSSetUmin(Mat,PetscReal); 1589014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDWPSetComputeNormU(Mat,PetscBool ); 1590e884886fSBarry Smith 1591014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *); 1592014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *); 15937dbadf16SMatthew Knepley 159497969023SHong Zhang /* 159597969023SHong Zhang PETSc interface to MUMPS 159697969023SHong Zhang */ 159797969023SHong Zhang #ifdef PETSC_HAVE_MUMPS 1598014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMumpsSetIcntl(Mat,PetscInt,PetscInt); 1599b9e9ea26SSatish Balay PETSC_EXTERN PetscErrorCode MatMumpsSetCntl(Mat,PetscInt,PetscReal); 160097969023SHong Zhang #endif 160123a5497aSJed Brown 1602d954db57SHong Zhang /* 1603d954db57SHong Zhang PETSc interface to SUPERLU 1604d954db57SHong Zhang */ 1605d954db57SHong Zhang #ifdef PETSC_HAVE_SUPERLU 1606014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSuperluSetILUDropTol(Mat,PetscReal); 1607d954db57SHong Zhang #endif 1608d954db57SHong Zhang 1609aa372e3fSPaul Mullowney #ifdef PETSC_HAVE_CUDA 16103f9c0db1SPaul Mullowney /*E 1611e057df02SPaul Mullowney MatCUSPARSEStorageFormat - indicates the storage format for CUSPARSE (GPU) 16122692e278SPaul Mullowney matrices. 1613e057df02SPaul Mullowney 1614e057df02SPaul Mullowney Not Collective 1615e057df02SPaul Mullowney 16168468deeeSKarl Rupp + MAT_CUSPARSE_CSR - Compressed Sparse Row 16172692e278SPaul Mullowney . MAT_CUSPARSE_ELL - Ellpack (requires CUDA 4.2 or later). 16182692e278SPaul Mullowney - MAT_CUSPARSE_HYB - Hybrid, a combination of Ellpack and Coordinate format (requires CUDA 4.2 or later). 1619e057df02SPaul Mullowney 1620e057df02SPaul Mullowney Level: intermediate 1621e057df02SPaul Mullowney 1622e057df02SPaul Mullowney Any additions/changes here MUST also be made in include/finclude/petscmat.h 1623e057df02SPaul Mullowney 1624e057df02SPaul Mullowney .seealso: MatCUSPARSESetFormat(), MatCUSPARSEFormatOperation 1625e057df02SPaul Mullowney E*/ 1626e057df02SPaul Mullowney 16273f9c0db1SPaul Mullowney typedef enum {MAT_CUSPARSE_CSR, MAT_CUSPARSE_ELL, MAT_CUSPARSE_HYB} MatCUSPARSEStorageFormat; 1628e057df02SPaul Mullowney 1629e057df02SPaul Mullowney /* these will be strings associated with enumerated type defined above */ 1630e057df02SPaul Mullowney PETSC_EXTERN const char *const MatCUSPARSEStorageFormats[]; 1631e057df02SPaul Mullowney 16323f9c0db1SPaul Mullowney /*E 1633e057df02SPaul Mullowney MatCUSPARSEFormatOperation - indicates the operation of CUSPARSE (GPU) 16342692e278SPaul Mullowney matrices whose operation should use a particular storage format. 1635e057df02SPaul Mullowney 1636e057df02SPaul Mullowney Not Collective 1637e057df02SPaul Mullowney 16388468deeeSKarl Rupp + MAT_CUSPARSE_MULT_DIAG - sets the storage format for the diagonal matrix in the parallel MatMult 16398468deeeSKarl Rupp . MAT_CUSPARSE_MULT_OFFDIAG - sets the storage format for the offdiagonal matrix in the parallel MatMult 16408468deeeSKarl Rupp . MAT_CUSPARSE_MULT - sets the storage format for the entire matrix in the serial (single GPU) MatMult 16418468deeeSKarl Rupp - MAT_CUSPARSE_ALL - sets the storage format for all CUSPARSE (GPU) matrices 1642e057df02SPaul Mullowney 1643e057df02SPaul Mullowney Level: intermediate 1644e057df02SPaul Mullowney 1645e057df02SPaul Mullowney .seealso: MatCUSPARSESetFormat(), MatCUSPARSEStorageFormat 1646e057df02SPaul Mullowney E*/ 164736d62e41SPaul Mullowney typedef enum {MAT_CUSPARSE_MULT_DIAG, MAT_CUSPARSE_MULT_OFFDIAG, MAT_CUSPARSE_MULT, MAT_CUSPARSE_ALL} MatCUSPARSEFormatOperation; 1648e057df02SPaul Mullowney 164910e9db80SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 165010e9db80SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCreateAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 1651e057df02SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCUSPARSESetFormat(Mat,MatCUSPARSEFormatOperation,MatCUSPARSEStorageFormat); 16529ae82921SPaul Mullowney #endif 16539ae82921SPaul Mullowney 165490c53902SBarry Smith #if defined(PETSC_HAVE_CUSP) 1655014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCUSP(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 1656014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateAIJCUSP(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 1657e057df02SPaul Mullowney 16583f9c0db1SPaul Mullowney /*E 1659e057df02SPaul Mullowney MatCUSPStorageFormat - indicates the storage format for CUSP (GPU) 166036d62e41SPaul Mullowney matrices. 1661e057df02SPaul Mullowney 1662e057df02SPaul Mullowney Not Collective 1663e057df02SPaul Mullowney 16648468deeeSKarl Rupp + MAT_CUSP_CSR - Compressed Sparse Row 16658468deeeSKarl Rupp . MAT_CUSP_DIA - Diagonal 16668468deeeSKarl Rupp - MAT_CUSP_ELL - Ellpack 1667e057df02SPaul Mullowney 1668e057df02SPaul Mullowney Level: intermediate 1669e057df02SPaul Mullowney 1670e057df02SPaul Mullowney Any additions/changes here MUST also be made in include/finclude/petscmat.h 1671e057df02SPaul Mullowney 1672e057df02SPaul Mullowney .seealso: MatCUSPSetFormat(), MatCUSPFormatOperation 1673e057df02SPaul Mullowney E*/ 16743f9c0db1SPaul Mullowney typedef enum {MAT_CUSP_CSR, MAT_CUSP_DIA, MAT_CUSP_ELL} MatCUSPStorageFormat; 1675e057df02SPaul Mullowney 1676e057df02SPaul Mullowney /* these will be strings associated with enumerated type defined above */ 1677e057df02SPaul Mullowney PETSC_EXTERN const char *const MatCUSPStorageFormats[]; 1678e057df02SPaul Mullowney 16793f9c0db1SPaul Mullowney /*E 1680e057df02SPaul Mullowney MatCUSPFormatOperation - indicates the operation of CUSP (GPU) 168136d62e41SPaul Mullowney matrices whose operation should use a particular storage format. 1682e057df02SPaul Mullowney 1683e057df02SPaul Mullowney Not Collective 1684e057df02SPaul Mullowney 16858468deeeSKarl Rupp + MAT_CUSP_MULT_DIAG - sets the storage format for the diagonal matrix in the parallel MatMult 16868468deeeSKarl Rupp . MAT_CUSP_MULT_OFFDIAG - sets the storage format for the offdiagonal matrix in the parallel MatMult 16878468deeeSKarl Rupp . MAT_CUSP_MULT - sets the storage format for the entire matrix in the serial (single GPU) MatMult 16888468deeeSKarl Rupp - MAT_CUSP_ALL - sets the storage format for all CUSP (GPU) matrices 1689e057df02SPaul Mullowney 1690e057df02SPaul Mullowney Level: intermediate 1691e057df02SPaul Mullowney 1692e057df02SPaul Mullowney Any additions/changes here MUST also be made in include/finclude/petscmat.h 1693e057df02SPaul Mullowney 1694e057df02SPaul Mullowney .seealso: MatCUSPSetFormat(), MatCUSPStorageFormat 1695e057df02SPaul Mullowney E*/ 169636d62e41SPaul Mullowney typedef enum {MAT_CUSP_MULT_DIAG, MAT_CUSP_MULT_OFFDIAG, MAT_CUSP_MULT, MAT_CUSP_ALL} MatCUSPFormatOperation; 1697e057df02SPaul Mullowney 1698e057df02SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCUSPSetFormat(Mat,MatCUSPFormatOperation,MatCUSPStorageFormat); 1699e057df02SPaul Mullowney #endif 170090c53902SBarry Smith 1701d67ff14aSKarl Rupp #if defined(PETSC_HAVE_VIENNACL) 1702d67ff14aSKarl Rupp PETSC_EXTERN PetscErrorCode MatCreateSeqAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 1703d67ff14aSKarl Rupp PETSC_EXTERN PetscErrorCode MatCreateAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 1704d67ff14aSKarl Rupp #endif 1705d67ff14aSKarl Rupp 170654efbe56SHong Zhang /* 170754efbe56SHong Zhang PETSc interface to FFTW 170854efbe56SHong Zhang */ 170954efbe56SHong Zhang #if defined(PETSC_HAVE_FFTW) 1710014dd563SJed Brown PETSC_EXTERN PetscErrorCode VecScatterPetscToFFTW(Mat,Vec,Vec); 1711014dd563SJed Brown PETSC_EXTERN PetscErrorCode VecScatterFFTWToPetsc(Mat,Vec,Vec); 1712014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetVecsFFTW(Mat,Vec*,Vec*,Vec*); 171354efbe56SHong Zhang #endif 171454efbe56SHong Zhang 1715db31f6deSJed Brown #if defined(PETSC_HAVE_ELEMENTAL) 1716607a6623SBarry Smith PETSC_EXTERN PetscErrorCode PetscElementalInitializePackage(void); 1717db31f6deSJed Brown PETSC_EXTERN PetscErrorCode PetscElementalFinalizePackage(void); 1718db31f6deSJed Brown #endif 1719db31f6deSJed Brown 1720014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateNest(MPI_Comm,PetscInt,const IS[],PetscInt,const IS[],const Mat[],Mat*); 1721014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSize(Mat,PetscInt*,PetscInt*); 1722014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetISs(Mat,IS[],IS[]); 1723014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetLocalISs(Mat,IS[],IS[]); 1724014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSubMats(Mat,PetscInt*,PetscInt*,Mat***); 1725014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSubMat(Mat,PetscInt,PetscInt,Mat*); 172619fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatNestSetVecType(Mat,VecType); 1727014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestSetSubMats(Mat,PetscInt,const IS[],PetscInt,const IS[],const Mat[]); 1728014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestSetSubMat(Mat,PetscInt,PetscInt,Mat); 1729d8588912SDave May 17304325cce7SMatthew G Knepley PETSC_EXTERN PetscErrorCode MatChop(Mat,PetscReal); 1731e0a58c10SMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatComputeBandwidth(Mat,PetscReal,PetscInt*); 17324325cce7SMatthew G Knepley 173323a5497aSJed Brown #endif 1734