12eac72dbSBarry Smith /* 22eac72dbSBarry Smith Include file for the matrix component of PETSc 32eac72dbSBarry Smith */ 40a835dfdSSatish Balay #ifndef __PETSCMAT_H 50a835dfdSSatish Balay #define __PETSCMAT_H 62c8e378dSBarry Smith #include <petscvec.h> 72eac72dbSBarry Smith 8d9274352SBarry Smith /*S 98f6c3df8SBarry Smith Mat - Abstract PETSc matrix object used to manage all linear operators in PETSc, even those without 108f6c3df8SBarry Smith an explicit sparse representation (such as matrix-free operators) 112eac72dbSBarry Smith 12d91e6319SBarry Smith Level: beginner 13d91e6319SBarry Smith 14d9274352SBarry Smith Concepts: matrix; linear operator 15d9274352SBarry Smith 168f6c3df8SBarry Smith .seealso: MatCreate(), MatType, MatSetType(), MatDestroy() 17d9274352SBarry Smith S*/ 18d9274352SBarry Smith typedef struct _p_Mat* Mat; 19d9274352SBarry Smith 2076bdecfbSBarry Smith /*J 218f6c3df8SBarry Smith MatType - String with the name of a PETSc matrix type 22d9274352SBarry Smith 23d9274352SBarry Smith Level: beginner 24d9274352SBarry Smith 258f6c3df8SBarry Smith .seealso: MatSetType(), Mat, MatSolverPackage, MatRegister() 2676bdecfbSBarry Smith J*/ 2719fd82e9SBarry Smith typedef const char* MatType; 28273d9f13SBarry Smith #define MATSAME "same" 295a11e1b2SBarry Smith #define MATMAIJ "maij" 30273d9f13SBarry Smith #define MATSEQMAIJ "seqmaij" 31273d9f13SBarry Smith #define MATMPIMAIJ "mpimaij" 32273d9f13SBarry Smith #define MATIS "is" 335a11e1b2SBarry Smith #define MATAIJ "aij" 34273d9f13SBarry Smith #define MATSEQAIJ "seqaij" 35273d9f13SBarry Smith #define MATMPIAIJ "mpiaij" 365a11e1b2SBarry Smith #define MATAIJCRL "aijcrl" 375a11e1b2SBarry Smith #define MATSEQAIJCRL "seqaijcrl" 385a11e1b2SBarry Smith #define MATMPIAIJCRL "mpiaijcrl" 398154be41SBarry Smith #define MATAIJCUSP "aijcusp" 408154be41SBarry Smith #define MATSEQAIJCUSP "seqaijcusp" 418154be41SBarry Smith #define MATMPIAIJCUSP "mpiaijcusp" 429ae82921SPaul Mullowney #define MATAIJCUSPARSE "aijcusparse" 439ae82921SPaul Mullowney #define MATSEQAIJCUSPARSE "seqaijcusparse" 449ae82921SPaul Mullowney #define MATMPIAIJCUSPARSE "mpiaijcusparse" 45d67ff14aSKarl Rupp #define MATAIJVIENNACL "aijviennacl" 46d67ff14aSKarl Rupp #define MATSEQAIJVIENNACL "seqaijviennacl" 47d67ff14aSKarl Rupp #define MATMPIAIJVIENNACL "mpiaijviennacl" 485a11e1b2SBarry Smith #define MATAIJPERM "aijperm" 495a11e1b2SBarry Smith #define MATSEQAIJPERM "seqaijperm" 505a11e1b2SBarry Smith #define MATMPIAIJPERM "mpiaijperm" 514a2a386eSRichard Tran Mills #define MATAIJMKL "aijmkl" 524a2a386eSRichard Tran Mills #define MATSEQAIJMKL "seqaijmkl" 534a2a386eSRichard Tran Mills #define MATMPIAIJMKL "mpiaijmkl" 54b5b72c8aSIrina Sokolova #define MATBAIJMKL "baijmkl" 55b5b72c8aSIrina Sokolova #define MATSEQBAIJMKL "seqbaijmkl" 56b5b72c8aSIrina Sokolova #define MATMPIBAIJMKL "mpibaijmkl" 57273d9f13SBarry Smith #define MATSHELL "shell" 585a11e1b2SBarry Smith #define MATDENSE "dense" 59209238afSKris Buschelman #define MATSEQDENSE "seqdense" 60273d9f13SBarry Smith #define MATMPIDENSE "mpidense" 61db31f6deSJed Brown #define MATELEMENTAL "elemental" 625a11e1b2SBarry Smith #define MATBAIJ "baij" 63273d9f13SBarry Smith #define MATSEQBAIJ "seqbaij" 64273d9f13SBarry Smith #define MATMPIBAIJ "mpibaij" 65273d9f13SBarry Smith #define MATMPIADJ "mpiadj" 665a11e1b2SBarry Smith #define MATSBAIJ "sbaij" 67273d9f13SBarry Smith #define MATSEQSBAIJ "seqsbaij" 68273d9f13SBarry Smith #define MATMPISBAIJ "mpisbaij" 69cebc7f6cSBarry Smith #define MATDAAD "daad" 70cebc7f6cSBarry Smith #define MATMFFD "mffd" 71c8a8475eSBarry Smith #define MATNORMAL "normal" 72c5e4d11fSDmitry Karpeev #define MATNORMALHERMITIAN "normalh" 73ab92ecdeSBarry Smith #define MATLRC "lrc" 742a6744ebSBarry Smith #define MATSCATTER "scatter" 75421e10b8SBarry Smith #define MATBLOCKMAT "blockmat" 76793850ffSBarry Smith #define MATCOMPOSITE "composite" 771f8c7532SHong Zhang #define MATFFT "fft" 781f8c7532SHong Zhang #define MATFFTW "fftw" 79e133240eSMatthew G Knepley #define MATSEQCUFFT "seqcufft" 80557cca28SSatish Balay #define MATTRANSPOSEMAT "transpose" 8172ca8751SBarry Smith #define MATSCHURCOMPLEMENT "schurcomplement" 821d6018f0SLisandro Dalcin #define MATPYTHON "python" 8363c07aadSStefano Zampini #define MATHYPRE "hypre" 84f91d8e95SBarry Smith #define MATHYPRESTRUCT "hyprestruct" 85a9e6138eSGlenn Hammond #define MATHYPRESSTRUCT "hypresstruct" 86ee1cef2cSJed Brown #define MATSUBMATRIX "submatrix" 872c0dbf93SJed Brown #define MATLOCALREF "localref" 88d8588912SDave May #define MATNEST "nest" 89c094ef40SMatthew G. Knepley #define MATPREALLOCATOR "preallocator" 90*d4002b98SHong Zhang #define MATSELL "sell" 91*d4002b98SHong Zhang #define MATSEQSELL "seqsell" 92*d4002b98SHong Zhang #define MATMPISELL "mpisell" 93a3b2e22bSHong Zhang #define MATDUMMY "dummy" 94773941d6SBarry Smith 9576bdecfbSBarry Smith /*J 96c7393fdbSBarry Smith MatSolverPackage - String with the name of a PETSc matrix solver type. 979989ab13SBarry Smith 98c2b89b5dSBarry Smith For example: "petsc" indicates what PETSc provides, "superlu_dist" the parallel SuperLU_DIST package etc 999989ab13SBarry Smith 1009989ab13SBarry Smith Level: beginner 1019989ab13SBarry Smith 1025c9eb25fSBarry Smith .seealso: MatGetFactor(), Mat, MatSetType(), MatType 10376bdecfbSBarry Smith J*/ 104c7393fdbSBarry Smith #define MatSolverPackage char* 1052692d6eeSBarry Smith #define MATSOLVERSUPERLU "superlu" 1062692d6eeSBarry Smith #define MATSOLVERSUPERLU_DIST "superlu_dist" 10708f5efcfSPieter Ghysels #define MATSOLVERSTRUMPACK "strumpack" 1082692d6eeSBarry Smith #define MATSOLVERUMFPACK "umfpack" 10920db9a53SJed Brown #define MATSOLVERCHOLMOD "cholmod" 110d89f5e7aSBarry Smith #define MATSOLVERKLU "klu" 111418810c4SBarry Smith #define MATSOLVERSPARSEELEMENTAL "sparseelemental" 112d89f5e7aSBarry Smith #define MATSOLVERELEMENTAL "elemental" 1132692d6eeSBarry Smith #define MATSOLVERESSL "essl" 1142692d6eeSBarry Smith #define MATSOLVERLUSOL "lusol" 1152692d6eeSBarry Smith #define MATSOLVERMUMPS "mumps" 116d615f992SSatish Balay #define MATSOLVERMKL_PARDISO "mkl_pardiso" 117d305a81bSVasiliy Kozyrev #define MATSOLVERMKL_CPARDISO "mkl_cpardiso" 1182692d6eeSBarry Smith #define MATSOLVERPASTIX "pastix" 1192692d6eeSBarry Smith #define MATSOLVERMATLAB "matlab" 1202692d6eeSBarry Smith #define MATSOLVERPETSC "petsc" 1212692d6eeSBarry Smith #define MATSOLVERBAS "bas" 1229ae82921SPaul Mullowney #define MATSOLVERCUSPARSE "cusparse" 123c0cdd4a1SDahai Guo 124b24902e0SBarry Smith /*E 125b24902e0SBarry Smith MatFactorType - indicates what type of factorization is requested 126b24902e0SBarry Smith 127b24902e0SBarry Smith Level: beginner 128b24902e0SBarry Smith 129af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 130b24902e0SBarry Smith 131c7393fdbSBarry Smith .seealso: MatSolverPackage, MatGetFactor() 132b24902e0SBarry Smith E*/ 133599ef60dSHong Zhang typedef enum {MAT_FACTOR_NONE, MAT_FACTOR_LU, MAT_FACTOR_CHOLESKY, MAT_FACTOR_ILU, MAT_FACTOR_ICC,MAT_FACTOR_ILUDT} MatFactorType; 134014dd563SJed Brown PETSC_EXTERN const char *const MatFactorTypes[]; 135e92e720dSBarry Smith 136014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactor(Mat,const MatSolverPackage,MatFactorType,Mat*); 137014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactorAvailable(Mat,const MatSolverPackage,MatFactorType,PetscBool *); 138014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFactorGetSolverPackage(Mat,const MatSolverPackage*); 139014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactorType(Mat,MatFactorType*); 14018713533SBarry Smith PETSC_EXTERN PetscErrorCode MatSolverPackageRegister(const MatSolverPackage,const MatType,MatFactorType,PetscErrorCode(*)(Mat,MatFactorType,Mat*)); 14142c9c57cSBarry Smith PETSC_EXTERN PetscErrorCode MatSolverPackageGet(const MatSolverPackage,const MatType,MatFactorType,PetscBool*,PetscBool*,PetscErrorCode (**)(Mat,MatFactorType,Mat*)); 1429989ab13SBarry Smith 143c06d978dSMatthew Knepley /* Logging support */ 1440700a824SBarry Smith #define MAT_FILE_CLASSID 1211216 /* used to indicate matrices in binary files */ 145014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_CLASSID; 146335efc43SPeter Brune PETSC_EXTERN PetscClassId MAT_COLORING_CLASSID; 147014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_FDCOLORING_CLASSID; 148014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_TRANSPOSECOLORING_CLASSID; 149014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_PARTITIONING_CLASSID; 150014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_COARSEN_CLASSID; 151014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_NULLSPACE_CLASSID; 152014dd563SJed Brown PETSC_EXTERN PetscClassId MATMFFD_CLASSID; 153c06d978dSMatthew Knepley 154ceb03754SKris Buschelman /*E 1557dae84e0SHong Zhang MatReuse - Indicates if matrices obtained from a previous call to MatCreateSubMatrices(), MatCreateSubMatrix(), MatConvert() or several other functions 156cf37664fSBarry Smith are to be reused to store the new matrix values. 157cf37664fSBarry Smith 158cf37664fSBarry Smith $ MAT_INITIAL_MATRIX - create a new matrix 159cf37664fSBarry Smith $ MAT_REUSE_MATRIX - reuse the matrix created with a previous call that used MAT_INITIAL_MATRIX 160cf37664fSBarry Smith $ MAT_INPLACE_MATRIX - replace the first input matrix with the new matrix (not applicable to all functions) 161cf37664fSBarry 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) 162ceb03754SKris Buschelman 163ceb03754SKris Buschelman Level: beginner 164ceb03754SKris Buschelman 165af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 166ceb03754SKris Buschelman 1677dae84e0SHong Zhang .seealso: MatCreateSubMatrices(), MatCreateSubMatrix(), MatDestroyMatrices(), MatConvert() 168ceb03754SKris Buschelman E*/ 169511c6705SHong Zhang typedef enum {MAT_INITIAL_MATRIX,MAT_REUSE_MATRIX,MAT_IGNORE_MATRIX,MAT_INPLACE_MATRIX} MatReuse; 170ceb03754SKris Buschelman 1715494a064SHong Zhang /*E 1727dae84e0SHong Zhang MatCreateSubMatrixOption - Indicates if matrices obtained from a call to MatCreateSubMatrices() 173829201f2SHong Zhang include the matrix values. Currently it is only used by MatGetSeqNonzerostructure(). 1745494a064SHong Zhang 1755494a064SHong Zhang Level: beginner 1765494a064SHong Zhang 177829201f2SHong Zhang .seealso: MatGetSeqNonzerostructure() 1785494a064SHong Zhang E*/ 1797dae84e0SHong Zhang typedef enum {MAT_DO_NOT_GET_VALUES,MAT_GET_VALUES} MatCreateSubMatrixOption; 1805494a064SHong Zhang 181607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatInitializePackage(void); 182c06d978dSMatthew Knepley 183014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreate(MPI_Comm,Mat*); 184014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetSizes(Mat,PetscInt,PetscInt,PetscInt,PetscInt); 18519fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatSetType(Mat,MatType); 186014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetFromOptions(Mat); 187685405a1SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatViewFromOptions(Mat A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);} 188bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatRegister(const char[],PetscErrorCode(*)(Mat)); 189014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRegisterBaseName(const char[],const char[],const char[]); 190014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetOptionsPrefix(Mat,const char[]); 191014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAppendOptionsPrefix(Mat,const char[]); 192014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOptionsPrefix(Mat,const char*[]); 19384d44b13SHong Zhang PETSC_EXTERN PetscErrorCode MatSetErrorIfFailure(Mat,PetscBool); 194f69a0ea3SMatthew Knepley 195140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatList; 196140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatColoringList; 197140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatPartitioningList; 19828988994SBarry Smith 1993b224e63SBarry Smith /*E 200aa6c7ce3SBarry Smith MatStructure - Indicates if two matrices have the same nonzero structure 2013b224e63SBarry Smith 2023b224e63SBarry Smith Level: beginner 2033b224e63SBarry Smith 204af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 2053b224e63SBarry Smith 206aa6c7ce3SBarry Smith .seealso: MatCopy(), MatAXPY() 2073b224e63SBarry Smith E*/ 208aa6c7ce3SBarry Smith typedef enum {DIFFERENT_NONZERO_PATTERN,SUBSET_NONZERO_PATTERN,SAME_NONZERO_PATTERN} MatStructure; 2093b224e63SBarry Smith 2109779e05dSSatish Balay #if defined PETSC_HAVE_MKL_SPARSE 21180095d54SIrina Sokolova PETSC_EXTERN PetscErrorCode MatCreateBAIJMKL(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 21280095d54SIrina Sokolova PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJMKL(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[],Mat *A); 21380095d54SIrina Sokolova #endif 214*d4002b98SHong Zhang PETSC_EXTERN PetscErrorCode MatCreateSeqSELL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 215*d4002b98SHong Zhang PETSC_EXTERN PetscErrorCode MatCreateSELL(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 216*d4002b98SHong Zhang PETSC_EXTERN PetscErrorCode MatSeqSELLSetPreallocation(Mat,PetscInt,const PetscInt[]); 217*d4002b98SHong Zhang PETSC_EXTERN PetscErrorCode MatMPISELLSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 21880095d54SIrina Sokolova 219014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqDense(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*); 220014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateDense(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*); 221014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 222014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 223014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *); 224014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithSplitArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],PetscInt[],PetscInt[],PetscScalar[],Mat*); 2253b00a383SHong Zhang PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithSeqAIJ(MPI_Comm,Mat,Mat,const PetscInt[],Mat*); 2263b00a383SHong Zhang 227014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 228014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 229014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat*); 230d21a29f3SJed Brown 231014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAdj(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscInt[],Mat*); 232014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 233d21a29f3SJed Brown 234014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 235014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPISBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *); 23638f409ebSLisandro Dalcin PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 237014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 2384cce697cSJed Brown PETSC_EXTERN PetscErrorCode MatXAIJSetPreallocation(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscInt[],const PetscInt[]); 239dfb205c3SBarry Smith 240014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateShell(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,void *,Mat*); 241014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateNormal(Mat,Mat*); 242c5e4d11fSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatCreateNormalHermitian(Mat,Mat*); 243986c4d72SJose E. Roman PETSC_EXTERN PetscErrorCode MatCreateLRC(Mat,Mat,Vec,Mat,Mat*); 244267ca84cSJose E. Roman PETSC_EXTERN PetscErrorCode MatLRCGetMats(Mat,Mat*,Mat*,Vec*,Mat*); 245c5e4d11fSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatCreateIS(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,ISLocalToGlobalMapping,ISLocalToGlobalMapping,Mat*); 246014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 247014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 248c0cdd4a1SDahai Guo 249014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateScatter(MPI_Comm,VecScatter,Mat*); 250014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScatterSetVecScatter(Mat,VecScatter); 251014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScatterGetVecScatter(Mat,VecScatter*); 252014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateBlockMat(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,Mat*); 253014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeAddMat(Mat,Mat); 254014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeMerge(Mat); 255014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateComposite(MPI_Comm,PetscInt,const Mat*,Mat*); 2566d7c1e57SBarry Smith typedef enum {MAT_COMPOSITE_ADDITIVE,MAT_COMPOSITE_MULTIPLICATIVE} MatCompositeType; 257014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeSetType(Mat,MatCompositeType); 2586d7c1e57SBarry Smith 25919fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatCreateFFT(MPI_Comm,PetscInt,const PetscInt[],MatType,Mat*); 260014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqCUFFT(MPI_Comm,PetscInt,const PetscInt[],Mat*); 261dedccee8SHong Zhang 262014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateTranspose(Mat,Mat*); 2638060fb66Sstefano_zampini PETSC_EXTERN PetscErrorCode MatTransposeGetMat(Mat,Mat*); 264d0de2241SAndrew Spott PETSC_EXTERN PetscErrorCode MatCreateHermitianTranspose(Mat,Mat*); 26554e05e6cSHong Zhang PETSC_EXTERN PetscErrorCode MatCreateSubMatrixVirtual(Mat,IS,IS,Mat*); 26654e05e6cSHong Zhang PETSC_EXTERN PetscErrorCode MatSubMatrixVirtualUpdate(Mat,Mat,IS,IS); 267014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateLocalRef(Mat,IS,IS,Mat*); 2681472f72bSBarry Smith 269978814f1SStefano Zampini #if defined(PETSC_HAVE_HYPRE) 270d975228cSstefano_zampini PETSC_EXTERN PetscErrorCode MatHYPRESetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 271978814f1SStefano Zampini #endif 272978814f1SStefano Zampini 273014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPythonSetType(Mat,const char[]); 2741d6018f0SLisandro Dalcin 275014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetUp(Mat); 276014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroy(Mat*); 277e56f5c9eSBarry Smith PETSC_EXTERN PetscErrorCode MatGetNonzeroState(Mat,PetscObjectState*); 27821c89e3eSBarry Smith 279014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatConjugate(Mat); 280014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRealPart(Mat); 281014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatImaginaryPart(Mat); 282014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetDiagonalBlock(Mat,Mat*); 283014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetTrace(Mat,PetscScalar*); 284713ccfa9SJed Brown PETSC_EXTERN PetscErrorCode MatInvertBlockDiagonal(Mat,const PetscScalar **); 28599cafbc1SBarry Smith 2868ed539a5SBarry Smith /* ------------------------------------------------------------*/ 287014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 288014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 289014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesRow(Mat,PetscInt,const PetscScalar[]); 290014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesRowLocal(Mat,PetscInt,const PetscScalar[]); 291014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBatch(Mat,PetscInt,PetscInt,PetscInt[],const PetscScalar[]); 29273a71a0fSBarry Smith PETSC_EXTERN PetscErrorCode MatSetRandom(Mat,PetscRandom); 29384cb2905SBarry Smith 2942ef4de8bSBarry Smith /*S 2952ef4de8bSBarry Smith MatStencil - Data structure (C struct) for storing information about a single row or 29662ef0227SBarry Smith column of a matrix as indexed on an associated grid. These are arguments to MatSetStencil() and MatSetBlockStencil() 29762ef0227SBarry Smith 29862ef0227SBarry 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). 29962ef0227SBarry 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. 30062ef0227SBarry Smith 30162ef0227SBarry Smith For stencil access to vectors see DMDAVecGetArray(), DMDAVecGetArrayF90(). 302be479b30SJed Brown 303be479b30SJed Brown Fortran usage is different, see MatSetValuesStencil() for details. 3042ef4de8bSBarry Smith 3052ef4de8bSBarry Smith Level: beginner 3062ef4de8bSBarry Smith 3072ef4de8bSBarry Smith Concepts: matrix; linear operator 3082ef4de8bSBarry Smith 30962ef0227SBarry Smith .seealso: MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockedStencil(), DMDAVecGetArray(), DMDAVecGetArrayF90() 3102ef4de8bSBarry Smith S*/ 311435da068SBarry Smith typedef struct { 312c1ac3661SBarry Smith PetscInt k,j,i,c; 313435da068SBarry Smith } MatStencil; 3142ef4de8bSBarry Smith 315014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode); 316014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode); 317014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt); 318435da068SBarry Smith 319d91e6319SBarry Smith /*E 320d91e6319SBarry Smith MatAssemblyType - Indicates if the matrix is now to be used, or if you plan 32162ef0227SBarry Smith to continue to add or insert values to it 322d91e6319SBarry Smith 323d91e6319SBarry Smith Level: beginner 324d91e6319SBarry Smith 325d91e6319SBarry Smith .seealso: MatAssemblyBegin(), MatAssemblyEnd() 326d91e6319SBarry Smith E*/ 3276d4a8577SBarry Smith typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType; 328014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssemblyBegin(Mat,MatAssemblyType); 329014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssemblyEnd(Mat,MatAssemblyType); 330014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssembled(Mat,PetscBool *); 3314f9c727eSBarry Smith 3321d73ed98SBarry Smith 33330de9b25SBarry Smith 334d91e6319SBarry Smith /*E 335d91e6319SBarry Smith MatOption - Options that may be set for a matrix and its behavior or storage 336d91e6319SBarry Smith 337d91e6319SBarry Smith Level: beginner 338d91e6319SBarry Smith 339af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 3400f8fb01aSBarry Smith Any additions/changes here must also be made in src/mat/interface/dlregismat.c in MatOptions[] 341d91e6319SBarry Smith 342382164b0SBarry Smith Developer Notes: 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); 397014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBlockSize(Mat,PetscInt *); 398014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSize(Mat,PetscInt); 399014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBlockSizes(Mat,PetscInt *,PetscInt *); 400014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSizes(Mat,PetscInt,PetscInt); 40133d57670SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSizesFromMats(Mat,Mat,Mat); 4021620fd73SBarry Smith 403014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMult(Mat,Vec,Vec); 404014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultDiagonalBlock(Mat,Vec,Vec); 405014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultAdd(Mat,Vec,Vec,Vec); 406014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTranspose(Mat,Vec,Vec); 407014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultHermitianTranspose(Mat,Vec,Vec); 408014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsTranspose(Mat,Mat,PetscReal,PetscBool *); 409014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitianTranspose(Mat,Mat,PetscReal,PetscBool *); 410014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeAdd(Mat,Vec,Vec,Vec); 411014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultHermitianTransposeAdd(Mat,Vec,Vec,Vec); 412014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultConstrained(Mat,Vec,Vec); 413014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeConstrained(Mat,Vec,Vec); 414014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatSolve(Mat,Mat,Mat); 415bdc285e1SStefano Zampini PETSC_EXTERN PetscErrorCode MatMatSolveTranspose(Mat,Mat,Mat); 416f9426fe0SMark Adams PETSC_EXTERN PetscErrorCode MatResidual(Mat,Vec,Vec,Vec); 417f5edf698SHong Zhang 418d91e6319SBarry Smith /*E 419d91e6319SBarry Smith MatDuplicateOption - Indicates if a duplicated sparse matrix should have 420d91e6319SBarry Smith its numerical values copied over or just its nonzero structure. 421d91e6319SBarry Smith 422d91e6319SBarry Smith Level: beginner 423d91e6319SBarry Smith 424af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 425d91e6319SBarry Smith 4260018da39SRichard Tran Mills $ MAT_DO_NOT_COPY_VALUES - Create a matrix using the same nonzero pattern as the original matrix, 4270018da39SRichard Tran Mills $ with zeros for the numerical values. 4280018da39SRichard Tran Mills $ MAT_COPY_VALUES - Create a matrix with the same nonzero pattern as the original matrix 4290018da39SRichard Tran Mills $ and with the same numerical values. 4300018da39SRichard Tran Mills $ MAT_SHARE_NONZERO_PATTERN - Create a matrix that shares the nonzero structure with the previous matrix 4310018da39SRichard Tran Mills $ and does not copy it, using zeros for the numerical values. The parent and 4320018da39SRichard Tran Mills $ child matrices will share their index (i and j) arrays, and you cannot 4330018da39SRichard Tran Mills $ insert new nonzero entries into either matrix. 4340018da39SRichard Tran Mills 4350018da39SRichard Tran Mills Notes: Many matrix types (including SeqAIJ) do not support the MAT_SHARE_NONZERO_PATTERN optimization; in 4360018da39SRichard Tran Mills this case the behavior is as if MAT_DO_NOT_COPY_VALUES has been specified. 43770dcbbb9SBarry Smith 438d91e6319SBarry Smith .seealso: MatDuplicate() 439d91e6319SBarry Smith E*/ 44070dcbbb9SBarry Smith typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES,MAT_SHARE_NONZERO_PATTERN} MatDuplicateOption; 4412e8a6d31SBarry Smith 44219fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatConvert(Mat,MatType,MatReuse,Mat*); 443014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDuplicate(Mat,MatDuplicateOption,Mat*); 44494a9d846SBarry Smith 445cb5b572fSBarry Smith 446014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCopy(Mat,Mat,MatStructure); 447014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatView(Mat,PetscViewer); 448014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsSymmetric(Mat,PetscReal,PetscBool *); 449014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsStructurallySymmetric(Mat,PetscBool *); 450014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitian(Mat,PetscReal,PetscBool *); 451014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsSymmetricKnown(Mat,PetscBool *,PetscBool *); 452014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitianKnown(Mat,PetscBool *,PetscBool *); 453014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMissingDiagonal(Mat,PetscBool *,PetscInt *); 454014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLoad(Mat, PetscViewer); 4557b80b807SBarry Smith 4561a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool *); 4571a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool *); 4581a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool *); 4591a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool *); 460d4fbbf0eSBarry Smith 461d91e6319SBarry Smith /*S 462d91e6319SBarry Smith MatInfo - Context of matrix information, used with MatGetInfo() 463d91e6319SBarry Smith 464d91e6319SBarry Smith In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE 465d91e6319SBarry Smith 466d91e6319SBarry Smith Level: intermediate 467d91e6319SBarry Smith 468d91e6319SBarry Smith Concepts: matrix^nonzero information 469d91e6319SBarry Smith 470d9274352SBarry Smith .seealso: MatGetInfo(), MatInfoType 471d91e6319SBarry Smith S*/ 4724e220ebcSLois Curfman McInnes typedef struct { 473b0a32e0cSBarry Smith PetscLogDouble block_size; /* block size */ 474b0a32e0cSBarry Smith PetscLogDouble nz_allocated,nz_used,nz_unneeded; /* number of nonzeros */ 475b0a32e0cSBarry Smith PetscLogDouble memory; /* memory allocated */ 476b0a32e0cSBarry Smith PetscLogDouble assemblies; /* number of matrix assemblies called */ 477b0a32e0cSBarry Smith PetscLogDouble mallocs; /* number of mallocs during MatSetValues() */ 478b0a32e0cSBarry Smith PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */ 479b0a32e0cSBarry Smith PetscLogDouble factor_mallocs; /* number of mallocs during factorization */ 4804e220ebcSLois Curfman McInnes } MatInfo; 4814e220ebcSLois Curfman McInnes 482d9274352SBarry Smith /*E 483d9274352SBarry Smith MatInfoType - Indicates if you want information about the local part of the matrix, 484d9274352SBarry Smith the entire parallel matrix or the maximum over all the local parts. 485d9274352SBarry Smith 486d9274352SBarry Smith Level: beginner 487d9274352SBarry Smith 488af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 489d9274352SBarry Smith 490d9274352SBarry Smith .seealso: MatGetInfo(), MatInfo 491d9274352SBarry Smith E*/ 4927b80b807SBarry Smith typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType; 493014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetInfo(Mat,MatInfoType,MatInfo*); 494014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetDiagonal(Mat,Vec); 495014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMax(Mat,Vec,PetscInt[]); 496014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMin(Mat,Vec,PetscInt[]); 497014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMaxAbs(Mat,Vec,PetscInt[]); 498014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMinAbs(Mat,Vec,PetscInt[]); 499014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowSum(Mat,Vec); 500014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTranspose(Mat,MatReuse,Mat*); 501014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatHermitianTranspose(Mat,MatReuse,Mat*); 502014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPermute(Mat,IS,IS,Mat*); 503014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalScale(Mat,Vec,Vec); 504014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalSet(Mat,Vec,InsertMode); 505a52676ecSHong Zhang 506014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatEqual(Mat,Mat,PetscBool*); 507014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultEqual(Mat,Mat,PetscInt,PetscBool*); 508014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultAddEqual(Mat,Mat,PetscInt,PetscBool*); 509014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeEqual(Mat,Mat,PetscInt,PetscBool*); 510014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscBool*); 511a52676ecSHong Zhang PETSC_EXTERN PetscErrorCode MatMatMultEqual(Mat,Mat,Mat,PetscInt,PetscBool*); 512a52676ecSHong Zhang PETSC_EXTERN PetscErrorCode MatTransposeMatMultEqual(Mat,Mat,Mat,PetscInt,PetscBool*); 5137b80b807SBarry Smith 514014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNorm(Mat,NormType,PetscReal*); 515014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnNorms(Mat,NormType,PetscReal*); 516014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroEntries(Mat); 517014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 518014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsIS(Mat,IS,PetscScalar,Vec,Vec); 519014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsStencil(Mat,PetscInt,const MatStencil [],PetscScalar,Vec,Vec); 520014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsStencil(Mat,PetscInt,const MatStencil[],PetscScalar,Vec,Vec); 521014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumns(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 522014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsIS(Mat,IS,PetscScalar,Vec,Vec); 5237b80b807SBarry Smith 524014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSize(Mat,PetscInt*,PetscInt*); 525014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalSize(Mat,PetscInt*,PetscInt*); 526014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRange(Mat,PetscInt*,PetscInt*); 527014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRanges(Mat,const PetscInt**); 528014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRangeColumn(Mat,PetscInt*,PetscInt*); 529014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRangesColumn(Mat,const PetscInt**); 530566876eaSJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipIS(Mat,IS*,IS*); 5315ef9f2a5SBarry Smith 5327dae84e0SHong Zhang PETSC_EXTERN PetscErrorCode MatCreateSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]); 533cd2534ddSHong 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);} 5347dae84e0SHong Zhang PETSC_EXTERN PetscErrorCode MatCreateSubMatricesMPI(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]); 535cd2534ddSHong 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);} 536014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroyMatrices(PetscInt,Mat *[]); 537df750dc8SHong Zhang PETSC_EXTERN PetscErrorCode MatDestroySubMatrices(PetscInt,Mat *[]); 5387dae84e0SHong Zhang PETSC_EXTERN PetscErrorCode MatCreateSubMatrix(Mat,IS,IS,MatReuse,Mat *); 539cd2534ddSHong 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);} 540014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalSubMatrix(Mat,IS,IS,Mat*); 541014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreLocalSubMatrix(Mat,IS,IS,Mat*); 542014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSeqNonzeroStructure(Mat,Mat*); 543014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroySeqNonzeroStructure(Mat*); 5447b80b807SBarry Smith 545014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJ(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*); 546014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJSymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*); 547014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJNumeric(Mat,Mat); 548014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMat(Mat,MatReuse,Mat*); 549014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*); 550014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,Mat*); 551014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetGhosts(Mat, PetscInt *,const PetscInt *[]); 5528efafbd8SBarry Smith 553014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt); 554aa1e27eaSFande Kong PETSC_EXTERN PetscErrorCode MatIncreaseOverlapSplit(Mat mat,PetscInt n,IS is[],PetscInt ov); 555d2b2a242SBarry Smith PETSC_EXTERN PetscErrorCode MatMPIAIJSetUseScalableIncreaseOverlap(Mat,PetscBool); 5567b80b807SBarry Smith 557014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*); 558014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMultSymbolic(Mat,Mat,PetscReal,Mat*); 559014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMultNumeric(Mat,Mat,Mat); 56022440eb1SKris Buschelman 5617bab7c10SHong Zhang PETSC_EXTERN PetscErrorCode MatMatMatMult(Mat,Mat,Mat,MatReuse,PetscReal,Mat*); 5622df6c5c3SPatrick Farrell PETSC_EXTERN PetscErrorCode MatGalerkin(Mat,Mat,Mat,MatReuse,PetscReal,Mat*); 5637bab7c10SHong Zhang 564014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*); 565014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAPSymbolic(Mat,Mat,PetscReal,Mat*); 566014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAPNumeric(Mat,Mat,Mat); 567014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARt(Mat,Mat,MatReuse,PetscReal,Mat*); 568014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARtSymbolic(Mat,Mat,PetscReal,Mat*); 569014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARtNumeric(Mat,Mat,Mat); 57022440eb1SKris Buschelman 571014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeMatMult(Mat,Mat,MatReuse,PetscReal,Mat*); 572014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposetMatMultSymbolic(Mat,Mat,PetscReal,Mat*); 573014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposetMatMultNumeric(Mat,Mat,Mat); 574014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMult(Mat,Mat,MatReuse,PetscReal,Mat*); 575014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMultSymbolic(Mat,Mat,PetscReal,Mat*); 576014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMultNumeric(Mat,Mat,Mat); 577bc011b1eSHong Zhang 578014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAXPY(Mat,PetscScalar,Mat,MatStructure); 579014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAYPX(Mat,PetscScalar,Mat,MatStructure); 5801c741599SHong Zhang 581014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScale(Mat,PetscScalar); 582014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShift(Mat,PetscScalar); 5837b80b807SBarry Smith 584014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping); 585014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*); 586a93ff8c4SPeter Brune PETSC_EXTERN PetscErrorCode MatGetLayouts(Mat,PetscLayout*,PetscLayout*); 587014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 588014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsLocalIS(Mat,IS,PetscScalar,Vec,Vec); 589014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 590014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocalIS(Mat,IS,PetscScalar,Vec,Vec); 591014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 592014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 593052efed2SBarry Smith 594014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStashSetInitialSize(Mat,PetscInt,PetscInt); 595014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*); 59690f02eecSBarry Smith 597014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInterpolate(Mat,Vec,Vec); 598014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInterpolateAdd(Mat,Vec,Vec,Vec); 599014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestrict(Mat,Vec,Vec); 6002a7a6963SBarry Smith PETSC_EXTERN PetscErrorCode MatCreateVecs(Mat,Vec*,Vec*); 60184cc2a5aSBarry Smith PETSC_DEPRECATED("Use MatCreateVecs()") PETSC_STATIC_INLINE PetscErrorCode MatGetVecs(Mat mat,Vec *x,Vec *y) {return MatCreateVecs(mat,x,y);} 60253cd1579SHong Zhang PETSC_EXTERN PetscErrorCode MatCreateRedundantMatrix(Mat,PetscInt,MPI_Comm,MatReuse,Mat*); 603014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetMultiProcBlock(Mat,MPI_Comm,MatReuse,Mat*); 604014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFindZeroDiagonals(Mat,IS*); 6053a062f41SBarry Smith PETSC_EXTERN PetscErrorCode MatFindOffBlockDiagonalEntries(Mat,IS*); 6069c8f2541SHong Zhang PETSC_EXTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat(MPI_Comm,Mat,PetscInt,MatReuse,Mat*); 607bd481603SBarry Smith 608bd481603SBarry Smith /*MC 6091620fd73SBarry Smith MatSetValue - Set a single entry into a matrix. 6101620fd73SBarry Smith 6111620fd73SBarry Smith Not collective 6121620fd73SBarry Smith 613a9834a7dSSatish Balay Synopsis: 614a9834a7dSSatish Balay #include <petscmat.h> 615a9834a7dSSatish Balay PetscErrorCode MatSetValue(Mat m,PetscInt row,PetscInt col,PetscScalar value,InsertMode mode) 616a9834a7dSSatish Balay 6171620fd73SBarry Smith Input Parameters: 6181620fd73SBarry Smith + m - the matrix 6191620fd73SBarry Smith . row - the row location of the entry 6201620fd73SBarry Smith . col - the column location of the entry 6211620fd73SBarry Smith . value - the value to insert 6221620fd73SBarry Smith - mode - either INSERT_VALUES or ADD_VALUES 6231620fd73SBarry Smith 6241620fd73SBarry Smith Notes: 6251620fd73SBarry Smith For efficiency one should use MatSetValues() and set several or many 6261620fd73SBarry Smith values simultaneously if possible. 6271620fd73SBarry Smith 6281620fd73SBarry Smith Level: beginner 6291620fd73SBarry Smith 6301620fd73SBarry Smith .seealso: MatSetValues(), MatSetValueLocal() 6311620fd73SBarry Smith M*/ 6321620fd73SBarry 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);} 6331620fd73SBarry Smith 6341620fd73SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatGetValue(Mat v,PetscInt i,PetscInt j,PetscScalar *va) {return MatGetValues(v,1,&i,1,&j,va);} 6351620fd73SBarry Smith 6361620fd73SBarry 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);} 6371620fd73SBarry Smith 6381620fd73SBarry Smith /*MC 6390d2a7ff2SSatish Balay MatPreallocateInitialize - Begins the block of code that will count the number of nonzeros per 640bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 641bd481603SBarry Smith 642bd481603SBarry Smith Synopsis: 643aaa7dc30SBarry Smith #include <petscmat.h> 644c1ac3661SBarry Smith PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz) 645bd481603SBarry Smith 646bd481603SBarry Smith Collective on MPI_Comm 647bd481603SBarry Smith 648bd481603SBarry Smith Input Parameters: 649bd481603SBarry Smith + comm - the communicator that will share the eventually allocated matrix 650859fcb39SBarry Smith . nrows - the number of LOCAL rows in the matrix 651859fcb39SBarry Smith - ncols - the number of LOCAL columns in the matrix 652bd481603SBarry Smith 653bd481603SBarry Smith Output Parameters: 654bd481603SBarry Smith + dnz - the array that will be passed to the matrix preallocation routines 655bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 656bd481603SBarry Smith 657bd481603SBarry Smith Level: intermediate 658bd481603SBarry Smith 659bd481603SBarry Smith Notes: 660a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 661bd481603SBarry Smith 6621d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 663bd481603SBarry Smith 664bd481603SBarry Smith Use MatPreallocateInitializeSymmetric() for symmetric matrices (MPISBAIJ matrices) 665bd481603SBarry Smith 6661620fd73SBarry Smith This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize(). 6671620fd73SBarry Smith 668bd481603SBarry Smith Concepts: preallocation^Matrix 669bd481603SBarry Smith 670d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(), 671d6e23781SBarry Smith MatPreallocateInitializeSymmetric(), MatPreallocateSymmetricSetLocalBlock() 672bd481603SBarry Smith M*/ 673c4f061fbSSatish Balay #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \ 6747c922b88SBarry Smith { \ 675a89bc211SBarry Smith PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__start,__end; \ 6767282cfc1SLisandro Dalcin _4_ierr = PetscCalloc2((size_t)__nrows,&dnz,(size_t)__nrows,&onz);CHKERRQ(_4_ierr); \ 6771795a4d1SJed Brown __start = 0; __end = __start; \ 678c1ac3661SBarry Smith _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ctmp;\ 679a89bc211SBarry Smith _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows; 6807c922b88SBarry Smith 681bd481603SBarry Smith /*MC 682bd481603SBarry Smith MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be 683bd481603SBarry Smith inserted using a local number of the rows and columns 684bd481603SBarry Smith 685bd481603SBarry Smith Synopsis: 686aaa7dc30SBarry Smith #include <petscmat.h> 687c1ac3661SBarry Smith PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 688bd481603SBarry Smith 689bd481603SBarry Smith Not Collective 690bd481603SBarry Smith 691bd481603SBarry Smith Input Parameters: 692784ac674SJed Brown + map - the row mapping from local numbering to global numbering 693bd481603SBarry Smith . nrows - the number of rows indicated 6941d73ed98SBarry Smith . rows - the indices of the rows 695784ac674SJed Brown . cmap - the column mapping from local to global numbering 696bd481603SBarry Smith . ncols - the number of columns in the matrix 697bd481603SBarry Smith . cols - the columns indicated 698bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 699bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 700bd481603SBarry Smith 701bd481603SBarry Smith Level: intermediate 702bd481603SBarry Smith 703bd481603SBarry Smith Notes: 704a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 705bd481603SBarry Smith 7061d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 707bd481603SBarry Smith 708bd481603SBarry Smith Concepts: preallocation^Matrix 709bd481603SBarry Smith 710d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(), 711c1154cd5SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocalRemoveDups() 712bd481603SBarry Smith M*/ 713784ac674SJed Brown #define MatPreallocateSetLocal(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \ 714c4f061fbSSatish Balay {\ 715c1ac3661SBarry Smith PetscInt __l;\ 716784ac674SJed Brown _4_ierr = ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\ 717784ac674SJed Brown _4_ierr = ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\ 718c4f061fbSSatish Balay for (__l=0;__l<nrows;__l++) {\ 719ef66eb69SBarry Smith _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 720c4f061fbSSatish Balay }\ 721c4f061fbSSatish Balay } 722c4f061fbSSatish Balay 723bd481603SBarry Smith /*MC 724c1154cd5SBarry Smith MatPreallocateSetLocalRemoveDups - Indicates the locations (rows and columns) in the matrix where nonzeros will be 725c1154cd5SBarry Smith inserted using a local number of the rows and columns. This version removes any duplicate columns in cols 726c1154cd5SBarry Smith 727c1154cd5SBarry Smith Synopsis: 728c1154cd5SBarry Smith #include <petscmat.h> 729c1154cd5SBarry Smith PetscErrorCode MatPreallocateSetLocalRemoveDups(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 730c1154cd5SBarry Smith 731c1154cd5SBarry Smith Not Collective 732c1154cd5SBarry Smith 733c1154cd5SBarry Smith Input Parameters: 734c1154cd5SBarry Smith + map - the row mapping from local numbering to global numbering 735c1154cd5SBarry Smith . nrows - the number of rows indicated 736c1154cd5SBarry Smith . rows - the indices of the rows (these values are mapped to the global values) 737c1154cd5SBarry Smith . cmap - the column mapping from local to global numbering 738c1154cd5SBarry Smith . ncols - the number of columns in the matrix (this value will be changed if duplicate columns are found) 739c1154cd5SBarry Smith . cols - the columns indicated (these values are mapped to the global values, they are then sorted and duplicates removed) 740c1154cd5SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 741c1154cd5SBarry Smith - ozn - the other array passed to the matrix preallocation routines 742c1154cd5SBarry Smith 743c1154cd5SBarry Smith Level: intermediate 744c1154cd5SBarry Smith 745c1154cd5SBarry Smith Notes: 746c1154cd5SBarry Smith See Users-Manual: ch_performance for more details. 747c1154cd5SBarry Smith 748c1154cd5SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 749c1154cd5SBarry Smith 750c1154cd5SBarry Smith Concepts: preallocation^Matrix 751c1154cd5SBarry Smith 752c1154cd5SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(), 753c1154cd5SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocal() 754c1154cd5SBarry Smith M*/ 755c1154cd5SBarry Smith #define MatPreallocateSetLocalRemoveDups(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \ 756c1154cd5SBarry Smith {\ 757c1154cd5SBarry Smith PetscInt __l;\ 758c1154cd5SBarry Smith _4_ierr = ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\ 759c1154cd5SBarry Smith _4_ierr = ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\ 760c1154cd5SBarry Smith _4_ierr = PetscSortRemoveDupsInt(&ncols,cols);CHKERRQ(_4_ierr);\ 761c1154cd5SBarry Smith for (__l=0;__l<nrows;__l++) {\ 762c1154cd5SBarry Smith _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 763c1154cd5SBarry Smith }\ 764c1154cd5SBarry Smith } 765c1154cd5SBarry Smith 766c1154cd5SBarry Smith /*MC 767d6e23781SBarry Smith MatPreallocateSetLocalBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be 768bd481603SBarry Smith inserted using a local number of the rows and columns 769bd481603SBarry Smith 770bd481603SBarry Smith Synopsis: 771aaa7dc30SBarry Smith #include <petscmat.h> 772d6e23781SBarry Smith PetscErrorCode MatPreallocateSetLocalBlock(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 773d6e23781SBarry Smith 774d6e23781SBarry Smith Not Collective 775d6e23781SBarry Smith 776d6e23781SBarry Smith Input Parameters: 777d6e23781SBarry Smith + map - the row mapping from local numbering to global numbering 778d6e23781SBarry Smith . nrows - the number of rows indicated 779d6e23781SBarry Smith . rows - the indices of the rows 780d6e23781SBarry Smith . cmap - the column mapping from local to global numbering 781d6e23781SBarry Smith . ncols - the number of columns in the matrix 782d6e23781SBarry Smith . cols - the columns indicated 783d6e23781SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 784d6e23781SBarry Smith - ozn - the other array passed to the matrix preallocation routines 785d6e23781SBarry Smith 786d6e23781SBarry Smith Level: intermediate 787d6e23781SBarry Smith 788d6e23781SBarry Smith Notes: 789d6e23781SBarry Smith See Users-Manual: ch_performance for more details. 790d6e23781SBarry Smith 791d6e23781SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 792d6e23781SBarry Smith 793d6e23781SBarry Smith Concepts: preallocation^Matrix 794d6e23781SBarry Smith 795d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(), 796d6e23781SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock() 797d6e23781SBarry Smith M*/ 798d6e23781SBarry Smith #define MatPreallocateSetLocalBlock(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \ 799d6e23781SBarry Smith {\ 800d6e23781SBarry Smith PetscInt __l;\ 801d6e23781SBarry Smith _4_ierr = ISLocalToGlobalMappingApplyBlock(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\ 802d6e23781SBarry Smith _4_ierr = ISLocalToGlobalMappingApplyBlock(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\ 803d6e23781SBarry Smith for (__l=0;__l<nrows;__l++) {\ 804d6e23781SBarry Smith _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 805d6e23781SBarry Smith }\ 806d6e23781SBarry Smith } 807d6e23781SBarry Smith 808d6e23781SBarry Smith /*MC 809d6e23781SBarry Smith MatPreallocateSymmetricSetLocalBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be 810d6e23781SBarry Smith inserted using a local number of the rows and columns 811d6e23781SBarry Smith 812d6e23781SBarry Smith Synopsis: 813d6e23781SBarry Smith #include <petscmat.h> 814d6e23781SBarry Smith PetscErrorCode MatPreallocateSymmetricSetLocalBlock(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 815bd481603SBarry Smith 816bd481603SBarry Smith Not Collective 817bd481603SBarry Smith 818bd481603SBarry Smith Input Parameters: 819bd481603SBarry Smith + map - the mapping between local numbering and global numbering 820bd481603SBarry Smith . nrows - the number of rows indicated 8211d73ed98SBarry Smith . rows - the indices of the rows 822bd481603SBarry Smith . ncols - the number of columns in the matrix 823bd481603SBarry Smith . cols - the columns indicated 824bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 825bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 826bd481603SBarry Smith 827bd481603SBarry Smith Level: intermediate 828bd481603SBarry Smith 829bd481603SBarry Smith Notes: 830a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 831bd481603SBarry Smith 832bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 833bd481603SBarry Smith 834bd481603SBarry Smith Concepts: preallocation^Matrix 835bd481603SBarry Smith 836d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateInitialize(), 837d6e23781SBarry Smith MatPreallocateInitialize(), MatPreallocateSetLocal() 838bd481603SBarry Smith M*/ 839d6e23781SBarry Smith #define MatPreallocateSymmetricSetLocalBlock(map,nrows,rows,ncols,cols,dnz,onz) 0;\ 840d3d32019SBarry Smith {\ 841c1ac3661SBarry Smith PetscInt __l;\ 842d6e23781SBarry Smith _4_ierr = ISLocalToGlobalMappingApplyBlock(map,nrows,rows,rows);CHKERRQ(_4_ierr);\ 843d6e23781SBarry Smith _4_ierr = ISLocalToGlobalMappingApplyBlock(map,ncols,cols,cols);CHKERRQ(_4_ierr);\ 844d3d32019SBarry Smith for (__l=0;__l<nrows;__l++) {\ 845d6e23781SBarry Smith _4_ierr = MatPreallocateSymmetricSetBlock((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 846d3d32019SBarry Smith }\ 847d3d32019SBarry Smith } 848bd481603SBarry Smith /*MC 849bd481603SBarry Smith MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be 850bd481603SBarry Smith inserted using a local number of the rows and columns 851bd481603SBarry Smith 852bd481603SBarry Smith Synopsis: 853aaa7dc30SBarry Smith #include <petscmat.h> 854c1ac3661SBarry Smith PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 855bd481603SBarry Smith 856bd481603SBarry Smith Not Collective 857bd481603SBarry Smith 858bd481603SBarry Smith Input Parameters: 85964075487SBarry Smith + row - the row 860bd481603SBarry Smith . ncols - the number of columns in the matrix 861a50f8bf6SHong Zhang - cols - the columns indicated 862a50f8bf6SHong Zhang 863a50f8bf6SHong Zhang Output Parameters: 864a50f8bf6SHong Zhang + dnz - the array that will be passed to the matrix preallocation routines 865bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 866bd481603SBarry Smith 867bd481603SBarry Smith Level: intermediate 868bd481603SBarry Smith 869bd481603SBarry Smith Notes: 870a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 871bd481603SBarry Smith 872bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 873bd481603SBarry Smith 8741620fd73SBarry Smith This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize(). 8751620fd73SBarry Smith 876bd481603SBarry Smith Concepts: preallocation^Matrix 877bd481603SBarry Smith 878d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(), 879d6e23781SBarry Smith MatPreallocateInitialize(), MatPreallocateSetLocal() 880bd481603SBarry Smith M*/ 881c4f061fbSSatish Balay #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\ 882c1ac3661SBarry Smith { PetscInt __i; \ 883e32f2f54SBarry 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);\ 884e32f2f54SBarry 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);\ 8857c922b88SBarry Smith for (__i=0; __i<nc; __i++) {\ 88664075487SBarry Smith if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \ 8877cc688d7SBarry Smith else dnz[row - __rstart]++;\ 8887c922b88SBarry Smith }\ 8897c922b88SBarry Smith } 8907c922b88SBarry Smith 891bd481603SBarry Smith /*MC 892d6e23781SBarry Smith MatPreallocateSymmetricSetBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be 893bd481603SBarry Smith inserted using a local number of the rows and columns 894bd481603SBarry Smith 895bd481603SBarry Smith Synopsis: 896aaa7dc30SBarry Smith #include <petscmat.h> 897d6e23781SBarry Smith PetscErrorCode MatPreallocateSymmetricSetBlock(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 898bd481603SBarry Smith 899bd481603SBarry Smith Not Collective 900bd481603SBarry Smith 901bd481603SBarry Smith Input Parameters: 902bd481603SBarry Smith + nrows - the number of rows indicated 9031d73ed98SBarry Smith . rows - the indices of the rows 904bd481603SBarry Smith . ncols - the number of columns in the matrix 905bd481603SBarry Smith . cols - the columns indicated 906bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 907bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 908bd481603SBarry Smith 909bd481603SBarry Smith Level: intermediate 910bd481603SBarry Smith 911bd481603SBarry Smith Notes: 912a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 913bd481603SBarry Smith 914bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 915bd481603SBarry Smith 9161620fd73SBarry Smith This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize(). 9171620fd73SBarry Smith 918bd481603SBarry Smith Concepts: preallocation^Matrix 919bd481603SBarry Smith 920d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateInitialize(), 921d6e23781SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocal() 922bd481603SBarry Smith M*/ 923d6e23781SBarry Smith #define MatPreallocateSymmetricSetBlock(row,nc,cols,dnz,onz) 0;\ 924c1ac3661SBarry Smith { PetscInt __i; \ 925d3d32019SBarry Smith for (__i=0; __i<nc; __i++) {\ 926d3d32019SBarry Smith if (cols[__i] >= __end) onz[row - __rstart]++; \ 927d3d32019SBarry Smith else if (cols[__i] >= row) dnz[row - __rstart]++;\ 928d3d32019SBarry Smith }\ 929d3d32019SBarry Smith } 930d3d32019SBarry Smith 931bd481603SBarry Smith /*MC 93216371a99SBarry Smith MatPreallocateLocation - An alternative to MatPreallocationSet() that puts the nonzero locations into the matrix if it exists 93316371a99SBarry Smith 93416371a99SBarry Smith Synopsis: 935aaa7dc30SBarry Smith #include <petscmat.h> 93616371a99SBarry Smith PetscErrorCode MatPreallocateLocations(Mat A,PetscInt row,PetscInt ncols,PetscInt *cols,PetscInt *dnz,PetscInt *onz) 93716371a99SBarry Smith 93816371a99SBarry Smith Not Collective 93916371a99SBarry Smith 94016371a99SBarry Smith Input Parameters: 94116371a99SBarry Smith . A - matrix 94216371a99SBarry Smith . row - row where values exist (must be local to this process) 94316371a99SBarry Smith . ncols - number of columns 94416371a99SBarry Smith . cols - columns with nonzeros 94516371a99SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 94616371a99SBarry Smith - ozn - the other array passed to the matrix preallocation routines 94716371a99SBarry Smith 94816371a99SBarry Smith Level: intermediate 94916371a99SBarry Smith 95016371a99SBarry Smith Notes: 951a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 95216371a99SBarry Smith 95316371a99SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 95416371a99SBarry Smith 95516371a99SBarry Smith This is a MACRO not a function because it uses a bunch of variables private to the MatPreallocation.... routines. 95616371a99SBarry Smith 95716371a99SBarry Smith Concepts: preallocation^Matrix 95816371a99SBarry Smith 959d6e23781SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(), 960d6e23781SBarry Smith MatPreallocateSymmetricSetLocalBlock() 96116371a99SBarry Smith M*/ 9620298fd71SBarry 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);} 96316371a99SBarry Smith 96416371a99SBarry Smith 96516371a99SBarry Smith /*MC 9660d2a7ff2SSatish Balay MatPreallocateFinalize - Ends the block of code that will count the number of nonzeros per 967bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 968bd481603SBarry Smith 969bd481603SBarry Smith Synopsis: 970aaa7dc30SBarry Smith #include <petscmat.h> 971c1ac3661SBarry Smith PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz) 972bd481603SBarry Smith 973bd481603SBarry Smith Collective on MPI_Comm 974bd481603SBarry Smith 975bd481603SBarry Smith Input Parameters: 97616371a99SBarry Smith + dnz - the array that was be passed to the matrix preallocation routines 977bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 978bd481603SBarry Smith 979bd481603SBarry Smith Level: intermediate 980bd481603SBarry Smith 981bd481603SBarry Smith Notes: 982a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 983bd481603SBarry Smith 984bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 985bd481603SBarry Smith 9861620fd73SBarry Smith This is a MACRO not a function because it closes the { started in MatPreallocateInitialize(). 9871620fd73SBarry Smith 988bd481603SBarry Smith Concepts: preallocation^Matrix 989bd481603SBarry Smith 990d6e23781SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(), 991d6e23781SBarry Smith MatPreallocateSymmetricSetLocalBlock() 992bd481603SBarry Smith M*/ 993a89bc211SBarry Smith #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree2(dnz,onz);CHKERRQ(_4_ierr);} 9947c922b88SBarry Smith 9957b80b807SBarry Smith /* Routines unique to particular data structures */ 996014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellGetContext(Mat,void *); 997698d4c6aSKris Buschelman 998014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInodeAdjustForInodes(Mat,IS*,IS*); 999014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *); 1000698d4c6aSKris Buschelman 1001014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetColumnIndices(Mat,PetscInt[]); 1002014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetColumnIndices(Mat,PetscInt[]); 1003014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 1004014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 1005014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 1006014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJFromTriple(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*,PetscInt,PetscBool); 10077b80b807SBarry Smith 1008a96a251dSBarry Smith #define MAT_SKIP_ALLOCATION -4 1009a96a251dSBarry Smith 1010014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]); 1011014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]); 1012014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]); 1013273d9f13SBarry Smith 1014014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 1015014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 1016014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 1017014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []); 1018014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 1019014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]); 1020014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 1021014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]); 1022014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIDenseSetPreallocation(Mat,PetscScalar[]); 1023014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqDenseSetPreallocation(Mat,PetscScalar[]); 10249230625dSJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,const PetscInt*[]); 10259230625dSJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,const PetscInt*[]); 1026014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAdjCreateNonemptySubcommMat(Mat,Mat*); 1027273d9f13SBarry Smith 10282e1947a5SStefano Zampini PETSC_EXTERN PetscErrorCode MatISSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 1029cf0a3239SStefano Zampini PETSC_EXTERN PetscErrorCode MatISSetUpSF(Mat); 1030014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqDenseSetLDA(Mat,PetscInt); 1031014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDenseGetLocalMatrix(Mat,Mat*); 10321b807ce4Svictorle 1033014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStoreValues(Mat); 1034014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRetrieveValues(Mat); 10352e8a6d31SBarry Smith 1036014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDAADSetCtx(Mat,void*); 10373a7fca6bSBarry Smith 1038014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFindNonzeroRows(Mat,IS*); 10397cfe41e4SPatrick Farrell PETSC_EXTERN PetscErrorCode MatFindZeroRows(Mat,IS*); 10407b80b807SBarry Smith /* 10417b80b807SBarry Smith These routines are not usually accessed directly, rather solving is 104294b7f48cSBarry Smith done through the KSP and PC interfaces. 10437b80b807SBarry Smith */ 10447b80b807SBarry Smith 104576bdecfbSBarry Smith /*J 10468f6c3df8SBarry Smith MatOrderingType - String with the name of a PETSc matrix ordering 1047d9274352SBarry Smith 1048d9274352SBarry Smith Level: beginner 1049d9274352SBarry Smith 1050d9274352SBarry Smith .seealso: MatGetOrdering() 105176bdecfbSBarry Smith J*/ 105219fd82e9SBarry Smith typedef const char* MatOrderingType; 10532692d6eeSBarry Smith #define MATORDERINGNATURAL "natural" 10542692d6eeSBarry Smith #define MATORDERINGND "nd" 10552692d6eeSBarry Smith #define MATORDERING1WD "1wd" 10562692d6eeSBarry Smith #define MATORDERINGRCM "rcm" 10572692d6eeSBarry Smith #define MATORDERINGQMD "qmd" 10582692d6eeSBarry Smith #define MATORDERINGROWLENGTH "rowlength" 1059510184a7SJed Brown #define MATORDERINGWBM "wbm" 1060fc1bef75SJed Brown #define MATORDERINGSPECTRAL "spectral" 1061312e372aSBarry Smith #define MATORDERINGAMD "amd" /* only works if UMFPACK is installed with PETSc */ 1062b12f92e5SBarry Smith 106319fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatGetOrdering(Mat,MatOrderingType,IS*,IS*); 1064140e18c1SBarry Smith PETSC_EXTERN PetscErrorCode MatGetOrderingList(PetscFunctionList*); 1065bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatOrderingRegister(const char[],PetscErrorCode(*)(Mat,MatOrderingType,IS*,IS*)); 1066140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatOrderingList; 1067d4fbbf0eSBarry Smith 1068014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS); 1069fc1bef75SJed Brown PETSC_EXTERN PetscErrorCode MatCreateLaplacian(Mat,PetscReal,PetscBool,Mat*); 1070a2ce50c7SBarry Smith 1071d91e6319SBarry Smith /*S 1072d90ac83dSHong Zhang MatFactorShiftType - Numeric Shift. 1073d90ac83dSHong Zhang 1074d90ac83dSHong Zhang Level: beginner 1075d90ac83dSHong Zhang 1076d90ac83dSHong Zhang S*/ 1077d90ac83dSHong Zhang typedef enum {MAT_SHIFT_NONE,MAT_SHIFT_NONZERO,MAT_SHIFT_POSITIVE_DEFINITE,MAT_SHIFT_INBLOCKS} MatFactorShiftType; 10786a6fc655SJed Brown PETSC_EXTERN const char *const MatFactorShiftTypes[]; 10795e9742b9SJed Brown PETSC_EXTERN const char *const MatFactorShiftTypesDetail[]; 1080d90ac83dSHong Zhang 1081d90ac83dSHong Zhang /*S 10829aa7eafdSHong Zhang MatFactorError - indicates what type of error in matrix factor 10839aa7eafdSHong Zhang 10849aa7eafdSHong Zhang Level: beginner 10859aa7eafdSHong Zhang 108600e125f8SBarry Smith Developer Notes: Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 108700e125f8SBarry Smith 108800e125f8SBarry Smith .seealso: MatGetFactor() 10899aa7eafdSHong Zhang S*/ 10905cd7cf9dSHong Zhang typedef enum {MAT_FACTOR_NOERROR,MAT_FACTOR_STRUCT_ZEROPIVOT,MAT_FACTOR_NUMERIC_ZEROPIVOT,MAT_FACTOR_OUTMEMORY,MAT_FACTOR_OTHER} MatFactorError; 10919aa7eafdSHong Zhang 109200e125f8SBarry Smith PETSC_EXTERN PetscErrorCode MatFactorGetError(Mat,MatFactorError*); 1093b8b68cfdSBarry Smith PETSC_EXTERN PetscErrorCode MatFactorClearError(Mat); 10947b6c816cSBarry Smith PETSC_EXTERN PetscErrorCode MatFactorGetErrorZeroPivot(Mat,PetscReal*,PetscInt*); 109500e125f8SBarry Smith 10969aa7eafdSHong Zhang /*S 1097422a814eSBarry Smith MatFactorInfo - Data passed into the matrix factorization routines, and information about the resulting factorization 1098b00f7748SHong Zhang 109961cad860SBarry Smith In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE, that is use 110061cad860SBarry Smith $ MatFactorInfo info(MAT_FACTORINFO_SIZE) 1101b00f7748SHong Zhang 110215e8a5b3SHong Zhang Notes: These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC. 1103b00f7748SHong Zhang 110461cad860SBarry Smith You can use MatFactorInfoInitialize() to set default values. 110561cad860SBarry Smith 1106b00f7748SHong Zhang Level: developer 1107b00f7748SHong Zhang 1108d7d82daaSBarry Smith .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(), 1109d7d82daaSBarry Smith MatFactorInfoInitialize() 1110b00f7748SHong Zhang 1111b00f7748SHong Zhang S*/ 1112b00f7748SHong Zhang typedef struct { 111315e8a5b3SHong Zhang PetscReal diagonal_fill; /* force diagonal to fill in if initially not filled */ 111485317021SBarry Smith PetscReal usedt; 111515e8a5b3SHong Zhang PetscReal dt; /* drop tolerance */ 1116b00f7748SHong Zhang PetscReal dtcol; /* tolerance for pivoting */ 111715e8a5b3SHong Zhang PetscReal dtcount; /* maximum nonzeros to be allowed per row */ 111867e28bfeSBarry Smith PetscReal fill; /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */ 1119348344bbSBarry Smith PetscReal levels; /* ICC/ILU(levels) */ 1120bcd9e38bSBarry Smith PetscReal pivotinblocks; /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0 1121bcd9e38bSBarry Smith factorization may be faster if do not pivot */ 112215e8a5b3SHong Zhang PetscReal zeropivot; /* pivot is called zero if less than this */ 1123f4db908eSBarry Smith PetscReal shifttype; /* type of shift added to matrix factor to prevent zero pivots */ 1124f4db908eSBarry Smith PetscReal shiftamount; /* how large the shift is */ 112515e8a5b3SHong Zhang } MatFactorInfo; 1126ffa6d0a5SLois Curfman McInnes 1127014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFactorInfoInitialize(MatFactorInfo*); 11289a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatCholeskyFactor(Mat,IS,const MatFactorInfo*); 11299a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatCholeskyFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*); 11309a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatCholeskyFactorNumeric(Mat,Mat,const MatFactorInfo*); 11319a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatLUFactor(Mat,IS,IS,const MatFactorInfo*); 11329a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatILUFactor(Mat,IS,IS,const MatFactorInfo*); 11339a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatLUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*); 11349a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatILUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*); 11359a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatICCFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*); 11369a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatICCFactor(Mat,IS,const MatFactorInfo*); 11379a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatLUFactorNumeric(Mat,Mat,const MatFactorInfo*); 1138014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*); 1139014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolve(Mat,Vec,Vec); 1140014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatForwardSolve(Mat,Vec,Vec); 1141014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatBackwardSolve(Mat,Vec,Vec); 1142014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveAdd(Mat,Vec,Vec,Vec); 1143014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveTranspose(Mat,Vec,Vec); 1144014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveTransposeAdd(Mat,Vec,Vec,Vec); 1145014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolves(Mat,Vecs,Vecs); 11467c2f51b8SStefano Zampini PETSC_EXTERN PetscErrorCode MatSetUnfactored(Mat); 11477c2f51b8SStefano Zampini 11487c2f51b8SStefano Zampini typedef enum {MAT_FACTOR_SCHUR_UNFACTORED, MAT_FACTOR_SCHUR_FACTORED, MAT_FACTOR_SCHUR_INVERTED} MatFactorSchurStatus; 11498e7ba810SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorSetSchurIS(Mat,IS); 11507c2f51b8SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorGetSchurComplement(Mat,Mat*,MatFactorSchurStatus*); 11517c2f51b8SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorRestoreSchurComplement(Mat,Mat*,MatFactorSchurStatus); 11525a05ddb0SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorInvertSchurComplement(Mat); 11537c2f51b8SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorCreateSchurComplement(Mat,Mat*,MatFactorSchurStatus*); 11545a05ddb0SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorSolveSchurComplement(Mat,Vec,Vec); 11555a05ddb0SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorSolveSchurComplementTranspose(Mat,Vec,Vec); 11566dba178dSStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorFactorizeSchurComplement(Mat); 1157bb5a7306SBarry Smith 1158d91e6319SBarry Smith /*E 1159d91e6319SBarry Smith MatSORType - What type of (S)SOR to perform 1160bb1eb677SSatish Balay 1161d91e6319SBarry Smith Level: beginner 1162d91e6319SBarry Smith 1163d9274352SBarry Smith May be bitwise ORd together 1164d9274352SBarry Smith 1165af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 1166d91e6319SBarry Smith 11674e7234bfSBarry Smith MatSORType may be bitwise ORd together, so do not change the numbers 11684e7234bfSBarry Smith 116941f059aeSBarry Smith .seealso: MatSOR() 1170d91e6319SBarry Smith E*/ 1171ee50ffe9SBarry Smith typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3, 1172ee50ffe9SBarry Smith SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8, 1173ee50ffe9SBarry Smith SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16, 117484cb2905SBarry Smith SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType; 1175014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSOR(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec); 11768ed539a5SBarry Smith 1177d4fbbf0eSBarry Smith /* 1178639f9d9dSBarry Smith These routines are for efficiently computing Jacobians via finite differences. 1179639f9d9dSBarry Smith */ 1180b12f92e5SBarry Smith 11817086a01eSPeter Brune /*S 11827086a01eSPeter Brune MatColoring - Object for managing the coloring of matrices. 11837086a01eSPeter Brune 11847086a01eSPeter Brune Level: beginner 11857086a01eSPeter Brune 11867086a01eSPeter Brune Concepts: matrix, coloring 11877086a01eSPeter Brune 11887086a01eSPeter Brune .seealso: MatFDColoringCreate() ISColoring MatFDColoring 11897086a01eSPeter Brune S*/ 11907086a01eSPeter Brune typedef struct _p_MatColoring* MatColoring; 119176bdecfbSBarry Smith /*J 11928f6c3df8SBarry Smith MatColoringType - String with the name of a PETSc matrix coloring 1193d9274352SBarry Smith 1194d9274352SBarry Smith Level: beginner 1195d9274352SBarry Smith 11967086a01eSPeter Brune .seealso: MatColoringSetType(), MatColoring 119776bdecfbSBarry Smith J*/ 11987086a01eSPeter Brune 119919fd82e9SBarry Smith typedef const char* MatColoringType; 12005567d71aSPeter Brune #define MATCOLORINGJP "jp" 1201b9acb617SPeter Brune #define MATCOLORINGPOWER "power" 12022692d6eeSBarry Smith #define MATCOLORINGNATURAL "natural" 12032692d6eeSBarry Smith #define MATCOLORINGSL "sl" 12042692d6eeSBarry Smith #define MATCOLORINGLF "lf" 12052692d6eeSBarry Smith #define MATCOLORINGID "id" 12068121bdceSPeter Brune #define MATCOLORINGGREEDY "greedy" 1207b12f92e5SBarry Smith 12088ac6da0aSPeter Brune /*E 12098ac6da0aSPeter Brune MatColoringWeightType - Type of weight scheme 12108ac6da0aSPeter Brune 12118ac6da0aSPeter Brune Not Collective 12128ac6da0aSPeter Brune 12138ac6da0aSPeter Brune + MAT_COLORING_RANDOM - Random weights 12148ac6da0aSPeter Brune . MAT_COLORING_LEXICAL - Lexical weighting based upon global numbering. 12158ac6da0aSPeter Brune - MAT_COLORING_LF - Last-first weighting. 12168ac6da0aSPeter Brune 12178ac6da0aSPeter Brune Level: intermediate 12188ac6da0aSPeter Brune 1219af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 12208ac6da0aSPeter Brune 12218ac6da0aSPeter Brune .seealso: MatCUSPSetFormat(), MatCUSPFormatOperation 12228ac6da0aSPeter Brune E*/ 1223301c140bSPeter Brune typedef enum {MAT_COLORING_WEIGHT_RANDOM,MAT_COLORING_WEIGHT_LEXICAL,MAT_COLORING_WEIGHT_LF,MAT_COLORING_WEIGHT_SL} MatColoringWeightType; 12248ac6da0aSPeter Brune 1225335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringCreate(Mat,MatColoring*); 1226d7f2a307SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringGetDegrees(Mat,PetscInt,PetscInt*); 1227335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringDestroy(MatColoring*); 1228335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringView(MatColoring,PetscViewer); 1229335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetType(MatColoring,MatColoringType); 1230335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetFromOptions(MatColoring); 1231335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetDistance(MatColoring,PetscInt); 1232335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringGetDistance(MatColoring,PetscInt*); 1233335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetMaxColors(MatColoring,PetscInt); 1234335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringGetMaxColors(MatColoring,PetscInt*); 1235335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringApply(MatColoring,ISColoring*); 1236335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringRegister(const char[],PetscErrorCode(*)(MatColoring)); 1237014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*); 12388ac6da0aSPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetWeightType(MatColoring,MatColoringWeightType); 1239c29971fcSPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetWeights(MatColoring,PetscReal*,PetscInt*); 12408ac6da0aSPeter Brune PETSC_EXTERN PetscErrorCode MatColoringCreateWeights(MatColoring,PetscReal **,PetscInt **lperm); 1241cb394dc5SBarry Smith PETSC_EXTERN PetscErrorCode MatColoringTestValid(MatColoring,ISColoring); 1242639f9d9dSBarry Smith 1243d9274352SBarry Smith /*S 1244d9274352SBarry Smith MatFDColoring - Object for computing a sparse Jacobian via finite differences 1245d9274352SBarry Smith and coloring 1246639f9d9dSBarry Smith 1247d9274352SBarry Smith Level: beginner 1248d9274352SBarry Smith 1249d9274352SBarry Smith Concepts: coloring, sparse Jacobian, finite differences 1250d9274352SBarry Smith 1251d9274352SBarry Smith .seealso: MatFDColoringCreate() 1252d9274352SBarry Smith S*/ 1253e2a1c21fSSatish Balay typedef struct _p_MatFDColoring* MatFDColoring; 1254639f9d9dSBarry Smith 1255014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringCreate(Mat,ISColoring,MatFDColoring *); 1256014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringDestroy(MatFDColoring*); 1257014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringView(MatFDColoring,PetscViewer); 1258014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*); 1259014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**); 1260014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal); 1261014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetFromOptions(MatFDColoring); 1262d1e9a80fSBarry Smith PETSC_EXTERN PetscErrorCode MatFDColoringApply(Mat,MatFDColoring,Vec,void *); 1263014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetF(MatFDColoring,Vec); 1264edaa7c33SBarry Smith PETSC_EXTERN PetscErrorCode MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,const PetscInt*[]); 1265f86b9fbaSHong Zhang PETSC_EXTERN PetscErrorCode MatFDColoringSetUp(Mat,ISColoring,MatFDColoring); 1266f86b9fbaSHong Zhang PETSC_EXTERN PetscErrorCode MatFDColoringSetBlockSize(MatFDColoring,PetscInt,PetscInt); 1267f86b9fbaSHong Zhang 1268b1683b59SHong Zhang 1269b1683b59SHong Zhang /*S 1270b9af6bddSHong Zhang MatTransposeColoring - Object for computing a sparse matrix product C=A*B^T via coloring 1271b1683b59SHong Zhang 1272b1683b59SHong Zhang Level: beginner 1273b1683b59SHong Zhang 1274b1683b59SHong Zhang Concepts: coloring, sparse matrix product 1275b1683b59SHong Zhang 1276b9af6bddSHong Zhang .seealso: MatTransposeColoringCreate() 1277b1683b59SHong Zhang S*/ 1278b9af6bddSHong Zhang typedef struct _p_MatTransposeColoring* MatTransposeColoring; 1279b1683b59SHong Zhang 1280014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeColoringCreate(Mat,ISColoring,MatTransposeColoring *); 1281014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransColoringApplySpToDen(MatTransposeColoring,Mat,Mat); 1282014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransColoringApplyDenToSp(MatTransposeColoring,Mat,Mat); 1283014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeColoringDestroy(MatTransposeColoring*); 1284b1683b59SHong Zhang 1285639f9d9dSBarry Smith /* 12860752156aSBarry Smith These routines are for partitioning matrices: currently used only 12873eda8832SBarry Smith for adjacency matrix, MatCreateMPIAdj(). 12880752156aSBarry Smith */ 1289ca161407SBarry Smith 1290d9274352SBarry Smith /*S 1291d9274352SBarry Smith MatPartitioning - Object for managing the partitioning of a matrix or graph 1292d9274352SBarry Smith 1293d9274352SBarry Smith Level: beginner 1294d9274352SBarry Smith 1295d9274352SBarry Smith Concepts: partitioning 1296d9274352SBarry Smith 1297743a1026Svictorle .seealso: MatPartitioningCreate(), MatPartitioningType 1298d9274352SBarry Smith S*/ 129991e9ee9fSBarry Smith typedef struct _p_MatPartitioning* MatPartitioning; 1300d9274352SBarry Smith 130176bdecfbSBarry Smith /*J 13028f6c3df8SBarry Smith MatPartitioningType - String with the name of a PETSc matrix partitioning 1303d9274352SBarry Smith 1304d9274352SBarry Smith Level: beginner 13050d04baf8SBarry Smith dm 1306b547a13bSvictorle .seealso: MatPartitioningCreate(), MatPartitioning 130776bdecfbSBarry Smith J*/ 130819fd82e9SBarry Smith typedef const char* MatPartitioningType; 13092692d6eeSBarry Smith #define MATPARTITIONINGCURRENT "current" 1310aa1e27eaSFande Kong #define MATPARTITIONINGAVERAGE "average" 13112692d6eeSBarry Smith #define MATPARTITIONINGSQUARE "square" 13122692d6eeSBarry Smith #define MATPARTITIONINGPARMETIS "parmetis" 13132692d6eeSBarry Smith #define MATPARTITIONINGCHACO "chaco" 13142692d6eeSBarry Smith #define MATPARTITIONINGPARTY "party" 131550d91057SBarry Smith #define MATPARTITIONINGPTSCOTCH "ptscotch" 131688d2ac2bSFande Kong #define MATPARTITIONINGHIERARCH "hierarch" 131771306c60Sjroman 1318ca161407SBarry Smith 1319014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningCreate(MPI_Comm,MatPartitioning*); 132019fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningSetType(MatPartitioning,MatPartitioningType); 1321014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetNParts(MatPartitioning,PetscInt); 1322014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetAdjacency(MatPartitioning,Mat); 1323014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]); 1324014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []); 1325014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningApply(MatPartitioning,IS*); 1326014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningDestroy(MatPartitioning*); 13272aabb6bbSBarry Smith 1328bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningRegister(const char[],PetscErrorCode (*)(MatPartitioning)); 132930de9b25SBarry Smith 1330f1144a30SSatish Balay 13312bad1931SBarry Smith 1332014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningView(MatPartitioning,PetscViewer); 1333014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetFromOptions(MatPartitioning); 133419fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningGetType(MatPartitioning,MatPartitioningType*); 1335ca161407SBarry Smith 133627538973SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetRepartition(MatPartitioning part); 1337014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetCoarseSequential(MatPartitioning); 1338014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningParmetisGetEdgeCut(MatPartitioning, PetscInt *); 13390752156aSBarry Smith 1340b6956c12SJose Roman typedef enum { MP_CHACO_MULTILEVEL=1,MP_CHACO_SPECTRAL=2,MP_CHACO_LINEAR=4,MP_CHACO_RANDOM=5,MP_CHACO_SCATTERED=6 } MPChacoGlobalType; 13416a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoGlobalTypes[]; 1342b6956c12SJose Roman typedef enum { MP_CHACO_KERNIGHAN=1,MP_CHACO_NONE=2 } MPChacoLocalType; 13436a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoLocalTypes[]; 1344b6956c12SJose Roman typedef enum { MP_CHACO_LANCZOS=0,MP_CHACO_RQI=1 } MPChacoEigenType; 13456a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoEigenTypes[]; 1346b6956c12SJose Roman 1347014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetGlobal(MatPartitioning,MPChacoGlobalType); 1348014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetGlobal(MatPartitioning,MPChacoGlobalType*); 1349014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetLocal(MatPartitioning,MPChacoLocalType); 1350014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetLocal(MatPartitioning,MPChacoLocalType*); 1351014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal); 1352014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType); 1353014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenSolver(MatPartitioning,MPChacoEigenType*); 1354014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenTol(MatPartitioning,PetscReal); 1355014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenTol(MatPartitioning,PetscReal*); 1356014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenNumber(MatPartitioning,PetscInt); 1357014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenNumber(MatPartitioning,PetscInt*); 135871306c60Sjroman 135971306c60Sjroman #define MP_PARTY_OPT "opt" 136071306c60Sjroman #define MP_PARTY_LIN "lin" 136171306c60Sjroman #define MP_PARTY_SCA "sca" 136271306c60Sjroman #define MP_PARTY_RAN "ran" 136371306c60Sjroman #define MP_PARTY_GBF "gbf" 136471306c60Sjroman #define MP_PARTY_GCF "gcf" 136571306c60Sjroman #define MP_PARTY_BUB "bub" 136671306c60Sjroman #define MP_PARTY_DEF "def" 1367014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetGlobal(MatPartitioning,const char*); 136871306c60Sjroman #define MP_PARTY_HELPFUL_SETS "hs" 136971306c60Sjroman #define MP_PARTY_KERNIGHAN_LIN "kl" 137071306c60Sjroman #define MP_PARTY_NONE "no" 1371014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetLocal(MatPartitioning,const char*); 1372014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal); 1373014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetBipart(MatPartitioning,PetscBool); 1374014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscBool); 137571306c60Sjroman 137650d91057SBarry Smith typedef enum { MP_PTSCOTCH_QUALITY,MP_PTSCOTCH_SPEED,MP_PTSCOTCH_BALANCE,MP_PTSCOTCH_SAFETY,MP_PTSCOTCH_SCALABILITY } MPPTScotchStrategyType; 13776a6fc655SJed Brown PETSC_EXTERN const char *const MPPTScotchStrategyTypes[]; 1378e0f1cffaSJose Roman 1379014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetImbalance(MatPartitioning,PetscReal); 1380014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetImbalance(MatPartitioning,PetscReal*); 1381014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetStrategy(MatPartitioning,MPPTScotchStrategyType); 1382014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetStrategy(MatPartitioning,MPPTScotchStrategyType*); 138371306c60Sjroman 1384b43b03e9SMark F. Adams /* 13856294fa2bSFande Kong * hierarchical partitioning 13866294fa2bSFande Kong */ 13874e713f55SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalGetFineparts(MatPartitioning,IS*); 13884e713f55SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalGetCoarseparts(MatPartitioning,IS*); 13894e713f55SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalSetNcoarseparts(MatPartitioning,PetscInt); 13904e713f55SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalSetNfineparts(MatPartitioning, PetscInt); 13916294fa2bSFande Kong 1392014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMeshToVertexGraph(Mat,PetscInt,Mat*); 1393014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMeshToCellGraph(Mat,PetscInt,Mat*); 1394591294e4SBarry Smith 13950752156aSBarry Smith /* 1396af0996ceSBarry Smith If you add entries here you must also add them to petsc/finclude/petscmat.h 1397d4fbbf0eSBarry Smith */ 13981c1c02c0SLois Curfman McInnes typedef enum { MATOP_SET_VALUES=0, 13991c1c02c0SLois Curfman McInnes MATOP_GET_ROW=1, 14001c1c02c0SLois Curfman McInnes MATOP_RESTORE_ROW=2, 14011c1c02c0SLois Curfman McInnes MATOP_MULT=3, 14021c1c02c0SLois Curfman McInnes MATOP_MULT_ADD=4, 14037c922b88SBarry Smith MATOP_MULT_TRANSPOSE=5, 14047c922b88SBarry Smith MATOP_MULT_TRANSPOSE_ADD=6, 14051c1c02c0SLois Curfman McInnes MATOP_SOLVE=7, 14061c1c02c0SLois Curfman McInnes MATOP_SOLVE_ADD=8, 14077c922b88SBarry Smith MATOP_SOLVE_TRANSPOSE=9, 14087c922b88SBarry Smith MATOP_SOLVE_TRANSPOSE_ADD=10, 14091c1c02c0SLois Curfman McInnes MATOP_LUFACTOR=11, 14101c1c02c0SLois Curfman McInnes MATOP_CHOLESKYFACTOR=12, 1411a714c06dSBarry Smith MATOP_SOR=13, 14121c1c02c0SLois Curfman McInnes MATOP_TRANSPOSE=14, 14131c1c02c0SLois Curfman McInnes MATOP_GETINFO=15, 14141c1c02c0SLois Curfman McInnes MATOP_EQUAL=16, 14151c1c02c0SLois Curfman McInnes MATOP_GET_DIAGONAL=17, 14161c1c02c0SLois Curfman McInnes MATOP_DIAGONAL_SCALE=18, 14171c1c02c0SLois Curfman McInnes MATOP_NORM=19, 14181c1c02c0SLois Curfman McInnes MATOP_ASSEMBLY_BEGIN=20, 14191c1c02c0SLois Curfman McInnes MATOP_ASSEMBLY_END=21, 1420d519adbfSMatthew Knepley MATOP_SET_OPTION=22, 1421d519adbfSMatthew Knepley MATOP_ZERO_ENTRIES=23, 1422d519adbfSMatthew Knepley MATOP_ZERO_ROWS=24, 1423d519adbfSMatthew Knepley MATOP_LUFACTOR_SYMBOLIC=25, 1424d519adbfSMatthew Knepley MATOP_LUFACTOR_NUMERIC=26, 1425d519adbfSMatthew Knepley MATOP_CHOLESKY_FACTOR_SYMBOLIC=27, 1426d519adbfSMatthew Knepley MATOP_CHOLESKY_FACTOR_NUMERIC=28, 1427d519adbfSMatthew Knepley MATOP_SETUP_PREALLOCATION=29, 1428d519adbfSMatthew Knepley MATOP_ILUFACTOR_SYMBOLIC=30, 1429d519adbfSMatthew Knepley MATOP_ICCFACTOR_SYMBOLIC=31, 1430a5b7ff6bSBarry Smith MATOP_GET_DIAGONAL_BLOCK=32, 14319d39f69dSJed Brown /* MATOP_PLACEHOLDER_33=33, */ 1432d519adbfSMatthew Knepley MATOP_DUPLICATE=34, 1433d519adbfSMatthew Knepley MATOP_FORWARD_SOLVE=35, 1434d519adbfSMatthew Knepley MATOP_BACKWARD_SOLVE=36, 1435d519adbfSMatthew Knepley MATOP_ILUFACTOR=37, 1436d519adbfSMatthew Knepley MATOP_ICCFACTOR=38, 1437d519adbfSMatthew Knepley MATOP_AXPY=39, 14387dae84e0SHong Zhang MATOP_CREATE_SUBMATRICES=40, 1439d519adbfSMatthew Knepley MATOP_INCREASE_OVERLAP=41, 1440d519adbfSMatthew Knepley MATOP_GET_VALUES=42, 1441d519adbfSMatthew Knepley MATOP_COPY=43, 1442d519adbfSMatthew Knepley MATOP_GET_ROW_MAX=44, 1443d519adbfSMatthew Knepley MATOP_SCALE=45, 1444d519adbfSMatthew Knepley MATOP_SHIFT=46, 144535153367SBarry Smith MATOP_DIAGONAL_SET=47, 14469d39f69dSJed Brown MATOP_ZERO_ROWS_COLUMNS=48, 14479d39f69dSJed Brown MATOP_SET_RANDOM=49, 1448d519adbfSMatthew Knepley MATOP_GET_ROW_IJ=50, 1449d519adbfSMatthew Knepley MATOP_RESTORE_ROW_IJ=51, 1450d519adbfSMatthew Knepley MATOP_GET_COLUMN_IJ=52, 1451d519adbfSMatthew Knepley MATOP_RESTORE_COLUMN_IJ=53, 1452d519adbfSMatthew Knepley MATOP_FDCOLORING_CREATE=54, 1453d519adbfSMatthew Knepley MATOP_COLORING_PATCH=55, 1454d519adbfSMatthew Knepley MATOP_SET_UNFACTORED=56, 1455d519adbfSMatthew Knepley MATOP_PERMUTE=57, 1456d519adbfSMatthew Knepley MATOP_SET_VALUES_BLOCKED=58, 14577dae84e0SHong Zhang MATOP_CREATE_SUBMATRIX=59, 1458d519adbfSMatthew Knepley MATOP_DESTROY=60, 1459d519adbfSMatthew Knepley MATOP_VIEW=61, 1460d519adbfSMatthew Knepley MATOP_CONVERT_FROM=62, 14617bab7c10SHong Zhang MATOP_MATMAT_MULT=63, 14627bab7c10SHong Zhang MATOP_MATMAT_MULT_SYMBOLIC=64, 14637bab7c10SHong Zhang MATOP_MATMAT_MULT_NUMERIC=65, 1464d519adbfSMatthew Knepley MATOP_SET_LOCAL_TO_GLOBAL_MAP=66, 1465d519adbfSMatthew Knepley MATOP_SET_VALUES_LOCAL=67, 1466d519adbfSMatthew Knepley MATOP_ZERO_ROWS_LOCAL=68, 1467d519adbfSMatthew Knepley MATOP_GET_ROW_MAX_ABS=69, 1468d519adbfSMatthew Knepley MATOP_GET_ROW_MIN_ABS=70, 1469d519adbfSMatthew Knepley MATOP_CONVERT=71, 1470d519adbfSMatthew Knepley MATOP_SET_COLORING=72, 14719d39f69dSJed Brown /* MATOP_PLACEHOLDER_73=73, */ 1472d519adbfSMatthew Knepley MATOP_SET_VALUES_ADIFOR=74, 1473d519adbfSMatthew Knepley MATOP_FD_COLORING_APPLY=75, 1474d519adbfSMatthew Knepley MATOP_SET_FROM_OPTIONS=76, 1475bda6c4cbSBarry Smith MATOP_MULT_CONSTRAINED=77, 1476bda6c4cbSBarry Smith MATOP_MULT_TRANSPOSE_CONSTRAIN=78, 14779d39f69dSJed Brown MATOP_FIND_ZERO_DIAGONALS=79, 1478d519adbfSMatthew Knepley MATOP_MULT_MULTIPLE=80, 1479d519adbfSMatthew Knepley MATOP_SOLVE_MULTIPLE=81, 1480d519adbfSMatthew Knepley MATOP_GET_INERTIA=82, 1481d519adbfSMatthew Knepley MATOP_LOAD=83, 1482d519adbfSMatthew Knepley MATOP_IS_SYMMETRIC=84, 1483d519adbfSMatthew Knepley MATOP_IS_HERMITIAN=85, 1484d519adbfSMatthew Knepley MATOP_IS_STRUCTURALLY_SYMMETRIC=86, 1485820469bcSHong Zhang MATOP_SET_VALUES_BLOCKEDLOCAL=87, 1486b41ce5d5SBarry Smith MATOP_CREATE_VECS=88, 1487d519adbfSMatthew Knepley MATOP_MAT_MULT=89, 1488d519adbfSMatthew Knepley MATOP_MAT_MULT_SYMBOLIC=90, 1489d519adbfSMatthew Knepley MATOP_MAT_MULT_NUMERIC=91, 1490d519adbfSMatthew Knepley MATOP_PTAP=92, 1491d519adbfSMatthew Knepley MATOP_PTAP_SYMBOLIC=93, 1492d519adbfSMatthew Knepley MATOP_PTAP_NUMERIC=94, 1493bda6c4cbSBarry Smith MATOP_MAT_TRANSPOSE_MULT=95, 1494bda6c4cbSBarry Smith MATOP_MAT_TRANSPOSE_MULT_SYMBO=96, 1495bda6c4cbSBarry Smith MATOP_MAT_TRANSPOSE_MULT_NUMER=97, 14969d39f69dSJed Brown /* MATOP_PLACEHOLDER_98=98, */ 14979d39f69dSJed Brown /* MATOP_PLACEHOLDER_99=99, */ 14989d39f69dSJed Brown /* MATOP_PLACEHOLDER_100=100, */ 14999d39f69dSJed Brown /* MATOP_PLACEHOLDER_101=101, */ 1500d519adbfSMatthew Knepley MATOP_CONJUGATE=102, 15019d39f69dSJed Brown /* MATOP_PLACEHOLDER_103=103, */ 1502d519adbfSMatthew Knepley MATOP_SET_VALUES_ROW=104, 1503d519adbfSMatthew Knepley MATOP_REAL_PART=105, 1504bda6c4cbSBarry Smith MATOP_IMAGINARY_PART=106, 1505bda6c4cbSBarry Smith MATOP_GET_ROW_UPPER_TRIANGULAR=107, 1506bda6c4cbSBarry Smith MATOP_RESTORE_ROW_UPPER_TRIANG=108, 1507bda6c4cbSBarry Smith MATOP_MAT_SOLVE=109, 150829eadf9eSStefano Zampini MATOP_MAT_SOLVE_TRANSPOSE=110, 1509d519adbfSMatthew Knepley MATOP_GET_ROW_MIN=111, 15100e8d04f7SBarry Smith MATOP_GET_COLUMN_VECTOR=112, 1511d519adbfSMatthew Knepley MATOP_MISSING_DIAGONAL=113, 15120e8d04f7SBarry Smith MATOP_GET_SEQ_NONZERO_STRUCTUR=114, 151389c6957cSBarry Smith MATOP_CREATE=115, 151489c6957cSBarry Smith MATOP_GET_GHOSTS=116, 15150e8d04f7SBarry Smith MATOP_GET_LOCAL_SUB_MATRIX=117, 15160e8d04f7SBarry Smith MATOP_RESTORE_LOCALSUB_MATRIX=118, 1517eeffb40dSHong Zhang MATOP_MULT_DIAGONAL_BLOCK=119, 15180e8d04f7SBarry Smith MATOP_HERMITIAN_TRANSPOSE=120, 15190e8d04f7SBarry Smith MATOP_MULT_HERMITIAN_TRANSPOSE=121, 15200e8d04f7SBarry Smith MATOP_MULT_HERMITIAN_TRANS_ADD=122, 15210e8d04f7SBarry Smith MATOP_GET_MULTI_PROC_BLOCK=123, 15229d39f69dSJed Brown MATOP_FIND_NONZERO_ROWS=124, 15230e8d04f7SBarry Smith MATOP_GET_COLUMN_NORMS=125, 15249d39f69dSJed Brown MATOP_INVERT_BLOCK_DIAGONAL=126, 15259d39f69dSJed Brown /* MATOP_PLACEHOLDER_127=127, */ 1526b41ce5d5SBarry Smith MATOP_CREATE_SUB_MATRICES_MPI=128, 15272b8ad9a3SHong Zhang MATOP_SET_VALUES_BATCH=129, 15280e8d04f7SBarry Smith MATOP_TRANSPOSE_MAT_MULT=130, 15290e8d04f7SBarry Smith MATOP_TRANSPOSE_MAT_MULT_SYMBO=131, 1530e9b602ebSSatish Balay MATOP_TRANSPOSE_MAT_MULT_NUMER=132, 15310e8d04f7SBarry Smith MATOP_TRANSPOSE_COLORING_CREAT=133, 15320e8d04f7SBarry Smith MATOP_TRANS_COLORING_APPLY_SPT=134, 15330e8d04f7SBarry Smith MATOP_TRANS_COLORING_APPLY_DEN=135, 15340e8d04f7SBarry Smith MATOP_RART=136, 15350e8d04f7SBarry Smith MATOP_RART_SYMBOLIC=137, 15360e8d04f7SBarry Smith MATOP_RART_NUMERIC=138, 1537e09a3074SHong Zhang MATOP_SET_BLOCK_SIZES=139, 1538f9426fe0SMark Adams MATOP_AYPX=140, 15391919a2e2SJed Brown MATOP_RESIDUAL=141, 15409c8f2541SHong Zhang MATOP_FDCOLORING_SETUP=142, 15412d033e1fSHong Zhang MATOP_MPICONCATENATESEQ=144, 15422d033e1fSHong Zhang MATOP_DESTROYSUBMATRICES=145 1543fae171e0SBarry Smith } MatOperation; 1544014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatHasOperation(Mat,MatOperation,PetscBool *); 1545014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellSetOperation(Mat,MatOperation,void(*)(void)); 1546014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellGetOperation(Mat,MatOperation,void(**)(void)); 1547014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellSetContext(Mat,void*); 1548112a2221SBarry Smith 154990ace30eSBarry Smith /* 155090ace30eSBarry Smith Codes for matrices stored on disk. By default they are 155190ace30eSBarry Smith stored in a universal format. By changing the format with 15526a9046bcSBarry Smith PetscViewerPushFormat(viewer,PETSC_VIEWER_NATIVE); the matrices will 155390ace30eSBarry Smith be stored in a way natural for the matrix, for example dense matrices 155490ace30eSBarry Smith would be stored as dense. Matrices stored this way may only be 1555f4403165SShri Abhyankar read into matrices of the same type. 155690ace30eSBarry Smith */ 155790ace30eSBarry Smith #define MATRIX_BINARY_FORMAT_DENSE -1 155890ace30eSBarry Smith 1559014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetHashTableFactor(Mat,PetscReal); 1560014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatISGetLocalMat(Mat,Mat*); 15613b3b1effSJed Brown PETSC_EXTERN PetscErrorCode MatISRestoreLocalMat(Mat,Mat*); 1562014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatISSetLocalMat(Mat,Mat); 1563b7ce53b6SStefano Zampini PETSC_EXTERN PetscErrorCode MatISGetMPIXAIJ(Mat,MatReuse,Mat*); 15641f4e1ec7SBarry Smith 1565d9274352SBarry Smith /*S 1566d9274352SBarry Smith MatNullSpace - Object that removes a null space from a vector, i.e. 1567d9274352SBarry Smith orthogonalizes the vector to a subsapce 1568d9274352SBarry Smith 1569f7a9e4ceSBarry Smith Level: advanced 1570d9274352SBarry Smith 1571d9274352SBarry Smith Concepts: matrix; linear operator, null space 1572d9274352SBarry Smith 15736e1639daSBarry Smith Users manual sections: 15746e1639daSBarry Smith . sec_singular 15756e1639daSBarry Smith 1576d9274352SBarry Smith .seealso: MatNullSpaceCreate() 1577d9274352SBarry Smith S*/ 157874637425SBarry Smith typedef struct _p_MatNullSpace* MatNullSpace; 1579d9274352SBarry Smith 1580014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceCreate(MPI_Comm,PetscBool ,PetscInt,const Vec[],MatNullSpace*); 1581014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(MatNullSpace,Vec,void*),void*); 1582014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceDestroy(MatNullSpace*); 1583d0195637SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceRemove(MatNullSpace,Vec); 1584014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNullSpace(Mat, MatNullSpace *); 15855fa7ec2dSBarry Smith PETSC_EXTERN PetscErrorCode MatGetTransposeNullSpace(Mat, MatNullSpace *); 15865fa7ec2dSBarry Smith PETSC_EXTERN PetscErrorCode MatSetTransposeNullSpace(Mat,MatNullSpace); 1587014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNullSpace(Mat,MatNullSpace); 1588014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNearNullSpace(Mat,MatNullSpace); 1589014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNearNullSpace(Mat,MatNullSpace*); 1590014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceTest(MatNullSpace,Mat,PetscBool *); 1591014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceView(MatNullSpace,PetscViewer); 1592014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceGetVecs(MatNullSpace,PetscBool*,PetscInt*,const Vec**); 1593014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceCreateRigidBody(Vec,MatNullSpace*); 159474637425SBarry Smith 1595014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatReorderingSeqSBAIJ(Mat,IS); 1596014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetHashTableFactor(Mat,PetscReal); 1597014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetColumnIndices(Mat,PetscInt *); 1598014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJInvertBlockDiagonal(Mat); 15993f1d51d7SBarry Smith 1600014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMAIJ(Mat,PetscInt,Mat*); 1601014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMAIJRedimension(Mat,PetscInt,Mat*); 1602014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMAIJGetAIJ(Mat,Mat*); 1603c4f061fbSSatish Balay 1604014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatComputeExplicitOperator(Mat,Mat*); 1605b0a32e0cSBarry Smith 1606014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalScaleLocal(Mat,Vec); 160704f1ad80SBarry Smith 1608014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*); 1609014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetBase(Mat,Vec,Vec); 1610014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*); 1611014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*)); 1612014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec)); 1613014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt); 1614014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDResetHHistory(Mat); 1615014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctionError(Mat,PetscReal); 1616014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetPeriod(Mat,PetscInt); 1617014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDGetH(Mat,PetscScalar *); 1618014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetOptionsPrefix(Mat,const char[]); 1619014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*); 1620014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*); 1621e884886fSBarry Smith 16226370053bSBarry Smith /*S 16236370053bSBarry Smith MatMFFD - A data structured used to manage the computation of the h differencing parameter for matrix-free 16246370053bSBarry Smith Jacobian vector products 1625e884886fSBarry Smith 16266370053bSBarry Smith Notes: MATMFFD is a specific MatType which uses the MatMFFD data structure 16276370053bSBarry Smith 16286370053bSBarry Smith MatMFFD*() methods actually take the Mat as their first argument. Not a MatMFFD data structure 16296370053bSBarry Smith 16306370053bSBarry Smith Level: developer 16316370053bSBarry Smith 16326370053bSBarry Smith .seealso: MATMFFD, MatCreateMFFD(), MatMFFDSetFuction(), MatMFFDSetType(), MatMFFDRegister() 16336370053bSBarry Smith S*/ 1634e884886fSBarry Smith typedef struct _p_MatMFFD* MatMFFD; 1635e884886fSBarry Smith 163676bdecfbSBarry Smith /*J 1637e884886fSBarry Smith MatMFFDType - algorithm used to compute the h used in computing matrix-vector products via differencing of the function 1638e884886fSBarry Smith 1639e884886fSBarry Smith Level: beginner 1640e884886fSBarry Smith 1641e884886fSBarry Smith .seealso: MatMFFDSetType(), MatMFFDRegister() 164276bdecfbSBarry Smith J*/ 164319fd82e9SBarry Smith typedef const char* MatMFFDType; 1644e884886fSBarry Smith #define MATMFFD_DS "ds" 1645e884886fSBarry Smith #define MATMFFD_WP "wp" 1646e884886fSBarry Smith 164719fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDSetType(Mat,MatMFFDType); 1648bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDRegister(const char[],PetscErrorCode (*)(MatMFFD)); 1649e884886fSBarry Smith 1650014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDDSSetUmin(Mat,PetscReal); 1651014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDWPSetComputeNormU(Mat,PetscBool ); 1652e884886fSBarry Smith 165361ab5df0SBarry Smith PETSC_EXTERN PetscErrorCode MatFDColoringSetType(MatFDColoring,MatMFFDType); 165461ab5df0SBarry Smith 1655014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *); 1656014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *); 16577dbadf16SMatthew Knepley 165897969023SHong Zhang /* 165997969023SHong Zhang PETSc interface to MUMPS 166097969023SHong Zhang */ 166197969023SHong Zhang #ifdef PETSC_HAVE_MUMPS 1662014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMumpsSetIcntl(Mat,PetscInt,PetscInt); 1663bc6112feSHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetIcntl(Mat,PetscInt,PetscInt*); 1664b9e9ea26SSatish Balay PETSC_EXTERN PetscErrorCode MatMumpsSetCntl(Mat,PetscInt,PetscReal); 1665bc6112feSHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetCntl(Mat,PetscInt,PetscReal*); 1666bc6112feSHong Zhang 1667ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetInfo(Mat,PetscInt,PetscInt*); 1668ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetInfog(Mat,PetscInt,PetscInt*); 1669ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetRinfo(Mat,PetscInt,PetscReal*); 1670ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetRinfog(Mat,PetscInt,PetscReal*); 167197969023SHong Zhang #endif 167223a5497aSJed Brown 1673d954db57SHong Zhang /* 1674d615f992SSatish Balay PETSc interface to Mkl_Pardiso 1675b500a6b7SJose David Bermeol */ 1676b500a6b7SJose David Bermeol #ifdef PETSC_HAVE_MKL_PARDISO 1677d615f992SSatish Balay PETSC_EXTERN PetscErrorCode MatMkl_PardisoSetCntl(Mat,PetscInt,PetscInt); 1678b500a6b7SJose David Bermeol #endif 1679b500a6b7SJose David Bermeol 1680b500a6b7SJose David Bermeol /* 1681d305a81bSVasiliy Kozyrev PETSc interface to Mkl_CPardiso 1682d305a81bSVasiliy Kozyrev */ 1683d305a81bSVasiliy Kozyrev #ifdef PETSC_HAVE_MKL_CPARDISO 1684d305a81bSVasiliy Kozyrev PETSC_EXTERN PetscErrorCode MatMkl_CPardisoSetCntl(Mat,PetscInt,PetscInt); 1685d305a81bSVasiliy Kozyrev #endif 1686d305a81bSVasiliy Kozyrev 1687d305a81bSVasiliy Kozyrev /* 1688d954db57SHong Zhang PETSc interface to SUPERLU 1689d954db57SHong Zhang */ 1690d954db57SHong Zhang #ifdef PETSC_HAVE_SUPERLU 1691014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSuperluSetILUDropTol(Mat,PetscReal); 1692d954db57SHong Zhang #endif 1693d954db57SHong Zhang 1694fb785984SHong Zhang /* 1695fb785984SHong Zhang PETSc interface to SUPERLU_DIST 1696fb785984SHong Zhang */ 1697fb785984SHong Zhang #ifdef PETSC_HAVE_SUPERLU_DIST 1698fb785984SHong Zhang PETSC_EXTERN PetscErrorCode MatSuperluDistGetDiagU(Mat,PetscScalar*); 1699fb785984SHong Zhang #endif 1700fb785984SHong Zhang 1701575704cbSPieter Ghysels /* 1702575704cbSPieter Ghysels PETSc interface to STRUMPACK 1703575704cbSPieter Ghysels */ 1704575704cbSPieter Ghysels #ifdef PETSC_HAVE_STRUMPACK 1705575704cbSPieter Ghysels PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSRelCompTol(Mat,PetscReal); 1706575704cbSPieter Ghysels PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSMinSize(Mat,PetscInt); 1707575704cbSPieter Ghysels PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetColPerm(Mat,PetscBool); 1708575704cbSPieter Ghysels #endif 1709575704cbSPieter Ghysels 1710575704cbSPieter Ghysels 1711aa372e3fSPaul Mullowney #ifdef PETSC_HAVE_CUDA 17123f9c0db1SPaul Mullowney /*E 1713e057df02SPaul Mullowney MatCUSPARSEStorageFormat - indicates the storage format for CUSPARSE (GPU) 17142692e278SPaul Mullowney matrices. 1715e057df02SPaul Mullowney 1716e057df02SPaul Mullowney Not Collective 1717e057df02SPaul Mullowney 17188468deeeSKarl Rupp + MAT_CUSPARSE_CSR - Compressed Sparse Row 17192692e278SPaul Mullowney . MAT_CUSPARSE_ELL - Ellpack (requires CUDA 4.2 or later). 17202692e278SPaul Mullowney - MAT_CUSPARSE_HYB - Hybrid, a combination of Ellpack and Coordinate format (requires CUDA 4.2 or later). 1721e057df02SPaul Mullowney 1722e057df02SPaul Mullowney Level: intermediate 1723e057df02SPaul Mullowney 1724af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 1725e057df02SPaul Mullowney 1726e057df02SPaul Mullowney .seealso: MatCUSPARSESetFormat(), MatCUSPARSEFormatOperation 1727e057df02SPaul Mullowney E*/ 1728e057df02SPaul Mullowney 17293f9c0db1SPaul Mullowney typedef enum {MAT_CUSPARSE_CSR, MAT_CUSPARSE_ELL, MAT_CUSPARSE_HYB} MatCUSPARSEStorageFormat; 1730e057df02SPaul Mullowney 1731e057df02SPaul Mullowney /* these will be strings associated with enumerated type defined above */ 1732e057df02SPaul Mullowney PETSC_EXTERN const char *const MatCUSPARSEStorageFormats[]; 1733e057df02SPaul Mullowney 17343f9c0db1SPaul Mullowney /*E 1735e057df02SPaul Mullowney MatCUSPARSEFormatOperation - indicates the operation of CUSPARSE (GPU) 17362692e278SPaul Mullowney matrices whose operation should use a particular storage format. 1737e057df02SPaul Mullowney 1738e057df02SPaul Mullowney Not Collective 1739e057df02SPaul Mullowney 17408468deeeSKarl Rupp + MAT_CUSPARSE_MULT_DIAG - sets the storage format for the diagonal matrix in the parallel MatMult 17418468deeeSKarl Rupp . MAT_CUSPARSE_MULT_OFFDIAG - sets the storage format for the offdiagonal matrix in the parallel MatMult 17428468deeeSKarl Rupp . MAT_CUSPARSE_MULT - sets the storage format for the entire matrix in the serial (single GPU) MatMult 17438468deeeSKarl Rupp - MAT_CUSPARSE_ALL - sets the storage format for all CUSPARSE (GPU) matrices 1744e057df02SPaul Mullowney 1745e057df02SPaul Mullowney Level: intermediate 1746e057df02SPaul Mullowney 1747e057df02SPaul Mullowney .seealso: MatCUSPARSESetFormat(), MatCUSPARSEStorageFormat 1748e057df02SPaul Mullowney E*/ 174936d62e41SPaul Mullowney typedef enum {MAT_CUSPARSE_MULT_DIAG, MAT_CUSPARSE_MULT_OFFDIAG, MAT_CUSPARSE_MULT, MAT_CUSPARSE_ALL} MatCUSPARSEFormatOperation; 1750e057df02SPaul Mullowney 175110e9db80SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 175210e9db80SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCreateAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 1753e057df02SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCUSPARSESetFormat(Mat,MatCUSPARSEFormatOperation,MatCUSPARSEStorageFormat); 17549ae82921SPaul Mullowney #endif 17559ae82921SPaul Mullowney 175690c53902SBarry Smith #if defined(PETSC_HAVE_CUSP) 1757014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCUSP(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 1758014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateAIJCUSP(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 1759e057df02SPaul Mullowney 17603f9c0db1SPaul Mullowney /*E 1761e057df02SPaul Mullowney MatCUSPStorageFormat - indicates the storage format for CUSP (GPU) 176236d62e41SPaul Mullowney matrices. 1763e057df02SPaul Mullowney 1764e057df02SPaul Mullowney Not Collective 1765e057df02SPaul Mullowney 17668468deeeSKarl Rupp + MAT_CUSP_CSR - Compressed Sparse Row 17678468deeeSKarl Rupp . MAT_CUSP_DIA - Diagonal 17688468deeeSKarl Rupp - MAT_CUSP_ELL - Ellpack 1769e057df02SPaul Mullowney 1770e057df02SPaul Mullowney Level: intermediate 1771e057df02SPaul Mullowney 1772af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 1773e057df02SPaul Mullowney 1774e057df02SPaul Mullowney .seealso: MatCUSPSetFormat(), MatCUSPFormatOperation 1775e057df02SPaul Mullowney E*/ 17763f9c0db1SPaul Mullowney typedef enum {MAT_CUSP_CSR, MAT_CUSP_DIA, MAT_CUSP_ELL} MatCUSPStorageFormat; 1777e057df02SPaul Mullowney 1778e057df02SPaul Mullowney /* these will be strings associated with enumerated type defined above */ 1779e057df02SPaul Mullowney PETSC_EXTERN const char *const MatCUSPStorageFormats[]; 1780e057df02SPaul Mullowney 17813f9c0db1SPaul Mullowney /*E 1782e057df02SPaul Mullowney MatCUSPFormatOperation - indicates the operation of CUSP (GPU) 178336d62e41SPaul Mullowney matrices whose operation should use a particular storage format. 1784e057df02SPaul Mullowney 1785e057df02SPaul Mullowney Not Collective 1786e057df02SPaul Mullowney 17878468deeeSKarl Rupp + MAT_CUSP_MULT_DIAG - sets the storage format for the diagonal matrix in the parallel MatMult 17888468deeeSKarl Rupp . MAT_CUSP_MULT_OFFDIAG - sets the storage format for the offdiagonal matrix in the parallel MatMult 17898468deeeSKarl Rupp . MAT_CUSP_MULT - sets the storage format for the entire matrix in the serial (single GPU) MatMult 17908468deeeSKarl Rupp - MAT_CUSP_ALL - sets the storage format for all CUSP (GPU) matrices 1791e057df02SPaul Mullowney 1792e057df02SPaul Mullowney Level: intermediate 1793e057df02SPaul Mullowney 1794af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 1795e057df02SPaul Mullowney 1796e057df02SPaul Mullowney .seealso: MatCUSPSetFormat(), MatCUSPStorageFormat 1797e057df02SPaul Mullowney E*/ 179836d62e41SPaul Mullowney typedef enum {MAT_CUSP_MULT_DIAG, MAT_CUSP_MULT_OFFDIAG, MAT_CUSP_MULT, MAT_CUSP_ALL} MatCUSPFormatOperation; 1799e057df02SPaul Mullowney 1800e057df02SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCUSPSetFormat(Mat,MatCUSPFormatOperation,MatCUSPStorageFormat); 1801e057df02SPaul Mullowney #endif 180290c53902SBarry Smith 1803d67ff14aSKarl Rupp #if defined(PETSC_HAVE_VIENNACL) 1804d67ff14aSKarl Rupp PETSC_EXTERN PetscErrorCode MatCreateSeqAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 1805d67ff14aSKarl Rupp PETSC_EXTERN PetscErrorCode MatCreateAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 1806d67ff14aSKarl Rupp #endif 1807d67ff14aSKarl Rupp 180854efbe56SHong Zhang /* 180954efbe56SHong Zhang PETSc interface to FFTW 181054efbe56SHong Zhang */ 181154efbe56SHong Zhang #if defined(PETSC_HAVE_FFTW) 1812014dd563SJed Brown PETSC_EXTERN PetscErrorCode VecScatterPetscToFFTW(Mat,Vec,Vec); 1813014dd563SJed Brown PETSC_EXTERN PetscErrorCode VecScatterFFTWToPetsc(Mat,Vec,Vec); 18142a7a6963SBarry Smith PETSC_EXTERN PetscErrorCode MatCreateVecsFFTW(Mat,Vec*,Vec*,Vec*); 181554efbe56SHong Zhang #endif 181654efbe56SHong Zhang 1817014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateNest(MPI_Comm,PetscInt,const IS[],PetscInt,const IS[],const Mat[],Mat*); 1818014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSize(Mat,PetscInt*,PetscInt*); 1819014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetISs(Mat,IS[],IS[]); 1820014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetLocalISs(Mat,IS[],IS[]); 1821014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSubMats(Mat,PetscInt*,PetscInt*,Mat***); 1822014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSubMat(Mat,PetscInt,PetscInt,Mat*); 182319fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatNestSetVecType(Mat,VecType); 1824014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestSetSubMats(Mat,PetscInt,const IS[],PetscInt,const IS[],const Mat[]); 1825014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestSetSubMat(Mat,PetscInt,PetscInt,Mat); 1826d8588912SDave May 18274325cce7SMatthew G Knepley PETSC_EXTERN PetscErrorCode MatChop(Mat,PetscReal); 1828e0a58c10SMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatComputeBandwidth(Mat,PetscReal,PetscInt*); 18294325cce7SMatthew G Knepley 1830b541e6a4SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSubdomainsCreateCoalesce(Mat,PetscInt,PetscInt*,IS**); 183153dd7562SDmitry Karpeev 1832c094ef40SMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatPreallocatorPreallocate(Mat,PetscBool,Mat); 1833c094ef40SMatthew G. Knepley 1834539c167fSBarry Smith PETSC_INTERN PetscErrorCode MatHeaderMerge(Mat,Mat*); 1835539c167fSBarry Smith PETSC_EXTERN PetscErrorCode MatHeaderReplace(Mat,Mat*); 1836539c167fSBarry Smith 183723a5497aSJed Brown #endif 1838