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" 109d615f992SSatish Balay #define MATSOLVERMKL_PARDISO "mkl_pardiso" 1102692d6eeSBarry Smith #define MATSOLVERPASTIX "pastix" 1112692d6eeSBarry Smith #define MATSOLVERMATLAB "matlab" 1122692d6eeSBarry Smith #define MATSOLVERPETSC "petsc" 1132692d6eeSBarry Smith #define MATSOLVERBAS "bas" 1149ae82921SPaul Mullowney #define MATSOLVERCUSPARSE "cusparse" 115aa5a9175SDahai Guo #define MATSOLVERBSTRM "bstrm" 116aa5a9175SDahai Guo #define MATSOLVERSBSTRM "sbstrm" 11715767789SHong Zhang #define MATSOLVERELEMENTAL "elemental" 11817f1a0eaSHong Zhang #define MATSOLVERCLIQUE "clique" 11969e15a41SShri Abhyankar #define MATSOLVERKLU "klu" 120c0cdd4a1SDahai Guo 121b24902e0SBarry Smith /*E 122b24902e0SBarry Smith MatFactorType - indicates what type of factorization is requested 123b24902e0SBarry Smith 124b24902e0SBarry Smith Level: beginner 125b24902e0SBarry Smith 126b24902e0SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 127b24902e0SBarry Smith 128c7393fdbSBarry Smith .seealso: MatSolverPackage, MatGetFactor() 129b24902e0SBarry Smith E*/ 130599ef60dSHong Zhang typedef enum {MAT_FACTOR_NONE, MAT_FACTOR_LU, MAT_FACTOR_CHOLESKY, MAT_FACTOR_ILU, MAT_FACTOR_ICC,MAT_FACTOR_ILUDT} MatFactorType; 131014dd563SJed Brown PETSC_EXTERN const char *const MatFactorTypes[]; 132e92e720dSBarry Smith 133014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactor(Mat,const MatSolverPackage,MatFactorType,Mat*); 134014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactorAvailable(Mat,const MatSolverPackage,MatFactorType,PetscBool *); 135014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFactorGetSolverPackage(Mat,const MatSolverPackage*); 136014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactorType(Mat,MatFactorType*); 1379989ab13SBarry Smith 138c06d978dSMatthew Knepley /* Logging support */ 1390700a824SBarry Smith #define MAT_FILE_CLASSID 1211216 /* used to indicate matrices in binary files */ 140014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_CLASSID; 141335efc43SPeter Brune PETSC_EXTERN PetscClassId MAT_COLORING_CLASSID; 142014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_FDCOLORING_CLASSID; 143014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_TRANSPOSECOLORING_CLASSID; 144014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_PARTITIONING_CLASSID; 145014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_COARSEN_CLASSID; 146014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_NULLSPACE_CLASSID; 147014dd563SJed Brown PETSC_EXTERN PetscClassId MATMFFD_CLASSID; 148c06d978dSMatthew Knepley 149ceb03754SKris Buschelman /*E 150ceb03754SKris Buschelman MatReuse - Indicates if matrices obtained from a previous call to MatGetSubMatrices() 151d6eff37eSBarry Smith or MatGetSubMatrix() are to be reused to store the new matrix values. For MatConvert() is used to indicate 152d6eff37eSBarry Smith that the input matrix is to be replaced with the converted matrix. 153ceb03754SKris Buschelman 154ceb03754SKris Buschelman Level: beginner 155ceb03754SKris Buschelman 156ceb03754SKris Buschelman Any additions/changes here MUST also be made in include/finclude/petscmat.h 157ceb03754SKris Buschelman 1580c451bc4SBarry Smith .seealso: MatGetSubMatrices(), MatGetSubMatrix(), MatDestroyMatrices(), MatConvert() 159ceb03754SKris Buschelman E*/ 160dfe085dbSJed Brown typedef enum {MAT_INITIAL_MATRIX,MAT_REUSE_MATRIX,MAT_IGNORE_MATRIX} MatReuse; 161ceb03754SKris Buschelman 1625494a064SHong Zhang /*E 1635494a064SHong Zhang MatGetSubMatrixOption - Indicates if matrices obtained from a call to MatGetSubMatrices() 164829201f2SHong Zhang include the matrix values. Currently it is only used by MatGetSeqNonzerostructure(). 1655494a064SHong Zhang 1665494a064SHong Zhang Level: beginner 1675494a064SHong Zhang 168829201f2SHong Zhang .seealso: MatGetSeqNonzerostructure() 1695494a064SHong Zhang E*/ 1705494a064SHong Zhang typedef enum {MAT_DO_NOT_GET_VALUES,MAT_GET_VALUES} MatGetSubMatrixOption; 1715494a064SHong Zhang 172607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatInitializePackage(void); 173c06d978dSMatthew Knepley 174014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreate(MPI_Comm,Mat*); 175014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetSizes(Mat,PetscInt,PetscInt,PetscInt,PetscInt); 17619fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatSetType(Mat,MatType); 177014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetFromOptions(Mat); 178ce1779c8SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatViewFromOptions(Mat A,const char prefix[],const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,prefix,name);} 179607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatRegisterAll(void); 180bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatRegister(const char[],PetscErrorCode(*)(Mat)); 181014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRegisterBaseName(const char[],const char[],const char[]); 182014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetOptionsPrefix(Mat,const char[]); 183014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAppendOptionsPrefix(Mat,const char[]); 184014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOptionsPrefix(Mat,const char*[]); 185f69a0ea3SMatthew Knepley 186014dd563SJed Brown PETSC_EXTERN PetscBool MatRegisterAllCalled; 187140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatList; 188140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatColoringList; 189140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatPartitioningList; 190140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatCoarsenList; 19128988994SBarry Smith 1923b224e63SBarry Smith /*E 193aa6c7ce3SBarry Smith MatStructure - Indicates if two matrices have the same nonzero structure 1943b224e63SBarry Smith 1953b224e63SBarry Smith Level: beginner 1963b224e63SBarry Smith 1973b224e63SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 1983b224e63SBarry Smith 199aa6c7ce3SBarry Smith .seealso: MatCopy(), MatAXPY() 2003b224e63SBarry Smith E*/ 201aa6c7ce3SBarry Smith typedef enum {DIFFERENT_NONZERO_PATTERN,SUBSET_NONZERO_PATTERN,SAME_NONZERO_PATTERN} MatStructure; 2023b224e63SBarry Smith 203014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqDense(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*); 204014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateDense(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*); 205014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 206014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 207014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *); 208014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithSplitArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],PetscInt[],PetscInt[],PetscScalar[],Mat*); 2098d7a6e47SBarry Smith 210014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 211014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 212*e561ad89SHong Zhang PETSC_EXTERN PetscErrorCode MatCreateMPIBAIJConcatenateSeqBAIJ(MPI_Comm,Mat,PetscInt,MatReuse,Mat*); 213*e561ad89SHong Zhang PETSC_EXTERN PetscErrorCode MatCreateMPIBAIJConcatenateSeqBAIJSymbolic(MPI_Comm,Mat,PetscInt,Mat*); 214*e561ad89SHong Zhang PETSC_EXTERN PetscErrorCode MatCreateMPIBAIJConcatenateSeqBAIJNumeric(MPI_Comm,Mat,PetscInt,Mat); 215014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat*); 216d21a29f3SJed Brown 217014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAdj(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscInt[],Mat*); 218014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 219d21a29f3SJed Brown 220014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 221014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPISBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *); 22238f409ebSLisandro Dalcin PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 223014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 2244cce697cSJed Brown PETSC_EXTERN PetscErrorCode MatXAIJSetPreallocation(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscInt[],const PetscInt[]); 225dfb205c3SBarry Smith 226014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateShell(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,void *,Mat*); 227014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateNormal(Mat,Mat*); 228014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateLRC(Mat,Mat,Mat,Mat*); 229014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateIS(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,ISLocalToGlobalMapping,Mat*); 230014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 231014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 232c0cdd4a1SDahai Guo 233014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 234014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 235014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqSBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 236014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPISBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 237c0cdd4a1SDahai Guo 238014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateScatter(MPI_Comm,VecScatter,Mat*); 239014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScatterSetVecScatter(Mat,VecScatter); 240014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScatterGetVecScatter(Mat,VecScatter*); 241014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateBlockMat(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,Mat*); 242014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeAddMat(Mat,Mat); 243014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeMerge(Mat); 244014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateComposite(MPI_Comm,PetscInt,const Mat*,Mat*); 2456d7c1e57SBarry Smith typedef enum {MAT_COMPOSITE_ADDITIVE,MAT_COMPOSITE_MULTIPLICATIVE} MatCompositeType; 246014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeSetType(Mat,MatCompositeType); 2476d7c1e57SBarry Smith 24819fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatCreateFFT(MPI_Comm,PetscInt,const PetscInt[],MatType,Mat*); 249014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqCUFFT(MPI_Comm,PetscInt,const PetscInt[],Mat*); 250dedccee8SHong Zhang 251014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateTranspose(Mat,Mat*); 252014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSubMatrix(Mat,IS,IS,Mat*); 253014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSubMatrixUpdate(Mat,Mat,IS,IS); 254014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateLocalRef(Mat,IS,IS,Mat*); 2551472f72bSBarry Smith 256014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPythonSetType(Mat,const char[]); 2571d6018f0SLisandro Dalcin 258014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetUp(Mat); 259014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroy(Mat*); 260e56f5c9eSBarry Smith PETSC_EXTERN PetscErrorCode MatGetNonzeroState(Mat,PetscObjectState*); 26121c89e3eSBarry Smith 262014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatConjugate(Mat); 263014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRealPart(Mat); 264014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatImaginaryPart(Mat); 265014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetDiagonalBlock(Mat,Mat*); 266014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetTrace(Mat,PetscScalar*); 267713ccfa9SJed Brown PETSC_EXTERN PetscErrorCode MatInvertBlockDiagonal(Mat,const PetscScalar **); 26899cafbc1SBarry Smith 2698ed539a5SBarry Smith /* ------------------------------------------------------------*/ 270014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 271014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 272014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesRow(Mat,PetscInt,const PetscScalar[]); 273014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesRowLocal(Mat,PetscInt,const PetscScalar[]); 274014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBatch(Mat,PetscInt,PetscInt,PetscInt[],const PetscScalar[]); 27573a71a0fSBarry Smith PETSC_EXTERN PetscErrorCode MatSetRandom(Mat,PetscRandom); 27684cb2905SBarry Smith 2772ef4de8bSBarry Smith /*S 2782ef4de8bSBarry Smith MatStencil - Data structure (C struct) for storing information about a single row or 279be479b30SJed Brown column of a matrix as indexed on an associated grid. 280be479b30SJed Brown 281be479b30SJed Brown Fortran usage is different, see MatSetValuesStencil() for details. 2822ef4de8bSBarry Smith 2832ef4de8bSBarry Smith Level: beginner 2842ef4de8bSBarry Smith 2852ef4de8bSBarry Smith Concepts: matrix; linear operator 2862ef4de8bSBarry Smith 287be479b30SJed Brown .seealso: MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockedStencil() 2882ef4de8bSBarry Smith S*/ 289435da068SBarry Smith typedef struct { 290c1ac3661SBarry Smith PetscInt k,j,i,c; 291435da068SBarry Smith } MatStencil; 2922ef4de8bSBarry Smith 293014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode); 294014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode); 295014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt); 296435da068SBarry Smith 297014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetColoring(Mat,ISColoring); 298014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesAdifor(Mat,PetscInt,void*); 2993a7fca6bSBarry Smith 300d91e6319SBarry Smith /*E 301d91e6319SBarry Smith MatAssemblyType - Indicates if the matrix is now to be used, or if you plan 302d91e6319SBarry Smith to continue to add values to it 303d91e6319SBarry Smith 304d91e6319SBarry Smith Level: beginner 305d91e6319SBarry Smith 306d91e6319SBarry Smith .seealso: MatAssemblyBegin(), MatAssemblyEnd() 307d91e6319SBarry Smith E*/ 3086d4a8577SBarry Smith typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType; 309014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssemblyBegin(Mat,MatAssemblyType); 310014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssemblyEnd(Mat,MatAssemblyType); 311014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssembled(Mat,PetscBool *); 3124f9c727eSBarry Smith 3131d73ed98SBarry Smith 31430de9b25SBarry Smith 315d91e6319SBarry Smith /*E 316d91e6319SBarry Smith MatOption - Options that may be set for a matrix and its behavior or storage 317d91e6319SBarry Smith 318d91e6319SBarry Smith Level: beginner 319d91e6319SBarry Smith 3200a835dfdSSatish Balay Any additions/changes here MUST also be made in include/finclude/petscmat.h 321d91e6319SBarry Smith 322382164b0SBarry Smith Developer Notes: Entries that are negative need not be called collectively by all processes. 323382164b0SBarry Smith 324d91e6319SBarry Smith .seealso: MatSetOption() 325d91e6319SBarry Smith E*/ 32692d4d306SBarry Smith typedef enum {MAT_OPTION_MIN = -5, 32792d4d306SBarry Smith MAT_NEW_NONZERO_LOCATION_ERR = -4, 32892d4d306SBarry Smith MAT_UNUSED_NONZERO_LOCATION_ERR = -3, 32992d4d306SBarry Smith MAT_NEW_NONZERO_ALLOCATION_ERR = -2, 33092d4d306SBarry Smith MAT_ROW_ORIENTED = -1, 331382164b0SBarry Smith MAT_SYMMETRIC = 1, 332382164b0SBarry Smith MAT_STRUCTURALLY_SYMMETRIC = 2, 333382164b0SBarry Smith MAT_NEW_DIAGONALS = 3, 334382164b0SBarry Smith MAT_IGNORE_OFF_PROC_ENTRIES = 4, 335382164b0SBarry Smith MAT_USE_HASH_TABLE = 5, 336382164b0SBarry Smith MAT_KEEP_NONZERO_PATTERN = 6, 337382164b0SBarry Smith MAT_IGNORE_ZERO_ENTRIES = 7, 338382164b0SBarry Smith MAT_USE_INODES = 8, 339382164b0SBarry Smith MAT_HERMITIAN = 9, 340382164b0SBarry Smith MAT_SYMMETRY_ETERNAL = 10, 34111e456e1SBarry Smith MAT_DUMMY = 11, 342382164b0SBarry Smith MAT_IGNORE_LOWER_TRIANGULAR = 12, 343382164b0SBarry Smith MAT_ERROR_LOWER_TRIANGULAR = 13, 344382164b0SBarry Smith MAT_GETROW_UPPERTRIANGULAR = 14, 345382164b0SBarry Smith MAT_SPD = 15, 34692d4d306SBarry Smith MAT_NO_OFF_PROC_ZERO_ROWS = 16, 34792d4d306SBarry Smith MAT_NO_OFF_PROC_ENTRIES = 17, 34892d4d306SBarry Smith MAT_NEW_NONZERO_LOCATIONS = 18, 34992d4d306SBarry Smith MAT_OPTION_MAX = 19} MatOption; 350e057df02SPaul Mullowney 351014dd563SJed Brown PETSC_EXTERN const char *MatOptions[]; 352014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetOption(Mat,MatOption,PetscBool ); 35319fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatGetType(Mat,MatType*); 35484cb2905SBarry Smith 355014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]); 356014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 357014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 358014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowUpperTriangular(Mat); 359014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRowUpperTriangular(Mat); 360014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 361014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 362014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnVector(Mat,Vec,PetscInt); 3638c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJGetArray(Mat,PetscScalar *[]); 3648c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJRestoreArray(Mat,PetscScalar *[]); 36521e72a00SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJGetMaxRowNonzeros(Mat,PetscInt*); 366bd04181cSBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJSetValuesLocalFast(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 3678c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatDenseGetArray(Mat,PetscScalar *[]); 3688c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatDenseRestoreArray(Mat,PetscScalar *[]); 369014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBlockSize(Mat,PetscInt *); 370014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSize(Mat,PetscInt); 371014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBlockSizes(Mat,PetscInt *,PetscInt *); 372014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSizes(Mat,PetscInt,PetscInt); 37333d57670SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSizesFromMats(Mat,Mat,Mat); 374014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNThreads(Mat,PetscInt); 375014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNThreads(Mat,PetscInt*); 3761620fd73SBarry Smith 377014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMult(Mat,Vec,Vec); 378014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultDiagonalBlock(Mat,Vec,Vec); 379014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultAdd(Mat,Vec,Vec,Vec); 380014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTranspose(Mat,Vec,Vec); 381014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultHermitianTranspose(Mat,Vec,Vec); 382014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsTranspose(Mat,Mat,PetscReal,PetscBool *); 383014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitianTranspose(Mat,Mat,PetscReal,PetscBool *); 384014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeAdd(Mat,Vec,Vec,Vec); 385014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultHermitianTransposeAdd(Mat,Vec,Vec,Vec); 386014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultConstrained(Mat,Vec,Vec); 387014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeConstrained(Mat,Vec,Vec); 388014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatSolve(Mat,Mat,Mat); 389f9426fe0SMark Adams PETSC_EXTERN PetscErrorCode MatResidual(Mat,Vec,Vec,Vec); 390f5edf698SHong Zhang 391d91e6319SBarry Smith /*E 392d91e6319SBarry Smith MatDuplicateOption - Indicates if a duplicated sparse matrix should have 393d91e6319SBarry Smith its numerical values copied over or just its nonzero structure. 394d91e6319SBarry Smith 395d91e6319SBarry Smith Level: beginner 396d91e6319SBarry Smith 397d91e6319SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 398d91e6319SBarry Smith 39970dcbbb9SBarry Smith $ MAT_SHARE_NONZERO_PATTERN - the i and j arrays in the new matrix will be shared with the original matrix 40070dcbbb9SBarry Smith $ this also triggers the MAT_DO_NOT_COPY_VALUES option. This is used when you 40170dcbbb9SBarry Smith $ have several matrices with the same nonzero pattern. 40270dcbbb9SBarry Smith 403d91e6319SBarry Smith .seealso: MatDuplicate() 404d91e6319SBarry Smith E*/ 40570dcbbb9SBarry Smith typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES,MAT_SHARE_NONZERO_PATTERN} MatDuplicateOption; 4062e8a6d31SBarry Smith 40719fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatConvert(Mat,MatType,MatReuse,Mat*); 408014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDuplicate(Mat,MatDuplicateOption,Mat*); 40994a9d846SBarry Smith 410cb5b572fSBarry Smith 411014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCopy(Mat,Mat,MatStructure); 412014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatView(Mat,PetscViewer); 413014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsSymmetric(Mat,PetscReal,PetscBool *); 414014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsStructurallySymmetric(Mat,PetscBool *); 415014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitian(Mat,PetscReal,PetscBool *); 416014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsSymmetricKnown(Mat,PetscBool *,PetscBool *); 417014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitianKnown(Mat,PetscBool *,PetscBool *); 418014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMissingDiagonal(Mat,PetscBool *,PetscInt *); 419014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLoad(Mat, PetscViewer); 4207b80b807SBarry Smith 4211a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool *); 4221a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool *); 4231a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool *); 4241a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool *); 425d4fbbf0eSBarry Smith 426d91e6319SBarry Smith /*S 427d91e6319SBarry Smith MatInfo - Context of matrix information, used with MatGetInfo() 428d91e6319SBarry Smith 429d91e6319SBarry Smith In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE 430d91e6319SBarry Smith 431d91e6319SBarry Smith Level: intermediate 432d91e6319SBarry Smith 433d91e6319SBarry Smith Concepts: matrix^nonzero information 434d91e6319SBarry Smith 435d9274352SBarry Smith .seealso: MatGetInfo(), MatInfoType 436d91e6319SBarry Smith S*/ 4374e220ebcSLois Curfman McInnes typedef struct { 438b0a32e0cSBarry Smith PetscLogDouble block_size; /* block size */ 439b0a32e0cSBarry Smith PetscLogDouble nz_allocated,nz_used,nz_unneeded; /* number of nonzeros */ 440b0a32e0cSBarry Smith PetscLogDouble memory; /* memory allocated */ 441b0a32e0cSBarry Smith PetscLogDouble assemblies; /* number of matrix assemblies called */ 442b0a32e0cSBarry Smith PetscLogDouble mallocs; /* number of mallocs during MatSetValues() */ 443b0a32e0cSBarry Smith PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */ 444b0a32e0cSBarry Smith PetscLogDouble factor_mallocs; /* number of mallocs during factorization */ 4454e220ebcSLois Curfman McInnes } MatInfo; 4464e220ebcSLois Curfman McInnes 447d9274352SBarry Smith /*E 448d9274352SBarry Smith MatInfoType - Indicates if you want information about the local part of the matrix, 449d9274352SBarry Smith the entire parallel matrix or the maximum over all the local parts. 450d9274352SBarry Smith 451d9274352SBarry Smith Level: beginner 452d9274352SBarry Smith 453d9274352SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 454d9274352SBarry Smith 455d9274352SBarry Smith .seealso: MatGetInfo(), MatInfo 456d9274352SBarry Smith E*/ 4577b80b807SBarry Smith typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType; 458014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetInfo(Mat,MatInfoType,MatInfo*); 459014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetDiagonal(Mat,Vec); 460014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMax(Mat,Vec,PetscInt[]); 461014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMin(Mat,Vec,PetscInt[]); 462014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMaxAbs(Mat,Vec,PetscInt[]); 463014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMinAbs(Mat,Vec,PetscInt[]); 464014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowSum(Mat,Vec); 465014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTranspose(Mat,MatReuse,Mat*); 466014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatHermitianTranspose(Mat,MatReuse,Mat*); 467014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPermute(Mat,IS,IS,Mat *); 468014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalScale(Mat,Vec,Vec); 469014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalSet(Mat,Vec,InsertMode); 470014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatEqual(Mat,Mat,PetscBool *); 471014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultEqual(Mat,Mat,PetscInt,PetscBool *); 472014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultAddEqual(Mat,Mat,PetscInt,PetscBool *); 473014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeEqual(Mat,Mat,PetscInt,PetscBool *); 474014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscBool *); 4757b80b807SBarry Smith 476014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNorm(Mat,NormType,PetscReal *); 477014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnNorms(Mat,NormType,PetscReal *); 478014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroEntries(Mat); 479014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 480014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsIS(Mat,IS,PetscScalar,Vec,Vec); 481014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsStencil(Mat,PetscInt,const MatStencil [],PetscScalar,Vec,Vec); 482014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsStencil(Mat,PetscInt,const MatStencil[],PetscScalar,Vec,Vec); 483014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumns(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 484014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsIS(Mat,IS,PetscScalar,Vec,Vec); 4857b80b807SBarry Smith 486014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSize(Mat,PetscInt*,PetscInt*); 487014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalSize(Mat,PetscInt*,PetscInt*); 488014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRange(Mat,PetscInt*,PetscInt*); 489014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRanges(Mat,const PetscInt**); 490014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRangeColumn(Mat,PetscInt*,PetscInt*); 491014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRangesColumn(Mat,const PetscInt**); 492566876eaSJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipIS(Mat,IS*,IS*); 4935ef9f2a5SBarry Smith 494014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]); 495014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSubMatricesParallel(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]); 496014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroyMatrices(PetscInt,Mat *[]); 497014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSubMatrix(Mat,IS,IS,MatReuse,Mat *); 498014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalSubMatrix(Mat,IS,IS,Mat*); 499014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreLocalSubMatrix(Mat,IS,IS,Mat*); 500014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSeqNonzeroStructure(Mat,Mat*); 501014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroySeqNonzeroStructure(Mat*); 5027b80b807SBarry Smith 503014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJConcatenateSeqAIJ(MPI_Comm,Mat,PetscInt,MatReuse,Mat*); 504014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJConcatenateSeqAIJSymbolic(MPI_Comm,Mat,PetscInt,Mat*); 505014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJConcatenateSeqAIJNumeric(MPI_Comm,Mat,PetscInt,Mat); 506014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJ(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*); 507014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJSymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*); 508014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJNumeric(Mat,Mat); 509014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMat(Mat,MatReuse,Mat*); 510014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*); 511014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,Mat*); 512014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetGhosts(Mat, PetscInt *,const PetscInt *[]); 5138efafbd8SBarry Smith 514014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt); 5157b80b807SBarry Smith 516014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*); 517014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMultSymbolic(Mat,Mat,PetscReal,Mat*); 518014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMultNumeric(Mat,Mat,Mat); 51922440eb1SKris Buschelman 5207bab7c10SHong Zhang PETSC_EXTERN PetscErrorCode MatMatMatMult(Mat,Mat,Mat,MatReuse,PetscReal,Mat*); 5217bab7c10SHong Zhang PETSC_EXTERN PetscErrorCode MatMatMatMultSymbolic(Mat,Mat,Mat,PetscReal,Mat*); 5227bab7c10SHong Zhang PETSC_EXTERN PetscErrorCode MatMatMatMultNumeric(Mat,Mat,Mat,Mat); 5237bab7c10SHong Zhang 524014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*); 525014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAPSymbolic(Mat,Mat,PetscReal,Mat*); 526014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAPNumeric(Mat,Mat,Mat); 527014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARt(Mat,Mat,MatReuse,PetscReal,Mat*); 528014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARtSymbolic(Mat,Mat,PetscReal,Mat*); 529014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARtNumeric(Mat,Mat,Mat); 53022440eb1SKris Buschelman 531014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeMatMult(Mat,Mat,MatReuse,PetscReal,Mat*); 532014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposetMatMultSymbolic(Mat,Mat,PetscReal,Mat*); 533014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposetMatMultNumeric(Mat,Mat,Mat); 534014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMult(Mat,Mat,MatReuse,PetscReal,Mat*); 535014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMultSymbolic(Mat,Mat,PetscReal,Mat*); 536014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMultNumeric(Mat,Mat,Mat); 537bc011b1eSHong Zhang 538014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAXPY(Mat,PetscScalar,Mat,MatStructure); 539014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAYPX(Mat,PetscScalar,Mat,MatStructure); 5401c741599SHong Zhang 541014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScale(Mat,PetscScalar); 542014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShift(Mat,PetscScalar); 5437b80b807SBarry Smith 544014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping); 545014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*); 546a93ff8c4SPeter Brune PETSC_EXTERN PetscErrorCode MatGetLayouts(Mat,PetscLayout*,PetscLayout*); 547014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 548014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsLocalIS(Mat,IS,PetscScalar,Vec,Vec); 549014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 550014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocalIS(Mat,IS,PetscScalar,Vec,Vec); 551014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 552014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 553052efed2SBarry Smith 554014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStashSetInitialSize(Mat,PetscInt,PetscInt); 555014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*); 55690f02eecSBarry Smith 557014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInterpolate(Mat,Vec,Vec); 558014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInterpolateAdd(Mat,Vec,Vec,Vec); 559014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestrict(Mat,Vec,Vec); 5602a7a6963SBarry Smith PETSC_EXTERN PetscErrorCode MatCreateVecs(Mat,Vec*,Vec*); 56184cc2a5aSBarry Smith PETSC_DEPRECATED("Use MatCreateVecs()") PETSC_STATIC_INLINE PetscErrorCode MatGetVecs(Mat mat,Vec *x,Vec *y) {return MatCreateVecs(mat,x,y);} 562b2bf6370SHong Zhang PETSC_EXTERN PetscErrorCode MatGetRedundantMatrix(Mat,PetscInt,MPI_Comm,MatReuse,Mat*); 563014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetMultiProcBlock(Mat,MPI_Comm,MatReuse,Mat*); 564014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFindZeroDiagonals(Mat,IS*); 5653a062f41SBarry Smith PETSC_EXTERN PetscErrorCode MatFindOffBlockDiagonalEntries(Mat,IS*); 566bd481603SBarry Smith 567bd481603SBarry Smith /*MC 5681620fd73SBarry Smith MatSetValue - Set a single entry into a matrix. 5691620fd73SBarry Smith 5701620fd73SBarry Smith Not collective 5711620fd73SBarry Smith 572a9834a7dSSatish Balay Synopsis: 573a9834a7dSSatish Balay #include <petscmat.h> 574a9834a7dSSatish Balay PetscErrorCode MatSetValue(Mat m,PetscInt row,PetscInt col,PetscScalar value,InsertMode mode) 575a9834a7dSSatish Balay 5761620fd73SBarry Smith Input Parameters: 5771620fd73SBarry Smith + m - the matrix 5781620fd73SBarry Smith . row - the row location of the entry 5791620fd73SBarry Smith . col - the column location of the entry 5801620fd73SBarry Smith . value - the value to insert 5811620fd73SBarry Smith - mode - either INSERT_VALUES or ADD_VALUES 5821620fd73SBarry Smith 5831620fd73SBarry Smith Notes: 5841620fd73SBarry Smith For efficiency one should use MatSetValues() and set several or many 5851620fd73SBarry Smith values simultaneously if possible. 5861620fd73SBarry Smith 5871620fd73SBarry Smith Level: beginner 5881620fd73SBarry Smith 5891620fd73SBarry Smith .seealso: MatSetValues(), MatSetValueLocal() 5901620fd73SBarry Smith M*/ 5911620fd73SBarry 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);} 5921620fd73SBarry Smith 5931620fd73SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatGetValue(Mat v,PetscInt i,PetscInt j,PetscScalar *va) {return MatGetValues(v,1,&i,1,&j,va);} 5941620fd73SBarry Smith 5951620fd73SBarry 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);} 5961620fd73SBarry Smith 5971620fd73SBarry Smith /*MC 5980d2a7ff2SSatish Balay MatPreallocateInitialize - Begins the block of code that will count the number of nonzeros per 599bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 600bd481603SBarry Smith 601bd481603SBarry Smith Synopsis: 602aaa7dc30SBarry Smith #include <petscmat.h> 603c1ac3661SBarry Smith PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz) 604bd481603SBarry Smith 605bd481603SBarry Smith Collective on MPI_Comm 606bd481603SBarry Smith 607bd481603SBarry Smith Input Parameters: 608bd481603SBarry Smith + comm - the communicator that will share the eventually allocated matrix 609859fcb39SBarry Smith . nrows - the number of LOCAL rows in the matrix 610859fcb39SBarry Smith - ncols - the number of LOCAL columns in the matrix 611bd481603SBarry Smith 612bd481603SBarry Smith Output Parameters: 613bd481603SBarry Smith + dnz - the array that will be passed to the matrix preallocation routines 614bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 615bd481603SBarry Smith 616bd481603SBarry Smith Level: intermediate 617bd481603SBarry Smith 618bd481603SBarry Smith Notes: 619a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 620bd481603SBarry Smith 6211d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 622bd481603SBarry Smith 623bd481603SBarry Smith Use MatPreallocateInitializeSymmetric() for symmetric matrices (MPISBAIJ matrices) 624bd481603SBarry Smith 6251620fd73SBarry Smith This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize(). 6261620fd73SBarry Smith 627bd481603SBarry Smith Concepts: preallocation^Matrix 628bd481603SBarry Smith 629d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(), 630d6e23781SBarry Smith MatPreallocateInitializeSymmetric(), MatPreallocateSymmetricSetLocalBlock() 631bd481603SBarry Smith M*/ 632c4f061fbSSatish Balay #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \ 6337c922b88SBarry Smith { \ 634a89bc211SBarry Smith PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__start,__end; \ 6351795a4d1SJed Brown _4_ierr = PetscCalloc2(__nrows,&dnz,__nrows,&onz);CHKERRQ(_4_ierr); \ 6361795a4d1SJed Brown __start = 0; __end = __start; \ 637c1ac3661SBarry Smith _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ctmp;\ 638a89bc211SBarry Smith _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows; 6397c922b88SBarry Smith 640bd481603SBarry Smith /*MC 641bd481603SBarry Smith MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be 642bd481603SBarry Smith inserted using a local number of the rows and columns 643bd481603SBarry Smith 644bd481603SBarry Smith Synopsis: 645aaa7dc30SBarry Smith #include <petscmat.h> 646c1ac3661SBarry Smith PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 647bd481603SBarry Smith 648bd481603SBarry Smith Not Collective 649bd481603SBarry Smith 650bd481603SBarry Smith Input Parameters: 651784ac674SJed Brown + map - the row mapping from local numbering to global numbering 652bd481603SBarry Smith . nrows - the number of rows indicated 6531d73ed98SBarry Smith . rows - the indices of the rows 654784ac674SJed Brown . cmap - the column mapping from local to global numbering 655bd481603SBarry Smith . ncols - the number of columns in the matrix 656bd481603SBarry Smith . cols - the columns indicated 657bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 658bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 659bd481603SBarry Smith 660bd481603SBarry Smith Level: intermediate 661bd481603SBarry Smith 662bd481603SBarry Smith Notes: 663a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 664bd481603SBarry Smith 6651d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 666bd481603SBarry Smith 667bd481603SBarry Smith Concepts: preallocation^Matrix 668bd481603SBarry Smith 669d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(), 670d6e23781SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock() 671bd481603SBarry Smith M*/ 672784ac674SJed Brown #define MatPreallocateSetLocal(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \ 673c4f061fbSSatish Balay {\ 674c1ac3661SBarry Smith PetscInt __l;\ 675784ac674SJed Brown _4_ierr = ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\ 676784ac674SJed Brown _4_ierr = ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\ 677c4f061fbSSatish Balay for (__l=0;__l<nrows;__l++) {\ 678ef66eb69SBarry Smith _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 679c4f061fbSSatish Balay }\ 680c4f061fbSSatish Balay } 681c4f061fbSSatish Balay 682bd481603SBarry Smith /*MC 683d6e23781SBarry Smith MatPreallocateSetLocalBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be 684bd481603SBarry Smith inserted using a local number of the rows and columns 685bd481603SBarry Smith 686bd481603SBarry Smith Synopsis: 687aaa7dc30SBarry Smith #include <petscmat.h> 688d6e23781SBarry Smith PetscErrorCode MatPreallocateSetLocalBlock(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 689d6e23781SBarry Smith 690d6e23781SBarry Smith Not Collective 691d6e23781SBarry Smith 692d6e23781SBarry Smith Input Parameters: 693d6e23781SBarry Smith + map - the row mapping from local numbering to global numbering 694d6e23781SBarry Smith . nrows - the number of rows indicated 695d6e23781SBarry Smith . rows - the indices of the rows 696d6e23781SBarry Smith . cmap - the column mapping from local to global numbering 697d6e23781SBarry Smith . ncols - the number of columns in the matrix 698d6e23781SBarry Smith . cols - the columns indicated 699d6e23781SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 700d6e23781SBarry Smith - ozn - the other array passed to the matrix preallocation routines 701d6e23781SBarry Smith 702d6e23781SBarry Smith Level: intermediate 703d6e23781SBarry Smith 704d6e23781SBarry Smith Notes: 705d6e23781SBarry Smith See Users-Manual: ch_performance for more details. 706d6e23781SBarry Smith 707d6e23781SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 708d6e23781SBarry Smith 709d6e23781SBarry Smith Concepts: preallocation^Matrix 710d6e23781SBarry Smith 711d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(), 712d6e23781SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock() 713d6e23781SBarry Smith M*/ 714d6e23781SBarry Smith #define MatPreallocateSetLocalBlock(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \ 715d6e23781SBarry Smith {\ 716d6e23781SBarry Smith PetscInt __l;\ 717d6e23781SBarry Smith _4_ierr = ISLocalToGlobalMappingApplyBlock(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\ 718d6e23781SBarry Smith _4_ierr = ISLocalToGlobalMappingApplyBlock(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\ 719d6e23781SBarry Smith for (__l=0;__l<nrows;__l++) {\ 720d6e23781SBarry Smith _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 721d6e23781SBarry Smith }\ 722d6e23781SBarry Smith } 723d6e23781SBarry Smith 724d6e23781SBarry Smith /*MC 725d6e23781SBarry Smith MatPreallocateSymmetricSetLocalBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be 726d6e23781SBarry Smith inserted using a local number of the rows and columns 727d6e23781SBarry Smith 728d6e23781SBarry Smith Synopsis: 729d6e23781SBarry Smith #include <petscmat.h> 730d6e23781SBarry Smith PetscErrorCode MatPreallocateSymmetricSetLocalBlock(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 731bd481603SBarry Smith 732bd481603SBarry Smith Not Collective 733bd481603SBarry Smith 734bd481603SBarry Smith Input Parameters: 735bd481603SBarry Smith + map - the mapping between local numbering and global numbering 736bd481603SBarry Smith . nrows - the number of rows indicated 7371d73ed98SBarry Smith . rows - the indices of the rows 738bd481603SBarry Smith . ncols - the number of columns in the matrix 739bd481603SBarry Smith . cols - the columns indicated 740bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 741bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 742bd481603SBarry Smith 743bd481603SBarry Smith Level: intermediate 744bd481603SBarry Smith 745bd481603SBarry Smith Notes: 746a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 747bd481603SBarry Smith 748bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 749bd481603SBarry Smith 750bd481603SBarry Smith Concepts: preallocation^Matrix 751bd481603SBarry Smith 752d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateInitialize(), 753d6e23781SBarry Smith MatPreallocateInitialize(), MatPreallocateSetLocal() 754bd481603SBarry Smith M*/ 755d6e23781SBarry Smith #define MatPreallocateSymmetricSetLocalBlock(map,nrows,rows,ncols,cols,dnz,onz) 0;\ 756d3d32019SBarry Smith {\ 757c1ac3661SBarry Smith PetscInt __l;\ 758d6e23781SBarry Smith _4_ierr = ISLocalToGlobalMappingApplyBlock(map,nrows,rows,rows);CHKERRQ(_4_ierr);\ 759d6e23781SBarry Smith _4_ierr = ISLocalToGlobalMappingApplyBlock(map,ncols,cols,cols);CHKERRQ(_4_ierr);\ 760d3d32019SBarry Smith for (__l=0;__l<nrows;__l++) {\ 761d6e23781SBarry Smith _4_ierr = MatPreallocateSymmetricSetBlock((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 762d3d32019SBarry Smith }\ 763d3d32019SBarry Smith } 764bd481603SBarry Smith /*MC 765bd481603SBarry Smith MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be 766bd481603SBarry Smith inserted using a local number of the rows and columns 767bd481603SBarry Smith 768bd481603SBarry Smith Synopsis: 769aaa7dc30SBarry Smith #include <petscmat.h> 770c1ac3661SBarry Smith PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 771bd481603SBarry Smith 772bd481603SBarry Smith Not Collective 773bd481603SBarry Smith 774bd481603SBarry Smith Input Parameters: 77564075487SBarry Smith + row - the row 776bd481603SBarry Smith . ncols - the number of columns in the matrix 777a50f8bf6SHong Zhang - cols - the columns indicated 778a50f8bf6SHong Zhang 779a50f8bf6SHong Zhang Output Parameters: 780a50f8bf6SHong Zhang + dnz - the array that will be passed to the matrix preallocation routines 781bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 782bd481603SBarry Smith 783bd481603SBarry Smith Level: intermediate 784bd481603SBarry Smith 785bd481603SBarry Smith Notes: 786a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 787bd481603SBarry Smith 788bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 789bd481603SBarry Smith 7901620fd73SBarry Smith This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize(). 7911620fd73SBarry Smith 792bd481603SBarry Smith Concepts: preallocation^Matrix 793bd481603SBarry Smith 794d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(), 795d6e23781SBarry Smith MatPreallocateInitialize(), MatPreallocateSetLocal() 796bd481603SBarry Smith M*/ 797c4f061fbSSatish Balay #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\ 798c1ac3661SBarry Smith { PetscInt __i; \ 799e32f2f54SBarry 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);\ 800e32f2f54SBarry 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);\ 8017c922b88SBarry Smith for (__i=0; __i<nc; __i++) {\ 80264075487SBarry Smith if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \ 8037cc688d7SBarry Smith else dnz[row - __rstart]++;\ 8047c922b88SBarry Smith }\ 8057c922b88SBarry Smith } 8067c922b88SBarry Smith 807bd481603SBarry Smith /*MC 808d6e23781SBarry Smith MatPreallocateSymmetricSetBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be 809bd481603SBarry Smith inserted using a local number of the rows and columns 810bd481603SBarry Smith 811bd481603SBarry Smith Synopsis: 812aaa7dc30SBarry Smith #include <petscmat.h> 813d6e23781SBarry Smith PetscErrorCode MatPreallocateSymmetricSetBlock(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 814bd481603SBarry Smith 815bd481603SBarry Smith Not Collective 816bd481603SBarry Smith 817bd481603SBarry Smith Input Parameters: 818bd481603SBarry Smith + nrows - the number of rows indicated 8191d73ed98SBarry Smith . rows - the indices of the rows 820bd481603SBarry Smith . ncols - the number of columns in the matrix 821bd481603SBarry Smith . cols - the columns indicated 822bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 823bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 824bd481603SBarry Smith 825bd481603SBarry Smith Level: intermediate 826bd481603SBarry Smith 827bd481603SBarry Smith Notes: 828a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 829bd481603SBarry Smith 830bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 831bd481603SBarry Smith 8321620fd73SBarry Smith This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize(). 8331620fd73SBarry Smith 834bd481603SBarry Smith Concepts: preallocation^Matrix 835bd481603SBarry Smith 836d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateInitialize(), 837d6e23781SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocal() 838bd481603SBarry Smith M*/ 839d6e23781SBarry Smith #define MatPreallocateSymmetricSetBlock(row,nc,cols,dnz,onz) 0;\ 840c1ac3661SBarry Smith { PetscInt __i; \ 841d3d32019SBarry Smith for (__i=0; __i<nc; __i++) {\ 842d3d32019SBarry Smith if (cols[__i] >= __end) onz[row - __rstart]++; \ 843d3d32019SBarry Smith else if (cols[__i] >= row) dnz[row - __rstart]++;\ 844d3d32019SBarry Smith }\ 845d3d32019SBarry Smith } 846d3d32019SBarry Smith 847bd481603SBarry Smith /*MC 84816371a99SBarry Smith MatPreallocateLocation - An alternative to MatPreallocationSet() that puts the nonzero locations into the matrix if it exists 84916371a99SBarry Smith 85016371a99SBarry Smith Synopsis: 851aaa7dc30SBarry Smith #include <petscmat.h> 85216371a99SBarry Smith PetscErrorCode MatPreallocateLocations(Mat A,PetscInt row,PetscInt ncols,PetscInt *cols,PetscInt *dnz,PetscInt *onz) 85316371a99SBarry Smith 85416371a99SBarry Smith Not Collective 85516371a99SBarry Smith 85616371a99SBarry Smith Input Parameters: 85716371a99SBarry Smith . A - matrix 85816371a99SBarry Smith . row - row where values exist (must be local to this process) 85916371a99SBarry Smith . ncols - number of columns 86016371a99SBarry Smith . cols - columns with nonzeros 86116371a99SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 86216371a99SBarry Smith - ozn - the other array passed to the matrix preallocation routines 86316371a99SBarry Smith 86416371a99SBarry Smith Level: intermediate 86516371a99SBarry Smith 86616371a99SBarry Smith Notes: 867a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 86816371a99SBarry Smith 86916371a99SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 87016371a99SBarry Smith 87116371a99SBarry Smith This is a MACRO not a function because it uses a bunch of variables private to the MatPreallocation.... routines. 87216371a99SBarry Smith 87316371a99SBarry Smith Concepts: preallocation^Matrix 87416371a99SBarry Smith 875d6e23781SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(), 876d6e23781SBarry Smith MatPreallocateSymmetricSetLocalBlock() 87716371a99SBarry Smith M*/ 8780298fd71SBarry 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);} 87916371a99SBarry Smith 88016371a99SBarry Smith 88116371a99SBarry Smith /*MC 8820d2a7ff2SSatish Balay MatPreallocateFinalize - Ends the block of code that will count the number of nonzeros per 883bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 884bd481603SBarry Smith 885bd481603SBarry Smith Synopsis: 886aaa7dc30SBarry Smith #include <petscmat.h> 887c1ac3661SBarry Smith PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz) 888bd481603SBarry Smith 889bd481603SBarry Smith Collective on MPI_Comm 890bd481603SBarry Smith 891bd481603SBarry Smith Input Parameters: 89216371a99SBarry Smith + dnz - the array that was be passed to the matrix preallocation routines 893bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 894bd481603SBarry Smith 895bd481603SBarry Smith Level: intermediate 896bd481603SBarry Smith 897bd481603SBarry Smith Notes: 898a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 899bd481603SBarry Smith 900bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 901bd481603SBarry Smith 9021620fd73SBarry Smith This is a MACRO not a function because it closes the { started in MatPreallocateInitialize(). 9031620fd73SBarry Smith 904bd481603SBarry Smith Concepts: preallocation^Matrix 905bd481603SBarry Smith 906d6e23781SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(), 907d6e23781SBarry Smith MatPreallocateSymmetricSetLocalBlock() 908bd481603SBarry Smith M*/ 909a89bc211SBarry Smith #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree2(dnz,onz);CHKERRQ(_4_ierr);} 9107c922b88SBarry Smith 9117b80b807SBarry Smith /* Routines unique to particular data structures */ 912014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellGetContext(Mat,void *); 913698d4c6aSKris Buschelman 914014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInodeAdjustForInodes(Mat,IS*,IS*); 915014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *); 916698d4c6aSKris Buschelman 917014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetColumnIndices(Mat,PetscInt[]); 918014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetColumnIndices(Mat,PetscInt[]); 919014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 920014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 921014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 922014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJFromTriple(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*,PetscInt,PetscBool); 9237b80b807SBarry Smith 924a96a251dSBarry Smith #define MAT_SKIP_ALLOCATION -4 925a96a251dSBarry Smith 926014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]); 927014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]); 928014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]); 929273d9f13SBarry Smith 930014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 931014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 932014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 933014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []); 934014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 935014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]); 936014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 937014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]); 938014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIDenseSetPreallocation(Mat,PetscScalar[]); 939014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqDenseSetPreallocation(Mat,PetscScalar[]); 9409230625dSJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,const PetscInt*[]); 9419230625dSJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,const PetscInt*[]); 942014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAdjCreateNonemptySubcommMat(Mat,Mat*); 943273d9f13SBarry Smith 944014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqDenseSetLDA(Mat,PetscInt); 945014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDenseGetLocalMatrix(Mat,Mat*); 9461b807ce4Svictorle 947014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStoreValues(Mat); 948014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRetrieveValues(Mat); 9492e8a6d31SBarry Smith 950014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDAADSetCtx(Mat,void*); 9513a7fca6bSBarry Smith 952014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFindNonzeroRows(Mat,IS*); 9537b80b807SBarry Smith /* 9547b80b807SBarry Smith These routines are not usually accessed directly, rather solving is 95594b7f48cSBarry Smith done through the KSP and PC interfaces. 9567b80b807SBarry Smith */ 9577b80b807SBarry Smith 95876bdecfbSBarry Smith /*J 9598f6c3df8SBarry Smith MatOrderingType - String with the name of a PETSc matrix ordering 960d9274352SBarry Smith 961d9274352SBarry Smith Level: beginner 962d9274352SBarry Smith 963d9274352SBarry Smith .seealso: MatGetOrdering() 96476bdecfbSBarry Smith J*/ 96519fd82e9SBarry Smith typedef const char* MatOrderingType; 9662692d6eeSBarry Smith #define MATORDERINGNATURAL "natural" 9672692d6eeSBarry Smith #define MATORDERINGND "nd" 9682692d6eeSBarry Smith #define MATORDERING1WD "1wd" 9692692d6eeSBarry Smith #define MATORDERINGRCM "rcm" 9702692d6eeSBarry Smith #define MATORDERINGQMD "qmd" 9712692d6eeSBarry Smith #define MATORDERINGROWLENGTH "rowlength" 972510184a7SJed Brown #define MATORDERINGWBM "wbm" 973fc1bef75SJed Brown #define MATORDERINGSPECTRAL "spectral" 974312e372aSBarry Smith #define MATORDERINGAMD "amd" /* only works if UMFPACK is installed with PETSc */ 975b12f92e5SBarry Smith 97619fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatGetOrdering(Mat,MatOrderingType,IS*,IS*); 977140e18c1SBarry Smith PETSC_EXTERN PetscErrorCode MatGetOrderingList(PetscFunctionList*); 978bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatOrderingRegister(const char[],PetscErrorCode(*)(Mat,MatOrderingType,IS*,IS*)); 979607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatOrderingRegisterAll(void); 980014dd563SJed Brown PETSC_EXTERN PetscBool MatOrderingRegisterAllCalled; 981140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatOrderingList; 982d4fbbf0eSBarry Smith 983014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS); 984fc1bef75SJed Brown PETSC_EXTERN PetscErrorCode MatCreateLaplacian(Mat,PetscReal,PetscBool,Mat*); 985a2ce50c7SBarry Smith 986d91e6319SBarry Smith /*S 987d90ac83dSHong Zhang MatFactorShiftType - Numeric Shift. 988d90ac83dSHong Zhang 989d90ac83dSHong Zhang Level: beginner 990d90ac83dSHong Zhang 991d90ac83dSHong Zhang S*/ 992d90ac83dSHong Zhang typedef enum {MAT_SHIFT_NONE,MAT_SHIFT_NONZERO,MAT_SHIFT_POSITIVE_DEFINITE,MAT_SHIFT_INBLOCKS} MatFactorShiftType; 9936a6fc655SJed Brown PETSC_EXTERN const char *const MatFactorShiftTypes[]; 9945e9742b9SJed Brown PETSC_EXTERN const char *const MatFactorShiftTypesDetail[]; 995d90ac83dSHong Zhang 996d90ac83dSHong Zhang /*S 9972401956bSBarry Smith MatFactorInfo - Data passed into the matrix factorization routines 998b00f7748SHong Zhang 99961cad860SBarry Smith In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE, that is use 100061cad860SBarry Smith $ MatFactorInfo info(MAT_FACTORINFO_SIZE) 1001b00f7748SHong Zhang 100215e8a5b3SHong Zhang Notes: These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC. 1003b00f7748SHong Zhang 100461cad860SBarry Smith You can use MatFactorInfoInitialize() to set default values. 100561cad860SBarry Smith 1006b00f7748SHong Zhang Level: developer 1007b00f7748SHong Zhang 1008d7d82daaSBarry Smith .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(), 1009d7d82daaSBarry Smith MatFactorInfoInitialize() 1010b00f7748SHong Zhang 1011b00f7748SHong Zhang S*/ 1012b00f7748SHong Zhang typedef struct { 101315e8a5b3SHong Zhang PetscReal diagonal_fill; /* force diagonal to fill in if initially not filled */ 101485317021SBarry Smith PetscReal usedt; 101515e8a5b3SHong Zhang PetscReal dt; /* drop tolerance */ 1016b00f7748SHong Zhang PetscReal dtcol; /* tolerance for pivoting */ 101715e8a5b3SHong Zhang PetscReal dtcount; /* maximum nonzeros to be allowed per row */ 101867e28bfeSBarry Smith PetscReal fill; /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */ 1019348344bbSBarry Smith PetscReal levels; /* ICC/ILU(levels) */ 1020bcd9e38bSBarry Smith PetscReal pivotinblocks; /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0 1021bcd9e38bSBarry Smith factorization may be faster if do not pivot */ 102215e8a5b3SHong Zhang PetscReal zeropivot; /* pivot is called zero if less than this */ 1023f4db908eSBarry Smith PetscReal shifttype; /* type of shift added to matrix factor to prevent zero pivots */ 1024f4db908eSBarry Smith PetscReal shiftamount; /* how large the shift is */ 102515e8a5b3SHong Zhang } MatFactorInfo; 1026ffa6d0a5SLois Curfman McInnes 1027014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFactorInfoInitialize(MatFactorInfo*); 1028014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCholeskyFactor(Mat,IS,const MatFactorInfo*); 1029014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCholeskyFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*); 1030014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCholeskyFactorNumeric(Mat,Mat,const MatFactorInfo*); 1031014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLUFactor(Mat,IS,IS,const MatFactorInfo*); 1032014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatILUFactor(Mat,IS,IS,const MatFactorInfo*); 1033014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*); 1034014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatILUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*); 1035014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatICCFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*); 1036014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatICCFactor(Mat,IS,const MatFactorInfo*); 1037014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLUFactorNumeric(Mat,Mat,const MatFactorInfo*); 1038014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*); 1039014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolve(Mat,Vec,Vec); 1040014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatForwardSolve(Mat,Vec,Vec); 1041014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatBackwardSolve(Mat,Vec,Vec); 1042014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveAdd(Mat,Vec,Vec,Vec); 1043014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveTranspose(Mat,Vec,Vec); 1044014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveTransposeAdd(Mat,Vec,Vec,Vec); 1045014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolves(Mat,Vecs,Vecs); 10468ed539a5SBarry Smith 1047014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetUnfactored(Mat); 1048bb5a7306SBarry Smith 1049d91e6319SBarry Smith /*E 1050d91e6319SBarry Smith MatSORType - What type of (S)SOR to perform 1051bb1eb677SSatish Balay 1052d91e6319SBarry Smith Level: beginner 1053d91e6319SBarry Smith 1054d9274352SBarry Smith May be bitwise ORd together 1055d9274352SBarry Smith 1056d91e6319SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 1057d91e6319SBarry Smith 10584e7234bfSBarry Smith MatSORType may be bitwise ORd together, so do not change the numbers 10594e7234bfSBarry Smith 106041f059aeSBarry Smith .seealso: MatSOR() 1061d91e6319SBarry Smith E*/ 1062ee50ffe9SBarry Smith typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3, 1063ee50ffe9SBarry Smith SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8, 1064ee50ffe9SBarry Smith SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16, 106584cb2905SBarry Smith SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType; 1066014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSOR(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec); 10678ed539a5SBarry Smith 1068d4fbbf0eSBarry Smith /* 1069639f9d9dSBarry Smith These routines are for efficiently computing Jacobians via finite differences. 1070639f9d9dSBarry Smith */ 1071b12f92e5SBarry Smith 10727086a01eSPeter Brune /*S 10737086a01eSPeter Brune MatColoring - Object for managing the coloring of matrices. 10747086a01eSPeter Brune 10757086a01eSPeter Brune Level: beginner 10767086a01eSPeter Brune 10777086a01eSPeter Brune Concepts: matrix, coloring 10787086a01eSPeter Brune 10797086a01eSPeter Brune .seealso: MatFDColoringCreate() ISColoring MatFDColoring 10807086a01eSPeter Brune S*/ 10817086a01eSPeter Brune typedef struct _p_MatColoring* MatColoring; 108276bdecfbSBarry Smith /*J 10838f6c3df8SBarry Smith MatColoringType - String with the name of a PETSc matrix coloring 1084d9274352SBarry Smith 1085d9274352SBarry Smith Level: beginner 1086d9274352SBarry Smith 10877086a01eSPeter Brune .seealso: MatColoringSetType(), MatColoring 108876bdecfbSBarry Smith J*/ 10897086a01eSPeter Brune 109019fd82e9SBarry Smith typedef const char* MatColoringType; 10915567d71aSPeter Brune #define MATCOLORINGJP "jp" 1092b9acb617SPeter Brune #define MATCOLORINGPOWER "power" 10932692d6eeSBarry Smith #define MATCOLORINGNATURAL "natural" 10942692d6eeSBarry Smith #define MATCOLORINGSL "sl" 10952692d6eeSBarry Smith #define MATCOLORINGLF "lf" 10962692d6eeSBarry Smith #define MATCOLORINGID "id" 10978121bdceSPeter Brune #define MATCOLORINGGREEDY "greedy" 1098b12f92e5SBarry Smith 10998ac6da0aSPeter Brune /*E 11008ac6da0aSPeter Brune MatColoringWeightType - Type of weight scheme 11018ac6da0aSPeter Brune 11028ac6da0aSPeter Brune Not Collective 11038ac6da0aSPeter Brune 11048ac6da0aSPeter Brune + MAT_COLORING_RANDOM - Random weights 11058ac6da0aSPeter Brune . MAT_COLORING_LEXICAL - Lexical weighting based upon global numbering. 11068ac6da0aSPeter Brune - MAT_COLORING_LF - Last-first weighting. 11078ac6da0aSPeter Brune 11088ac6da0aSPeter Brune Level: intermediate 11098ac6da0aSPeter Brune 11108ac6da0aSPeter Brune Any additions/changes here MUST also be made in include/finclude/petscmat.h 11118ac6da0aSPeter Brune 11128ac6da0aSPeter Brune .seealso: MatCUSPSetFormat(), MatCUSPFormatOperation 11138ac6da0aSPeter Brune E*/ 1114301c140bSPeter Brune typedef enum {MAT_COLORING_WEIGHT_RANDOM,MAT_COLORING_WEIGHT_LEXICAL,MAT_COLORING_WEIGHT_LF,MAT_COLORING_WEIGHT_SL} MatColoringWeightType; 11158ac6da0aSPeter Brune 1116335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringCreate(Mat,MatColoring*); 1117d7f2a307SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringGetDegrees(Mat,PetscInt,PetscInt*); 1118335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringDestroy(MatColoring*); 1119335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringView(MatColoring,PetscViewer); 1120335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetType(MatColoring,MatColoringType); 1121335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetFromOptions(MatColoring); 1122335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetDistance(MatColoring,PetscInt); 1123335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringGetDistance(MatColoring,PetscInt*); 1124335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetMaxColors(MatColoring,PetscInt); 1125335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringGetMaxColors(MatColoring,PetscInt*); 1126335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringApply(MatColoring,ISColoring*); 1127607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatColoringRegisterAll(void); 1128335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringRegister(const char[],PetscErrorCode(*)(MatColoring)); 1129335efc43SPeter Brune PETSC_EXTERN PetscBool MatColoringRegisterAllCalled; 1130014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*); 11318ac6da0aSPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetWeightType(MatColoring,MatColoringWeightType); 1132c29971fcSPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetWeights(MatColoring,PetscReal*,PetscInt*); 11338ac6da0aSPeter Brune PETSC_EXTERN PetscErrorCode MatColoringCreateWeights(MatColoring,PetscReal **,PetscInt **lperm); 1134639f9d9dSBarry Smith 1135d9274352SBarry Smith /*S 1136d9274352SBarry Smith MatFDColoring - Object for computing a sparse Jacobian via finite differences 1137d9274352SBarry Smith and coloring 1138639f9d9dSBarry Smith 1139d9274352SBarry Smith Level: beginner 1140d9274352SBarry Smith 1141d9274352SBarry Smith Concepts: coloring, sparse Jacobian, finite differences 1142d9274352SBarry Smith 1143d9274352SBarry Smith .seealso: MatFDColoringCreate() 1144d9274352SBarry Smith S*/ 1145e2a1c21fSSatish Balay typedef struct _p_MatFDColoring* MatFDColoring; 1146639f9d9dSBarry Smith 1147014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringCreate(Mat,ISColoring,MatFDColoring *); 1148014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringDestroy(MatFDColoring*); 1149014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringView(MatFDColoring,PetscViewer); 1150014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*); 1151014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**); 1152014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal); 1153014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetFromOptions(MatFDColoring); 1154d1e9a80fSBarry Smith PETSC_EXTERN PetscErrorCode MatFDColoringApply(Mat,MatFDColoring,Vec,void *); 1155014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetF(MatFDColoring,Vec); 1156014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,PetscInt*[]); 1157f86b9fbaSHong Zhang PETSC_EXTERN PetscErrorCode MatFDColoringSetUp(Mat,ISColoring,MatFDColoring); 1158f86b9fbaSHong Zhang PETSC_EXTERN PetscErrorCode MatFDColoringSetBlockSize(MatFDColoring,PetscInt,PetscInt); 1159f86b9fbaSHong Zhang 1160b1683b59SHong Zhang 1161b1683b59SHong Zhang /*S 1162b9af6bddSHong Zhang MatTransposeColoring - Object for computing a sparse matrix product C=A*B^T via coloring 1163b1683b59SHong Zhang 1164b1683b59SHong Zhang Level: beginner 1165b1683b59SHong Zhang 1166b1683b59SHong Zhang Concepts: coloring, sparse matrix product 1167b1683b59SHong Zhang 1168b9af6bddSHong Zhang .seealso: MatTransposeColoringCreate() 1169b1683b59SHong Zhang S*/ 1170b9af6bddSHong Zhang typedef struct _p_MatTransposeColoring* MatTransposeColoring; 1171b1683b59SHong Zhang 1172014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeColoringCreate(Mat,ISColoring,MatTransposeColoring *); 1173014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransColoringApplySpToDen(MatTransposeColoring,Mat,Mat); 1174014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransColoringApplyDenToSp(MatTransposeColoring,Mat,Mat); 1175014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeColoringDestroy(MatTransposeColoring*); 1176b1683b59SHong Zhang 1177639f9d9dSBarry Smith /* 11780752156aSBarry Smith These routines are for partitioning matrices: currently used only 11793eda8832SBarry Smith for adjacency matrix, MatCreateMPIAdj(). 11800752156aSBarry Smith */ 1181ca161407SBarry Smith 1182d9274352SBarry Smith /*S 1183d9274352SBarry Smith MatPartitioning - Object for managing the partitioning of a matrix or graph 1184d9274352SBarry Smith 1185d9274352SBarry Smith Level: beginner 1186d9274352SBarry Smith 1187d9274352SBarry Smith Concepts: partitioning 1188d9274352SBarry Smith 1189743a1026Svictorle .seealso: MatPartitioningCreate(), MatPartitioningType 1190d9274352SBarry Smith S*/ 119191e9ee9fSBarry Smith typedef struct _p_MatPartitioning* MatPartitioning; 1192d9274352SBarry Smith 119376bdecfbSBarry Smith /*J 11948f6c3df8SBarry Smith MatPartitioningType - String with the name of a PETSc matrix partitioning 1195d9274352SBarry Smith 1196d9274352SBarry Smith Level: beginner 11970d04baf8SBarry Smith dm 1198b547a13bSvictorle .seealso: MatPartitioningCreate(), MatPartitioning 119976bdecfbSBarry Smith J*/ 120019fd82e9SBarry Smith typedef const char* MatPartitioningType; 12012692d6eeSBarry Smith #define MATPARTITIONINGCURRENT "current" 12022692d6eeSBarry Smith #define MATPARTITIONINGSQUARE "square" 12032692d6eeSBarry Smith #define MATPARTITIONINGPARMETIS "parmetis" 12042692d6eeSBarry Smith #define MATPARTITIONINGCHACO "chaco" 12052692d6eeSBarry Smith #define MATPARTITIONINGPARTY "party" 120650d91057SBarry Smith #define MATPARTITIONINGPTSCOTCH "ptscotch" 120771306c60Sjroman 1208ca161407SBarry Smith 1209014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningCreate(MPI_Comm,MatPartitioning*); 121019fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningSetType(MatPartitioning,MatPartitioningType); 1211014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetNParts(MatPartitioning,PetscInt); 1212014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetAdjacency(MatPartitioning,Mat); 1213014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]); 1214014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []); 1215014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningApply(MatPartitioning,IS*); 1216014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningDestroy(MatPartitioning*); 12172aabb6bbSBarry Smith 1218bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningRegister(const char[],PetscErrorCode (*)(MatPartitioning)); 121930de9b25SBarry Smith 1220014dd563SJed Brown PETSC_EXTERN PetscBool MatPartitioningRegisterAllCalled; 1221f1144a30SSatish Balay 1222607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningRegisterAll(void); 12232bad1931SBarry Smith 1224014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningView(MatPartitioning,PetscViewer); 1225014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetFromOptions(MatPartitioning); 122619fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningGetType(MatPartitioning,MatPartitioningType*); 1227ca161407SBarry Smith 1228014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetCoarseSequential(MatPartitioning); 1229014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningParmetisGetEdgeCut(MatPartitioning, PetscInt *); 12300752156aSBarry Smith 1231b6956c12SJose Roman typedef enum { MP_CHACO_MULTILEVEL=1,MP_CHACO_SPECTRAL=2,MP_CHACO_LINEAR=4,MP_CHACO_RANDOM=5,MP_CHACO_SCATTERED=6 } MPChacoGlobalType; 12326a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoGlobalTypes[]; 1233b6956c12SJose Roman typedef enum { MP_CHACO_KERNIGHAN=1,MP_CHACO_NONE=2 } MPChacoLocalType; 12346a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoLocalTypes[]; 1235b6956c12SJose Roman typedef enum { MP_CHACO_LANCZOS=0,MP_CHACO_RQI=1 } MPChacoEigenType; 12366a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoEigenTypes[]; 1237b6956c12SJose Roman 1238014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetGlobal(MatPartitioning,MPChacoGlobalType); 1239014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetGlobal(MatPartitioning,MPChacoGlobalType*); 1240014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetLocal(MatPartitioning,MPChacoLocalType); 1241014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetLocal(MatPartitioning,MPChacoLocalType*); 1242014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal); 1243014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType); 1244014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenSolver(MatPartitioning,MPChacoEigenType*); 1245014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenTol(MatPartitioning,PetscReal); 1246014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenTol(MatPartitioning,PetscReal*); 1247014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenNumber(MatPartitioning,PetscInt); 1248014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenNumber(MatPartitioning,PetscInt*); 124971306c60Sjroman 125071306c60Sjroman #define MP_PARTY_OPT "opt" 125171306c60Sjroman #define MP_PARTY_LIN "lin" 125271306c60Sjroman #define MP_PARTY_SCA "sca" 125371306c60Sjroman #define MP_PARTY_RAN "ran" 125471306c60Sjroman #define MP_PARTY_GBF "gbf" 125571306c60Sjroman #define MP_PARTY_GCF "gcf" 125671306c60Sjroman #define MP_PARTY_BUB "bub" 125771306c60Sjroman #define MP_PARTY_DEF "def" 1258014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetGlobal(MatPartitioning,const char*); 125971306c60Sjroman #define MP_PARTY_HELPFUL_SETS "hs" 126071306c60Sjroman #define MP_PARTY_KERNIGHAN_LIN "kl" 126171306c60Sjroman #define MP_PARTY_NONE "no" 1262014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetLocal(MatPartitioning,const char*); 1263014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal); 1264014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetBipart(MatPartitioning,PetscBool); 1265014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscBool); 126671306c60Sjroman 126750d91057SBarry Smith typedef enum { MP_PTSCOTCH_QUALITY,MP_PTSCOTCH_SPEED,MP_PTSCOTCH_BALANCE,MP_PTSCOTCH_SAFETY,MP_PTSCOTCH_SCALABILITY } MPPTScotchStrategyType; 12686a6fc655SJed Brown PETSC_EXTERN const char *const MPPTScotchStrategyTypes[]; 1269e0f1cffaSJose Roman 1270014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetImbalance(MatPartitioning,PetscReal); 1271014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetImbalance(MatPartitioning,PetscReal*); 1272014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetStrategy(MatPartitioning,MPPTScotchStrategyType); 1273014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetStrategy(MatPartitioning,MPPTScotchStrategyType*); 127471306c60Sjroman 1275b43b03e9SMark F. Adams /* 1276b43b03e9SMark F. Adams These routines are for coarsening matrices: 1277b43b03e9SMark F. Adams */ 1278b43b03e9SMark F. Adams 1279b43b03e9SMark F. Adams /*S 1280b43b03e9SMark F. Adams MatCoarsen - Object for managing the coarsening of a graph (symmetric matrix) 1281b43b03e9SMark F. Adams 1282b43b03e9SMark F. Adams Level: beginner 1283b43b03e9SMark F. Adams 1284b43b03e9SMark F. Adams Concepts: coarsen 1285b43b03e9SMark F. Adams 1286b43b03e9SMark F. Adams .seealso: MatCoarsenCreate), MatCoarsenType 1287b43b03e9SMark F. Adams S*/ 1288b43b03e9SMark F. Adams typedef struct _p_MatCoarsen* MatCoarsen; 1289b43b03e9SMark F. Adams 1290b43b03e9SMark F. Adams /*J 12918f6c3df8SBarry Smith MatCoarsenType - String with the name of a PETSc matrix coarsen 1292b43b03e9SMark F. Adams 1293b43b03e9SMark F. Adams Level: beginner 12948f6c3df8SBarry Smith 1295b43b03e9SMark F. Adams .seealso: MatCoarsenCreate(), MatCoarsen 1296b43b03e9SMark F. Adams J*/ 129719fd82e9SBarry Smith typedef const char* MatCoarsenType; 1298b43b03e9SMark F. Adams #define MATCOARSENMIS "mis" 12999057884aSMark F. Adams #define MATCOARSENHEM "hem" 1300b43b03e9SMark F. Adams 13010cbbd2e1SMark F. Adams /* linked list for aggregates */ 130241b27cdeSMark F. Adams typedef struct _PetscCDIntNd{ 130341b27cdeSMark F. Adams struct _PetscCDIntNd *next; 13040cbbd2e1SMark F. Adams PetscInt gid; 130541b27cdeSMark F. Adams }PetscCDIntNd; 13060cbbd2e1SMark F. Adams 13070cbbd2e1SMark F. Adams /* only used by node pool */ 130841b27cdeSMark F. Adams typedef struct _PetscCDArrNd{ 130941b27cdeSMark F. Adams struct _PetscCDArrNd *next; 131041b27cdeSMark F. Adams struct _PetscCDIntNd *array; 131141b27cdeSMark F. Adams }PetscCDArrNd; 13120cbbd2e1SMark F. Adams 13130cbbd2e1SMark F. Adams typedef struct _PetscCoarsenData{ 131482c86c8fSBarry Smith PetscCDArrNd pool_list; /* node pool */ 131541b27cdeSMark F. Adams PetscCDIntNd *new_node; 13160cbbd2e1SMark F. Adams PetscInt new_left; 13170cbbd2e1SMark F. Adams PetscInt chk_sz; 131841b27cdeSMark F. Adams PetscCDIntNd *extra_nodes; 131982c86c8fSBarry Smith PetscCDIntNd **array; /* Array of lists */ 13200cbbd2e1SMark F. Adams PetscInt size; 132182c86c8fSBarry Smith Mat mat; /* cache a Mat for communication data */ 13220cbbd2e1SMark F. Adams }PetscCoarsenData; 13230cbbd2e1SMark F. Adams 1324014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenCreate(MPI_Comm,MatCoarsen*); 132519fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenSetType(MatCoarsen,MatCoarsenType); 1326014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetAdjacency(MatCoarsen,Mat); 1327014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetGreedyOrdering(MatCoarsen,const IS); 1328014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetStrictAggs(MatCoarsen,PetscBool); 1329014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetVerbose(MatCoarsen,PetscInt); 1330014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenGetData( MatCoarsen, PetscCoarsenData ** ); 1331014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenApply(MatCoarsen); 1332014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenDestroy(MatCoarsen*); 1333b43b03e9SMark F. Adams 1334bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenRegister(const char[],PetscErrorCode (*)(MatCoarsen)); 1335b43b03e9SMark F. Adams 1336014dd563SJed Brown PETSC_EXTERN PetscBool MatCoarsenRegisterAllCalled; 1337b43b03e9SMark F. Adams 1338607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenRegisterAll(void); 1339b43b03e9SMark F. Adams 1340014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenView(MatCoarsen,PetscViewer); 1341014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetFromOptions(MatCoarsen); 134219fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenGetType(MatCoarsen,MatCoarsenType*); 1343ce1779c8SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatCoarsenViewFromOptions(MatCoarsen A,const char prefix[],const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,prefix,name);} 1344b43b03e9SMark F. Adams 1345014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMeshToVertexGraph(Mat,PetscInt,Mat*); 1346014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMeshToCellGraph(Mat,PetscInt,Mat*); 1347591294e4SBarry Smith 13480752156aSBarry Smith /* 13490a835dfdSSatish Balay If you add entries here you must also add them to finclude/petscmat.h 1350d4fbbf0eSBarry Smith */ 13511c1c02c0SLois Curfman McInnes typedef enum { MATOP_SET_VALUES=0, 13521c1c02c0SLois Curfman McInnes MATOP_GET_ROW=1, 13531c1c02c0SLois Curfman McInnes MATOP_RESTORE_ROW=2, 13541c1c02c0SLois Curfman McInnes MATOP_MULT=3, 13551c1c02c0SLois Curfman McInnes MATOP_MULT_ADD=4, 13567c922b88SBarry Smith MATOP_MULT_TRANSPOSE=5, 13577c922b88SBarry Smith MATOP_MULT_TRANSPOSE_ADD=6, 13581c1c02c0SLois Curfman McInnes MATOP_SOLVE=7, 13591c1c02c0SLois Curfman McInnes MATOP_SOLVE_ADD=8, 13607c922b88SBarry Smith MATOP_SOLVE_TRANSPOSE=9, 13617c922b88SBarry Smith MATOP_SOLVE_TRANSPOSE_ADD=10, 13621c1c02c0SLois Curfman McInnes MATOP_LUFACTOR=11, 13631c1c02c0SLois Curfman McInnes MATOP_CHOLESKYFACTOR=12, 1364a714c06dSBarry Smith MATOP_SOR=13, 13651c1c02c0SLois Curfman McInnes MATOP_TRANSPOSE=14, 13661c1c02c0SLois Curfman McInnes MATOP_GETINFO=15, 13671c1c02c0SLois Curfman McInnes MATOP_EQUAL=16, 13681c1c02c0SLois Curfman McInnes MATOP_GET_DIAGONAL=17, 13691c1c02c0SLois Curfman McInnes MATOP_DIAGONAL_SCALE=18, 13701c1c02c0SLois Curfman McInnes MATOP_NORM=19, 13711c1c02c0SLois Curfman McInnes MATOP_ASSEMBLY_BEGIN=20, 13721c1c02c0SLois Curfman McInnes MATOP_ASSEMBLY_END=21, 1373d519adbfSMatthew Knepley MATOP_SET_OPTION=22, 1374d519adbfSMatthew Knepley MATOP_ZERO_ENTRIES=23, 1375d519adbfSMatthew Knepley MATOP_ZERO_ROWS=24, 1376d519adbfSMatthew Knepley MATOP_LUFACTOR_SYMBOLIC=25, 1377d519adbfSMatthew Knepley MATOP_LUFACTOR_NUMERIC=26, 1378d519adbfSMatthew Knepley MATOP_CHOLESKY_FACTOR_SYMBOLIC=27, 1379d519adbfSMatthew Knepley MATOP_CHOLESKY_FACTOR_NUMERIC=28, 1380d519adbfSMatthew Knepley MATOP_SETUP_PREALLOCATION=29, 1381d519adbfSMatthew Knepley MATOP_ILUFACTOR_SYMBOLIC=30, 1382d519adbfSMatthew Knepley MATOP_ICCFACTOR_SYMBOLIC=31, 13839d39f69dSJed Brown /* MATOP_PLACEHOLDER_32=32, */ 13849d39f69dSJed Brown /* MATOP_PLACEHOLDER_33=33, */ 1385d519adbfSMatthew Knepley MATOP_DUPLICATE=34, 1386d519adbfSMatthew Knepley MATOP_FORWARD_SOLVE=35, 1387d519adbfSMatthew Knepley MATOP_BACKWARD_SOLVE=36, 1388d519adbfSMatthew Knepley MATOP_ILUFACTOR=37, 1389d519adbfSMatthew Knepley MATOP_ICCFACTOR=38, 1390d519adbfSMatthew Knepley MATOP_AXPY=39, 1391d519adbfSMatthew Knepley MATOP_GET_SUBMATRICES=40, 1392d519adbfSMatthew Knepley MATOP_INCREASE_OVERLAP=41, 1393d519adbfSMatthew Knepley MATOP_GET_VALUES=42, 1394d519adbfSMatthew Knepley MATOP_COPY=43, 1395d519adbfSMatthew Knepley MATOP_GET_ROW_MAX=44, 1396d519adbfSMatthew Knepley MATOP_SCALE=45, 1397d519adbfSMatthew Knepley MATOP_SHIFT=46, 139835153367SBarry Smith MATOP_DIAGONAL_SET=47, 13999d39f69dSJed Brown MATOP_ZERO_ROWS_COLUMNS=48, 14009d39f69dSJed Brown MATOP_SET_RANDOM=49, 1401d519adbfSMatthew Knepley MATOP_GET_ROW_IJ=50, 1402d519adbfSMatthew Knepley MATOP_RESTORE_ROW_IJ=51, 1403d519adbfSMatthew Knepley MATOP_GET_COLUMN_IJ=52, 1404d519adbfSMatthew Knepley MATOP_RESTORE_COLUMN_IJ=53, 1405d519adbfSMatthew Knepley MATOP_FDCOLORING_CREATE=54, 1406d519adbfSMatthew Knepley MATOP_COLORING_PATCH=55, 1407d519adbfSMatthew Knepley MATOP_SET_UNFACTORED=56, 1408d519adbfSMatthew Knepley MATOP_PERMUTE=57, 1409d519adbfSMatthew Knepley MATOP_SET_VALUES_BLOCKED=58, 1410d519adbfSMatthew Knepley MATOP_GET_SUBMATRIX=59, 1411d519adbfSMatthew Knepley MATOP_DESTROY=60, 1412d519adbfSMatthew Knepley MATOP_VIEW=61, 1413d519adbfSMatthew Knepley MATOP_CONVERT_FROM=62, 14147bab7c10SHong Zhang MATOP_MATMAT_MULT=63, 14157bab7c10SHong Zhang MATOP_MATMAT_MULT_SYMBOLIC=64, 14167bab7c10SHong Zhang MATOP_MATMAT_MULT_NUMERIC=65, 1417d519adbfSMatthew Knepley MATOP_SET_LOCAL_TO_GLOBAL_MAP=66, 1418d519adbfSMatthew Knepley MATOP_SET_VALUES_LOCAL=67, 1419d519adbfSMatthew Knepley MATOP_ZERO_ROWS_LOCAL=68, 1420d519adbfSMatthew Knepley MATOP_GET_ROW_MAX_ABS=69, 1421d519adbfSMatthew Knepley MATOP_GET_ROW_MIN_ABS=70, 1422d519adbfSMatthew Knepley MATOP_CONVERT=71, 1423d519adbfSMatthew Knepley MATOP_SET_COLORING=72, 14249d39f69dSJed Brown /* MATOP_PLACEHOLDER_73=73, */ 1425d519adbfSMatthew Knepley MATOP_SET_VALUES_ADIFOR=74, 1426d519adbfSMatthew Knepley MATOP_FD_COLORING_APPLY=75, 1427d519adbfSMatthew Knepley MATOP_SET_FROM_OPTIONS=76, 1428bda6c4cbSBarry Smith MATOP_MULT_CONSTRAINED=77, 1429bda6c4cbSBarry Smith MATOP_MULT_TRANSPOSE_CONSTRAIN=78, 14309d39f69dSJed Brown MATOP_FIND_ZERO_DIAGONALS=79, 1431d519adbfSMatthew Knepley MATOP_MULT_MULTIPLE=80, 1432d519adbfSMatthew Knepley MATOP_SOLVE_MULTIPLE=81, 1433d519adbfSMatthew Knepley MATOP_GET_INERTIA=82, 1434d519adbfSMatthew Knepley MATOP_LOAD=83, 1435d519adbfSMatthew Knepley MATOP_IS_SYMMETRIC=84, 1436d519adbfSMatthew Knepley MATOP_IS_HERMITIAN=85, 1437d519adbfSMatthew Knepley MATOP_IS_STRUCTURALLY_SYMMETRIC=86, 1438820469bcSHong Zhang MATOP_SET_VALUES_BLOCKEDLOCAL=87, 1439d519adbfSMatthew Knepley MATOP_GET_VECS=88, 1440d519adbfSMatthew Knepley MATOP_MAT_MULT=89, 1441d519adbfSMatthew Knepley MATOP_MAT_MULT_SYMBOLIC=90, 1442d519adbfSMatthew Knepley MATOP_MAT_MULT_NUMERIC=91, 1443d519adbfSMatthew Knepley MATOP_PTAP=92, 1444d519adbfSMatthew Knepley MATOP_PTAP_SYMBOLIC=93, 1445d519adbfSMatthew Knepley MATOP_PTAP_NUMERIC=94, 1446bda6c4cbSBarry Smith MATOP_MAT_TRANSPOSE_MULT=95, 1447bda6c4cbSBarry Smith MATOP_MAT_TRANSPOSE_MULT_SYMBO=96, 1448bda6c4cbSBarry Smith MATOP_MAT_TRANSPOSE_MULT_NUMER=97, 14499d39f69dSJed Brown /* MATOP_PLACEHOLDER_98=98, */ 14509d39f69dSJed Brown /* MATOP_PLACEHOLDER_99=99, */ 14519d39f69dSJed Brown /* MATOP_PLACEHOLDER_100=100, */ 14529d39f69dSJed Brown /* MATOP_PLACEHOLDER_101=101, */ 1453d519adbfSMatthew Knepley MATOP_CONJUGATE=102, 14549d39f69dSJed Brown /* MATOP_PLACEHOLDER_103=103, */ 1455d519adbfSMatthew Knepley MATOP_SET_VALUES_ROW=104, 1456d519adbfSMatthew Knepley MATOP_REAL_PART=105, 1457bda6c4cbSBarry Smith MATOP_IMAGINARY_PART=106, 1458bda6c4cbSBarry Smith MATOP_GET_ROW_UPPER_TRIANGULAR=107, 1459bda6c4cbSBarry Smith MATOP_RESTORE_ROW_UPPER_TRIANG=108, 1460bda6c4cbSBarry Smith MATOP_MAT_SOLVE=109, 14610e8d04f7SBarry Smith MATOP_GET_REDUNDANT_MATRIX=110, 1462d519adbfSMatthew Knepley MATOP_GET_ROW_MIN=111, 14630e8d04f7SBarry Smith MATOP_GET_COLUMN_VECTOR=112, 1464d519adbfSMatthew Knepley MATOP_MISSING_DIAGONAL=113, 14650e8d04f7SBarry Smith MATOP_GET_SEQ_NONZERO_STRUCTUR=114, 146689c6957cSBarry Smith MATOP_CREATE=115, 146789c6957cSBarry Smith MATOP_GET_GHOSTS=116, 14680e8d04f7SBarry Smith MATOP_GET_LOCAL_SUB_MATRIX=117, 14690e8d04f7SBarry Smith MATOP_RESTORE_LOCALSUB_MATRIX=118, 1470eeffb40dSHong Zhang MATOP_MULT_DIAGONAL_BLOCK=119, 14710e8d04f7SBarry Smith MATOP_HERMITIAN_TRANSPOSE=120, 14720e8d04f7SBarry Smith MATOP_MULT_HERMITIAN_TRANSPOSE=121, 14730e8d04f7SBarry Smith MATOP_MULT_HERMITIAN_TRANS_ADD=122, 14740e8d04f7SBarry Smith MATOP_GET_MULTI_PROC_BLOCK=123, 14759d39f69dSJed Brown MATOP_FIND_NONZERO_ROWS=124, 14760e8d04f7SBarry Smith MATOP_GET_COLUMN_NORMS=125, 14779d39f69dSJed Brown MATOP_INVERT_BLOCK_DIAGONAL=126, 14789d39f69dSJed Brown /* MATOP_PLACEHOLDER_127=127, */ 14790e8d04f7SBarry Smith MATOP_GET_SUB_MATRICES_PARALLE=128, 14802b8ad9a3SHong Zhang MATOP_SET_VALUES_BATCH=129, 14810e8d04f7SBarry Smith MATOP_TRANSPOSE_MAT_MULT=130, 14820e8d04f7SBarry Smith MATOP_TRANSPOSE_MAT_MULT_SYMBO=131, 14830e8d04f7SBarry Smith MATOP_TRANSPOSE_MAT_MULT_NUMER=132, 14840e8d04f7SBarry Smith MATOP_TRANSPOSE_COLORING_CREAT=133, 14850e8d04f7SBarry Smith MATOP_TRANS_COLORING_APPLY_SPT=134, 14860e8d04f7SBarry Smith MATOP_TRANS_COLORING_APPLY_DEN=135, 14870e8d04f7SBarry Smith MATOP_RART=136, 14880e8d04f7SBarry Smith MATOP_RART_SYMBOLIC=137, 14890e8d04f7SBarry Smith MATOP_RART_NUMERIC=138, 1490e09a3074SHong Zhang MATOP_SET_BLOCK_SIZES=139, 1491f9426fe0SMark Adams MATOP_AYPX=140, 14921919a2e2SJed Brown MATOP_RESIDUAL=141, 14931919a2e2SJed Brown MATOP_FDCOLORING_SETUP= 142 1494fae171e0SBarry Smith } MatOperation; 1495014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatHasOperation(Mat,MatOperation,PetscBool *); 1496014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellSetOperation(Mat,MatOperation,void(*)(void)); 1497014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellGetOperation(Mat,MatOperation,void(**)(void)); 1498014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellSetContext(Mat,void*); 1499112a2221SBarry Smith 150090ace30eSBarry Smith /* 150190ace30eSBarry Smith Codes for matrices stored on disk. By default they are 150290ace30eSBarry Smith stored in a universal format. By changing the format with 15037973ad04SBarry Smith PetscViewerSetFormat(viewer,PETSC_VIEWER_NATIVE); the matrices will 150490ace30eSBarry Smith be stored in a way natural for the matrix, for example dense matrices 150590ace30eSBarry Smith would be stored as dense. Matrices stored this way may only be 1506f4403165SShri Abhyankar read into matrices of the same type. 150790ace30eSBarry Smith */ 150890ace30eSBarry Smith #define MATRIX_BINARY_FORMAT_DENSE -1 150990ace30eSBarry Smith 1510014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetHashTableFactor(Mat,PetscReal); 1511014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatISGetLocalMat(Mat,Mat*); 1512014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatISSetLocalMat(Mat,Mat); 1513b7ce53b6SStefano Zampini PETSC_EXTERN PetscErrorCode MatISGetMPIXAIJ(Mat,MatReuse,Mat*); 15141f4e1ec7SBarry Smith 1515d9274352SBarry Smith /*S 1516d9274352SBarry Smith MatNullSpace - Object that removes a null space from a vector, i.e. 1517d9274352SBarry Smith orthogonalizes the vector to a subsapce 1518d9274352SBarry Smith 1519f7a9e4ceSBarry Smith Level: advanced 1520d9274352SBarry Smith 1521d9274352SBarry Smith Concepts: matrix; linear operator, null space 1522d9274352SBarry Smith 15236e1639daSBarry Smith Users manual sections: 15246e1639daSBarry Smith . sec_singular 15256e1639daSBarry Smith 1526d9274352SBarry Smith .seealso: MatNullSpaceCreate() 1527d9274352SBarry Smith S*/ 152874637425SBarry Smith typedef struct _p_MatNullSpace* MatNullSpace; 1529d9274352SBarry Smith 1530014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceCreate(MPI_Comm,PetscBool ,PetscInt,const Vec[],MatNullSpace*); 1531014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(MatNullSpace,Vec,void*),void*); 1532014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceDestroy(MatNullSpace*); 1533d0195637SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceRemove(MatNullSpace,Vec); 1534014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNullSpace(Mat, MatNullSpace *); 1535014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNullSpace(Mat,MatNullSpace); 1536014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNearNullSpace(Mat,MatNullSpace); 1537014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNearNullSpace(Mat,MatNullSpace*); 1538014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceTest(MatNullSpace,Mat,PetscBool *); 1539014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceView(MatNullSpace,PetscViewer); 1540014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceGetVecs(MatNullSpace,PetscBool*,PetscInt*,const Vec**); 1541014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceCreateRigidBody(Vec,MatNullSpace*); 154274637425SBarry Smith 1543014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatReorderingSeqSBAIJ(Mat,IS); 1544014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetHashTableFactor(Mat,PetscReal); 1545014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetColumnIndices(Mat,PetscInt *); 1546014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJInvertBlockDiagonal(Mat); 15473f1d51d7SBarry Smith 1548014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMAIJ(Mat,PetscInt,Mat*); 1549014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMAIJRedimension(Mat,PetscInt,Mat*); 1550014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMAIJGetAIJ(Mat,Mat*); 1551c4f061fbSSatish Balay 1552014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatComputeExplicitOperator(Mat,Mat*); 1553b0a32e0cSBarry Smith 1554014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalScaleLocal(Mat,Vec); 155504f1ad80SBarry Smith 1556014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*); 1557014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetBase(Mat,Vec,Vec); 1558014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*); 1559014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*)); 1560014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec)); 1561014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDAddNullSpace(Mat,MatNullSpace); 1562014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt); 1563014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDResetHHistory(Mat); 1564014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctionError(Mat,PetscReal); 1565014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetPeriod(Mat,PetscInt); 1566014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDGetH(Mat,PetscScalar *); 1567014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetOptionsPrefix(Mat,const char[]); 1568014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*); 1569014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*); 1570e884886fSBarry Smith 15716370053bSBarry Smith /*S 15726370053bSBarry Smith MatMFFD - A data structured used to manage the computation of the h differencing parameter for matrix-free 15736370053bSBarry Smith Jacobian vector products 1574e884886fSBarry Smith 15756370053bSBarry Smith Notes: MATMFFD is a specific MatType which uses the MatMFFD data structure 15766370053bSBarry Smith 15776370053bSBarry Smith MatMFFD*() methods actually take the Mat as their first argument. Not a MatMFFD data structure 15786370053bSBarry Smith 15796370053bSBarry Smith Level: developer 15806370053bSBarry Smith 15816370053bSBarry Smith .seealso: MATMFFD, MatCreateMFFD(), MatMFFDSetFuction(), MatMFFDSetType(), MatMFFDRegister() 15826370053bSBarry Smith S*/ 1583e884886fSBarry Smith typedef struct _p_MatMFFD* MatMFFD; 1584e884886fSBarry Smith 158576bdecfbSBarry Smith /*J 1586e884886fSBarry Smith MatMFFDType - algorithm used to compute the h used in computing matrix-vector products via differencing of the function 1587e884886fSBarry Smith 1588e884886fSBarry Smith Level: beginner 1589e884886fSBarry Smith 1590e884886fSBarry Smith .seealso: MatMFFDSetType(), MatMFFDRegister() 159176bdecfbSBarry Smith J*/ 159219fd82e9SBarry Smith typedef const char* MatMFFDType; 1593e884886fSBarry Smith #define MATMFFD_DS "ds" 1594e884886fSBarry Smith #define MATMFFD_WP "wp" 1595e884886fSBarry Smith 159619fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDSetType(Mat,MatMFFDType); 1597bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDRegister(const char[],PetscErrorCode (*)(MatMFFD)); 1598e884886fSBarry Smith 1599607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDRegisterAll(void); 1600014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDDSSetUmin(Mat,PetscReal); 1601014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDWPSetComputeNormU(Mat,PetscBool ); 1602e884886fSBarry Smith 1603014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *); 1604014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *); 16057dbadf16SMatthew Knepley 160697969023SHong Zhang /* 160797969023SHong Zhang PETSc interface to MUMPS 160897969023SHong Zhang */ 160997969023SHong Zhang #ifdef PETSC_HAVE_MUMPS 1610014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMumpsSetIcntl(Mat,PetscInt,PetscInt); 1611bc6112feSHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetIcntl(Mat,PetscInt,PetscInt*); 1612b9e9ea26SSatish Balay PETSC_EXTERN PetscErrorCode MatMumpsSetCntl(Mat,PetscInt,PetscReal); 1613bc6112feSHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetCntl(Mat,PetscInt,PetscReal*); 1614bc6112feSHong Zhang 1615ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetInfo(Mat,PetscInt,PetscInt*); 1616ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetInfog(Mat,PetscInt,PetscInt*); 1617ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetRinfo(Mat,PetscInt,PetscReal*); 1618ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetRinfog(Mat,PetscInt,PetscReal*); 161997969023SHong Zhang #endif 162023a5497aSJed Brown 1621d954db57SHong Zhang /* 1622d615f992SSatish Balay PETSc interface to Mkl_Pardiso 1623b500a6b7SJose David Bermeol */ 1624b500a6b7SJose David Bermeol #ifdef PETSC_HAVE_MKL_PARDISO 1625d615f992SSatish Balay PETSC_EXTERN PetscErrorCode MatMkl_PardisoSetCntl(Mat,PetscInt,PetscInt); 1626b500a6b7SJose David Bermeol #endif 1627b500a6b7SJose David Bermeol 1628b500a6b7SJose David Bermeol /* 1629d954db57SHong Zhang PETSc interface to SUPERLU 1630d954db57SHong Zhang */ 1631d954db57SHong Zhang #ifdef PETSC_HAVE_SUPERLU 1632014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSuperluSetILUDropTol(Mat,PetscReal); 1633d954db57SHong Zhang #endif 1634d954db57SHong Zhang 1635aa372e3fSPaul Mullowney #ifdef PETSC_HAVE_CUDA 16363f9c0db1SPaul Mullowney /*E 1637e057df02SPaul Mullowney MatCUSPARSEStorageFormat - indicates the storage format for CUSPARSE (GPU) 16382692e278SPaul Mullowney matrices. 1639e057df02SPaul Mullowney 1640e057df02SPaul Mullowney Not Collective 1641e057df02SPaul Mullowney 16428468deeeSKarl Rupp + MAT_CUSPARSE_CSR - Compressed Sparse Row 16432692e278SPaul Mullowney . MAT_CUSPARSE_ELL - Ellpack (requires CUDA 4.2 or later). 16442692e278SPaul Mullowney - MAT_CUSPARSE_HYB - Hybrid, a combination of Ellpack and Coordinate format (requires CUDA 4.2 or later). 1645e057df02SPaul Mullowney 1646e057df02SPaul Mullowney Level: intermediate 1647e057df02SPaul Mullowney 1648e057df02SPaul Mullowney Any additions/changes here MUST also be made in include/finclude/petscmat.h 1649e057df02SPaul Mullowney 1650e057df02SPaul Mullowney .seealso: MatCUSPARSESetFormat(), MatCUSPARSEFormatOperation 1651e057df02SPaul Mullowney E*/ 1652e057df02SPaul Mullowney 16533f9c0db1SPaul Mullowney typedef enum {MAT_CUSPARSE_CSR, MAT_CUSPARSE_ELL, MAT_CUSPARSE_HYB} MatCUSPARSEStorageFormat; 1654e057df02SPaul Mullowney 1655e057df02SPaul Mullowney /* these will be strings associated with enumerated type defined above */ 1656e057df02SPaul Mullowney PETSC_EXTERN const char *const MatCUSPARSEStorageFormats[]; 1657e057df02SPaul Mullowney 16583f9c0db1SPaul Mullowney /*E 1659e057df02SPaul Mullowney MatCUSPARSEFormatOperation - indicates the operation of CUSPARSE (GPU) 16602692e278SPaul Mullowney matrices whose operation should use a particular storage format. 1661e057df02SPaul Mullowney 1662e057df02SPaul Mullowney Not Collective 1663e057df02SPaul Mullowney 16648468deeeSKarl Rupp + MAT_CUSPARSE_MULT_DIAG - sets the storage format for the diagonal matrix in the parallel MatMult 16658468deeeSKarl Rupp . MAT_CUSPARSE_MULT_OFFDIAG - sets the storage format for the offdiagonal matrix in the parallel MatMult 16668468deeeSKarl Rupp . MAT_CUSPARSE_MULT - sets the storage format for the entire matrix in the serial (single GPU) MatMult 16678468deeeSKarl Rupp - MAT_CUSPARSE_ALL - sets the storage format for all CUSPARSE (GPU) matrices 1668e057df02SPaul Mullowney 1669e057df02SPaul Mullowney Level: intermediate 1670e057df02SPaul Mullowney 1671e057df02SPaul Mullowney .seealso: MatCUSPARSESetFormat(), MatCUSPARSEStorageFormat 1672e057df02SPaul Mullowney E*/ 167336d62e41SPaul Mullowney typedef enum {MAT_CUSPARSE_MULT_DIAG, MAT_CUSPARSE_MULT_OFFDIAG, MAT_CUSPARSE_MULT, MAT_CUSPARSE_ALL} MatCUSPARSEFormatOperation; 1674e057df02SPaul Mullowney 167510e9db80SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 167610e9db80SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCreateAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 1677e057df02SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCUSPARSESetFormat(Mat,MatCUSPARSEFormatOperation,MatCUSPARSEStorageFormat); 16789ae82921SPaul Mullowney #endif 16799ae82921SPaul Mullowney 168090c53902SBarry Smith #if defined(PETSC_HAVE_CUSP) 1681014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCUSP(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 1682014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateAIJCUSP(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 1683e057df02SPaul Mullowney 16843f9c0db1SPaul Mullowney /*E 1685e057df02SPaul Mullowney MatCUSPStorageFormat - indicates the storage format for CUSP (GPU) 168636d62e41SPaul Mullowney matrices. 1687e057df02SPaul Mullowney 1688e057df02SPaul Mullowney Not Collective 1689e057df02SPaul Mullowney 16908468deeeSKarl Rupp + MAT_CUSP_CSR - Compressed Sparse Row 16918468deeeSKarl Rupp . MAT_CUSP_DIA - Diagonal 16928468deeeSKarl Rupp - MAT_CUSP_ELL - Ellpack 1693e057df02SPaul Mullowney 1694e057df02SPaul Mullowney Level: intermediate 1695e057df02SPaul Mullowney 1696e057df02SPaul Mullowney Any additions/changes here MUST also be made in include/finclude/petscmat.h 1697e057df02SPaul Mullowney 1698e057df02SPaul Mullowney .seealso: MatCUSPSetFormat(), MatCUSPFormatOperation 1699e057df02SPaul Mullowney E*/ 17003f9c0db1SPaul Mullowney typedef enum {MAT_CUSP_CSR, MAT_CUSP_DIA, MAT_CUSP_ELL} MatCUSPStorageFormat; 1701e057df02SPaul Mullowney 1702e057df02SPaul Mullowney /* these will be strings associated with enumerated type defined above */ 1703e057df02SPaul Mullowney PETSC_EXTERN const char *const MatCUSPStorageFormats[]; 1704e057df02SPaul Mullowney 17053f9c0db1SPaul Mullowney /*E 1706e057df02SPaul Mullowney MatCUSPFormatOperation - indicates the operation of CUSP (GPU) 170736d62e41SPaul Mullowney matrices whose operation should use a particular storage format. 1708e057df02SPaul Mullowney 1709e057df02SPaul Mullowney Not Collective 1710e057df02SPaul Mullowney 17118468deeeSKarl Rupp + MAT_CUSP_MULT_DIAG - sets the storage format for the diagonal matrix in the parallel MatMult 17128468deeeSKarl Rupp . MAT_CUSP_MULT_OFFDIAG - sets the storage format for the offdiagonal matrix in the parallel MatMult 17138468deeeSKarl Rupp . MAT_CUSP_MULT - sets the storage format for the entire matrix in the serial (single GPU) MatMult 17148468deeeSKarl Rupp - MAT_CUSP_ALL - sets the storage format for all CUSP (GPU) matrices 1715e057df02SPaul Mullowney 1716e057df02SPaul Mullowney Level: intermediate 1717e057df02SPaul Mullowney 1718e057df02SPaul Mullowney Any additions/changes here MUST also be made in include/finclude/petscmat.h 1719e057df02SPaul Mullowney 1720e057df02SPaul Mullowney .seealso: MatCUSPSetFormat(), MatCUSPStorageFormat 1721e057df02SPaul Mullowney E*/ 172236d62e41SPaul Mullowney typedef enum {MAT_CUSP_MULT_DIAG, MAT_CUSP_MULT_OFFDIAG, MAT_CUSP_MULT, MAT_CUSP_ALL} MatCUSPFormatOperation; 1723e057df02SPaul Mullowney 1724e057df02SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCUSPSetFormat(Mat,MatCUSPFormatOperation,MatCUSPStorageFormat); 1725e057df02SPaul Mullowney #endif 172690c53902SBarry Smith 1727d67ff14aSKarl Rupp #if defined(PETSC_HAVE_VIENNACL) 1728d67ff14aSKarl Rupp PETSC_EXTERN PetscErrorCode MatCreateSeqAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 1729d67ff14aSKarl Rupp PETSC_EXTERN PetscErrorCode MatCreateAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 1730d67ff14aSKarl Rupp #endif 1731d67ff14aSKarl Rupp 173254efbe56SHong Zhang /* 173354efbe56SHong Zhang PETSc interface to FFTW 173454efbe56SHong Zhang */ 173554efbe56SHong Zhang #if defined(PETSC_HAVE_FFTW) 1736014dd563SJed Brown PETSC_EXTERN PetscErrorCode VecScatterPetscToFFTW(Mat,Vec,Vec); 1737014dd563SJed Brown PETSC_EXTERN PetscErrorCode VecScatterFFTWToPetsc(Mat,Vec,Vec); 17382a7a6963SBarry Smith PETSC_EXTERN PetscErrorCode MatCreateVecsFFTW(Mat,Vec*,Vec*,Vec*); 173954efbe56SHong Zhang #endif 174054efbe56SHong Zhang 1741db31f6deSJed Brown #if defined(PETSC_HAVE_ELEMENTAL) 1742607a6623SBarry Smith PETSC_EXTERN PetscErrorCode PetscElementalInitializePackage(void); 1743db31f6deSJed Brown PETSC_EXTERN PetscErrorCode PetscElementalFinalizePackage(void); 1744db31f6deSJed Brown #endif 1745db31f6deSJed Brown 1746014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateNest(MPI_Comm,PetscInt,const IS[],PetscInt,const IS[],const Mat[],Mat*); 1747014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSize(Mat,PetscInt*,PetscInt*); 1748014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetISs(Mat,IS[],IS[]); 1749014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetLocalISs(Mat,IS[],IS[]); 1750014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSubMats(Mat,PetscInt*,PetscInt*,Mat***); 1751014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSubMat(Mat,PetscInt,PetscInt,Mat*); 175219fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatNestSetVecType(Mat,VecType); 1753014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestSetSubMats(Mat,PetscInt,const IS[],PetscInt,const IS[],const Mat[]); 1754014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestSetSubMat(Mat,PetscInt,PetscInt,Mat); 1755d8588912SDave May 17564325cce7SMatthew G Knepley PETSC_EXTERN PetscErrorCode MatChop(Mat,PetscReal); 1757e0a58c10SMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatComputeBandwidth(Mat,PetscReal,PetscInt*); 17584325cce7SMatthew G Knepley 175923a5497aSJed Brown #endif 1760