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 253ca39a21SBarry Smith .seealso: MatSetType(), Mat, MatSolverType, MatRegister() 2676bdecfbSBarry Smith J*/ 2719fd82e9SBarry Smith typedef const char* MatType; 28273d9f13SBarry Smith #define MATSAME "same" 295a11e1b2SBarry Smith #define MATMAIJ "maij" 30273d9f13SBarry Smith #define MATSEQMAIJ "seqmaij" 31273d9f13SBarry Smith #define MATMPIMAIJ "mpimaij" 32273d9f13SBarry Smith #define MATIS "is" 335a11e1b2SBarry Smith #define MATAIJ "aij" 34273d9f13SBarry Smith #define MATSEQAIJ "seqaij" 35273d9f13SBarry Smith #define MATMPIAIJ "mpiaij" 365a11e1b2SBarry Smith #define MATAIJCRL "aijcrl" 375a11e1b2SBarry Smith #define MATSEQAIJCRL "seqaijcrl" 385a11e1b2SBarry Smith #define MATMPIAIJCRL "mpiaijcrl" 399ae82921SPaul Mullowney #define MATAIJCUSPARSE "aijcusparse" 409ae82921SPaul Mullowney #define MATSEQAIJCUSPARSE "seqaijcusparse" 419ae82921SPaul Mullowney #define MATMPIAIJCUSPARSE "mpiaijcusparse" 42d67ff14aSKarl Rupp #define MATAIJVIENNACL "aijviennacl" 43d67ff14aSKarl Rupp #define MATSEQAIJVIENNACL "seqaijviennacl" 44d67ff14aSKarl Rupp #define MATMPIAIJVIENNACL "mpiaijviennacl" 455a11e1b2SBarry Smith #define MATAIJPERM "aijperm" 465a11e1b2SBarry Smith #define MATSEQAIJPERM "seqaijperm" 475a11e1b2SBarry Smith #define MATMPIAIJPERM "mpiaijperm" 484a2a386eSRichard Tran Mills #define MATAIJMKL "aijmkl" 494a2a386eSRichard Tran Mills #define MATSEQAIJMKL "seqaijmkl" 504a2a386eSRichard Tran Mills #define MATMPIAIJMKL "mpiaijmkl" 51b5b72c8aSIrina Sokolova #define MATBAIJMKL "baijmkl" 52b5b72c8aSIrina Sokolova #define MATSEQBAIJMKL "seqbaijmkl" 53b5b72c8aSIrina Sokolova #define MATMPIBAIJMKL "mpibaijmkl" 54273d9f13SBarry Smith #define MATSHELL "shell" 555a11e1b2SBarry Smith #define MATDENSE "dense" 56209238afSKris Buschelman #define MATSEQDENSE "seqdense" 57273d9f13SBarry Smith #define MATMPIDENSE "mpidense" 58db31f6deSJed Brown #define MATELEMENTAL "elemental" 595a11e1b2SBarry Smith #define MATBAIJ "baij" 60273d9f13SBarry Smith #define MATSEQBAIJ "seqbaij" 61273d9f13SBarry Smith #define MATMPIBAIJ "mpibaij" 62273d9f13SBarry Smith #define MATMPIADJ "mpiadj" 635a11e1b2SBarry Smith #define MATSBAIJ "sbaij" 64273d9f13SBarry Smith #define MATSEQSBAIJ "seqsbaij" 65273d9f13SBarry Smith #define MATMPISBAIJ "mpisbaij" 66cebc7f6cSBarry Smith #define MATDAAD "daad" 67cebc7f6cSBarry Smith #define MATMFFD "mffd" 68c8a8475eSBarry Smith #define MATNORMAL "normal" 69c5e4d11fSDmitry Karpeev #define MATNORMALHERMITIAN "normalh" 70ab92ecdeSBarry Smith #define MATLRC "lrc" 712a6744ebSBarry Smith #define MATSCATTER "scatter" 72421e10b8SBarry Smith #define MATBLOCKMAT "blockmat" 73793850ffSBarry Smith #define MATCOMPOSITE "composite" 741f8c7532SHong Zhang #define MATFFT "fft" 751f8c7532SHong Zhang #define MATFFTW "fftw" 76e133240eSMatthew G Knepley #define MATSEQCUFFT "seqcufft" 77557cca28SSatish Balay #define MATTRANSPOSEMAT "transpose" 7872ca8751SBarry Smith #define MATSCHURCOMPLEMENT "schurcomplement" 791d6018f0SLisandro Dalcin #define MATPYTHON "python" 8063c07aadSStefano Zampini #define MATHYPRE "hypre" 81f91d8e95SBarry Smith #define MATHYPRESTRUCT "hyprestruct" 82a9e6138eSGlenn Hammond #define MATHYPRESSTRUCT "hypresstruct" 83ee1cef2cSJed Brown #define MATSUBMATRIX "submatrix" 842c0dbf93SJed Brown #define MATLOCALREF "localref" 85d8588912SDave May #define MATNEST "nest" 86c094ef40SMatthew G. Knepley #define MATPREALLOCATOR "preallocator" 87d4002b98SHong Zhang #define MATSELL "sell" 88d4002b98SHong Zhang #define MATSEQSELL "seqsell" 89d4002b98SHong Zhang #define MATMPISELL "mpisell" 90a3b2e22bSHong Zhang #define MATDUMMY "dummy" 91773941d6SBarry Smith 9276bdecfbSBarry Smith /*J 933ca39a21SBarry Smith MatSolverType - String with the name of a PETSc matrix solver type. 949989ab13SBarry Smith 95c2b89b5dSBarry Smith For example: "petsc" indicates what PETSc provides, "superlu_dist" the parallel SuperLU_DIST package etc 969989ab13SBarry Smith 979989ab13SBarry Smith Level: beginner 989989ab13SBarry Smith 995c9eb25fSBarry Smith .seealso: MatGetFactor(), Mat, MatSetType(), MatType 10076bdecfbSBarry Smith J*/ 101ea799195SBarry Smith typedef const char* MatSolverType; 1022692d6eeSBarry Smith #define MATSOLVERSUPERLU "superlu" 1032692d6eeSBarry Smith #define MATSOLVERSUPERLU_DIST "superlu_dist" 10408f5efcfSPieter Ghysels #define MATSOLVERSTRUMPACK "strumpack" 1052692d6eeSBarry Smith #define MATSOLVERUMFPACK "umfpack" 10620db9a53SJed Brown #define MATSOLVERCHOLMOD "cholmod" 107d89f5e7aSBarry Smith #define MATSOLVERKLU "klu" 108418810c4SBarry Smith #define MATSOLVERSPARSEELEMENTAL "sparseelemental" 109d89f5e7aSBarry Smith #define MATSOLVERELEMENTAL "elemental" 1102692d6eeSBarry Smith #define MATSOLVERESSL "essl" 1112692d6eeSBarry Smith #define MATSOLVERLUSOL "lusol" 1122692d6eeSBarry Smith #define MATSOLVERMUMPS "mumps" 113d615f992SSatish Balay #define MATSOLVERMKL_PARDISO "mkl_pardiso" 114d305a81bSVasiliy Kozyrev #define MATSOLVERMKL_CPARDISO "mkl_cpardiso" 1152692d6eeSBarry Smith #define MATSOLVERPASTIX "pastix" 1162692d6eeSBarry Smith #define MATSOLVERMATLAB "matlab" 1172692d6eeSBarry Smith #define MATSOLVERPETSC "petsc" 1182692d6eeSBarry Smith #define MATSOLVERBAS "bas" 1199ae82921SPaul Mullowney #define MATSOLVERCUSPARSE "cusparse" 120c0cdd4a1SDahai Guo 121b24902e0SBarry Smith /*E 122b24902e0SBarry Smith MatFactorType - indicates what type of factorization is requested 123b24902e0SBarry Smith 124b24902e0SBarry Smith Level: beginner 125b24902e0SBarry Smith 126af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 127b24902e0SBarry Smith 1283ca39a21SBarry Smith .seealso: MatSolverType, 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 133ea799195SBarry Smith PETSC_EXTERN PetscErrorCode MatGetFactor(Mat,MatSolverType,MatFactorType,Mat*); 134ea799195SBarry Smith PETSC_EXTERN PetscErrorCode MatGetFactorAvailable(Mat,MatSolverType,MatFactorType,PetscBool *); 135ea799195SBarry Smith PETSC_EXTERN PetscErrorCode MatFactorGetSolverType(Mat,MatSolverType*); 136014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactorType(Mat,MatFactorType*); 137ea799195SBarry Smith PETSC_EXTERN PetscErrorCode MatSolverTypeRegister(MatSolverType,const MatType,MatFactorType,PetscErrorCode(*)(Mat,MatFactorType,Mat*)); 138ea799195SBarry Smith PETSC_EXTERN PetscErrorCode MatSolverTypeGet(MatSolverType,const MatType,MatFactorType,PetscBool*,PetscBool*,PetscErrorCode (**)(Mat,MatFactorType,Mat*)); 1399989ab13SBarry Smith 140c06d978dSMatthew Knepley /* Logging support */ 1410700a824SBarry Smith #define MAT_FILE_CLASSID 1211216 /* used to indicate matrices in binary files */ 142014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_CLASSID; 143335efc43SPeter Brune PETSC_EXTERN PetscClassId MAT_COLORING_CLASSID; 144014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_FDCOLORING_CLASSID; 145014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_TRANSPOSECOLORING_CLASSID; 146014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_PARTITIONING_CLASSID; 147014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_COARSEN_CLASSID; 148014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_NULLSPACE_CLASSID; 149014dd563SJed Brown PETSC_EXTERN PetscClassId MATMFFD_CLASSID; 150c06d978dSMatthew Knepley 151ceb03754SKris Buschelman /*E 1527dae84e0SHong Zhang MatReuse - Indicates if matrices obtained from a previous call to MatCreateSubMatrices(), MatCreateSubMatrix(), MatConvert() or several other functions 153cf37664fSBarry Smith are to be reused to store the new matrix values. 154cf37664fSBarry Smith 155cf37664fSBarry Smith $ MAT_INITIAL_MATRIX - create a new matrix 156cf37664fSBarry Smith $ MAT_REUSE_MATRIX - reuse the matrix created with a previous call that used MAT_INITIAL_MATRIX 157cf37664fSBarry Smith $ MAT_INPLACE_MATRIX - replace the first input matrix with the new matrix (not applicable to all functions) 158cf37664fSBarry Smith $ MAT_IGNORE_MATRIX - do not create a new matrix or reuse a give matrix, just ignore that matrix argument (not applicable to all functions) 159ceb03754SKris Buschelman 160ceb03754SKris Buschelman Level: beginner 161ceb03754SKris Buschelman 162af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 163ceb03754SKris Buschelman 1647dae84e0SHong Zhang .seealso: MatCreateSubMatrices(), MatCreateSubMatrix(), MatDestroyMatrices(), MatConvert() 165ceb03754SKris Buschelman E*/ 166511c6705SHong Zhang typedef enum {MAT_INITIAL_MATRIX,MAT_REUSE_MATRIX,MAT_IGNORE_MATRIX,MAT_INPLACE_MATRIX} MatReuse; 167ceb03754SKris Buschelman 1685494a064SHong Zhang /*E 1697dae84e0SHong Zhang MatCreateSubMatrixOption - Indicates if matrices obtained from a call to MatCreateSubMatrices() 170829201f2SHong Zhang include the matrix values. Currently it is only used by MatGetSeqNonzerostructure(). 1715494a064SHong Zhang 1725494a064SHong Zhang Level: beginner 1735494a064SHong Zhang 174829201f2SHong Zhang .seealso: MatGetSeqNonzerostructure() 1755494a064SHong Zhang E*/ 1767dae84e0SHong Zhang typedef enum {MAT_DO_NOT_GET_VALUES,MAT_GET_VALUES} MatCreateSubMatrixOption; 1775494a064SHong Zhang 178607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatInitializePackage(void); 179c06d978dSMatthew Knepley 180014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreate(MPI_Comm,Mat*); 181014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetSizes(Mat,PetscInt,PetscInt,PetscInt,PetscInt); 18219fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatSetType(Mat,MatType); 183014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetFromOptions(Mat); 184685405a1SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatViewFromOptions(Mat A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);} 185bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatRegister(const char[],PetscErrorCode(*)(Mat)); 186014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRegisterBaseName(const char[],const char[],const char[]); 187014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetOptionsPrefix(Mat,const char[]); 188014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAppendOptionsPrefix(Mat,const char[]); 189014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOptionsPrefix(Mat,const char*[]); 19084d44b13SHong Zhang PETSC_EXTERN PetscErrorCode MatSetErrorIfFailure(Mat,PetscBool); 191f69a0ea3SMatthew Knepley 192140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatList; 193140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatColoringList; 194140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatPartitioningList; 19528988994SBarry Smith 1963b224e63SBarry Smith /*E 197aa6c7ce3SBarry Smith MatStructure - Indicates if two matrices have the same nonzero structure 1983b224e63SBarry Smith 1993b224e63SBarry Smith Level: beginner 2003b224e63SBarry Smith 201af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 2023b224e63SBarry Smith 203aa6c7ce3SBarry Smith .seealso: MatCopy(), MatAXPY() 2043b224e63SBarry Smith E*/ 205aa6c7ce3SBarry Smith typedef enum {DIFFERENT_NONZERO_PATTERN,SUBSET_NONZERO_PATTERN,SAME_NONZERO_PATTERN} MatStructure; 2063b224e63SBarry Smith 2079779e05dSSatish Balay #if defined PETSC_HAVE_MKL_SPARSE 20880095d54SIrina Sokolova PETSC_EXTERN PetscErrorCode MatCreateBAIJMKL(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 20980095d54SIrina Sokolova PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJMKL(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[],Mat *A); 21080095d54SIrina Sokolova #endif 211d4002b98SHong Zhang PETSC_EXTERN PetscErrorCode MatCreateSeqSELL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 212d4002b98SHong Zhang PETSC_EXTERN PetscErrorCode MatCreateSELL(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 213d4002b98SHong Zhang PETSC_EXTERN PetscErrorCode MatSeqSELLSetPreallocation(Mat,PetscInt,const PetscInt[]); 214d4002b98SHong Zhang PETSC_EXTERN PetscErrorCode MatMPISELLSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 21580095d54SIrina Sokolova 216014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqDense(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*); 217014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateDense(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*); 218014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 219014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 220014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *); 221014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithSplitArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],PetscInt[],PetscInt[],PetscScalar[],Mat*); 2223b00a383SHong Zhang PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithSeqAIJ(MPI_Comm,Mat,Mat,const PetscInt[],Mat*); 2233b00a383SHong Zhang 224014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 225014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 226014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat*); 227d21a29f3SJed Brown 228014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAdj(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscInt[],Mat*); 229014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 230d21a29f3SJed Brown 231014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 232014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPISBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *); 23338f409ebSLisandro Dalcin PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 234014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 2354cce697cSJed Brown PETSC_EXTERN PetscErrorCode MatXAIJSetPreallocation(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscInt[],const PetscInt[]); 236dfb205c3SBarry Smith 237014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateShell(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,void *,Mat*); 238014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateNormal(Mat,Mat*); 239c5e4d11fSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatCreateNormalHermitian(Mat,Mat*); 240986c4d72SJose E. Roman PETSC_EXTERN PetscErrorCode MatCreateLRC(Mat,Mat,Vec,Mat,Mat*); 241267ca84cSJose E. Roman PETSC_EXTERN PetscErrorCode MatLRCGetMats(Mat,Mat*,Mat*,Vec*,Mat*); 242c5e4d11fSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatCreateIS(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,ISLocalToGlobalMapping,ISLocalToGlobalMapping,Mat*); 243014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 244014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 245c0cdd4a1SDahai Guo 246014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateScatter(MPI_Comm,VecScatter,Mat*); 247014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScatterSetVecScatter(Mat,VecScatter); 248014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScatterGetVecScatter(Mat,VecScatter*); 249014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateBlockMat(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,Mat*); 250014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeAddMat(Mat,Mat); 251014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeMerge(Mat); 252014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateComposite(MPI_Comm,PetscInt,const Mat*,Mat*); 2536d7c1e57SBarry Smith typedef enum {MAT_COMPOSITE_ADDITIVE,MAT_COMPOSITE_MULTIPLICATIVE} MatCompositeType; 254014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeSetType(Mat,MatCompositeType); 2556d7c1e57SBarry Smith 25619fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatCreateFFT(MPI_Comm,PetscInt,const PetscInt[],MatType,Mat*); 257014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqCUFFT(MPI_Comm,PetscInt,const PetscInt[],Mat*); 258dedccee8SHong Zhang 259014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateTranspose(Mat,Mat*); 2608060fb66Sstefano_zampini PETSC_EXTERN PetscErrorCode MatTransposeGetMat(Mat,Mat*); 261d0de2241SAndrew Spott PETSC_EXTERN PetscErrorCode MatCreateHermitianTranspose(Mat,Mat*); 26254e05e6cSHong Zhang PETSC_EXTERN PetscErrorCode MatCreateSubMatrixVirtual(Mat,IS,IS,Mat*); 26354e05e6cSHong Zhang PETSC_EXTERN PetscErrorCode MatSubMatrixVirtualUpdate(Mat,Mat,IS,IS); 264014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateLocalRef(Mat,IS,IS,Mat*); 2651472f72bSBarry Smith 266978814f1SStefano Zampini #if defined(PETSC_HAVE_HYPRE) 267d975228cSstefano_zampini PETSC_EXTERN PetscErrorCode MatHYPRESetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 268978814f1SStefano Zampini #endif 269978814f1SStefano Zampini 270014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPythonSetType(Mat,const char[]); 2711d6018f0SLisandro Dalcin 272846b4da1SFande Kong PETSC_EXTERN PetscErrorCode MatResetPreallocation(Mat); 273014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetUp(Mat); 274014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroy(Mat*); 275e56f5c9eSBarry Smith PETSC_EXTERN PetscErrorCode MatGetNonzeroState(Mat,PetscObjectState*); 27621c89e3eSBarry Smith 277014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatConjugate(Mat); 278014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRealPart(Mat); 279014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatImaginaryPart(Mat); 280014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetDiagonalBlock(Mat,Mat*); 281014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetTrace(Mat,PetscScalar*); 282713ccfa9SJed Brown PETSC_EXTERN PetscErrorCode MatInvertBlockDiagonal(Mat,const PetscScalar **); 283479a70c0SJed Brown PETSC_EXTERN PetscErrorCode MatInvertBlockDiagonalMat(Mat,Mat); 28499cafbc1SBarry Smith 2858ed539a5SBarry Smith /* ------------------------------------------------------------*/ 286014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 287014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 288014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesRow(Mat,PetscInt,const PetscScalar[]); 289014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesRowLocal(Mat,PetscInt,const PetscScalar[]); 290014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBatch(Mat,PetscInt,PetscInt,PetscInt[],const PetscScalar[]); 29173a71a0fSBarry Smith PETSC_EXTERN PetscErrorCode MatSetRandom(Mat,PetscRandom); 29284cb2905SBarry Smith 2932ef4de8bSBarry Smith /*S 2942ef4de8bSBarry Smith MatStencil - Data structure (C struct) for storing information about a single row or 29562ef0227SBarry Smith column of a matrix as indexed on an associated grid. These are arguments to MatSetStencil() and MatSetBlockStencil() 29662ef0227SBarry Smith 29762ef0227SBarry Smith The i,j, and k represent the logical coordinates over the entire grid (for 2 and 1 dimensional problems the k and j entries are ignored). 29862ef0227SBarry Smith The c represents the the degrees of freedom at each grid point (the dof argument to DMDASetDOF()). If dof is 1 then this entry is ignored. 29962ef0227SBarry Smith 30062ef0227SBarry Smith For stencil access to vectors see DMDAVecGetArray(), DMDAVecGetArrayF90(). 301be479b30SJed Brown 302be479b30SJed Brown Fortran usage is different, see MatSetValuesStencil() for details. 3032ef4de8bSBarry Smith 3042ef4de8bSBarry Smith Level: beginner 3052ef4de8bSBarry Smith 3062ef4de8bSBarry Smith Concepts: matrix; linear operator 3072ef4de8bSBarry Smith 30862ef0227SBarry Smith .seealso: MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockedStencil(), DMDAVecGetArray(), DMDAVecGetArrayF90() 3092ef4de8bSBarry Smith S*/ 310435da068SBarry Smith typedef struct { 311c1ac3661SBarry Smith PetscInt k,j,i,c; 312435da068SBarry Smith } MatStencil; 3132ef4de8bSBarry Smith 314014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode); 315014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode); 316014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt); 317435da068SBarry Smith 318d91e6319SBarry Smith /*E 319d91e6319SBarry Smith MatAssemblyType - Indicates if the matrix is now to be used, or if you plan 32062ef0227SBarry Smith to continue to add or insert values to it 321d91e6319SBarry Smith 322d91e6319SBarry Smith Level: beginner 323d91e6319SBarry Smith 324d91e6319SBarry Smith .seealso: MatAssemblyBegin(), MatAssemblyEnd() 325d91e6319SBarry Smith E*/ 3266d4a8577SBarry Smith typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType; 327014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssemblyBegin(Mat,MatAssemblyType); 328014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssemblyEnd(Mat,MatAssemblyType); 329014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssembled(Mat,PetscBool *); 3304f9c727eSBarry Smith 3311d73ed98SBarry Smith 33230de9b25SBarry Smith 333d91e6319SBarry Smith /*E 334d91e6319SBarry Smith MatOption - Options that may be set for a matrix and its behavior or storage 335d91e6319SBarry Smith 336d91e6319SBarry Smith Level: beginner 337d91e6319SBarry Smith 338af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 3390f8fb01aSBarry Smith Any additions/changes here must also be made in src/mat/interface/dlregismat.c in MatOptions[] 340d91e6319SBarry Smith 34195452b02SPatrick Sanan Developer Notes: 34295452b02SPatrick Sanan Entries that are negative need not be called collectively by all processes. 343382164b0SBarry Smith 344d91e6319SBarry Smith .seealso: MatSetOption() 345d91e6319SBarry Smith E*/ 346c5e4d11fSDmitry Karpeev typedef enum {MAT_OPTION_MIN = -3, 347c5e4d11fSDmitry Karpeev MAT_UNUSED_NONZERO_LOCATION_ERR = -2, 34892d4d306SBarry Smith MAT_ROW_ORIENTED = -1, 349382164b0SBarry Smith MAT_SYMMETRIC = 1, 350382164b0SBarry Smith MAT_STRUCTURALLY_SYMMETRIC = 2, 351382164b0SBarry Smith MAT_NEW_DIAGONALS = 3, 352382164b0SBarry Smith MAT_IGNORE_OFF_PROC_ENTRIES = 4, 353382164b0SBarry Smith MAT_USE_HASH_TABLE = 5, 354382164b0SBarry Smith MAT_KEEP_NONZERO_PATTERN = 6, 355382164b0SBarry Smith MAT_IGNORE_ZERO_ENTRIES = 7, 356382164b0SBarry Smith MAT_USE_INODES = 8, 357382164b0SBarry Smith MAT_HERMITIAN = 9, 358382164b0SBarry Smith MAT_SYMMETRY_ETERNAL = 10, 359c5e4d11fSDmitry Karpeev MAT_NEW_NONZERO_LOCATION_ERR = 11, 360382164b0SBarry Smith MAT_IGNORE_LOWER_TRIANGULAR = 12, 361382164b0SBarry Smith MAT_ERROR_LOWER_TRIANGULAR = 13, 362382164b0SBarry Smith MAT_GETROW_UPPERTRIANGULAR = 14, 363382164b0SBarry Smith MAT_SPD = 15, 36492d4d306SBarry Smith MAT_NO_OFF_PROC_ZERO_ROWS = 16, 36592d4d306SBarry Smith MAT_NO_OFF_PROC_ENTRIES = 17, 36692d4d306SBarry Smith MAT_NEW_NONZERO_LOCATIONS = 18, 367c5e4d11fSDmitry Karpeev MAT_NEW_NONZERO_ALLOCATION_ERR = 19, 368c5e4d11fSDmitry Karpeev MAT_SUBSET_OFF_PROC_ENTRIES = 20, 369c10200c1SHong Zhang MAT_SUBMAT_SINGLEIS = 21, 370957cac9fSHong Zhang MAT_STRUCTURE_ONLY = 22, 371957cac9fSHong Zhang MAT_OPTION_MAX = 23} MatOption; 372e057df02SPaul Mullowney 3730f8fb01aSBarry Smith PETSC_EXTERN const char *const *MatOptions; 374014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetOption(Mat,MatOption,PetscBool); 3757d68702bSBarry Smith PETSC_EXTERN PetscErrorCode MatGetOption(Mat,MatOption,PetscBool*); 37619fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatGetType(Mat,MatType*); 37784cb2905SBarry Smith 378014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]); 379014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 380014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 381014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowUpperTriangular(Mat); 382014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRowUpperTriangular(Mat); 383014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnVector(Mat,Vec,PetscInt); 3848c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJGetArray(Mat,PetscScalar *[]); 3858c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJRestoreArray(Mat,PetscScalar *[]); 38621e72a00SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJGetMaxRowNonzeros(Mat,PetscInt*); 387bd04181cSBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJSetValuesLocalFast(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 3884099cc6bSBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJSetType(Mat,MatType); 389388d47a6SSatish Balay PETSC_EXTERN PetscErrorCode MatSeqAIJRegister(const char[],PetscErrorCode (*)(Mat,MatType,MatReuse,Mat *)); 3904099cc6bSBarry Smith PETSC_EXTERN PetscFunctionList MatSeqAIJList; 3918397e458SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqSBAIJGetArray(Mat,PetscScalar *[]); 3928397e458SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqSBAIJRestoreArray(Mat,PetscScalar *[]); 3938c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatDenseGetArray(Mat,PetscScalar *[]); 3948c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatDenseRestoreArray(Mat,PetscScalar *[]); 395d3042a70SBarry Smith PETSC_EXTERN PetscErrorCode MatDensePlaceArray(Mat,const PetscScalar[]); 396d3042a70SBarry Smith PETSC_EXTERN PetscErrorCode MatDenseResetArray(Mat); 3978572280aSBarry Smith PETSC_EXTERN PetscErrorCode MatDenseGetArrayRead(Mat,const PetscScalar *[]); 3988572280aSBarry Smith PETSC_EXTERN PetscErrorCode MatDenseRestoreArrayRead(Mat,const PetscScalar *[]); 399014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBlockSize(Mat,PetscInt *); 400014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSize(Mat,PetscInt); 401014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBlockSizes(Mat,PetscInt *,PetscInt *); 402014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSizes(Mat,PetscInt,PetscInt); 40333d57670SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSizesFromMats(Mat,Mat,Mat); 40486aefd0dSHong Zhang PETSC_EXTERN PetscErrorCode MatDenseGetColumn(Mat,PetscInt,PetscScalar *[]); 40586aefd0dSHong Zhang PETSC_EXTERN PetscErrorCode MatDenseRestoreColumn(Mat,PetscScalar *[]); 4061620fd73SBarry Smith 407014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMult(Mat,Vec,Vec); 408014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultDiagonalBlock(Mat,Vec,Vec); 409014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultAdd(Mat,Vec,Vec,Vec); 410014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTranspose(Mat,Vec,Vec); 411014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultHermitianTranspose(Mat,Vec,Vec); 412014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsTranspose(Mat,Mat,PetscReal,PetscBool *); 413014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitianTranspose(Mat,Mat,PetscReal,PetscBool *); 414014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeAdd(Mat,Vec,Vec,Vec); 415014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultHermitianTransposeAdd(Mat,Vec,Vec,Vec); 416014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultConstrained(Mat,Vec,Vec); 417014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeConstrained(Mat,Vec,Vec); 418014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatSolve(Mat,Mat,Mat); 419bdc285e1SStefano Zampini PETSC_EXTERN PetscErrorCode MatMatSolveTranspose(Mat,Mat,Mat); 420f9426fe0SMark Adams PETSC_EXTERN PetscErrorCode MatResidual(Mat,Vec,Vec,Vec); 421f5edf698SHong Zhang 422d91e6319SBarry Smith /*E 423d91e6319SBarry Smith MatDuplicateOption - Indicates if a duplicated sparse matrix should have 424d91e6319SBarry Smith its numerical values copied over or just its nonzero structure. 425d91e6319SBarry Smith 426d91e6319SBarry Smith Level: beginner 427d91e6319SBarry Smith 428af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 429d91e6319SBarry Smith 4300018da39SRichard Tran Mills $ MAT_DO_NOT_COPY_VALUES - Create a matrix using the same nonzero pattern as the original matrix, 4310018da39SRichard Tran Mills $ with zeros for the numerical values. 4320018da39SRichard Tran Mills $ MAT_COPY_VALUES - Create a matrix with the same nonzero pattern as the original matrix 4330018da39SRichard Tran Mills $ and with the same numerical values. 4340018da39SRichard Tran Mills $ MAT_SHARE_NONZERO_PATTERN - Create a matrix that shares the nonzero structure with the previous matrix 4350018da39SRichard Tran Mills $ and does not copy it, using zeros for the numerical values. The parent and 4360018da39SRichard Tran Mills $ child matrices will share their index (i and j) arrays, and you cannot 4370018da39SRichard Tran Mills $ insert new nonzero entries into either matrix. 4380018da39SRichard Tran Mills 43995452b02SPatrick Sanan Notes: 44095452b02SPatrick Sanan Many matrix types (including SeqAIJ) do not support the MAT_SHARE_NONZERO_PATTERN optimization; in 4410018da39SRichard Tran Mills this case the behavior is as if MAT_DO_NOT_COPY_VALUES has been specified. 44270dcbbb9SBarry Smith 443d91e6319SBarry Smith .seealso: MatDuplicate() 444d91e6319SBarry Smith E*/ 44570dcbbb9SBarry Smith typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES,MAT_SHARE_NONZERO_PATTERN} MatDuplicateOption; 4462e8a6d31SBarry Smith 44719fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatConvert(Mat,MatType,MatReuse,Mat*); 448014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDuplicate(Mat,MatDuplicateOption,Mat*); 44994a9d846SBarry Smith 450cb5b572fSBarry Smith 451014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCopy(Mat,Mat,MatStructure); 452014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatView(Mat,PetscViewer); 453014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsSymmetric(Mat,PetscReal,PetscBool *); 454014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsStructurallySymmetric(Mat,PetscBool *); 455014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitian(Mat,PetscReal,PetscBool *); 456014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsSymmetricKnown(Mat,PetscBool *,PetscBool *); 457014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitianKnown(Mat,PetscBool *,PetscBool *); 458014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMissingDiagonal(Mat,PetscBool *,PetscInt *); 459014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLoad(Mat, PetscViewer); 4607b80b807SBarry Smith 4611a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool *); 4621a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool *); 4631a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool *); 4641a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool *); 465d4fbbf0eSBarry Smith 466d91e6319SBarry Smith /*S 467d91e6319SBarry Smith MatInfo - Context of matrix information, used with MatGetInfo() 468d91e6319SBarry Smith 469d91e6319SBarry Smith In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE 470d91e6319SBarry Smith 471d91e6319SBarry Smith Level: intermediate 472d91e6319SBarry Smith 473d91e6319SBarry Smith Concepts: matrix^nonzero information 474d91e6319SBarry Smith 475d9274352SBarry Smith .seealso: MatGetInfo(), MatInfoType 476d91e6319SBarry Smith S*/ 4774e220ebcSLois Curfman McInnes typedef struct { 478b0a32e0cSBarry Smith PetscLogDouble block_size; /* block size */ 479b0a32e0cSBarry Smith PetscLogDouble nz_allocated,nz_used,nz_unneeded; /* number of nonzeros */ 480b0a32e0cSBarry Smith PetscLogDouble memory; /* memory allocated */ 481b0a32e0cSBarry Smith PetscLogDouble assemblies; /* number of matrix assemblies called */ 482b0a32e0cSBarry Smith PetscLogDouble mallocs; /* number of mallocs during MatSetValues() */ 483b0a32e0cSBarry Smith PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */ 484b0a32e0cSBarry Smith PetscLogDouble factor_mallocs; /* number of mallocs during factorization */ 4854e220ebcSLois Curfman McInnes } MatInfo; 4864e220ebcSLois Curfman McInnes 487d9274352SBarry Smith /*E 488d9274352SBarry Smith MatInfoType - Indicates if you want information about the local part of the matrix, 489d9274352SBarry Smith the entire parallel matrix or the maximum over all the local parts. 490d9274352SBarry Smith 491d9274352SBarry Smith Level: beginner 492d9274352SBarry Smith 493af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 494d9274352SBarry Smith 495d9274352SBarry Smith .seealso: MatGetInfo(), MatInfo 496d9274352SBarry Smith E*/ 4977b80b807SBarry Smith typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType; 498014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetInfo(Mat,MatInfoType,MatInfo*); 499014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetDiagonal(Mat,Vec); 500014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMax(Mat,Vec,PetscInt[]); 501014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMin(Mat,Vec,PetscInt[]); 502014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMaxAbs(Mat,Vec,PetscInt[]); 503014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMinAbs(Mat,Vec,PetscInt[]); 504014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowSum(Mat,Vec); 505014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTranspose(Mat,MatReuse,Mat*); 506014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatHermitianTranspose(Mat,MatReuse,Mat*); 507014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPermute(Mat,IS,IS,Mat*); 508014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalScale(Mat,Vec,Vec); 509014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalSet(Mat,Vec,InsertMode); 510a52676ecSHong Zhang 511014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatEqual(Mat,Mat,PetscBool*); 512014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultEqual(Mat,Mat,PetscInt,PetscBool*); 513014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultAddEqual(Mat,Mat,PetscInt,PetscBool*); 514014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeEqual(Mat,Mat,PetscInt,PetscBool*); 515014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscBool*); 516a52676ecSHong Zhang PETSC_EXTERN PetscErrorCode MatMatMultEqual(Mat,Mat,Mat,PetscInt,PetscBool*); 517a52676ecSHong Zhang PETSC_EXTERN PetscErrorCode MatTransposeMatMultEqual(Mat,Mat,Mat,PetscInt,PetscBool*); 5187b80b807SBarry Smith 519014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNorm(Mat,NormType,PetscReal*); 520014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnNorms(Mat,NormType,PetscReal*); 521014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroEntries(Mat); 522014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 523014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsIS(Mat,IS,PetscScalar,Vec,Vec); 524014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsStencil(Mat,PetscInt,const MatStencil [],PetscScalar,Vec,Vec); 525014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsStencil(Mat,PetscInt,const MatStencil[],PetscScalar,Vec,Vec); 526014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumns(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 527014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsIS(Mat,IS,PetscScalar,Vec,Vec); 5287b80b807SBarry Smith 529014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSize(Mat,PetscInt*,PetscInt*); 530014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalSize(Mat,PetscInt*,PetscInt*); 531014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRange(Mat,PetscInt*,PetscInt*); 532014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRanges(Mat,const PetscInt**); 533014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRangeColumn(Mat,PetscInt*,PetscInt*); 534014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRangesColumn(Mat,const PetscInt**); 535566876eaSJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipIS(Mat,IS*,IS*); 5365ef9f2a5SBarry Smith 5377dae84e0SHong Zhang PETSC_EXTERN PetscErrorCode MatCreateSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]); 538cd2534ddSHong Zhang PETSC_DEPRECATED("Use MatCreateSubMatrices()") PETSC_STATIC_INLINE PetscErrorCode MatGetSubMatrices(Mat mat,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *submat[]) {return MatCreateSubMatrices(mat,n,irow,icol,scall,submat);} 5397dae84e0SHong Zhang PETSC_EXTERN PetscErrorCode MatCreateSubMatricesMPI(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]); 540cd2534ddSHong Zhang PETSC_DEPRECATED("Use MatCreateSubMatricesMPI()") PETSC_STATIC_INLINE PetscErrorCode MatGetSubMatricesMPI(Mat mat,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *submat[]) {return MatCreateSubMatricesMPI(mat,n,irow,icol,scall,submat);} 541014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroyMatrices(PetscInt,Mat *[]); 542df750dc8SHong Zhang PETSC_EXTERN PetscErrorCode MatDestroySubMatrices(PetscInt,Mat *[]); 5437dae84e0SHong Zhang PETSC_EXTERN PetscErrorCode MatCreateSubMatrix(Mat,IS,IS,MatReuse,Mat *); 544cd2534ddSHong Zhang PETSC_DEPRECATED("Use MatCreateSubMatrix()") PETSC_STATIC_INLINE PetscErrorCode MatGetSubMatrix(Mat mat,IS isrow,IS iscol,MatReuse cll,Mat *newmat) {return MatCreateSubMatrix(mat,isrow,iscol,cll,newmat);} 545014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalSubMatrix(Mat,IS,IS,Mat*); 546014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreLocalSubMatrix(Mat,IS,IS,Mat*); 547014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSeqNonzeroStructure(Mat,Mat*); 548014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroySeqNonzeroStructure(Mat*); 5497b80b807SBarry Smith 550014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJ(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*); 551014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJSymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*); 552014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJNumeric(Mat,Mat); 553014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMat(Mat,MatReuse,Mat*); 554014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*); 555014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,Mat*); 556014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetGhosts(Mat, PetscInt *,const PetscInt *[]); 5578efafbd8SBarry Smith 558014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt); 559aa1e27eaSFande Kong PETSC_EXTERN PetscErrorCode MatIncreaseOverlapSplit(Mat mat,PetscInt n,IS is[],PetscInt ov); 560d2b2a242SBarry Smith PETSC_EXTERN PetscErrorCode MatMPIAIJSetUseScalableIncreaseOverlap(Mat,PetscBool); 5617b80b807SBarry Smith 562014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*); 563014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMultSymbolic(Mat,Mat,PetscReal,Mat*); 564014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMultNumeric(Mat,Mat,Mat); 56522440eb1SKris Buschelman 5667bab7c10SHong Zhang PETSC_EXTERN PetscErrorCode MatMatMatMult(Mat,Mat,Mat,MatReuse,PetscReal,Mat*); 5672df6c5c3SPatrick Farrell PETSC_EXTERN PetscErrorCode MatGalerkin(Mat,Mat,Mat,MatReuse,PetscReal,Mat*); 5687bab7c10SHong Zhang 569014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*); 570014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAPSymbolic(Mat,Mat,PetscReal,Mat*); 571014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAPNumeric(Mat,Mat,Mat); 572014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARt(Mat,Mat,MatReuse,PetscReal,Mat*); 573014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARtSymbolic(Mat,Mat,PetscReal,Mat*); 574014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARtNumeric(Mat,Mat,Mat); 57522440eb1SKris Buschelman 576014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeMatMult(Mat,Mat,MatReuse,PetscReal,Mat*); 577014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposetMatMultSymbolic(Mat,Mat,PetscReal,Mat*); 578014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposetMatMultNumeric(Mat,Mat,Mat); 579014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMult(Mat,Mat,MatReuse,PetscReal,Mat*); 580014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMultSymbolic(Mat,Mat,PetscReal,Mat*); 581014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMultNumeric(Mat,Mat,Mat); 582bc011b1eSHong Zhang 583014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAXPY(Mat,PetscScalar,Mat,MatStructure); 584014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAYPX(Mat,PetscScalar,Mat,MatStructure); 5851c741599SHong Zhang 586014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScale(Mat,PetscScalar); 587014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShift(Mat,PetscScalar); 5887b80b807SBarry Smith 589014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping); 590014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*); 591a93ff8c4SPeter Brune PETSC_EXTERN PetscErrorCode MatGetLayouts(Mat,PetscLayout*,PetscLayout*); 592014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 593014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsLocalIS(Mat,IS,PetscScalar,Vec,Vec); 594014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 595014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocalIS(Mat,IS,PetscScalar,Vec,Vec); 596014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 597014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 598052efed2SBarry Smith 599014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStashSetInitialSize(Mat,PetscInt,PetscInt); 600014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*); 60190f02eecSBarry Smith 602014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInterpolate(Mat,Vec,Vec); 603014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInterpolateAdd(Mat,Vec,Vec,Vec); 604014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestrict(Mat,Vec,Vec); 6052a7a6963SBarry Smith PETSC_EXTERN PetscErrorCode MatCreateVecs(Mat,Vec*,Vec*); 60684cc2a5aSBarry Smith PETSC_DEPRECATED("Use MatCreateVecs()") PETSC_STATIC_INLINE PetscErrorCode MatGetVecs(Mat mat,Vec *x,Vec *y) {return MatCreateVecs(mat,x,y);} 60753cd1579SHong Zhang PETSC_EXTERN PetscErrorCode MatCreateRedundantMatrix(Mat,PetscInt,MPI_Comm,MatReuse,Mat*); 608014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetMultiProcBlock(Mat,MPI_Comm,MatReuse,Mat*); 609014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFindZeroDiagonals(Mat,IS*); 6103a062f41SBarry Smith PETSC_EXTERN PetscErrorCode MatFindOffBlockDiagonalEntries(Mat,IS*); 6119c8f2541SHong Zhang PETSC_EXTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat(MPI_Comm,Mat,PetscInt,MatReuse,Mat*); 612bd481603SBarry Smith 613bd481603SBarry Smith /*MC 6141620fd73SBarry Smith MatSetValue - Set a single entry into a matrix. 6151620fd73SBarry Smith 6161620fd73SBarry Smith Not collective 6171620fd73SBarry Smith 618a9834a7dSSatish Balay Synopsis: 619a9834a7dSSatish Balay #include <petscmat.h> 620a9834a7dSSatish Balay PetscErrorCode MatSetValue(Mat m,PetscInt row,PetscInt col,PetscScalar value,InsertMode mode) 621a9834a7dSSatish Balay 6221620fd73SBarry Smith Input Parameters: 6231620fd73SBarry Smith + m - the matrix 6241620fd73SBarry Smith . row - the row location of the entry 6251620fd73SBarry Smith . col - the column location of the entry 6261620fd73SBarry Smith . value - the value to insert 6271620fd73SBarry Smith - mode - either INSERT_VALUES or ADD_VALUES 6281620fd73SBarry Smith 6291620fd73SBarry Smith Notes: 6301620fd73SBarry Smith For efficiency one should use MatSetValues() and set several or many 6311620fd73SBarry Smith values simultaneously if possible. 6321620fd73SBarry Smith 6331620fd73SBarry Smith Level: beginner 6341620fd73SBarry Smith 6351620fd73SBarry Smith .seealso: MatSetValues(), MatSetValueLocal() 6361620fd73SBarry Smith M*/ 6371620fd73SBarry 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);} 6381620fd73SBarry Smith 6391620fd73SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatGetValue(Mat v,PetscInt i,PetscInt j,PetscScalar *va) {return MatGetValues(v,1,&i,1,&j,va);} 6401620fd73SBarry Smith 6411620fd73SBarry 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);} 6421620fd73SBarry Smith 6431620fd73SBarry Smith /*MC 6440d2a7ff2SSatish Balay MatPreallocateInitialize - Begins the block of code that will count the number of nonzeros per 645bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 646bd481603SBarry Smith 647bd481603SBarry Smith Synopsis: 648aaa7dc30SBarry Smith #include <petscmat.h> 649c1ac3661SBarry Smith PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz) 650bd481603SBarry Smith 651bd481603SBarry Smith Collective on MPI_Comm 652bd481603SBarry Smith 653bd481603SBarry Smith Input Parameters: 654bd481603SBarry Smith + comm - the communicator that will share the eventually allocated matrix 655859fcb39SBarry Smith . nrows - the number of LOCAL rows in the matrix 656859fcb39SBarry Smith - ncols - the number of LOCAL columns in the matrix 657bd481603SBarry Smith 658bd481603SBarry Smith Output Parameters: 659bd481603SBarry Smith + dnz - the array that will be passed to the matrix preallocation routines 660bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 661bd481603SBarry Smith 662bd481603SBarry Smith Level: intermediate 663bd481603SBarry Smith 664bd481603SBarry Smith Notes: 665a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 666bd481603SBarry Smith 6671d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 668bd481603SBarry Smith 669bd481603SBarry Smith Use MatPreallocateInitializeSymmetric() for symmetric matrices (MPISBAIJ matrices) 670bd481603SBarry Smith 6711620fd73SBarry Smith This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize(). 6721620fd73SBarry Smith 673bd481603SBarry Smith Concepts: preallocation^Matrix 674bd481603SBarry Smith 675d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(), 676d6e23781SBarry Smith MatPreallocateInitializeSymmetric(), MatPreallocateSymmetricSetLocalBlock() 677bd481603SBarry Smith M*/ 678c4f061fbSSatish Balay #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \ 6797c922b88SBarry Smith { \ 680a89bc211SBarry Smith PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__start,__end; \ 6817282cfc1SLisandro Dalcin _4_ierr = PetscCalloc2((size_t)__nrows,&dnz,(size_t)__nrows,&onz);CHKERRQ(_4_ierr); \ 6821795a4d1SJed Brown __start = 0; __end = __start; \ 683c1ac3661SBarry Smith _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ctmp;\ 684a89bc211SBarry Smith _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows; 6857c922b88SBarry Smith 686bd481603SBarry Smith /*MC 687bd481603SBarry Smith MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be 688bd481603SBarry Smith inserted using a local number of the rows and columns 689bd481603SBarry Smith 690bd481603SBarry Smith Synopsis: 691aaa7dc30SBarry Smith #include <petscmat.h> 692c1ac3661SBarry Smith PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 693bd481603SBarry Smith 694bd481603SBarry Smith Not Collective 695bd481603SBarry Smith 696bd481603SBarry Smith Input Parameters: 697784ac674SJed Brown + map - the row mapping from local numbering to global numbering 698bd481603SBarry Smith . nrows - the number of rows indicated 6991d73ed98SBarry Smith . rows - the indices of the rows 700784ac674SJed Brown . cmap - the column mapping from local to global numbering 701bd481603SBarry Smith . ncols - the number of columns in the matrix 702bd481603SBarry Smith . cols - the columns indicated 703bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 704bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 705bd481603SBarry Smith 706bd481603SBarry Smith Level: intermediate 707bd481603SBarry Smith 708bd481603SBarry Smith Notes: 709a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 710bd481603SBarry Smith 7111d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 712bd481603SBarry Smith 713bd481603SBarry Smith Concepts: preallocation^Matrix 714bd481603SBarry Smith 715d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(), 716c1154cd5SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocalRemoveDups() 717bd481603SBarry Smith M*/ 718784ac674SJed Brown #define MatPreallocateSetLocal(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \ 719c4f061fbSSatish Balay {\ 720c1ac3661SBarry Smith PetscInt __l;\ 721784ac674SJed Brown _4_ierr = ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\ 722784ac674SJed Brown _4_ierr = ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\ 723c4f061fbSSatish Balay for (__l=0;__l<nrows;__l++) {\ 724ef66eb69SBarry Smith _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 725c4f061fbSSatish Balay }\ 726c4f061fbSSatish Balay } 727c4f061fbSSatish Balay 728bd481603SBarry Smith /*MC 729c1154cd5SBarry Smith MatPreallocateSetLocalRemoveDups - Indicates the locations (rows and columns) in the matrix where nonzeros will be 730c1154cd5SBarry Smith inserted using a local number of the rows and columns. This version removes any duplicate columns in cols 731c1154cd5SBarry Smith 732c1154cd5SBarry Smith Synopsis: 733c1154cd5SBarry Smith #include <petscmat.h> 734c1154cd5SBarry Smith PetscErrorCode MatPreallocateSetLocalRemoveDups(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 735c1154cd5SBarry Smith 736c1154cd5SBarry Smith Not Collective 737c1154cd5SBarry Smith 738c1154cd5SBarry Smith Input Parameters: 739c1154cd5SBarry Smith + map - the row mapping from local numbering to global numbering 740c1154cd5SBarry Smith . nrows - the number of rows indicated 741c1154cd5SBarry Smith . rows - the indices of the rows (these values are mapped to the global values) 742c1154cd5SBarry Smith . cmap - the column mapping from local to global numbering 743c1154cd5SBarry Smith . ncols - the number of columns in the matrix (this value will be changed if duplicate columns are found) 744c1154cd5SBarry Smith . cols - the columns indicated (these values are mapped to the global values, they are then sorted and duplicates removed) 745c1154cd5SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 746c1154cd5SBarry Smith - ozn - the other array passed to the matrix preallocation routines 747c1154cd5SBarry Smith 748c1154cd5SBarry Smith Level: intermediate 749c1154cd5SBarry Smith 750c1154cd5SBarry Smith Notes: 751c1154cd5SBarry Smith See Users-Manual: ch_performance for more details. 752c1154cd5SBarry Smith 753c1154cd5SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 754c1154cd5SBarry Smith 755c1154cd5SBarry Smith Concepts: preallocation^Matrix 756c1154cd5SBarry Smith 757c1154cd5SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(), 758c1154cd5SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocal() 759c1154cd5SBarry Smith M*/ 760c1154cd5SBarry Smith #define MatPreallocateSetLocalRemoveDups(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \ 761c1154cd5SBarry Smith {\ 762c1154cd5SBarry Smith PetscInt __l;\ 763c1154cd5SBarry Smith _4_ierr = ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\ 764c1154cd5SBarry Smith _4_ierr = ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\ 765c1154cd5SBarry Smith _4_ierr = PetscSortRemoveDupsInt(&ncols,cols);CHKERRQ(_4_ierr);\ 766c1154cd5SBarry Smith for (__l=0;__l<nrows;__l++) {\ 767c1154cd5SBarry Smith _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 768c1154cd5SBarry Smith }\ 769c1154cd5SBarry Smith } 770c1154cd5SBarry Smith 771c1154cd5SBarry Smith /*MC 772d6e23781SBarry Smith MatPreallocateSetLocalBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be 773bd481603SBarry Smith inserted using a local number of the rows and columns 774bd481603SBarry Smith 775bd481603SBarry Smith Synopsis: 776aaa7dc30SBarry Smith #include <petscmat.h> 777d6e23781SBarry Smith PetscErrorCode MatPreallocateSetLocalBlock(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 778d6e23781SBarry Smith 779d6e23781SBarry Smith Not Collective 780d6e23781SBarry Smith 781d6e23781SBarry Smith Input Parameters: 782d6e23781SBarry Smith + map - the row mapping from local numbering to global numbering 783d6e23781SBarry Smith . nrows - the number of rows indicated 784d6e23781SBarry Smith . rows - the indices of the rows 785d6e23781SBarry Smith . cmap - the column mapping from local to global numbering 786d6e23781SBarry Smith . ncols - the number of columns in the matrix 787d6e23781SBarry Smith . cols - the columns indicated 788d6e23781SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 789d6e23781SBarry Smith - ozn - the other array passed to the matrix preallocation routines 790d6e23781SBarry Smith 791d6e23781SBarry Smith Level: intermediate 792d6e23781SBarry Smith 793d6e23781SBarry Smith Notes: 794d6e23781SBarry Smith See Users-Manual: ch_performance for more details. 795d6e23781SBarry Smith 796d6e23781SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 797d6e23781SBarry Smith 798d6e23781SBarry Smith Concepts: preallocation^Matrix 799d6e23781SBarry Smith 800d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(), 801d6e23781SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock() 802d6e23781SBarry Smith M*/ 803d6e23781SBarry Smith #define MatPreallocateSetLocalBlock(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \ 804d6e23781SBarry Smith {\ 805d6e23781SBarry Smith PetscInt __l;\ 806d6e23781SBarry Smith _4_ierr = ISLocalToGlobalMappingApplyBlock(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\ 807d6e23781SBarry Smith _4_ierr = ISLocalToGlobalMappingApplyBlock(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\ 808d6e23781SBarry Smith for (__l=0;__l<nrows;__l++) {\ 809d6e23781SBarry Smith _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 810d6e23781SBarry Smith }\ 811d6e23781SBarry Smith } 812d6e23781SBarry Smith 813d6e23781SBarry Smith /*MC 814d6e23781SBarry Smith MatPreallocateSymmetricSetLocalBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be 815d6e23781SBarry Smith inserted using a local number of the rows and columns 816d6e23781SBarry Smith 817d6e23781SBarry Smith Synopsis: 818d6e23781SBarry Smith #include <petscmat.h> 819d6e23781SBarry Smith PetscErrorCode MatPreallocateSymmetricSetLocalBlock(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 820bd481603SBarry Smith 821bd481603SBarry Smith Not Collective 822bd481603SBarry Smith 823bd481603SBarry Smith Input Parameters: 824bd481603SBarry Smith + map - the mapping between local numbering and global numbering 825bd481603SBarry Smith . nrows - the number of rows indicated 8261d73ed98SBarry Smith . rows - the indices of the rows 827bd481603SBarry Smith . ncols - the number of columns in the matrix 828bd481603SBarry Smith . cols - the columns indicated 829bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 830bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 831bd481603SBarry Smith 832bd481603SBarry Smith Level: intermediate 833bd481603SBarry Smith 834bd481603SBarry Smith Notes: 835a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 836bd481603SBarry Smith 837bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 838bd481603SBarry Smith 839bd481603SBarry Smith Concepts: preallocation^Matrix 840bd481603SBarry Smith 841d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateInitialize(), 842d6e23781SBarry Smith MatPreallocateInitialize(), MatPreallocateSetLocal() 843bd481603SBarry Smith M*/ 844d6e23781SBarry Smith #define MatPreallocateSymmetricSetLocalBlock(map,nrows,rows,ncols,cols,dnz,onz) 0;\ 845d3d32019SBarry Smith {\ 846c1ac3661SBarry Smith PetscInt __l;\ 847d6e23781SBarry Smith _4_ierr = ISLocalToGlobalMappingApplyBlock(map,nrows,rows,rows);CHKERRQ(_4_ierr);\ 848d6e23781SBarry Smith _4_ierr = ISLocalToGlobalMappingApplyBlock(map,ncols,cols,cols);CHKERRQ(_4_ierr);\ 849d3d32019SBarry Smith for (__l=0;__l<nrows;__l++) {\ 850d6e23781SBarry Smith _4_ierr = MatPreallocateSymmetricSetBlock((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 851d3d32019SBarry Smith }\ 852d3d32019SBarry Smith } 853bd481603SBarry Smith /*MC 854bd481603SBarry Smith MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be 855bd481603SBarry Smith inserted using a local number of the rows and columns 856bd481603SBarry Smith 857bd481603SBarry Smith Synopsis: 858aaa7dc30SBarry Smith #include <petscmat.h> 859c1ac3661SBarry Smith PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 860bd481603SBarry Smith 861bd481603SBarry Smith Not Collective 862bd481603SBarry Smith 863bd481603SBarry Smith Input Parameters: 86464075487SBarry Smith + row - the row 865bd481603SBarry Smith . ncols - the number of columns in the matrix 866a50f8bf6SHong Zhang - cols - the columns indicated 867a50f8bf6SHong Zhang 868a50f8bf6SHong Zhang Output Parameters: 869a50f8bf6SHong Zhang + dnz - the array that will be passed to the matrix preallocation routines 870bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 871bd481603SBarry Smith 872bd481603SBarry Smith Level: intermediate 873bd481603SBarry Smith 874bd481603SBarry Smith Notes: 875a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 876bd481603SBarry Smith 877bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 878bd481603SBarry Smith 8791620fd73SBarry Smith This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize(). 8801620fd73SBarry Smith 881bd481603SBarry Smith Concepts: preallocation^Matrix 882bd481603SBarry Smith 883d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(), 884d6e23781SBarry Smith MatPreallocateInitialize(), MatPreallocateSetLocal() 885bd481603SBarry Smith M*/ 886c4f061fbSSatish Balay #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\ 887c1ac3661SBarry Smith { PetscInt __i; \ 888e32f2f54SBarry 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);\ 889e32f2f54SBarry 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);\ 8907c922b88SBarry Smith for (__i=0; __i<nc; __i++) {\ 89164075487SBarry Smith if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \ 8927cc688d7SBarry Smith else dnz[row - __rstart]++;\ 8937c922b88SBarry Smith }\ 8947c922b88SBarry Smith } 8957c922b88SBarry Smith 896bd481603SBarry Smith /*MC 897d6e23781SBarry Smith MatPreallocateSymmetricSetBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be 898bd481603SBarry Smith inserted using a local number of the rows and columns 899bd481603SBarry Smith 900bd481603SBarry Smith Synopsis: 901aaa7dc30SBarry Smith #include <petscmat.h> 902d6e23781SBarry Smith PetscErrorCode MatPreallocateSymmetricSetBlock(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 903bd481603SBarry Smith 904bd481603SBarry Smith Not Collective 905bd481603SBarry Smith 906bd481603SBarry Smith Input Parameters: 907bd481603SBarry Smith + nrows - the number of rows indicated 9081d73ed98SBarry Smith . rows - the indices of the rows 909bd481603SBarry Smith . ncols - the number of columns in the matrix 910bd481603SBarry Smith . cols - the columns indicated 911bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 912bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 913bd481603SBarry Smith 914bd481603SBarry Smith Level: intermediate 915bd481603SBarry Smith 916bd481603SBarry Smith Notes: 917a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 918bd481603SBarry Smith 919bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 920bd481603SBarry Smith 9211620fd73SBarry Smith This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize(). 9221620fd73SBarry Smith 923bd481603SBarry Smith Concepts: preallocation^Matrix 924bd481603SBarry Smith 925d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateInitialize(), 926d6e23781SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocal() 927bd481603SBarry Smith M*/ 928d6e23781SBarry Smith #define MatPreallocateSymmetricSetBlock(row,nc,cols,dnz,onz) 0;\ 929c1ac3661SBarry Smith { PetscInt __i; \ 930d3d32019SBarry Smith for (__i=0; __i<nc; __i++) {\ 931d3d32019SBarry Smith if (cols[__i] >= __end) onz[row - __rstart]++; \ 932d3d32019SBarry Smith else if (cols[__i] >= row) dnz[row - __rstart]++;\ 933d3d32019SBarry Smith }\ 934d3d32019SBarry Smith } 935d3d32019SBarry Smith 936bd481603SBarry Smith /*MC 93716371a99SBarry Smith MatPreallocateLocation - An alternative to MatPreallocationSet() that puts the nonzero locations into the matrix if it exists 93816371a99SBarry Smith 93916371a99SBarry Smith Synopsis: 940aaa7dc30SBarry Smith #include <petscmat.h> 94116371a99SBarry Smith PetscErrorCode MatPreallocateLocations(Mat A,PetscInt row,PetscInt ncols,PetscInt *cols,PetscInt *dnz,PetscInt *onz) 94216371a99SBarry Smith 94316371a99SBarry Smith Not Collective 94416371a99SBarry Smith 94516371a99SBarry Smith Input Parameters: 94616371a99SBarry Smith . A - matrix 94716371a99SBarry Smith . row - row where values exist (must be local to this process) 94816371a99SBarry Smith . ncols - number of columns 94916371a99SBarry Smith . cols - columns with nonzeros 95016371a99SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 95116371a99SBarry Smith - ozn - the other array passed to the matrix preallocation routines 95216371a99SBarry Smith 95316371a99SBarry Smith Level: intermediate 95416371a99SBarry Smith 95516371a99SBarry Smith Notes: 956a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 95716371a99SBarry Smith 95816371a99SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 95916371a99SBarry Smith 96016371a99SBarry Smith This is a MACRO not a function because it uses a bunch of variables private to the MatPreallocation.... routines. 96116371a99SBarry Smith 96216371a99SBarry Smith Concepts: preallocation^Matrix 96316371a99SBarry Smith 964d6e23781SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(), 965d6e23781SBarry Smith MatPreallocateSymmetricSetLocalBlock() 96616371a99SBarry Smith M*/ 9670298fd71SBarry 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);} 96816371a99SBarry Smith 96916371a99SBarry Smith 97016371a99SBarry Smith /*MC 9710d2a7ff2SSatish Balay MatPreallocateFinalize - Ends the block of code that will count the number of nonzeros per 972bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 973bd481603SBarry Smith 974bd481603SBarry Smith Synopsis: 975aaa7dc30SBarry Smith #include <petscmat.h> 976c1ac3661SBarry Smith PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz) 977bd481603SBarry Smith 978bd481603SBarry Smith Collective on MPI_Comm 979bd481603SBarry Smith 980bd481603SBarry Smith Input Parameters: 98116371a99SBarry Smith + dnz - the array that was be passed to the matrix preallocation routines 982bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 983bd481603SBarry Smith 984bd481603SBarry Smith Level: intermediate 985bd481603SBarry Smith 986bd481603SBarry Smith Notes: 987a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 988bd481603SBarry Smith 989bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 990bd481603SBarry Smith 9911620fd73SBarry Smith This is a MACRO not a function because it closes the { started in MatPreallocateInitialize(). 9921620fd73SBarry Smith 993bd481603SBarry Smith Concepts: preallocation^Matrix 994bd481603SBarry Smith 995d6e23781SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(), 996d6e23781SBarry Smith MatPreallocateSymmetricSetLocalBlock() 997bd481603SBarry Smith M*/ 998a89bc211SBarry Smith #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree2(dnz,onz);CHKERRQ(_4_ierr);} 9997c922b88SBarry Smith 10007b80b807SBarry Smith /* Routines unique to particular data structures */ 1001014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellGetContext(Mat,void *); 1002698d4c6aSKris Buschelman 1003014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInodeAdjustForInodes(Mat,IS*,IS*); 1004014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *); 1005698d4c6aSKris Buschelman 1006014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetColumnIndices(Mat,PetscInt[]); 1007014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetColumnIndices(Mat,PetscInt[]); 1008014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 1009014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 1010014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 1011014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJFromTriple(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*,PetscInt,PetscBool); 10127b80b807SBarry Smith 1013a96a251dSBarry Smith #define MAT_SKIP_ALLOCATION -4 1014a96a251dSBarry Smith 1015014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]); 1016014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]); 1017014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]); 1018273d9f13SBarry Smith 1019014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 1020014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 1021014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 1022014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []); 1023014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 1024014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]); 1025014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 1026014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]); 1027b44e7856SBarry Smith PETSC_EXTERN PetscErrorCode MatMPIAdjToSeq(Mat,Mat*); 1028014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIDenseSetPreallocation(Mat,PetscScalar[]); 1029014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqDenseSetPreallocation(Mat,PetscScalar[]); 10309230625dSJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,const PetscInt*[]); 10319230625dSJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,const PetscInt*[]); 1032014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAdjCreateNonemptySubcommMat(Mat,Mat*); 1033273d9f13SBarry Smith 1034014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqDenseSetLDA(Mat,PetscInt); 1035014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDenseGetLocalMatrix(Mat,Mat*); 10361b807ce4Svictorle 1037014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStoreValues(Mat); 1038014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRetrieveValues(Mat); 10392e8a6d31SBarry Smith 1040014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDAADSetCtx(Mat,void*); 10413a7fca6bSBarry Smith 1042014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFindNonzeroRows(Mat,IS*); 10437cfe41e4SPatrick Farrell PETSC_EXTERN PetscErrorCode MatFindZeroRows(Mat,IS*); 10447b80b807SBarry Smith /* 10457b80b807SBarry Smith These routines are not usually accessed directly, rather solving is 104694b7f48cSBarry Smith done through the KSP and PC interfaces. 10477b80b807SBarry Smith */ 10487b80b807SBarry Smith 104976bdecfbSBarry Smith /*J 10508f6c3df8SBarry Smith MatOrderingType - String with the name of a PETSc matrix ordering 1051d9274352SBarry Smith 1052d9274352SBarry Smith Level: beginner 1053d9274352SBarry Smith 1054d9274352SBarry Smith .seealso: MatGetOrdering() 105576bdecfbSBarry Smith J*/ 105619fd82e9SBarry Smith typedef const char* MatOrderingType; 10572692d6eeSBarry Smith #define MATORDERINGNATURAL "natural" 10582692d6eeSBarry Smith #define MATORDERINGND "nd" 10592692d6eeSBarry Smith #define MATORDERING1WD "1wd" 10602692d6eeSBarry Smith #define MATORDERINGRCM "rcm" 10612692d6eeSBarry Smith #define MATORDERINGQMD "qmd" 10622692d6eeSBarry Smith #define MATORDERINGROWLENGTH "rowlength" 1063510184a7SJed Brown #define MATORDERINGWBM "wbm" 1064fc1bef75SJed Brown #define MATORDERINGSPECTRAL "spectral" 1065312e372aSBarry Smith #define MATORDERINGAMD "amd" /* only works if UMFPACK is installed with PETSc */ 1066b12f92e5SBarry Smith 106719fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatGetOrdering(Mat,MatOrderingType,IS*,IS*); 1068140e18c1SBarry Smith PETSC_EXTERN PetscErrorCode MatGetOrderingList(PetscFunctionList*); 1069bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatOrderingRegister(const char[],PetscErrorCode(*)(Mat,MatOrderingType,IS*,IS*)); 1070140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatOrderingList; 1071d4fbbf0eSBarry Smith 1072014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS); 1073fc1bef75SJed Brown PETSC_EXTERN PetscErrorCode MatCreateLaplacian(Mat,PetscReal,PetscBool,Mat*); 1074a2ce50c7SBarry Smith 1075d91e6319SBarry Smith /*S 1076d90ac83dSHong Zhang MatFactorShiftType - Numeric Shift. 1077d90ac83dSHong Zhang 1078d90ac83dSHong Zhang Level: beginner 1079d90ac83dSHong Zhang 1080d90ac83dSHong Zhang S*/ 1081d90ac83dSHong Zhang typedef enum {MAT_SHIFT_NONE,MAT_SHIFT_NONZERO,MAT_SHIFT_POSITIVE_DEFINITE,MAT_SHIFT_INBLOCKS} MatFactorShiftType; 10826a6fc655SJed Brown PETSC_EXTERN const char *const MatFactorShiftTypes[]; 10835e9742b9SJed Brown PETSC_EXTERN const char *const MatFactorShiftTypesDetail[]; 1084d90ac83dSHong Zhang 1085d90ac83dSHong Zhang /*S 10869aa7eafdSHong Zhang MatFactorError - indicates what type of error in matrix factor 10879aa7eafdSHong Zhang 10889aa7eafdSHong Zhang Level: beginner 10899aa7eafdSHong Zhang 109095452b02SPatrick Sanan Developer Notes: 109195452b02SPatrick Sanan Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 109200e125f8SBarry Smith 109300e125f8SBarry Smith .seealso: MatGetFactor() 10949aa7eafdSHong Zhang S*/ 10955cd7cf9dSHong Zhang typedef enum {MAT_FACTOR_NOERROR,MAT_FACTOR_STRUCT_ZEROPIVOT,MAT_FACTOR_NUMERIC_ZEROPIVOT,MAT_FACTOR_OUTMEMORY,MAT_FACTOR_OTHER} MatFactorError; 10969aa7eafdSHong Zhang 109700e125f8SBarry Smith PETSC_EXTERN PetscErrorCode MatFactorGetError(Mat,MatFactorError*); 1098b8b68cfdSBarry Smith PETSC_EXTERN PetscErrorCode MatFactorClearError(Mat); 10997b6c816cSBarry Smith PETSC_EXTERN PetscErrorCode MatFactorGetErrorZeroPivot(Mat,PetscReal*,PetscInt*); 110000e125f8SBarry Smith 11019aa7eafdSHong Zhang /*S 1102422a814eSBarry Smith MatFactorInfo - Data passed into the matrix factorization routines, and information about the resulting factorization 1103b00f7748SHong Zhang 110461cad860SBarry Smith In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE, that is use 110561cad860SBarry Smith $ MatFactorInfo info(MAT_FACTORINFO_SIZE) 1106b00f7748SHong Zhang 110795452b02SPatrick Sanan Notes: 110895452b02SPatrick Sanan These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC. 1109b00f7748SHong Zhang 111061cad860SBarry Smith You can use MatFactorInfoInitialize() to set default values. 111161cad860SBarry Smith 1112b00f7748SHong Zhang Level: developer 1113b00f7748SHong Zhang 1114d7d82daaSBarry Smith .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(), 1115d7d82daaSBarry Smith MatFactorInfoInitialize() 1116b00f7748SHong Zhang 1117b00f7748SHong Zhang S*/ 1118b00f7748SHong Zhang typedef struct { 111915e8a5b3SHong Zhang PetscReal diagonal_fill; /* force diagonal to fill in if initially not filled */ 112085317021SBarry Smith PetscReal usedt; 112115e8a5b3SHong Zhang PetscReal dt; /* drop tolerance */ 1122b00f7748SHong Zhang PetscReal dtcol; /* tolerance for pivoting */ 112315e8a5b3SHong Zhang PetscReal dtcount; /* maximum nonzeros to be allowed per row */ 112467e28bfeSBarry Smith PetscReal fill; /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */ 1125348344bbSBarry Smith PetscReal levels; /* ICC/ILU(levels) */ 1126bcd9e38bSBarry Smith PetscReal pivotinblocks; /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0 1127bcd9e38bSBarry Smith factorization may be faster if do not pivot */ 112815e8a5b3SHong Zhang PetscReal zeropivot; /* pivot is called zero if less than this */ 1129f4db908eSBarry Smith PetscReal shifttype; /* type of shift added to matrix factor to prevent zero pivots */ 1130f4db908eSBarry Smith PetscReal shiftamount; /* how large the shift is */ 113115e8a5b3SHong Zhang } MatFactorInfo; 1132ffa6d0a5SLois Curfman McInnes 1133014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFactorInfoInitialize(MatFactorInfo*); 11349a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatCholeskyFactor(Mat,IS,const MatFactorInfo*); 11359a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatCholeskyFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*); 11369a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatCholeskyFactorNumeric(Mat,Mat,const MatFactorInfo*); 11379a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatLUFactor(Mat,IS,IS,const MatFactorInfo*); 11389a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatILUFactor(Mat,IS,IS,const MatFactorInfo*); 11399a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatLUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*); 11409a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatILUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*); 11419a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatICCFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*); 11429a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatICCFactor(Mat,IS,const MatFactorInfo*); 11439a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatLUFactorNumeric(Mat,Mat,const MatFactorInfo*); 1144014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*); 1145014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolve(Mat,Vec,Vec); 1146014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatForwardSolve(Mat,Vec,Vec); 1147014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatBackwardSolve(Mat,Vec,Vec); 1148014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveAdd(Mat,Vec,Vec,Vec); 1149014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveTranspose(Mat,Vec,Vec); 1150014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveTransposeAdd(Mat,Vec,Vec,Vec); 1151014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolves(Mat,Vecs,Vecs); 11527c2f51b8SStefano Zampini PETSC_EXTERN PetscErrorCode MatSetUnfactored(Mat); 11537c2f51b8SStefano Zampini 11547c2f51b8SStefano Zampini typedef enum {MAT_FACTOR_SCHUR_UNFACTORED, MAT_FACTOR_SCHUR_FACTORED, MAT_FACTOR_SCHUR_INVERTED} MatFactorSchurStatus; 11558e7ba810SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorSetSchurIS(Mat,IS); 11567c2f51b8SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorGetSchurComplement(Mat,Mat*,MatFactorSchurStatus*); 11577c2f51b8SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorRestoreSchurComplement(Mat,Mat*,MatFactorSchurStatus); 11585a05ddb0SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorInvertSchurComplement(Mat); 11597c2f51b8SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorCreateSchurComplement(Mat,Mat*,MatFactorSchurStatus*); 11605a05ddb0SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorSolveSchurComplement(Mat,Vec,Vec); 11615a05ddb0SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorSolveSchurComplementTranspose(Mat,Vec,Vec); 11626dba178dSStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorFactorizeSchurComplement(Mat); 1163bb5a7306SBarry Smith 1164d91e6319SBarry Smith /*E 1165d91e6319SBarry Smith MatSORType - What type of (S)SOR to perform 1166bb1eb677SSatish Balay 1167d91e6319SBarry Smith Level: beginner 1168d91e6319SBarry Smith 1169d9274352SBarry Smith May be bitwise ORd together 1170d9274352SBarry Smith 1171af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 1172d91e6319SBarry Smith 11734e7234bfSBarry Smith MatSORType may be bitwise ORd together, so do not change the numbers 11744e7234bfSBarry Smith 117541f059aeSBarry Smith .seealso: MatSOR() 1176d91e6319SBarry Smith E*/ 1177ee50ffe9SBarry Smith typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3, 1178ee50ffe9SBarry Smith SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8, 1179ee50ffe9SBarry Smith SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16, 118084cb2905SBarry Smith SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType; 1181014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSOR(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec); 11828ed539a5SBarry Smith 1183d4fbbf0eSBarry Smith /* 1184639f9d9dSBarry Smith These routines are for efficiently computing Jacobians via finite differences. 1185639f9d9dSBarry Smith */ 1186b12f92e5SBarry Smith 11877086a01eSPeter Brune /*S 11887086a01eSPeter Brune MatColoring - Object for managing the coloring of matrices. 11897086a01eSPeter Brune 11907086a01eSPeter Brune Level: beginner 11917086a01eSPeter Brune 11927086a01eSPeter Brune Concepts: matrix, coloring 11937086a01eSPeter Brune 11947086a01eSPeter Brune .seealso: MatFDColoringCreate() ISColoring MatFDColoring 11957086a01eSPeter Brune S*/ 11967086a01eSPeter Brune typedef struct _p_MatColoring* MatColoring; 119776bdecfbSBarry Smith /*J 11988f6c3df8SBarry Smith MatColoringType - String with the name of a PETSc matrix coloring 1199d9274352SBarry Smith 1200d9274352SBarry Smith Level: beginner 1201d9274352SBarry Smith 12027086a01eSPeter Brune .seealso: MatColoringSetType(), MatColoring 120376bdecfbSBarry Smith J*/ 12047086a01eSPeter Brune 120519fd82e9SBarry Smith typedef const char* MatColoringType; 12065567d71aSPeter Brune #define MATCOLORINGJP "jp" 1207b9acb617SPeter Brune #define MATCOLORINGPOWER "power" 12082692d6eeSBarry Smith #define MATCOLORINGNATURAL "natural" 12092692d6eeSBarry Smith #define MATCOLORINGSL "sl" 12102692d6eeSBarry Smith #define MATCOLORINGLF "lf" 12112692d6eeSBarry Smith #define MATCOLORINGID "id" 12128121bdceSPeter Brune #define MATCOLORINGGREEDY "greedy" 1213b12f92e5SBarry Smith 12148ac6da0aSPeter Brune /*E 12158ac6da0aSPeter Brune MatColoringWeightType - Type of weight scheme 12168ac6da0aSPeter Brune 12178ac6da0aSPeter Brune Not Collective 12188ac6da0aSPeter Brune 12198ac6da0aSPeter Brune + MAT_COLORING_RANDOM - Random weights 12208ac6da0aSPeter Brune . MAT_COLORING_LEXICAL - Lexical weighting based upon global numbering. 12218ac6da0aSPeter Brune - MAT_COLORING_LF - Last-first weighting. 12228ac6da0aSPeter Brune 12238ac6da0aSPeter Brune Level: intermediate 12248ac6da0aSPeter Brune 1225af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 12268ac6da0aSPeter Brune E*/ 1227301c140bSPeter Brune typedef enum {MAT_COLORING_WEIGHT_RANDOM,MAT_COLORING_WEIGHT_LEXICAL,MAT_COLORING_WEIGHT_LF,MAT_COLORING_WEIGHT_SL} MatColoringWeightType; 12288ac6da0aSPeter Brune 1229335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringCreate(Mat,MatColoring*); 1230d7f2a307SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringGetDegrees(Mat,PetscInt,PetscInt*); 1231335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringDestroy(MatColoring*); 1232335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringView(MatColoring,PetscViewer); 1233335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetType(MatColoring,MatColoringType); 1234335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetFromOptions(MatColoring); 1235335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetDistance(MatColoring,PetscInt); 1236335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringGetDistance(MatColoring,PetscInt*); 1237335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetMaxColors(MatColoring,PetscInt); 1238335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringGetMaxColors(MatColoring,PetscInt*); 1239335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringApply(MatColoring,ISColoring*); 1240335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringRegister(const char[],PetscErrorCode(*)(MatColoring)); 1241014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*); 12428ac6da0aSPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetWeightType(MatColoring,MatColoringWeightType); 1243c29971fcSPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetWeights(MatColoring,PetscReal*,PetscInt*); 12448ac6da0aSPeter Brune PETSC_EXTERN PetscErrorCode MatColoringCreateWeights(MatColoring,PetscReal **,PetscInt **lperm); 1245bc76d17bSHong Zhang PETSC_EXTERN PetscErrorCode MatColoringTest(MatColoring,ISColoring); 1246bc76d17bSHong Zhang PETSC_EXTERN PetscErrorCode MatISColoringTest(Mat,ISColoring); 1247639f9d9dSBarry Smith 1248d9274352SBarry Smith /*S 1249d9274352SBarry Smith MatFDColoring - Object for computing a sparse Jacobian via finite differences 1250d9274352SBarry Smith and coloring 1251639f9d9dSBarry Smith 1252d9274352SBarry Smith Level: beginner 1253d9274352SBarry Smith 1254d9274352SBarry Smith Concepts: coloring, sparse Jacobian, finite differences 1255d9274352SBarry Smith 1256d9274352SBarry Smith .seealso: MatFDColoringCreate() 1257d9274352SBarry Smith S*/ 1258e2a1c21fSSatish Balay typedef struct _p_MatFDColoring* MatFDColoring; 1259639f9d9dSBarry Smith 1260014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringCreate(Mat,ISColoring,MatFDColoring *); 1261014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringDestroy(MatFDColoring*); 1262014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringView(MatFDColoring,PetscViewer); 1263014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*); 1264014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**); 1265014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal); 1266014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetFromOptions(MatFDColoring); 1267d1e9a80fSBarry Smith PETSC_EXTERN PetscErrorCode MatFDColoringApply(Mat,MatFDColoring,Vec,void *); 1268014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetF(MatFDColoring,Vec); 1269edaa7c33SBarry Smith PETSC_EXTERN PetscErrorCode MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,const PetscInt*[]); 1270f86b9fbaSHong Zhang PETSC_EXTERN PetscErrorCode MatFDColoringSetUp(Mat,ISColoring,MatFDColoring); 1271f86b9fbaSHong Zhang PETSC_EXTERN PetscErrorCode MatFDColoringSetBlockSize(MatFDColoring,PetscInt,PetscInt); 1272f86b9fbaSHong Zhang 1273b1683b59SHong Zhang 1274b1683b59SHong Zhang /*S 1275b9af6bddSHong Zhang MatTransposeColoring - Object for computing a sparse matrix product C=A*B^T via coloring 1276b1683b59SHong Zhang 1277b1683b59SHong Zhang Level: beginner 1278b1683b59SHong Zhang 1279b1683b59SHong Zhang Concepts: coloring, sparse matrix product 1280b1683b59SHong Zhang 1281b9af6bddSHong Zhang .seealso: MatTransposeColoringCreate() 1282b1683b59SHong Zhang S*/ 1283b9af6bddSHong Zhang typedef struct _p_MatTransposeColoring* MatTransposeColoring; 1284b1683b59SHong Zhang 1285014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeColoringCreate(Mat,ISColoring,MatTransposeColoring *); 1286014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransColoringApplySpToDen(MatTransposeColoring,Mat,Mat); 1287014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransColoringApplyDenToSp(MatTransposeColoring,Mat,Mat); 1288014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeColoringDestroy(MatTransposeColoring*); 1289b1683b59SHong Zhang 1290639f9d9dSBarry Smith /* 12910752156aSBarry Smith These routines are for partitioning matrices: currently used only 12923eda8832SBarry Smith for adjacency matrix, MatCreateMPIAdj(). 12930752156aSBarry Smith */ 1294ca161407SBarry Smith 1295d9274352SBarry Smith /*S 1296d9274352SBarry Smith MatPartitioning - Object for managing the partitioning of a matrix or graph 1297d9274352SBarry Smith 1298d9274352SBarry Smith Level: beginner 1299d9274352SBarry Smith 1300d9274352SBarry Smith Concepts: partitioning 1301d9274352SBarry Smith 1302743a1026Svictorle .seealso: MatPartitioningCreate(), MatPartitioningType 1303d9274352SBarry Smith S*/ 130491e9ee9fSBarry Smith typedef struct _p_MatPartitioning* MatPartitioning; 1305d9274352SBarry Smith 130676bdecfbSBarry Smith /*J 13078f6c3df8SBarry Smith MatPartitioningType - String with the name of a PETSc matrix partitioning 1308d9274352SBarry Smith 1309d9274352SBarry Smith Level: beginner 13100d04baf8SBarry Smith dm 1311b547a13bSvictorle .seealso: MatPartitioningCreate(), MatPartitioning 131276bdecfbSBarry Smith J*/ 131319fd82e9SBarry Smith typedef const char* MatPartitioningType; 13142692d6eeSBarry Smith #define MATPARTITIONINGCURRENT "current" 1315aa1e27eaSFande Kong #define MATPARTITIONINGAVERAGE "average" 13162692d6eeSBarry Smith #define MATPARTITIONINGSQUARE "square" 13172692d6eeSBarry Smith #define MATPARTITIONINGPARMETIS "parmetis" 13182692d6eeSBarry Smith #define MATPARTITIONINGCHACO "chaco" 13192692d6eeSBarry Smith #define MATPARTITIONINGPARTY "party" 132050d91057SBarry Smith #define MATPARTITIONINGPTSCOTCH "ptscotch" 132188d2ac2bSFande Kong #define MATPARTITIONINGHIERARCH "hierarch" 132271306c60Sjroman 1323ca161407SBarry Smith 1324014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningCreate(MPI_Comm,MatPartitioning*); 132519fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningSetType(MatPartitioning,MatPartitioningType); 1326014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetNParts(MatPartitioning,PetscInt); 1327014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetAdjacency(MatPartitioning,Mat); 1328014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]); 1329014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []); 1330014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningApply(MatPartitioning,IS*); 1331014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningDestroy(MatPartitioning*); 13322aabb6bbSBarry Smith 1333bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningRegister(const char[],PetscErrorCode (*)(MatPartitioning)); 133430de9b25SBarry Smith 1335f1144a30SSatish Balay 13362bad1931SBarry Smith 1337014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningView(MatPartitioning,PetscViewer); 1338014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetFromOptions(MatPartitioning); 133919fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningGetType(MatPartitioning,MatPartitioningType*); 1340ca161407SBarry Smith 134127538973SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetRepartition(MatPartitioning part); 1342014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetCoarseSequential(MatPartitioning); 1343014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningParmetisGetEdgeCut(MatPartitioning, PetscInt *); 13440752156aSBarry Smith 1345b6956c12SJose Roman typedef enum { MP_CHACO_MULTILEVEL=1,MP_CHACO_SPECTRAL=2,MP_CHACO_LINEAR=4,MP_CHACO_RANDOM=5,MP_CHACO_SCATTERED=6 } MPChacoGlobalType; 13466a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoGlobalTypes[]; 1347b6956c12SJose Roman typedef enum { MP_CHACO_KERNIGHAN=1,MP_CHACO_NONE=2 } MPChacoLocalType; 13486a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoLocalTypes[]; 1349b6956c12SJose Roman typedef enum { MP_CHACO_LANCZOS=0,MP_CHACO_RQI=1 } MPChacoEigenType; 13506a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoEigenTypes[]; 1351b6956c12SJose Roman 1352014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetGlobal(MatPartitioning,MPChacoGlobalType); 1353014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetGlobal(MatPartitioning,MPChacoGlobalType*); 1354014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetLocal(MatPartitioning,MPChacoLocalType); 1355014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetLocal(MatPartitioning,MPChacoLocalType*); 1356014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal); 1357014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType); 1358014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenSolver(MatPartitioning,MPChacoEigenType*); 1359014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenTol(MatPartitioning,PetscReal); 1360014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenTol(MatPartitioning,PetscReal*); 1361014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenNumber(MatPartitioning,PetscInt); 1362014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenNumber(MatPartitioning,PetscInt*); 136371306c60Sjroman 136471306c60Sjroman #define MP_PARTY_OPT "opt" 136571306c60Sjroman #define MP_PARTY_LIN "lin" 136671306c60Sjroman #define MP_PARTY_SCA "sca" 136771306c60Sjroman #define MP_PARTY_RAN "ran" 136871306c60Sjroman #define MP_PARTY_GBF "gbf" 136971306c60Sjroman #define MP_PARTY_GCF "gcf" 137071306c60Sjroman #define MP_PARTY_BUB "bub" 137171306c60Sjroman #define MP_PARTY_DEF "def" 1372014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetGlobal(MatPartitioning,const char*); 137371306c60Sjroman #define MP_PARTY_HELPFUL_SETS "hs" 137471306c60Sjroman #define MP_PARTY_KERNIGHAN_LIN "kl" 137571306c60Sjroman #define MP_PARTY_NONE "no" 1376014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetLocal(MatPartitioning,const char*); 1377014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal); 1378014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetBipart(MatPartitioning,PetscBool); 1379014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscBool); 138071306c60Sjroman 1381b0ca845eSVaclav Hapla typedef enum { MP_PTSCOTCH_DEFAULT,MP_PTSCOTCH_QUALITY,MP_PTSCOTCH_SPEED,MP_PTSCOTCH_BALANCE,MP_PTSCOTCH_SAFETY,MP_PTSCOTCH_SCALABILITY } MPPTScotchStrategyType; 13826a6fc655SJed Brown PETSC_EXTERN const char *const MPPTScotchStrategyTypes[]; 1383e0f1cffaSJose Roman 1384014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetImbalance(MatPartitioning,PetscReal); 1385014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetImbalance(MatPartitioning,PetscReal*); 1386014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetStrategy(MatPartitioning,MPPTScotchStrategyType); 1387014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetStrategy(MatPartitioning,MPPTScotchStrategyType*); 138871306c60Sjroman 1389b43b03e9SMark F. Adams /* 13906294fa2bSFande Kong * hierarchical partitioning 13916294fa2bSFande Kong */ 13924e713f55SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalGetFineparts(MatPartitioning,IS*); 13934e713f55SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalGetCoarseparts(MatPartitioning,IS*); 13944e713f55SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalSetNcoarseparts(MatPartitioning,PetscInt); 13954e713f55SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalSetNfineparts(MatPartitioning, PetscInt); 13966294fa2bSFande Kong 1397014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMeshToVertexGraph(Mat,PetscInt,Mat*); 1398014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMeshToCellGraph(Mat,PetscInt,Mat*); 1399591294e4SBarry Smith 14000752156aSBarry Smith /* 1401af0996ceSBarry Smith If you add entries here you must also add them to petsc/finclude/petscmat.h 1402d4fbbf0eSBarry Smith */ 14031c1c02c0SLois Curfman McInnes typedef enum { MATOP_SET_VALUES=0, 14041c1c02c0SLois Curfman McInnes MATOP_GET_ROW=1, 14051c1c02c0SLois Curfman McInnes MATOP_RESTORE_ROW=2, 14061c1c02c0SLois Curfman McInnes MATOP_MULT=3, 14071c1c02c0SLois Curfman McInnes MATOP_MULT_ADD=4, 14087c922b88SBarry Smith MATOP_MULT_TRANSPOSE=5, 14097c922b88SBarry Smith MATOP_MULT_TRANSPOSE_ADD=6, 14101c1c02c0SLois Curfman McInnes MATOP_SOLVE=7, 14111c1c02c0SLois Curfman McInnes MATOP_SOLVE_ADD=8, 14127c922b88SBarry Smith MATOP_SOLVE_TRANSPOSE=9, 14137c922b88SBarry Smith MATOP_SOLVE_TRANSPOSE_ADD=10, 14141c1c02c0SLois Curfman McInnes MATOP_LUFACTOR=11, 14151c1c02c0SLois Curfman McInnes MATOP_CHOLESKYFACTOR=12, 1416a714c06dSBarry Smith MATOP_SOR=13, 14171c1c02c0SLois Curfman McInnes MATOP_TRANSPOSE=14, 14181c1c02c0SLois Curfman McInnes MATOP_GETINFO=15, 14191c1c02c0SLois Curfman McInnes MATOP_EQUAL=16, 14201c1c02c0SLois Curfman McInnes MATOP_GET_DIAGONAL=17, 14211c1c02c0SLois Curfman McInnes MATOP_DIAGONAL_SCALE=18, 14221c1c02c0SLois Curfman McInnes MATOP_NORM=19, 14231c1c02c0SLois Curfman McInnes MATOP_ASSEMBLY_BEGIN=20, 14241c1c02c0SLois Curfman McInnes MATOP_ASSEMBLY_END=21, 1425d519adbfSMatthew Knepley MATOP_SET_OPTION=22, 1426d519adbfSMatthew Knepley MATOP_ZERO_ENTRIES=23, 1427d519adbfSMatthew Knepley MATOP_ZERO_ROWS=24, 1428d519adbfSMatthew Knepley MATOP_LUFACTOR_SYMBOLIC=25, 1429d519adbfSMatthew Knepley MATOP_LUFACTOR_NUMERIC=26, 1430d519adbfSMatthew Knepley MATOP_CHOLESKY_FACTOR_SYMBOLIC=27, 1431d519adbfSMatthew Knepley MATOP_CHOLESKY_FACTOR_NUMERIC=28, 1432d519adbfSMatthew Knepley MATOP_SETUP_PREALLOCATION=29, 1433d519adbfSMatthew Knepley MATOP_ILUFACTOR_SYMBOLIC=30, 1434d519adbfSMatthew Knepley MATOP_ICCFACTOR_SYMBOLIC=31, 1435a5b7ff6bSBarry Smith MATOP_GET_DIAGONAL_BLOCK=32, 14369d39f69dSJed Brown /* MATOP_PLACEHOLDER_33=33, */ 1437d519adbfSMatthew Knepley MATOP_DUPLICATE=34, 1438d519adbfSMatthew Knepley MATOP_FORWARD_SOLVE=35, 1439d519adbfSMatthew Knepley MATOP_BACKWARD_SOLVE=36, 1440d519adbfSMatthew Knepley MATOP_ILUFACTOR=37, 1441d519adbfSMatthew Knepley MATOP_ICCFACTOR=38, 1442d519adbfSMatthew Knepley MATOP_AXPY=39, 14437dae84e0SHong Zhang MATOP_CREATE_SUBMATRICES=40, 1444d519adbfSMatthew Knepley MATOP_INCREASE_OVERLAP=41, 1445d519adbfSMatthew Knepley MATOP_GET_VALUES=42, 1446d519adbfSMatthew Knepley MATOP_COPY=43, 1447d519adbfSMatthew Knepley MATOP_GET_ROW_MAX=44, 1448d519adbfSMatthew Knepley MATOP_SCALE=45, 1449d519adbfSMatthew Knepley MATOP_SHIFT=46, 145035153367SBarry Smith MATOP_DIAGONAL_SET=47, 14519d39f69dSJed Brown MATOP_ZERO_ROWS_COLUMNS=48, 14529d39f69dSJed Brown MATOP_SET_RANDOM=49, 1453d519adbfSMatthew Knepley MATOP_GET_ROW_IJ=50, 1454d519adbfSMatthew Knepley MATOP_RESTORE_ROW_IJ=51, 1455d519adbfSMatthew Knepley MATOP_GET_COLUMN_IJ=52, 1456d519adbfSMatthew Knepley MATOP_RESTORE_COLUMN_IJ=53, 1457d519adbfSMatthew Knepley MATOP_FDCOLORING_CREATE=54, 1458d519adbfSMatthew Knepley MATOP_COLORING_PATCH=55, 1459d519adbfSMatthew Knepley MATOP_SET_UNFACTORED=56, 1460d519adbfSMatthew Knepley MATOP_PERMUTE=57, 1461d519adbfSMatthew Knepley MATOP_SET_VALUES_BLOCKED=58, 14627dae84e0SHong Zhang MATOP_CREATE_SUBMATRIX=59, 1463d519adbfSMatthew Knepley MATOP_DESTROY=60, 1464d519adbfSMatthew Knepley MATOP_VIEW=61, 1465d519adbfSMatthew Knepley MATOP_CONVERT_FROM=62, 14667bab7c10SHong Zhang MATOP_MATMAT_MULT=63, 14677bab7c10SHong Zhang MATOP_MATMAT_MULT_SYMBOLIC=64, 14687bab7c10SHong Zhang MATOP_MATMAT_MULT_NUMERIC=65, 1469d519adbfSMatthew Knepley MATOP_SET_LOCAL_TO_GLOBAL_MAP=66, 1470d519adbfSMatthew Knepley MATOP_SET_VALUES_LOCAL=67, 1471d519adbfSMatthew Knepley MATOP_ZERO_ROWS_LOCAL=68, 1472d519adbfSMatthew Knepley MATOP_GET_ROW_MAX_ABS=69, 1473d519adbfSMatthew Knepley MATOP_GET_ROW_MIN_ABS=70, 1474d519adbfSMatthew Knepley MATOP_CONVERT=71, 1475d519adbfSMatthew Knepley MATOP_SET_COLORING=72, 14769d39f69dSJed Brown /* MATOP_PLACEHOLDER_73=73, */ 1477d519adbfSMatthew Knepley MATOP_SET_VALUES_ADIFOR=74, 1478d519adbfSMatthew Knepley MATOP_FD_COLORING_APPLY=75, 1479d519adbfSMatthew Knepley MATOP_SET_FROM_OPTIONS=76, 1480bda6c4cbSBarry Smith MATOP_MULT_CONSTRAINED=77, 1481bda6c4cbSBarry Smith MATOP_MULT_TRANSPOSE_CONSTRAIN=78, 14829d39f69dSJed Brown MATOP_FIND_ZERO_DIAGONALS=79, 1483d519adbfSMatthew Knepley MATOP_MULT_MULTIPLE=80, 1484d519adbfSMatthew Knepley MATOP_SOLVE_MULTIPLE=81, 1485d519adbfSMatthew Knepley MATOP_GET_INERTIA=82, 1486d519adbfSMatthew Knepley MATOP_LOAD=83, 1487d519adbfSMatthew Knepley MATOP_IS_SYMMETRIC=84, 1488d519adbfSMatthew Knepley MATOP_IS_HERMITIAN=85, 1489d519adbfSMatthew Knepley MATOP_IS_STRUCTURALLY_SYMMETRIC=86, 1490820469bcSHong Zhang MATOP_SET_VALUES_BLOCKEDLOCAL=87, 1491b41ce5d5SBarry Smith MATOP_CREATE_VECS=88, 1492d519adbfSMatthew Knepley MATOP_MAT_MULT=89, 1493d519adbfSMatthew Knepley MATOP_MAT_MULT_SYMBOLIC=90, 1494d519adbfSMatthew Knepley MATOP_MAT_MULT_NUMERIC=91, 1495d519adbfSMatthew Knepley MATOP_PTAP=92, 1496d519adbfSMatthew Knepley MATOP_PTAP_SYMBOLIC=93, 1497d519adbfSMatthew Knepley MATOP_PTAP_NUMERIC=94, 1498bda6c4cbSBarry Smith MATOP_MAT_TRANSPOSE_MULT=95, 1499bda6c4cbSBarry Smith MATOP_MAT_TRANSPOSE_MULT_SYMBO=96, 1500bda6c4cbSBarry Smith MATOP_MAT_TRANSPOSE_MULT_NUMER=97, 15019d39f69dSJed Brown /* MATOP_PLACEHOLDER_98=98, */ 15029d39f69dSJed Brown /* MATOP_PLACEHOLDER_99=99, */ 15039d39f69dSJed Brown /* MATOP_PLACEHOLDER_100=100, */ 15049d39f69dSJed Brown /* MATOP_PLACEHOLDER_101=101, */ 1505d519adbfSMatthew Knepley MATOP_CONJUGATE=102, 15069d39f69dSJed Brown /* MATOP_PLACEHOLDER_103=103, */ 1507d519adbfSMatthew Knepley MATOP_SET_VALUES_ROW=104, 1508d519adbfSMatthew Knepley MATOP_REAL_PART=105, 1509bda6c4cbSBarry Smith MATOP_IMAGINARY_PART=106, 1510bda6c4cbSBarry Smith MATOP_GET_ROW_UPPER_TRIANGULAR=107, 1511bda6c4cbSBarry Smith MATOP_RESTORE_ROW_UPPER_TRIANG=108, 1512bda6c4cbSBarry Smith MATOP_MAT_SOLVE=109, 151329eadf9eSStefano Zampini MATOP_MAT_SOLVE_TRANSPOSE=110, 1514d519adbfSMatthew Knepley MATOP_GET_ROW_MIN=111, 15150e8d04f7SBarry Smith MATOP_GET_COLUMN_VECTOR=112, 1516d519adbfSMatthew Knepley MATOP_MISSING_DIAGONAL=113, 15170e8d04f7SBarry Smith MATOP_GET_SEQ_NONZERO_STRUCTUR=114, 151889c6957cSBarry Smith MATOP_CREATE=115, 151989c6957cSBarry Smith MATOP_GET_GHOSTS=116, 15200e8d04f7SBarry Smith MATOP_GET_LOCAL_SUB_MATRIX=117, 15210e8d04f7SBarry Smith MATOP_RESTORE_LOCALSUB_MATRIX=118, 1522eeffb40dSHong Zhang MATOP_MULT_DIAGONAL_BLOCK=119, 15230e8d04f7SBarry Smith MATOP_HERMITIAN_TRANSPOSE=120, 15240e8d04f7SBarry Smith MATOP_MULT_HERMITIAN_TRANSPOSE=121, 15250e8d04f7SBarry Smith MATOP_MULT_HERMITIAN_TRANS_ADD=122, 15260e8d04f7SBarry Smith MATOP_GET_MULTI_PROC_BLOCK=123, 15279d39f69dSJed Brown MATOP_FIND_NONZERO_ROWS=124, 15280e8d04f7SBarry Smith MATOP_GET_COLUMN_NORMS=125, 15299d39f69dSJed Brown MATOP_INVERT_BLOCK_DIAGONAL=126, 15309d39f69dSJed Brown /* MATOP_PLACEHOLDER_127=127, */ 1531b41ce5d5SBarry Smith MATOP_CREATE_SUB_MATRICES_MPI=128, 15322b8ad9a3SHong Zhang MATOP_SET_VALUES_BATCH=129, 15330e8d04f7SBarry Smith MATOP_TRANSPOSE_MAT_MULT=130, 15340e8d04f7SBarry Smith MATOP_TRANSPOSE_MAT_MULT_SYMBO=131, 1535e9b602ebSSatish Balay MATOP_TRANSPOSE_MAT_MULT_NUMER=132, 15360e8d04f7SBarry Smith MATOP_TRANSPOSE_COLORING_CREAT=133, 15370e8d04f7SBarry Smith MATOP_TRANS_COLORING_APPLY_SPT=134, 15380e8d04f7SBarry Smith MATOP_TRANS_COLORING_APPLY_DEN=135, 15390e8d04f7SBarry Smith MATOP_RART=136, 15400e8d04f7SBarry Smith MATOP_RART_SYMBOLIC=137, 15410e8d04f7SBarry Smith MATOP_RART_NUMERIC=138, 1542e09a3074SHong Zhang MATOP_SET_BLOCK_SIZES=139, 1543f9426fe0SMark Adams MATOP_AYPX=140, 15441919a2e2SJed Brown MATOP_RESIDUAL=141, 15459c8f2541SHong Zhang MATOP_FDCOLORING_SETUP=142, 15462d033e1fSHong Zhang MATOP_MPICONCATENATESEQ=144, 15472d033e1fSHong Zhang MATOP_DESTROYSUBMATRICES=145 1548fae171e0SBarry Smith } MatOperation; 15490c0fd78eSBarry Smith PETSC_EXTERN PetscErrorCode MatSetOperation(Mat,MatOperation,void(*)(void)); 1550976e8c5aSLisandro Dalcin PETSC_EXTERN PetscErrorCode MatGetOperation(Mat,MatOperation,void(**)(void)); 1551976e8c5aSLisandro Dalcin PETSC_EXTERN PetscErrorCode MatHasOperation(Mat,MatOperation,PetscBool *); 1552976e8c5aSLisandro Dalcin 1553014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellSetOperation(Mat,MatOperation,void(*)(void)); 1554014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellGetOperation(Mat,MatOperation,void(**)(void)); 1555014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellSetContext(Mat,void*); 1556f3b1f45cSBarry Smith PETSC_EXTERN PetscErrorCode MatShellTestMult(Mat,PetscErrorCode (*)(void*,Vec,Vec),Vec,void*,PetscBool*); 1557f3b1f45cSBarry Smith PETSC_EXTERN PetscErrorCode MatShellTestMultTranspose(Mat,PetscErrorCode (*)(void*,Vec,Vec),Vec,void*,PetscBool*); 15580c0fd78eSBarry Smith PETSC_EXTERN PetscErrorCode MatShellSetManageScalingShifts(Mat); 1559112a2221SBarry Smith 156090ace30eSBarry Smith /* 156190ace30eSBarry Smith Codes for matrices stored on disk. By default they are 156290ace30eSBarry Smith stored in a universal format. By changing the format with 15636a9046bcSBarry Smith PetscViewerPushFormat(viewer,PETSC_VIEWER_NATIVE); the matrices will 156490ace30eSBarry Smith be stored in a way natural for the matrix, for example dense matrices 156590ace30eSBarry Smith would be stored as dense. Matrices stored this way may only be 1566f4403165SShri Abhyankar read into matrices of the same type. 156790ace30eSBarry Smith */ 156890ace30eSBarry Smith #define MATRIX_BINARY_FORMAT_DENSE -1 156990ace30eSBarry Smith 1570014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetHashTableFactor(Mat,PetscReal); 1571*75d48cdbSStefano Zampini 1572*75d48cdbSStefano Zampini PETSC_EXTERN PetscErrorCode MatISSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 1573*75d48cdbSStefano Zampini PETSC_EXTERN PetscErrorCode MatISSetUpSF(Mat); 1574*75d48cdbSStefano Zampini PETSC_EXTERN PetscErrorCode MatISStoreL2L(Mat,PetscBool); 1575014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatISGetLocalMat(Mat,Mat*); 15763b3b1effSJed Brown PETSC_EXTERN PetscErrorCode MatISRestoreLocalMat(Mat,Mat*); 1577014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatISSetLocalMat(Mat,Mat); 1578b7ce53b6SStefano Zampini PETSC_EXTERN PetscErrorCode MatISGetMPIXAIJ(Mat,MatReuse,Mat*); 15791f4e1ec7SBarry Smith 1580d9274352SBarry Smith /*S 1581d9274352SBarry Smith MatNullSpace - Object that removes a null space from a vector, i.e. 1582d9274352SBarry Smith orthogonalizes the vector to a subsapce 1583d9274352SBarry Smith 1584f7a9e4ceSBarry Smith Level: advanced 1585d9274352SBarry Smith 1586d9274352SBarry Smith Concepts: matrix; linear operator, null space 1587d9274352SBarry Smith 15886e1639daSBarry Smith Users manual sections: 15896e1639daSBarry Smith . sec_singular 15906e1639daSBarry Smith 1591d9274352SBarry Smith .seealso: MatNullSpaceCreate() 1592d9274352SBarry Smith S*/ 159374637425SBarry Smith typedef struct _p_MatNullSpace* MatNullSpace; 1594d9274352SBarry Smith 1595014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceCreate(MPI_Comm,PetscBool ,PetscInt,const Vec[],MatNullSpace*); 1596014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(MatNullSpace,Vec,void*),void*); 1597014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceDestroy(MatNullSpace*); 1598d0195637SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceRemove(MatNullSpace,Vec); 1599014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNullSpace(Mat, MatNullSpace *); 16005fa7ec2dSBarry Smith PETSC_EXTERN PetscErrorCode MatGetTransposeNullSpace(Mat, MatNullSpace *); 16015fa7ec2dSBarry Smith PETSC_EXTERN PetscErrorCode MatSetTransposeNullSpace(Mat,MatNullSpace); 1602014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNullSpace(Mat,MatNullSpace); 1603014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNearNullSpace(Mat,MatNullSpace); 1604014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNearNullSpace(Mat,MatNullSpace*); 1605014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceTest(MatNullSpace,Mat,PetscBool *); 1606014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceView(MatNullSpace,PetscViewer); 1607014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceGetVecs(MatNullSpace,PetscBool*,PetscInt*,const Vec**); 1608014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceCreateRigidBody(Vec,MatNullSpace*); 160974637425SBarry Smith 1610014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatReorderingSeqSBAIJ(Mat,IS); 1611014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetHashTableFactor(Mat,PetscReal); 1612014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetColumnIndices(Mat,PetscInt *); 1613014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJInvertBlockDiagonal(Mat); 16143f1d51d7SBarry Smith 1615014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMAIJ(Mat,PetscInt,Mat*); 1616014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMAIJRedimension(Mat,PetscInt,Mat*); 1617014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMAIJGetAIJ(Mat,Mat*); 1618c4f061fbSSatish Balay 1619014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatComputeExplicitOperator(Mat,Mat*); 1620f3b1f45cSBarry Smith PETSC_EXTERN PetscErrorCode MatComputeExplicitOperatorTranspose(Mat,Mat*); 1621b0a32e0cSBarry Smith 1622014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalScaleLocal(Mat,Vec); 162304f1ad80SBarry Smith 1624014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*); 1625014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetBase(Mat,Vec,Vec); 1626014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*); 1627014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*)); 1628014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec)); 1629014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt); 1630014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDResetHHistory(Mat); 1631014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctionError(Mat,PetscReal); 1632014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetPeriod(Mat,PetscInt); 1633014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDGetH(Mat,PetscScalar *); 1634014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetOptionsPrefix(Mat,const char[]); 1635014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*); 1636014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*); 1637e884886fSBarry Smith 16386370053bSBarry Smith /*S 16396370053bSBarry Smith MatMFFD - A data structured used to manage the computation of the h differencing parameter for matrix-free 16406370053bSBarry Smith Jacobian vector products 1641e884886fSBarry Smith 164295452b02SPatrick Sanan Notes: 164395452b02SPatrick Sanan MATMFFD is a specific MatType which uses the MatMFFD data structure 16446370053bSBarry Smith 16456370053bSBarry Smith MatMFFD*() methods actually take the Mat as their first argument. Not a MatMFFD data structure 16466370053bSBarry Smith 16476370053bSBarry Smith Level: developer 16486370053bSBarry Smith 16496370053bSBarry Smith .seealso: MATMFFD, MatCreateMFFD(), MatMFFDSetFuction(), MatMFFDSetType(), MatMFFDRegister() 16506370053bSBarry Smith S*/ 1651e884886fSBarry Smith typedef struct _p_MatMFFD* MatMFFD; 1652e884886fSBarry Smith 165376bdecfbSBarry Smith /*J 1654e884886fSBarry Smith MatMFFDType - algorithm used to compute the h used in computing matrix-vector products via differencing of the function 1655e884886fSBarry Smith 1656e884886fSBarry Smith Level: beginner 1657e884886fSBarry Smith 1658e884886fSBarry Smith .seealso: MatMFFDSetType(), MatMFFDRegister() 165976bdecfbSBarry Smith J*/ 166019fd82e9SBarry Smith typedef const char* MatMFFDType; 1661e884886fSBarry Smith #define MATMFFD_DS "ds" 1662e884886fSBarry Smith #define MATMFFD_WP "wp" 1663e884886fSBarry Smith 166419fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDSetType(Mat,MatMFFDType); 1665bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDRegister(const char[],PetscErrorCode (*)(MatMFFD)); 1666e884886fSBarry Smith 1667014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDDSSetUmin(Mat,PetscReal); 1668014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDWPSetComputeNormU(Mat,PetscBool ); 1669e884886fSBarry Smith 167061ab5df0SBarry Smith PETSC_EXTERN PetscErrorCode MatFDColoringSetType(MatFDColoring,MatMFFDType); 167161ab5df0SBarry Smith 1672014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *); 1673014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *); 16747dbadf16SMatthew Knepley 167597969023SHong Zhang /* 167697969023SHong Zhang PETSc interface to MUMPS 167797969023SHong Zhang */ 167897969023SHong Zhang #ifdef PETSC_HAVE_MUMPS 1679014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMumpsSetIcntl(Mat,PetscInt,PetscInt); 1680bc6112feSHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetIcntl(Mat,PetscInt,PetscInt*); 1681b9e9ea26SSatish Balay PETSC_EXTERN PetscErrorCode MatMumpsSetCntl(Mat,PetscInt,PetscReal); 1682bc6112feSHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetCntl(Mat,PetscInt,PetscReal*); 1683bc6112feSHong Zhang 1684ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetInfo(Mat,PetscInt,PetscInt*); 1685ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetInfog(Mat,PetscInt,PetscInt*); 1686ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetRinfo(Mat,PetscInt,PetscReal*); 1687ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetRinfog(Mat,PetscInt,PetscReal*); 168889a9c03aSHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetInverse(Mat,Mat); 168997969023SHong Zhang #endif 169023a5497aSJed Brown 1691d954db57SHong Zhang /* 1692d615f992SSatish Balay PETSc interface to Mkl_Pardiso 1693b500a6b7SJose David Bermeol */ 1694b500a6b7SJose David Bermeol #ifdef PETSC_HAVE_MKL_PARDISO 1695d615f992SSatish Balay PETSC_EXTERN PetscErrorCode MatMkl_PardisoSetCntl(Mat,PetscInt,PetscInt); 1696b500a6b7SJose David Bermeol #endif 1697b500a6b7SJose David Bermeol 1698b500a6b7SJose David Bermeol /* 1699d305a81bSVasiliy Kozyrev PETSc interface to Mkl_CPardiso 1700d305a81bSVasiliy Kozyrev */ 1701d305a81bSVasiliy Kozyrev #ifdef PETSC_HAVE_MKL_CPARDISO 1702d305a81bSVasiliy Kozyrev PETSC_EXTERN PetscErrorCode MatMkl_CPardisoSetCntl(Mat,PetscInt,PetscInt); 1703d305a81bSVasiliy Kozyrev #endif 1704d305a81bSVasiliy Kozyrev 1705d305a81bSVasiliy Kozyrev /* 1706d954db57SHong Zhang PETSc interface to SUPERLU 1707d954db57SHong Zhang */ 1708d954db57SHong Zhang #ifdef PETSC_HAVE_SUPERLU 1709014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSuperluSetILUDropTol(Mat,PetscReal); 1710d954db57SHong Zhang #endif 1711d954db57SHong Zhang 1712fb785984SHong Zhang /* 1713fb785984SHong Zhang PETSc interface to SUPERLU_DIST 1714fb785984SHong Zhang */ 1715fb785984SHong Zhang #ifdef PETSC_HAVE_SUPERLU_DIST 1716fb785984SHong Zhang PETSC_EXTERN PetscErrorCode MatSuperluDistGetDiagU(Mat,PetscScalar*); 1717fb785984SHong Zhang #endif 1718fb785984SHong Zhang 1719575704cbSPieter Ghysels /* 1720575704cbSPieter Ghysels PETSc interface to STRUMPACK 1721575704cbSPieter Ghysels */ 1722575704cbSPieter Ghysels #ifdef PETSC_HAVE_STRUMPACK 1723542aee0fSPieter Ghysels /*E 1724542aee0fSPieter Ghysels MatSTRUMPACKReordering - sparsity reducing ordering to be used in STRUMPACK 1725542aee0fSPieter Ghysels 1726542aee0fSPieter Ghysels Level: intermediate 1727542aee0fSPieter Ghysels E*/ 172834f43fa5SPieter Ghysels typedef enum {MAT_STRUMPACK_NATURAL=0, 172934f43fa5SPieter Ghysels MAT_STRUMPACK_METIS=1, 173034f43fa5SPieter Ghysels MAT_STRUMPACK_PARMETIS=2, 173134f43fa5SPieter Ghysels MAT_STRUMPACK_SCOTCH=3, 173234f43fa5SPieter Ghysels MAT_STRUMPACK_PTSCOTCH=4, 173334f43fa5SPieter Ghysels MAT_STRUMPACK_RCM=5} MatSTRUMPACKReordering; 173434f43fa5SPieter Ghysels PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetReordering(Mat,MatSTRUMPACKReordering); 1735575704cbSPieter Ghysels PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetColPerm(Mat,PetscBool); 173634f43fa5SPieter Ghysels PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSRelTol(Mat,PetscReal); 173734f43fa5SPieter Ghysels PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSAbsTol(Mat,PetscReal); 173834f43fa5SPieter Ghysels PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSMinSepSize(Mat,PetscInt); 173934f43fa5SPieter Ghysels PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSMaxRank(Mat,PetscInt); 1740a36bf211SPieter Ghysels PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSLeafSize(Mat,PetscInt); 1741575704cbSPieter Ghysels #endif 1742575704cbSPieter Ghysels 1743575704cbSPieter Ghysels 1744aa372e3fSPaul Mullowney #ifdef PETSC_HAVE_CUDA 17453f9c0db1SPaul Mullowney /*E 1746e057df02SPaul Mullowney MatCUSPARSEStorageFormat - indicates the storage format for CUSPARSE (GPU) 17472692e278SPaul Mullowney matrices. 1748e057df02SPaul Mullowney 1749e057df02SPaul Mullowney Not Collective 1750e057df02SPaul Mullowney 17518468deeeSKarl Rupp + MAT_CUSPARSE_CSR - Compressed Sparse Row 17522692e278SPaul Mullowney . MAT_CUSPARSE_ELL - Ellpack (requires CUDA 4.2 or later). 17532692e278SPaul Mullowney - MAT_CUSPARSE_HYB - Hybrid, a combination of Ellpack and Coordinate format (requires CUDA 4.2 or later). 1754e057df02SPaul Mullowney 1755e057df02SPaul Mullowney Level: intermediate 1756e057df02SPaul Mullowney 1757af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 1758e057df02SPaul Mullowney 1759e057df02SPaul Mullowney .seealso: MatCUSPARSESetFormat(), MatCUSPARSEFormatOperation 1760e057df02SPaul Mullowney E*/ 1761e057df02SPaul Mullowney 17623f9c0db1SPaul Mullowney typedef enum {MAT_CUSPARSE_CSR, MAT_CUSPARSE_ELL, MAT_CUSPARSE_HYB} MatCUSPARSEStorageFormat; 1763e057df02SPaul Mullowney 1764e057df02SPaul Mullowney /* these will be strings associated with enumerated type defined above */ 1765e057df02SPaul Mullowney PETSC_EXTERN const char *const MatCUSPARSEStorageFormats[]; 1766e057df02SPaul Mullowney 17673f9c0db1SPaul Mullowney /*E 1768e057df02SPaul Mullowney MatCUSPARSEFormatOperation - indicates the operation of CUSPARSE (GPU) 17692692e278SPaul Mullowney matrices whose operation should use a particular storage format. 1770e057df02SPaul Mullowney 1771e057df02SPaul Mullowney Not Collective 1772e057df02SPaul Mullowney 17738468deeeSKarl Rupp + MAT_CUSPARSE_MULT_DIAG - sets the storage format for the diagonal matrix in the parallel MatMult 17748468deeeSKarl Rupp . MAT_CUSPARSE_MULT_OFFDIAG - sets the storage format for the offdiagonal matrix in the parallel MatMult 17758468deeeSKarl Rupp . MAT_CUSPARSE_MULT - sets the storage format for the entire matrix in the serial (single GPU) MatMult 17768468deeeSKarl Rupp - MAT_CUSPARSE_ALL - sets the storage format for all CUSPARSE (GPU) matrices 1777e057df02SPaul Mullowney 1778e057df02SPaul Mullowney Level: intermediate 1779e057df02SPaul Mullowney 1780e057df02SPaul Mullowney .seealso: MatCUSPARSESetFormat(), MatCUSPARSEStorageFormat 1781e057df02SPaul Mullowney E*/ 178236d62e41SPaul Mullowney typedef enum {MAT_CUSPARSE_MULT_DIAG, MAT_CUSPARSE_MULT_OFFDIAG, MAT_CUSPARSE_MULT, MAT_CUSPARSE_ALL} MatCUSPARSEFormatOperation; 1783e057df02SPaul Mullowney 178410e9db80SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 178510e9db80SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCreateAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 1786e057df02SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCUSPARSESetFormat(Mat,MatCUSPARSEFormatOperation,MatCUSPARSEStorageFormat); 17879ae82921SPaul Mullowney #endif 17889ae82921SPaul Mullowney 1789d67ff14aSKarl Rupp #if defined(PETSC_HAVE_VIENNACL) 1790d67ff14aSKarl Rupp PETSC_EXTERN PetscErrorCode MatCreateSeqAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 1791d67ff14aSKarl Rupp PETSC_EXTERN PetscErrorCode MatCreateAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 1792d67ff14aSKarl Rupp #endif 1793d67ff14aSKarl Rupp 179454efbe56SHong Zhang /* 179554efbe56SHong Zhang PETSc interface to FFTW 179654efbe56SHong Zhang */ 179754efbe56SHong Zhang #if defined(PETSC_HAVE_FFTW) 1798014dd563SJed Brown PETSC_EXTERN PetscErrorCode VecScatterPetscToFFTW(Mat,Vec,Vec); 1799014dd563SJed Brown PETSC_EXTERN PetscErrorCode VecScatterFFTWToPetsc(Mat,Vec,Vec); 18002a7a6963SBarry Smith PETSC_EXTERN PetscErrorCode MatCreateVecsFFTW(Mat,Vec*,Vec*,Vec*); 180154efbe56SHong Zhang #endif 180254efbe56SHong Zhang 1803014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateNest(MPI_Comm,PetscInt,const IS[],PetscInt,const IS[],const Mat[],Mat*); 1804014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSize(Mat,PetscInt*,PetscInt*); 1805014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetISs(Mat,IS[],IS[]); 1806014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetLocalISs(Mat,IS[],IS[]); 1807014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSubMats(Mat,PetscInt*,PetscInt*,Mat***); 1808014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSubMat(Mat,PetscInt,PetscInt,Mat*); 180919fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatNestSetVecType(Mat,VecType); 1810014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestSetSubMats(Mat,PetscInt,const IS[],PetscInt,const IS[],const Mat[]); 1811014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestSetSubMat(Mat,PetscInt,PetscInt,Mat); 1812d8588912SDave May 18134325cce7SMatthew G Knepley PETSC_EXTERN PetscErrorCode MatChop(Mat,PetscReal); 1814e0a58c10SMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatComputeBandwidth(Mat,PetscReal,PetscInt*); 18154325cce7SMatthew G Knepley 1816b541e6a4SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSubdomainsCreateCoalesce(Mat,PetscInt,PetscInt*,IS**); 181753dd7562SDmitry Karpeev 1818c094ef40SMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatPreallocatorPreallocate(Mat,PetscBool,Mat); 1819c094ef40SMatthew G. Knepley 1820539c167fSBarry Smith PETSC_INTERN PetscErrorCode MatHeaderMerge(Mat,Mat*); 1821539c167fSBarry Smith PETSC_EXTERN PetscErrorCode MatHeaderReplace(Mat,Mat*); 1822539c167fSBarry Smith 182323a5497aSJed Brown #endif 1824