12eac72dbSBarry Smith /* 22eac72dbSBarry Smith Include file for the matrix component of PETSc 32eac72dbSBarry Smith */ 40a835dfdSSatish Balay #ifndef __PETSCMAT_H 50a835dfdSSatish Balay #define __PETSCMAT_H 62c8e378dSBarry Smith #include <petscvec.h> 72eac72dbSBarry Smith 8d9274352SBarry Smith /*S 98f6c3df8SBarry Smith Mat - Abstract PETSc matrix object used to manage all linear operators in PETSc, even those without 108f6c3df8SBarry Smith an explicit sparse representation (such as matrix-free operators) 112eac72dbSBarry Smith 12d91e6319SBarry Smith Level: beginner 13d91e6319SBarry Smith 14d9274352SBarry Smith Concepts: matrix; linear operator 15d9274352SBarry Smith 168f6c3df8SBarry Smith .seealso: MatCreate(), MatType, MatSetType(), MatDestroy() 17d9274352SBarry Smith S*/ 18d9274352SBarry Smith typedef struct _p_Mat* Mat; 19d9274352SBarry Smith 2076bdecfbSBarry Smith /*J 218f6c3df8SBarry Smith MatType - String with the name of a PETSc matrix type 22d9274352SBarry Smith 23d9274352SBarry Smith Level: beginner 24d9274352SBarry Smith 258f6c3df8SBarry Smith .seealso: MatSetType(), Mat, MatSolverPackage, MatRegister() 2676bdecfbSBarry Smith J*/ 2719fd82e9SBarry Smith typedef const char* MatType; 28273d9f13SBarry Smith #define MATSAME "same" 295a11e1b2SBarry Smith #define MATMAIJ "maij" 30273d9f13SBarry Smith #define MATSEQMAIJ "seqmaij" 31273d9f13SBarry Smith #define MATMPIMAIJ "mpimaij" 32273d9f13SBarry Smith #define MATIS "is" 335a11e1b2SBarry Smith #define MATAIJ "aij" 34273d9f13SBarry Smith #define MATSEQAIJ "seqaij" 357d6a0e61SBarry Smith #define MATSEQAIJPTHREAD "seqaijpthread" 36bf2c1783SBarry Smith #define MATAIJPTHREAD "aijpthread" 37273d9f13SBarry Smith #define MATMPIAIJ "mpiaij" 385a11e1b2SBarry Smith #define MATAIJCRL "aijcrl" 395a11e1b2SBarry Smith #define MATSEQAIJCRL "seqaijcrl" 405a11e1b2SBarry Smith #define MATMPIAIJCRL "mpiaijcrl" 418154be41SBarry Smith #define MATAIJCUSP "aijcusp" 428154be41SBarry Smith #define MATSEQAIJCUSP "seqaijcusp" 438154be41SBarry Smith #define MATMPIAIJCUSP "mpiaijcusp" 449ae82921SPaul Mullowney #define MATAIJCUSPARSE "aijcusparse" 459ae82921SPaul Mullowney #define MATSEQAIJCUSPARSE "seqaijcusparse" 469ae82921SPaul Mullowney #define MATMPIAIJCUSPARSE "mpiaijcusparse" 47d67ff14aSKarl Rupp #define MATAIJVIENNACL "aijviennacl" 48d67ff14aSKarl Rupp #define MATSEQAIJVIENNACL "seqaijviennacl" 49d67ff14aSKarl Rupp #define MATMPIAIJVIENNACL "mpiaijviennacl" 505a11e1b2SBarry Smith #define MATAIJPERM "aijperm" 515a11e1b2SBarry Smith #define MATSEQAIJPERM "seqaijperm" 525a11e1b2SBarry Smith #define MATMPIAIJPERM "mpiaijperm" 53273d9f13SBarry Smith #define MATSHELL "shell" 545a11e1b2SBarry Smith #define MATDENSE "dense" 55209238afSKris Buschelman #define MATSEQDENSE "seqdense" 56273d9f13SBarry Smith #define MATMPIDENSE "mpidense" 57db31f6deSJed Brown #define MATELEMENTAL "elemental" 585a11e1b2SBarry Smith #define MATBAIJ "baij" 59273d9f13SBarry Smith #define MATSEQBAIJ "seqbaij" 60273d9f13SBarry Smith #define MATMPIBAIJ "mpibaij" 61273d9f13SBarry Smith #define MATMPIADJ "mpiadj" 625a11e1b2SBarry Smith #define MATSBAIJ "sbaij" 63273d9f13SBarry Smith #define MATSEQSBAIJ "seqsbaij" 64273d9f13SBarry Smith #define MATMPISBAIJ "mpisbaij" 65c0cdd4a1SDahai Guo #define MATSEQBSTRM "seqbstrm" 66c0cdd4a1SDahai Guo #define MATMPIBSTRM "mpibstrm" 67c0cdd4a1SDahai Guo #define MATBSTRM "bstrm" 68aa5a9175SDahai Guo #define MATSEQSBSTRM "seqsbstrm" 69aa5a9175SDahai Guo #define MATMPISBSTRM "mpisbstrm" 70aa5a9175SDahai Guo #define MATSBSTRM "sbstrm" 71cebc7f6cSBarry Smith #define MATDAAD "daad" 72cebc7f6cSBarry Smith #define MATMFFD "mffd" 73c8a8475eSBarry Smith #define MATNORMAL "normal" 74ab92ecdeSBarry Smith #define MATLRC "lrc" 752a6744ebSBarry Smith #define MATSCATTER "scatter" 76421e10b8SBarry Smith #define MATBLOCKMAT "blockmat" 77793850ffSBarry Smith #define MATCOMPOSITE "composite" 781f8c7532SHong Zhang #define MATFFT "fft" 791f8c7532SHong Zhang #define MATFFTW "fftw" 80e133240eSMatthew G Knepley #define MATSEQCUFFT "seqcufft" 81557cca28SSatish Balay #define MATTRANSPOSEMAT "transpose" 8272ca8751SBarry Smith #define MATSCHURCOMPLEMENT "schurcomplement" 831d6018f0SLisandro Dalcin #define MATPYTHON "python" 84f91d8e95SBarry Smith #define MATHYPRESTRUCT "hyprestruct" 85a9e6138eSGlenn Hammond #define MATHYPRESSTRUCT "hypresstruct" 86ee1cef2cSJed Brown #define MATSUBMATRIX "submatrix" 872c0dbf93SJed Brown #define MATLOCALREF "localref" 88d8588912SDave May #define MATNEST "nest" 89773941d6SBarry Smith 9076bdecfbSBarry Smith /*J 91c7393fdbSBarry Smith MatSolverPackage - String with the name of a PETSc matrix solver type. 929989ab13SBarry Smith 93c2b89b5dSBarry Smith For example: "petsc" indicates what PETSc provides, "superlu_dist" the parallel SuperLU_DIST package etc 949989ab13SBarry Smith 959989ab13SBarry Smith Level: beginner 969989ab13SBarry Smith 975c9eb25fSBarry Smith .seealso: MatGetFactor(), Mat, MatSetType(), MatType 9876bdecfbSBarry Smith J*/ 99c7393fdbSBarry Smith #define MatSolverPackage char* 1002692d6eeSBarry Smith #define MATSOLVERSUPERLU "superlu" 1012692d6eeSBarry Smith #define MATSOLVERSUPERLU_DIST "superlu_dist" 1022692d6eeSBarry Smith #define MATSOLVERUMFPACK "umfpack" 10320db9a53SJed Brown #define MATSOLVERCHOLMOD "cholmod" 1042692d6eeSBarry Smith #define MATSOLVERESSL "essl" 1052692d6eeSBarry Smith #define MATSOLVERLUSOL "lusol" 1062692d6eeSBarry Smith #define MATSOLVERMUMPS "mumps" 107d615f992SSatish Balay #define MATSOLVERMKL_PARDISO "mkl_pardiso" 108d305a81bSVasiliy Kozyrev #define MATSOLVERMKL_CPARDISO "mkl_cpardiso" 1092692d6eeSBarry Smith #define MATSOLVERPASTIX "pastix" 1102692d6eeSBarry Smith #define MATSOLVERMATLAB "matlab" 1112692d6eeSBarry Smith #define MATSOLVERPETSC "petsc" 1122692d6eeSBarry Smith #define MATSOLVERBAS "bas" 1139ae82921SPaul Mullowney #define MATSOLVERCUSPARSE "cusparse" 114aa5a9175SDahai Guo #define MATSOLVERBSTRM "bstrm" 115aa5a9175SDahai Guo #define MATSOLVERSBSTRM "sbstrm" 11615767789SHong Zhang #define MATSOLVERELEMENTAL "elemental" 11717f1a0eaSHong Zhang #define MATSOLVERCLIQUE "clique" 11869e15a41SShri Abhyankar #define MATSOLVERKLU "klu" 119c0cdd4a1SDahai Guo 120b24902e0SBarry Smith /*E 121b24902e0SBarry Smith MatFactorType - indicates what type of factorization is requested 122b24902e0SBarry Smith 123b24902e0SBarry Smith Level: beginner 124b24902e0SBarry Smith 125af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 126b24902e0SBarry Smith 127c7393fdbSBarry Smith .seealso: MatSolverPackage, MatGetFactor() 128b24902e0SBarry Smith E*/ 129599ef60dSHong Zhang typedef enum {MAT_FACTOR_NONE, MAT_FACTOR_LU, MAT_FACTOR_CHOLESKY, MAT_FACTOR_ILU, MAT_FACTOR_ICC,MAT_FACTOR_ILUDT} MatFactorType; 130014dd563SJed Brown PETSC_EXTERN const char *const MatFactorTypes[]; 131e92e720dSBarry Smith 132014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactor(Mat,const MatSolverPackage,MatFactorType,Mat*); 133014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactorAvailable(Mat,const MatSolverPackage,MatFactorType,PetscBool *); 134014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFactorGetSolverPackage(Mat,const MatSolverPackage*); 135014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactorType(Mat,MatFactorType*); 13618713533SBarry Smith PETSC_EXTERN PetscErrorCode MatSolverPackageRegister(const MatSolverPackage,const MatType,MatFactorType,PetscErrorCode(*)(Mat,MatFactorType,Mat*)); 13742c9c57cSBarry Smith PETSC_EXTERN PetscErrorCode MatSolverPackageGet(const MatSolverPackage,const MatType,MatFactorType,PetscBool*,PetscBool*,PetscErrorCode (**)(Mat,MatFactorType,Mat*)); 1389989ab13SBarry Smith 139c06d978dSMatthew Knepley /* Logging support */ 1400700a824SBarry Smith #define MAT_FILE_CLASSID 1211216 /* used to indicate matrices in binary files */ 141014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_CLASSID; 142335efc43SPeter Brune PETSC_EXTERN PetscClassId MAT_COLORING_CLASSID; 143014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_FDCOLORING_CLASSID; 144014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_TRANSPOSECOLORING_CLASSID; 145014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_PARTITIONING_CLASSID; 146014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_COARSEN_CLASSID; 147014dd563SJed Brown PETSC_EXTERN PetscClassId MAT_NULLSPACE_CLASSID; 148014dd563SJed Brown PETSC_EXTERN PetscClassId MATMFFD_CLASSID; 149c06d978dSMatthew Knepley 150ceb03754SKris Buschelman /*E 151ceb03754SKris Buschelman MatReuse - Indicates if matrices obtained from a previous call to MatGetSubMatrices() 152d6eff37eSBarry Smith or MatGetSubMatrix() are to be reused to store the new matrix values. For MatConvert() is used to indicate 153d6eff37eSBarry Smith that the input matrix is to be replaced with the converted matrix. 154ceb03754SKris Buschelman 155ceb03754SKris Buschelman Level: beginner 156ceb03754SKris Buschelman 157af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 158ceb03754SKris Buschelman 1590c451bc4SBarry Smith .seealso: MatGetSubMatrices(), MatGetSubMatrix(), MatDestroyMatrices(), MatConvert() 160ceb03754SKris Buschelman E*/ 161dfe085dbSJed Brown typedef enum {MAT_INITIAL_MATRIX,MAT_REUSE_MATRIX,MAT_IGNORE_MATRIX} MatReuse; 162ceb03754SKris Buschelman 1635494a064SHong Zhang /*E 1645494a064SHong Zhang MatGetSubMatrixOption - Indicates if matrices obtained from a call to MatGetSubMatrices() 165829201f2SHong Zhang include the matrix values. Currently it is only used by MatGetSeqNonzerostructure(). 1665494a064SHong Zhang 1675494a064SHong Zhang Level: beginner 1685494a064SHong Zhang 169829201f2SHong Zhang .seealso: MatGetSeqNonzerostructure() 1705494a064SHong Zhang E*/ 1715494a064SHong Zhang typedef enum {MAT_DO_NOT_GET_VALUES,MAT_GET_VALUES} MatGetSubMatrixOption; 1725494a064SHong Zhang 173607a6623SBarry Smith PETSC_EXTERN PetscErrorCode MatInitializePackage(void); 174c06d978dSMatthew Knepley 175014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreate(MPI_Comm,Mat*); 176014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetSizes(Mat,PetscInt,PetscInt,PetscInt,PetscInt); 17719fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatSetType(Mat,MatType); 178014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetFromOptions(Mat); 179685405a1SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatViewFromOptions(Mat A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);} 180bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatRegister(const char[],PetscErrorCode(*)(Mat)); 181014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRegisterBaseName(const char[],const char[],const char[]); 182014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetOptionsPrefix(Mat,const char[]); 183014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAppendOptionsPrefix(Mat,const char[]); 184014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOptionsPrefix(Mat,const char*[]); 185422a814eSBarry Smith PETSC_EXTERN PetscErrorCode MatSetErrorIfFPE(Mat,PetscBool); 186f69a0ea3SMatthew Knepley 187140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatList; 188140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatColoringList; 189140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatPartitioningList; 190140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatCoarsenList; 19128988994SBarry Smith 1923b224e63SBarry Smith /*E 193aa6c7ce3SBarry Smith MatStructure - Indicates if two matrices have the same nonzero structure 1943b224e63SBarry Smith 1953b224e63SBarry Smith Level: beginner 1963b224e63SBarry Smith 197af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 1983b224e63SBarry Smith 199aa6c7ce3SBarry Smith .seealso: MatCopy(), MatAXPY() 2003b224e63SBarry Smith E*/ 201aa6c7ce3SBarry Smith typedef enum {DIFFERENT_NONZERO_PATTERN,SUBSET_NONZERO_PATTERN,SAME_NONZERO_PATTERN} MatStructure; 2023b224e63SBarry Smith 203014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqDense(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*); 204014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateDense(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*); 205014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 206014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 207014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *); 208014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithSplitArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],PetscInt[],PetscInt[],PetscScalar[],Mat*); 2098d7a6e47SBarry Smith 210014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 211014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 212014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat*); 213d21a29f3SJed Brown 214014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAdj(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscInt[],Mat*); 215014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 216d21a29f3SJed Brown 217014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 218014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPISBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *); 21938f409ebSLisandro Dalcin PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 220014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 2214cce697cSJed Brown PETSC_EXTERN PetscErrorCode MatXAIJSetPreallocation(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscInt[],const PetscInt[]); 222dfb205c3SBarry Smith 223014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateShell(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,void *,Mat*); 224014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateNormal(Mat,Mat*); 225014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateLRC(Mat,Mat,Mat,Mat*); 226014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateIS(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,ISLocalToGlobalMapping,Mat*); 227014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 228014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 229c0cdd4a1SDahai Guo 230014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 231014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 232014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqSBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 233014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPISBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 234c0cdd4a1SDahai Guo 235014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateScatter(MPI_Comm,VecScatter,Mat*); 236014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScatterSetVecScatter(Mat,VecScatter); 237014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScatterGetVecScatter(Mat,VecScatter*); 238014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateBlockMat(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,Mat*); 239014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeAddMat(Mat,Mat); 240014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeMerge(Mat); 241014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateComposite(MPI_Comm,PetscInt,const Mat*,Mat*); 2426d7c1e57SBarry Smith typedef enum {MAT_COMPOSITE_ADDITIVE,MAT_COMPOSITE_MULTIPLICATIVE} MatCompositeType; 243014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCompositeSetType(Mat,MatCompositeType); 2446d7c1e57SBarry Smith 24519fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatCreateFFT(MPI_Comm,PetscInt,const PetscInt[],MatType,Mat*); 246014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqCUFFT(MPI_Comm,PetscInt,const PetscInt[],Mat*); 247dedccee8SHong Zhang 248014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateTranspose(Mat,Mat*); 249d0de2241SAndrew Spott PETSC_EXTERN PetscErrorCode MatCreateHermitianTranspose(Mat,Mat*); 250014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSubMatrix(Mat,IS,IS,Mat*); 251014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSubMatrixUpdate(Mat,Mat,IS,IS); 252014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateLocalRef(Mat,IS,IS,Mat*); 2531472f72bSBarry Smith 254014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPythonSetType(Mat,const char[]); 2551d6018f0SLisandro Dalcin 256014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetUp(Mat); 257014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroy(Mat*); 258e56f5c9eSBarry Smith PETSC_EXTERN PetscErrorCode MatGetNonzeroState(Mat,PetscObjectState*); 25921c89e3eSBarry Smith 260014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatConjugate(Mat); 261014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRealPart(Mat); 262014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatImaginaryPart(Mat); 263014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetDiagonalBlock(Mat,Mat*); 264014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetTrace(Mat,PetscScalar*); 265713ccfa9SJed Brown PETSC_EXTERN PetscErrorCode MatInvertBlockDiagonal(Mat,const PetscScalar **); 26699cafbc1SBarry Smith 2678ed539a5SBarry Smith /* ------------------------------------------------------------*/ 268014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 269014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 270014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesRow(Mat,PetscInt,const PetscScalar[]); 271014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesRowLocal(Mat,PetscInt,const PetscScalar[]); 272014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBatch(Mat,PetscInt,PetscInt,PetscInt[],const PetscScalar[]); 27373a71a0fSBarry Smith PETSC_EXTERN PetscErrorCode MatSetRandom(Mat,PetscRandom); 27484cb2905SBarry Smith 2752ef4de8bSBarry Smith /*S 2762ef4de8bSBarry Smith MatStencil - Data structure (C struct) for storing information about a single row or 277be479b30SJed Brown column of a matrix as indexed on an associated grid. 278be479b30SJed Brown 279be479b30SJed Brown Fortran usage is different, see MatSetValuesStencil() for details. 2802ef4de8bSBarry Smith 2812ef4de8bSBarry Smith Level: beginner 2822ef4de8bSBarry Smith 2832ef4de8bSBarry Smith Concepts: matrix; linear operator 2842ef4de8bSBarry Smith 285be479b30SJed Brown .seealso: MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockedStencil() 2862ef4de8bSBarry Smith S*/ 287435da068SBarry Smith typedef struct { 288c1ac3661SBarry Smith PetscInt k,j,i,c; 289435da068SBarry Smith } MatStencil; 2902ef4de8bSBarry Smith 291014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode); 292014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode); 293014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt); 294435da068SBarry Smith 295014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetColoring(Mat,ISColoring); 296014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesAdifor(Mat,PetscInt,void*); 2973a7fca6bSBarry Smith 298d91e6319SBarry Smith /*E 299d91e6319SBarry Smith MatAssemblyType - Indicates if the matrix is now to be used, or if you plan 300d91e6319SBarry Smith to continue to add values to it 301d91e6319SBarry Smith 302d91e6319SBarry Smith Level: beginner 303d91e6319SBarry Smith 304d91e6319SBarry Smith .seealso: MatAssemblyBegin(), MatAssemblyEnd() 305d91e6319SBarry Smith E*/ 3066d4a8577SBarry Smith typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType; 307014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssemblyBegin(Mat,MatAssemblyType); 308014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssemblyEnd(Mat,MatAssemblyType); 309014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAssembled(Mat,PetscBool *); 3104f9c727eSBarry Smith 3111d73ed98SBarry Smith 31230de9b25SBarry Smith 313d91e6319SBarry Smith /*E 314d91e6319SBarry Smith MatOption - Options that may be set for a matrix and its behavior or storage 315d91e6319SBarry Smith 316d91e6319SBarry Smith Level: beginner 317d91e6319SBarry Smith 318af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 319d91e6319SBarry Smith 320382164b0SBarry Smith Developer Notes: Entries that are negative need not be called collectively by all processes. 321382164b0SBarry Smith 322d91e6319SBarry Smith .seealso: MatSetOption() 323d91e6319SBarry Smith E*/ 32492d4d306SBarry Smith typedef enum {MAT_OPTION_MIN = -5, 32592d4d306SBarry Smith MAT_NEW_NONZERO_LOCATION_ERR = -4, 32692d4d306SBarry Smith MAT_UNUSED_NONZERO_LOCATION_ERR = -3, 32792d4d306SBarry Smith MAT_NEW_NONZERO_ALLOCATION_ERR = -2, 32892d4d306SBarry Smith MAT_ROW_ORIENTED = -1, 329382164b0SBarry Smith MAT_SYMMETRIC = 1, 330382164b0SBarry Smith MAT_STRUCTURALLY_SYMMETRIC = 2, 331382164b0SBarry Smith MAT_NEW_DIAGONALS = 3, 332382164b0SBarry Smith MAT_IGNORE_OFF_PROC_ENTRIES = 4, 333382164b0SBarry Smith MAT_USE_HASH_TABLE = 5, 334382164b0SBarry Smith MAT_KEEP_NONZERO_PATTERN = 6, 335382164b0SBarry Smith MAT_IGNORE_ZERO_ENTRIES = 7, 336382164b0SBarry Smith MAT_USE_INODES = 8, 337382164b0SBarry Smith MAT_HERMITIAN = 9, 338382164b0SBarry Smith MAT_SYMMETRY_ETERNAL = 10, 33911e456e1SBarry Smith MAT_DUMMY = 11, 340382164b0SBarry Smith MAT_IGNORE_LOWER_TRIANGULAR = 12, 341382164b0SBarry Smith MAT_ERROR_LOWER_TRIANGULAR = 13, 342382164b0SBarry Smith MAT_GETROW_UPPERTRIANGULAR = 14, 343382164b0SBarry Smith MAT_SPD = 15, 34492d4d306SBarry Smith MAT_NO_OFF_PROC_ZERO_ROWS = 16, 34592d4d306SBarry Smith MAT_NO_OFF_PROC_ENTRIES = 17, 34692d4d306SBarry Smith MAT_NEW_NONZERO_LOCATIONS = 18, 34792d4d306SBarry Smith MAT_OPTION_MAX = 19} MatOption; 348e057df02SPaul Mullowney 349014dd563SJed Brown PETSC_EXTERN const char *MatOptions[]; 350014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetOption(Mat,MatOption,PetscBool); 3517d68702bSBarry Smith PETSC_EXTERN PetscErrorCode MatGetOption(Mat,MatOption,PetscBool*); 35219fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatGetType(Mat,MatType*); 35384cb2905SBarry Smith 354014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]); 355014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 356014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 357014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowUpperTriangular(Mat); 358014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRowUpperTriangular(Mat); 359014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 360014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 361014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnVector(Mat,Vec,PetscInt); 3628c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJGetArray(Mat,PetscScalar *[]); 3638c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJRestoreArray(Mat,PetscScalar *[]); 36421e72a00SBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJGetMaxRowNonzeros(Mat,PetscInt*); 365bd04181cSBarry Smith PETSC_EXTERN PetscErrorCode MatSeqAIJSetValuesLocalFast(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 3668c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatDenseGetArray(Mat,PetscScalar *[]); 3678c778c55SBarry Smith PETSC_EXTERN PetscErrorCode MatDenseRestoreArray(Mat,PetscScalar *[]); 368014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBlockSize(Mat,PetscInt *); 369014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSize(Mat,PetscInt); 370014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBlockSizes(Mat,PetscInt *,PetscInt *); 371014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSizes(Mat,PetscInt,PetscInt); 37233d57670SJed Brown PETSC_EXTERN PetscErrorCode MatSetBlockSizesFromMats(Mat,Mat,Mat); 373014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNThreads(Mat,PetscInt); 374014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNThreads(Mat,PetscInt*); 3751620fd73SBarry Smith 376014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMult(Mat,Vec,Vec); 377014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultDiagonalBlock(Mat,Vec,Vec); 378014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultAdd(Mat,Vec,Vec,Vec); 379014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTranspose(Mat,Vec,Vec); 380014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultHermitianTranspose(Mat,Vec,Vec); 381014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsTranspose(Mat,Mat,PetscReal,PetscBool *); 382014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitianTranspose(Mat,Mat,PetscReal,PetscBool *); 383014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeAdd(Mat,Vec,Vec,Vec); 384014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultHermitianTransposeAdd(Mat,Vec,Vec,Vec); 385014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultConstrained(Mat,Vec,Vec); 386014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeConstrained(Mat,Vec,Vec); 387014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatSolve(Mat,Mat,Mat); 388f9426fe0SMark Adams PETSC_EXTERN PetscErrorCode MatResidual(Mat,Vec,Vec,Vec); 389f5edf698SHong Zhang 390d91e6319SBarry Smith /*E 391d91e6319SBarry Smith MatDuplicateOption - Indicates if a duplicated sparse matrix should have 392d91e6319SBarry Smith its numerical values copied over or just its nonzero structure. 393d91e6319SBarry Smith 394d91e6319SBarry Smith Level: beginner 395d91e6319SBarry Smith 396af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 397d91e6319SBarry Smith 39870dcbbb9SBarry Smith $ MAT_SHARE_NONZERO_PATTERN - the i and j arrays in the new matrix will be shared with the original matrix 39970dcbbb9SBarry Smith $ this also triggers the MAT_DO_NOT_COPY_VALUES option. This is used when you 40070dcbbb9SBarry Smith $ have several matrices with the same nonzero pattern. 40170dcbbb9SBarry Smith 402d91e6319SBarry Smith .seealso: MatDuplicate() 403d91e6319SBarry Smith E*/ 40470dcbbb9SBarry Smith typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES,MAT_SHARE_NONZERO_PATTERN} MatDuplicateOption; 4052e8a6d31SBarry Smith 40619fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatConvert(Mat,MatType,MatReuse,Mat*); 407014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDuplicate(Mat,MatDuplicateOption,Mat*); 40894a9d846SBarry Smith 409cb5b572fSBarry Smith 410014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCopy(Mat,Mat,MatStructure); 411014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatView(Mat,PetscViewer); 412014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsSymmetric(Mat,PetscReal,PetscBool *); 413014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsStructurallySymmetric(Mat,PetscBool *); 414014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitian(Mat,PetscReal,PetscBool *); 415014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsSymmetricKnown(Mat,PetscBool *,PetscBool *); 416014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIsHermitianKnown(Mat,PetscBool *,PetscBool *); 417014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMissingDiagonal(Mat,PetscBool *,PetscInt *); 418014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLoad(Mat, PetscViewer); 4197b80b807SBarry Smith 4201a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool *); 4211a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool *); 4221a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool *); 4231a83f524SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool *); 424d4fbbf0eSBarry Smith 425d91e6319SBarry Smith /*S 426d91e6319SBarry Smith MatInfo - Context of matrix information, used with MatGetInfo() 427d91e6319SBarry Smith 428d91e6319SBarry Smith In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE 429d91e6319SBarry Smith 430d91e6319SBarry Smith Level: intermediate 431d91e6319SBarry Smith 432d91e6319SBarry Smith Concepts: matrix^nonzero information 433d91e6319SBarry Smith 434d9274352SBarry Smith .seealso: MatGetInfo(), MatInfoType 435d91e6319SBarry Smith S*/ 4364e220ebcSLois Curfman McInnes typedef struct { 437b0a32e0cSBarry Smith PetscLogDouble block_size; /* block size */ 438b0a32e0cSBarry Smith PetscLogDouble nz_allocated,nz_used,nz_unneeded; /* number of nonzeros */ 439b0a32e0cSBarry Smith PetscLogDouble memory; /* memory allocated */ 440b0a32e0cSBarry Smith PetscLogDouble assemblies; /* number of matrix assemblies called */ 441b0a32e0cSBarry Smith PetscLogDouble mallocs; /* number of mallocs during MatSetValues() */ 442b0a32e0cSBarry Smith PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */ 443b0a32e0cSBarry Smith PetscLogDouble factor_mallocs; /* number of mallocs during factorization */ 4444e220ebcSLois Curfman McInnes } MatInfo; 4454e220ebcSLois Curfman McInnes 446d9274352SBarry Smith /*E 447d9274352SBarry Smith MatInfoType - Indicates if you want information about the local part of the matrix, 448d9274352SBarry Smith the entire parallel matrix or the maximum over all the local parts. 449d9274352SBarry Smith 450d9274352SBarry Smith Level: beginner 451d9274352SBarry Smith 452af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 453d9274352SBarry Smith 454d9274352SBarry Smith .seealso: MatGetInfo(), MatInfo 455d9274352SBarry Smith E*/ 4567b80b807SBarry Smith typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType; 457014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetInfo(Mat,MatInfoType,MatInfo*); 458014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetDiagonal(Mat,Vec); 459014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMax(Mat,Vec,PetscInt[]); 460014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMin(Mat,Vec,PetscInt[]); 461014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMaxAbs(Mat,Vec,PetscInt[]); 462014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowMinAbs(Mat,Vec,PetscInt[]); 463014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetRowSum(Mat,Vec); 464014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTranspose(Mat,MatReuse,Mat*); 465014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatHermitianTranspose(Mat,MatReuse,Mat*); 466014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPermute(Mat,IS,IS,Mat *); 467014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalScale(Mat,Vec,Vec); 468014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalSet(Mat,Vec,InsertMode); 469014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatEqual(Mat,Mat,PetscBool *); 470014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultEqual(Mat,Mat,PetscInt,PetscBool *); 471014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultAddEqual(Mat,Mat,PetscInt,PetscBool *); 472014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeEqual(Mat,Mat,PetscInt,PetscBool *); 473014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscBool *); 4747b80b807SBarry Smith 475014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNorm(Mat,NormType,PetscReal *); 476014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetColumnNorms(Mat,NormType,PetscReal *); 477014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroEntries(Mat); 478014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 479014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsIS(Mat,IS,PetscScalar,Vec,Vec); 480014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsStencil(Mat,PetscInt,const MatStencil [],PetscScalar,Vec,Vec); 481014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsStencil(Mat,PetscInt,const MatStencil[],PetscScalar,Vec,Vec); 482014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumns(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 483014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsIS(Mat,IS,PetscScalar,Vec,Vec); 4847b80b807SBarry Smith 485014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSize(Mat,PetscInt*,PetscInt*); 486014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalSize(Mat,PetscInt*,PetscInt*); 487014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRange(Mat,PetscInt*,PetscInt*); 488014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRanges(Mat,const PetscInt**); 489014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRangeColumn(Mat,PetscInt*,PetscInt*); 490014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipRangesColumn(Mat,const PetscInt**); 491566876eaSJed Brown PETSC_EXTERN PetscErrorCode MatGetOwnershipIS(Mat,IS*,IS*); 4925ef9f2a5SBarry Smith 493014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]); 49453dd7562SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatGetSubMatricesMPI(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]); 495014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroyMatrices(PetscInt,Mat *[]); 496014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSubMatrix(Mat,IS,IS,MatReuse,Mat *); 497014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalSubMatrix(Mat,IS,IS,Mat*); 498014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestoreLocalSubMatrix(Mat,IS,IS,Mat*); 499014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetSeqNonzeroStructure(Mat,Mat*); 500014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDestroySeqNonzeroStructure(Mat*); 5017b80b807SBarry Smith 502014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJ(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*); 503014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJSymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*); 504014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJNumeric(Mat,Mat); 505014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMat(Mat,MatReuse,Mat*); 506014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*); 507014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,Mat*); 508014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetGhosts(Mat, PetscInt *,const PetscInt *[]); 5098efafbd8SBarry Smith 510014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt); 511*aa1e27eaSFande Kong PETSC_EXTERN PetscErrorCode MatIncreaseOverlapSplit(Mat mat,PetscInt n,IS is[],PetscInt ov); 512d2b2a242SBarry Smith PETSC_EXTERN PetscErrorCode MatMPIAIJSetUseScalableIncreaseOverlap(Mat,PetscBool); 5137b80b807SBarry Smith 514014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*); 515014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMultSymbolic(Mat,Mat,PetscReal,Mat*); 516014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatMultNumeric(Mat,Mat,Mat); 51722440eb1SKris Buschelman 5187bab7c10SHong Zhang PETSC_EXTERN PetscErrorCode MatMatMatMult(Mat,Mat,Mat,MatReuse,PetscReal,Mat*); 5197bab7c10SHong Zhang PETSC_EXTERN PetscErrorCode MatMatMatMultSymbolic(Mat,Mat,Mat,PetscReal,Mat*); 5207bab7c10SHong Zhang PETSC_EXTERN PetscErrorCode MatMatMatMultNumeric(Mat,Mat,Mat,Mat); 5217bab7c10SHong Zhang 522014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*); 523014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAPSymbolic(Mat,Mat,PetscReal,Mat*); 524014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPtAPNumeric(Mat,Mat,Mat); 525014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARt(Mat,Mat,MatReuse,PetscReal,Mat*); 526014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARtSymbolic(Mat,Mat,PetscReal,Mat*); 527014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRARtNumeric(Mat,Mat,Mat); 52822440eb1SKris Buschelman 529014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeMatMult(Mat,Mat,MatReuse,PetscReal,Mat*); 530014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposetMatMultSymbolic(Mat,Mat,PetscReal,Mat*); 531014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposetMatMultNumeric(Mat,Mat,Mat); 532014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMult(Mat,Mat,MatReuse,PetscReal,Mat*); 533014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMultSymbolic(Mat,Mat,PetscReal,Mat*); 534014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMatTransposeMultNumeric(Mat,Mat,Mat); 535bc011b1eSHong Zhang 536014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAXPY(Mat,PetscScalar,Mat,MatStructure); 537014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatAYPX(Mat,PetscScalar,Mat,MatStructure); 5381c741599SHong Zhang 539014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatScale(Mat,PetscScalar); 540014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShift(Mat,PetscScalar); 5417b80b807SBarry Smith 542014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping); 543014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*); 544a93ff8c4SPeter Brune PETSC_EXTERN PetscErrorCode MatGetLayouts(Mat,PetscLayout*,PetscLayout*); 545014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 546014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsLocalIS(Mat,IS,PetscScalar,Vec,Vec); 547014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 548014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocalIS(Mat,IS,PetscScalar,Vec,Vec); 549014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 550014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 551052efed2SBarry Smith 552014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStashSetInitialSize(Mat,PetscInt,PetscInt); 553014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*); 55490f02eecSBarry Smith 555014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInterpolate(Mat,Vec,Vec); 556014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInterpolateAdd(Mat,Vec,Vec,Vec); 557014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRestrict(Mat,Vec,Vec); 5582a7a6963SBarry Smith PETSC_EXTERN PetscErrorCode MatCreateVecs(Mat,Vec*,Vec*); 55984cc2a5aSBarry Smith PETSC_DEPRECATED("Use MatCreateVecs()") PETSC_STATIC_INLINE PetscErrorCode MatGetVecs(Mat mat,Vec *x,Vec *y) {return MatCreateVecs(mat,x,y);} 56053cd1579SHong Zhang PETSC_EXTERN PetscErrorCode MatCreateRedundantMatrix(Mat,PetscInt,MPI_Comm,MatReuse,Mat*); 561014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetMultiProcBlock(Mat,MPI_Comm,MatReuse,Mat*); 562014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFindZeroDiagonals(Mat,IS*); 5633a062f41SBarry Smith PETSC_EXTERN PetscErrorCode MatFindOffBlockDiagonalEntries(Mat,IS*); 5649c8f2541SHong Zhang PETSC_EXTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat(MPI_Comm,Mat,PetscInt,MatReuse,Mat*); 565bd481603SBarry Smith 566bd481603SBarry Smith /*MC 5671620fd73SBarry Smith MatSetValue - Set a single entry into a matrix. 5681620fd73SBarry Smith 5691620fd73SBarry Smith Not collective 5701620fd73SBarry Smith 571a9834a7dSSatish Balay Synopsis: 572a9834a7dSSatish Balay #include <petscmat.h> 573a9834a7dSSatish Balay PetscErrorCode MatSetValue(Mat m,PetscInt row,PetscInt col,PetscScalar value,InsertMode mode) 574a9834a7dSSatish Balay 5751620fd73SBarry Smith Input Parameters: 5761620fd73SBarry Smith + m - the matrix 5771620fd73SBarry Smith . row - the row location of the entry 5781620fd73SBarry Smith . col - the column location of the entry 5791620fd73SBarry Smith . value - the value to insert 5801620fd73SBarry Smith - mode - either INSERT_VALUES or ADD_VALUES 5811620fd73SBarry Smith 5821620fd73SBarry Smith Notes: 5831620fd73SBarry Smith For efficiency one should use MatSetValues() and set several or many 5841620fd73SBarry Smith values simultaneously if possible. 5851620fd73SBarry Smith 5861620fd73SBarry Smith Level: beginner 5871620fd73SBarry Smith 5881620fd73SBarry Smith .seealso: MatSetValues(), MatSetValueLocal() 5891620fd73SBarry Smith M*/ 5901620fd73SBarry 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);} 5911620fd73SBarry Smith 5921620fd73SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatGetValue(Mat v,PetscInt i,PetscInt j,PetscScalar *va) {return MatGetValues(v,1,&i,1,&j,va);} 5931620fd73SBarry Smith 5941620fd73SBarry 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);} 5951620fd73SBarry Smith 5961620fd73SBarry Smith /*MC 5970d2a7ff2SSatish Balay MatPreallocateInitialize - Begins the block of code that will count the number of nonzeros per 598bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 599bd481603SBarry Smith 600bd481603SBarry Smith Synopsis: 601aaa7dc30SBarry Smith #include <petscmat.h> 602c1ac3661SBarry Smith PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz) 603bd481603SBarry Smith 604bd481603SBarry Smith Collective on MPI_Comm 605bd481603SBarry Smith 606bd481603SBarry Smith Input Parameters: 607bd481603SBarry Smith + comm - the communicator that will share the eventually allocated matrix 608859fcb39SBarry Smith . nrows - the number of LOCAL rows in the matrix 609859fcb39SBarry Smith - ncols - the number of LOCAL columns in the matrix 610bd481603SBarry Smith 611bd481603SBarry Smith Output Parameters: 612bd481603SBarry Smith + dnz - the array that will be passed to the matrix preallocation routines 613bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 614bd481603SBarry Smith 615bd481603SBarry Smith Level: intermediate 616bd481603SBarry Smith 617bd481603SBarry Smith Notes: 618a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 619bd481603SBarry Smith 6201d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 621bd481603SBarry Smith 622bd481603SBarry Smith Use MatPreallocateInitializeSymmetric() for symmetric matrices (MPISBAIJ matrices) 623bd481603SBarry Smith 6241620fd73SBarry Smith This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize(). 6251620fd73SBarry Smith 626bd481603SBarry Smith Concepts: preallocation^Matrix 627bd481603SBarry Smith 628d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(), 629d6e23781SBarry Smith MatPreallocateInitializeSymmetric(), MatPreallocateSymmetricSetLocalBlock() 630bd481603SBarry Smith M*/ 631c4f061fbSSatish Balay #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \ 6327c922b88SBarry Smith { \ 633a89bc211SBarry Smith PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__start,__end; \ 6341795a4d1SJed Brown _4_ierr = PetscCalloc2(__nrows,&dnz,__nrows,&onz);CHKERRQ(_4_ierr); \ 6351795a4d1SJed Brown __start = 0; __end = __start; \ 636c1ac3661SBarry Smith _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ctmp;\ 637a89bc211SBarry Smith _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows; 6387c922b88SBarry Smith 639bd481603SBarry Smith /*MC 640bd481603SBarry Smith MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be 641bd481603SBarry Smith inserted using a local number of the rows and columns 642bd481603SBarry Smith 643bd481603SBarry Smith Synopsis: 644aaa7dc30SBarry Smith #include <petscmat.h> 645c1ac3661SBarry Smith PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 646bd481603SBarry Smith 647bd481603SBarry Smith Not Collective 648bd481603SBarry Smith 649bd481603SBarry Smith Input Parameters: 650784ac674SJed Brown + map - the row mapping from local numbering to global numbering 651bd481603SBarry Smith . nrows - the number of rows indicated 6521d73ed98SBarry Smith . rows - the indices of the rows 653784ac674SJed Brown . cmap - the column mapping from local to global numbering 654bd481603SBarry Smith . ncols - the number of columns in the matrix 655bd481603SBarry Smith . cols - the columns indicated 656bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 657bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 658bd481603SBarry Smith 659bd481603SBarry Smith Level: intermediate 660bd481603SBarry Smith 661bd481603SBarry Smith Notes: 662a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 663bd481603SBarry Smith 6641d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 665bd481603SBarry Smith 666bd481603SBarry Smith Concepts: preallocation^Matrix 667bd481603SBarry Smith 668d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(), 669d6e23781SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock() 670bd481603SBarry Smith M*/ 671784ac674SJed Brown #define MatPreallocateSetLocal(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \ 672c4f061fbSSatish Balay {\ 673c1ac3661SBarry Smith PetscInt __l;\ 674784ac674SJed Brown _4_ierr = ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\ 675784ac674SJed Brown _4_ierr = ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\ 676c4f061fbSSatish Balay for (__l=0;__l<nrows;__l++) {\ 677ef66eb69SBarry Smith _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 678c4f061fbSSatish Balay }\ 679c4f061fbSSatish Balay } 680c4f061fbSSatish Balay 681bd481603SBarry Smith /*MC 682d6e23781SBarry Smith MatPreallocateSetLocalBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be 683bd481603SBarry Smith inserted using a local number of the rows and columns 684bd481603SBarry Smith 685bd481603SBarry Smith Synopsis: 686aaa7dc30SBarry Smith #include <petscmat.h> 687d6e23781SBarry Smith PetscErrorCode MatPreallocateSetLocalBlock(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 688d6e23781SBarry Smith 689d6e23781SBarry Smith Not Collective 690d6e23781SBarry Smith 691d6e23781SBarry Smith Input Parameters: 692d6e23781SBarry Smith + map - the row mapping from local numbering to global numbering 693d6e23781SBarry Smith . nrows - the number of rows indicated 694d6e23781SBarry Smith . rows - the indices of the rows 695d6e23781SBarry Smith . cmap - the column mapping from local to global numbering 696d6e23781SBarry Smith . ncols - the number of columns in the matrix 697d6e23781SBarry Smith . cols - the columns indicated 698d6e23781SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 699d6e23781SBarry Smith - ozn - the other array passed to the matrix preallocation routines 700d6e23781SBarry Smith 701d6e23781SBarry Smith Level: intermediate 702d6e23781SBarry Smith 703d6e23781SBarry Smith Notes: 704d6e23781SBarry Smith See Users-Manual: ch_performance for more details. 705d6e23781SBarry Smith 706d6e23781SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 707d6e23781SBarry Smith 708d6e23781SBarry Smith Concepts: preallocation^Matrix 709d6e23781SBarry Smith 710d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(), 711d6e23781SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock() 712d6e23781SBarry Smith M*/ 713d6e23781SBarry Smith #define MatPreallocateSetLocalBlock(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \ 714d6e23781SBarry Smith {\ 715d6e23781SBarry Smith PetscInt __l;\ 716d6e23781SBarry Smith _4_ierr = ISLocalToGlobalMappingApplyBlock(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\ 717d6e23781SBarry Smith _4_ierr = ISLocalToGlobalMappingApplyBlock(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\ 718d6e23781SBarry Smith for (__l=0;__l<nrows;__l++) {\ 719d6e23781SBarry Smith _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 720d6e23781SBarry Smith }\ 721d6e23781SBarry Smith } 722d6e23781SBarry Smith 723d6e23781SBarry Smith /*MC 724d6e23781SBarry Smith MatPreallocateSymmetricSetLocalBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be 725d6e23781SBarry Smith inserted using a local number of the rows and columns 726d6e23781SBarry Smith 727d6e23781SBarry Smith Synopsis: 728d6e23781SBarry Smith #include <petscmat.h> 729d6e23781SBarry Smith PetscErrorCode MatPreallocateSymmetricSetLocalBlock(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 730bd481603SBarry Smith 731bd481603SBarry Smith Not Collective 732bd481603SBarry Smith 733bd481603SBarry Smith Input Parameters: 734bd481603SBarry Smith + map - the mapping between local numbering and global numbering 735bd481603SBarry Smith . nrows - the number of rows indicated 7361d73ed98SBarry Smith . rows - the indices of the rows 737bd481603SBarry Smith . ncols - the number of columns in the matrix 738bd481603SBarry Smith . cols - the columns indicated 739bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 740bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 741bd481603SBarry Smith 742bd481603SBarry Smith Level: intermediate 743bd481603SBarry Smith 744bd481603SBarry Smith Notes: 745a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 746bd481603SBarry Smith 747bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 748bd481603SBarry Smith 749bd481603SBarry Smith Concepts: preallocation^Matrix 750bd481603SBarry Smith 751d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateInitialize(), 752d6e23781SBarry Smith MatPreallocateInitialize(), MatPreallocateSetLocal() 753bd481603SBarry Smith M*/ 754d6e23781SBarry Smith #define MatPreallocateSymmetricSetLocalBlock(map,nrows,rows,ncols,cols,dnz,onz) 0;\ 755d3d32019SBarry Smith {\ 756c1ac3661SBarry Smith PetscInt __l;\ 757d6e23781SBarry Smith _4_ierr = ISLocalToGlobalMappingApplyBlock(map,nrows,rows,rows);CHKERRQ(_4_ierr);\ 758d6e23781SBarry Smith _4_ierr = ISLocalToGlobalMappingApplyBlock(map,ncols,cols,cols);CHKERRQ(_4_ierr);\ 759d3d32019SBarry Smith for (__l=0;__l<nrows;__l++) {\ 760d6e23781SBarry Smith _4_ierr = MatPreallocateSymmetricSetBlock((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 761d3d32019SBarry Smith }\ 762d3d32019SBarry Smith } 763bd481603SBarry Smith /*MC 764bd481603SBarry Smith MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be 765bd481603SBarry Smith inserted using a local number of the rows and columns 766bd481603SBarry Smith 767bd481603SBarry Smith Synopsis: 768aaa7dc30SBarry Smith #include <petscmat.h> 769c1ac3661SBarry Smith PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 770bd481603SBarry Smith 771bd481603SBarry Smith Not Collective 772bd481603SBarry Smith 773bd481603SBarry Smith Input Parameters: 77464075487SBarry Smith + row - the row 775bd481603SBarry Smith . ncols - the number of columns in the matrix 776a50f8bf6SHong Zhang - cols - the columns indicated 777a50f8bf6SHong Zhang 778a50f8bf6SHong Zhang Output Parameters: 779a50f8bf6SHong Zhang + dnz - the array that will be passed to the matrix preallocation routines 780bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 781bd481603SBarry Smith 782bd481603SBarry Smith Level: intermediate 783bd481603SBarry Smith 784bd481603SBarry Smith Notes: 785a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 786bd481603SBarry Smith 787bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 788bd481603SBarry Smith 7891620fd73SBarry Smith This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize(). 7901620fd73SBarry Smith 791bd481603SBarry Smith Concepts: preallocation^Matrix 792bd481603SBarry Smith 793d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateInitialize(), 794d6e23781SBarry Smith MatPreallocateInitialize(), MatPreallocateSetLocal() 795bd481603SBarry Smith M*/ 796c4f061fbSSatish Balay #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\ 797c1ac3661SBarry Smith { PetscInt __i; \ 798e32f2f54SBarry 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);\ 799e32f2f54SBarry 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);\ 8007c922b88SBarry Smith for (__i=0; __i<nc; __i++) {\ 80164075487SBarry Smith if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \ 8027cc688d7SBarry Smith else dnz[row - __rstart]++;\ 8037c922b88SBarry Smith }\ 8047c922b88SBarry Smith } 8057c922b88SBarry Smith 806bd481603SBarry Smith /*MC 807d6e23781SBarry Smith MatPreallocateSymmetricSetBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be 808bd481603SBarry Smith inserted using a local number of the rows and columns 809bd481603SBarry Smith 810bd481603SBarry Smith Synopsis: 811aaa7dc30SBarry Smith #include <petscmat.h> 812d6e23781SBarry Smith PetscErrorCode MatPreallocateSymmetricSetBlock(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 813bd481603SBarry Smith 814bd481603SBarry Smith Not Collective 815bd481603SBarry Smith 816bd481603SBarry Smith Input Parameters: 817bd481603SBarry Smith + nrows - the number of rows indicated 8181d73ed98SBarry Smith . rows - the indices of the rows 819bd481603SBarry Smith . ncols - the number of columns in the matrix 820bd481603SBarry Smith . cols - the columns indicated 821bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 822bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 823bd481603SBarry Smith 824bd481603SBarry Smith Level: intermediate 825bd481603SBarry Smith 826bd481603SBarry Smith Notes: 827a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 828bd481603SBarry Smith 829bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 830bd481603SBarry Smith 8311620fd73SBarry Smith This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize(). 8321620fd73SBarry Smith 833bd481603SBarry Smith Concepts: preallocation^Matrix 834bd481603SBarry Smith 835d6e23781SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateInitialize(), 836d6e23781SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocal() 837bd481603SBarry Smith M*/ 838d6e23781SBarry Smith #define MatPreallocateSymmetricSetBlock(row,nc,cols,dnz,onz) 0;\ 839c1ac3661SBarry Smith { PetscInt __i; \ 840d3d32019SBarry Smith for (__i=0; __i<nc; __i++) {\ 841d3d32019SBarry Smith if (cols[__i] >= __end) onz[row - __rstart]++; \ 842d3d32019SBarry Smith else if (cols[__i] >= row) dnz[row - __rstart]++;\ 843d3d32019SBarry Smith }\ 844d3d32019SBarry Smith } 845d3d32019SBarry Smith 846bd481603SBarry Smith /*MC 84716371a99SBarry Smith MatPreallocateLocation - An alternative to MatPreallocationSet() that puts the nonzero locations into the matrix if it exists 84816371a99SBarry Smith 84916371a99SBarry Smith Synopsis: 850aaa7dc30SBarry Smith #include <petscmat.h> 85116371a99SBarry Smith PetscErrorCode MatPreallocateLocations(Mat A,PetscInt row,PetscInt ncols,PetscInt *cols,PetscInt *dnz,PetscInt *onz) 85216371a99SBarry Smith 85316371a99SBarry Smith Not Collective 85416371a99SBarry Smith 85516371a99SBarry Smith Input Parameters: 85616371a99SBarry Smith . A - matrix 85716371a99SBarry Smith . row - row where values exist (must be local to this process) 85816371a99SBarry Smith . ncols - number of columns 85916371a99SBarry Smith . cols - columns with nonzeros 86016371a99SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 86116371a99SBarry Smith - ozn - the other array passed to the matrix preallocation routines 86216371a99SBarry Smith 86316371a99SBarry Smith Level: intermediate 86416371a99SBarry Smith 86516371a99SBarry Smith Notes: 866a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 86716371a99SBarry Smith 86816371a99SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 86916371a99SBarry Smith 87016371a99SBarry Smith This is a MACRO not a function because it uses a bunch of variables private to the MatPreallocation.... routines. 87116371a99SBarry Smith 87216371a99SBarry Smith Concepts: preallocation^Matrix 87316371a99SBarry Smith 874d6e23781SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(), 875d6e23781SBarry Smith MatPreallocateSymmetricSetLocalBlock() 87616371a99SBarry Smith M*/ 8770298fd71SBarry 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);} 87816371a99SBarry Smith 87916371a99SBarry Smith 88016371a99SBarry Smith /*MC 8810d2a7ff2SSatish Balay MatPreallocateFinalize - Ends the block of code that will count the number of nonzeros per 882bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 883bd481603SBarry Smith 884bd481603SBarry Smith Synopsis: 885aaa7dc30SBarry Smith #include <petscmat.h> 886c1ac3661SBarry Smith PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz) 887bd481603SBarry Smith 888bd481603SBarry Smith Collective on MPI_Comm 889bd481603SBarry Smith 890bd481603SBarry Smith Input Parameters: 89116371a99SBarry Smith + dnz - the array that was be passed to the matrix preallocation routines 892bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 893bd481603SBarry Smith 894bd481603SBarry Smith Level: intermediate 895bd481603SBarry Smith 896bd481603SBarry Smith Notes: 897a7f22e61SSatish Balay See Users-Manual: ch_performance for more details. 898bd481603SBarry Smith 899bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 900bd481603SBarry Smith 9011620fd73SBarry Smith This is a MACRO not a function because it closes the { started in MatPreallocateInitialize(). 9021620fd73SBarry Smith 903bd481603SBarry Smith Concepts: preallocation^Matrix 904bd481603SBarry Smith 905d6e23781SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(), 906d6e23781SBarry Smith MatPreallocateSymmetricSetLocalBlock() 907bd481603SBarry Smith M*/ 908a89bc211SBarry Smith #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree2(dnz,onz);CHKERRQ(_4_ierr);} 9097c922b88SBarry Smith 9107b80b807SBarry Smith /* Routines unique to particular data structures */ 911014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellGetContext(Mat,void *); 912698d4c6aSKris Buschelman 913014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInodeAdjustForInodes(Mat,IS*,IS*); 914014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *); 915698d4c6aSKris Buschelman 916014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetColumnIndices(Mat,PetscInt[]); 917014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetColumnIndices(Mat,PetscInt[]); 918014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 919014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 920014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 921014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJFromTriple(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*,PetscInt,PetscBool); 9227b80b807SBarry Smith 923a96a251dSBarry Smith #define MAT_SKIP_ALLOCATION -4 924a96a251dSBarry Smith 925014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]); 926014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]); 927014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]); 928273d9f13SBarry Smith 929014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 930014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 931014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 932014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []); 933014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 934014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]); 935014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 936014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]); 937014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIDenseSetPreallocation(Mat,PetscScalar[]); 938014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqDenseSetPreallocation(Mat,PetscScalar[]); 9399230625dSJed Brown PETSC_EXTERN PetscErrorCode MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,const PetscInt*[]); 9409230625dSJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,const PetscInt*[]); 941014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIAdjCreateNonemptySubcommMat(Mat,Mat*); 942273d9f13SBarry Smith 9432e1947a5SStefano Zampini PETSC_EXTERN PetscErrorCode MatISSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 944014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqDenseSetLDA(Mat,PetscInt); 945014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDenseGetLocalMatrix(Mat,Mat*); 9461b807ce4Svictorle 947014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatStoreValues(Mat); 948014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatRetrieveValues(Mat); 9492e8a6d31SBarry Smith 950014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDAADSetCtx(Mat,void*); 9513a7fca6bSBarry Smith 952014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFindNonzeroRows(Mat,IS*); 9537b80b807SBarry Smith /* 9547b80b807SBarry Smith These routines are not usually accessed directly, rather solving is 95594b7f48cSBarry Smith done through the KSP and PC interfaces. 9567b80b807SBarry Smith */ 9577b80b807SBarry Smith 95876bdecfbSBarry Smith /*J 9598f6c3df8SBarry Smith MatOrderingType - String with the name of a PETSc matrix ordering 960d9274352SBarry Smith 961d9274352SBarry Smith Level: beginner 962d9274352SBarry Smith 963d9274352SBarry Smith .seealso: MatGetOrdering() 96476bdecfbSBarry Smith J*/ 96519fd82e9SBarry Smith typedef const char* MatOrderingType; 9662692d6eeSBarry Smith #define MATORDERINGNATURAL "natural" 9672692d6eeSBarry Smith #define MATORDERINGND "nd" 9682692d6eeSBarry Smith #define MATORDERING1WD "1wd" 9692692d6eeSBarry Smith #define MATORDERINGRCM "rcm" 9702692d6eeSBarry Smith #define MATORDERINGQMD "qmd" 9712692d6eeSBarry Smith #define MATORDERINGROWLENGTH "rowlength" 972510184a7SJed Brown #define MATORDERINGWBM "wbm" 973fc1bef75SJed Brown #define MATORDERINGSPECTRAL "spectral" 974312e372aSBarry Smith #define MATORDERINGAMD "amd" /* only works if UMFPACK is installed with PETSc */ 975b12f92e5SBarry Smith 97619fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatGetOrdering(Mat,MatOrderingType,IS*,IS*); 977140e18c1SBarry Smith PETSC_EXTERN PetscErrorCode MatGetOrderingList(PetscFunctionList*); 978bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatOrderingRegister(const char[],PetscErrorCode(*)(Mat,MatOrderingType,IS*,IS*)); 979140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList MatOrderingList; 980d4fbbf0eSBarry Smith 981014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS); 982fc1bef75SJed Brown PETSC_EXTERN PetscErrorCode MatCreateLaplacian(Mat,PetscReal,PetscBool,Mat*); 983a2ce50c7SBarry Smith 984d91e6319SBarry Smith /*S 985d90ac83dSHong Zhang MatFactorShiftType - Numeric Shift. 986d90ac83dSHong Zhang 987d90ac83dSHong Zhang Level: beginner 988d90ac83dSHong Zhang 989d90ac83dSHong Zhang S*/ 990d90ac83dSHong Zhang typedef enum {MAT_SHIFT_NONE,MAT_SHIFT_NONZERO,MAT_SHIFT_POSITIVE_DEFINITE,MAT_SHIFT_INBLOCKS} MatFactorShiftType; 9916a6fc655SJed Brown PETSC_EXTERN const char *const MatFactorShiftTypes[]; 9925e9742b9SJed Brown PETSC_EXTERN const char *const MatFactorShiftTypesDetail[]; 993d90ac83dSHong Zhang 994d90ac83dSHong Zhang /*S 995422a814eSBarry Smith MatFactorInfo - Data passed into the matrix factorization routines, and information about the resulting factorization 996b00f7748SHong Zhang 99761cad860SBarry Smith In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE, that is use 99861cad860SBarry Smith $ MatFactorInfo info(MAT_FACTORINFO_SIZE) 999b00f7748SHong Zhang 100015e8a5b3SHong Zhang Notes: These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC. 1001b00f7748SHong Zhang 100261cad860SBarry Smith You can use MatFactorInfoInitialize() to set default values. 100361cad860SBarry Smith 1004b00f7748SHong Zhang Level: developer 1005b00f7748SHong Zhang 1006d7d82daaSBarry Smith .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(), 1007d7d82daaSBarry Smith MatFactorInfoInitialize() 1008b00f7748SHong Zhang 1009b00f7748SHong Zhang S*/ 1010b00f7748SHong Zhang typedef struct { 101115e8a5b3SHong Zhang PetscReal diagonal_fill; /* force diagonal to fill in if initially not filled */ 101285317021SBarry Smith PetscReal usedt; 101315e8a5b3SHong Zhang PetscReal dt; /* drop tolerance */ 1014b00f7748SHong Zhang PetscReal dtcol; /* tolerance for pivoting */ 101515e8a5b3SHong Zhang PetscReal dtcount; /* maximum nonzeros to be allowed per row */ 101667e28bfeSBarry Smith PetscReal fill; /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */ 1017348344bbSBarry Smith PetscReal levels; /* ICC/ILU(levels) */ 1018bcd9e38bSBarry Smith PetscReal pivotinblocks; /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0 1019bcd9e38bSBarry Smith factorization may be faster if do not pivot */ 102015e8a5b3SHong Zhang PetscReal zeropivot; /* pivot is called zero if less than this */ 1021f4db908eSBarry Smith PetscReal shifttype; /* type of shift added to matrix factor to prevent zero pivots */ 1022f4db908eSBarry Smith PetscReal shiftamount; /* how large the shift is */ 102315e8a5b3SHong Zhang } MatFactorInfo; 1024ffa6d0a5SLois Curfman McInnes 1025014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFactorInfoInitialize(MatFactorInfo*); 1026014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCholeskyFactor(Mat,IS,const MatFactorInfo*); 1027014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCholeskyFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*); 1028014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCholeskyFactorNumeric(Mat,Mat,const MatFactorInfo*); 1029014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLUFactor(Mat,IS,IS,const MatFactorInfo*); 1030014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatILUFactor(Mat,IS,IS,const MatFactorInfo*); 1031014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*); 1032014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatILUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*); 1033014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatICCFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*); 1034014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatICCFactor(Mat,IS,const MatFactorInfo*); 1035014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatLUFactorNumeric(Mat,Mat,const MatFactorInfo*); 1036014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*); 1037014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolve(Mat,Vec,Vec); 1038014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatForwardSolve(Mat,Vec,Vec); 1039014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatBackwardSolve(Mat,Vec,Vec); 1040014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveAdd(Mat,Vec,Vec,Vec); 1041014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveTranspose(Mat,Vec,Vec); 1042014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolveTransposeAdd(Mat,Vec,Vec,Vec); 1043014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSolves(Mat,Vecs,Vecs); 10448ed539a5SBarry Smith 1045014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetUnfactored(Mat); 1046bb5a7306SBarry Smith 1047d91e6319SBarry Smith /*E 1048d91e6319SBarry Smith MatSORType - What type of (S)SOR to perform 1049bb1eb677SSatish Balay 1050d91e6319SBarry Smith Level: beginner 1051d91e6319SBarry Smith 1052d9274352SBarry Smith May be bitwise ORd together 1053d9274352SBarry Smith 1054af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 1055d91e6319SBarry Smith 10564e7234bfSBarry Smith MatSORType may be bitwise ORd together, so do not change the numbers 10574e7234bfSBarry Smith 105841f059aeSBarry Smith .seealso: MatSOR() 1059d91e6319SBarry Smith E*/ 1060ee50ffe9SBarry Smith typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3, 1061ee50ffe9SBarry Smith SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8, 1062ee50ffe9SBarry Smith SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16, 106384cb2905SBarry Smith SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType; 1064014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSOR(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec); 10658ed539a5SBarry Smith 1066d4fbbf0eSBarry Smith /* 1067639f9d9dSBarry Smith These routines are for efficiently computing Jacobians via finite differences. 1068639f9d9dSBarry Smith */ 1069b12f92e5SBarry Smith 10707086a01eSPeter Brune /*S 10717086a01eSPeter Brune MatColoring - Object for managing the coloring of matrices. 10727086a01eSPeter Brune 10737086a01eSPeter Brune Level: beginner 10747086a01eSPeter Brune 10757086a01eSPeter Brune Concepts: matrix, coloring 10767086a01eSPeter Brune 10777086a01eSPeter Brune .seealso: MatFDColoringCreate() ISColoring MatFDColoring 10787086a01eSPeter Brune S*/ 10797086a01eSPeter Brune typedef struct _p_MatColoring* MatColoring; 108076bdecfbSBarry Smith /*J 10818f6c3df8SBarry Smith MatColoringType - String with the name of a PETSc matrix coloring 1082d9274352SBarry Smith 1083d9274352SBarry Smith Level: beginner 1084d9274352SBarry Smith 10857086a01eSPeter Brune .seealso: MatColoringSetType(), MatColoring 108676bdecfbSBarry Smith J*/ 10877086a01eSPeter Brune 108819fd82e9SBarry Smith typedef const char* MatColoringType; 10895567d71aSPeter Brune #define MATCOLORINGJP "jp" 1090b9acb617SPeter Brune #define MATCOLORINGPOWER "power" 10912692d6eeSBarry Smith #define MATCOLORINGNATURAL "natural" 10922692d6eeSBarry Smith #define MATCOLORINGSL "sl" 10932692d6eeSBarry Smith #define MATCOLORINGLF "lf" 10942692d6eeSBarry Smith #define MATCOLORINGID "id" 10958121bdceSPeter Brune #define MATCOLORINGGREEDY "greedy" 1096b12f92e5SBarry Smith 10978ac6da0aSPeter Brune /*E 10988ac6da0aSPeter Brune MatColoringWeightType - Type of weight scheme 10998ac6da0aSPeter Brune 11008ac6da0aSPeter Brune Not Collective 11018ac6da0aSPeter Brune 11028ac6da0aSPeter Brune + MAT_COLORING_RANDOM - Random weights 11038ac6da0aSPeter Brune . MAT_COLORING_LEXICAL - Lexical weighting based upon global numbering. 11048ac6da0aSPeter Brune - MAT_COLORING_LF - Last-first weighting. 11058ac6da0aSPeter Brune 11068ac6da0aSPeter Brune Level: intermediate 11078ac6da0aSPeter Brune 1108af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 11098ac6da0aSPeter Brune 11108ac6da0aSPeter Brune .seealso: MatCUSPSetFormat(), MatCUSPFormatOperation 11118ac6da0aSPeter Brune E*/ 1112301c140bSPeter Brune typedef enum {MAT_COLORING_WEIGHT_RANDOM,MAT_COLORING_WEIGHT_LEXICAL,MAT_COLORING_WEIGHT_LF,MAT_COLORING_WEIGHT_SL} MatColoringWeightType; 11138ac6da0aSPeter Brune 1114335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringCreate(Mat,MatColoring*); 1115d7f2a307SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringGetDegrees(Mat,PetscInt,PetscInt*); 1116335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringDestroy(MatColoring*); 1117335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringView(MatColoring,PetscViewer); 1118335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetType(MatColoring,MatColoringType); 1119335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetFromOptions(MatColoring); 1120335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetDistance(MatColoring,PetscInt); 1121335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringGetDistance(MatColoring,PetscInt*); 1122335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetMaxColors(MatColoring,PetscInt); 1123335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringGetMaxColors(MatColoring,PetscInt*); 1124335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringApply(MatColoring,ISColoring*); 1125335efc43SPeter Brune PETSC_EXTERN PetscErrorCode MatColoringRegister(const char[],PetscErrorCode(*)(MatColoring)); 1126014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*); 11278ac6da0aSPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetWeightType(MatColoring,MatColoringWeightType); 1128c29971fcSPeter Brune PETSC_EXTERN PetscErrorCode MatColoringSetWeights(MatColoring,PetscReal*,PetscInt*); 11298ac6da0aSPeter Brune PETSC_EXTERN PetscErrorCode MatColoringCreateWeights(MatColoring,PetscReal **,PetscInt **lperm); 1130639f9d9dSBarry Smith 1131d9274352SBarry Smith /*S 1132d9274352SBarry Smith MatFDColoring - Object for computing a sparse Jacobian via finite differences 1133d9274352SBarry Smith and coloring 1134639f9d9dSBarry Smith 1135d9274352SBarry Smith Level: beginner 1136d9274352SBarry Smith 1137d9274352SBarry Smith Concepts: coloring, sparse Jacobian, finite differences 1138d9274352SBarry Smith 1139d9274352SBarry Smith .seealso: MatFDColoringCreate() 1140d9274352SBarry Smith S*/ 1141e2a1c21fSSatish Balay typedef struct _p_MatFDColoring* MatFDColoring; 1142639f9d9dSBarry Smith 1143014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringCreate(Mat,ISColoring,MatFDColoring *); 1144014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringDestroy(MatFDColoring*); 1145014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringView(MatFDColoring,PetscViewer); 1146014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*); 1147014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**); 1148014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal); 1149014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetFromOptions(MatFDColoring); 1150d1e9a80fSBarry Smith PETSC_EXTERN PetscErrorCode MatFDColoringApply(Mat,MatFDColoring,Vec,void *); 1151014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringSetF(MatFDColoring,Vec); 1152014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,PetscInt*[]); 1153f86b9fbaSHong Zhang PETSC_EXTERN PetscErrorCode MatFDColoringSetUp(Mat,ISColoring,MatFDColoring); 1154f86b9fbaSHong Zhang PETSC_EXTERN PetscErrorCode MatFDColoringSetBlockSize(MatFDColoring,PetscInt,PetscInt); 1155f86b9fbaSHong Zhang 1156b1683b59SHong Zhang 1157b1683b59SHong Zhang /*S 1158b9af6bddSHong Zhang MatTransposeColoring - Object for computing a sparse matrix product C=A*B^T via coloring 1159b1683b59SHong Zhang 1160b1683b59SHong Zhang Level: beginner 1161b1683b59SHong Zhang 1162b1683b59SHong Zhang Concepts: coloring, sparse matrix product 1163b1683b59SHong Zhang 1164b9af6bddSHong Zhang .seealso: MatTransposeColoringCreate() 1165b1683b59SHong Zhang S*/ 1166b9af6bddSHong Zhang typedef struct _p_MatTransposeColoring* MatTransposeColoring; 1167b1683b59SHong Zhang 1168014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeColoringCreate(Mat,ISColoring,MatTransposeColoring *); 1169014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransColoringApplySpToDen(MatTransposeColoring,Mat,Mat); 1170014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransColoringApplyDenToSp(MatTransposeColoring,Mat,Mat); 1171014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatTransposeColoringDestroy(MatTransposeColoring*); 1172b1683b59SHong Zhang 1173639f9d9dSBarry Smith /* 11740752156aSBarry Smith These routines are for partitioning matrices: currently used only 11753eda8832SBarry Smith for adjacency matrix, MatCreateMPIAdj(). 11760752156aSBarry Smith */ 1177ca161407SBarry Smith 1178d9274352SBarry Smith /*S 1179d9274352SBarry Smith MatPartitioning - Object for managing the partitioning of a matrix or graph 1180d9274352SBarry Smith 1181d9274352SBarry Smith Level: beginner 1182d9274352SBarry Smith 1183d9274352SBarry Smith Concepts: partitioning 1184d9274352SBarry Smith 1185743a1026Svictorle .seealso: MatPartitioningCreate(), MatPartitioningType 1186d9274352SBarry Smith S*/ 118791e9ee9fSBarry Smith typedef struct _p_MatPartitioning* MatPartitioning; 1188d9274352SBarry Smith 118976bdecfbSBarry Smith /*J 11908f6c3df8SBarry Smith MatPartitioningType - String with the name of a PETSc matrix partitioning 1191d9274352SBarry Smith 1192d9274352SBarry Smith Level: beginner 11930d04baf8SBarry Smith dm 1194b547a13bSvictorle .seealso: MatPartitioningCreate(), MatPartitioning 119576bdecfbSBarry Smith J*/ 119619fd82e9SBarry Smith typedef const char* MatPartitioningType; 11972692d6eeSBarry Smith #define MATPARTITIONINGCURRENT "current" 1198*aa1e27eaSFande Kong #define MATPARTITIONINGAVERAGE "average" 11992692d6eeSBarry Smith #define MATPARTITIONINGSQUARE "square" 12002692d6eeSBarry Smith #define MATPARTITIONINGPARMETIS "parmetis" 12012692d6eeSBarry Smith #define MATPARTITIONINGCHACO "chaco" 12022692d6eeSBarry Smith #define MATPARTITIONINGPARTY "party" 120350d91057SBarry Smith #define MATPARTITIONINGPTSCOTCH "ptscotch" 120488d2ac2bSFande Kong #define MATPARTITIONINGHIERARCH "hierarch" 120571306c60Sjroman 1206ca161407SBarry Smith 1207014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningCreate(MPI_Comm,MatPartitioning*); 120819fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningSetType(MatPartitioning,MatPartitioningType); 1209014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetNParts(MatPartitioning,PetscInt); 1210014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetAdjacency(MatPartitioning,Mat); 1211014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]); 1212014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []); 1213014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningApply(MatPartitioning,IS*); 1214014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningDestroy(MatPartitioning*); 12152aabb6bbSBarry Smith 1216bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningRegister(const char[],PetscErrorCode (*)(MatPartitioning)); 121730de9b25SBarry Smith 1218f1144a30SSatish Balay 12192bad1931SBarry Smith 1220014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningView(MatPartitioning,PetscViewer); 1221014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningSetFromOptions(MatPartitioning); 122219fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningGetType(MatPartitioning,MatPartitioningType*); 1223ca161407SBarry Smith 122427538973SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetRepartition(MatPartitioning part); 1225014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetCoarseSequential(MatPartitioning); 1226014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningParmetisGetEdgeCut(MatPartitioning, PetscInt *); 12270752156aSBarry Smith 1228b6956c12SJose Roman typedef enum { MP_CHACO_MULTILEVEL=1,MP_CHACO_SPECTRAL=2,MP_CHACO_LINEAR=4,MP_CHACO_RANDOM=5,MP_CHACO_SCATTERED=6 } MPChacoGlobalType; 12296a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoGlobalTypes[]; 1230b6956c12SJose Roman typedef enum { MP_CHACO_KERNIGHAN=1,MP_CHACO_NONE=2 } MPChacoLocalType; 12316a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoLocalTypes[]; 1232b6956c12SJose Roman typedef enum { MP_CHACO_LANCZOS=0,MP_CHACO_RQI=1 } MPChacoEigenType; 12336a6fc655SJed Brown PETSC_EXTERN const char *const MPChacoEigenTypes[]; 1234b6956c12SJose Roman 1235014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetGlobal(MatPartitioning,MPChacoGlobalType); 1236014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetGlobal(MatPartitioning,MPChacoGlobalType*); 1237014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetLocal(MatPartitioning,MPChacoLocalType); 1238014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetLocal(MatPartitioning,MPChacoLocalType*); 1239014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal); 1240014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType); 1241014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenSolver(MatPartitioning,MPChacoEigenType*); 1242014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenTol(MatPartitioning,PetscReal); 1243014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenTol(MatPartitioning,PetscReal*); 1244014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenNumber(MatPartitioning,PetscInt); 1245014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenNumber(MatPartitioning,PetscInt*); 124671306c60Sjroman 124771306c60Sjroman #define MP_PARTY_OPT "opt" 124871306c60Sjroman #define MP_PARTY_LIN "lin" 124971306c60Sjroman #define MP_PARTY_SCA "sca" 125071306c60Sjroman #define MP_PARTY_RAN "ran" 125171306c60Sjroman #define MP_PARTY_GBF "gbf" 125271306c60Sjroman #define MP_PARTY_GCF "gcf" 125371306c60Sjroman #define MP_PARTY_BUB "bub" 125471306c60Sjroman #define MP_PARTY_DEF "def" 1255014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetGlobal(MatPartitioning,const char*); 125671306c60Sjroman #define MP_PARTY_HELPFUL_SETS "hs" 125771306c60Sjroman #define MP_PARTY_KERNIGHAN_LIN "kl" 125871306c60Sjroman #define MP_PARTY_NONE "no" 1259014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetLocal(MatPartitioning,const char*); 1260014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal); 1261014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetBipart(MatPartitioning,PetscBool); 1262014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscBool); 126371306c60Sjroman 126450d91057SBarry Smith typedef enum { MP_PTSCOTCH_QUALITY,MP_PTSCOTCH_SPEED,MP_PTSCOTCH_BALANCE,MP_PTSCOTCH_SAFETY,MP_PTSCOTCH_SCALABILITY } MPPTScotchStrategyType; 12656a6fc655SJed Brown PETSC_EXTERN const char *const MPPTScotchStrategyTypes[]; 1266e0f1cffaSJose Roman 1267014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetImbalance(MatPartitioning,PetscReal); 1268014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetImbalance(MatPartitioning,PetscReal*); 1269014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetStrategy(MatPartitioning,MPPTScotchStrategyType); 1270014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetStrategy(MatPartitioning,MPPTScotchStrategyType*); 127171306c60Sjroman 1272b43b03e9SMark F. Adams /* 12736294fa2bSFande Kong * hierarchical partitioning 12746294fa2bSFande Kong */ 12754e713f55SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalGetFineparts(MatPartitioning,IS*); 12764e713f55SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalGetCoarseparts(MatPartitioning,IS*); 12774e713f55SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalSetNcoarseparts(MatPartitioning,PetscInt); 12784e713f55SFande Kong PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalSetNfineparts(MatPartitioning, PetscInt); 12796294fa2bSFande Kong 12806294fa2bSFande Kong /* 1281b43b03e9SMark F. Adams These routines are for coarsening matrices: 1282b43b03e9SMark F. Adams */ 1283b43b03e9SMark F. Adams 1284b43b03e9SMark F. Adams /*S 1285b43b03e9SMark F. Adams MatCoarsen - Object for managing the coarsening of a graph (symmetric matrix) 1286b43b03e9SMark F. Adams 1287b43b03e9SMark F. Adams Level: beginner 1288b43b03e9SMark F. Adams 1289b43b03e9SMark F. Adams Concepts: coarsen 1290b43b03e9SMark F. Adams 1291b43b03e9SMark F. Adams .seealso: MatCoarsenCreate), MatCoarsenType 1292b43b03e9SMark F. Adams S*/ 1293b43b03e9SMark F. Adams typedef struct _p_MatCoarsen* MatCoarsen; 1294b43b03e9SMark F. Adams 1295b43b03e9SMark F. Adams /*J 12968f6c3df8SBarry Smith MatCoarsenType - String with the name of a PETSc matrix coarsen 1297b43b03e9SMark F. Adams 1298b43b03e9SMark F. Adams Level: beginner 12998f6c3df8SBarry Smith 1300b43b03e9SMark F. Adams .seealso: MatCoarsenCreate(), MatCoarsen 1301b43b03e9SMark F. Adams J*/ 130219fd82e9SBarry Smith typedef const char* MatCoarsenType; 1303b43b03e9SMark F. Adams #define MATCOARSENMIS "mis" 13049057884aSMark F. Adams #define MATCOARSENHEM "hem" 1305b43b03e9SMark F. Adams 13060cbbd2e1SMark F. Adams /* linked list for aggregates */ 130741b27cdeSMark F. Adams typedef struct _PetscCDIntNd{ 130841b27cdeSMark F. Adams struct _PetscCDIntNd *next; 13090cbbd2e1SMark F. Adams PetscInt gid; 131041b27cdeSMark F. Adams }PetscCDIntNd; 13110cbbd2e1SMark F. Adams 13120cbbd2e1SMark F. Adams /* only used by node pool */ 131341b27cdeSMark F. Adams typedef struct _PetscCDArrNd{ 131441b27cdeSMark F. Adams struct _PetscCDArrNd *next; 131541b27cdeSMark F. Adams struct _PetscCDIntNd *array; 131641b27cdeSMark F. Adams }PetscCDArrNd; 13170cbbd2e1SMark F. Adams 13180cbbd2e1SMark F. Adams typedef struct _PetscCoarsenData{ 131982c86c8fSBarry Smith PetscCDArrNd pool_list; /* node pool */ 132041b27cdeSMark F. Adams PetscCDIntNd *new_node; 13210cbbd2e1SMark F. Adams PetscInt new_left; 13220cbbd2e1SMark F. Adams PetscInt chk_sz; 132341b27cdeSMark F. Adams PetscCDIntNd *extra_nodes; 132482c86c8fSBarry Smith PetscCDIntNd **array; /* Array of lists */ 13250cbbd2e1SMark F. Adams PetscInt size; 132682c86c8fSBarry Smith Mat mat; /* cache a Mat for communication data */ 13270cbbd2e1SMark F. Adams }PetscCoarsenData; 13280cbbd2e1SMark F. Adams 1329014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenCreate(MPI_Comm,MatCoarsen*); 133019fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenSetType(MatCoarsen,MatCoarsenType); 1331014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetAdjacency(MatCoarsen,Mat); 1332014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetGreedyOrdering(MatCoarsen,const IS); 1333014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetStrictAggs(MatCoarsen,PetscBool); 1334014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenGetData( MatCoarsen, PetscCoarsenData ** ); 1335014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenApply(MatCoarsen); 1336014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenDestroy(MatCoarsen*); 1337b43b03e9SMark F. Adams 1338bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenRegister(const char[],PetscErrorCode (*)(MatCoarsen)); 1339b43b03e9SMark F. Adams 1340b43b03e9SMark F. Adams 1341b43b03e9SMark F. Adams 1342014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenView(MatCoarsen,PetscViewer); 1343014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCoarsenSetFromOptions(MatCoarsen); 134419fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenGetType(MatCoarsen,MatCoarsenType*); 1345685405a1SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatCoarsenViewFromOptions(MatCoarsen A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);} 1346b43b03e9SMark F. Adams 1347014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMeshToVertexGraph(Mat,PetscInt,Mat*); 1348014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMeshToCellGraph(Mat,PetscInt,Mat*); 1349591294e4SBarry Smith 13500752156aSBarry Smith /* 1351af0996ceSBarry Smith If you add entries here you must also add them to petsc/finclude/petscmat.h 1352d4fbbf0eSBarry Smith */ 13531c1c02c0SLois Curfman McInnes typedef enum { MATOP_SET_VALUES=0, 13541c1c02c0SLois Curfman McInnes MATOP_GET_ROW=1, 13551c1c02c0SLois Curfman McInnes MATOP_RESTORE_ROW=2, 13561c1c02c0SLois Curfman McInnes MATOP_MULT=3, 13571c1c02c0SLois Curfman McInnes MATOP_MULT_ADD=4, 13587c922b88SBarry Smith MATOP_MULT_TRANSPOSE=5, 13597c922b88SBarry Smith MATOP_MULT_TRANSPOSE_ADD=6, 13601c1c02c0SLois Curfman McInnes MATOP_SOLVE=7, 13611c1c02c0SLois Curfman McInnes MATOP_SOLVE_ADD=8, 13627c922b88SBarry Smith MATOP_SOLVE_TRANSPOSE=9, 13637c922b88SBarry Smith MATOP_SOLVE_TRANSPOSE_ADD=10, 13641c1c02c0SLois Curfman McInnes MATOP_LUFACTOR=11, 13651c1c02c0SLois Curfman McInnes MATOP_CHOLESKYFACTOR=12, 1366a714c06dSBarry Smith MATOP_SOR=13, 13671c1c02c0SLois Curfman McInnes MATOP_TRANSPOSE=14, 13681c1c02c0SLois Curfman McInnes MATOP_GETINFO=15, 13691c1c02c0SLois Curfman McInnes MATOP_EQUAL=16, 13701c1c02c0SLois Curfman McInnes MATOP_GET_DIAGONAL=17, 13711c1c02c0SLois Curfman McInnes MATOP_DIAGONAL_SCALE=18, 13721c1c02c0SLois Curfman McInnes MATOP_NORM=19, 13731c1c02c0SLois Curfman McInnes MATOP_ASSEMBLY_BEGIN=20, 13741c1c02c0SLois Curfman McInnes MATOP_ASSEMBLY_END=21, 1375d519adbfSMatthew Knepley MATOP_SET_OPTION=22, 1376d519adbfSMatthew Knepley MATOP_ZERO_ENTRIES=23, 1377d519adbfSMatthew Knepley MATOP_ZERO_ROWS=24, 1378d519adbfSMatthew Knepley MATOP_LUFACTOR_SYMBOLIC=25, 1379d519adbfSMatthew Knepley MATOP_LUFACTOR_NUMERIC=26, 1380d519adbfSMatthew Knepley MATOP_CHOLESKY_FACTOR_SYMBOLIC=27, 1381d519adbfSMatthew Knepley MATOP_CHOLESKY_FACTOR_NUMERIC=28, 1382d519adbfSMatthew Knepley MATOP_SETUP_PREALLOCATION=29, 1383d519adbfSMatthew Knepley MATOP_ILUFACTOR_SYMBOLIC=30, 1384d519adbfSMatthew Knepley MATOP_ICCFACTOR_SYMBOLIC=31, 13859d39f69dSJed Brown /* MATOP_PLACEHOLDER_32=32, */ 13869d39f69dSJed Brown /* MATOP_PLACEHOLDER_33=33, */ 1387d519adbfSMatthew Knepley MATOP_DUPLICATE=34, 1388d519adbfSMatthew Knepley MATOP_FORWARD_SOLVE=35, 1389d519adbfSMatthew Knepley MATOP_BACKWARD_SOLVE=36, 1390d519adbfSMatthew Knepley MATOP_ILUFACTOR=37, 1391d519adbfSMatthew Knepley MATOP_ICCFACTOR=38, 1392d519adbfSMatthew Knepley MATOP_AXPY=39, 1393d519adbfSMatthew Knepley MATOP_GET_SUBMATRICES=40, 1394d519adbfSMatthew Knepley MATOP_INCREASE_OVERLAP=41, 1395d519adbfSMatthew Knepley MATOP_GET_VALUES=42, 1396d519adbfSMatthew Knepley MATOP_COPY=43, 1397d519adbfSMatthew Knepley MATOP_GET_ROW_MAX=44, 1398d519adbfSMatthew Knepley MATOP_SCALE=45, 1399d519adbfSMatthew Knepley MATOP_SHIFT=46, 140035153367SBarry Smith MATOP_DIAGONAL_SET=47, 14019d39f69dSJed Brown MATOP_ZERO_ROWS_COLUMNS=48, 14029d39f69dSJed Brown MATOP_SET_RANDOM=49, 1403d519adbfSMatthew Knepley MATOP_GET_ROW_IJ=50, 1404d519adbfSMatthew Knepley MATOP_RESTORE_ROW_IJ=51, 1405d519adbfSMatthew Knepley MATOP_GET_COLUMN_IJ=52, 1406d519adbfSMatthew Knepley MATOP_RESTORE_COLUMN_IJ=53, 1407d519adbfSMatthew Knepley MATOP_FDCOLORING_CREATE=54, 1408d519adbfSMatthew Knepley MATOP_COLORING_PATCH=55, 1409d519adbfSMatthew Knepley MATOP_SET_UNFACTORED=56, 1410d519adbfSMatthew Knepley MATOP_PERMUTE=57, 1411d519adbfSMatthew Knepley MATOP_SET_VALUES_BLOCKED=58, 1412d519adbfSMatthew Knepley MATOP_GET_SUBMATRIX=59, 1413d519adbfSMatthew Knepley MATOP_DESTROY=60, 1414d519adbfSMatthew Knepley MATOP_VIEW=61, 1415d519adbfSMatthew Knepley MATOP_CONVERT_FROM=62, 14167bab7c10SHong Zhang MATOP_MATMAT_MULT=63, 14177bab7c10SHong Zhang MATOP_MATMAT_MULT_SYMBOLIC=64, 14187bab7c10SHong Zhang MATOP_MATMAT_MULT_NUMERIC=65, 1419d519adbfSMatthew Knepley MATOP_SET_LOCAL_TO_GLOBAL_MAP=66, 1420d519adbfSMatthew Knepley MATOP_SET_VALUES_LOCAL=67, 1421d519adbfSMatthew Knepley MATOP_ZERO_ROWS_LOCAL=68, 1422d519adbfSMatthew Knepley MATOP_GET_ROW_MAX_ABS=69, 1423d519adbfSMatthew Knepley MATOP_GET_ROW_MIN_ABS=70, 1424d519adbfSMatthew Knepley MATOP_CONVERT=71, 1425d519adbfSMatthew Knepley MATOP_SET_COLORING=72, 14269d39f69dSJed Brown /* MATOP_PLACEHOLDER_73=73, */ 1427d519adbfSMatthew Knepley MATOP_SET_VALUES_ADIFOR=74, 1428d519adbfSMatthew Knepley MATOP_FD_COLORING_APPLY=75, 1429d519adbfSMatthew Knepley MATOP_SET_FROM_OPTIONS=76, 1430bda6c4cbSBarry Smith MATOP_MULT_CONSTRAINED=77, 1431bda6c4cbSBarry Smith MATOP_MULT_TRANSPOSE_CONSTRAIN=78, 14329d39f69dSJed Brown MATOP_FIND_ZERO_DIAGONALS=79, 1433d519adbfSMatthew Knepley MATOP_MULT_MULTIPLE=80, 1434d519adbfSMatthew Knepley MATOP_SOLVE_MULTIPLE=81, 1435d519adbfSMatthew Knepley MATOP_GET_INERTIA=82, 1436d519adbfSMatthew Knepley MATOP_LOAD=83, 1437d519adbfSMatthew Knepley MATOP_IS_SYMMETRIC=84, 1438d519adbfSMatthew Knepley MATOP_IS_HERMITIAN=85, 1439d519adbfSMatthew Knepley MATOP_IS_STRUCTURALLY_SYMMETRIC=86, 1440820469bcSHong Zhang MATOP_SET_VALUES_BLOCKEDLOCAL=87, 1441d519adbfSMatthew Knepley MATOP_GET_VECS=88, 1442d519adbfSMatthew Knepley MATOP_MAT_MULT=89, 1443d519adbfSMatthew Knepley MATOP_MAT_MULT_SYMBOLIC=90, 1444d519adbfSMatthew Knepley MATOP_MAT_MULT_NUMERIC=91, 1445d519adbfSMatthew Knepley MATOP_PTAP=92, 1446d519adbfSMatthew Knepley MATOP_PTAP_SYMBOLIC=93, 1447d519adbfSMatthew Knepley MATOP_PTAP_NUMERIC=94, 1448bda6c4cbSBarry Smith MATOP_MAT_TRANSPOSE_MULT=95, 1449bda6c4cbSBarry Smith MATOP_MAT_TRANSPOSE_MULT_SYMBO=96, 1450bda6c4cbSBarry Smith MATOP_MAT_TRANSPOSE_MULT_NUMER=97, 14519d39f69dSJed Brown /* MATOP_PLACEHOLDER_98=98, */ 14529d39f69dSJed Brown /* MATOP_PLACEHOLDER_99=99, */ 14539d39f69dSJed Brown /* MATOP_PLACEHOLDER_100=100, */ 14549d39f69dSJed Brown /* MATOP_PLACEHOLDER_101=101, */ 1455d519adbfSMatthew Knepley MATOP_CONJUGATE=102, 14569d39f69dSJed Brown /* MATOP_PLACEHOLDER_103=103, */ 1457d519adbfSMatthew Knepley MATOP_SET_VALUES_ROW=104, 1458d519adbfSMatthew Knepley MATOP_REAL_PART=105, 1459bda6c4cbSBarry Smith MATOP_IMAGINARY_PART=106, 1460bda6c4cbSBarry Smith MATOP_GET_ROW_UPPER_TRIANGULAR=107, 1461bda6c4cbSBarry Smith MATOP_RESTORE_ROW_UPPER_TRIANG=108, 1462bda6c4cbSBarry Smith MATOP_MAT_SOLVE=109, 14630e8d04f7SBarry Smith MATOP_GET_REDUNDANT_MATRIX=110, 1464d519adbfSMatthew Knepley MATOP_GET_ROW_MIN=111, 14650e8d04f7SBarry Smith MATOP_GET_COLUMN_VECTOR=112, 1466d519adbfSMatthew Knepley MATOP_MISSING_DIAGONAL=113, 14670e8d04f7SBarry Smith MATOP_GET_SEQ_NONZERO_STRUCTUR=114, 146889c6957cSBarry Smith MATOP_CREATE=115, 146989c6957cSBarry Smith MATOP_GET_GHOSTS=116, 14700e8d04f7SBarry Smith MATOP_GET_LOCAL_SUB_MATRIX=117, 14710e8d04f7SBarry Smith MATOP_RESTORE_LOCALSUB_MATRIX=118, 1472eeffb40dSHong Zhang MATOP_MULT_DIAGONAL_BLOCK=119, 14730e8d04f7SBarry Smith MATOP_HERMITIAN_TRANSPOSE=120, 14740e8d04f7SBarry Smith MATOP_MULT_HERMITIAN_TRANSPOSE=121, 14750e8d04f7SBarry Smith MATOP_MULT_HERMITIAN_TRANS_ADD=122, 14760e8d04f7SBarry Smith MATOP_GET_MULTI_PROC_BLOCK=123, 14779d39f69dSJed Brown MATOP_FIND_NONZERO_ROWS=124, 14780e8d04f7SBarry Smith MATOP_GET_COLUMN_NORMS=125, 14799d39f69dSJed Brown MATOP_INVERT_BLOCK_DIAGONAL=126, 14809d39f69dSJed Brown /* MATOP_PLACEHOLDER_127=127, */ 14810e8d04f7SBarry Smith MATOP_GET_SUB_MATRICES_PARALLE=128, 14822b8ad9a3SHong Zhang MATOP_SET_VALUES_BATCH=129, 14830e8d04f7SBarry Smith MATOP_TRANSPOSE_MAT_MULT=130, 14840e8d04f7SBarry Smith MATOP_TRANSPOSE_MAT_MULT_SYMBO=131, 1485e9b602ebSSatish Balay MATOP_TRANSPOSE_MAT_MULT_NUMER=132, 14860e8d04f7SBarry Smith MATOP_TRANSPOSE_COLORING_CREAT=133, 14870e8d04f7SBarry Smith MATOP_TRANS_COLORING_APPLY_SPT=134, 14880e8d04f7SBarry Smith MATOP_TRANS_COLORING_APPLY_DEN=135, 14890e8d04f7SBarry Smith MATOP_RART=136, 14900e8d04f7SBarry Smith MATOP_RART_SYMBOLIC=137, 14910e8d04f7SBarry Smith MATOP_RART_NUMERIC=138, 1492e09a3074SHong Zhang MATOP_SET_BLOCK_SIZES=139, 1493f9426fe0SMark Adams MATOP_AYPX=140, 14941919a2e2SJed Brown MATOP_RESIDUAL=141, 14959c8f2541SHong Zhang MATOP_FDCOLORING_SETUP=142, 14969c8f2541SHong Zhang MATOP_MPICONCATENATESEQ=144 1497fae171e0SBarry Smith } MatOperation; 1498014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatHasOperation(Mat,MatOperation,PetscBool *); 1499014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellSetOperation(Mat,MatOperation,void(*)(void)); 1500014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellGetOperation(Mat,MatOperation,void(**)(void)); 1501014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatShellSetContext(Mat,void*); 1502112a2221SBarry Smith 150390ace30eSBarry Smith /* 150490ace30eSBarry Smith Codes for matrices stored on disk. By default they are 150590ace30eSBarry Smith stored in a universal format. By changing the format with 15067973ad04SBarry Smith PetscViewerSetFormat(viewer,PETSC_VIEWER_NATIVE); the matrices will 150790ace30eSBarry Smith be stored in a way natural for the matrix, for example dense matrices 150890ace30eSBarry Smith would be stored as dense. Matrices stored this way may only be 1509f4403165SShri Abhyankar read into matrices of the same type. 151090ace30eSBarry Smith */ 151190ace30eSBarry Smith #define MATRIX_BINARY_FORMAT_DENSE -1 151290ace30eSBarry Smith 1513014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPIBAIJSetHashTableFactor(Mat,PetscReal); 1514014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatISGetLocalMat(Mat,Mat*); 1515014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatISSetLocalMat(Mat,Mat); 1516b7ce53b6SStefano Zampini PETSC_EXTERN PetscErrorCode MatISGetMPIXAIJ(Mat,MatReuse,Mat*); 15171f4e1ec7SBarry Smith 1518d9274352SBarry Smith /*S 1519d9274352SBarry Smith MatNullSpace - Object that removes a null space from a vector, i.e. 1520d9274352SBarry Smith orthogonalizes the vector to a subsapce 1521d9274352SBarry Smith 1522f7a9e4ceSBarry Smith Level: advanced 1523d9274352SBarry Smith 1524d9274352SBarry Smith Concepts: matrix; linear operator, null space 1525d9274352SBarry Smith 15266e1639daSBarry Smith Users manual sections: 15276e1639daSBarry Smith . sec_singular 15286e1639daSBarry Smith 1529d9274352SBarry Smith .seealso: MatNullSpaceCreate() 1530d9274352SBarry Smith S*/ 153174637425SBarry Smith typedef struct _p_MatNullSpace* MatNullSpace; 1532d9274352SBarry Smith 1533014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceCreate(MPI_Comm,PetscBool ,PetscInt,const Vec[],MatNullSpace*); 1534014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(MatNullSpace,Vec,void*),void*); 1535014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceDestroy(MatNullSpace*); 1536d0195637SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceRemove(MatNullSpace,Vec); 1537014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNullSpace(Mat, MatNullSpace *); 15385fa7ec2dSBarry Smith PETSC_EXTERN PetscErrorCode MatGetTransposeNullSpace(Mat, MatNullSpace *); 15395fa7ec2dSBarry Smith PETSC_EXTERN PetscErrorCode MatSetTransposeNullSpace(Mat,MatNullSpace); 1540014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNullSpace(Mat,MatNullSpace); 1541014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSetNearNullSpace(Mat,MatNullSpace); 1542014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatGetNearNullSpace(Mat,MatNullSpace*); 1543014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceTest(MatNullSpace,Mat,PetscBool *); 1544014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceView(MatNullSpace,PetscViewer); 1545014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceGetVecs(MatNullSpace,PetscBool*,PetscInt*,const Vec**); 1546014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNullSpaceCreateRigidBody(Vec,MatNullSpace*); 154774637425SBarry Smith 1548014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatReorderingSeqSBAIJ(Mat,IS); 1549014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMPISBAIJSetHashTableFactor(Mat,PetscReal); 1550014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetColumnIndices(Mat,PetscInt *); 1551014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSeqBAIJInvertBlockDiagonal(Mat); 15523f1d51d7SBarry Smith 1553014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMAIJ(Mat,PetscInt,Mat*); 1554014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMAIJRedimension(Mat,PetscInt,Mat*); 1555014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMAIJGetAIJ(Mat,Mat*); 1556c4f061fbSSatish Balay 1557014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatComputeExplicitOperator(Mat,Mat*); 1558b0a32e0cSBarry Smith 1559014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatDiagonalScaleLocal(Mat,Vec); 156004f1ad80SBarry Smith 1561014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*); 1562014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetBase(Mat,Vec,Vec); 1563014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*); 1564014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*)); 1565014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec)); 1566014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDAddNullSpace(Mat,MatNullSpace); 1567014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt); 1568014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDResetHHistory(Mat); 1569014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetFunctionError(Mat,PetscReal); 1570014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetPeriod(Mat,PetscInt); 1571014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDGetH(Mat,PetscScalar *); 1572014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetOptionsPrefix(Mat,const char[]); 1573014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*); 1574014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*); 1575e884886fSBarry Smith 15766370053bSBarry Smith /*S 15776370053bSBarry Smith MatMFFD - A data structured used to manage the computation of the h differencing parameter for matrix-free 15786370053bSBarry Smith Jacobian vector products 1579e884886fSBarry Smith 15806370053bSBarry Smith Notes: MATMFFD is a specific MatType which uses the MatMFFD data structure 15816370053bSBarry Smith 15826370053bSBarry Smith MatMFFD*() methods actually take the Mat as their first argument. Not a MatMFFD data structure 15836370053bSBarry Smith 15846370053bSBarry Smith Level: developer 15856370053bSBarry Smith 15866370053bSBarry Smith .seealso: MATMFFD, MatCreateMFFD(), MatMFFDSetFuction(), MatMFFDSetType(), MatMFFDRegister() 15876370053bSBarry Smith S*/ 1588e884886fSBarry Smith typedef struct _p_MatMFFD* MatMFFD; 1589e884886fSBarry Smith 159076bdecfbSBarry Smith /*J 1591e884886fSBarry Smith MatMFFDType - algorithm used to compute the h used in computing matrix-vector products via differencing of the function 1592e884886fSBarry Smith 1593e884886fSBarry Smith Level: beginner 1594e884886fSBarry Smith 1595e884886fSBarry Smith .seealso: MatMFFDSetType(), MatMFFDRegister() 159676bdecfbSBarry Smith J*/ 159719fd82e9SBarry Smith typedef const char* MatMFFDType; 1598e884886fSBarry Smith #define MATMFFD_DS "ds" 1599e884886fSBarry Smith #define MATMFFD_WP "wp" 1600e884886fSBarry Smith 160119fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDSetType(Mat,MatMFFDType); 1602bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDRegister(const char[],PetscErrorCode (*)(MatMFFD)); 1603e884886fSBarry Smith 1604014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDDSSetUmin(Mat,PetscReal); 1605014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMFFDWPSetComputeNormU(Mat,PetscBool ); 1606e884886fSBarry Smith 1607014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *); 1608014dd563SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *); 16097dbadf16SMatthew Knepley 161097969023SHong Zhang /* 161197969023SHong Zhang PETSc interface to MUMPS 161297969023SHong Zhang */ 161397969023SHong Zhang #ifdef PETSC_HAVE_MUMPS 1614014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatMumpsSetIcntl(Mat,PetscInt,PetscInt); 1615bc6112feSHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetIcntl(Mat,PetscInt,PetscInt*); 1616b9e9ea26SSatish Balay PETSC_EXTERN PetscErrorCode MatMumpsSetCntl(Mat,PetscInt,PetscReal); 1617bc6112feSHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetCntl(Mat,PetscInt,PetscReal*); 1618bc6112feSHong Zhang 1619ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetInfo(Mat,PetscInt,PetscInt*); 1620ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetInfog(Mat,PetscInt,PetscInt*); 1621ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetRinfo(Mat,PetscInt,PetscReal*); 1622ca810319SHong Zhang PETSC_EXTERN PetscErrorCode MatMumpsGetRinfog(Mat,PetscInt,PetscReal*); 16236444a565SStefano Zampini 16246444a565SStefano Zampini PETSC_EXTERN PetscErrorCode MatMumpsSetSchurIndices(Mat,PetscInt,PetscInt[]); 162559ac8732SStefano Zampini PETSC_EXTERN PetscErrorCode MatMumpsInvertSchurComplement(Mat); 162659ac8732SStefano Zampini PETSC_EXTERN PetscErrorCode MatMumpsCreateSchurComplement(Mat,Mat*); 16276444a565SStefano Zampini PETSC_EXTERN PetscErrorCode MatMumpsGetSchurComplement(Mat,Mat*); 162859ac8732SStefano Zampini PETSC_EXTERN PetscErrorCode MatMumpsRestoreSchurComplement(Mat,Mat*); 1629e807eca7SStefano Zampini PETSC_EXTERN PetscErrorCode MatMumpsSolveSchurComplement(Mat,Vec,Vec); 16307404bcfbSStefano Zampini PETSC_EXTERN PetscErrorCode MatMumpsSolveSchurComplementTranspose(Mat,Vec,Vec); 163197969023SHong Zhang #endif 163223a5497aSJed Brown 1633d954db57SHong Zhang /* 1634d615f992SSatish Balay PETSc interface to Mkl_Pardiso 1635b500a6b7SJose David Bermeol */ 1636b500a6b7SJose David Bermeol #ifdef PETSC_HAVE_MKL_PARDISO 1637d615f992SSatish Balay PETSC_EXTERN PetscErrorCode MatMkl_PardisoSetCntl(Mat,PetscInt,PetscInt); 1638b500a6b7SJose David Bermeol #endif 1639b500a6b7SJose David Bermeol 1640b500a6b7SJose David Bermeol /* 1641d305a81bSVasiliy Kozyrev PETSc interface to Mkl_CPardiso 1642d305a81bSVasiliy Kozyrev */ 1643d305a81bSVasiliy Kozyrev #ifdef PETSC_HAVE_MKL_CPARDISO 1644d305a81bSVasiliy Kozyrev PETSC_EXTERN PetscErrorCode MatMkl_CPardisoSetCntl(Mat,PetscInt,PetscInt); 1645d305a81bSVasiliy Kozyrev #endif 1646d305a81bSVasiliy Kozyrev 1647d305a81bSVasiliy Kozyrev /* 1648d954db57SHong Zhang PETSc interface to SUPERLU 1649d954db57SHong Zhang */ 1650d954db57SHong Zhang #ifdef PETSC_HAVE_SUPERLU 1651014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSuperluSetILUDropTol(Mat,PetscReal); 1652d954db57SHong Zhang #endif 1653d954db57SHong Zhang 1654aa372e3fSPaul Mullowney #ifdef PETSC_HAVE_CUDA 16553f9c0db1SPaul Mullowney /*E 1656e057df02SPaul Mullowney MatCUSPARSEStorageFormat - indicates the storage format for CUSPARSE (GPU) 16572692e278SPaul Mullowney matrices. 1658e057df02SPaul Mullowney 1659e057df02SPaul Mullowney Not Collective 1660e057df02SPaul Mullowney 16618468deeeSKarl Rupp + MAT_CUSPARSE_CSR - Compressed Sparse Row 16622692e278SPaul Mullowney . MAT_CUSPARSE_ELL - Ellpack (requires CUDA 4.2 or later). 16632692e278SPaul Mullowney - MAT_CUSPARSE_HYB - Hybrid, a combination of Ellpack and Coordinate format (requires CUDA 4.2 or later). 1664e057df02SPaul Mullowney 1665e057df02SPaul Mullowney Level: intermediate 1666e057df02SPaul Mullowney 1667af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 1668e057df02SPaul Mullowney 1669e057df02SPaul Mullowney .seealso: MatCUSPARSESetFormat(), MatCUSPARSEFormatOperation 1670e057df02SPaul Mullowney E*/ 1671e057df02SPaul Mullowney 16723f9c0db1SPaul Mullowney typedef enum {MAT_CUSPARSE_CSR, MAT_CUSPARSE_ELL, MAT_CUSPARSE_HYB} MatCUSPARSEStorageFormat; 1673e057df02SPaul Mullowney 1674e057df02SPaul Mullowney /* these will be strings associated with enumerated type defined above */ 1675e057df02SPaul Mullowney PETSC_EXTERN const char *const MatCUSPARSEStorageFormats[]; 1676e057df02SPaul Mullowney 16773f9c0db1SPaul Mullowney /*E 1678e057df02SPaul Mullowney MatCUSPARSEFormatOperation - indicates the operation of CUSPARSE (GPU) 16792692e278SPaul Mullowney matrices whose operation should use a particular storage format. 1680e057df02SPaul Mullowney 1681e057df02SPaul Mullowney Not Collective 1682e057df02SPaul Mullowney 16838468deeeSKarl Rupp + MAT_CUSPARSE_MULT_DIAG - sets the storage format for the diagonal matrix in the parallel MatMult 16848468deeeSKarl Rupp . MAT_CUSPARSE_MULT_OFFDIAG - sets the storage format for the offdiagonal matrix in the parallel MatMult 16858468deeeSKarl Rupp . MAT_CUSPARSE_MULT - sets the storage format for the entire matrix in the serial (single GPU) MatMult 16868468deeeSKarl Rupp - MAT_CUSPARSE_ALL - sets the storage format for all CUSPARSE (GPU) matrices 1687e057df02SPaul Mullowney 1688e057df02SPaul Mullowney Level: intermediate 1689e057df02SPaul Mullowney 1690e057df02SPaul Mullowney .seealso: MatCUSPARSESetFormat(), MatCUSPARSEStorageFormat 1691e057df02SPaul Mullowney E*/ 169236d62e41SPaul Mullowney typedef enum {MAT_CUSPARSE_MULT_DIAG, MAT_CUSPARSE_MULT_OFFDIAG, MAT_CUSPARSE_MULT, MAT_CUSPARSE_ALL} MatCUSPARSEFormatOperation; 1693e057df02SPaul Mullowney 169410e9db80SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 169510e9db80SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCreateAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 1696e057df02SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCUSPARSESetFormat(Mat,MatCUSPARSEFormatOperation,MatCUSPARSEStorageFormat); 16979ae82921SPaul Mullowney #endif 16989ae82921SPaul Mullowney 169990c53902SBarry Smith #if defined(PETSC_HAVE_CUSP) 1700014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCUSP(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 1701014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateAIJCUSP(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 1702e057df02SPaul Mullowney 17033f9c0db1SPaul Mullowney /*E 1704e057df02SPaul Mullowney MatCUSPStorageFormat - indicates the storage format for CUSP (GPU) 170536d62e41SPaul Mullowney matrices. 1706e057df02SPaul Mullowney 1707e057df02SPaul Mullowney Not Collective 1708e057df02SPaul Mullowney 17098468deeeSKarl Rupp + MAT_CUSP_CSR - Compressed Sparse Row 17108468deeeSKarl Rupp . MAT_CUSP_DIA - Diagonal 17118468deeeSKarl Rupp - MAT_CUSP_ELL - Ellpack 1712e057df02SPaul Mullowney 1713e057df02SPaul Mullowney Level: intermediate 1714e057df02SPaul Mullowney 1715af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 1716e057df02SPaul Mullowney 1717e057df02SPaul Mullowney .seealso: MatCUSPSetFormat(), MatCUSPFormatOperation 1718e057df02SPaul Mullowney E*/ 17193f9c0db1SPaul Mullowney typedef enum {MAT_CUSP_CSR, MAT_CUSP_DIA, MAT_CUSP_ELL} MatCUSPStorageFormat; 1720e057df02SPaul Mullowney 1721e057df02SPaul Mullowney /* these will be strings associated with enumerated type defined above */ 1722e057df02SPaul Mullowney PETSC_EXTERN const char *const MatCUSPStorageFormats[]; 1723e057df02SPaul Mullowney 17243f9c0db1SPaul Mullowney /*E 1725e057df02SPaul Mullowney MatCUSPFormatOperation - indicates the operation of CUSP (GPU) 172636d62e41SPaul Mullowney matrices whose operation should use a particular storage format. 1727e057df02SPaul Mullowney 1728e057df02SPaul Mullowney Not Collective 1729e057df02SPaul Mullowney 17308468deeeSKarl Rupp + MAT_CUSP_MULT_DIAG - sets the storage format for the diagonal matrix in the parallel MatMult 17318468deeeSKarl Rupp . MAT_CUSP_MULT_OFFDIAG - sets the storage format for the offdiagonal matrix in the parallel MatMult 17328468deeeSKarl Rupp . MAT_CUSP_MULT - sets the storage format for the entire matrix in the serial (single GPU) MatMult 17338468deeeSKarl Rupp - MAT_CUSP_ALL - sets the storage format for all CUSP (GPU) matrices 1734e057df02SPaul Mullowney 1735e057df02SPaul Mullowney Level: intermediate 1736e057df02SPaul Mullowney 1737af0996ceSBarry Smith Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 1738e057df02SPaul Mullowney 1739e057df02SPaul Mullowney .seealso: MatCUSPSetFormat(), MatCUSPStorageFormat 1740e057df02SPaul Mullowney E*/ 174136d62e41SPaul Mullowney typedef enum {MAT_CUSP_MULT_DIAG, MAT_CUSP_MULT_OFFDIAG, MAT_CUSP_MULT, MAT_CUSP_ALL} MatCUSPFormatOperation; 1742e057df02SPaul Mullowney 1743e057df02SPaul Mullowney PETSC_EXTERN PetscErrorCode MatCUSPSetFormat(Mat,MatCUSPFormatOperation,MatCUSPStorageFormat); 1744e057df02SPaul Mullowney #endif 174590c53902SBarry Smith 1746d67ff14aSKarl Rupp #if defined(PETSC_HAVE_VIENNACL) 1747d67ff14aSKarl Rupp PETSC_EXTERN PetscErrorCode MatCreateSeqAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 1748d67ff14aSKarl Rupp PETSC_EXTERN PetscErrorCode MatCreateAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 1749d67ff14aSKarl Rupp #endif 1750d67ff14aSKarl Rupp 175154efbe56SHong Zhang /* 175254efbe56SHong Zhang PETSc interface to FFTW 175354efbe56SHong Zhang */ 175454efbe56SHong Zhang #if defined(PETSC_HAVE_FFTW) 1755014dd563SJed Brown PETSC_EXTERN PetscErrorCode VecScatterPetscToFFTW(Mat,Vec,Vec); 1756014dd563SJed Brown PETSC_EXTERN PetscErrorCode VecScatterFFTWToPetsc(Mat,Vec,Vec); 17572a7a6963SBarry Smith PETSC_EXTERN PetscErrorCode MatCreateVecsFFTW(Mat,Vec*,Vec*,Vec*); 175854efbe56SHong Zhang #endif 175954efbe56SHong Zhang 1760382fd914SHong Zhang /* 1761382fd914SHong Zhang PETSc interface to ELEMENTAL 1762382fd914SHong Zhang */ 1763db31f6deSJed Brown #if defined(PETSC_HAVE_ELEMENTAL) 17642ef15aa2SHong Zhang #if defined(PETSC_USE_COMPLEX) 17652ef15aa2SHong Zhang typedef El::Complex<PetscReal> PetscElemScalar; 17662ef15aa2SHong Zhang #else 17672ef15aa2SHong Zhang typedef PetscScalar PetscElemScalar; 1768382fd914SHong Zhang #endif 1769607a6623SBarry Smith PETSC_EXTERN PetscErrorCode PetscElementalInitializePackage(void); 1770db31f6deSJed Brown PETSC_EXTERN PetscErrorCode PetscElementalFinalizePackage(void); 1771db31f6deSJed Brown #endif 1772db31f6deSJed Brown 1773014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateNest(MPI_Comm,PetscInt,const IS[],PetscInt,const IS[],const Mat[],Mat*); 1774014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSize(Mat,PetscInt*,PetscInt*); 1775014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetISs(Mat,IS[],IS[]); 1776014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetLocalISs(Mat,IS[],IS[]); 1777014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSubMats(Mat,PetscInt*,PetscInt*,Mat***); 1778014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestGetSubMat(Mat,PetscInt,PetscInt,Mat*); 177919fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode MatNestSetVecType(Mat,VecType); 1780014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestSetSubMats(Mat,PetscInt,const IS[],PetscInt,const IS[],const Mat[]); 1781014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatNestSetSubMat(Mat,PetscInt,PetscInt,Mat); 1782d8588912SDave May 17834325cce7SMatthew G Knepley PETSC_EXTERN PetscErrorCode MatChop(Mat,PetscReal); 1784e0a58c10SMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatComputeBandwidth(Mat,PetscReal,PetscInt*); 17854325cce7SMatthew G Knepley 1786b541e6a4SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSubdomainsCreateCoalesce(Mat,PetscInt,PetscInt*,IS**); 178753dd7562SDmitry Karpeev 178823a5497aSJed Brown #endif 1789