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*); 212014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat*); 213d21a29f3SJed Brown 214014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAdj(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscInt[],Mat*); 215014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 216d21a29f3SJed Brown 217014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 218014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPISBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *); 21938f409ebSLisandro Dalcin PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 220014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 2214cce697cSJed Brown PETSC_EXTERN PetscErrorCode MatXAIJSetPreallocation(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscInt[],const PetscInt[]); 222dfb205c3SBarry Smith 223014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateShell(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,void *,Mat*); 224014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateNormal(Mat,Mat*); 225014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateLRC(Mat,Mat,Mat,Mat*); 226014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateIS(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,ISLocalToGlobalMapping,Mat*); 227014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 228014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 229c0cdd4a1SDahai Guo 230014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 231014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 232014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqSBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 233014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPISBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 234c0cdd4a1SDahai Guo 235014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateScatter(MPI_Comm,VecScatter,Mat*); 236014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScatterSetVecScatter(Mat,VecScatter); 237014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScatterGetVecScatter(Mat,VecScatter*); 238014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateBlockMat(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,Mat*); 239014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeAddMat(Mat,Mat); 240014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeMerge(Mat); 241014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateComposite(MPI_Comm,PetscInt,const Mat*,Mat*); 2426d7c1e57SBarry Smith typedef enum {MAT_COMPOSITE_ADDITIVE,MAT_COMPOSITE_MULTIPLICATIVE} MatCompositeType; 243014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeSetType(Mat,MatCompositeType); 2446d7c1e57SBarry Smith 24519fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatCreateFFT(MPI_Comm,PetscInt,const PetscInt[],MatType,Mat*); 246014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqCUFFT(MPI_Comm,PetscInt,const PetscInt[],Mat*); 247dedccee8SHong Zhang 248014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateTranspose(Mat,Mat*); 249014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSubMatrix(Mat,IS,IS,Mat*); 250014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSubMatrixUpdate(Mat,Mat,IS,IS); 251014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateLocalRef(Mat,IS,IS,Mat*); 2521472f72bSBarry Smith 253014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPythonSetType(Mat,const char[]); 2541d6018f0SLisandro Dalcin 255014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetUp(Mat); 256014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroy(Mat*); 257e56f5c9eSBarry Smith PETSC_EXTERN PetscErrorCode MatGetNonzeroState(Mat,PetscObjectState*); 25821c89e3eSBarry Smith 259014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatConjugate(Mat); 260014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRealPart(Mat); 261014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatImaginaryPart(Mat); 262014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetDiagonalBlock(Mat,Mat*); 263014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetTrace(Mat,PetscScalar*); 264713ccfa9SJed Brown PETSC_EXTERN PetscErrorCode MatInvertBlockDiagonal(Mat,const PetscScalar **); 26599cafbc1SBarry Smith 2668ed539a5SBarry Smith /* ------------------------------------------------------------*/ 267014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 268014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 269014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesRow(Mat,PetscInt,const PetscScalar[]); 270014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesRowLocal(Mat,PetscInt,const PetscScalar[]); 271014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBatch(Mat,PetscInt,PetscInt,PetscInt[],const PetscScalar[]); 27273a71a0fSBarry Smith PETSC_EXTERN PetscErrorCode MatSetRandom(Mat,PetscRandom); 27384cb2905SBarry Smith 2742ef4de8bSBarry Smith /*S 2752ef4de8bSBarry Smith MatStencil - Data structure (C struct) for storing information about a single row or 276be479b30SJed Brown column of a matrix as indexed on an associated grid. 277be479b30SJed Brown 278be479b30SJed Brown Fortran usage is different, see MatSetValuesStencil() for details. 2792ef4de8bSBarry Smith 2802ef4de8bSBarry Smith Level: beginner 2812ef4de8bSBarry Smith 2822ef4de8bSBarry Smith Concepts: matrix; linear operator 2832ef4de8bSBarry Smith 284be479b30SJed Brown .seealso: MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockedStencil() 2852ef4de8bSBarry Smith S*/ 286435da068SBarry Smith typedef struct { 287c1ac3661SBarry Smith PetscInt k,j,i,c; 288435da068SBarry Smith } MatStencil; 2892ef4de8bSBarry Smith 290014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode); 291014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode); 292014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt); 293435da068SBarry Smith 294014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetColoring(Mat,ISColoring); 295014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesAdifor(Mat,PetscInt,void*); 2963a7fca6bSBarry Smith 297d91e6319SBarry Smith /*E 298d91e6319SBarry Smith MatAssemblyType - Indicates if the matrix is now to be used, or if you plan 299d91e6319SBarry Smith to continue to add values to it 300d91e6319SBarry Smith 301d91e6319SBarry Smith Level: beginner 302d91e6319SBarry Smith 303d91e6319SBarry Smith .seealso: MatAssemblyBegin(), MatAssemblyEnd() 304d91e6319SBarry Smith E*/ 3056d4a8577SBarry Smith typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType; 306014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssemblyBegin(Mat,MatAssemblyType); 307014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssemblyEnd(Mat,MatAssemblyType); 308014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssembled(Mat,PetscBool *); 3094f9c727eSBarry Smith 3101d73ed98SBarry Smith 31130de9b25SBarry Smith 312d91e6319SBarry Smith /*E 313d91e6319SBarry Smith MatOption - Options that may be set for a matrix and its behavior or storage 314d91e6319SBarry Smith 315d91e6319SBarry Smith Level: beginner 316d91e6319SBarry Smith 3170a835dfdSSatish Balay Any additions/changes here MUST also be made in include/finclude/petscmat.h 318d91e6319SBarry Smith 319382164b0SBarry Smith Developer Notes: Entries that are negative need not be called collectively by all processes. 320382164b0SBarry Smith 321d91e6319SBarry Smith .seealso: MatSetOption() 322d91e6319SBarry Smith E*/ 32392d4d306SBarry Smith typedef enum {MAT_OPTION_MIN = -5, 32492d4d306SBarry Smith MAT_NEW_NONZERO_LOCATION_ERR = -4, 32592d4d306SBarry Smith MAT_UNUSED_NONZERO_LOCATION_ERR = -3, 32692d4d306SBarry Smith MAT_NEW_NONZERO_ALLOCATION_ERR = -2, 32792d4d306SBarry Smith MAT_ROW_ORIENTED = -1, 328382164b0SBarry Smith MAT_SYMMETRIC = 1, 329382164b0SBarry Smith MAT_STRUCTURALLY_SYMMETRIC = 2, 330382164b0SBarry Smith MAT_NEW_DIAGONALS = 3, 331382164b0SBarry Smith MAT_IGNORE_OFF_PROC_ENTRIES = 4, 332382164b0SBarry Smith MAT_USE_HASH_TABLE = 5, 333382164b0SBarry Smith MAT_KEEP_NONZERO_PATTERN = 6, 334382164b0SBarry Smith MAT_IGNORE_ZERO_ENTRIES = 7, 335382164b0SBarry Smith MAT_USE_INODES = 8, 336382164b0SBarry Smith MAT_HERMITIAN = 9, 337382164b0SBarry Smith MAT_SYMMETRY_ETERNAL = 10, 33811e456e1SBarry Smith MAT_DUMMY = 11, 339382164b0SBarry Smith MAT_IGNORE_LOWER_TRIANGULAR = 12, 340382164b0SBarry Smith MAT_ERROR_LOWER_TRIANGULAR = 13, 341382164b0SBarry Smith MAT_GETROW_UPPERTRIANGULAR = 14, 342382164b0SBarry Smith MAT_SPD = 15, 34392d4d306SBarry Smith MAT_NO_OFF_PROC_ZERO_ROWS = 16, 34492d4d306SBarry Smith MAT_NO_OFF_PROC_ENTRIES = 17, 34592d4d306SBarry Smith MAT_NEW_NONZERO_LOCATIONS = 18, 34692d4d306SBarry Smith MAT_OPTION_MAX = 19} MatOption; 347e057df02SPaul Mullowney 348014dd563SJed Brown PETSC_EXTERN const char *MatOptions[]; 349014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetOption(Mat,MatOption,PetscBool ); 35019fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatGetType(Mat,MatType*); 35184cb2905SBarry Smith 352014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]); 353014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 354014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 355014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowUpperTriangular(Mat); 356014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRowUpperTriangular(Mat); 357014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 358014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 359014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnVector(Mat,Vec,PetscInt); 3608c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJGetArray(Mat,PetscScalar *[]); 3618c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJRestoreArray(Mat,PetscScalar *[]); 362*21e72a00SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJGetMaxRowNonzeros(Mat,PetscInt*); 363bd04181cSBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJSetValuesLocalFast(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 3648c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatDenseGetArray(Mat,PetscScalar *[]); 3658c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatDenseRestoreArray(Mat,PetscScalar *[]); 366014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBlockSize(Mat,PetscInt *); 367014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSize(Mat,PetscInt); 368014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBlockSizes(Mat,PetscInt *,PetscInt *); 369014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSizes(Mat,PetscInt,PetscInt); 37033d57670SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSizesFromMats(Mat,Mat,Mat); 371014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNThreads(Mat,PetscInt); 372014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNThreads(Mat,PetscInt*); 3731620fd73SBarry Smith 374014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMult(Mat,Vec,Vec); 375014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultDiagonalBlock(Mat,Vec,Vec); 376014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultAdd(Mat,Vec,Vec,Vec); 377014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTranspose(Mat,Vec,Vec); 378014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultHermitianTranspose(Mat,Vec,Vec); 379014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsTranspose(Mat,Mat,PetscReal,PetscBool *); 380014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitianTranspose(Mat,Mat,PetscReal,PetscBool *); 381014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeAdd(Mat,Vec,Vec,Vec); 382014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultHermitianTransposeAdd(Mat,Vec,Vec,Vec); 383014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultConstrained(Mat,Vec,Vec); 384014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeConstrained(Mat,Vec,Vec); 385014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatSolve(Mat,Mat,Mat); 386f9426fe0SMark Adams PETSC_EXTERN PetscErrorCode MatResidual(Mat,Vec,Vec,Vec); 387f5edf698SHong Zhang 388d91e6319SBarry Smith /*E 389d91e6319SBarry Smith MatDuplicateOption - Indicates if a duplicated sparse matrix should have 390d91e6319SBarry Smith its numerical values copied over or just its nonzero structure. 391d91e6319SBarry Smith 392d91e6319SBarry Smith Level: beginner 393d91e6319SBarry Smith 394d91e6319SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 395d91e6319SBarry Smith 39670dcbbb9SBarry Smith $ MAT_SHARE_NONZERO_PATTERN - the i and j arrays in the new matrix will be shared with the original matrix 39770dcbbb9SBarry Smith $ this also triggers the MAT_DO_NOT_COPY_VALUES option. This is used when you 39870dcbbb9SBarry Smith $ have several matrices with the same nonzero pattern. 39970dcbbb9SBarry Smith 400d91e6319SBarry Smith .seealso: MatDuplicate() 401d91e6319SBarry Smith E*/ 40270dcbbb9SBarry Smith typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES,MAT_SHARE_NONZERO_PATTERN} MatDuplicateOption; 4032e8a6d31SBarry Smith 40419fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatConvert(Mat,MatType,MatReuse,Mat*); 405014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDuplicate(Mat,MatDuplicateOption,Mat*); 40694a9d846SBarry Smith 407cb5b572fSBarry Smith 408014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCopy(Mat,Mat,MatStructure); 409014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatView(Mat,PetscViewer); 410014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsSymmetric(Mat,PetscReal,PetscBool *); 411014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsStructurallySymmetric(Mat,PetscBool *); 412014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitian(Mat,PetscReal,PetscBool *); 413014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsSymmetricKnown(Mat,PetscBool *,PetscBool *); 414014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitianKnown(Mat,PetscBool *,PetscBool *); 415014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMissingDiagonal(Mat,PetscBool *,PetscInt *); 416014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLoad(Mat, PetscViewer); 4177b80b807SBarry Smith 4181a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool *); 4191a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool *); 4201a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool *); 4211a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool *); 422d4fbbf0eSBarry Smith 423d91e6319SBarry Smith /*S 424d91e6319SBarry Smith MatInfo - Context of matrix information, used with MatGetInfo() 425d91e6319SBarry Smith 426d91e6319SBarry Smith In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE 427d91e6319SBarry Smith 428d91e6319SBarry Smith Level: intermediate 429d91e6319SBarry Smith 430d91e6319SBarry Smith Concepts: matrix^nonzero information 431d91e6319SBarry Smith 432d9274352SBarry Smith .seealso: MatGetInfo(), MatInfoType 433d91e6319SBarry Smith S*/ 4344e220ebcSLois Curfman McInnes typedef struct { 435b0a32e0cSBarry Smith PetscLogDouble block_size; /* block size */ 436b0a32e0cSBarry Smith PetscLogDouble nz_allocated,nz_used,nz_unneeded; /* number of nonzeros */ 437b0a32e0cSBarry Smith PetscLogDouble memory; /* memory allocated */ 438b0a32e0cSBarry Smith PetscLogDouble assemblies; /* number of matrix assemblies called */ 439b0a32e0cSBarry Smith PetscLogDouble mallocs; /* number of mallocs during MatSetValues() */ 440b0a32e0cSBarry Smith PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */ 441b0a32e0cSBarry Smith PetscLogDouble factor_mallocs; /* number of mallocs during factorization */ 4424e220ebcSLois Curfman McInnes } MatInfo; 4434e220ebcSLois Curfman McInnes 444d9274352SBarry Smith /*E 445d9274352SBarry Smith MatInfoType - Indicates if you want information about the local part of the matrix, 446d9274352SBarry Smith the entire parallel matrix or the maximum over all the local parts. 447d9274352SBarry Smith 448d9274352SBarry Smith Level: beginner 449d9274352SBarry Smith 450d9274352SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 451d9274352SBarry Smith 452d9274352SBarry Smith .seealso: MatGetInfo(), MatInfo 453d9274352SBarry Smith E*/ 4547b80b807SBarry Smith typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType; 455014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetInfo(Mat,MatInfoType,MatInfo*); 456014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetDiagonal(Mat,Vec); 457014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMax(Mat,Vec,PetscInt[]); 458014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMin(Mat,Vec,PetscInt[]); 459014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMaxAbs(Mat,Vec,PetscInt[]); 460014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMinAbs(Mat,Vec,PetscInt[]); 461014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowSum(Mat,Vec); 462014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTranspose(Mat,MatReuse,Mat*); 463014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatHermitianTranspose(Mat,MatReuse,Mat*); 464014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPermute(Mat,IS,IS,Mat *); 465014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalScale(Mat,Vec,Vec); 466014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalSet(Mat,Vec,InsertMode); 467014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatEqual(Mat,Mat,PetscBool *); 468014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultEqual(Mat,Mat,PetscInt,PetscBool *); 469014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultAddEqual(Mat,Mat,PetscInt,PetscBool *); 470014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeEqual(Mat,Mat,PetscInt,PetscBool *); 471014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscBool *); 4727b80b807SBarry Smith 473014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNorm(Mat,NormType,PetscReal *); 474014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnNorms(Mat,NormType,PetscReal *); 475014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroEntries(Mat); 476014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 477014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsIS(Mat,IS,PetscScalar,Vec,Vec); 478014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsStencil(Mat,PetscInt,const MatStencil [],PetscScalar,Vec,Vec); 479014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsStencil(Mat,PetscInt,const MatStencil[],PetscScalar,Vec,Vec); 480014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumns(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 481014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsIS(Mat,IS,PetscScalar,Vec,Vec); 4827b80b807SBarry Smith 483014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSize(Mat,PetscInt*,PetscInt*); 484014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalSize(Mat,PetscInt*,PetscInt*); 485014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRange(Mat,PetscInt*,PetscInt*); 486014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRanges(Mat,const PetscInt**); 487014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRangeColumn(Mat,PetscInt*,PetscInt*); 488014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRangesColumn(Mat,const PetscInt**); 489566876eaSJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipIS(Mat,IS*,IS*); 4905ef9f2a5SBarry Smith 491014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]); 492014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSubMatricesParallel(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]); 493014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroyMatrices(PetscInt,Mat *[]); 494014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSubMatrix(Mat,IS,IS,MatReuse,Mat *); 495014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalSubMatrix(Mat,IS,IS,Mat*); 496014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreLocalSubMatrix(Mat,IS,IS,Mat*); 497014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSeqNonzeroStructure(Mat,Mat*); 498014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroySeqNonzeroStructure(Mat*); 4997b80b807SBarry Smith 500014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJConcatenateSeqAIJ(MPI_Comm,Mat,PetscInt,MatReuse,Mat*); 501014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJConcatenateSeqAIJSymbolic(MPI_Comm,Mat,PetscInt,Mat*); 502014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJConcatenateSeqAIJNumeric(MPI_Comm,Mat,PetscInt,Mat); 503014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJ(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*); 504014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJSymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*); 505014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJNumeric(Mat,Mat); 506014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMat(Mat,MatReuse,Mat*); 507014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*); 508014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,Mat*); 509014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetGhosts(Mat, PetscInt *,const PetscInt *[]); 5108efafbd8SBarry Smith 511014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt); 5127b80b807SBarry Smith 513014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*); 514014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMultSymbolic(Mat,Mat,PetscReal,Mat*); 515014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMultNumeric(Mat,Mat,Mat); 51622440eb1SKris Buschelman 5177bab7c10SHong Zhang PETSC_EXTERN PetscErrorCode MatMatMatMult(Mat,Mat,Mat,MatReuse,PetscReal,Mat*); 5187bab7c10SHong Zhang PETSC_EXTERN PetscErrorCode MatMatMatMultSymbolic(Mat,Mat,Mat,PetscReal,Mat*); 5197bab7c10SHong Zhang PETSC_EXTERN PetscErrorCode MatMatMatMultNumeric(Mat,Mat,Mat,Mat); 5207bab7c10SHong Zhang 521014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*); 522014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAPSymbolic(Mat,Mat,PetscReal,Mat*); 523014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAPNumeric(Mat,Mat,Mat); 524014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARt(Mat,Mat,MatReuse,PetscReal,Mat*); 525014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARtSymbolic(Mat,Mat,PetscReal,Mat*); 526014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARtNumeric(Mat,Mat,Mat); 52722440eb1SKris Buschelman 528014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeMatMult(Mat,Mat,MatReuse,PetscReal,Mat*); 529014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposetMatMultSymbolic(Mat,Mat,PetscReal,Mat*); 530014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposetMatMultNumeric(Mat,Mat,Mat); 531014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMult(Mat,Mat,MatReuse,PetscReal,Mat*); 532014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMultSymbolic(Mat,Mat,PetscReal,Mat*); 533014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMultNumeric(Mat,Mat,Mat); 534bc011b1eSHong Zhang 535014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAXPY(Mat,PetscScalar,Mat,MatStructure); 536014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAYPX(Mat,PetscScalar,Mat,MatStructure); 5371c741599SHong Zhang 538014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScale(Mat,PetscScalar); 539014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShift(Mat,PetscScalar); 5407b80b807SBarry Smith 541014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping); 542014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*); 543a93ff8c4SPeter Brune PETSC_EXTERN PetscErrorCode MatGetLayouts(Mat,PetscLayout*,PetscLayout*); 544014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 545014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsLocalIS(Mat,IS,PetscScalar,Vec,Vec); 546014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 547014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocalIS(Mat,IS,PetscScalar,Vec,Vec); 548014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 549014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 550052efed2SBarry Smith 551014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStashSetInitialSize(Mat,PetscInt,PetscInt); 552014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*); 55390f02eecSBarry Smith 554014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInterpolate(Mat,Vec,Vec); 555014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInterpolateAdd(Mat,Vec,Vec,Vec); 556014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestrict(Mat,Vec,Vec); 557014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetVecs(Mat,Vec*,Vec*); 558b2bf6370SHong Zhang PETSC_EXTERN PetscErrorCode MatGetRedundantMatrix(Mat,PetscInt,MPI_Comm,MatReuse,Mat*); 559014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetMultiProcBlock(Mat,MPI_Comm,MatReuse,Mat*); 560014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFindZeroDiagonals(Mat,IS*); 5613a062f41SBarry Smith PETSC_EXTERN PetscErrorCode MatFindOffBlockDiagonalEntries(Mat,IS*); 562bd481603SBarry Smith 563bd481603SBarry Smith /*MC 5641620fd73SBarry Smith MatSetValue - Set a single entry into a matrix. 5651620fd73SBarry Smith 5661620fd73SBarry Smith Not collective 5671620fd73SBarry Smith 568a9834a7dSSatish Balay Synopsis: 569a9834a7dSSatish Balay #include <petscmat.h> 570a9834a7dSSatish Balay PetscErrorCode MatSetValue(Mat m,PetscInt row,PetscInt col,PetscScalar value,InsertMode mode) 571a9834a7dSSatish Balay 5721620fd73SBarry Smith Input Parameters: 5731620fd73SBarry Smith + m - the matrix 5741620fd73SBarry Smith . row - the row location of the entry 5751620fd73SBarry Smith . col - the column location of the entry 5761620fd73SBarry Smith . value - the value to insert 5771620fd73SBarry Smith - mode - either INSERT_VALUES or ADD_VALUES 5781620fd73SBarry Smith 5791620fd73SBarry Smith Notes: 5801620fd73SBarry Smith For efficiency one should use MatSetValues() and set several or many 5811620fd73SBarry Smith values simultaneously if possible. 5821620fd73SBarry Smith 5831620fd73SBarry Smith Level: beginner 5841620fd73SBarry Smith 5851620fd73SBarry Smith .seealso: MatSetValues(), MatSetValueLocal() 5861620fd73SBarry Smith M*/ 5871620fd73SBarry 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);} 5881620fd73SBarry Smith 5891620fd73SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatGetValue(Mat v,PetscInt i,PetscInt j,PetscScalar *va) {return MatGetValues(v,1,&i,1,&j,va);} 5901620fd73SBarry Smith 5911620fd73SBarry 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);} 5921620fd73SBarry Smith 5931620fd73SBarry Smith /*MC 5940d2a7ff2SSatish Balay MatPreallocateInitialize - Begins the block of code that will count the number of nonzeros per 595bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 596bd481603SBarry Smith 597bd481603SBarry Smith Synopsis: 598aaa7dc30SBarry Smith #include <petscmat.h> 599c1ac3661SBarry Smith PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz) 600bd481603SBarry Smith 601bd481603SBarry Smith Collective on MPI_Comm 602bd481603SBarry Smith 603bd481603SBarry Smith Input Parameters: 604bd481603SBarry Smith + comm - the communicator that will share the eventually allocated matrix 605859fcb39SBarry Smith . nrows - the number of LOCAL rows in the matrix 606859fcb39SBarry Smith - ncols - the number of LOCAL columns in the matrix 607bd481603SBarry Smith 608bd481603SBarry Smith Output Parameters: 609bd481603SBarry Smith + dnz - the array that will be passed to the matrix preallocation routines 610bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 611bd481603SBarry Smith 612bd481603SBarry Smith Level: intermediate 613bd481603SBarry Smith 614bd481603SBarry Smith Notes: 615a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 616bd481603SBarry Smith 6171d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 618bd481603SBarry Smith 619bd481603SBarry Smith Use MatPreallocateInitializeSymmetric() for symmetric matrices (MPISBAIJ matrices) 620bd481603SBarry Smith 6211620fd73SBarry Smith This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize(). 6221620fd73SBarry Smith 623bd481603SBarry Smith Concepts: preallocation^Matrix 624bd481603SBarry Smith 625d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(), 626d6e23781SBarry Smith MatPreallocateInitializeSymmetric(), MatPreallocateSymmetricSetLocalBlock() 627bd481603SBarry Smith M*/ 628c4f061fbSSatish Balay #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \ 6297c922b88SBarry Smith { \ 630a89bc211SBarry Smith PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__start,__end; \ 6311795a4d1SJed Brown _4_ierr = PetscCalloc2(__nrows,&dnz,__nrows,&onz);CHKERRQ(_4_ierr); \ 6321795a4d1SJed Brown __start = 0; __end = __start; \ 633c1ac3661SBarry Smith _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ctmp;\ 634a89bc211SBarry Smith _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows; 6357c922b88SBarry Smith 636bd481603SBarry Smith /*MC 637bd481603SBarry Smith MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be 638bd481603SBarry Smith inserted using a local number of the rows and columns 639bd481603SBarry Smith 640bd481603SBarry Smith Synopsis: 641aaa7dc30SBarry Smith #include <petscmat.h> 642c1ac3661SBarry Smith PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 643bd481603SBarry Smith 644bd481603SBarry Smith Not Collective 645bd481603SBarry Smith 646bd481603SBarry Smith Input Parameters: 647784ac674SJed Brown + map - the row mapping from local numbering to global numbering 648bd481603SBarry Smith . nrows - the number of rows indicated 6491d73ed98SBarry Smith . rows - the indices of the rows 650784ac674SJed Brown . cmap - the column mapping from local to global numbering 651bd481603SBarry Smith . ncols - the number of columns in the matrix 652bd481603SBarry Smith . cols - the columns indicated 653bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 654bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 655bd481603SBarry Smith 656bd481603SBarry Smith Level: intermediate 657bd481603SBarry Smith 658bd481603SBarry Smith Notes: 659a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 660bd481603SBarry Smith 6611d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 662bd481603SBarry Smith 663bd481603SBarry Smith Concepts: preallocation^Matrix 664bd481603SBarry Smith 665d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(), 666d6e23781SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock() 667bd481603SBarry Smith M*/ 668784ac674SJed Brown #define MatPreallocateSetLocal(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \ 669c4f061fbSSatish Balay {\ 670c1ac3661SBarry Smith PetscInt __l;\ 671784ac674SJed Brown _4_ierr = ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\ 672784ac674SJed Brown _4_ierr = ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\ 673c4f061fbSSatish Balay for (__l=0;__l<nrows;__l++) {\ 674ef66eb69SBarry Smith _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 675c4f061fbSSatish Balay }\ 676c4f061fbSSatish Balay } 677c4f061fbSSatish Balay 678bd481603SBarry Smith /*MC 679d6e23781SBarry Smith MatPreallocateSetLocalBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be 680bd481603SBarry Smith inserted using a local number of the rows and columns 681bd481603SBarry Smith 682bd481603SBarry Smith Synopsis: 683aaa7dc30SBarry Smith #include <petscmat.h> 684d6e23781SBarry Smith PetscErrorCode MatPreallocateSetLocalBlock(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 685d6e23781SBarry Smith 686d6e23781SBarry Smith Not Collective 687d6e23781SBarry Smith 688d6e23781SBarry Smith Input Parameters: 689d6e23781SBarry Smith + map - the row mapping from local numbering to global numbering 690d6e23781SBarry Smith . nrows - the number of rows indicated 691d6e23781SBarry Smith . rows - the indices of the rows 692d6e23781SBarry Smith . cmap - the column mapping from local to global numbering 693d6e23781SBarry Smith . ncols - the number of columns in the matrix 694d6e23781SBarry Smith . cols - the columns indicated 695d6e23781SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 696d6e23781SBarry Smith - ozn - the other array passed to the matrix preallocation routines 697d6e23781SBarry Smith 698d6e23781SBarry Smith Level: intermediate 699d6e23781SBarry Smith 700d6e23781SBarry Smith Notes: 701d6e23781SBarry Smith See Users-Manual: ch_performance for more details. 702d6e23781SBarry Smith 703d6e23781SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 704d6e23781SBarry Smith 705d6e23781SBarry Smith Concepts: preallocation^Matrix 706d6e23781SBarry Smith 707d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(), 708d6e23781SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock() 709d6e23781SBarry Smith M*/ 710d6e23781SBarry Smith #define MatPreallocateSetLocalBlock(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \ 711d6e23781SBarry Smith {\ 712d6e23781SBarry Smith PetscInt __l;\ 713d6e23781SBarry Smith _4_ierr = ISLocalToGlobalMappingApplyBlock(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\ 714d6e23781SBarry Smith _4_ierr = ISLocalToGlobalMappingApplyBlock(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\ 715d6e23781SBarry Smith for (__l=0;__l<nrows;__l++) {\ 716d6e23781SBarry Smith _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 717d6e23781SBarry Smith }\ 718d6e23781SBarry Smith } 719d6e23781SBarry Smith 720d6e23781SBarry Smith /*MC 721d6e23781SBarry Smith MatPreallocateSymmetricSetLocalBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be 722d6e23781SBarry Smith inserted using a local number of the rows and columns 723d6e23781SBarry Smith 724d6e23781SBarry Smith Synopsis: 725d6e23781SBarry Smith #include <petscmat.h> 726d6e23781SBarry Smith PetscErrorCode MatPreallocateSymmetricSetLocalBlock(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 727bd481603SBarry Smith 728bd481603SBarry Smith Not Collective 729bd481603SBarry Smith 730bd481603SBarry Smith Input Parameters: 731bd481603SBarry Smith + map - the mapping between local numbering and global numbering 732bd481603SBarry Smith . nrows - the number of rows indicated 7331d73ed98SBarry Smith . rows - the indices of the rows 734bd481603SBarry Smith . ncols - the number of columns in the matrix 735bd481603SBarry Smith . cols - the columns indicated 736bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 737bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 738bd481603SBarry Smith 739bd481603SBarry Smith Level: intermediate 740bd481603SBarry Smith 741bd481603SBarry Smith Notes: 742a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 743bd481603SBarry Smith 744bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 745bd481603SBarry Smith 746bd481603SBarry Smith Concepts: preallocation^Matrix 747bd481603SBarry Smith 748d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateInitialize(), 749d6e23781SBarry Smith MatPreallocateInitialize(), MatPreallocateSetLocal() 750bd481603SBarry Smith M*/ 751d6e23781SBarry Smith #define MatPreallocateSymmetricSetLocalBlock(map,nrows,rows,ncols,cols,dnz,onz) 0;\ 752d3d32019SBarry Smith {\ 753c1ac3661SBarry Smith PetscInt __l;\ 754d6e23781SBarry Smith _4_ierr = ISLocalToGlobalMappingApplyBlock(map,nrows,rows,rows);CHKERRQ(_4_ierr);\ 755d6e23781SBarry Smith _4_ierr = ISLocalToGlobalMappingApplyBlock(map,ncols,cols,cols);CHKERRQ(_4_ierr);\ 756d3d32019SBarry Smith for (__l=0;__l<nrows;__l++) {\ 757d6e23781SBarry Smith _4_ierr = MatPreallocateSymmetricSetBlock((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 758d3d32019SBarry Smith }\ 759d3d32019SBarry Smith } 760bd481603SBarry Smith /*MC 761bd481603SBarry Smith MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be 762bd481603SBarry Smith inserted using a local number of the rows and columns 763bd481603SBarry Smith 764bd481603SBarry Smith Synopsis: 765aaa7dc30SBarry Smith #include <petscmat.h> 766c1ac3661SBarry Smith PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 767bd481603SBarry Smith 768bd481603SBarry Smith Not Collective 769bd481603SBarry Smith 770bd481603SBarry Smith Input Parameters: 77164075487SBarry Smith + row - the row 772bd481603SBarry Smith . ncols - the number of columns in the matrix 773a50f8bf6SHong Zhang - cols - the columns indicated 774a50f8bf6SHong Zhang 775a50f8bf6SHong Zhang Output Parameters: 776a50f8bf6SHong Zhang + dnz - the array that will be passed to the matrix preallocation routines 777bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 778bd481603SBarry Smith 779bd481603SBarry Smith Level: intermediate 780bd481603SBarry Smith 781bd481603SBarry Smith Notes: 782a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 783bd481603SBarry Smith 784bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 785bd481603SBarry Smith 7861620fd73SBarry Smith This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize(). 7871620fd73SBarry Smith 788bd481603SBarry Smith Concepts: preallocation^Matrix 789bd481603SBarry Smith 790d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(), 791d6e23781SBarry Smith MatPreallocateInitialize(), MatPreallocateSetLocal() 792bd481603SBarry Smith M*/ 793c4f061fbSSatish Balay #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\ 794c1ac3661SBarry Smith { PetscInt __i; \ 795e32f2f54SBarry 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);\ 796e32f2f54SBarry 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);\ 7977c922b88SBarry Smith for (__i=0; __i<nc; __i++) {\ 79864075487SBarry Smith if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \ 7997cc688d7SBarry Smith else dnz[row - __rstart]++;\ 8007c922b88SBarry Smith }\ 8017c922b88SBarry Smith } 8027c922b88SBarry Smith 803bd481603SBarry Smith /*MC 804d6e23781SBarry Smith MatPreallocateSymmetricSetBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be 805bd481603SBarry Smith inserted using a local number of the rows and columns 806bd481603SBarry Smith 807bd481603SBarry Smith Synopsis: 808aaa7dc30SBarry Smith #include <petscmat.h> 809d6e23781SBarry Smith PetscErrorCode MatPreallocateSymmetricSetBlock(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 810bd481603SBarry Smith 811bd481603SBarry Smith Not Collective 812bd481603SBarry Smith 813bd481603SBarry Smith Input Parameters: 814bd481603SBarry Smith + nrows - the number of rows indicated 8151d73ed98SBarry Smith . rows - the indices of the rows 816bd481603SBarry Smith . ncols - the number of columns in the matrix 817bd481603SBarry Smith . cols - the columns indicated 818bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 819bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 820bd481603SBarry Smith 821bd481603SBarry Smith Level: intermediate 822bd481603SBarry Smith 823bd481603SBarry Smith Notes: 824a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 825bd481603SBarry Smith 826bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 827bd481603SBarry Smith 8281620fd73SBarry Smith This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize(). 8291620fd73SBarry Smith 830bd481603SBarry Smith Concepts: preallocation^Matrix 831bd481603SBarry Smith 832d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateInitialize(), 833d6e23781SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocal() 834bd481603SBarry Smith M*/ 835d6e23781SBarry Smith #define MatPreallocateSymmetricSetBlock(row,nc,cols,dnz,onz) 0;\ 836c1ac3661SBarry Smith { PetscInt __i; \ 837d3d32019SBarry Smith for (__i=0; __i<nc; __i++) {\ 838d3d32019SBarry Smith if (cols[__i] >= __end) onz[row - __rstart]++; \ 839d3d32019SBarry Smith else if (cols[__i] >= row) dnz[row - __rstart]++;\ 840d3d32019SBarry Smith }\ 841d3d32019SBarry Smith } 842d3d32019SBarry Smith 843bd481603SBarry Smith /*MC 84416371a99SBarry Smith MatPreallocateLocation - An alternative to MatPreallocationSet() that puts the nonzero locations into the matrix if it exists 84516371a99SBarry Smith 84616371a99SBarry Smith Synopsis: 847aaa7dc30SBarry Smith #include <petscmat.h> 84816371a99SBarry Smith PetscErrorCode MatPreallocateLocations(Mat A,PetscInt row,PetscInt ncols,PetscInt *cols,PetscInt *dnz,PetscInt *onz) 84916371a99SBarry Smith 85016371a99SBarry Smith Not Collective 85116371a99SBarry Smith 85216371a99SBarry Smith Input Parameters: 85316371a99SBarry Smith . A - matrix 85416371a99SBarry Smith . row - row where values exist (must be local to this process) 85516371a99SBarry Smith . ncols - number of columns 85616371a99SBarry Smith . cols - columns with nonzeros 85716371a99SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 85816371a99SBarry Smith - ozn - the other array passed to the matrix preallocation routines 85916371a99SBarry Smith 86016371a99SBarry Smith Level: intermediate 86116371a99SBarry Smith 86216371a99SBarry Smith Notes: 863a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 86416371a99SBarry Smith 86516371a99SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 86616371a99SBarry Smith 86716371a99SBarry Smith This is a MACRO not a function because it uses a bunch of variables private to the MatPreallocation.... routines. 86816371a99SBarry Smith 86916371a99SBarry Smith Concepts: preallocation^Matrix 87016371a99SBarry Smith 871d6e23781SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(), 872d6e23781SBarry Smith MatPreallocateSymmetricSetLocalBlock() 87316371a99SBarry Smith M*/ 8740298fd71SBarry 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);} 87516371a99SBarry Smith 87616371a99SBarry Smith 87716371a99SBarry Smith /*MC 8780d2a7ff2SSatish Balay MatPreallocateFinalize - Ends the block of code that will count the number of nonzeros per 879bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 880bd481603SBarry Smith 881bd481603SBarry Smith Synopsis: 882aaa7dc30SBarry Smith #include <petscmat.h> 883c1ac3661SBarry Smith PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz) 884bd481603SBarry Smith 885bd481603SBarry Smith Collective on MPI_Comm 886bd481603SBarry Smith 887bd481603SBarry Smith Input Parameters: 88816371a99SBarry Smith + dnz - the array that was be passed to the matrix preallocation routines 889bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 890bd481603SBarry Smith 891bd481603SBarry Smith Level: intermediate 892bd481603SBarry Smith 893bd481603SBarry Smith Notes: 894a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 895bd481603SBarry Smith 896bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 897bd481603SBarry Smith 8981620fd73SBarry Smith This is a MACRO not a function because it closes the { started in MatPreallocateInitialize(). 8991620fd73SBarry Smith 900bd481603SBarry Smith Concepts: preallocation^Matrix 901bd481603SBarry Smith 902d6e23781SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(), 903d6e23781SBarry Smith MatPreallocateSymmetricSetLocalBlock() 904bd481603SBarry Smith M*/ 905a89bc211SBarry Smith #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree2(dnz,onz);CHKERRQ(_4_ierr);} 9067c922b88SBarry Smith 9077b80b807SBarry Smith /* Routines unique to particular data structures */ 908014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellGetContext(Mat,void *); 909698d4c6aSKris Buschelman 910014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInodeAdjustForInodes(Mat,IS*,IS*); 911014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *); 912698d4c6aSKris Buschelman 913014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetColumnIndices(Mat,PetscInt[]); 914014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetColumnIndices(Mat,PetscInt[]); 915014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 916014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 917014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 918014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJFromTriple(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*,PetscInt,PetscBool); 9197b80b807SBarry Smith 920a96a251dSBarry Smith #define MAT_SKIP_ALLOCATION -4 921a96a251dSBarry Smith 922014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]); 923014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]); 924014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]); 925273d9f13SBarry Smith 926014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 927014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 928014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 929014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []); 930014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 931014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]); 932014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 933014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]); 934014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIDenseSetPreallocation(Mat,PetscScalar[]); 935014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqDenseSetPreallocation(Mat,PetscScalar[]); 9369230625dSJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,const PetscInt*[]); 9379230625dSJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,const PetscInt*[]); 938014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAdjCreateNonemptySubcommMat(Mat,Mat*); 939273d9f13SBarry Smith 940014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqDenseSetLDA(Mat,PetscInt); 941014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDenseGetLocalMatrix(Mat,Mat*); 9421b807ce4Svictorle 943014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStoreValues(Mat); 944014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRetrieveValues(Mat); 9452e8a6d31SBarry Smith 946014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDAADSetCtx(Mat,void*); 9473a7fca6bSBarry Smith 948014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFindNonzeroRows(Mat,IS*); 9497b80b807SBarry Smith /* 9507b80b807SBarry Smith These routines are not usually accessed directly, rather solving is 95194b7f48cSBarry Smith done through the KSP and PC interfaces. 9527b80b807SBarry Smith */ 9537b80b807SBarry Smith 95476bdecfbSBarry Smith /*J 9558f6c3df8SBarry Smith MatOrderingType - String with the name of a PETSc matrix ordering 956d9274352SBarry Smith 957d9274352SBarry Smith Level: beginner 958d9274352SBarry Smith 959d9274352SBarry Smith .seealso: MatGetOrdering() 96076bdecfbSBarry Smith J*/ 96119fd82e9SBarry Smith typedef const char* MatOrderingType; 9622692d6eeSBarry Smith #define MATORDERINGNATURAL "natural" 9632692d6eeSBarry Smith #define MATORDERINGND "nd" 9642692d6eeSBarry Smith #define MATORDERING1WD "1wd" 9652692d6eeSBarry Smith #define MATORDERINGRCM "rcm" 9662692d6eeSBarry Smith #define MATORDERINGQMD "qmd" 9672692d6eeSBarry Smith #define MATORDERINGROWLENGTH "rowlength" 968510184a7SJed Brown #define MATORDERINGWBM "wbm" 969fc1bef75SJed Brown #define MATORDERINGSPECTRAL "spectral" 970312e372aSBarry Smith #define MATORDERINGAMD "amd" /* only works if UMFPACK is installed with PETSc */ 971b12f92e5SBarry Smith 97219fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatGetOrdering(Mat,MatOrderingType,IS*,IS*); 973140e18c1SBarry Smith PETSC_EXTERN PetscErrorCode MatGetOrderingList(PetscFunctionList*); 974bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatOrderingRegister(const char[],PetscErrorCode(*)(Mat,MatOrderingType,IS*,IS*)); 975607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatOrderingRegisterAll(void); 976014dd563SJed Brown PETSC_EXTERN PetscBool MatOrderingRegisterAllCalled; 977140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatOrderingList; 978d4fbbf0eSBarry Smith 979014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS); 980fc1bef75SJed Brown PETSC_EXTERN PetscErrorCode MatCreateLaplacian(Mat,PetscReal,PetscBool,Mat*); 981a2ce50c7SBarry Smith 982d91e6319SBarry Smith /*S 983d90ac83dSHong Zhang MatFactorShiftType - Numeric Shift. 984d90ac83dSHong Zhang 985d90ac83dSHong Zhang Level: beginner 986d90ac83dSHong Zhang 987d90ac83dSHong Zhang S*/ 988d90ac83dSHong Zhang typedef enum {MAT_SHIFT_NONE,MAT_SHIFT_NONZERO,MAT_SHIFT_POSITIVE_DEFINITE,MAT_SHIFT_INBLOCKS} MatFactorShiftType; 9896a6fc655SJed Brown PETSC_EXTERN const char *const MatFactorShiftTypes[]; 9905e9742b9SJed Brown PETSC_EXTERN const char *const MatFactorShiftTypesDetail[]; 991d90ac83dSHong Zhang 992d90ac83dSHong Zhang /*S 9932401956bSBarry Smith MatFactorInfo - Data passed into the matrix factorization routines 994b00f7748SHong Zhang 99561cad860SBarry Smith In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE, that is use 99661cad860SBarry Smith $ MatFactorInfo info(MAT_FACTORINFO_SIZE) 997b00f7748SHong Zhang 99815e8a5b3SHong Zhang Notes: These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC. 999b00f7748SHong Zhang 100061cad860SBarry Smith You can use MatFactorInfoInitialize() to set default values. 100161cad860SBarry Smith 1002b00f7748SHong Zhang Level: developer 1003b00f7748SHong Zhang 1004d7d82daaSBarry Smith .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(), 1005d7d82daaSBarry Smith MatFactorInfoInitialize() 1006b00f7748SHong Zhang 1007b00f7748SHong Zhang S*/ 1008b00f7748SHong Zhang typedef struct { 100915e8a5b3SHong Zhang PetscReal diagonal_fill; /* force diagonal to fill in if initially not filled */ 101085317021SBarry Smith PetscReal usedt; 101115e8a5b3SHong Zhang PetscReal dt; /* drop tolerance */ 1012b00f7748SHong Zhang PetscReal dtcol; /* tolerance for pivoting */ 101315e8a5b3SHong Zhang PetscReal dtcount; /* maximum nonzeros to be allowed per row */ 101467e28bfeSBarry Smith PetscReal fill; /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */ 1015348344bbSBarry Smith PetscReal levels; /* ICC/ILU(levels) */ 1016bcd9e38bSBarry Smith PetscReal pivotinblocks; /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0 1017bcd9e38bSBarry Smith factorization may be faster if do not pivot */ 101815e8a5b3SHong Zhang PetscReal zeropivot; /* pivot is called zero if less than this */ 1019f4db908eSBarry Smith PetscReal shifttype; /* type of shift added to matrix factor to prevent zero pivots */ 1020f4db908eSBarry Smith PetscReal shiftamount; /* how large the shift is */ 102115e8a5b3SHong Zhang } MatFactorInfo; 1022ffa6d0a5SLois Curfman McInnes 1023014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFactorInfoInitialize(MatFactorInfo*); 1024014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCholeskyFactor(Mat,IS,const MatFactorInfo*); 1025014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCholeskyFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*); 1026014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCholeskyFactorNumeric(Mat,Mat,const MatFactorInfo*); 1027014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLUFactor(Mat,IS,IS,const MatFactorInfo*); 1028014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatILUFactor(Mat,IS,IS,const MatFactorInfo*); 1029014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*); 1030014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatILUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*); 1031014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatICCFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*); 1032014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatICCFactor(Mat,IS,const MatFactorInfo*); 1033014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLUFactorNumeric(Mat,Mat,const MatFactorInfo*); 1034014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*); 1035014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolve(Mat,Vec,Vec); 1036014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatForwardSolve(Mat,Vec,Vec); 1037014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatBackwardSolve(Mat,Vec,Vec); 1038014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveAdd(Mat,Vec,Vec,Vec); 1039014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveTranspose(Mat,Vec,Vec); 1040014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveTransposeAdd(Mat,Vec,Vec,Vec); 1041014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolves(Mat,Vecs,Vecs); 10428ed539a5SBarry Smith 1043014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetUnfactored(Mat); 1044bb5a7306SBarry Smith 1045d91e6319SBarry Smith /*E 1046d91e6319SBarry Smith MatSORType - What type of (S)SOR to perform 1047bb1eb677SSatish Balay 1048d91e6319SBarry Smith Level: beginner 1049d91e6319SBarry Smith 1050d9274352SBarry Smith May be bitwise ORd together 1051d9274352SBarry Smith 1052d91e6319SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 1053d91e6319SBarry Smith 10544e7234bfSBarry Smith MatSORType may be bitwise ORd together, so do not change the numbers 10554e7234bfSBarry Smith 105641f059aeSBarry Smith .seealso: MatSOR() 1057d91e6319SBarry Smith E*/ 1058ee50ffe9SBarry Smith typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3, 1059ee50ffe9SBarry Smith SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8, 1060ee50ffe9SBarry Smith SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16, 106184cb2905SBarry Smith SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType; 1062014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSOR(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec); 10638ed539a5SBarry Smith 1064d4fbbf0eSBarry Smith /* 1065639f9d9dSBarry Smith These routines are for efficiently computing Jacobians via finite differences. 1066639f9d9dSBarry Smith */ 1067b12f92e5SBarry Smith 10687086a01eSPeter Brune /*S 10697086a01eSPeter Brune MatColoring - Object for managing the coloring of matrices. 10707086a01eSPeter Brune 10717086a01eSPeter Brune Level: beginner 10727086a01eSPeter Brune 10737086a01eSPeter Brune Concepts: matrix, coloring 10747086a01eSPeter Brune 10757086a01eSPeter Brune .seealso: MatFDColoringCreate() ISColoring MatFDColoring 10767086a01eSPeter Brune S*/ 10777086a01eSPeter Brune typedef struct _p_MatColoring* MatColoring; 107876bdecfbSBarry Smith /*J 10798f6c3df8SBarry Smith MatColoringType - String with the name of a PETSc matrix coloring 1080d9274352SBarry Smith 1081d9274352SBarry Smith Level: beginner 1082d9274352SBarry Smith 10837086a01eSPeter Brune .seealso: MatColoringSetType(), MatColoring 108476bdecfbSBarry Smith J*/ 10857086a01eSPeter Brune 108619fd82e9SBarry Smith typedef const char* MatColoringType; 10875567d71aSPeter Brune #define MATCOLORINGJP "jp" 1088b9acb617SPeter Brune #define MATCOLORINGPOWER "power" 10892692d6eeSBarry Smith #define MATCOLORINGNATURAL "natural" 10902692d6eeSBarry Smith #define MATCOLORINGSL "sl" 10912692d6eeSBarry Smith #define MATCOLORINGLF "lf" 10922692d6eeSBarry Smith #define MATCOLORINGID "id" 10938121bdceSPeter Brune #define MATCOLORINGGREEDY "greedy" 1094b12f92e5SBarry Smith 10958ac6da0aSPeter Brune /*E 10968ac6da0aSPeter Brune MatColoringWeightType - Type of weight scheme 10978ac6da0aSPeter Brune 10988ac6da0aSPeter Brune Not Collective 10998ac6da0aSPeter Brune 11008ac6da0aSPeter Brune + MAT_COLORING_RANDOM - Random weights 11018ac6da0aSPeter Brune . MAT_COLORING_LEXICAL - Lexical weighting based upon global numbering. 11028ac6da0aSPeter Brune - MAT_COLORING_LF - Last-first weighting. 11038ac6da0aSPeter Brune 11048ac6da0aSPeter Brune Level: intermediate 11058ac6da0aSPeter Brune 11068ac6da0aSPeter Brune Any additions/changes here MUST also be made in include/finclude/petscmat.h 11078ac6da0aSPeter Brune 11088ac6da0aSPeter Brune .seealso: MatCUSPSetFormat(), MatCUSPFormatOperation 11098ac6da0aSPeter Brune E*/ 1110301c140bSPeter Brune typedef enum {MAT_COLORING_WEIGHT_RANDOM,MAT_COLORING_WEIGHT_LEXICAL,MAT_COLORING_WEIGHT_LF,MAT_COLORING_WEIGHT_SL} MatColoringWeightType; 11118ac6da0aSPeter Brune 1112335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringCreate(Mat,MatColoring*); 1113d7f2a307SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringGetDegrees(Mat,PetscInt,PetscInt*); 1114335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringDestroy(MatColoring*); 1115335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringView(MatColoring,PetscViewer); 1116335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetType(MatColoring,MatColoringType); 1117335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetFromOptions(MatColoring); 1118335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetDistance(MatColoring,PetscInt); 1119335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringGetDistance(MatColoring,PetscInt*); 1120335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetMaxColors(MatColoring,PetscInt); 1121335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringGetMaxColors(MatColoring,PetscInt*); 1122335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringApply(MatColoring,ISColoring*); 1123607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatColoringRegisterAll(void); 1124335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringRegister(const char[],PetscErrorCode(*)(MatColoring)); 1125335efc43SPeter Brune PETSC_EXTERN PetscBool MatColoringRegisterAllCalled; 1126014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*); 11278ac6da0aSPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetWeightType(MatColoring,MatColoringWeightType); 1128c29971fcSPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetWeights(MatColoring,PetscReal*,PetscInt*); 11298ac6da0aSPeter Brune PETSC_EXTERN PetscErrorCode MatColoringCreateWeights(MatColoring,PetscReal **,PetscInt **lperm); 1130639f9d9dSBarry Smith 1131d9274352SBarry Smith /*S 1132d9274352SBarry Smith MatFDColoring - Object for computing a sparse Jacobian via finite differences 1133d9274352SBarry Smith and coloring 1134639f9d9dSBarry Smith 1135d9274352SBarry Smith Level: beginner 1136d9274352SBarry Smith 1137d9274352SBarry Smith Concepts: coloring, sparse Jacobian, finite differences 1138d9274352SBarry Smith 1139d9274352SBarry Smith .seealso: MatFDColoringCreate() 1140d9274352SBarry Smith S*/ 1141e2a1c21fSSatish Balay typedef struct _p_MatFDColoring* MatFDColoring; 1142639f9d9dSBarry Smith 1143014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringCreate(Mat,ISColoring,MatFDColoring *); 1144014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringDestroy(MatFDColoring*); 1145014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringView(MatFDColoring,PetscViewer); 1146014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*); 1147014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**); 1148014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal); 1149014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetFromOptions(MatFDColoring); 1150d1e9a80fSBarry Smith PETSC_EXTERN PetscErrorCode MatFDColoringApply(Mat,MatFDColoring,Vec,void *); 1151014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetF(MatFDColoring,Vec); 1152014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,PetscInt*[]); 1153f86b9fbaSHong Zhang PETSC_EXTERN PetscErrorCode MatFDColoringSetUp(Mat,ISColoring,MatFDColoring); 1154f86b9fbaSHong Zhang PETSC_EXTERN PetscErrorCode MatFDColoringSetBlockSize(MatFDColoring,PetscInt,PetscInt); 1155f86b9fbaSHong Zhang 1156b1683b59SHong Zhang 1157b1683b59SHong Zhang /*S 1158b9af6bddSHong Zhang MatTransposeColoring - Object for computing a sparse matrix product C=A*B^T via coloring 1159b1683b59SHong Zhang 1160b1683b59SHong Zhang Level: beginner 1161b1683b59SHong Zhang 1162b1683b59SHong Zhang Concepts: coloring, sparse matrix product 1163b1683b59SHong Zhang 1164b9af6bddSHong Zhang .seealso: MatTransposeColoringCreate() 1165b1683b59SHong Zhang S*/ 1166b9af6bddSHong Zhang typedef struct _p_MatTransposeColoring* MatTransposeColoring; 1167b1683b59SHong Zhang 1168014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeColoringCreate(Mat,ISColoring,MatTransposeColoring *); 1169014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransColoringApplySpToDen(MatTransposeColoring,Mat,Mat); 1170014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransColoringApplyDenToSp(MatTransposeColoring,Mat,Mat); 1171014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeColoringDestroy(MatTransposeColoring*); 1172b1683b59SHong Zhang 1173639f9d9dSBarry Smith /* 11740752156aSBarry Smith These routines are for partitioning matrices: currently used only 11753eda8832SBarry Smith for adjacency matrix, MatCreateMPIAdj(). 11760752156aSBarry Smith */ 1177ca161407SBarry Smith 1178d9274352SBarry Smith /*S 1179d9274352SBarry Smith MatPartitioning - Object for managing the partitioning of a matrix or graph 1180d9274352SBarry Smith 1181d9274352SBarry Smith Level: beginner 1182d9274352SBarry Smith 1183d9274352SBarry Smith Concepts: partitioning 1184d9274352SBarry Smith 1185743a1026Svictorle .seealso: MatPartitioningCreate(), MatPartitioningType 1186d9274352SBarry Smith S*/ 118791e9ee9fSBarry Smith typedef struct _p_MatPartitioning* MatPartitioning; 1188d9274352SBarry Smith 118976bdecfbSBarry Smith /*J 11908f6c3df8SBarry Smith MatPartitioningType - String with the name of a PETSc matrix partitioning 1191d9274352SBarry Smith 1192d9274352SBarry Smith Level: beginner 11930d04baf8SBarry Smith dm 1194b547a13bSvictorle .seealso: MatPartitioningCreate(), MatPartitioning 119576bdecfbSBarry Smith J*/ 119619fd82e9SBarry Smith typedef const char* MatPartitioningType; 11972692d6eeSBarry Smith #define MATPARTITIONINGCURRENT "current" 11982692d6eeSBarry Smith #define MATPARTITIONINGSQUARE "square" 11992692d6eeSBarry Smith #define MATPARTITIONINGPARMETIS "parmetis" 12002692d6eeSBarry Smith #define MATPARTITIONINGCHACO "chaco" 12012692d6eeSBarry Smith #define MATPARTITIONINGPARTY "party" 120250d91057SBarry Smith #define MATPARTITIONINGPTSCOTCH "ptscotch" 120371306c60Sjroman 1204ca161407SBarry Smith 1205014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningCreate(MPI_Comm,MatPartitioning*); 120619fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningSetType(MatPartitioning,MatPartitioningType); 1207014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetNParts(MatPartitioning,PetscInt); 1208014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetAdjacency(MatPartitioning,Mat); 1209014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]); 1210014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []); 1211014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningApply(MatPartitioning,IS*); 1212014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningDestroy(MatPartitioning*); 12132aabb6bbSBarry Smith 1214bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningRegister(const char[],PetscErrorCode (*)(MatPartitioning)); 121530de9b25SBarry Smith 1216014dd563SJed Brown PETSC_EXTERN PetscBool MatPartitioningRegisterAllCalled; 1217f1144a30SSatish Balay 1218607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningRegisterAll(void); 12192bad1931SBarry Smith 1220014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningView(MatPartitioning,PetscViewer); 1221014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetFromOptions(MatPartitioning); 122219fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningGetType(MatPartitioning,MatPartitioningType*); 1223ca161407SBarry Smith 1224014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetCoarseSequential(MatPartitioning); 1225014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningParmetisGetEdgeCut(MatPartitioning, PetscInt *); 12260752156aSBarry Smith 1227b6956c12SJose Roman typedef enum { MP_CHACO_MULTILEVEL=1,MP_CHACO_SPECTRAL=2,MP_CHACO_LINEAR=4,MP_CHACO_RANDOM=5,MP_CHACO_SCATTERED=6 } MPChacoGlobalType; 12286a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoGlobalTypes[]; 1229b6956c12SJose Roman typedef enum { MP_CHACO_KERNIGHAN=1,MP_CHACO_NONE=2 } MPChacoLocalType; 12306a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoLocalTypes[]; 1231b6956c12SJose Roman typedef enum { MP_CHACO_LANCZOS=0,MP_CHACO_RQI=1 } MPChacoEigenType; 12326a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoEigenTypes[]; 1233b6956c12SJose Roman 1234014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetGlobal(MatPartitioning,MPChacoGlobalType); 1235014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetGlobal(MatPartitioning,MPChacoGlobalType*); 1236014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetLocal(MatPartitioning,MPChacoLocalType); 1237014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetLocal(MatPartitioning,MPChacoLocalType*); 1238014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal); 1239014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType); 1240014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenSolver(MatPartitioning,MPChacoEigenType*); 1241014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenTol(MatPartitioning,PetscReal); 1242014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenTol(MatPartitioning,PetscReal*); 1243014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenNumber(MatPartitioning,PetscInt); 1244014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenNumber(MatPartitioning,PetscInt*); 124571306c60Sjroman 124671306c60Sjroman #define MP_PARTY_OPT "opt" 124771306c60Sjroman #define MP_PARTY_LIN "lin" 124871306c60Sjroman #define MP_PARTY_SCA "sca" 124971306c60Sjroman #define MP_PARTY_RAN "ran" 125071306c60Sjroman #define MP_PARTY_GBF "gbf" 125171306c60Sjroman #define MP_PARTY_GCF "gcf" 125271306c60Sjroman #define MP_PARTY_BUB "bub" 125371306c60Sjroman #define MP_PARTY_DEF "def" 1254014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetGlobal(MatPartitioning,const char*); 125571306c60Sjroman #define MP_PARTY_HELPFUL_SETS "hs" 125671306c60Sjroman #define MP_PARTY_KERNIGHAN_LIN "kl" 125771306c60Sjroman #define MP_PARTY_NONE "no" 1258014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetLocal(MatPartitioning,const char*); 1259014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal); 1260014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetBipart(MatPartitioning,PetscBool); 1261014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscBool); 126271306c60Sjroman 126350d91057SBarry Smith typedef enum { MP_PTSCOTCH_QUALITY,MP_PTSCOTCH_SPEED,MP_PTSCOTCH_BALANCE,MP_PTSCOTCH_SAFETY,MP_PTSCOTCH_SCALABILITY } MPPTScotchStrategyType; 12646a6fc655SJed Brown PETSC_EXTERN const char *const MPPTScotchStrategyTypes[]; 1265e0f1cffaSJose Roman 1266014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetImbalance(MatPartitioning,PetscReal); 1267014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetImbalance(MatPartitioning,PetscReal*); 1268014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetStrategy(MatPartitioning,MPPTScotchStrategyType); 1269014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetStrategy(MatPartitioning,MPPTScotchStrategyType*); 127071306c60Sjroman 1271b43b03e9SMark F. Adams /* 1272b43b03e9SMark F. Adams These routines are for coarsening matrices: 1273b43b03e9SMark F. Adams */ 1274b43b03e9SMark F. Adams 1275b43b03e9SMark F. Adams /*S 1276b43b03e9SMark F. Adams MatCoarsen - Object for managing the coarsening of a graph (symmetric matrix) 1277b43b03e9SMark F. Adams 1278b43b03e9SMark F. Adams Level: beginner 1279b43b03e9SMark F. Adams 1280b43b03e9SMark F. Adams Concepts: coarsen 1281b43b03e9SMark F. Adams 1282b43b03e9SMark F. Adams .seealso: MatCoarsenCreate), MatCoarsenType 1283b43b03e9SMark F. Adams S*/ 1284b43b03e9SMark F. Adams typedef struct _p_MatCoarsen* MatCoarsen; 1285b43b03e9SMark F. Adams 1286b43b03e9SMark F. Adams /*J 12878f6c3df8SBarry Smith MatCoarsenType - String with the name of a PETSc matrix coarsen 1288b43b03e9SMark F. Adams 1289b43b03e9SMark F. Adams Level: beginner 12908f6c3df8SBarry Smith 1291b43b03e9SMark F. Adams .seealso: MatCoarsenCreate(), MatCoarsen 1292b43b03e9SMark F. Adams J*/ 129319fd82e9SBarry Smith typedef const char* MatCoarsenType; 1294b43b03e9SMark F. Adams #define MATCOARSENMIS "mis" 12959057884aSMark F. Adams #define MATCOARSENHEM "hem" 1296b43b03e9SMark F. Adams 12970cbbd2e1SMark F. Adams /* linked list for aggregates */ 129841b27cdeSMark F. Adams typedef struct _PetscCDIntNd{ 129941b27cdeSMark F. Adams struct _PetscCDIntNd *next; 13000cbbd2e1SMark F. Adams PetscInt gid; 130141b27cdeSMark F. Adams }PetscCDIntNd; 13020cbbd2e1SMark F. Adams 13030cbbd2e1SMark F. Adams /* only used by node pool */ 130441b27cdeSMark F. Adams typedef struct _PetscCDArrNd{ 130541b27cdeSMark F. Adams struct _PetscCDArrNd *next; 130641b27cdeSMark F. Adams struct _PetscCDIntNd *array; 130741b27cdeSMark F. Adams }PetscCDArrNd; 13080cbbd2e1SMark F. Adams 13090cbbd2e1SMark F. Adams typedef struct _PetscCoarsenData{ 131082c86c8fSBarry Smith PetscCDArrNd pool_list; /* node pool */ 131141b27cdeSMark F. Adams PetscCDIntNd *new_node; 13120cbbd2e1SMark F. Adams PetscInt new_left; 13130cbbd2e1SMark F. Adams PetscInt chk_sz; 131441b27cdeSMark F. Adams PetscCDIntNd *extra_nodes; 131582c86c8fSBarry Smith PetscCDIntNd **array; /* Array of lists */ 13160cbbd2e1SMark F. Adams PetscInt size; 131782c86c8fSBarry Smith Mat mat; /* cache a Mat for communication data */ 13180cbbd2e1SMark F. Adams }PetscCoarsenData; 13190cbbd2e1SMark F. Adams 1320014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenCreate(MPI_Comm,MatCoarsen*); 132119fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenSetType(MatCoarsen,MatCoarsenType); 1322014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetAdjacency(MatCoarsen,Mat); 1323014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetGreedyOrdering(MatCoarsen,const IS); 1324014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetStrictAggs(MatCoarsen,PetscBool); 1325014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetVerbose(MatCoarsen,PetscInt); 1326014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenGetData( MatCoarsen, PetscCoarsenData ** ); 1327014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenApply(MatCoarsen); 1328014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenDestroy(MatCoarsen*); 1329b43b03e9SMark F. Adams 1330bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenRegister(const char[],PetscErrorCode (*)(MatCoarsen)); 1331b43b03e9SMark F. Adams 1332014dd563SJed Brown PETSC_EXTERN PetscBool MatCoarsenRegisterAllCalled; 1333b43b03e9SMark F. Adams 1334607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenRegisterAll(void); 1335b43b03e9SMark F. Adams 1336014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenView(MatCoarsen,PetscViewer); 1337014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetFromOptions(MatCoarsen); 133819fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenGetType(MatCoarsen,MatCoarsenType*); 1339ce1779c8SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatCoarsenViewFromOptions(MatCoarsen A,const char prefix[],const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,prefix,name);} 1340b43b03e9SMark F. Adams 1341014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMeshToVertexGraph(Mat,PetscInt,Mat*); 1342014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMeshToCellGraph(Mat,PetscInt,Mat*); 1343591294e4SBarry Smith 13440752156aSBarry Smith /* 13450a835dfdSSatish Balay If you add entries here you must also add them to finclude/petscmat.h 1346d4fbbf0eSBarry Smith */ 13471c1c02c0SLois Curfman McInnes typedef enum { MATOP_SET_VALUES=0, 13481c1c02c0SLois Curfman McInnes MATOP_GET_ROW=1, 13491c1c02c0SLois Curfman McInnes MATOP_RESTORE_ROW=2, 13501c1c02c0SLois Curfman McInnes MATOP_MULT=3, 13511c1c02c0SLois Curfman McInnes MATOP_MULT_ADD=4, 13527c922b88SBarry Smith MATOP_MULT_TRANSPOSE=5, 13537c922b88SBarry Smith MATOP_MULT_TRANSPOSE_ADD=6, 13541c1c02c0SLois Curfman McInnes MATOP_SOLVE=7, 13551c1c02c0SLois Curfman McInnes MATOP_SOLVE_ADD=8, 13567c922b88SBarry Smith MATOP_SOLVE_TRANSPOSE=9, 13577c922b88SBarry Smith MATOP_SOLVE_TRANSPOSE_ADD=10, 13581c1c02c0SLois Curfman McInnes MATOP_LUFACTOR=11, 13591c1c02c0SLois Curfman McInnes MATOP_CHOLESKYFACTOR=12, 1360a714c06dSBarry Smith MATOP_SOR=13, 13611c1c02c0SLois Curfman McInnes MATOP_TRANSPOSE=14, 13621c1c02c0SLois Curfman McInnes MATOP_GETINFO=15, 13631c1c02c0SLois Curfman McInnes MATOP_EQUAL=16, 13641c1c02c0SLois Curfman McInnes MATOP_GET_DIAGONAL=17, 13651c1c02c0SLois Curfman McInnes MATOP_DIAGONAL_SCALE=18, 13661c1c02c0SLois Curfman McInnes MATOP_NORM=19, 13671c1c02c0SLois Curfman McInnes MATOP_ASSEMBLY_BEGIN=20, 13681c1c02c0SLois Curfman McInnes MATOP_ASSEMBLY_END=21, 1369d519adbfSMatthew Knepley MATOP_SET_OPTION=22, 1370d519adbfSMatthew Knepley MATOP_ZERO_ENTRIES=23, 1371d519adbfSMatthew Knepley MATOP_ZERO_ROWS=24, 1372d519adbfSMatthew Knepley MATOP_LUFACTOR_SYMBOLIC=25, 1373d519adbfSMatthew Knepley MATOP_LUFACTOR_NUMERIC=26, 1374d519adbfSMatthew Knepley MATOP_CHOLESKY_FACTOR_SYMBOLIC=27, 1375d519adbfSMatthew Knepley MATOP_CHOLESKY_FACTOR_NUMERIC=28, 1376d519adbfSMatthew Knepley MATOP_SETUP_PREALLOCATION=29, 1377d519adbfSMatthew Knepley MATOP_ILUFACTOR_SYMBOLIC=30, 1378d519adbfSMatthew Knepley MATOP_ICCFACTOR_SYMBOLIC=31, 13799d39f69dSJed Brown /* MATOP_PLACEHOLDER_32=32, */ 13809d39f69dSJed Brown /* MATOP_PLACEHOLDER_33=33, */ 1381d519adbfSMatthew Knepley MATOP_DUPLICATE=34, 1382d519adbfSMatthew Knepley MATOP_FORWARD_SOLVE=35, 1383d519adbfSMatthew Knepley MATOP_BACKWARD_SOLVE=36, 1384d519adbfSMatthew Knepley MATOP_ILUFACTOR=37, 1385d519adbfSMatthew Knepley MATOP_ICCFACTOR=38, 1386d519adbfSMatthew Knepley MATOP_AXPY=39, 1387d519adbfSMatthew Knepley MATOP_GET_SUBMATRICES=40, 1388d519adbfSMatthew Knepley MATOP_INCREASE_OVERLAP=41, 1389d519adbfSMatthew Knepley MATOP_GET_VALUES=42, 1390d519adbfSMatthew Knepley MATOP_COPY=43, 1391d519adbfSMatthew Knepley MATOP_GET_ROW_MAX=44, 1392d519adbfSMatthew Knepley MATOP_SCALE=45, 1393d519adbfSMatthew Knepley MATOP_SHIFT=46, 139435153367SBarry Smith MATOP_DIAGONAL_SET=47, 13959d39f69dSJed Brown MATOP_ZERO_ROWS_COLUMNS=48, 13969d39f69dSJed Brown MATOP_SET_RANDOM=49, 1397d519adbfSMatthew Knepley MATOP_GET_ROW_IJ=50, 1398d519adbfSMatthew Knepley MATOP_RESTORE_ROW_IJ=51, 1399d519adbfSMatthew Knepley MATOP_GET_COLUMN_IJ=52, 1400d519adbfSMatthew Knepley MATOP_RESTORE_COLUMN_IJ=53, 1401d519adbfSMatthew Knepley MATOP_FDCOLORING_CREATE=54, 1402d519adbfSMatthew Knepley MATOP_COLORING_PATCH=55, 1403d519adbfSMatthew Knepley MATOP_SET_UNFACTORED=56, 1404d519adbfSMatthew Knepley MATOP_PERMUTE=57, 1405d519adbfSMatthew Knepley MATOP_SET_VALUES_BLOCKED=58, 1406d519adbfSMatthew Knepley MATOP_GET_SUBMATRIX=59, 1407d519adbfSMatthew Knepley MATOP_DESTROY=60, 1408d519adbfSMatthew Knepley MATOP_VIEW=61, 1409d519adbfSMatthew Knepley MATOP_CONVERT_FROM=62, 14107bab7c10SHong Zhang MATOP_MATMAT_MULT=63, 14117bab7c10SHong Zhang MATOP_MATMAT_MULT_SYMBOLIC=64, 14127bab7c10SHong Zhang MATOP_MATMAT_MULT_NUMERIC=65, 1413d519adbfSMatthew Knepley MATOP_SET_LOCAL_TO_GLOBAL_MAP=66, 1414d519adbfSMatthew Knepley MATOP_SET_VALUES_LOCAL=67, 1415d519adbfSMatthew Knepley MATOP_ZERO_ROWS_LOCAL=68, 1416d519adbfSMatthew Knepley MATOP_GET_ROW_MAX_ABS=69, 1417d519adbfSMatthew Knepley MATOP_GET_ROW_MIN_ABS=70, 1418d519adbfSMatthew Knepley MATOP_CONVERT=71, 1419d519adbfSMatthew Knepley MATOP_SET_COLORING=72, 14209d39f69dSJed Brown /* MATOP_PLACEHOLDER_73=73, */ 1421d519adbfSMatthew Knepley MATOP_SET_VALUES_ADIFOR=74, 1422d519adbfSMatthew Knepley MATOP_FD_COLORING_APPLY=75, 1423d519adbfSMatthew Knepley MATOP_SET_FROM_OPTIONS=76, 1424bda6c4cbSBarry Smith MATOP_MULT_CONSTRAINED=77, 1425bda6c4cbSBarry Smith MATOP_MULT_TRANSPOSE_CONSTRAIN=78, 14269d39f69dSJed Brown MATOP_FIND_ZERO_DIAGONALS=79, 1427d519adbfSMatthew Knepley MATOP_MULT_MULTIPLE=80, 1428d519adbfSMatthew Knepley MATOP_SOLVE_MULTIPLE=81, 1429d519adbfSMatthew Knepley MATOP_GET_INERTIA=82, 1430d519adbfSMatthew Knepley MATOP_LOAD=83, 1431d519adbfSMatthew Knepley MATOP_IS_SYMMETRIC=84, 1432d519adbfSMatthew Knepley MATOP_IS_HERMITIAN=85, 1433d519adbfSMatthew Knepley MATOP_IS_STRUCTURALLY_SYMMETRIC=86, 1434820469bcSHong Zhang MATOP_SET_VALUES_BLOCKEDLOCAL=87, 1435d519adbfSMatthew Knepley MATOP_GET_VECS=88, 1436d519adbfSMatthew Knepley MATOP_MAT_MULT=89, 1437d519adbfSMatthew Knepley MATOP_MAT_MULT_SYMBOLIC=90, 1438d519adbfSMatthew Knepley MATOP_MAT_MULT_NUMERIC=91, 1439d519adbfSMatthew Knepley MATOP_PTAP=92, 1440d519adbfSMatthew Knepley MATOP_PTAP_SYMBOLIC=93, 1441d519adbfSMatthew Knepley MATOP_PTAP_NUMERIC=94, 1442bda6c4cbSBarry Smith MATOP_MAT_TRANSPOSE_MULT=95, 1443bda6c4cbSBarry Smith MATOP_MAT_TRANSPOSE_MULT_SYMBO=96, 1444bda6c4cbSBarry Smith MATOP_MAT_TRANSPOSE_MULT_NUMER=97, 14459d39f69dSJed Brown /* MATOP_PLACEHOLDER_98=98, */ 14469d39f69dSJed Brown /* MATOP_PLACEHOLDER_99=99, */ 14479d39f69dSJed Brown /* MATOP_PLACEHOLDER_100=100, */ 14489d39f69dSJed Brown /* MATOP_PLACEHOLDER_101=101, */ 1449d519adbfSMatthew Knepley MATOP_CONJUGATE=102, 14509d39f69dSJed Brown /* MATOP_PLACEHOLDER_103=103, */ 1451d519adbfSMatthew Knepley MATOP_SET_VALUES_ROW=104, 1452d519adbfSMatthew Knepley MATOP_REAL_PART=105, 1453bda6c4cbSBarry Smith MATOP_IMAGINARY_PART=106, 1454bda6c4cbSBarry Smith MATOP_GET_ROW_UPPER_TRIANGULAR=107, 1455bda6c4cbSBarry Smith MATOP_RESTORE_ROW_UPPER_TRIANG=108, 1456bda6c4cbSBarry Smith MATOP_MAT_SOLVE=109, 14570e8d04f7SBarry Smith MATOP_GET_REDUNDANT_MATRIX=110, 1458d519adbfSMatthew Knepley MATOP_GET_ROW_MIN=111, 14590e8d04f7SBarry Smith MATOP_GET_COLUMN_VECTOR=112, 1460d519adbfSMatthew Knepley MATOP_MISSING_DIAGONAL=113, 14610e8d04f7SBarry Smith MATOP_GET_SEQ_NONZERO_STRUCTUR=114, 146289c6957cSBarry Smith MATOP_CREATE=115, 146389c6957cSBarry Smith MATOP_GET_GHOSTS=116, 14640e8d04f7SBarry Smith MATOP_GET_LOCAL_SUB_MATRIX=117, 14650e8d04f7SBarry Smith MATOP_RESTORE_LOCALSUB_MATRIX=118, 1466eeffb40dSHong Zhang MATOP_MULT_DIAGONAL_BLOCK=119, 14670e8d04f7SBarry Smith MATOP_HERMITIAN_TRANSPOSE=120, 14680e8d04f7SBarry Smith MATOP_MULT_HERMITIAN_TRANSPOSE=121, 14690e8d04f7SBarry Smith MATOP_MULT_HERMITIAN_TRANS_ADD=122, 14700e8d04f7SBarry Smith MATOP_GET_MULTI_PROC_BLOCK=123, 14719d39f69dSJed Brown MATOP_FIND_NONZERO_ROWS=124, 14720e8d04f7SBarry Smith MATOP_GET_COLUMN_NORMS=125, 14739d39f69dSJed Brown MATOP_INVERT_BLOCK_DIAGONAL=126, 14749d39f69dSJed Brown /* MATOP_PLACEHOLDER_127=127, */ 14750e8d04f7SBarry Smith MATOP_GET_SUB_MATRICES_PARALLE=128, 14762b8ad9a3SHong Zhang MATOP_SET_VALUES_BATCH=129, 14770e8d04f7SBarry Smith MATOP_TRANSPOSE_MAT_MULT=130, 14780e8d04f7SBarry Smith MATOP_TRANSPOSE_MAT_MULT_SYMBO=131, 14790e8d04f7SBarry Smith MATOP_TRANSPOSE_MAT_MULT_NUMER=132, 14800e8d04f7SBarry Smith MATOP_TRANSPOSE_COLORING_CREAT=133, 14810e8d04f7SBarry Smith MATOP_TRANS_COLORING_APPLY_SPT=134, 14820e8d04f7SBarry Smith MATOP_TRANS_COLORING_APPLY_DEN=135, 14830e8d04f7SBarry Smith MATOP_RART=136, 14840e8d04f7SBarry Smith MATOP_RART_SYMBOLIC=137, 14850e8d04f7SBarry Smith MATOP_RART_NUMERIC=138, 1486e09a3074SHong Zhang MATOP_SET_BLOCK_SIZES=139, 1487f9426fe0SMark Adams MATOP_AYPX=140, 14881919a2e2SJed Brown MATOP_RESIDUAL=141, 14891919a2e2SJed Brown MATOP_FDCOLORING_SETUP= 142 1490fae171e0SBarry Smith } MatOperation; 1491014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatHasOperation(Mat,MatOperation,PetscBool *); 1492014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellSetOperation(Mat,MatOperation,void(*)(void)); 1493014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellGetOperation(Mat,MatOperation,void(**)(void)); 1494014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellSetContext(Mat,void*); 1495112a2221SBarry Smith 149690ace30eSBarry Smith /* 149790ace30eSBarry Smith Codes for matrices stored on disk. By default they are 149890ace30eSBarry Smith stored in a universal format. By changing the format with 14997973ad04SBarry Smith PetscViewerSetFormat(viewer,PETSC_VIEWER_NATIVE); the matrices will 150090ace30eSBarry Smith be stored in a way natural for the matrix, for example dense matrices 150190ace30eSBarry Smith would be stored as dense. Matrices stored this way may only be 1502f4403165SShri Abhyankar read into matrices of the same type. 150390ace30eSBarry Smith */ 150490ace30eSBarry Smith #define MATRIX_BINARY_FORMAT_DENSE -1 150590ace30eSBarry Smith 1506014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetHashTableFactor(Mat,PetscReal); 1507014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatISGetLocalMat(Mat,Mat*); 1508014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatISSetLocalMat(Mat,Mat); 15091f4e1ec7SBarry Smith 1510d9274352SBarry Smith /*S 1511d9274352SBarry Smith MatNullSpace - Object that removes a null space from a vector, i.e. 1512d9274352SBarry Smith orthogonalizes the vector to a subsapce 1513d9274352SBarry Smith 1514f7a9e4ceSBarry Smith Level: advanced 1515d9274352SBarry Smith 1516d9274352SBarry Smith Concepts: matrix; linear operator, null space 1517d9274352SBarry Smith 15186e1639daSBarry Smith Users manual sections: 15196e1639daSBarry Smith . sec_singular 15206e1639daSBarry Smith 1521d9274352SBarry Smith .seealso: MatNullSpaceCreate() 1522d9274352SBarry Smith S*/ 152374637425SBarry Smith typedef struct _p_MatNullSpace* MatNullSpace; 1524d9274352SBarry Smith 1525014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceCreate(MPI_Comm,PetscBool ,PetscInt,const Vec[],MatNullSpace*); 1526014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(MatNullSpace,Vec,void*),void*); 1527014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceDestroy(MatNullSpace*); 1528d0195637SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceRemove(MatNullSpace,Vec); 1529014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNullSpace(Mat, MatNullSpace *); 1530014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNullSpace(Mat,MatNullSpace); 1531014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNearNullSpace(Mat,MatNullSpace); 1532014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNearNullSpace(Mat,MatNullSpace*); 1533014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceTest(MatNullSpace,Mat,PetscBool *); 1534014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceView(MatNullSpace,PetscViewer); 1535014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceGetVecs(MatNullSpace,PetscBool*,PetscInt*,const Vec**); 1536014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceCreateRigidBody(Vec,MatNullSpace*); 153774637425SBarry Smith 1538014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatReorderingSeqSBAIJ(Mat,IS); 1539014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetHashTableFactor(Mat,PetscReal); 1540014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetColumnIndices(Mat,PetscInt *); 1541014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJInvertBlockDiagonal(Mat); 15423f1d51d7SBarry Smith 1543014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMAIJ(Mat,PetscInt,Mat*); 1544014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMAIJRedimension(Mat,PetscInt,Mat*); 1545014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMAIJGetAIJ(Mat,Mat*); 1546c4f061fbSSatish Balay 1547014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatComputeExplicitOperator(Mat,Mat*); 1548b0a32e0cSBarry Smith 1549014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalScaleLocal(Mat,Vec); 155004f1ad80SBarry Smith 1551014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*); 1552014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetBase(Mat,Vec,Vec); 1553014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*); 1554014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*)); 1555014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec)); 1556014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDAddNullSpace(Mat,MatNullSpace); 1557014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt); 1558014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDResetHHistory(Mat); 1559014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctionError(Mat,PetscReal); 1560014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetPeriod(Mat,PetscInt); 1561014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDGetH(Mat,PetscScalar *); 1562014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetOptionsPrefix(Mat,const char[]); 1563014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*); 1564014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*); 1565e884886fSBarry Smith 15666370053bSBarry Smith /*S 15676370053bSBarry Smith MatMFFD - A data structured used to manage the computation of the h differencing parameter for matrix-free 15686370053bSBarry Smith Jacobian vector products 1569e884886fSBarry Smith 15706370053bSBarry Smith Notes: MATMFFD is a specific MatType which uses the MatMFFD data structure 15716370053bSBarry Smith 15726370053bSBarry Smith MatMFFD*() methods actually take the Mat as their first argument. Not a MatMFFD data structure 15736370053bSBarry Smith 15746370053bSBarry Smith Level: developer 15756370053bSBarry Smith 15766370053bSBarry Smith .seealso: MATMFFD, MatCreateMFFD(), MatMFFDSetFuction(), MatMFFDSetType(), MatMFFDRegister() 15776370053bSBarry Smith S*/ 1578e884886fSBarry Smith typedef struct _p_MatMFFD* MatMFFD; 1579e884886fSBarry Smith 158076bdecfbSBarry Smith /*J 1581e884886fSBarry Smith MatMFFDType - algorithm used to compute the h used in computing matrix-vector products via differencing of the function 1582e884886fSBarry Smith 1583e884886fSBarry Smith Level: beginner 1584e884886fSBarry Smith 1585e884886fSBarry Smith .seealso: MatMFFDSetType(), MatMFFDRegister() 158676bdecfbSBarry Smith J*/ 158719fd82e9SBarry Smith typedef const char* MatMFFDType; 1588e884886fSBarry Smith #define MATMFFD_DS "ds" 1589e884886fSBarry Smith #define MATMFFD_WP "wp" 1590e884886fSBarry Smith 159119fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDSetType(Mat,MatMFFDType); 1592bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDRegister(const char[],PetscErrorCode (*)(MatMFFD)); 1593e884886fSBarry Smith 1594607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDRegisterAll(void); 1595014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDDSSetUmin(Mat,PetscReal); 1596014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDWPSetComputeNormU(Mat,PetscBool ); 1597e884886fSBarry Smith 1598014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *); 1599014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *); 16007dbadf16SMatthew Knepley 160197969023SHong Zhang /* 160297969023SHong Zhang PETSc interface to MUMPS 160397969023SHong Zhang */ 160497969023SHong Zhang #ifdef PETSC_HAVE_MUMPS 1605014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMumpsSetIcntl(Mat,PetscInt,PetscInt); 1606bc6112feSHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetIcntl(Mat,PetscInt,PetscInt*); 1607b9e9ea26SSatish Balay PETSC_EXTERN PetscErrorCode MatMumpsSetCntl(Mat,PetscInt,PetscReal); 1608bc6112feSHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetCntl(Mat,PetscInt,PetscReal*); 1609bc6112feSHong Zhang 1610ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetInfo(Mat,PetscInt,PetscInt*); 1611ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetInfog(Mat,PetscInt,PetscInt*); 1612ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetRinfo(Mat,PetscInt,PetscReal*); 1613ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetRinfog(Mat,PetscInt,PetscReal*); 161497969023SHong Zhang #endif 161523a5497aSJed Brown 1616d954db57SHong Zhang /* 1617d615f992SSatish Balay PETSc interface to Mkl_Pardiso 1618b500a6b7SJose David Bermeol */ 1619b500a6b7SJose David Bermeol #ifdef PETSC_HAVE_MKL_PARDISO 1620d615f992SSatish Balay PETSC_EXTERN PetscErrorCode MatMkl_PardisoSetCntl(Mat,PetscInt,PetscInt); 1621b500a6b7SJose David Bermeol #endif 1622b500a6b7SJose David Bermeol 1623b500a6b7SJose David Bermeol /* 1624d954db57SHong Zhang PETSc interface to SUPERLU 1625d954db57SHong Zhang */ 1626d954db57SHong Zhang #ifdef PETSC_HAVE_SUPERLU 1627014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSuperluSetILUDropTol(Mat,PetscReal); 1628d954db57SHong Zhang #endif 1629d954db57SHong Zhang 1630aa372e3fSPaul Mullowney #ifdef PETSC_HAVE_CUDA 16313f9c0db1SPaul Mullowney /*E 1632e057df02SPaul Mullowney MatCUSPARSEStorageFormat - indicates the storage format for CUSPARSE (GPU) 16332692e278SPaul Mullowney matrices. 1634e057df02SPaul Mullowney 1635e057df02SPaul Mullowney Not Collective 1636e057df02SPaul Mullowney 16378468deeeSKarl Rupp + MAT_CUSPARSE_CSR - Compressed Sparse Row 16382692e278SPaul Mullowney . MAT_CUSPARSE_ELL - Ellpack (requires CUDA 4.2 or later). 16392692e278SPaul Mullowney - MAT_CUSPARSE_HYB - Hybrid, a combination of Ellpack and Coordinate format (requires CUDA 4.2 or later). 1640e057df02SPaul Mullowney 1641e057df02SPaul Mullowney Level: intermediate 1642e057df02SPaul Mullowney 1643e057df02SPaul Mullowney Any additions/changes here MUST also be made in include/finclude/petscmat.h 1644e057df02SPaul Mullowney 1645e057df02SPaul Mullowney .seealso: MatCUSPARSESetFormat(), MatCUSPARSEFormatOperation 1646e057df02SPaul Mullowney E*/ 1647e057df02SPaul Mullowney 16483f9c0db1SPaul Mullowney typedef enum {MAT_CUSPARSE_CSR, MAT_CUSPARSE_ELL, MAT_CUSPARSE_HYB} MatCUSPARSEStorageFormat; 1649e057df02SPaul Mullowney 1650e057df02SPaul Mullowney /* these will be strings associated with enumerated type defined above */ 1651e057df02SPaul Mullowney PETSC_EXTERN const char *const MatCUSPARSEStorageFormats[]; 1652e057df02SPaul Mullowney 16533f9c0db1SPaul Mullowney /*E 1654e057df02SPaul Mullowney MatCUSPARSEFormatOperation - indicates the operation of CUSPARSE (GPU) 16552692e278SPaul Mullowney matrices whose operation should use a particular storage format. 1656e057df02SPaul Mullowney 1657e057df02SPaul Mullowney Not Collective 1658e057df02SPaul Mullowney 16598468deeeSKarl Rupp + MAT_CUSPARSE_MULT_DIAG - sets the storage format for the diagonal matrix in the parallel MatMult 16608468deeeSKarl Rupp . MAT_CUSPARSE_MULT_OFFDIAG - sets the storage format for the offdiagonal matrix in the parallel MatMult 16618468deeeSKarl Rupp . MAT_CUSPARSE_MULT - sets the storage format for the entire matrix in the serial (single GPU) MatMult 16628468deeeSKarl Rupp - MAT_CUSPARSE_ALL - sets the storage format for all CUSPARSE (GPU) matrices 1663e057df02SPaul Mullowney 1664e057df02SPaul Mullowney Level: intermediate 1665e057df02SPaul Mullowney 1666e057df02SPaul Mullowney .seealso: MatCUSPARSESetFormat(), MatCUSPARSEStorageFormat 1667e057df02SPaul Mullowney E*/ 166836d62e41SPaul Mullowney typedef enum {MAT_CUSPARSE_MULT_DIAG, MAT_CUSPARSE_MULT_OFFDIAG, MAT_CUSPARSE_MULT, MAT_CUSPARSE_ALL} MatCUSPARSEFormatOperation; 1669e057df02SPaul Mullowney 167010e9db80SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 167110e9db80SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCreateAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 1672e057df02SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCUSPARSESetFormat(Mat,MatCUSPARSEFormatOperation,MatCUSPARSEStorageFormat); 16739ae82921SPaul Mullowney #endif 16749ae82921SPaul Mullowney 167590c53902SBarry Smith #if defined(PETSC_HAVE_CUSP) 1676014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCUSP(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 1677014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateAIJCUSP(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 1678e057df02SPaul Mullowney 16793f9c0db1SPaul Mullowney /*E 1680e057df02SPaul Mullowney MatCUSPStorageFormat - indicates the storage format for CUSP (GPU) 168136d62e41SPaul Mullowney matrices. 1682e057df02SPaul Mullowney 1683e057df02SPaul Mullowney Not Collective 1684e057df02SPaul Mullowney 16858468deeeSKarl Rupp + MAT_CUSP_CSR - Compressed Sparse Row 16868468deeeSKarl Rupp . MAT_CUSP_DIA - Diagonal 16878468deeeSKarl Rupp - MAT_CUSP_ELL - Ellpack 1688e057df02SPaul Mullowney 1689e057df02SPaul Mullowney Level: intermediate 1690e057df02SPaul Mullowney 1691e057df02SPaul Mullowney Any additions/changes here MUST also be made in include/finclude/petscmat.h 1692e057df02SPaul Mullowney 1693e057df02SPaul Mullowney .seealso: MatCUSPSetFormat(), MatCUSPFormatOperation 1694e057df02SPaul Mullowney E*/ 16953f9c0db1SPaul Mullowney typedef enum {MAT_CUSP_CSR, MAT_CUSP_DIA, MAT_CUSP_ELL} MatCUSPStorageFormat; 1696e057df02SPaul Mullowney 1697e057df02SPaul Mullowney /* these will be strings associated with enumerated type defined above */ 1698e057df02SPaul Mullowney PETSC_EXTERN const char *const MatCUSPStorageFormats[]; 1699e057df02SPaul Mullowney 17003f9c0db1SPaul Mullowney /*E 1701e057df02SPaul Mullowney MatCUSPFormatOperation - indicates the operation of CUSP (GPU) 170236d62e41SPaul Mullowney matrices whose operation should use a particular storage format. 1703e057df02SPaul Mullowney 1704e057df02SPaul Mullowney Not Collective 1705e057df02SPaul Mullowney 17068468deeeSKarl Rupp + MAT_CUSP_MULT_DIAG - sets the storage format for the diagonal matrix in the parallel MatMult 17078468deeeSKarl Rupp . MAT_CUSP_MULT_OFFDIAG - sets the storage format for the offdiagonal matrix in the parallel MatMult 17088468deeeSKarl Rupp . MAT_CUSP_MULT - sets the storage format for the entire matrix in the serial (single GPU) MatMult 17098468deeeSKarl Rupp - MAT_CUSP_ALL - sets the storage format for all CUSP (GPU) matrices 1710e057df02SPaul Mullowney 1711e057df02SPaul Mullowney Level: intermediate 1712e057df02SPaul Mullowney 1713e057df02SPaul Mullowney Any additions/changes here MUST also be made in include/finclude/petscmat.h 1714e057df02SPaul Mullowney 1715e057df02SPaul Mullowney .seealso: MatCUSPSetFormat(), MatCUSPStorageFormat 1716e057df02SPaul Mullowney E*/ 171736d62e41SPaul Mullowney typedef enum {MAT_CUSP_MULT_DIAG, MAT_CUSP_MULT_OFFDIAG, MAT_CUSP_MULT, MAT_CUSP_ALL} MatCUSPFormatOperation; 1718e057df02SPaul Mullowney 1719e057df02SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCUSPSetFormat(Mat,MatCUSPFormatOperation,MatCUSPStorageFormat); 1720e057df02SPaul Mullowney #endif 172190c53902SBarry Smith 1722d67ff14aSKarl Rupp #if defined(PETSC_HAVE_VIENNACL) 1723d67ff14aSKarl Rupp PETSC_EXTERN PetscErrorCode MatCreateSeqAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 1724d67ff14aSKarl Rupp PETSC_EXTERN PetscErrorCode MatCreateAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 1725d67ff14aSKarl Rupp #endif 1726d67ff14aSKarl Rupp 172754efbe56SHong Zhang /* 172854efbe56SHong Zhang PETSc interface to FFTW 172954efbe56SHong Zhang */ 173054efbe56SHong Zhang #if defined(PETSC_HAVE_FFTW) 1731014dd563SJed Brown PETSC_EXTERN PetscErrorCode VecScatterPetscToFFTW(Mat,Vec,Vec); 1732014dd563SJed Brown PETSC_EXTERN PetscErrorCode VecScatterFFTWToPetsc(Mat,Vec,Vec); 1733014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetVecsFFTW(Mat,Vec*,Vec*,Vec*); 173454efbe56SHong Zhang #endif 173554efbe56SHong Zhang 1736db31f6deSJed Brown #if defined(PETSC_HAVE_ELEMENTAL) 1737607a6623SBarry Smith PETSC_EXTERN PetscErrorCode PetscElementalInitializePackage(void); 1738db31f6deSJed Brown PETSC_EXTERN PetscErrorCode PetscElementalFinalizePackage(void); 1739db31f6deSJed Brown #endif 1740db31f6deSJed Brown 1741014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateNest(MPI_Comm,PetscInt,const IS[],PetscInt,const IS[],const Mat[],Mat*); 1742014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSize(Mat,PetscInt*,PetscInt*); 1743014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetISs(Mat,IS[],IS[]); 1744014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetLocalISs(Mat,IS[],IS[]); 1745014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSubMats(Mat,PetscInt*,PetscInt*,Mat***); 1746014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSubMat(Mat,PetscInt,PetscInt,Mat*); 174719fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatNestSetVecType(Mat,VecType); 1748014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestSetSubMats(Mat,PetscInt,const IS[],PetscInt,const IS[],const Mat[]); 1749014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestSetSubMat(Mat,PetscInt,PetscInt,Mat); 1750d8588912SDave May 17514325cce7SMatthew G Knepley PETSC_EXTERN PetscErrorCode MatChop(Mat,PetscReal); 1752e0a58c10SMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatComputeBandwidth(Mat,PetscReal,PetscInt*); 17534325cce7SMatthew G Knepley 175423a5497aSJed Brown #endif 1755