12eac72dbSBarry Smith /* 22eac72dbSBarry Smith Include file for the matrix component of PETSc 32eac72dbSBarry Smith */ 40a835dfdSSatish Balay #ifndef __PETSCMAT_H 50a835dfdSSatish Balay #define __PETSCMAT_H 62c8e378dSBarry Smith #include <petscvec.h> 72eac72dbSBarry Smith 8d9274352SBarry Smith /*S 98f6c3df8SBarry Smith Mat - Abstract PETSc matrix object used to manage all linear operators in PETSc, even those without 108f6c3df8SBarry Smith an explicit sparse representation (such as matrix-free operators) 112eac72dbSBarry Smith 12d91e6319SBarry Smith Level: beginner 13d91e6319SBarry Smith 14d9274352SBarry Smith Concepts: matrix; linear operator 15d9274352SBarry Smith 168f6c3df8SBarry Smith .seealso: MatCreate(), MatType, MatSetType(), MatDestroy() 17d9274352SBarry Smith S*/ 18d9274352SBarry Smith typedef struct _p_Mat* Mat; 19d9274352SBarry Smith 2076bdecfbSBarry Smith /*J 218f6c3df8SBarry Smith MatType - String with the name of a PETSc matrix type 22d9274352SBarry Smith 23d9274352SBarry Smith Level: beginner 24d9274352SBarry Smith 253ca39a21SBarry Smith .seealso: MatSetType(), Mat, MatSolverType, MatRegister() 2676bdecfbSBarry Smith J*/ 2719fd82e9SBarry Smith typedef const char* MatType; 28273d9f13SBarry Smith #define MATSAME "same" 295a11e1b2SBarry Smith #define MATMAIJ "maij" 30273d9f13SBarry Smith #define MATSEQMAIJ "seqmaij" 31273d9f13SBarry Smith #define MATMPIMAIJ "mpimaij" 32273d9f13SBarry Smith #define MATIS "is" 335a11e1b2SBarry Smith #define MATAIJ "aij" 34273d9f13SBarry Smith #define MATSEQAIJ "seqaij" 35273d9f13SBarry Smith #define MATMPIAIJ "mpiaij" 365a11e1b2SBarry Smith #define MATAIJCRL "aijcrl" 375a11e1b2SBarry Smith #define MATSEQAIJCRL "seqaijcrl" 385a11e1b2SBarry Smith #define MATMPIAIJCRL "mpiaijcrl" 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" 90d4002b98SHong Zhang #define MATSELL "sell" 91d4002b98SHong Zhang #define MATSEQSELL "seqsell" 92d4002b98SHong Zhang #define MATMPISELL "mpisell" 93a3b2e22bSHong Zhang #define MATDUMMY "dummy" 94773941d6SBarry Smith 9576bdecfbSBarry Smith /*J 963ca39a21SBarry Smith MatSolverType - 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*/ 104ea799195SBarry Smith typedef const char* MatSolverType; 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 1313ca39a21SBarry Smith .seealso: MatSolverType, 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 136ea799195SBarry Smith PETSC_EXTERN PetscErrorCode MatGetFactor(Mat,MatSolverType,MatFactorType,Mat*); 137ea799195SBarry Smith PETSC_EXTERN PetscErrorCode MatGetFactorAvailable(Mat,MatSolverType,MatFactorType,PetscBool *); 138ea799195SBarry Smith PETSC_EXTERN PetscErrorCode MatFactorGetSolverType(Mat,MatSolverType*); 139014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactorType(Mat,MatFactorType*); 140ea799195SBarry Smith PETSC_EXTERN PetscErrorCode MatSolverTypeRegister(MatSolverType,const MatType,MatFactorType,PetscErrorCode(*)(Mat,MatFactorType,Mat*)); 141ea799195SBarry Smith PETSC_EXTERN PetscErrorCode MatSolverTypeGet(MatSolverType,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 214d4002b98SHong Zhang PETSC_EXTERN PetscErrorCode MatCreateSeqSELL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 215d4002b98SHong Zhang PETSC_EXTERN PetscErrorCode MatCreateSELL(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 216d4002b98SHong Zhang PETSC_EXTERN PetscErrorCode MatSeqSELLSetPreallocation(Mat,PetscInt,const PetscInt[]); 217d4002b98SHong 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 275846b4da1SFande Kong PETSC_EXTERN PetscErrorCode MatResetPreallocation(Mat); 276014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetUp(Mat); 277014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroy(Mat*); 278e56f5c9eSBarry Smith PETSC_EXTERN PetscErrorCode MatGetNonzeroState(Mat,PetscObjectState*); 27921c89e3eSBarry Smith 280014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatConjugate(Mat); 281014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRealPart(Mat); 282014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatImaginaryPart(Mat); 283014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetDiagonalBlock(Mat,Mat*); 284014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetTrace(Mat,PetscScalar*); 285713ccfa9SJed Brown PETSC_EXTERN PetscErrorCode MatInvertBlockDiagonal(Mat,const PetscScalar **); 286479a70c0SJed Brown PETSC_EXTERN PetscErrorCode MatInvertBlockDiagonalMat(Mat,Mat); 28799cafbc1SBarry Smith 2888ed539a5SBarry Smith /* ------------------------------------------------------------*/ 289014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 290014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 291014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesRow(Mat,PetscInt,const PetscScalar[]); 292014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesRowLocal(Mat,PetscInt,const PetscScalar[]); 293014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBatch(Mat,PetscInt,PetscInt,PetscInt[],const PetscScalar[]); 29473a71a0fSBarry Smith PETSC_EXTERN PetscErrorCode MatSetRandom(Mat,PetscRandom); 29584cb2905SBarry Smith 2962ef4de8bSBarry Smith /*S 2972ef4de8bSBarry Smith MatStencil - Data structure (C struct) for storing information about a single row or 29862ef0227SBarry Smith column of a matrix as indexed on an associated grid. These are arguments to MatSetStencil() and MatSetBlockStencil() 29962ef0227SBarry Smith 30062ef0227SBarry 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). 30162ef0227SBarry 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. 30262ef0227SBarry Smith 30362ef0227SBarry Smith For stencil access to vectors see DMDAVecGetArray(), DMDAVecGetArrayF90(). 304be479b30SJed Brown 305be479b30SJed Brown Fortran usage is different, see MatSetValuesStencil() for details. 3062ef4de8bSBarry Smith 3072ef4de8bSBarry Smith Level: beginner 3082ef4de8bSBarry Smith 3092ef4de8bSBarry Smith Concepts: matrix; linear operator 3102ef4de8bSBarry Smith 31162ef0227SBarry Smith .seealso: MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockedStencil(), DMDAVecGetArray(), DMDAVecGetArrayF90() 3122ef4de8bSBarry Smith S*/ 313435da068SBarry Smith typedef struct { 314c1ac3661SBarry Smith PetscInt k,j,i,c; 315435da068SBarry Smith } MatStencil; 3162ef4de8bSBarry Smith 317014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode); 318014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode); 319014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt); 320435da068SBarry Smith 321d91e6319SBarry Smith /*E 322d91e6319SBarry Smith MatAssemblyType - Indicates if the matrix is now to be used, or if you plan 32362ef0227SBarry Smith to continue to add or insert values to it 324d91e6319SBarry Smith 325d91e6319SBarry Smith Level: beginner 326d91e6319SBarry Smith 327d91e6319SBarry Smith .seealso: MatAssemblyBegin(), MatAssemblyEnd() 328d91e6319SBarry Smith E*/ 3296d4a8577SBarry Smith typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType; 330014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssemblyBegin(Mat,MatAssemblyType); 331014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssemblyEnd(Mat,MatAssemblyType); 332014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssembled(Mat,PetscBool *); 3334f9c727eSBarry Smith 3341d73ed98SBarry Smith 33530de9b25SBarry Smith 336d91e6319SBarry Smith /*E 337d91e6319SBarry Smith MatOption - Options that may be set for a matrix and its behavior or storage 338d91e6319SBarry Smith 339d91e6319SBarry Smith Level: beginner 340d91e6319SBarry Smith 341af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 3420f8fb01aSBarry Smith Any additions/changes here must also be made in src/mat/interface/dlregismat.c in MatOptions[] 343d91e6319SBarry Smith 344382164b0SBarry Smith Developer Notes: Entries that are negative need not be called collectively by all processes. 345382164b0SBarry Smith 346d91e6319SBarry Smith .seealso: MatSetOption() 347d91e6319SBarry Smith E*/ 348c5e4d11fSDmitry Karpeev typedef enum {MAT_OPTION_MIN = -3, 349c5e4d11fSDmitry Karpeev MAT_UNUSED_NONZERO_LOCATION_ERR = -2, 35092d4d306SBarry Smith MAT_ROW_ORIENTED = -1, 351382164b0SBarry Smith MAT_SYMMETRIC = 1, 352382164b0SBarry Smith MAT_STRUCTURALLY_SYMMETRIC = 2, 353382164b0SBarry Smith MAT_NEW_DIAGONALS = 3, 354382164b0SBarry Smith MAT_IGNORE_OFF_PROC_ENTRIES = 4, 355382164b0SBarry Smith MAT_USE_HASH_TABLE = 5, 356382164b0SBarry Smith MAT_KEEP_NONZERO_PATTERN = 6, 357382164b0SBarry Smith MAT_IGNORE_ZERO_ENTRIES = 7, 358382164b0SBarry Smith MAT_USE_INODES = 8, 359382164b0SBarry Smith MAT_HERMITIAN = 9, 360382164b0SBarry Smith MAT_SYMMETRY_ETERNAL = 10, 361c5e4d11fSDmitry Karpeev MAT_NEW_NONZERO_LOCATION_ERR = 11, 362382164b0SBarry Smith MAT_IGNORE_LOWER_TRIANGULAR = 12, 363382164b0SBarry Smith MAT_ERROR_LOWER_TRIANGULAR = 13, 364382164b0SBarry Smith MAT_GETROW_UPPERTRIANGULAR = 14, 365382164b0SBarry Smith MAT_SPD = 15, 36692d4d306SBarry Smith MAT_NO_OFF_PROC_ZERO_ROWS = 16, 36792d4d306SBarry Smith MAT_NO_OFF_PROC_ENTRIES = 17, 36892d4d306SBarry Smith MAT_NEW_NONZERO_LOCATIONS = 18, 369c5e4d11fSDmitry Karpeev MAT_NEW_NONZERO_ALLOCATION_ERR = 19, 370c5e4d11fSDmitry Karpeev MAT_SUBSET_OFF_PROC_ENTRIES = 20, 371c10200c1SHong Zhang MAT_SUBMAT_SINGLEIS = 21, 372957cac9fSHong Zhang MAT_STRUCTURE_ONLY = 22, 373957cac9fSHong Zhang MAT_OPTION_MAX = 23} MatOption; 374e057df02SPaul Mullowney 3750f8fb01aSBarry Smith PETSC_EXTERN const char *const *MatOptions; 376014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetOption(Mat,MatOption,PetscBool); 3777d68702bSBarry Smith PETSC_EXTERN PetscErrorCode MatGetOption(Mat,MatOption,PetscBool*); 37819fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatGetType(Mat,MatType*); 37984cb2905SBarry Smith 380014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]); 381014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 382014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 383014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowUpperTriangular(Mat); 384014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRowUpperTriangular(Mat); 385014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnVector(Mat,Vec,PetscInt); 3868c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJGetArray(Mat,PetscScalar *[]); 3878c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJRestoreArray(Mat,PetscScalar *[]); 38821e72a00SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJGetMaxRowNonzeros(Mat,PetscInt*); 389bd04181cSBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJSetValuesLocalFast(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 3904099cc6bSBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJSetType(Mat,MatType); 391388d47a6SSatish Balay PETSC_EXTERN PetscErrorCode MatSeqAIJRegister(const char[],PetscErrorCode (*)(Mat,MatType,MatReuse,Mat *)); 3924099cc6bSBarry Smith PETSC_EXTERN PetscFunctionList MatSeqAIJList; 3938397e458SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqSBAIJGetArray(Mat,PetscScalar *[]); 3948397e458SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqSBAIJRestoreArray(Mat,PetscScalar *[]); 3958c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatDenseGetArray(Mat,PetscScalar *[]); 3968c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatDenseRestoreArray(Mat,PetscScalar *[]); 397d3042a70SBarry Smith PETSC_EXTERN PetscErrorCode MatDensePlaceArray(Mat,const PetscScalar[]); 398d3042a70SBarry Smith PETSC_EXTERN PetscErrorCode MatDenseResetArray(Mat); 3998572280aSBarry Smith PETSC_EXTERN PetscErrorCode MatDenseGetArrayRead(Mat,const PetscScalar *[]); 4008572280aSBarry Smith PETSC_EXTERN PetscErrorCode MatDenseRestoreArrayRead(Mat,const PetscScalar *[]); 401014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBlockSize(Mat,PetscInt *); 402014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSize(Mat,PetscInt); 403014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBlockSizes(Mat,PetscInt *,PetscInt *); 404014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSizes(Mat,PetscInt,PetscInt); 40533d57670SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSizesFromMats(Mat,Mat,Mat); 4061620fd73SBarry Smith 407014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMult(Mat,Vec,Vec); 408014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultDiagonalBlock(Mat,Vec,Vec); 409014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultAdd(Mat,Vec,Vec,Vec); 410014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTranspose(Mat,Vec,Vec); 411014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultHermitianTranspose(Mat,Vec,Vec); 412014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsTranspose(Mat,Mat,PetscReal,PetscBool *); 413014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitianTranspose(Mat,Mat,PetscReal,PetscBool *); 414014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeAdd(Mat,Vec,Vec,Vec); 415014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultHermitianTransposeAdd(Mat,Vec,Vec,Vec); 416014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultConstrained(Mat,Vec,Vec); 417014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeConstrained(Mat,Vec,Vec); 418014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatSolve(Mat,Mat,Mat); 419bdc285e1SStefano Zampini PETSC_EXTERN PetscErrorCode MatMatSolveTranspose(Mat,Mat,Mat); 420f9426fe0SMark Adams PETSC_EXTERN PetscErrorCode MatResidual(Mat,Vec,Vec,Vec); 421f5edf698SHong Zhang 422d91e6319SBarry Smith /*E 423d91e6319SBarry Smith MatDuplicateOption - Indicates if a duplicated sparse matrix should have 424d91e6319SBarry Smith its numerical values copied over or just its nonzero structure. 425d91e6319SBarry Smith 426d91e6319SBarry Smith Level: beginner 427d91e6319SBarry Smith 428af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 429d91e6319SBarry Smith 4300018da39SRichard Tran Mills $ MAT_DO_NOT_COPY_VALUES - Create a matrix using the same nonzero pattern as the original matrix, 4310018da39SRichard Tran Mills $ with zeros for the numerical values. 4320018da39SRichard Tran Mills $ MAT_COPY_VALUES - Create a matrix with the same nonzero pattern as the original matrix 4330018da39SRichard Tran Mills $ and with the same numerical values. 4340018da39SRichard Tran Mills $ MAT_SHARE_NONZERO_PATTERN - Create a matrix that shares the nonzero structure with the previous matrix 4350018da39SRichard Tran Mills $ and does not copy it, using zeros for the numerical values. The parent and 4360018da39SRichard Tran Mills $ child matrices will share their index (i and j) arrays, and you cannot 4370018da39SRichard Tran Mills $ insert new nonzero entries into either matrix. 4380018da39SRichard Tran Mills 4390018da39SRichard Tran Mills Notes: Many matrix types (including SeqAIJ) do not support the MAT_SHARE_NONZERO_PATTERN optimization; in 4400018da39SRichard Tran Mills this case the behavior is as if MAT_DO_NOT_COPY_VALUES has been specified. 44170dcbbb9SBarry Smith 442d91e6319SBarry Smith .seealso: MatDuplicate() 443d91e6319SBarry Smith E*/ 44470dcbbb9SBarry Smith typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES,MAT_SHARE_NONZERO_PATTERN} MatDuplicateOption; 4452e8a6d31SBarry Smith 44619fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatConvert(Mat,MatType,MatReuse,Mat*); 447014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDuplicate(Mat,MatDuplicateOption,Mat*); 44894a9d846SBarry Smith 449cb5b572fSBarry Smith 450014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCopy(Mat,Mat,MatStructure); 451014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatView(Mat,PetscViewer); 452014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsSymmetric(Mat,PetscReal,PetscBool *); 453014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsStructurallySymmetric(Mat,PetscBool *); 454014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitian(Mat,PetscReal,PetscBool *); 455014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsSymmetricKnown(Mat,PetscBool *,PetscBool *); 456014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitianKnown(Mat,PetscBool *,PetscBool *); 457014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMissingDiagonal(Mat,PetscBool *,PetscInt *); 458014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLoad(Mat, PetscViewer); 4597b80b807SBarry Smith 4601a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool *); 4611a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool *); 4621a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool *); 4631a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool *); 464d4fbbf0eSBarry Smith 465d91e6319SBarry Smith /*S 466d91e6319SBarry Smith MatInfo - Context of matrix information, used with MatGetInfo() 467d91e6319SBarry Smith 468d91e6319SBarry Smith In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE 469d91e6319SBarry Smith 470d91e6319SBarry Smith Level: intermediate 471d91e6319SBarry Smith 472d91e6319SBarry Smith Concepts: matrix^nonzero information 473d91e6319SBarry Smith 474d9274352SBarry Smith .seealso: MatGetInfo(), MatInfoType 475d91e6319SBarry Smith S*/ 4764e220ebcSLois Curfman McInnes typedef struct { 477b0a32e0cSBarry Smith PetscLogDouble block_size; /* block size */ 478b0a32e0cSBarry Smith PetscLogDouble nz_allocated,nz_used,nz_unneeded; /* number of nonzeros */ 479b0a32e0cSBarry Smith PetscLogDouble memory; /* memory allocated */ 480b0a32e0cSBarry Smith PetscLogDouble assemblies; /* number of matrix assemblies called */ 481b0a32e0cSBarry Smith PetscLogDouble mallocs; /* number of mallocs during MatSetValues() */ 482b0a32e0cSBarry Smith PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */ 483b0a32e0cSBarry Smith PetscLogDouble factor_mallocs; /* number of mallocs during factorization */ 4844e220ebcSLois Curfman McInnes } MatInfo; 4854e220ebcSLois Curfman McInnes 486d9274352SBarry Smith /*E 487d9274352SBarry Smith MatInfoType - Indicates if you want information about the local part of the matrix, 488d9274352SBarry Smith the entire parallel matrix or the maximum over all the local parts. 489d9274352SBarry Smith 490d9274352SBarry Smith Level: beginner 491d9274352SBarry Smith 492af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 493d9274352SBarry Smith 494d9274352SBarry Smith .seealso: MatGetInfo(), MatInfo 495d9274352SBarry Smith E*/ 4967b80b807SBarry Smith typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType; 497014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetInfo(Mat,MatInfoType,MatInfo*); 498014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetDiagonal(Mat,Vec); 499014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMax(Mat,Vec,PetscInt[]); 500014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMin(Mat,Vec,PetscInt[]); 501014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMaxAbs(Mat,Vec,PetscInt[]); 502014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMinAbs(Mat,Vec,PetscInt[]); 503014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowSum(Mat,Vec); 504014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTranspose(Mat,MatReuse,Mat*); 505014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatHermitianTranspose(Mat,MatReuse,Mat*); 506014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPermute(Mat,IS,IS,Mat*); 507014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalScale(Mat,Vec,Vec); 508014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalSet(Mat,Vec,InsertMode); 509a52676ecSHong Zhang 510014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatEqual(Mat,Mat,PetscBool*); 511014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultEqual(Mat,Mat,PetscInt,PetscBool*); 512014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultAddEqual(Mat,Mat,PetscInt,PetscBool*); 513014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeEqual(Mat,Mat,PetscInt,PetscBool*); 514014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscBool*); 515a52676ecSHong Zhang PETSC_EXTERN PetscErrorCode MatMatMultEqual(Mat,Mat,Mat,PetscInt,PetscBool*); 516a52676ecSHong Zhang PETSC_EXTERN PetscErrorCode MatTransposeMatMultEqual(Mat,Mat,Mat,PetscInt,PetscBool*); 5177b80b807SBarry Smith 518014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNorm(Mat,NormType,PetscReal*); 519014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnNorms(Mat,NormType,PetscReal*); 520014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroEntries(Mat); 521014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 522014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsIS(Mat,IS,PetscScalar,Vec,Vec); 523014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsStencil(Mat,PetscInt,const MatStencil [],PetscScalar,Vec,Vec); 524014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsStencil(Mat,PetscInt,const MatStencil[],PetscScalar,Vec,Vec); 525014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumns(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 526014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsIS(Mat,IS,PetscScalar,Vec,Vec); 5277b80b807SBarry Smith 528014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSize(Mat,PetscInt*,PetscInt*); 529014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalSize(Mat,PetscInt*,PetscInt*); 530014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRange(Mat,PetscInt*,PetscInt*); 531014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRanges(Mat,const PetscInt**); 532014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRangeColumn(Mat,PetscInt*,PetscInt*); 533014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRangesColumn(Mat,const PetscInt**); 534566876eaSJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipIS(Mat,IS*,IS*); 5355ef9f2a5SBarry Smith 5367dae84e0SHong Zhang PETSC_EXTERN PetscErrorCode MatCreateSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]); 537cd2534ddSHong 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);} 5387dae84e0SHong Zhang PETSC_EXTERN PetscErrorCode MatCreateSubMatricesMPI(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]); 539cd2534ddSHong 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);} 540014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroyMatrices(PetscInt,Mat *[]); 541df750dc8SHong Zhang PETSC_EXTERN PetscErrorCode MatDestroySubMatrices(PetscInt,Mat *[]); 5427dae84e0SHong Zhang PETSC_EXTERN PetscErrorCode MatCreateSubMatrix(Mat,IS,IS,MatReuse,Mat *); 543cd2534ddSHong 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);} 544014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalSubMatrix(Mat,IS,IS,Mat*); 545014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreLocalSubMatrix(Mat,IS,IS,Mat*); 546014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSeqNonzeroStructure(Mat,Mat*); 547014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroySeqNonzeroStructure(Mat*); 5487b80b807SBarry Smith 549014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJ(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*); 550014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJSymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*); 551014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJNumeric(Mat,Mat); 552014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMat(Mat,MatReuse,Mat*); 553014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*); 554014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,Mat*); 555014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetGhosts(Mat, PetscInt *,const PetscInt *[]); 5568efafbd8SBarry Smith 557014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt); 558aa1e27eaSFande Kong PETSC_EXTERN PetscErrorCode MatIncreaseOverlapSplit(Mat mat,PetscInt n,IS is[],PetscInt ov); 559d2b2a242SBarry Smith PETSC_EXTERN PetscErrorCode MatMPIAIJSetUseScalableIncreaseOverlap(Mat,PetscBool); 5607b80b807SBarry Smith 561014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*); 562014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMultSymbolic(Mat,Mat,PetscReal,Mat*); 563014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMultNumeric(Mat,Mat,Mat); 56422440eb1SKris Buschelman 5657bab7c10SHong Zhang PETSC_EXTERN PetscErrorCode MatMatMatMult(Mat,Mat,Mat,MatReuse,PetscReal,Mat*); 5662df6c5c3SPatrick Farrell PETSC_EXTERN PetscErrorCode MatGalerkin(Mat,Mat,Mat,MatReuse,PetscReal,Mat*); 5677bab7c10SHong Zhang 568014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*); 569014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAPSymbolic(Mat,Mat,PetscReal,Mat*); 570014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAPNumeric(Mat,Mat,Mat); 571014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARt(Mat,Mat,MatReuse,PetscReal,Mat*); 572014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARtSymbolic(Mat,Mat,PetscReal,Mat*); 573014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARtNumeric(Mat,Mat,Mat); 57422440eb1SKris Buschelman 575014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeMatMult(Mat,Mat,MatReuse,PetscReal,Mat*); 576014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposetMatMultSymbolic(Mat,Mat,PetscReal,Mat*); 577014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposetMatMultNumeric(Mat,Mat,Mat); 578014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMult(Mat,Mat,MatReuse,PetscReal,Mat*); 579014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMultSymbolic(Mat,Mat,PetscReal,Mat*); 580014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMultNumeric(Mat,Mat,Mat); 581bc011b1eSHong Zhang 582014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAXPY(Mat,PetscScalar,Mat,MatStructure); 583014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAYPX(Mat,PetscScalar,Mat,MatStructure); 5841c741599SHong Zhang 585014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScale(Mat,PetscScalar); 586014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShift(Mat,PetscScalar); 5877b80b807SBarry Smith 588014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping); 589014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*); 590a93ff8c4SPeter Brune PETSC_EXTERN PetscErrorCode MatGetLayouts(Mat,PetscLayout*,PetscLayout*); 591014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 592014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsLocalIS(Mat,IS,PetscScalar,Vec,Vec); 593014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 594014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocalIS(Mat,IS,PetscScalar,Vec,Vec); 595014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 596014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 597052efed2SBarry Smith 598014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStashSetInitialSize(Mat,PetscInt,PetscInt); 599014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*); 60090f02eecSBarry Smith 601014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInterpolate(Mat,Vec,Vec); 602014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInterpolateAdd(Mat,Vec,Vec,Vec); 603014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestrict(Mat,Vec,Vec); 6042a7a6963SBarry Smith PETSC_EXTERN PetscErrorCode MatCreateVecs(Mat,Vec*,Vec*); 60584cc2a5aSBarry Smith PETSC_DEPRECATED("Use MatCreateVecs()") PETSC_STATIC_INLINE PetscErrorCode MatGetVecs(Mat mat,Vec *x,Vec *y) {return MatCreateVecs(mat,x,y);} 60653cd1579SHong Zhang PETSC_EXTERN PetscErrorCode MatCreateRedundantMatrix(Mat,PetscInt,MPI_Comm,MatReuse,Mat*); 607014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetMultiProcBlock(Mat,MPI_Comm,MatReuse,Mat*); 608014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFindZeroDiagonals(Mat,IS*); 6093a062f41SBarry Smith PETSC_EXTERN PetscErrorCode MatFindOffBlockDiagonalEntries(Mat,IS*); 6109c8f2541SHong Zhang PETSC_EXTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat(MPI_Comm,Mat,PetscInt,MatReuse,Mat*); 611bd481603SBarry Smith 612bd481603SBarry Smith /*MC 6131620fd73SBarry Smith MatSetValue - Set a single entry into a matrix. 6141620fd73SBarry Smith 6151620fd73SBarry Smith Not collective 6161620fd73SBarry Smith 617a9834a7dSSatish Balay Synopsis: 618a9834a7dSSatish Balay #include <petscmat.h> 619a9834a7dSSatish Balay PetscErrorCode MatSetValue(Mat m,PetscInt row,PetscInt col,PetscScalar value,InsertMode mode) 620a9834a7dSSatish Balay 6211620fd73SBarry Smith Input Parameters: 6221620fd73SBarry Smith + m - the matrix 6231620fd73SBarry Smith . row - the row location of the entry 6241620fd73SBarry Smith . col - the column location of the entry 6251620fd73SBarry Smith . value - the value to insert 6261620fd73SBarry Smith - mode - either INSERT_VALUES or ADD_VALUES 6271620fd73SBarry Smith 6281620fd73SBarry Smith Notes: 6291620fd73SBarry Smith For efficiency one should use MatSetValues() and set several or many 6301620fd73SBarry Smith values simultaneously if possible. 6311620fd73SBarry Smith 6321620fd73SBarry Smith Level: beginner 6331620fd73SBarry Smith 6341620fd73SBarry Smith .seealso: MatSetValues(), MatSetValueLocal() 6351620fd73SBarry Smith M*/ 6361620fd73SBarry 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);} 6371620fd73SBarry Smith 6381620fd73SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatGetValue(Mat v,PetscInt i,PetscInt j,PetscScalar *va) {return MatGetValues(v,1,&i,1,&j,va);} 6391620fd73SBarry Smith 6401620fd73SBarry 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);} 6411620fd73SBarry Smith 6421620fd73SBarry Smith /*MC 6430d2a7ff2SSatish Balay MatPreallocateInitialize - Begins the block of code that will count the number of nonzeros per 644bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 645bd481603SBarry Smith 646bd481603SBarry Smith Synopsis: 647aaa7dc30SBarry Smith #include <petscmat.h> 648c1ac3661SBarry Smith PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz) 649bd481603SBarry Smith 650bd481603SBarry Smith Collective on MPI_Comm 651bd481603SBarry Smith 652bd481603SBarry Smith Input Parameters: 653bd481603SBarry Smith + comm - the communicator that will share the eventually allocated matrix 654859fcb39SBarry Smith . nrows - the number of LOCAL rows in the matrix 655859fcb39SBarry Smith - ncols - the number of LOCAL columns in the matrix 656bd481603SBarry Smith 657bd481603SBarry Smith Output Parameters: 658bd481603SBarry Smith + dnz - the array that will be passed to the matrix preallocation routines 659bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 660bd481603SBarry Smith 661bd481603SBarry Smith Level: intermediate 662bd481603SBarry Smith 663bd481603SBarry Smith Notes: 664a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 665bd481603SBarry Smith 6661d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 667bd481603SBarry Smith 668bd481603SBarry Smith Use MatPreallocateInitializeSymmetric() for symmetric matrices (MPISBAIJ matrices) 669bd481603SBarry Smith 6701620fd73SBarry Smith This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize(). 6711620fd73SBarry Smith 672bd481603SBarry Smith Concepts: preallocation^Matrix 673bd481603SBarry Smith 674d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(), 675d6e23781SBarry Smith MatPreallocateInitializeSymmetric(), MatPreallocateSymmetricSetLocalBlock() 676bd481603SBarry Smith M*/ 677c4f061fbSSatish Balay #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \ 6787c922b88SBarry Smith { \ 679a89bc211SBarry Smith PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__start,__end; \ 6807282cfc1SLisandro Dalcin _4_ierr = PetscCalloc2((size_t)__nrows,&dnz,(size_t)__nrows,&onz);CHKERRQ(_4_ierr); \ 6811795a4d1SJed Brown __start = 0; __end = __start; \ 682c1ac3661SBarry Smith _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ctmp;\ 683a89bc211SBarry Smith _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows; 6847c922b88SBarry Smith 685bd481603SBarry Smith /*MC 686bd481603SBarry Smith MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be 687bd481603SBarry Smith inserted using a local number of the rows and columns 688bd481603SBarry Smith 689bd481603SBarry Smith Synopsis: 690aaa7dc30SBarry Smith #include <petscmat.h> 691c1ac3661SBarry Smith PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 692bd481603SBarry Smith 693bd481603SBarry Smith Not Collective 694bd481603SBarry Smith 695bd481603SBarry Smith Input Parameters: 696784ac674SJed Brown + map - the row mapping from local numbering to global numbering 697bd481603SBarry Smith . nrows - the number of rows indicated 6981d73ed98SBarry Smith . rows - the indices of the rows 699784ac674SJed Brown . cmap - the column mapping from local to global numbering 700bd481603SBarry Smith . ncols - the number of columns in the matrix 701bd481603SBarry Smith . cols - the columns indicated 702bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 703bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 704bd481603SBarry Smith 705bd481603SBarry Smith Level: intermediate 706bd481603SBarry Smith 707bd481603SBarry Smith Notes: 708a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 709bd481603SBarry Smith 7101d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 711bd481603SBarry Smith 712bd481603SBarry Smith Concepts: preallocation^Matrix 713bd481603SBarry Smith 714d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(), 715c1154cd5SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocalRemoveDups() 716bd481603SBarry Smith M*/ 717784ac674SJed Brown #define MatPreallocateSetLocal(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \ 718c4f061fbSSatish Balay {\ 719c1ac3661SBarry Smith PetscInt __l;\ 720784ac674SJed Brown _4_ierr = ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\ 721784ac674SJed Brown _4_ierr = ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\ 722c4f061fbSSatish Balay for (__l=0;__l<nrows;__l++) {\ 723ef66eb69SBarry Smith _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 724c4f061fbSSatish Balay }\ 725c4f061fbSSatish Balay } 726c4f061fbSSatish Balay 727bd481603SBarry Smith /*MC 728c1154cd5SBarry Smith MatPreallocateSetLocalRemoveDups - Indicates the locations (rows and columns) in the matrix where nonzeros will be 729c1154cd5SBarry Smith inserted using a local number of the rows and columns. This version removes any duplicate columns in cols 730c1154cd5SBarry Smith 731c1154cd5SBarry Smith Synopsis: 732c1154cd5SBarry Smith #include <petscmat.h> 733c1154cd5SBarry Smith PetscErrorCode MatPreallocateSetLocalRemoveDups(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 734c1154cd5SBarry Smith 735c1154cd5SBarry Smith Not Collective 736c1154cd5SBarry Smith 737c1154cd5SBarry Smith Input Parameters: 738c1154cd5SBarry Smith + map - the row mapping from local numbering to global numbering 739c1154cd5SBarry Smith . nrows - the number of rows indicated 740c1154cd5SBarry Smith . rows - the indices of the rows (these values are mapped to the global values) 741c1154cd5SBarry Smith . cmap - the column mapping from local to global numbering 742c1154cd5SBarry Smith . ncols - the number of columns in the matrix (this value will be changed if duplicate columns are found) 743c1154cd5SBarry Smith . cols - the columns indicated (these values are mapped to the global values, they are then sorted and duplicates removed) 744c1154cd5SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 745c1154cd5SBarry Smith - ozn - the other array passed to the matrix preallocation routines 746c1154cd5SBarry Smith 747c1154cd5SBarry Smith Level: intermediate 748c1154cd5SBarry Smith 749c1154cd5SBarry Smith Notes: 750c1154cd5SBarry Smith See Users-Manual: ch_performance for more details. 751c1154cd5SBarry Smith 752c1154cd5SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 753c1154cd5SBarry Smith 754c1154cd5SBarry Smith Concepts: preallocation^Matrix 755c1154cd5SBarry Smith 756c1154cd5SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(), 757c1154cd5SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocal() 758c1154cd5SBarry Smith M*/ 759c1154cd5SBarry Smith #define MatPreallocateSetLocalRemoveDups(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \ 760c1154cd5SBarry Smith {\ 761c1154cd5SBarry Smith PetscInt __l;\ 762c1154cd5SBarry Smith _4_ierr = ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\ 763c1154cd5SBarry Smith _4_ierr = ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\ 764c1154cd5SBarry Smith _4_ierr = PetscSortRemoveDupsInt(&ncols,cols);CHKERRQ(_4_ierr);\ 765c1154cd5SBarry Smith for (__l=0;__l<nrows;__l++) {\ 766c1154cd5SBarry Smith _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 767c1154cd5SBarry Smith }\ 768c1154cd5SBarry Smith } 769c1154cd5SBarry Smith 770c1154cd5SBarry Smith /*MC 771d6e23781SBarry Smith MatPreallocateSetLocalBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be 772bd481603SBarry Smith inserted using a local number of the rows and columns 773bd481603SBarry Smith 774bd481603SBarry Smith Synopsis: 775aaa7dc30SBarry Smith #include <petscmat.h> 776d6e23781SBarry Smith PetscErrorCode MatPreallocateSetLocalBlock(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 777d6e23781SBarry Smith 778d6e23781SBarry Smith Not Collective 779d6e23781SBarry Smith 780d6e23781SBarry Smith Input Parameters: 781d6e23781SBarry Smith + map - the row mapping from local numbering to global numbering 782d6e23781SBarry Smith . nrows - the number of rows indicated 783d6e23781SBarry Smith . rows - the indices of the rows 784d6e23781SBarry Smith . cmap - the column mapping from local to global numbering 785d6e23781SBarry Smith . ncols - the number of columns in the matrix 786d6e23781SBarry Smith . cols - the columns indicated 787d6e23781SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 788d6e23781SBarry Smith - ozn - the other array passed to the matrix preallocation routines 789d6e23781SBarry Smith 790d6e23781SBarry Smith Level: intermediate 791d6e23781SBarry Smith 792d6e23781SBarry Smith Notes: 793d6e23781SBarry Smith See Users-Manual: ch_performance for more details. 794d6e23781SBarry Smith 795d6e23781SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 796d6e23781SBarry Smith 797d6e23781SBarry Smith Concepts: preallocation^Matrix 798d6e23781SBarry Smith 799d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(), 800d6e23781SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock() 801d6e23781SBarry Smith M*/ 802d6e23781SBarry Smith #define MatPreallocateSetLocalBlock(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \ 803d6e23781SBarry Smith {\ 804d6e23781SBarry Smith PetscInt __l;\ 805d6e23781SBarry Smith _4_ierr = ISLocalToGlobalMappingApplyBlock(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\ 806d6e23781SBarry Smith _4_ierr = ISLocalToGlobalMappingApplyBlock(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\ 807d6e23781SBarry Smith for (__l=0;__l<nrows;__l++) {\ 808d6e23781SBarry Smith _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 809d6e23781SBarry Smith }\ 810d6e23781SBarry Smith } 811d6e23781SBarry Smith 812d6e23781SBarry Smith /*MC 813d6e23781SBarry Smith MatPreallocateSymmetricSetLocalBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be 814d6e23781SBarry Smith inserted using a local number of the rows and columns 815d6e23781SBarry Smith 816d6e23781SBarry Smith Synopsis: 817d6e23781SBarry Smith #include <petscmat.h> 818d6e23781SBarry Smith PetscErrorCode MatPreallocateSymmetricSetLocalBlock(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 819bd481603SBarry Smith 820bd481603SBarry Smith Not Collective 821bd481603SBarry Smith 822bd481603SBarry Smith Input Parameters: 823bd481603SBarry Smith + map - the mapping between local numbering and global numbering 824bd481603SBarry Smith . nrows - the number of rows indicated 8251d73ed98SBarry Smith . rows - the indices of the rows 826bd481603SBarry Smith . ncols - the number of columns in the matrix 827bd481603SBarry Smith . cols - the columns indicated 828bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 829bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 830bd481603SBarry Smith 831bd481603SBarry Smith Level: intermediate 832bd481603SBarry Smith 833bd481603SBarry Smith Notes: 834a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 835bd481603SBarry Smith 836bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 837bd481603SBarry Smith 838bd481603SBarry Smith Concepts: preallocation^Matrix 839bd481603SBarry Smith 840d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateInitialize(), 841d6e23781SBarry Smith MatPreallocateInitialize(), MatPreallocateSetLocal() 842bd481603SBarry Smith M*/ 843d6e23781SBarry Smith #define MatPreallocateSymmetricSetLocalBlock(map,nrows,rows,ncols,cols,dnz,onz) 0;\ 844d3d32019SBarry Smith {\ 845c1ac3661SBarry Smith PetscInt __l;\ 846d6e23781SBarry Smith _4_ierr = ISLocalToGlobalMappingApplyBlock(map,nrows,rows,rows);CHKERRQ(_4_ierr);\ 847d6e23781SBarry Smith _4_ierr = ISLocalToGlobalMappingApplyBlock(map,ncols,cols,cols);CHKERRQ(_4_ierr);\ 848d3d32019SBarry Smith for (__l=0;__l<nrows;__l++) {\ 849d6e23781SBarry Smith _4_ierr = MatPreallocateSymmetricSetBlock((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 850d3d32019SBarry Smith }\ 851d3d32019SBarry Smith } 852bd481603SBarry Smith /*MC 853bd481603SBarry Smith MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be 854bd481603SBarry Smith inserted using a local number of the rows and columns 855bd481603SBarry Smith 856bd481603SBarry Smith Synopsis: 857aaa7dc30SBarry Smith #include <petscmat.h> 858c1ac3661SBarry Smith PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 859bd481603SBarry Smith 860bd481603SBarry Smith Not Collective 861bd481603SBarry Smith 862bd481603SBarry Smith Input Parameters: 86364075487SBarry Smith + row - the row 864bd481603SBarry Smith . ncols - the number of columns in the matrix 865a50f8bf6SHong Zhang - cols - the columns indicated 866a50f8bf6SHong Zhang 867a50f8bf6SHong Zhang Output Parameters: 868a50f8bf6SHong Zhang + dnz - the array that will be passed to the matrix preallocation routines 869bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 870bd481603SBarry Smith 871bd481603SBarry Smith Level: intermediate 872bd481603SBarry Smith 873bd481603SBarry Smith Notes: 874a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 875bd481603SBarry Smith 876bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 877bd481603SBarry Smith 8781620fd73SBarry Smith This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize(). 8791620fd73SBarry Smith 880bd481603SBarry Smith Concepts: preallocation^Matrix 881bd481603SBarry Smith 882d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(), 883d6e23781SBarry Smith MatPreallocateInitialize(), MatPreallocateSetLocal() 884bd481603SBarry Smith M*/ 885c4f061fbSSatish Balay #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\ 886c1ac3661SBarry Smith { PetscInt __i; \ 887e32f2f54SBarry 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);\ 888e32f2f54SBarry 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);\ 8897c922b88SBarry Smith for (__i=0; __i<nc; __i++) {\ 89064075487SBarry Smith if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \ 8917cc688d7SBarry Smith else dnz[row - __rstart]++;\ 8927c922b88SBarry Smith }\ 8937c922b88SBarry Smith } 8947c922b88SBarry Smith 895bd481603SBarry Smith /*MC 896d6e23781SBarry Smith MatPreallocateSymmetricSetBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be 897bd481603SBarry Smith inserted using a local number of the rows and columns 898bd481603SBarry Smith 899bd481603SBarry Smith Synopsis: 900aaa7dc30SBarry Smith #include <petscmat.h> 901d6e23781SBarry Smith PetscErrorCode MatPreallocateSymmetricSetBlock(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 902bd481603SBarry Smith 903bd481603SBarry Smith Not Collective 904bd481603SBarry Smith 905bd481603SBarry Smith Input Parameters: 906bd481603SBarry Smith + nrows - the number of rows indicated 9071d73ed98SBarry Smith . rows - the indices of the rows 908bd481603SBarry Smith . ncols - the number of columns in the matrix 909bd481603SBarry Smith . cols - the columns indicated 910bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 911bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 912bd481603SBarry Smith 913bd481603SBarry Smith Level: intermediate 914bd481603SBarry Smith 915bd481603SBarry Smith Notes: 916a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 917bd481603SBarry Smith 918bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 919bd481603SBarry Smith 9201620fd73SBarry Smith This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize(). 9211620fd73SBarry Smith 922bd481603SBarry Smith Concepts: preallocation^Matrix 923bd481603SBarry Smith 924d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateInitialize(), 925d6e23781SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocal() 926bd481603SBarry Smith M*/ 927d6e23781SBarry Smith #define MatPreallocateSymmetricSetBlock(row,nc,cols,dnz,onz) 0;\ 928c1ac3661SBarry Smith { PetscInt __i; \ 929d3d32019SBarry Smith for (__i=0; __i<nc; __i++) {\ 930d3d32019SBarry Smith if (cols[__i] >= __end) onz[row - __rstart]++; \ 931d3d32019SBarry Smith else if (cols[__i] >= row) dnz[row - __rstart]++;\ 932d3d32019SBarry Smith }\ 933d3d32019SBarry Smith } 934d3d32019SBarry Smith 935bd481603SBarry Smith /*MC 93616371a99SBarry Smith MatPreallocateLocation - An alternative to MatPreallocationSet() that puts the nonzero locations into the matrix if it exists 93716371a99SBarry Smith 93816371a99SBarry Smith Synopsis: 939aaa7dc30SBarry Smith #include <petscmat.h> 94016371a99SBarry Smith PetscErrorCode MatPreallocateLocations(Mat A,PetscInt row,PetscInt ncols,PetscInt *cols,PetscInt *dnz,PetscInt *onz) 94116371a99SBarry Smith 94216371a99SBarry Smith Not Collective 94316371a99SBarry Smith 94416371a99SBarry Smith Input Parameters: 94516371a99SBarry Smith . A - matrix 94616371a99SBarry Smith . row - row where values exist (must be local to this process) 94716371a99SBarry Smith . ncols - number of columns 94816371a99SBarry Smith . cols - columns with nonzeros 94916371a99SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 95016371a99SBarry Smith - ozn - the other array passed to the matrix preallocation routines 95116371a99SBarry Smith 95216371a99SBarry Smith Level: intermediate 95316371a99SBarry Smith 95416371a99SBarry Smith Notes: 955a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 95616371a99SBarry Smith 95716371a99SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 95816371a99SBarry Smith 95916371a99SBarry Smith This is a MACRO not a function because it uses a bunch of variables private to the MatPreallocation.... routines. 96016371a99SBarry Smith 96116371a99SBarry Smith Concepts: preallocation^Matrix 96216371a99SBarry Smith 963d6e23781SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(), 964d6e23781SBarry Smith MatPreallocateSymmetricSetLocalBlock() 96516371a99SBarry Smith M*/ 9660298fd71SBarry 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);} 96716371a99SBarry Smith 96816371a99SBarry Smith 96916371a99SBarry Smith /*MC 9700d2a7ff2SSatish Balay MatPreallocateFinalize - Ends the block of code that will count the number of nonzeros per 971bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 972bd481603SBarry Smith 973bd481603SBarry Smith Synopsis: 974aaa7dc30SBarry Smith #include <petscmat.h> 975c1ac3661SBarry Smith PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz) 976bd481603SBarry Smith 977bd481603SBarry Smith Collective on MPI_Comm 978bd481603SBarry Smith 979bd481603SBarry Smith Input Parameters: 98016371a99SBarry Smith + dnz - the array that was be passed to the matrix preallocation routines 981bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 982bd481603SBarry Smith 983bd481603SBarry Smith Level: intermediate 984bd481603SBarry Smith 985bd481603SBarry Smith Notes: 986a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 987bd481603SBarry Smith 988bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 989bd481603SBarry Smith 9901620fd73SBarry Smith This is a MACRO not a function because it closes the { started in MatPreallocateInitialize(). 9911620fd73SBarry Smith 992bd481603SBarry Smith Concepts: preallocation^Matrix 993bd481603SBarry Smith 994d6e23781SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(), 995d6e23781SBarry Smith MatPreallocateSymmetricSetLocalBlock() 996bd481603SBarry Smith M*/ 997a89bc211SBarry Smith #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree2(dnz,onz);CHKERRQ(_4_ierr);} 9987c922b88SBarry Smith 9997b80b807SBarry Smith /* Routines unique to particular data structures */ 1000014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellGetContext(Mat,void *); 1001698d4c6aSKris Buschelman 1002014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInodeAdjustForInodes(Mat,IS*,IS*); 1003014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *); 1004698d4c6aSKris Buschelman 1005014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetColumnIndices(Mat,PetscInt[]); 1006014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetColumnIndices(Mat,PetscInt[]); 1007014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 1008014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 1009014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 1010014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJFromTriple(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*,PetscInt,PetscBool); 10117b80b807SBarry Smith 1012a96a251dSBarry Smith #define MAT_SKIP_ALLOCATION -4 1013a96a251dSBarry Smith 1014014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]); 1015014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]); 1016014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]); 1017273d9f13SBarry Smith 1018014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 1019014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 1020014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 1021014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []); 1022014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 1023014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]); 1024014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 1025014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]); 1026014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIDenseSetPreallocation(Mat,PetscScalar[]); 1027014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqDenseSetPreallocation(Mat,PetscScalar[]); 10289230625dSJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,const PetscInt*[]); 10299230625dSJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,const PetscInt*[]); 1030014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAdjCreateNonemptySubcommMat(Mat,Mat*); 1031273d9f13SBarry Smith 10322e1947a5SStefano Zampini PETSC_EXTERN PetscErrorCode MatISSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 1033cf0a3239SStefano Zampini PETSC_EXTERN PetscErrorCode MatISSetUpSF(Mat); 1034014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqDenseSetLDA(Mat,PetscInt); 1035014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDenseGetLocalMatrix(Mat,Mat*); 10361b807ce4Svictorle 1037014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStoreValues(Mat); 1038014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRetrieveValues(Mat); 10392e8a6d31SBarry Smith 1040014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDAADSetCtx(Mat,void*); 10413a7fca6bSBarry Smith 1042014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFindNonzeroRows(Mat,IS*); 10437cfe41e4SPatrick Farrell PETSC_EXTERN PetscErrorCode MatFindZeroRows(Mat,IS*); 10447b80b807SBarry Smith /* 10457b80b807SBarry Smith These routines are not usually accessed directly, rather solving is 104694b7f48cSBarry Smith done through the KSP and PC interfaces. 10477b80b807SBarry Smith */ 10487b80b807SBarry Smith 104976bdecfbSBarry Smith /*J 10508f6c3df8SBarry Smith MatOrderingType - String with the name of a PETSc matrix ordering 1051d9274352SBarry Smith 1052d9274352SBarry Smith Level: beginner 1053d9274352SBarry Smith 1054d9274352SBarry Smith .seealso: MatGetOrdering() 105576bdecfbSBarry Smith J*/ 105619fd82e9SBarry Smith typedef const char* MatOrderingType; 10572692d6eeSBarry Smith #define MATORDERINGNATURAL "natural" 10582692d6eeSBarry Smith #define MATORDERINGND "nd" 10592692d6eeSBarry Smith #define MATORDERING1WD "1wd" 10602692d6eeSBarry Smith #define MATORDERINGRCM "rcm" 10612692d6eeSBarry Smith #define MATORDERINGQMD "qmd" 10622692d6eeSBarry Smith #define MATORDERINGROWLENGTH "rowlength" 1063510184a7SJed Brown #define MATORDERINGWBM "wbm" 1064fc1bef75SJed Brown #define MATORDERINGSPECTRAL "spectral" 1065312e372aSBarry Smith #define MATORDERINGAMD "amd" /* only works if UMFPACK is installed with PETSc */ 1066b12f92e5SBarry Smith 106719fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatGetOrdering(Mat,MatOrderingType,IS*,IS*); 1068140e18c1SBarry Smith PETSC_EXTERN PetscErrorCode MatGetOrderingList(PetscFunctionList*); 1069bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatOrderingRegister(const char[],PetscErrorCode(*)(Mat,MatOrderingType,IS*,IS*)); 1070140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatOrderingList; 1071d4fbbf0eSBarry Smith 1072014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS); 1073fc1bef75SJed Brown PETSC_EXTERN PetscErrorCode MatCreateLaplacian(Mat,PetscReal,PetscBool,Mat*); 1074a2ce50c7SBarry Smith 1075d91e6319SBarry Smith /*S 1076d90ac83dSHong Zhang MatFactorShiftType - Numeric Shift. 1077d90ac83dSHong Zhang 1078d90ac83dSHong Zhang Level: beginner 1079d90ac83dSHong Zhang 1080d90ac83dSHong Zhang S*/ 1081d90ac83dSHong Zhang typedef enum {MAT_SHIFT_NONE,MAT_SHIFT_NONZERO,MAT_SHIFT_POSITIVE_DEFINITE,MAT_SHIFT_INBLOCKS} MatFactorShiftType; 10826a6fc655SJed Brown PETSC_EXTERN const char *const MatFactorShiftTypes[]; 10835e9742b9SJed Brown PETSC_EXTERN const char *const MatFactorShiftTypesDetail[]; 1084d90ac83dSHong Zhang 1085d90ac83dSHong Zhang /*S 10869aa7eafdSHong Zhang MatFactorError - indicates what type of error in matrix factor 10879aa7eafdSHong Zhang 10889aa7eafdSHong Zhang Level: beginner 10899aa7eafdSHong Zhang 109000e125f8SBarry Smith Developer Notes: Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 109100e125f8SBarry Smith 109200e125f8SBarry Smith .seealso: MatGetFactor() 10939aa7eafdSHong Zhang S*/ 10945cd7cf9dSHong Zhang typedef enum {MAT_FACTOR_NOERROR,MAT_FACTOR_STRUCT_ZEROPIVOT,MAT_FACTOR_NUMERIC_ZEROPIVOT,MAT_FACTOR_OUTMEMORY,MAT_FACTOR_OTHER} MatFactorError; 10959aa7eafdSHong Zhang 109600e125f8SBarry Smith PETSC_EXTERN PetscErrorCode MatFactorGetError(Mat,MatFactorError*); 1097b8b68cfdSBarry Smith PETSC_EXTERN PetscErrorCode MatFactorClearError(Mat); 10987b6c816cSBarry Smith PETSC_EXTERN PetscErrorCode MatFactorGetErrorZeroPivot(Mat,PetscReal*,PetscInt*); 109900e125f8SBarry Smith 11009aa7eafdSHong Zhang /*S 1101422a814eSBarry Smith MatFactorInfo - Data passed into the matrix factorization routines, and information about the resulting factorization 1102b00f7748SHong Zhang 110361cad860SBarry Smith In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE, that is use 110461cad860SBarry Smith $ MatFactorInfo info(MAT_FACTORINFO_SIZE) 1105b00f7748SHong Zhang 110615e8a5b3SHong Zhang Notes: These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC. 1107b00f7748SHong Zhang 110861cad860SBarry Smith You can use MatFactorInfoInitialize() to set default values. 110961cad860SBarry Smith 1110b00f7748SHong Zhang Level: developer 1111b00f7748SHong Zhang 1112d7d82daaSBarry Smith .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(), 1113d7d82daaSBarry Smith MatFactorInfoInitialize() 1114b00f7748SHong Zhang 1115b00f7748SHong Zhang S*/ 1116b00f7748SHong Zhang typedef struct { 111715e8a5b3SHong Zhang PetscReal diagonal_fill; /* force diagonal to fill in if initially not filled */ 111885317021SBarry Smith PetscReal usedt; 111915e8a5b3SHong Zhang PetscReal dt; /* drop tolerance */ 1120b00f7748SHong Zhang PetscReal dtcol; /* tolerance for pivoting */ 112115e8a5b3SHong Zhang PetscReal dtcount; /* maximum nonzeros to be allowed per row */ 112267e28bfeSBarry Smith PetscReal fill; /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */ 1123348344bbSBarry Smith PetscReal levels; /* ICC/ILU(levels) */ 1124bcd9e38bSBarry Smith PetscReal pivotinblocks; /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0 1125bcd9e38bSBarry Smith factorization may be faster if do not pivot */ 112615e8a5b3SHong Zhang PetscReal zeropivot; /* pivot is called zero if less than this */ 1127f4db908eSBarry Smith PetscReal shifttype; /* type of shift added to matrix factor to prevent zero pivots */ 1128f4db908eSBarry Smith PetscReal shiftamount; /* how large the shift is */ 112915e8a5b3SHong Zhang } MatFactorInfo; 1130ffa6d0a5SLois Curfman McInnes 1131014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFactorInfoInitialize(MatFactorInfo*); 11329a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatCholeskyFactor(Mat,IS,const MatFactorInfo*); 11339a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatCholeskyFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*); 11349a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatCholeskyFactorNumeric(Mat,Mat,const MatFactorInfo*); 11359a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatLUFactor(Mat,IS,IS,const MatFactorInfo*); 11369a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatILUFactor(Mat,IS,IS,const MatFactorInfo*); 11379a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatLUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*); 11389a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatILUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*); 11399a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatICCFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*); 11409a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatICCFactor(Mat,IS,const MatFactorInfo*); 11419a625307SHong Zhang PETSC_EXTERN PetscErrorCode MatLUFactorNumeric(Mat,Mat,const MatFactorInfo*); 1142014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*); 1143014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolve(Mat,Vec,Vec); 1144014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatForwardSolve(Mat,Vec,Vec); 1145014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatBackwardSolve(Mat,Vec,Vec); 1146014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveAdd(Mat,Vec,Vec,Vec); 1147014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveTranspose(Mat,Vec,Vec); 1148014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveTransposeAdd(Mat,Vec,Vec,Vec); 1149014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolves(Mat,Vecs,Vecs); 11507c2f51b8SStefano Zampini PETSC_EXTERN PetscErrorCode MatSetUnfactored(Mat); 11517c2f51b8SStefano Zampini 11527c2f51b8SStefano Zampini typedef enum {MAT_FACTOR_SCHUR_UNFACTORED, MAT_FACTOR_SCHUR_FACTORED, MAT_FACTOR_SCHUR_INVERTED} MatFactorSchurStatus; 11538e7ba810SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorSetSchurIS(Mat,IS); 11547c2f51b8SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorGetSchurComplement(Mat,Mat*,MatFactorSchurStatus*); 11557c2f51b8SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorRestoreSchurComplement(Mat,Mat*,MatFactorSchurStatus); 11565a05ddb0SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorInvertSchurComplement(Mat); 11577c2f51b8SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorCreateSchurComplement(Mat,Mat*,MatFactorSchurStatus*); 11585a05ddb0SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorSolveSchurComplement(Mat,Vec,Vec); 11595a05ddb0SStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorSolveSchurComplementTranspose(Mat,Vec,Vec); 11606dba178dSStefano Zampini PETSC_EXTERN PetscErrorCode MatFactorFactorizeSchurComplement(Mat); 1161bb5a7306SBarry Smith 1162d91e6319SBarry Smith /*E 1163d91e6319SBarry Smith MatSORType - What type of (S)SOR to perform 1164bb1eb677SSatish Balay 1165d91e6319SBarry Smith Level: beginner 1166d91e6319SBarry Smith 1167d9274352SBarry Smith May be bitwise ORd together 1168d9274352SBarry Smith 1169af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 1170d91e6319SBarry Smith 11714e7234bfSBarry Smith MatSORType may be bitwise ORd together, so do not change the numbers 11724e7234bfSBarry Smith 117341f059aeSBarry Smith .seealso: MatSOR() 1174d91e6319SBarry Smith E*/ 1175ee50ffe9SBarry Smith typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3, 1176ee50ffe9SBarry Smith SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8, 1177ee50ffe9SBarry Smith SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16, 117884cb2905SBarry Smith SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType; 1179014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSOR(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec); 11808ed539a5SBarry Smith 1181d4fbbf0eSBarry Smith /* 1182639f9d9dSBarry Smith These routines are for efficiently computing Jacobians via finite differences. 1183639f9d9dSBarry Smith */ 1184b12f92e5SBarry Smith 11857086a01eSPeter Brune /*S 11867086a01eSPeter Brune MatColoring - Object for managing the coloring of matrices. 11877086a01eSPeter Brune 11887086a01eSPeter Brune Level: beginner 11897086a01eSPeter Brune 11907086a01eSPeter Brune Concepts: matrix, coloring 11917086a01eSPeter Brune 11927086a01eSPeter Brune .seealso: MatFDColoringCreate() ISColoring MatFDColoring 11937086a01eSPeter Brune S*/ 11947086a01eSPeter Brune typedef struct _p_MatColoring* MatColoring; 119576bdecfbSBarry Smith /*J 11968f6c3df8SBarry Smith MatColoringType - String with the name of a PETSc matrix coloring 1197d9274352SBarry Smith 1198d9274352SBarry Smith Level: beginner 1199d9274352SBarry Smith 12007086a01eSPeter Brune .seealso: MatColoringSetType(), MatColoring 120176bdecfbSBarry Smith J*/ 12027086a01eSPeter Brune 120319fd82e9SBarry Smith typedef const char* MatColoringType; 12045567d71aSPeter Brune #define MATCOLORINGJP "jp" 1205b9acb617SPeter Brune #define MATCOLORINGPOWER "power" 12062692d6eeSBarry Smith #define MATCOLORINGNATURAL "natural" 12072692d6eeSBarry Smith #define MATCOLORINGSL "sl" 12082692d6eeSBarry Smith #define MATCOLORINGLF "lf" 12092692d6eeSBarry Smith #define MATCOLORINGID "id" 12108121bdceSPeter Brune #define MATCOLORINGGREEDY "greedy" 1211b12f92e5SBarry Smith 12128ac6da0aSPeter Brune /*E 12138ac6da0aSPeter Brune MatColoringWeightType - Type of weight scheme 12148ac6da0aSPeter Brune 12158ac6da0aSPeter Brune Not Collective 12168ac6da0aSPeter Brune 12178ac6da0aSPeter Brune + MAT_COLORING_RANDOM - Random weights 12188ac6da0aSPeter Brune . MAT_COLORING_LEXICAL - Lexical weighting based upon global numbering. 12198ac6da0aSPeter Brune - MAT_COLORING_LF - Last-first weighting. 12208ac6da0aSPeter Brune 12218ac6da0aSPeter Brune Level: intermediate 12228ac6da0aSPeter Brune 1223af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 12248ac6da0aSPeter Brune 12258ac6da0aSPeter Brune .seealso: MatCUSPSetFormat(), MatCUSPFormatOperation 12268ac6da0aSPeter Brune E*/ 1227301c140bSPeter Brune typedef enum {MAT_COLORING_WEIGHT_RANDOM,MAT_COLORING_WEIGHT_LEXICAL,MAT_COLORING_WEIGHT_LF,MAT_COLORING_WEIGHT_SL} MatColoringWeightType; 12288ac6da0aSPeter Brune 1229335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringCreate(Mat,MatColoring*); 1230d7f2a307SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringGetDegrees(Mat,PetscInt,PetscInt*); 1231335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringDestroy(MatColoring*); 1232335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringView(MatColoring,PetscViewer); 1233335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetType(MatColoring,MatColoringType); 1234335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetFromOptions(MatColoring); 1235335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetDistance(MatColoring,PetscInt); 1236335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringGetDistance(MatColoring,PetscInt*); 1237335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetMaxColors(MatColoring,PetscInt); 1238335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringGetMaxColors(MatColoring,PetscInt*); 1239335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringApply(MatColoring,ISColoring*); 1240335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringRegister(const char[],PetscErrorCode(*)(MatColoring)); 1241014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*); 12428ac6da0aSPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetWeightType(MatColoring,MatColoringWeightType); 1243c29971fcSPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetWeights(MatColoring,PetscReal*,PetscInt*); 12448ac6da0aSPeter Brune PETSC_EXTERN PetscErrorCode MatColoringCreateWeights(MatColoring,PetscReal **,PetscInt **lperm); 1245cb394dc5SBarry Smith PETSC_EXTERN PetscErrorCode MatColoringTestValid(MatColoring,ISColoring); 1246639f9d9dSBarry Smith 1247d9274352SBarry Smith /*S 1248d9274352SBarry Smith MatFDColoring - Object for computing a sparse Jacobian via finite differences 1249d9274352SBarry Smith and coloring 1250639f9d9dSBarry Smith 1251d9274352SBarry Smith Level: beginner 1252d9274352SBarry Smith 1253d9274352SBarry Smith Concepts: coloring, sparse Jacobian, finite differences 1254d9274352SBarry Smith 1255d9274352SBarry Smith .seealso: MatFDColoringCreate() 1256d9274352SBarry Smith S*/ 1257e2a1c21fSSatish Balay typedef struct _p_MatFDColoring* MatFDColoring; 1258639f9d9dSBarry Smith 1259014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringCreate(Mat,ISColoring,MatFDColoring *); 1260014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringDestroy(MatFDColoring*); 1261014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringView(MatFDColoring,PetscViewer); 1262014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*); 1263014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**); 1264014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal); 1265014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetFromOptions(MatFDColoring); 1266d1e9a80fSBarry Smith PETSC_EXTERN PetscErrorCode MatFDColoringApply(Mat,MatFDColoring,Vec,void *); 1267014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetF(MatFDColoring,Vec); 1268edaa7c33SBarry Smith PETSC_EXTERN PetscErrorCode MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,const PetscInt*[]); 1269f86b9fbaSHong Zhang PETSC_EXTERN PetscErrorCode MatFDColoringSetUp(Mat,ISColoring,MatFDColoring); 1270f86b9fbaSHong Zhang PETSC_EXTERN PetscErrorCode MatFDColoringSetBlockSize(MatFDColoring,PetscInt,PetscInt); 1271f86b9fbaSHong Zhang 1272b1683b59SHong Zhang 1273b1683b59SHong Zhang /*S 1274b9af6bddSHong Zhang MatTransposeColoring - Object for computing a sparse matrix product C=A*B^T via coloring 1275b1683b59SHong Zhang 1276b1683b59SHong Zhang Level: beginner 1277b1683b59SHong Zhang 1278b1683b59SHong Zhang Concepts: coloring, sparse matrix product 1279b1683b59SHong Zhang 1280b9af6bddSHong Zhang .seealso: MatTransposeColoringCreate() 1281b1683b59SHong Zhang S*/ 1282b9af6bddSHong Zhang typedef struct _p_MatTransposeColoring* MatTransposeColoring; 1283b1683b59SHong Zhang 1284014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeColoringCreate(Mat,ISColoring,MatTransposeColoring *); 1285014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransColoringApplySpToDen(MatTransposeColoring,Mat,Mat); 1286014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransColoringApplyDenToSp(MatTransposeColoring,Mat,Mat); 1287014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeColoringDestroy(MatTransposeColoring*); 1288b1683b59SHong Zhang 1289639f9d9dSBarry Smith /* 12900752156aSBarry Smith These routines are for partitioning matrices: currently used only 12913eda8832SBarry Smith for adjacency matrix, MatCreateMPIAdj(). 12920752156aSBarry Smith */ 1293ca161407SBarry Smith 1294d9274352SBarry Smith /*S 1295d9274352SBarry Smith MatPartitioning - Object for managing the partitioning of a matrix or graph 1296d9274352SBarry Smith 1297d9274352SBarry Smith Level: beginner 1298d9274352SBarry Smith 1299d9274352SBarry Smith Concepts: partitioning 1300d9274352SBarry Smith 1301743a1026Svictorle .seealso: MatPartitioningCreate(), MatPartitioningType 1302d9274352SBarry Smith S*/ 130391e9ee9fSBarry Smith typedef struct _p_MatPartitioning* MatPartitioning; 1304d9274352SBarry Smith 130576bdecfbSBarry Smith /*J 13068f6c3df8SBarry Smith MatPartitioningType - String with the name of a PETSc matrix partitioning 1307d9274352SBarry Smith 1308d9274352SBarry Smith Level: beginner 13090d04baf8SBarry Smith dm 1310b547a13bSvictorle .seealso: MatPartitioningCreate(), MatPartitioning 131176bdecfbSBarry Smith J*/ 131219fd82e9SBarry Smith typedef const char* MatPartitioningType; 13132692d6eeSBarry Smith #define MATPARTITIONINGCURRENT "current" 1314aa1e27eaSFande Kong #define MATPARTITIONINGAVERAGE "average" 13152692d6eeSBarry Smith #define MATPARTITIONINGSQUARE "square" 13162692d6eeSBarry Smith #define MATPARTITIONINGPARMETIS "parmetis" 13172692d6eeSBarry Smith #define MATPARTITIONINGCHACO "chaco" 13182692d6eeSBarry Smith #define MATPARTITIONINGPARTY "party" 131950d91057SBarry Smith #define MATPARTITIONINGPTSCOTCH "ptscotch" 132088d2ac2bSFande Kong #define MATPARTITIONINGHIERARCH "hierarch" 132171306c60Sjroman 1322ca161407SBarry Smith 1323014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningCreate(MPI_Comm,MatPartitioning*); 132419fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningSetType(MatPartitioning,MatPartitioningType); 1325014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetNParts(MatPartitioning,PetscInt); 1326014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetAdjacency(MatPartitioning,Mat); 1327014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]); 1328014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []); 1329014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningApply(MatPartitioning,IS*); 1330014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningDestroy(MatPartitioning*); 13312aabb6bbSBarry Smith 1332bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningRegister(const char[],PetscErrorCode (*)(MatPartitioning)); 133330de9b25SBarry Smith 1334f1144a30SSatish Balay 13352bad1931SBarry Smith 1336014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningView(MatPartitioning,PetscViewer); 1337014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetFromOptions(MatPartitioning); 133819fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningGetType(MatPartitioning,MatPartitioningType*); 1339ca161407SBarry Smith 134027538973SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetRepartition(MatPartitioning part); 1341014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetCoarseSequential(MatPartitioning); 1342014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningParmetisGetEdgeCut(MatPartitioning, PetscInt *); 13430752156aSBarry Smith 1344b6956c12SJose Roman typedef enum { MP_CHACO_MULTILEVEL=1,MP_CHACO_SPECTRAL=2,MP_CHACO_LINEAR=4,MP_CHACO_RANDOM=5,MP_CHACO_SCATTERED=6 } MPChacoGlobalType; 13456a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoGlobalTypes[]; 1346b6956c12SJose Roman typedef enum { MP_CHACO_KERNIGHAN=1,MP_CHACO_NONE=2 } MPChacoLocalType; 13476a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoLocalTypes[]; 1348b6956c12SJose Roman typedef enum { MP_CHACO_LANCZOS=0,MP_CHACO_RQI=1 } MPChacoEigenType; 13496a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoEigenTypes[]; 1350b6956c12SJose Roman 1351014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetGlobal(MatPartitioning,MPChacoGlobalType); 1352014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetGlobal(MatPartitioning,MPChacoGlobalType*); 1353014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetLocal(MatPartitioning,MPChacoLocalType); 1354014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetLocal(MatPartitioning,MPChacoLocalType*); 1355014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal); 1356014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType); 1357014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenSolver(MatPartitioning,MPChacoEigenType*); 1358014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenTol(MatPartitioning,PetscReal); 1359014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenTol(MatPartitioning,PetscReal*); 1360014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenNumber(MatPartitioning,PetscInt); 1361014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenNumber(MatPartitioning,PetscInt*); 136271306c60Sjroman 136371306c60Sjroman #define MP_PARTY_OPT "opt" 136471306c60Sjroman #define MP_PARTY_LIN "lin" 136571306c60Sjroman #define MP_PARTY_SCA "sca" 136671306c60Sjroman #define MP_PARTY_RAN "ran" 136771306c60Sjroman #define MP_PARTY_GBF "gbf" 136871306c60Sjroman #define MP_PARTY_GCF "gcf" 136971306c60Sjroman #define MP_PARTY_BUB "bub" 137071306c60Sjroman #define MP_PARTY_DEF "def" 1371014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetGlobal(MatPartitioning,const char*); 137271306c60Sjroman #define MP_PARTY_HELPFUL_SETS "hs" 137371306c60Sjroman #define MP_PARTY_KERNIGHAN_LIN "kl" 137471306c60Sjroman #define MP_PARTY_NONE "no" 1375014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetLocal(MatPartitioning,const char*); 1376014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal); 1377014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetBipart(MatPartitioning,PetscBool); 1378014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscBool); 137971306c60Sjroman 138050d91057SBarry Smith typedef enum { MP_PTSCOTCH_QUALITY,MP_PTSCOTCH_SPEED,MP_PTSCOTCH_BALANCE,MP_PTSCOTCH_SAFETY,MP_PTSCOTCH_SCALABILITY } MPPTScotchStrategyType; 13816a6fc655SJed Brown PETSC_EXTERN const char *const MPPTScotchStrategyTypes[]; 1382e0f1cffaSJose Roman 1383014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetImbalance(MatPartitioning,PetscReal); 1384014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetImbalance(MatPartitioning,PetscReal*); 1385014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetStrategy(MatPartitioning,MPPTScotchStrategyType); 1386014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetStrategy(MatPartitioning,MPPTScotchStrategyType*); 138771306c60Sjroman 1388b43b03e9SMark F. Adams /* 13896294fa2bSFande Kong * hierarchical partitioning 13906294fa2bSFande Kong */ 13914e713f55SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalGetFineparts(MatPartitioning,IS*); 13924e713f55SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalGetCoarseparts(MatPartitioning,IS*); 13934e713f55SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalSetNcoarseparts(MatPartitioning,PetscInt); 13944e713f55SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalSetNfineparts(MatPartitioning, PetscInt); 13956294fa2bSFande Kong 1396014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMeshToVertexGraph(Mat,PetscInt,Mat*); 1397014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMeshToCellGraph(Mat,PetscInt,Mat*); 1398591294e4SBarry Smith 13990752156aSBarry Smith /* 1400af0996ceSBarry Smith If you add entries here you must also add them to petsc/finclude/petscmat.h 1401d4fbbf0eSBarry Smith */ 14021c1c02c0SLois Curfman McInnes typedef enum { MATOP_SET_VALUES=0, 14031c1c02c0SLois Curfman McInnes MATOP_GET_ROW=1, 14041c1c02c0SLois Curfman McInnes MATOP_RESTORE_ROW=2, 14051c1c02c0SLois Curfman McInnes MATOP_MULT=3, 14061c1c02c0SLois Curfman McInnes MATOP_MULT_ADD=4, 14077c922b88SBarry Smith MATOP_MULT_TRANSPOSE=5, 14087c922b88SBarry Smith MATOP_MULT_TRANSPOSE_ADD=6, 14091c1c02c0SLois Curfman McInnes MATOP_SOLVE=7, 14101c1c02c0SLois Curfman McInnes MATOP_SOLVE_ADD=8, 14117c922b88SBarry Smith MATOP_SOLVE_TRANSPOSE=9, 14127c922b88SBarry Smith MATOP_SOLVE_TRANSPOSE_ADD=10, 14131c1c02c0SLois Curfman McInnes MATOP_LUFACTOR=11, 14141c1c02c0SLois Curfman McInnes MATOP_CHOLESKYFACTOR=12, 1415a714c06dSBarry Smith MATOP_SOR=13, 14161c1c02c0SLois Curfman McInnes MATOP_TRANSPOSE=14, 14171c1c02c0SLois Curfman McInnes MATOP_GETINFO=15, 14181c1c02c0SLois Curfman McInnes MATOP_EQUAL=16, 14191c1c02c0SLois Curfman McInnes MATOP_GET_DIAGONAL=17, 14201c1c02c0SLois Curfman McInnes MATOP_DIAGONAL_SCALE=18, 14211c1c02c0SLois Curfman McInnes MATOP_NORM=19, 14221c1c02c0SLois Curfman McInnes MATOP_ASSEMBLY_BEGIN=20, 14231c1c02c0SLois Curfman McInnes MATOP_ASSEMBLY_END=21, 1424d519adbfSMatthew Knepley MATOP_SET_OPTION=22, 1425d519adbfSMatthew Knepley MATOP_ZERO_ENTRIES=23, 1426d519adbfSMatthew Knepley MATOP_ZERO_ROWS=24, 1427d519adbfSMatthew Knepley MATOP_LUFACTOR_SYMBOLIC=25, 1428d519adbfSMatthew Knepley MATOP_LUFACTOR_NUMERIC=26, 1429d519adbfSMatthew Knepley MATOP_CHOLESKY_FACTOR_SYMBOLIC=27, 1430d519adbfSMatthew Knepley MATOP_CHOLESKY_FACTOR_NUMERIC=28, 1431d519adbfSMatthew Knepley MATOP_SETUP_PREALLOCATION=29, 1432d519adbfSMatthew Knepley MATOP_ILUFACTOR_SYMBOLIC=30, 1433d519adbfSMatthew Knepley MATOP_ICCFACTOR_SYMBOLIC=31, 1434a5b7ff6bSBarry Smith MATOP_GET_DIAGONAL_BLOCK=32, 14359d39f69dSJed Brown /* MATOP_PLACEHOLDER_33=33, */ 1436d519adbfSMatthew Knepley MATOP_DUPLICATE=34, 1437d519adbfSMatthew Knepley MATOP_FORWARD_SOLVE=35, 1438d519adbfSMatthew Knepley MATOP_BACKWARD_SOLVE=36, 1439d519adbfSMatthew Knepley MATOP_ILUFACTOR=37, 1440d519adbfSMatthew Knepley MATOP_ICCFACTOR=38, 1441d519adbfSMatthew Knepley MATOP_AXPY=39, 14427dae84e0SHong Zhang MATOP_CREATE_SUBMATRICES=40, 1443d519adbfSMatthew Knepley MATOP_INCREASE_OVERLAP=41, 1444d519adbfSMatthew Knepley MATOP_GET_VALUES=42, 1445d519adbfSMatthew Knepley MATOP_COPY=43, 1446d519adbfSMatthew Knepley MATOP_GET_ROW_MAX=44, 1447d519adbfSMatthew Knepley MATOP_SCALE=45, 1448d519adbfSMatthew Knepley MATOP_SHIFT=46, 144935153367SBarry Smith MATOP_DIAGONAL_SET=47, 14509d39f69dSJed Brown MATOP_ZERO_ROWS_COLUMNS=48, 14519d39f69dSJed Brown MATOP_SET_RANDOM=49, 1452d519adbfSMatthew Knepley MATOP_GET_ROW_IJ=50, 1453d519adbfSMatthew Knepley MATOP_RESTORE_ROW_IJ=51, 1454d519adbfSMatthew Knepley MATOP_GET_COLUMN_IJ=52, 1455d519adbfSMatthew Knepley MATOP_RESTORE_COLUMN_IJ=53, 1456d519adbfSMatthew Knepley MATOP_FDCOLORING_CREATE=54, 1457d519adbfSMatthew Knepley MATOP_COLORING_PATCH=55, 1458d519adbfSMatthew Knepley MATOP_SET_UNFACTORED=56, 1459d519adbfSMatthew Knepley MATOP_PERMUTE=57, 1460d519adbfSMatthew Knepley MATOP_SET_VALUES_BLOCKED=58, 14617dae84e0SHong Zhang MATOP_CREATE_SUBMATRIX=59, 1462d519adbfSMatthew Knepley MATOP_DESTROY=60, 1463d519adbfSMatthew Knepley MATOP_VIEW=61, 1464d519adbfSMatthew Knepley MATOP_CONVERT_FROM=62, 14657bab7c10SHong Zhang MATOP_MATMAT_MULT=63, 14667bab7c10SHong Zhang MATOP_MATMAT_MULT_SYMBOLIC=64, 14677bab7c10SHong Zhang MATOP_MATMAT_MULT_NUMERIC=65, 1468d519adbfSMatthew Knepley MATOP_SET_LOCAL_TO_GLOBAL_MAP=66, 1469d519adbfSMatthew Knepley MATOP_SET_VALUES_LOCAL=67, 1470d519adbfSMatthew Knepley MATOP_ZERO_ROWS_LOCAL=68, 1471d519adbfSMatthew Knepley MATOP_GET_ROW_MAX_ABS=69, 1472d519adbfSMatthew Knepley MATOP_GET_ROW_MIN_ABS=70, 1473d519adbfSMatthew Knepley MATOP_CONVERT=71, 1474d519adbfSMatthew Knepley MATOP_SET_COLORING=72, 14759d39f69dSJed Brown /* MATOP_PLACEHOLDER_73=73, */ 1476d519adbfSMatthew Knepley MATOP_SET_VALUES_ADIFOR=74, 1477d519adbfSMatthew Knepley MATOP_FD_COLORING_APPLY=75, 1478d519adbfSMatthew Knepley MATOP_SET_FROM_OPTIONS=76, 1479bda6c4cbSBarry Smith MATOP_MULT_CONSTRAINED=77, 1480bda6c4cbSBarry Smith MATOP_MULT_TRANSPOSE_CONSTRAIN=78, 14819d39f69dSJed Brown MATOP_FIND_ZERO_DIAGONALS=79, 1482d519adbfSMatthew Knepley MATOP_MULT_MULTIPLE=80, 1483d519adbfSMatthew Knepley MATOP_SOLVE_MULTIPLE=81, 1484d519adbfSMatthew Knepley MATOP_GET_INERTIA=82, 1485d519adbfSMatthew Knepley MATOP_LOAD=83, 1486d519adbfSMatthew Knepley MATOP_IS_SYMMETRIC=84, 1487d519adbfSMatthew Knepley MATOP_IS_HERMITIAN=85, 1488d519adbfSMatthew Knepley MATOP_IS_STRUCTURALLY_SYMMETRIC=86, 1489820469bcSHong Zhang MATOP_SET_VALUES_BLOCKEDLOCAL=87, 1490b41ce5d5SBarry Smith MATOP_CREATE_VECS=88, 1491d519adbfSMatthew Knepley MATOP_MAT_MULT=89, 1492d519adbfSMatthew Knepley MATOP_MAT_MULT_SYMBOLIC=90, 1493d519adbfSMatthew Knepley MATOP_MAT_MULT_NUMERIC=91, 1494d519adbfSMatthew Knepley MATOP_PTAP=92, 1495d519adbfSMatthew Knepley MATOP_PTAP_SYMBOLIC=93, 1496d519adbfSMatthew Knepley MATOP_PTAP_NUMERIC=94, 1497bda6c4cbSBarry Smith MATOP_MAT_TRANSPOSE_MULT=95, 1498bda6c4cbSBarry Smith MATOP_MAT_TRANSPOSE_MULT_SYMBO=96, 1499bda6c4cbSBarry Smith MATOP_MAT_TRANSPOSE_MULT_NUMER=97, 15009d39f69dSJed Brown /* MATOP_PLACEHOLDER_98=98, */ 15019d39f69dSJed Brown /* MATOP_PLACEHOLDER_99=99, */ 15029d39f69dSJed Brown /* MATOP_PLACEHOLDER_100=100, */ 15039d39f69dSJed Brown /* MATOP_PLACEHOLDER_101=101, */ 1504d519adbfSMatthew Knepley MATOP_CONJUGATE=102, 15059d39f69dSJed Brown /* MATOP_PLACEHOLDER_103=103, */ 1506d519adbfSMatthew Knepley MATOP_SET_VALUES_ROW=104, 1507d519adbfSMatthew Knepley MATOP_REAL_PART=105, 1508bda6c4cbSBarry Smith MATOP_IMAGINARY_PART=106, 1509bda6c4cbSBarry Smith MATOP_GET_ROW_UPPER_TRIANGULAR=107, 1510bda6c4cbSBarry Smith MATOP_RESTORE_ROW_UPPER_TRIANG=108, 1511bda6c4cbSBarry Smith MATOP_MAT_SOLVE=109, 151229eadf9eSStefano Zampini MATOP_MAT_SOLVE_TRANSPOSE=110, 1513d519adbfSMatthew Knepley MATOP_GET_ROW_MIN=111, 15140e8d04f7SBarry Smith MATOP_GET_COLUMN_VECTOR=112, 1515d519adbfSMatthew Knepley MATOP_MISSING_DIAGONAL=113, 15160e8d04f7SBarry Smith MATOP_GET_SEQ_NONZERO_STRUCTUR=114, 151789c6957cSBarry Smith MATOP_CREATE=115, 151889c6957cSBarry Smith MATOP_GET_GHOSTS=116, 15190e8d04f7SBarry Smith MATOP_GET_LOCAL_SUB_MATRIX=117, 15200e8d04f7SBarry Smith MATOP_RESTORE_LOCALSUB_MATRIX=118, 1521eeffb40dSHong Zhang MATOP_MULT_DIAGONAL_BLOCK=119, 15220e8d04f7SBarry Smith MATOP_HERMITIAN_TRANSPOSE=120, 15230e8d04f7SBarry Smith MATOP_MULT_HERMITIAN_TRANSPOSE=121, 15240e8d04f7SBarry Smith MATOP_MULT_HERMITIAN_TRANS_ADD=122, 15250e8d04f7SBarry Smith MATOP_GET_MULTI_PROC_BLOCK=123, 15269d39f69dSJed Brown MATOP_FIND_NONZERO_ROWS=124, 15270e8d04f7SBarry Smith MATOP_GET_COLUMN_NORMS=125, 15289d39f69dSJed Brown MATOP_INVERT_BLOCK_DIAGONAL=126, 15299d39f69dSJed Brown /* MATOP_PLACEHOLDER_127=127, */ 1530b41ce5d5SBarry Smith MATOP_CREATE_SUB_MATRICES_MPI=128, 15312b8ad9a3SHong Zhang MATOP_SET_VALUES_BATCH=129, 15320e8d04f7SBarry Smith MATOP_TRANSPOSE_MAT_MULT=130, 15330e8d04f7SBarry Smith MATOP_TRANSPOSE_MAT_MULT_SYMBO=131, 1534e9b602ebSSatish Balay MATOP_TRANSPOSE_MAT_MULT_NUMER=132, 15350e8d04f7SBarry Smith MATOP_TRANSPOSE_COLORING_CREAT=133, 15360e8d04f7SBarry Smith MATOP_TRANS_COLORING_APPLY_SPT=134, 15370e8d04f7SBarry Smith MATOP_TRANS_COLORING_APPLY_DEN=135, 15380e8d04f7SBarry Smith MATOP_RART=136, 15390e8d04f7SBarry Smith MATOP_RART_SYMBOLIC=137, 15400e8d04f7SBarry Smith MATOP_RART_NUMERIC=138, 1541e09a3074SHong Zhang MATOP_SET_BLOCK_SIZES=139, 1542f9426fe0SMark Adams MATOP_AYPX=140, 15431919a2e2SJed Brown MATOP_RESIDUAL=141, 15449c8f2541SHong Zhang MATOP_FDCOLORING_SETUP=142, 15452d033e1fSHong Zhang MATOP_MPICONCATENATESEQ=144, 15462d033e1fSHong Zhang MATOP_DESTROYSUBMATRICES=145 1547fae171e0SBarry Smith } MatOperation; 1548014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatHasOperation(Mat,MatOperation,PetscBool *); 1549014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellSetOperation(Mat,MatOperation,void(*)(void)); 1550014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellGetOperation(Mat,MatOperation,void(**)(void)); 1551014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellSetContext(Mat,void*); 1552*f3b1f45cSBarry Smith PETSC_EXTERN PetscErrorCode MatShellTestMult(Mat,PetscErrorCode (*)(void*,Vec,Vec),Vec,void*,PetscBool*); 1553*f3b1f45cSBarry Smith PETSC_EXTERN PetscErrorCode MatShellTestMultTranspose(Mat,PetscErrorCode (*)(void*,Vec,Vec),Vec,void*,PetscBool*); 1554112a2221SBarry Smith 155590ace30eSBarry Smith /* 155690ace30eSBarry Smith Codes for matrices stored on disk. By default they are 155790ace30eSBarry Smith stored in a universal format. By changing the format with 15586a9046bcSBarry Smith PetscViewerPushFormat(viewer,PETSC_VIEWER_NATIVE); the matrices will 155990ace30eSBarry Smith be stored in a way natural for the matrix, for example dense matrices 156090ace30eSBarry Smith would be stored as dense. Matrices stored this way may only be 1561f4403165SShri Abhyankar read into matrices of the same type. 156290ace30eSBarry Smith */ 156390ace30eSBarry Smith #define MATRIX_BINARY_FORMAT_DENSE -1 156490ace30eSBarry Smith 1565014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetHashTableFactor(Mat,PetscReal); 1566014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatISGetLocalMat(Mat,Mat*); 15673b3b1effSJed Brown PETSC_EXTERN PetscErrorCode MatISRestoreLocalMat(Mat,Mat*); 1568014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatISSetLocalMat(Mat,Mat); 1569b7ce53b6SStefano Zampini PETSC_EXTERN PetscErrorCode MatISGetMPIXAIJ(Mat,MatReuse,Mat*); 15701f4e1ec7SBarry Smith 1571d9274352SBarry Smith /*S 1572d9274352SBarry Smith MatNullSpace - Object that removes a null space from a vector, i.e. 1573d9274352SBarry Smith orthogonalizes the vector to a subsapce 1574d9274352SBarry Smith 1575f7a9e4ceSBarry Smith Level: advanced 1576d9274352SBarry Smith 1577d9274352SBarry Smith Concepts: matrix; linear operator, null space 1578d9274352SBarry Smith 15796e1639daSBarry Smith Users manual sections: 15806e1639daSBarry Smith . sec_singular 15816e1639daSBarry Smith 1582d9274352SBarry Smith .seealso: MatNullSpaceCreate() 1583d9274352SBarry Smith S*/ 158474637425SBarry Smith typedef struct _p_MatNullSpace* MatNullSpace; 1585d9274352SBarry Smith 1586014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceCreate(MPI_Comm,PetscBool ,PetscInt,const Vec[],MatNullSpace*); 1587014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(MatNullSpace,Vec,void*),void*); 1588014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceDestroy(MatNullSpace*); 1589d0195637SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceRemove(MatNullSpace,Vec); 1590014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNullSpace(Mat, MatNullSpace *); 15915fa7ec2dSBarry Smith PETSC_EXTERN PetscErrorCode MatGetTransposeNullSpace(Mat, MatNullSpace *); 15925fa7ec2dSBarry Smith PETSC_EXTERN PetscErrorCode MatSetTransposeNullSpace(Mat,MatNullSpace); 1593014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNullSpace(Mat,MatNullSpace); 1594014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNearNullSpace(Mat,MatNullSpace); 1595014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNearNullSpace(Mat,MatNullSpace*); 1596014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceTest(MatNullSpace,Mat,PetscBool *); 1597014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceView(MatNullSpace,PetscViewer); 1598014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceGetVecs(MatNullSpace,PetscBool*,PetscInt*,const Vec**); 1599014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceCreateRigidBody(Vec,MatNullSpace*); 160074637425SBarry Smith 1601014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatReorderingSeqSBAIJ(Mat,IS); 1602014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetHashTableFactor(Mat,PetscReal); 1603014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetColumnIndices(Mat,PetscInt *); 1604014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJInvertBlockDiagonal(Mat); 16053f1d51d7SBarry Smith 1606014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMAIJ(Mat,PetscInt,Mat*); 1607014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMAIJRedimension(Mat,PetscInt,Mat*); 1608014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMAIJGetAIJ(Mat,Mat*); 1609c4f061fbSSatish Balay 1610014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatComputeExplicitOperator(Mat,Mat*); 1611*f3b1f45cSBarry Smith PETSC_EXTERN PetscErrorCode MatComputeExplicitOperatorTranspose(Mat,Mat*); 1612b0a32e0cSBarry Smith 1613014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalScaleLocal(Mat,Vec); 161404f1ad80SBarry Smith 1615014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*); 1616014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetBase(Mat,Vec,Vec); 1617014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*); 1618014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*)); 1619014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec)); 1620014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt); 1621014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDResetHHistory(Mat); 1622014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctionError(Mat,PetscReal); 1623014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetPeriod(Mat,PetscInt); 1624014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDGetH(Mat,PetscScalar *); 1625014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetOptionsPrefix(Mat,const char[]); 1626014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*); 1627014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*); 1628e884886fSBarry Smith 16296370053bSBarry Smith /*S 16306370053bSBarry Smith MatMFFD - A data structured used to manage the computation of the h differencing parameter for matrix-free 16316370053bSBarry Smith Jacobian vector products 1632e884886fSBarry Smith 16336370053bSBarry Smith Notes: MATMFFD is a specific MatType which uses the MatMFFD data structure 16346370053bSBarry Smith 16356370053bSBarry Smith MatMFFD*() methods actually take the Mat as their first argument. Not a MatMFFD data structure 16366370053bSBarry Smith 16376370053bSBarry Smith Level: developer 16386370053bSBarry Smith 16396370053bSBarry Smith .seealso: MATMFFD, MatCreateMFFD(), MatMFFDSetFuction(), MatMFFDSetType(), MatMFFDRegister() 16406370053bSBarry Smith S*/ 1641e884886fSBarry Smith typedef struct _p_MatMFFD* MatMFFD; 1642e884886fSBarry Smith 164376bdecfbSBarry Smith /*J 1644e884886fSBarry Smith MatMFFDType - algorithm used to compute the h used in computing matrix-vector products via differencing of the function 1645e884886fSBarry Smith 1646e884886fSBarry Smith Level: beginner 1647e884886fSBarry Smith 1648e884886fSBarry Smith .seealso: MatMFFDSetType(), MatMFFDRegister() 164976bdecfbSBarry Smith J*/ 165019fd82e9SBarry Smith typedef const char* MatMFFDType; 1651e884886fSBarry Smith #define MATMFFD_DS "ds" 1652e884886fSBarry Smith #define MATMFFD_WP "wp" 1653e884886fSBarry Smith 165419fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDSetType(Mat,MatMFFDType); 1655bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDRegister(const char[],PetscErrorCode (*)(MatMFFD)); 1656e884886fSBarry Smith 1657014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDDSSetUmin(Mat,PetscReal); 1658014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDWPSetComputeNormU(Mat,PetscBool ); 1659e884886fSBarry Smith 166061ab5df0SBarry Smith PETSC_EXTERN PetscErrorCode MatFDColoringSetType(MatFDColoring,MatMFFDType); 166161ab5df0SBarry Smith 1662014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *); 1663014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *); 16647dbadf16SMatthew Knepley 166597969023SHong Zhang /* 166697969023SHong Zhang PETSc interface to MUMPS 166797969023SHong Zhang */ 166897969023SHong Zhang #ifdef PETSC_HAVE_MUMPS 1669014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMumpsSetIcntl(Mat,PetscInt,PetscInt); 1670bc6112feSHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetIcntl(Mat,PetscInt,PetscInt*); 1671b9e9ea26SSatish Balay PETSC_EXTERN PetscErrorCode MatMumpsSetCntl(Mat,PetscInt,PetscReal); 1672bc6112feSHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetCntl(Mat,PetscInt,PetscReal*); 1673bc6112feSHong Zhang 1674ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetInfo(Mat,PetscInt,PetscInt*); 1675ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetInfog(Mat,PetscInt,PetscInt*); 1676ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetRinfo(Mat,PetscInt,PetscReal*); 1677ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetRinfog(Mat,PetscInt,PetscReal*); 167897969023SHong Zhang #endif 167923a5497aSJed Brown 1680d954db57SHong Zhang /* 1681d615f992SSatish Balay PETSc interface to Mkl_Pardiso 1682b500a6b7SJose David Bermeol */ 1683b500a6b7SJose David Bermeol #ifdef PETSC_HAVE_MKL_PARDISO 1684d615f992SSatish Balay PETSC_EXTERN PetscErrorCode MatMkl_PardisoSetCntl(Mat,PetscInt,PetscInt); 1685b500a6b7SJose David Bermeol #endif 1686b500a6b7SJose David Bermeol 1687b500a6b7SJose David Bermeol /* 1688d305a81bSVasiliy Kozyrev PETSc interface to Mkl_CPardiso 1689d305a81bSVasiliy Kozyrev */ 1690d305a81bSVasiliy Kozyrev #ifdef PETSC_HAVE_MKL_CPARDISO 1691d305a81bSVasiliy Kozyrev PETSC_EXTERN PetscErrorCode MatMkl_CPardisoSetCntl(Mat,PetscInt,PetscInt); 1692d305a81bSVasiliy Kozyrev #endif 1693d305a81bSVasiliy Kozyrev 1694d305a81bSVasiliy Kozyrev /* 1695d954db57SHong Zhang PETSc interface to SUPERLU 1696d954db57SHong Zhang */ 1697d954db57SHong Zhang #ifdef PETSC_HAVE_SUPERLU 1698014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSuperluSetILUDropTol(Mat,PetscReal); 1699d954db57SHong Zhang #endif 1700d954db57SHong Zhang 1701fb785984SHong Zhang /* 1702fb785984SHong Zhang PETSc interface to SUPERLU_DIST 1703fb785984SHong Zhang */ 1704fb785984SHong Zhang #ifdef PETSC_HAVE_SUPERLU_DIST 1705fb785984SHong Zhang PETSC_EXTERN PetscErrorCode MatSuperluDistGetDiagU(Mat,PetscScalar*); 1706fb785984SHong Zhang #endif 1707fb785984SHong Zhang 1708575704cbSPieter Ghysels /* 1709575704cbSPieter Ghysels PETSc interface to STRUMPACK 1710575704cbSPieter Ghysels */ 1711575704cbSPieter Ghysels #ifdef PETSC_HAVE_STRUMPACK 1712542aee0fSPieter Ghysels /*E 1713542aee0fSPieter Ghysels MatSTRUMPACKReordering - sparsity reducing ordering to be used in STRUMPACK 1714542aee0fSPieter Ghysels 1715542aee0fSPieter Ghysels Level: intermediate 1716542aee0fSPieter Ghysels E*/ 171734f43fa5SPieter Ghysels typedef enum {MAT_STRUMPACK_NATURAL=0, 171834f43fa5SPieter Ghysels MAT_STRUMPACK_METIS=1, 171934f43fa5SPieter Ghysels MAT_STRUMPACK_PARMETIS=2, 172034f43fa5SPieter Ghysels MAT_STRUMPACK_SCOTCH=3, 172134f43fa5SPieter Ghysels MAT_STRUMPACK_PTSCOTCH=4, 172234f43fa5SPieter Ghysels MAT_STRUMPACK_RCM=5} MatSTRUMPACKReordering; 172334f43fa5SPieter Ghysels PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetReordering(Mat,MatSTRUMPACKReordering); 1724575704cbSPieter Ghysels PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetColPerm(Mat,PetscBool); 172534f43fa5SPieter Ghysels PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSRelTol(Mat,PetscReal); 172634f43fa5SPieter Ghysels PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSAbsTol(Mat,PetscReal); 172734f43fa5SPieter Ghysels PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSMinSepSize(Mat,PetscInt); 172834f43fa5SPieter Ghysels PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSMaxRank(Mat,PetscInt); 1729a36bf211SPieter Ghysels PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSLeafSize(Mat,PetscInt); 1730575704cbSPieter Ghysels #endif 1731575704cbSPieter Ghysels 1732575704cbSPieter Ghysels 1733aa372e3fSPaul Mullowney #ifdef PETSC_HAVE_CUDA 17343f9c0db1SPaul Mullowney /*E 1735e057df02SPaul Mullowney MatCUSPARSEStorageFormat - indicates the storage format for CUSPARSE (GPU) 17362692e278SPaul Mullowney matrices. 1737e057df02SPaul Mullowney 1738e057df02SPaul Mullowney Not Collective 1739e057df02SPaul Mullowney 17408468deeeSKarl Rupp + MAT_CUSPARSE_CSR - Compressed Sparse Row 17412692e278SPaul Mullowney . MAT_CUSPARSE_ELL - Ellpack (requires CUDA 4.2 or later). 17422692e278SPaul Mullowney - MAT_CUSPARSE_HYB - Hybrid, a combination of Ellpack and Coordinate format (requires CUDA 4.2 or later). 1743e057df02SPaul Mullowney 1744e057df02SPaul Mullowney Level: intermediate 1745e057df02SPaul Mullowney 1746af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 1747e057df02SPaul Mullowney 1748e057df02SPaul Mullowney .seealso: MatCUSPARSESetFormat(), MatCUSPARSEFormatOperation 1749e057df02SPaul Mullowney E*/ 1750e057df02SPaul Mullowney 17513f9c0db1SPaul Mullowney typedef enum {MAT_CUSPARSE_CSR, MAT_CUSPARSE_ELL, MAT_CUSPARSE_HYB} MatCUSPARSEStorageFormat; 1752e057df02SPaul Mullowney 1753e057df02SPaul Mullowney /* these will be strings associated with enumerated type defined above */ 1754e057df02SPaul Mullowney PETSC_EXTERN const char *const MatCUSPARSEStorageFormats[]; 1755e057df02SPaul Mullowney 17563f9c0db1SPaul Mullowney /*E 1757e057df02SPaul Mullowney MatCUSPARSEFormatOperation - indicates the operation of CUSPARSE (GPU) 17582692e278SPaul Mullowney matrices whose operation should use a particular storage format. 1759e057df02SPaul Mullowney 1760e057df02SPaul Mullowney Not Collective 1761e057df02SPaul Mullowney 17628468deeeSKarl Rupp + MAT_CUSPARSE_MULT_DIAG - sets the storage format for the diagonal matrix in the parallel MatMult 17638468deeeSKarl Rupp . MAT_CUSPARSE_MULT_OFFDIAG - sets the storage format for the offdiagonal matrix in the parallel MatMult 17648468deeeSKarl Rupp . MAT_CUSPARSE_MULT - sets the storage format for the entire matrix in the serial (single GPU) MatMult 17658468deeeSKarl Rupp - MAT_CUSPARSE_ALL - sets the storage format for all CUSPARSE (GPU) matrices 1766e057df02SPaul Mullowney 1767e057df02SPaul Mullowney Level: intermediate 1768e057df02SPaul Mullowney 1769e057df02SPaul Mullowney .seealso: MatCUSPARSESetFormat(), MatCUSPARSEStorageFormat 1770e057df02SPaul Mullowney E*/ 177136d62e41SPaul Mullowney typedef enum {MAT_CUSPARSE_MULT_DIAG, MAT_CUSPARSE_MULT_OFFDIAG, MAT_CUSPARSE_MULT, MAT_CUSPARSE_ALL} MatCUSPARSEFormatOperation; 1772e057df02SPaul Mullowney 177310e9db80SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 177410e9db80SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCreateAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 1775e057df02SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCUSPARSESetFormat(Mat,MatCUSPARSEFormatOperation,MatCUSPARSEStorageFormat); 17769ae82921SPaul Mullowney #endif 17779ae82921SPaul Mullowney 177890c53902SBarry Smith #if defined(PETSC_HAVE_CUSP) 1779014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCUSP(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 1780014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateAIJCUSP(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 1781e057df02SPaul Mullowney 17823f9c0db1SPaul Mullowney /*E 1783e057df02SPaul Mullowney MatCUSPStorageFormat - indicates the storage format for CUSP (GPU) 178436d62e41SPaul Mullowney matrices. 1785e057df02SPaul Mullowney 1786e057df02SPaul Mullowney Not Collective 1787e057df02SPaul Mullowney 17888468deeeSKarl Rupp + MAT_CUSP_CSR - Compressed Sparse Row 17898468deeeSKarl Rupp . MAT_CUSP_DIA - Diagonal 17908468deeeSKarl Rupp - MAT_CUSP_ELL - Ellpack 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(), MatCUSPFormatOperation 1797e057df02SPaul Mullowney E*/ 17983f9c0db1SPaul Mullowney typedef enum {MAT_CUSP_CSR, MAT_CUSP_DIA, MAT_CUSP_ELL} MatCUSPStorageFormat; 1799e057df02SPaul Mullowney 1800e057df02SPaul Mullowney /* these will be strings associated with enumerated type defined above */ 1801e057df02SPaul Mullowney PETSC_EXTERN const char *const MatCUSPStorageFormats[]; 1802e057df02SPaul Mullowney 18033f9c0db1SPaul Mullowney /*E 1804e057df02SPaul Mullowney MatCUSPFormatOperation - indicates the operation of CUSP (GPU) 180536d62e41SPaul Mullowney matrices whose operation should use a particular storage format. 1806e057df02SPaul Mullowney 1807e057df02SPaul Mullowney Not Collective 1808e057df02SPaul Mullowney 18098468deeeSKarl Rupp + MAT_CUSP_MULT_DIAG - sets the storage format for the diagonal matrix in the parallel MatMult 18108468deeeSKarl Rupp . MAT_CUSP_MULT_OFFDIAG - sets the storage format for the offdiagonal matrix in the parallel MatMult 18118468deeeSKarl Rupp . MAT_CUSP_MULT - sets the storage format for the entire matrix in the serial (single GPU) MatMult 18128468deeeSKarl Rupp - MAT_CUSP_ALL - sets the storage format for all CUSP (GPU) matrices 1813e057df02SPaul Mullowney 1814e057df02SPaul Mullowney Level: intermediate 1815e057df02SPaul Mullowney 1816af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 1817e057df02SPaul Mullowney 1818e057df02SPaul Mullowney .seealso: MatCUSPSetFormat(), MatCUSPStorageFormat 1819e057df02SPaul Mullowney E*/ 182036d62e41SPaul Mullowney typedef enum {MAT_CUSP_MULT_DIAG, MAT_CUSP_MULT_OFFDIAG, MAT_CUSP_MULT, MAT_CUSP_ALL} MatCUSPFormatOperation; 1821e057df02SPaul Mullowney 1822e057df02SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCUSPSetFormat(Mat,MatCUSPFormatOperation,MatCUSPStorageFormat); 1823e057df02SPaul Mullowney #endif 182490c53902SBarry Smith 1825d67ff14aSKarl Rupp #if defined(PETSC_HAVE_VIENNACL) 1826d67ff14aSKarl Rupp PETSC_EXTERN PetscErrorCode MatCreateSeqAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 1827d67ff14aSKarl Rupp PETSC_EXTERN PetscErrorCode MatCreateAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 1828d67ff14aSKarl Rupp #endif 1829d67ff14aSKarl Rupp 183054efbe56SHong Zhang /* 183154efbe56SHong Zhang PETSc interface to FFTW 183254efbe56SHong Zhang */ 183354efbe56SHong Zhang #if defined(PETSC_HAVE_FFTW) 1834014dd563SJed Brown PETSC_EXTERN PetscErrorCode VecScatterPetscToFFTW(Mat,Vec,Vec); 1835014dd563SJed Brown PETSC_EXTERN PetscErrorCode VecScatterFFTWToPetsc(Mat,Vec,Vec); 18362a7a6963SBarry Smith PETSC_EXTERN PetscErrorCode MatCreateVecsFFTW(Mat,Vec*,Vec*,Vec*); 183754efbe56SHong Zhang #endif 183854efbe56SHong Zhang 1839014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateNest(MPI_Comm,PetscInt,const IS[],PetscInt,const IS[],const Mat[],Mat*); 1840014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSize(Mat,PetscInt*,PetscInt*); 1841014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetISs(Mat,IS[],IS[]); 1842014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetLocalISs(Mat,IS[],IS[]); 1843014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSubMats(Mat,PetscInt*,PetscInt*,Mat***); 1844014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSubMat(Mat,PetscInt,PetscInt,Mat*); 184519fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatNestSetVecType(Mat,VecType); 1846014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestSetSubMats(Mat,PetscInt,const IS[],PetscInt,const IS[],const Mat[]); 1847014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestSetSubMat(Mat,PetscInt,PetscInt,Mat); 1848d8588912SDave May 18494325cce7SMatthew G Knepley PETSC_EXTERN PetscErrorCode MatChop(Mat,PetscReal); 1850e0a58c10SMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatComputeBandwidth(Mat,PetscReal,PetscInt*); 18514325cce7SMatthew G Knepley 1852b541e6a4SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSubdomainsCreateCoalesce(Mat,PetscInt,PetscInt*,IS**); 185353dd7562SDmitry Karpeev 1854c094ef40SMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatPreallocatorPreallocate(Mat,PetscBool,Mat); 1855c094ef40SMatthew G. Knepley 1856539c167fSBarry Smith PETSC_INTERN PetscErrorCode MatHeaderMerge(Mat,Mat*); 1857539c167fSBarry Smith PETSC_EXTERN PetscErrorCode MatHeaderReplace(Mat,Mat*); 1858539c167fSBarry Smith 185923a5497aSJed Brown #endif 1860