12eac72dbSBarry Smith /* 22eac72dbSBarry Smith Include file for the matrix component of PETSc 32eac72dbSBarry Smith */ 40a835dfdSSatish Balay #ifndef __PETSCMAT_H 50a835dfdSSatish Balay #define __PETSCMAT_H 60a835dfdSSatish Balay #include "petscvec.h" 7e9fa29b7SSatish Balay PETSC_EXTERN_CXX_BEGIN 82eac72dbSBarry Smith 9d9274352SBarry Smith /*S 10d9274352SBarry Smith Mat - Abstract PETSc matrix object 112eac72dbSBarry Smith 12d91e6319SBarry Smith Level: beginner 13d91e6319SBarry Smith 14d9274352SBarry Smith Concepts: matrix; linear operator 15d9274352SBarry Smith 16d9274352SBarry Smith .seealso: MatCreate(), MatType, MatSetType() 17d9274352SBarry Smith S*/ 18d9274352SBarry Smith typedef struct _p_Mat* Mat; 19d9274352SBarry Smith 20d9274352SBarry Smith /*E 21d9274352SBarry Smith MatType - String with the name of a PETSc matrix or the creation function 22d9274352SBarry Smith with an optional dynamic library name, for example 23d9274352SBarry Smith http://www.mcs.anl.gov/petsc/lib.a:mymatcreate() 24d9274352SBarry Smith 25d9274352SBarry Smith Level: beginner 26d9274352SBarry Smith 27c7393fdbSBarry Smith .seealso: MatSetType(), Mat, MatSolverPackage 28d91e6319SBarry Smith E*/ 29a313700dSBarry Smith #define MatType char* 30273d9f13SBarry Smith #define MATSAME "same" 315a11e1b2SBarry Smith #define MATMAIJ "maij" 32273d9f13SBarry Smith #define MATSEQMAIJ "seqmaij" 33273d9f13SBarry Smith #define MATMPIMAIJ "mpimaij" 34273d9f13SBarry Smith #define MATIS "is" 355a11e1b2SBarry Smith #define MATAIJ "aij" 36273d9f13SBarry Smith #define MATSEQAIJ "seqaij" 3751d315f7SKerry Stevens #define MATSEQPTHREADAIJ "seqpthreadaij" 38273d9f13SBarry Smith #define MATMPIAIJ "mpiaij" 395a11e1b2SBarry Smith #define MATAIJCRL "aijcrl" 405a11e1b2SBarry Smith #define MATSEQAIJCRL "seqaijcrl" 415a11e1b2SBarry Smith #define MATMPIAIJCRL "mpiaijcrl" 428154be41SBarry Smith #define MATAIJCUSP "aijcusp" 438154be41SBarry Smith #define MATSEQAIJCUSP "seqaijcusp" 448154be41SBarry Smith #define MATMPIAIJCUSP "mpiaijcusp" 455a11e1b2SBarry Smith #define MATAIJPERM "aijperm" 465a11e1b2SBarry Smith #define MATSEQAIJPERM "seqaijperm" 475a11e1b2SBarry Smith #define MATMPIAIJPERM "mpiaijperm" 48273d9f13SBarry Smith #define MATSHELL "shell" 495a11e1b2SBarry Smith #define MATDENSE "dense" 50209238afSKris Buschelman #define MATSEQDENSE "seqdense" 51273d9f13SBarry Smith #define MATMPIDENSE "mpidense" 525a11e1b2SBarry Smith #define MATBAIJ "baij" 53273d9f13SBarry Smith #define MATSEQBAIJ "seqbaij" 54273d9f13SBarry Smith #define MATMPIBAIJ "mpibaij" 55273d9f13SBarry Smith #define MATMPIADJ "mpiadj" 565a11e1b2SBarry Smith #define MATSBAIJ "sbaij" 57273d9f13SBarry Smith #define MATSEQSBAIJ "seqsbaij" 58273d9f13SBarry Smith #define MATMPISBAIJ "mpisbaij" 59c0cdd4a1SDahai Guo 60c0cdd4a1SDahai Guo #define MATSEQBSTRM "seqbstrm" 61c0cdd4a1SDahai Guo #define MATMPIBSTRM "mpibstrm" 62c0cdd4a1SDahai Guo #define MATBSTRM "bstrm" 63aa5a9175SDahai Guo #define MATSEQSBSTRM "seqsbstrm" 64aa5a9175SDahai Guo #define MATMPISBSTRM "mpisbstrm" 65aa5a9175SDahai Guo #define MATSBSTRM "sbstrm" 66c0cdd4a1SDahai Guo 67cebc7f6cSBarry Smith #define MATDAAD "daad" 68cebc7f6cSBarry Smith #define MATMFFD "mffd" 69c8a8475eSBarry Smith #define MATNORMAL "normal" 70ab92ecdeSBarry Smith #define MATLRC "lrc" 712a6744ebSBarry Smith #define MATSCATTER "scatter" 72421e10b8SBarry Smith #define MATBLOCKMAT "blockmat" 73793850ffSBarry Smith #define MATCOMPOSITE "composite" 741f8c7532SHong Zhang #define MATFFT "fft" 751f8c7532SHong Zhang #define MATFFTW "fftw" 76e133240eSMatthew G Knepley #define MATSEQCUFFT "seqcufft" 77557cca28SSatish Balay #define MATTRANSPOSEMAT "transpose" 7872ca8751SBarry Smith #define MATSCHURCOMPLEMENT "schurcomplement" 791d6018f0SLisandro Dalcin #define MATPYTHON "python" 80f91d8e95SBarry Smith #define MATHYPRESTRUCT "hyprestruct" 81a9e6138eSGlenn Hammond #define MATHYPRESSTRUCT "hypresstruct" 82ee1cef2cSJed Brown #define MATSUBMATRIX "submatrix" 832c0dbf93SJed Brown #define MATLOCALREF "localref" 84d8588912SDave May #define MATNEST "nest" 85773941d6SBarry Smith 869989ab13SBarry Smith /*E 87c7393fdbSBarry Smith MatSolverPackage - String with the name of a PETSc matrix solver type. 889989ab13SBarry Smith 899989ab13SBarry Smith For example: "petsc" indicates what PETSc provides, "superlu" indicates either 909989ab13SBarry Smith SuperLU or SuperLU_Dist etc. 919989ab13SBarry Smith 929989ab13SBarry Smith 939989ab13SBarry Smith Level: beginner 949989ab13SBarry Smith 955c9eb25fSBarry Smith .seealso: MatGetFactor(), Mat, MatSetType(), MatType 969989ab13SBarry Smith E*/ 97c7393fdbSBarry Smith #define MatSolverPackage char* 982692d6eeSBarry Smith #define MATSOLVERSPOOLES "spooles" 992692d6eeSBarry Smith #define MATSOLVERSUPERLU "superlu" 1002692d6eeSBarry Smith #define MATSOLVERSUPERLU_DIST "superlu_dist" 1012692d6eeSBarry Smith #define MATSOLVERUMFPACK "umfpack" 10220db9a53SJed Brown #define MATSOLVERCHOLMOD "cholmod" 1032692d6eeSBarry Smith #define MATSOLVERESSL "essl" 1042692d6eeSBarry Smith #define MATSOLVERLUSOL "lusol" 1052692d6eeSBarry Smith #define MATSOLVERMUMPS "mumps" 1062692d6eeSBarry Smith #define MATSOLVERPASTIX "pastix" 1072692d6eeSBarry Smith #define MATSOLVERDSCPACK "dscpack" 1082692d6eeSBarry Smith #define MATSOLVERMATLAB "matlab" 1092692d6eeSBarry Smith #define MATSOLVERPETSC "petsc" 1102692d6eeSBarry Smith #define MATSOLVERPLAPACK "plapack" 1112692d6eeSBarry Smith #define MATSOLVERBAS "bas" 112773941d6SBarry Smith 113aa5a9175SDahai Guo #define MATSOLVERBSTRM "bstrm" 114aa5a9175SDahai Guo #define MATSOLVERSBSTRM "sbstrm" 115c0cdd4a1SDahai Guo 116b24902e0SBarry Smith /*E 117b24902e0SBarry Smith MatFactorType - indicates what type of factorization is requested 118b24902e0SBarry Smith 119b24902e0SBarry Smith Level: beginner 120b24902e0SBarry Smith 121b24902e0SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 122b24902e0SBarry Smith 123c7393fdbSBarry Smith .seealso: MatSolverPackage, MatGetFactor() 124b24902e0SBarry Smith E*/ 125599ef60dSHong Zhang typedef enum {MAT_FACTOR_NONE, MAT_FACTOR_LU, MAT_FACTOR_CHOLESKY, MAT_FACTOR_ILU, MAT_FACTOR_ICC,MAT_FACTOR_ILUDT} MatFactorType; 12634918c53SJed Brown extern const char *const MatFactorTypes[]; 127e92e720dSBarry Smith 1287087cfbeSBarry Smith extern PetscErrorCode MatGetFactor(Mat,const MatSolverPackage,MatFactorType,Mat*); 1297087cfbeSBarry Smith extern PetscErrorCode MatGetFactorAvailable(Mat,const MatSolverPackage,MatFactorType,PetscBool *); 1307087cfbeSBarry Smith extern PetscErrorCode MatFactorGetSolverPackage(Mat,const MatSolverPackage*); 1317087cfbeSBarry Smith extern PetscErrorCode MatGetFactorType(Mat,MatFactorType*); 1329989ab13SBarry Smith 133c06d978dSMatthew Knepley /* Logging support */ 1340700a824SBarry Smith #define MAT_FILE_CLASSID 1211216 /* used to indicate matrices in binary files */ 1357087cfbeSBarry Smith extern PetscClassId MAT_CLASSID; 1367087cfbeSBarry Smith extern PetscClassId MAT_FDCOLORING_CLASSID; 1377087cfbeSBarry Smith extern PetscClassId MAT_PARTITIONING_CLASSID; 1387087cfbeSBarry Smith extern PetscClassId MAT_NULLSPACE_CLASSID; 1397087cfbeSBarry Smith extern PetscClassId MATMFFD_CLASSID; 140c06d978dSMatthew Knepley 141ceb03754SKris Buschelman /*E 142ceb03754SKris Buschelman MatReuse - Indicates if matrices obtained from a previous call to MatGetSubMatrices() 143d6eff37eSBarry Smith or MatGetSubMatrix() are to be reused to store the new matrix values. For MatConvert() is used to indicate 144d6eff37eSBarry Smith that the input matrix is to be replaced with the converted matrix. 145ceb03754SKris Buschelman 146ceb03754SKris Buschelman Level: beginner 147ceb03754SKris Buschelman 148ceb03754SKris Buschelman Any additions/changes here MUST also be made in include/finclude/petscmat.h 149ceb03754SKris Buschelman 1500c451bc4SBarry Smith .seealso: MatGetSubMatrices(), MatGetSubMatrix(), MatDestroyMatrices(), MatConvert() 151ceb03754SKris Buschelman E*/ 152dfe085dbSJed Brown typedef enum {MAT_INITIAL_MATRIX,MAT_REUSE_MATRIX,MAT_IGNORE_MATRIX} MatReuse; 153ceb03754SKris Buschelman 1545494a064SHong Zhang /*E 1555494a064SHong Zhang MatGetSubMatrixOption - Indicates if matrices obtained from a call to MatGetSubMatrices() 156829201f2SHong Zhang include the matrix values. Currently it is only used by MatGetSeqNonzerostructure(). 1575494a064SHong Zhang 1585494a064SHong Zhang Level: beginner 1595494a064SHong Zhang 160829201f2SHong Zhang .seealso: MatGetSeqNonzerostructure() 1615494a064SHong Zhang E*/ 1625494a064SHong Zhang typedef enum {MAT_DO_NOT_GET_VALUES,MAT_GET_VALUES} MatGetSubMatrixOption; 1635494a064SHong Zhang 1647087cfbeSBarry Smith extern PetscErrorCode MatInitializePackage(const char[]); 165c06d978dSMatthew Knepley 1667087cfbeSBarry Smith extern PetscErrorCode MatCreate(MPI_Comm,Mat*); 167a69a960dSMatthew Knepley PetscPolymorphicFunction(MatCreate,(MPI_Comm comm),(comm,&A),Mat,A) 168a69a960dSMatthew Knepley PetscPolymorphicFunction(MatCreate,(),(PETSC_COMM_WORLD,&A),Mat,A) 1697087cfbeSBarry Smith extern PetscErrorCode MatSetSizes(Mat,PetscInt,PetscInt,PetscInt,PetscInt); 1707087cfbeSBarry Smith extern PetscErrorCode MatSetType(Mat,const MatType); 1717087cfbeSBarry Smith extern PetscErrorCode MatSetFromOptions(Mat); 1727087cfbeSBarry Smith extern PetscErrorCode MatSetUpPreallocation(Mat); 1737087cfbeSBarry Smith extern PetscErrorCode MatRegisterAll(const char[]); 1747087cfbeSBarry Smith extern PetscErrorCode MatRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat)); 17501bebe75SBarry Smith extern PetscErrorCode MatRegisterBaseName(const char[],const char[],const char[]); 1767087cfbeSBarry Smith extern PetscErrorCode MatSetOptionsPrefix(Mat,const char[]); 1777087cfbeSBarry Smith extern PetscErrorCode MatAppendOptionsPrefix(Mat,const char[]); 1787087cfbeSBarry Smith extern PetscErrorCode MatGetOptionsPrefix(Mat,const char*[]); 179f69a0ea3SMatthew Knepley 18030de9b25SBarry Smith /*MC 18130de9b25SBarry Smith MatRegisterDynamic - Adds a new matrix type 18230de9b25SBarry Smith 18330de9b25SBarry Smith Synopsis: 1841890ba74SBarry Smith PetscErrorCode MatRegisterDynamic(const char *name,const char *path,const char *name_create,PetscErrorCode (*routine_create)(Mat)) 18530de9b25SBarry Smith 18630de9b25SBarry Smith Not Collective 18730de9b25SBarry Smith 18830de9b25SBarry Smith Input Parameters: 18930de9b25SBarry Smith + name - name of a new user-defined matrix type 19030de9b25SBarry Smith . path - path (either absolute or relative) the library containing this solver 19130de9b25SBarry Smith . name_create - name of routine to create method context 19230de9b25SBarry Smith - routine_create - routine to create method context 19330de9b25SBarry Smith 19430de9b25SBarry Smith Notes: 19530de9b25SBarry Smith MatRegisterDynamic() may be called multiple times to add several user-defined solvers. 19630de9b25SBarry Smith 19730de9b25SBarry Smith If dynamic libraries are used, then the fourth input argument (routine_create) 19830de9b25SBarry Smith is ignored. 19930de9b25SBarry Smith 20030de9b25SBarry Smith Sample usage: 20130de9b25SBarry Smith .vb 20230de9b25SBarry Smith MatRegisterDynamic("my_mat",/home/username/my_lib/lib/libO/solaris/mylib.a, 20330de9b25SBarry Smith "MyMatCreate",MyMatCreate); 20430de9b25SBarry Smith .ve 20530de9b25SBarry Smith 20630de9b25SBarry Smith Then, your solver can be chosen with the procedural interface via 20730de9b25SBarry Smith $ MatSetType(Mat,"my_mat") 20830de9b25SBarry Smith or at runtime via the option 20930de9b25SBarry Smith $ -mat_type my_mat 21030de9b25SBarry Smith 21130de9b25SBarry Smith Level: advanced 21230de9b25SBarry Smith 213ab901514SBarry Smith Notes: ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values. 21430de9b25SBarry Smith If your function is not being put into a shared library then use VecRegister() instead 21530de9b25SBarry Smith 21630de9b25SBarry Smith .keywords: Mat, register 21730de9b25SBarry Smith 21830de9b25SBarry Smith .seealso: MatRegisterAll(), MatRegisterDestroy() 21930de9b25SBarry Smith 22030de9b25SBarry Smith M*/ 221273d9f13SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 222273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,0) 223273d9f13SBarry Smith #else 224273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,d) 22530de9b25SBarry Smith #endif 22630de9b25SBarry Smith 227ace3abfcSBarry Smith extern PetscBool MatRegisterAllCalled; 228b0a32e0cSBarry Smith extern PetscFList MatList; 229b022a5c1SBarry Smith extern PetscFList MatColoringList; 230b022a5c1SBarry Smith extern PetscFList MatPartitioningList; 23128988994SBarry Smith 2323b224e63SBarry Smith /*E 2333b224e63SBarry Smith MatStructure - Indicates if the matrix has the same nonzero structure 2343b224e63SBarry Smith 2353b224e63SBarry Smith Level: beginner 2363b224e63SBarry Smith 2373b224e63SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 2383b224e63SBarry Smith 2393b224e63SBarry Smith .seealso: MatCopy(), KSPSetOperators(), PCSetOperators() 2403b224e63SBarry Smith E*/ 241214bc6a2SJed Brown typedef enum {DIFFERENT_NONZERO_PATTERN,SUBSET_NONZERO_PATTERN,SAME_NONZERO_PATTERN,SAME_PRECONDITIONER} MatStructure; 2423b224e63SBarry Smith 2437087cfbeSBarry Smith extern PetscErrorCode MatCreateSeqDense(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*); 2447087cfbeSBarry Smith extern PetscErrorCode MatCreateMPIDense(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*); 2457087cfbeSBarry Smith extern PetscErrorCode MatCreateSeqAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 246ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,nz,nnz,&A),Mat,A) 247ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,&A),Mat,A) 248ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,0,nnz,&A),Mat,A) 249ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,&A),Mat,A) 250ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,A)) 251ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,m,n,0,nnz,A)) 252ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,A)) 2537087cfbeSBarry Smith extern PetscErrorCode MatCreateMPIAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 254ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPIAIJ,(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,const PetscInt nnz[],PetscInt onz,const PetscInt onnz[]),(comm,m,n,M,N,nz,nnz,onz,onnz,&A),Mat,A) 255ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPIAIJ,(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,PetscInt nnz),(comm,m,n,M,N,nz,PETSC_NULL,nnz,PETSC_NULL,&A),Mat,A) 256ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPIAIJ,(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt nnz[],const PetscInt onz[]),(comm,m,n,M,N,0,nnz,0,onz,&A),Mat,A) 257ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPIAIJ,(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N),(comm,m,n,M,N,0,PETSC_NULL,0,PETSC_NULL,&A),Mat,A) 258ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPIAIJ,(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,PetscInt nnz,Mat *A),(comm,m,n,M,N,nz,PETSC_NULL,nnz,PETSC_NULL,A)) 259ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPIAIJ,(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt nnz[],const PetscInt onz[],Mat *A),(comm,m,n,M,N,0,nnz,0,onz,A)) 260ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPIAIJ,(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,Mat *A),(comm,m,n,M,N,0,PETSC_NULL,0,PETSC_NULL,A)) 261ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPIAIJ,(PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,const PetscInt nnz[],PetscInt onz,const PetscInt onnz[]),(PETSC_COMM_WORLD,m,n,M,N,nz,nnz,onz,onnz,&A),Mat,A) 262ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPIAIJ,(PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,PetscInt nnz),(PETSC_COMM_WORLD,m,n,M,N,nz,PETSC_NULL,nnz,PETSC_NULL,&A),Mat,A) 263ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPIAIJ,(PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt nnz[],const PetscInt onz[]),(PETSC_COMM_WORLD,m,n,M,N,0,nnz,0,onz,&A),Mat,A) 264ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPIAIJ,(PetscInt m,PetscInt n,PetscInt M,PetscInt N),(PETSC_COMM_WORLD,m,n,M,N,0,PETSC_NULL,0,PETSC_NULL,&A),Mat,A) 265ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPIAIJ,(PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,PetscInt nnz,Mat *A),(PETSC_COMM_WORLD,m,n,M,N,nz,PETSC_NULL,nnz,PETSC_NULL,A)) 266ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPIAIJ,(PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt nnz[],const PetscInt onz[],Mat *A),(PETSC_COMM_WORLD,m,n,M,N,0,nnz,0,onz,A)) 267ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPIAIJ,(PetscInt m,PetscInt n,PetscInt M,PetscInt N,Mat *A),(PETSC_COMM_WORLD,m,n,M,N,0,PETSC_NULL,0,PETSC_NULL,A)) 2687087cfbeSBarry Smith extern PetscErrorCode MatCreateMPIAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *); 2697087cfbeSBarry Smith extern PetscErrorCode MatCreateMPIAIJWithSplitArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],PetscInt[],PetscInt[],PetscScalar[],Mat*); 2708d7a6e47SBarry Smith 2717087cfbeSBarry Smith extern PetscErrorCode MatCreateSeqBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 272ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,nz,nnz,&A),Mat,A) 273ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A) 274ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A) 275ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A) 276ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A)) 277ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A)) 278ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A)) 2797087cfbeSBarry Smith extern PetscErrorCode MatCreateMPIBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 280ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPIBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,const PetscInt nnz[],PetscInt onz,const PetscInt onnz[]),(comm,bs,m,n,M,N,nz,nnz,onz,onnz,&A),Mat,A) 281ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPIBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,PetscInt nnz),(comm,bs,m,n,M,N,nz,PETSC_NULL,nnz,PETSC_NULL,&A),Mat,A) 282ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPIBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt nnz[],const PetscInt onz[]),(comm,bs,m,n,M,N,0,nnz,0,onz,&A),Mat,A) 283ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPIBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N),(comm,bs,m,n,M,N,0,PETSC_NULL,0,PETSC_NULL,&A),Mat,A) 284ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPIBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,PetscInt nnz,Mat *A),(comm,bs,m,n,M,N,nz,PETSC_NULL,nnz,PETSC_NULL,A)) 285ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPIBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt nnz[],const PetscInt onz[],Mat *A),(comm,bs,m,n,M,N,0,nnz,0,onz,A)) 286ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPIBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,Mat *A),(comm,bs,m,n,M,N,0,PETSC_NULL,0,PETSC_NULL,A)) 287ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPIBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,const PetscInt nnz[],PetscInt onz,const PetscInt onnz[]),(PETSC_COMM_WORLD,bs,m,n,M,N,nz,nnz,onz,onnz,&A),Mat,A) 288ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPIBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,PetscInt nnz),(PETSC_COMM_WORLD,bs,m,n,M,N,nz,PETSC_NULL,nnz,PETSC_NULL,&A),Mat,A) 289ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPIBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt nnz[],const PetscInt onz[]),(PETSC_COMM_WORLD,bs,m,n,M,N,0,nnz,0,onz,&A),Mat,A) 290ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPIBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N),(PETSC_COMM_WORLD,bs,m,n,M,N,0,PETSC_NULL,0,PETSC_NULL,&A),Mat,A) 291ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPIBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,PetscInt nnz,Mat *A),(PETSC_COMM_WORLD,bs,m,n,M,N,nz,PETSC_NULL,nnz,PETSC_NULL,A)) 292ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPIBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt nnz[],const PetscInt onz[],Mat *A),(PETSC_COMM_WORLD,bs,m,n,M,N,0,nnz,0,onz,A)) 293ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPIBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,Mat *A),(PETSC_COMM_WORLD,bs,m,n,M,N,0,PETSC_NULL,0,PETSC_NULL,A)) 2947087cfbeSBarry Smith extern PetscErrorCode MatCreateMPIBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat*); 295d21a29f3SJed Brown 2967087cfbeSBarry Smith extern PetscErrorCode MatCreateMPIAdj(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscInt[],Mat*); 2977087cfbeSBarry Smith extern PetscErrorCode MatCreateSeqSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 298ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,nz,nnz,&A),Mat,A) 299ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A) 300ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A) 301ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A) 302ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A)) 303ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A)) 304ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A)) 305d21a29f3SJed Brown 3067087cfbeSBarry Smith extern PetscErrorCode MatCreateMPISBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 307ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPISBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,const PetscInt nnz[],PetscInt onz,const PetscInt onnz[]),(comm,bs,m,n,M,N,nz,nnz,onz,onnz,&A),Mat,A) 308ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPISBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,PetscInt nnz),(comm,bs,m,n,M,N,nz,PETSC_NULL,nnz,PETSC_NULL,&A),Mat,A) 309ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPISBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt nnz[],const PetscInt onz[]),(comm,bs,m,n,M,N,0,nnz,0,onz,&A),Mat,A) 310ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPISBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N),(comm,bs,m,n,M,N,0,PETSC_NULL,0,PETSC_NULL,&A),Mat,A) 311ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPISBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,PetscInt nnz,Mat *A),(comm,bs,m,n,M,N,nz,PETSC_NULL,nnz,PETSC_NULL,A)) 312ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPISBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt nnz[],const PetscInt onz[],Mat *A),(comm,bs,m,n,M,N,0,nnz,0,onz,A)) 313ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPISBAIJ,(MPI_Comm comm,PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,Mat *A),(comm,bs,m,n,M,N,0,PETSC_NULL,0,PETSC_NULL,A)) 314ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPISBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,const PetscInt nnz[],PetscInt onz,const PetscInt onnz[]),(PETSC_COMM_WORLD,bs,m,n,M,N,nz,nnz,onz,onnz,&A),Mat,A) 315ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPISBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,PetscInt nnz),(PETSC_COMM_WORLD,bs,m,n,M,N,nz,PETSC_NULL,nnz,PETSC_NULL,&A),Mat,A) 316ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPISBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt nnz[],const PetscInt onz[]),(PETSC_COMM_WORLD,bs,m,n,M,N,0,nnz,0,onz,&A),Mat,A) 317ba966639SSatish Balay PetscPolymorphicFunction(MatCreateMPISBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N),(PETSC_COMM_WORLD,bs,m,n,M,N,0,PETSC_NULL,0,PETSC_NULL,&A),Mat,A) 318ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPISBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt nz,PetscInt nnz,Mat *A),(PETSC_COMM_WORLD,bs,m,n,M,N,nz,PETSC_NULL,nnz,PETSC_NULL,A)) 319ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPISBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,const PetscInt nnz[],const PetscInt onz[],Mat *A),(PETSC_COMM_WORLD,bs,m,n,M,N,0,nnz,0,onz,A)) 320ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateMPISBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt M,PetscInt N,Mat *A),(PETSC_COMM_WORLD,bs,m,n,M,N,0,PETSC_NULL,0,PETSC_NULL,A)) 3217087cfbeSBarry Smith extern PetscErrorCode MatCreateMPISBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *); 3227087cfbeSBarry Smith extern PetscErrorCode MatMPISBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 323dfb205c3SBarry Smith 3247087cfbeSBarry Smith extern PetscErrorCode MatCreateShell(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,void *,Mat*); 325ba966639SSatish Balay PetscPolymorphicFunction(MatCreateShell,(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx),(comm,m,n,M,N,ctx,&A),Mat,A) 326ba966639SSatish Balay PetscPolymorphicFunction(MatCreateShell,(PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx),(PETSC_COMM_WORLD,m,n,M,N,ctx,&A),Mat,A) 3277087cfbeSBarry Smith extern PetscErrorCode MatCreateAdic(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,void (*)(void),Mat*); 3287087cfbeSBarry Smith extern PetscErrorCode MatCreateNormal(Mat,Mat*); 329ba966639SSatish Balay PetscPolymorphicFunction(MatCreateNormal,(Mat mat),(mat,&A),Mat,A) 3307087cfbeSBarry Smith extern PetscErrorCode MatCreateLRC(Mat,Mat,Mat,Mat*); 3317087cfbeSBarry Smith extern PetscErrorCode MatCreateIS(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,ISLocalToGlobalMapping,Mat*); 3327087cfbeSBarry Smith extern PetscErrorCode MatCreateSeqAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 3337087cfbeSBarry Smith extern PetscErrorCode MatCreateMPIAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 334c0cdd4a1SDahai Guo 335c0cdd4a1SDahai Guo extern PetscErrorCode MatCreateSeqBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 336c0cdd4a1SDahai Guo extern PetscErrorCode MatCreateMPIBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 337aa5a9175SDahai Guo extern PetscErrorCode MatCreateSeqSBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 338aa5a9175SDahai Guo extern PetscErrorCode MatCreateMPISBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 339c0cdd4a1SDahai Guo 3407087cfbeSBarry Smith extern PetscErrorCode MatCreateScatter(MPI_Comm,VecScatter,Mat*); 3417087cfbeSBarry Smith extern PetscErrorCode MatScatterSetVecScatter(Mat,VecScatter); 3427087cfbeSBarry Smith extern PetscErrorCode MatScatterGetVecScatter(Mat,VecScatter*); 3437087cfbeSBarry Smith extern PetscErrorCode MatCreateBlockMat(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,Mat*); 3447087cfbeSBarry Smith extern PetscErrorCode MatCompositeAddMat(Mat,Mat); 3457087cfbeSBarry Smith extern PetscErrorCode MatCompositeMerge(Mat); 3467087cfbeSBarry Smith extern PetscErrorCode MatCreateComposite(MPI_Comm,PetscInt,const Mat*,Mat*); 3476d7c1e57SBarry Smith typedef enum {MAT_COMPOSITE_ADDITIVE,MAT_COMPOSITE_MULTIPLICATIVE} MatCompositeType; 3487087cfbeSBarry Smith extern PetscErrorCode MatCompositeSetType(Mat,MatCompositeType); 3496d7c1e57SBarry Smith 350dedccee8SHong Zhang extern PetscErrorCode MatCreateFFT(MPI_Comm,PetscInt,const PetscInt[],const MatType,Mat*); 3517087cfbeSBarry Smith extern PetscErrorCode MatCreateSeqCUFFT(MPI_Comm,PetscInt,const PetscInt[],Mat*); 352dedccee8SHong Zhang 3537087cfbeSBarry Smith extern PetscErrorCode MatCreateTranspose(Mat,Mat*); 3547087cfbeSBarry Smith extern PetscErrorCode MatCreateSubMatrix(Mat,IS,IS,Mat*); 3557087cfbeSBarry Smith extern PetscErrorCode MatSubMatrixUpdate(Mat,Mat,IS,IS); 3567087cfbeSBarry Smith extern PetscErrorCode MatCreateLocalRef(Mat,IS,IS,Mat*); 3571472f72bSBarry Smith 3587087cfbeSBarry Smith extern PetscErrorCode MatPythonSetType(Mat,const char[]); 3591d6018f0SLisandro Dalcin 3607087cfbeSBarry Smith extern PetscErrorCode MatSetUp(Mat); 361fcfd50ebSBarry Smith extern PetscErrorCode MatDestroy(Mat*); 36221c89e3eSBarry Smith 3637087cfbeSBarry Smith extern PetscErrorCode MatConjugate(Mat); 3647087cfbeSBarry Smith extern PetscErrorCode MatRealPart(Mat); 3657087cfbeSBarry Smith extern PetscErrorCode MatImaginaryPart(Mat); 36611bd1e4dSLisandro Dalcin extern PetscErrorCode MatGetDiagonalBlock(Mat,Mat*); 3677087cfbeSBarry Smith extern PetscErrorCode MatGetTrace(Mat,PetscScalar*); 368*bbead8a2SBarry Smith extern PetscErrorCode MatInvertBlockDiagonal(Mat,PetscScalar **); 36999cafbc1SBarry Smith 3708ed539a5SBarry Smith /* ------------------------------------------------------------*/ 3717087cfbeSBarry Smith extern PetscErrorCode MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 3727087cfbeSBarry Smith extern PetscErrorCode MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 3737087cfbeSBarry Smith extern PetscErrorCode MatSetValuesRow(Mat,PetscInt,const PetscScalar[]); 3747087cfbeSBarry Smith extern PetscErrorCode MatSetValuesRowLocal(Mat,PetscInt,const PetscScalar[]); 37584cb2905SBarry Smith 3762ef4de8bSBarry Smith /*S 3772ef4de8bSBarry Smith MatStencil - Data structure (C struct) for storing information about a single row or 3782ef4de8bSBarry Smith column of a matrix as index on an associated grid. 3792ef4de8bSBarry Smith 3802ef4de8bSBarry Smith Level: beginner 3812ef4de8bSBarry Smith 3822ef4de8bSBarry Smith Concepts: matrix; linear operator 3832ef4de8bSBarry Smith 384d7bd93baSBarry Smith .seealso: MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockStencil() 3852ef4de8bSBarry Smith S*/ 386435da068SBarry Smith typedef struct { 387c1ac3661SBarry Smith PetscInt k,j,i,c; 388435da068SBarry Smith } MatStencil; 3892ef4de8bSBarry Smith 3907087cfbeSBarry Smith extern PetscErrorCode MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode); 3917087cfbeSBarry Smith extern PetscErrorCode MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode); 3927087cfbeSBarry Smith extern PetscErrorCode MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt); 393435da068SBarry Smith 3947087cfbeSBarry Smith extern PetscErrorCode MatSetColoring(Mat,ISColoring); 3957087cfbeSBarry Smith extern PetscErrorCode MatSetValuesAdic(Mat,void*); 3967087cfbeSBarry Smith extern PetscErrorCode MatSetValuesAdifor(Mat,PetscInt,void*); 3973a7fca6bSBarry Smith 398d91e6319SBarry Smith /*E 399d91e6319SBarry Smith MatAssemblyType - Indicates if the matrix is now to be used, or if you plan 400d91e6319SBarry Smith to continue to add values to it 401d91e6319SBarry Smith 402d91e6319SBarry Smith Level: beginner 403d91e6319SBarry Smith 404d91e6319SBarry Smith .seealso: MatAssemblyBegin(), MatAssemblyEnd() 405d91e6319SBarry Smith E*/ 4066d4a8577SBarry Smith typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType; 4077087cfbeSBarry Smith extern PetscErrorCode MatAssemblyBegin(Mat,MatAssemblyType); 4087087cfbeSBarry Smith extern PetscErrorCode MatAssemblyEnd(Mat,MatAssemblyType); 4097087cfbeSBarry Smith extern PetscErrorCode MatAssembled(Mat,PetscBool *); 4104f9c727eSBarry Smith 4111d73ed98SBarry Smith 41230de9b25SBarry Smith 413d91e6319SBarry Smith /*E 414d91e6319SBarry Smith MatOption - Options that may be set for a matrix and its behavior or storage 415d91e6319SBarry Smith 416d91e6319SBarry Smith Level: beginner 417d91e6319SBarry Smith 4180a835dfdSSatish Balay Any additions/changes here MUST also be made in include/finclude/petscmat.h 419d91e6319SBarry Smith 420d91e6319SBarry Smith .seealso: MatSetOption() 421d91e6319SBarry Smith E*/ 422512a5fc5SBarry Smith typedef enum {MAT_ROW_ORIENTED,MAT_NEW_NONZERO_LOCATIONS, 4234e0d8c25SBarry Smith MAT_SYMMETRIC, 4244e0d8c25SBarry Smith MAT_STRUCTURALLY_SYMMETRIC, 4258047e77bSBarry Smith MAT_NEW_DIAGONALS,MAT_IGNORE_OFF_PROC_ENTRIES, 4264e0d8c25SBarry Smith MAT_NEW_NONZERO_LOCATION_ERR, 4274e0d8c25SBarry Smith MAT_NEW_NONZERO_ALLOCATION_ERR,MAT_USE_HASH_TABLE, 428a9817697SBarry Smith MAT_KEEP_NONZERO_PATTERN,MAT_IGNORE_ZERO_ENTRIES, 4294e0d8c25SBarry Smith MAT_USE_INODES, 4304e0d8c25SBarry Smith MAT_HERMITIAN, 4314e0d8c25SBarry Smith MAT_SYMMETRY_ETERNAL, 432cd6b891eSBarry Smith MAT_CHECK_COMPRESSED_ROW, 4334e0d8c25SBarry Smith MAT_IGNORE_LOWER_TRIANGULAR,MAT_ERROR_LOWER_TRIANGULAR, 43428b2fa4aSMatthew Knepley MAT_GETROW_UPPERTRIANGULAR,MAT_UNUSED_NONZERO_LOCATION_ERR, 4354cb17eb5SBarry Smith MAT_SPD,MAT_NO_OFF_PROC_ENTRIES,MAT_NO_OFF_PROC_ZERO_ROWS, 43628b2fa4aSMatthew Knepley NUM_MAT_OPTIONS} MatOption; 437290bbb0aSBarry Smith extern const char *MatOptions[]; 4387087cfbeSBarry Smith extern PetscErrorCode MatSetOption(Mat,MatOption,PetscBool ); 4397087cfbeSBarry Smith extern PetscErrorCode MatGetType(Mat,const MatType*); 440a313700dSBarry Smith PetscPolymorphicFunction(MatGetType,(Mat mat),(mat,&t),const MatType,t) 44184cb2905SBarry Smith 4427087cfbeSBarry Smith extern PetscErrorCode MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]); 4437087cfbeSBarry Smith extern PetscErrorCode MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 4447087cfbeSBarry Smith extern PetscErrorCode MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 4457087cfbeSBarry Smith extern PetscErrorCode MatGetRowUpperTriangular(Mat); 4467087cfbeSBarry Smith extern PetscErrorCode MatRestoreRowUpperTriangular(Mat); 4477087cfbeSBarry Smith extern PetscErrorCode MatGetColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 4487087cfbeSBarry Smith extern PetscErrorCode MatRestoreColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 4497087cfbeSBarry Smith extern PetscErrorCode MatGetColumnVector(Mat,Vec,PetscInt); 4507087cfbeSBarry Smith extern PetscErrorCode MatGetArray(Mat,PetscScalar *[]); 451ba966639SSatish Balay PetscPolymorphicFunction(MatGetArray,(Mat mat),(mat,&a),PetscScalar*,a) 4527087cfbeSBarry Smith extern PetscErrorCode MatRestoreArray(Mat,PetscScalar *[]); 4537087cfbeSBarry Smith extern PetscErrorCode MatGetBlockSize(Mat,PetscInt *); 454ba966639SSatish Balay PetscPolymorphicFunction(MatGetBlockSize,(Mat mat),(mat,&a),PetscInt,a) 4557087cfbeSBarry Smith extern PetscErrorCode MatSetBlockSize(Mat,PetscInt); 4567b80b807SBarry Smith 4571620fd73SBarry Smith 4587087cfbeSBarry Smith extern PetscErrorCode MatMult(Mat,Vec,Vec); 4597087cfbeSBarry Smith extern PetscErrorCode MatMultDiagonalBlock(Mat,Vec,Vec); 4607087cfbeSBarry Smith extern PetscErrorCode MatMultAdd(Mat,Vec,Vec,Vec); 461ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultAdd,(Mat A,Vec x,Vec y),(A,x,y,y)) 4627087cfbeSBarry Smith extern PetscErrorCode MatMultTranspose(Mat,Vec,Vec); 4637087cfbeSBarry Smith extern PetscErrorCode MatMultHermitianTranspose(Mat,Vec,Vec); 4647087cfbeSBarry Smith extern PetscErrorCode MatIsTranspose(Mat,Mat,PetscReal,PetscBool *); 465ace3abfcSBarry Smith PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B,PetscReal tol),(A,B,tol,&t),PetscBool ,t) 466ace3abfcSBarry Smith PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B),(A,B,0,&t),PetscBool ,t) 4677087cfbeSBarry Smith extern PetscErrorCode MatIsHermitianTranspose(Mat,Mat,PetscReal,PetscBool *); 4687087cfbeSBarry Smith extern PetscErrorCode MatMultTransposeAdd(Mat,Vec,Vec,Vec); 4697087cfbeSBarry Smith extern PetscErrorCode MatMultHermitianTransposeAdd(Mat,Vec,Vec,Vec); 470ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultTransposeAdd,(Mat A,Vec x,Vec y),(A,x,y,y)) 4717087cfbeSBarry Smith extern PetscErrorCode MatMultConstrained(Mat,Vec,Vec); 4727087cfbeSBarry Smith extern PetscErrorCode MatMultTransposeConstrained(Mat,Vec,Vec); 4737087cfbeSBarry Smith extern PetscErrorCode MatMatSolve(Mat,Mat,Mat); 474f5edf698SHong Zhang 475d91e6319SBarry Smith /*E 476d91e6319SBarry Smith MatDuplicateOption - Indicates if a duplicated sparse matrix should have 477d91e6319SBarry Smith its numerical values copied over or just its nonzero structure. 478d91e6319SBarry Smith 479d91e6319SBarry Smith Level: beginner 480d91e6319SBarry Smith 481d91e6319SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 482d91e6319SBarry Smith 48370dcbbb9SBarry Smith $ MAT_SHARE_NONZERO_PATTERN - the i and j arrays in the new matrix will be shared with the original matrix 48470dcbbb9SBarry Smith $ this also triggers the MAT_DO_NOT_COPY_VALUES option. This is used when you 48570dcbbb9SBarry Smith $ have several matrices with the same nonzero pattern. 48670dcbbb9SBarry Smith 487d91e6319SBarry Smith .seealso: MatDuplicate() 488d91e6319SBarry Smith E*/ 48970dcbbb9SBarry Smith typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES,MAT_SHARE_NONZERO_PATTERN} MatDuplicateOption; 4902e8a6d31SBarry Smith 4917087cfbeSBarry Smith extern PetscErrorCode MatConvert(Mat,const MatType,MatReuse,Mat*); 492a313700dSBarry Smith PetscPolymorphicFunction(MatConvert,(Mat A,const MatType t),(A,t,MAT_INITIAL_MATRIX,&a),Mat,a) 4937087cfbeSBarry Smith extern PetscErrorCode MatDuplicate(Mat,MatDuplicateOption,Mat*); 494ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A,MatDuplicateOption o),(A,o,&a),Mat,a) 495ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A),(A,MAT_COPY_VALUES,&a),Mat,a) 49694a9d846SBarry Smith 497cb5b572fSBarry Smith 4987087cfbeSBarry Smith extern PetscErrorCode MatCopy(Mat,Mat,MatStructure); 4997087cfbeSBarry Smith extern PetscErrorCode MatView(Mat,PetscViewer); 5007087cfbeSBarry Smith extern PetscErrorCode MatIsSymmetric(Mat,PetscReal,PetscBool *); 501ace3abfcSBarry Smith PetscPolymorphicFunction(MatIsSymmetric,(Mat A,PetscReal tol),(A,tol,&t),PetscBool ,t) 502ace3abfcSBarry Smith PetscPolymorphicFunction(MatIsSymmetric,(Mat A),(A,0,&t),PetscBool ,t) 5037087cfbeSBarry Smith extern PetscErrorCode MatIsStructurallySymmetric(Mat,PetscBool *); 504ace3abfcSBarry Smith PetscPolymorphicFunction(MatIsStructurallySymmetric,(Mat A),(A,&t),PetscBool ,t) 5057087cfbeSBarry Smith extern PetscErrorCode MatIsHermitian(Mat,PetscReal,PetscBool *); 506ace3abfcSBarry Smith PetscPolymorphicFunction(MatIsHermitian,(Mat A),(A,0,&t),PetscBool ,t) 5077087cfbeSBarry Smith extern PetscErrorCode MatIsSymmetricKnown(Mat,PetscBool *,PetscBool *); 5087087cfbeSBarry Smith extern PetscErrorCode MatIsHermitianKnown(Mat,PetscBool *,PetscBool *); 5097087cfbeSBarry Smith extern PetscErrorCode MatMissingDiagonal(Mat,PetscBool *,PetscInt *); 5107087cfbeSBarry Smith extern PetscErrorCode MatLoad(Mat, PetscViewer); 5117b80b807SBarry Smith 5127087cfbeSBarry Smith extern PetscErrorCode MatGetRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,PetscInt *[],PetscInt *[],PetscBool *); 5137087cfbeSBarry Smith extern PetscErrorCode MatRestoreRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,PetscInt *[],PetscInt *[],PetscBool *); 5147087cfbeSBarry Smith extern PetscErrorCode MatGetColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,PetscInt *[],PetscInt *[],PetscBool *); 5157087cfbeSBarry Smith extern PetscErrorCode MatRestoreColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,PetscInt *[],PetscInt *[],PetscBool *); 516d4fbbf0eSBarry Smith 517d91e6319SBarry Smith /*S 518d91e6319SBarry Smith MatInfo - Context of matrix information, used with MatGetInfo() 519d91e6319SBarry Smith 520d91e6319SBarry Smith In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE 521d91e6319SBarry Smith 522d91e6319SBarry Smith Level: intermediate 523d91e6319SBarry Smith 524d91e6319SBarry Smith Concepts: matrix^nonzero information 525d91e6319SBarry Smith 526d9274352SBarry Smith .seealso: MatGetInfo(), MatInfoType 527d91e6319SBarry Smith S*/ 5284e220ebcSLois Curfman McInnes typedef struct { 529b0a32e0cSBarry Smith PetscLogDouble block_size; /* block size */ 530b0a32e0cSBarry Smith PetscLogDouble nz_allocated,nz_used,nz_unneeded; /* number of nonzeros */ 531b0a32e0cSBarry Smith PetscLogDouble memory; /* memory allocated */ 532b0a32e0cSBarry Smith PetscLogDouble assemblies; /* number of matrix assemblies called */ 533b0a32e0cSBarry Smith PetscLogDouble mallocs; /* number of mallocs during MatSetValues() */ 534b0a32e0cSBarry Smith PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */ 535b0a32e0cSBarry Smith PetscLogDouble factor_mallocs; /* number of mallocs during factorization */ 5364e220ebcSLois Curfman McInnes } MatInfo; 5374e220ebcSLois Curfman McInnes 538d9274352SBarry Smith /*E 539d9274352SBarry Smith MatInfoType - Indicates if you want information about the local part of the matrix, 540d9274352SBarry Smith the entire parallel matrix or the maximum over all the local parts. 541d9274352SBarry Smith 542d9274352SBarry Smith Level: beginner 543d9274352SBarry Smith 544d9274352SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 545d9274352SBarry Smith 546d9274352SBarry Smith .seealso: MatGetInfo(), MatInfo 547d9274352SBarry Smith E*/ 5487b80b807SBarry Smith typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType; 5497087cfbeSBarry Smith extern PetscErrorCode MatGetInfo(Mat,MatInfoType,MatInfo*); 5507087cfbeSBarry Smith extern PetscErrorCode MatGetDiagonal(Mat,Vec); 5517087cfbeSBarry Smith extern PetscErrorCode MatGetRowMax(Mat,Vec,PetscInt[]); 5527087cfbeSBarry Smith extern PetscErrorCode MatGetRowMin(Mat,Vec,PetscInt[]); 5537087cfbeSBarry Smith extern PetscErrorCode MatGetRowMaxAbs(Mat,Vec,PetscInt[]); 5547087cfbeSBarry Smith extern PetscErrorCode MatGetRowMinAbs(Mat,Vec,PetscInt[]); 5557087cfbeSBarry Smith extern PetscErrorCode MatGetRowSum(Mat,Vec); 5567087cfbeSBarry Smith extern PetscErrorCode MatTranspose(Mat,MatReuse,Mat*); 557fc4dec0aSBarry Smith PetscPolymorphicFunction(MatTranspose,(Mat A),(A,MAT_INITIAL_MATRIX,&t),Mat,t) 5587087cfbeSBarry Smith extern PetscErrorCode MatHermitianTranspose(Mat,MatReuse,Mat*); 5597087cfbeSBarry Smith extern PetscErrorCode MatPermute(Mat,IS,IS,Mat *); 560ba966639SSatish Balay PetscPolymorphicFunction(MatPermute,(Mat A,IS is1,IS is2),(A,is1,is2,&t),Mat,t) 5617087cfbeSBarry Smith extern PetscErrorCode MatDiagonalScale(Mat,Vec,Vec); 5627087cfbeSBarry Smith extern PetscErrorCode MatDiagonalSet(Mat,Vec,InsertMode); 5637087cfbeSBarry Smith extern PetscErrorCode MatEqual(Mat,Mat,PetscBool *); 564ace3abfcSBarry Smith PetscPolymorphicFunction(MatEqual,(Mat A,Mat B),(A,B,&t),PetscBool ,t) 5657087cfbeSBarry Smith extern PetscErrorCode MatMultEqual(Mat,Mat,PetscInt,PetscBool *); 5667087cfbeSBarry Smith extern PetscErrorCode MatMultAddEqual(Mat,Mat,PetscInt,PetscBool *); 5677087cfbeSBarry Smith extern PetscErrorCode MatMultTransposeEqual(Mat,Mat,PetscInt,PetscBool *); 5687087cfbeSBarry Smith extern PetscErrorCode MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscBool *); 5697b80b807SBarry Smith 5707087cfbeSBarry Smith extern PetscErrorCode MatNorm(Mat,NormType,PetscReal *); 571ba966639SSatish Balay PetscPolymorphicFunction(MatNorm,(Mat A,NormType t),(A,t,&n),PetscReal,n) 57209573ac7SBarry Smith extern PetscErrorCode MatGetColumnNorms(Mat,NormType,PetscReal *); 5737087cfbeSBarry Smith extern PetscErrorCode MatZeroEntries(Mat); 5747087cfbeSBarry Smith extern PetscErrorCode MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 5757087cfbeSBarry Smith extern PetscErrorCode MatZeroRowsIS(Mat,IS,PetscScalar,Vec,Vec); 5767087cfbeSBarry Smith extern PetscErrorCode MatZeroRowsStencil(Mat,PetscInt,const MatStencil [],PetscScalar,Vec,Vec); 5777087cfbeSBarry Smith extern PetscErrorCode MatZeroRowsColumns(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 5787087cfbeSBarry Smith extern PetscErrorCode MatZeroRowsColumnsIS(Mat,IS,PetscScalar,Vec,Vec); 5797b80b807SBarry Smith 5807087cfbeSBarry Smith extern PetscErrorCode MatUseScaledForm(Mat,PetscBool ); 5817087cfbeSBarry Smith extern PetscErrorCode MatScaleSystem(Mat,Vec,Vec); 5827087cfbeSBarry Smith extern PetscErrorCode MatUnScaleSystem(Mat,Vec,Vec); 5835ef9f2a5SBarry Smith 5847087cfbeSBarry Smith extern PetscErrorCode MatGetSize(Mat,PetscInt*,PetscInt*); 5857087cfbeSBarry Smith extern PetscErrorCode MatGetLocalSize(Mat,PetscInt*,PetscInt*); 5867087cfbeSBarry Smith extern PetscErrorCode MatGetOwnershipRange(Mat,PetscInt*,PetscInt*); 5877087cfbeSBarry Smith extern PetscErrorCode MatGetOwnershipRanges(Mat,const PetscInt**); 5887087cfbeSBarry Smith extern PetscErrorCode MatGetOwnershipRangeColumn(Mat,PetscInt*,PetscInt*); 5897087cfbeSBarry Smith extern PetscErrorCode MatGetOwnershipRangesColumn(Mat,const PetscInt**); 5907b80b807SBarry Smith 5917087cfbeSBarry Smith extern PetscErrorCode MatGetSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]); 5927087cfbeSBarry Smith extern PetscErrorCode MatGetSubMatricesParallel(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]); 5937087cfbeSBarry Smith extern PetscErrorCode MatDestroyMatrices(PetscInt,Mat *[]); 5947087cfbeSBarry Smith extern PetscErrorCode MatGetSubMatrix(Mat,IS,IS,MatReuse,Mat *); 5957087cfbeSBarry Smith extern PetscErrorCode MatGetLocalSubMatrix(Mat,IS,IS,Mat*); 5967087cfbeSBarry Smith extern PetscErrorCode MatRestoreLocalSubMatrix(Mat,IS,IS,Mat*); 5977087cfbeSBarry Smith extern PetscErrorCode MatGetSeqNonzeroStructure(Mat,Mat*); 5987087cfbeSBarry Smith extern PetscErrorCode MatDestroySeqNonzeroStructure(Mat*); 5995494a064SHong Zhang 6007087cfbeSBarry Smith extern PetscErrorCode MatMerge(MPI_Comm,Mat,PetscInt,MatReuse,Mat*); 6017087cfbeSBarry Smith extern PetscErrorCode MatMerge_SeqsToMPI(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*); 6027087cfbeSBarry Smith extern PetscErrorCode MatMerge_SeqsToMPISymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*); 6037087cfbeSBarry Smith extern PetscErrorCode MatMerge_SeqsToMPINumeric(Mat,Mat); 6044a2b5492SBarry Smith extern PetscErrorCode MatMPIAIJGetLocalMat(Mat,MatReuse,Mat*); 6054a2b5492SBarry Smith extern PetscErrorCode MatMPIAIJGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*); 6067087cfbeSBarry Smith extern PetscErrorCode MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,PetscInt*,Mat*); 6077087cfbeSBarry Smith extern PetscErrorCode MatGetBrowsOfAoCols(Mat,Mat,MatReuse,PetscInt**,PetscInt**,MatScalar**,Mat*); 60843eb5e2fSMatthew Knepley #if defined (PETSC_USE_CTABLE) 609cd116e26SSatish Balay #include "petscctable.h" 6107087cfbeSBarry Smith extern PetscErrorCode MatGetCommunicationStructs(Mat, Vec *, PetscTable *, VecScatter *); 61143eb5e2fSMatthew Knepley #else 6127087cfbeSBarry Smith extern PetscErrorCode MatGetCommunicationStructs(Mat, Vec *, PetscInt *[], VecScatter *); 61343eb5e2fSMatthew Knepley #endif 6147087cfbeSBarry Smith extern PetscErrorCode MatGetGhosts(Mat, PetscInt *,const PetscInt *[]); 6158efafbd8SBarry Smith 6167087cfbeSBarry Smith extern PetscErrorCode MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt); 6177b80b807SBarry Smith 6187087cfbeSBarry Smith extern PetscErrorCode MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*); 6197087cfbeSBarry Smith extern PetscErrorCode MatMatMultSymbolic(Mat,Mat,PetscReal,Mat*); 6207087cfbeSBarry Smith extern PetscErrorCode MatMatMultNumeric(Mat,Mat,Mat); 62122440eb1SKris Buschelman 6227087cfbeSBarry Smith extern PetscErrorCode MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*); 6237087cfbeSBarry Smith extern PetscErrorCode MatPtAPSymbolic(Mat,Mat,PetscReal,Mat*); 6247087cfbeSBarry Smith extern PetscErrorCode MatPtAPNumeric(Mat,Mat,Mat); 62522440eb1SKris Buschelman 6267087cfbeSBarry Smith extern PetscErrorCode MatMatMultTranspose(Mat,Mat,MatReuse,PetscReal,Mat*); 6277087cfbeSBarry Smith extern PetscErrorCode MatMatMultTransposeSymbolic(Mat,Mat,PetscReal,Mat*); 6287087cfbeSBarry Smith extern PetscErrorCode MatMatMultTransposeNumeric(Mat,Mat,Mat); 629bc011b1eSHong Zhang 6307087cfbeSBarry Smith extern PetscErrorCode MatAXPY(Mat,PetscScalar,Mat,MatStructure); 6317087cfbeSBarry Smith extern PetscErrorCode MatAYPX(Mat,PetscScalar,Mat,MatStructure); 6321c741599SHong Zhang 6337087cfbeSBarry Smith extern PetscErrorCode MatScale(Mat,PetscScalar); 6347087cfbeSBarry Smith extern PetscErrorCode MatShift(Mat,PetscScalar); 6357b80b807SBarry Smith 6367087cfbeSBarry Smith extern PetscErrorCode MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping); 6377087cfbeSBarry Smith extern PetscErrorCode MatSetLocalToGlobalMappingBlock(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping); 6387087cfbeSBarry Smith extern PetscErrorCode MatGetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*); 6397087cfbeSBarry Smith extern PetscErrorCode MatGetLocalToGlobalMappingBlock(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*); 6407087cfbeSBarry Smith extern PetscErrorCode MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 6417087cfbeSBarry Smith extern PetscErrorCode MatZeroRowsLocalIS(Mat,IS,PetscScalar,Vec,Vec); 6427087cfbeSBarry Smith extern PetscErrorCode MatZeroRowsColumnsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 6437087cfbeSBarry Smith extern PetscErrorCode MatZeroRowsColumnsLocalIS(Mat,IS,PetscScalar,Vec,Vec); 6447087cfbeSBarry Smith extern PetscErrorCode MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 6457087cfbeSBarry Smith extern PetscErrorCode MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 646052efed2SBarry Smith 6477087cfbeSBarry Smith extern PetscErrorCode MatStashSetInitialSize(Mat,PetscInt,PetscInt); 6487087cfbeSBarry Smith extern PetscErrorCode MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*); 64990f02eecSBarry Smith 6507087cfbeSBarry Smith extern PetscErrorCode MatInterpolate(Mat,Vec,Vec); 6517087cfbeSBarry Smith extern PetscErrorCode MatInterpolateAdd(Mat,Vec,Vec,Vec); 652ba966639SSatish Balay PetscPolymorphicSubroutine(MatInterpolateAdd,(Mat A,Vec x,Vec y),(A,x,y,y)) 6537087cfbeSBarry Smith extern PetscErrorCode MatRestrict(Mat,Vec,Vec); 6547087cfbeSBarry Smith extern PetscErrorCode MatGetVecs(Mat,Vec*,Vec*); 6557087cfbeSBarry Smith extern PetscErrorCode MatGetRedundantMatrix(Mat,PetscInt,MPI_Comm,PetscInt,MatReuse,Mat*); 6567087cfbeSBarry Smith extern PetscErrorCode MatGetMultiProcBlock(Mat,MPI_Comm,Mat*); 657fc9bc008SSatish Balay extern PetscErrorCode MatFindZeroDiagonals(Mat,IS*); 658bd481603SBarry Smith 659bd481603SBarry Smith /*MC 6601620fd73SBarry Smith MatSetValue - Set a single entry into a matrix. 6611620fd73SBarry Smith 6621620fd73SBarry Smith Not collective 6631620fd73SBarry Smith 6641620fd73SBarry Smith Input Parameters: 6651620fd73SBarry Smith + m - the matrix 6661620fd73SBarry Smith . row - the row location of the entry 6671620fd73SBarry Smith . col - the column location of the entry 6681620fd73SBarry Smith . value - the value to insert 6691620fd73SBarry Smith - mode - either INSERT_VALUES or ADD_VALUES 6701620fd73SBarry Smith 6711620fd73SBarry Smith Notes: 6721620fd73SBarry Smith For efficiency one should use MatSetValues() and set several or many 6731620fd73SBarry Smith values simultaneously if possible. 6741620fd73SBarry Smith 6751620fd73SBarry Smith Level: beginner 6761620fd73SBarry Smith 6771620fd73SBarry Smith .seealso: MatSetValues(), MatSetValueLocal() 6781620fd73SBarry Smith M*/ 6791620fd73SBarry 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);} 6801620fd73SBarry Smith 6811620fd73SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatGetValue(Mat v,PetscInt i,PetscInt j,PetscScalar *va) {return MatGetValues(v,1,&i,1,&j,va);} 6821620fd73SBarry Smith 6831620fd73SBarry 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);} 6841620fd73SBarry Smith 68519180038SMatthew G Knepley extern PetscErrorCode MatSeqAIJSetValuesBatch(Mat, PetscInt, PetscInt, PetscInt *, PetscScalar *); 68619180038SMatthew G Knepley extern PetscErrorCode MatMPIAIJSetValuesBatch(Mat, PetscInt, PetscInt, PetscInt *, PetscScalar *); 6876acf6a18SMatthew G Knepley 6881620fd73SBarry Smith /*MC 6890d2a7ff2SSatish Balay MatPreallocateInitialize - Begins the block of code that will count the number of nonzeros per 690bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 691bd481603SBarry Smith 692bd481603SBarry Smith Synopsis: 693c1ac3661SBarry Smith PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz) 694bd481603SBarry Smith 695bd481603SBarry Smith Collective on MPI_Comm 696bd481603SBarry Smith 697bd481603SBarry Smith Input Parameters: 698bd481603SBarry Smith + comm - the communicator that will share the eventually allocated matrix 699859fcb39SBarry Smith . nrows - the number of LOCAL rows in the matrix 700859fcb39SBarry Smith - ncols - the number of LOCAL columns in the matrix 701bd481603SBarry Smith 702bd481603SBarry Smith Output Parameters: 703bd481603SBarry Smith + dnz - the array that will be passed to the matrix preallocation routines 704bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 705bd481603SBarry Smith 706bd481603SBarry Smith 707bd481603SBarry Smith Level: intermediate 708bd481603SBarry Smith 709bd481603SBarry Smith Notes: 7100598bfebSBarry Smith See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details. 711bd481603SBarry Smith 7121d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 713bd481603SBarry Smith 714bd481603SBarry Smith Use MatPreallocateInitializeSymmetric() for symmetric matrices (MPISBAIJ matrices) 715bd481603SBarry Smith 7161620fd73SBarry Smith This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize(). 7171620fd73SBarry Smith 718bd481603SBarry Smith Concepts: preallocation^Matrix 719bd481603SBarry Smith 720bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(), 721bd481603SBarry Smith MatPreallocateInitializeSymmetric(), MatPreallocateSymmetricSetLocal() 722bd481603SBarry Smith M*/ 723c4f061fbSSatish Balay #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \ 7247c922b88SBarry Smith { \ 725a89bc211SBarry Smith PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__start,__end; \ 726a89bc211SBarry Smith _4_ierr = PetscMalloc2(__nrows,PetscInt,&dnz,__nrows,PetscInt,&onz);CHKERRQ(_4_ierr); \ 727a89bc211SBarry Smith _4_ierr = PetscMemzero(dnz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\ 728a89bc211SBarry Smith _4_ierr = PetscMemzero(onz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\ 729c1ac3661SBarry Smith _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ctmp;\ 730a89bc211SBarry Smith _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows; 7317c922b88SBarry Smith 732bd481603SBarry Smith /*MC 7330d2a7ff2SSatish Balay MatPreallocateSymmetricInitialize - Begins the block of code that will count the number of nonzeros per 734bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 735bd481603SBarry Smith 736bd481603SBarry Smith Synopsis: 737c1ac3661SBarry Smith PetscErrorCode MatPreallocateSymmetricInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz) 738bd481603SBarry Smith 739bd481603SBarry Smith Collective on MPI_Comm 740bd481603SBarry Smith 741bd481603SBarry Smith Input Parameters: 742bd481603SBarry Smith + comm - the communicator that will share the eventually allocated matrix 743859fcb39SBarry Smith . nrows - the number of LOCAL rows in the matrix 744859fcb39SBarry Smith - ncols - the number of LOCAL columns in the matrix 745bd481603SBarry Smith 746bd481603SBarry Smith Output Parameters: 747bd481603SBarry Smith + dnz - the array that will be passed to the matrix preallocation routines 748bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 749bd481603SBarry Smith 750bd481603SBarry Smith 751bd481603SBarry Smith Level: intermediate 752bd481603SBarry Smith 753bd481603SBarry Smith Notes: 7540598bfebSBarry Smith See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details. 755bd481603SBarry Smith 7561d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 757bd481603SBarry Smith 7581620fd73SBarry Smith This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize(). 7591620fd73SBarry Smith 760bd481603SBarry Smith Concepts: preallocation^Matrix 761bd481603SBarry Smith 762bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(), 763bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal() 764bd481603SBarry Smith M*/ 765222b16d4SBarry Smith #define MatPreallocateSymmetricInitialize(comm,nrows,ncols,dnz,onz) 0; \ 766222b16d4SBarry Smith { \ 767a89bc211SBarry Smith PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__end; \ 768a89bc211SBarry Smith _4_ierr = PetscMalloc2(__nrows,PetscInt,&dnz,__nrows,PetscInt,&onz);CHKERRQ(_4_ierr); \ 769a89bc211SBarry Smith _4_ierr = PetscMemzero(dnz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\ 770a89bc211SBarry Smith _4_ierr = PetscMemzero(onz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\ 771c1ac3661SBarry Smith _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr);\ 772a89bc211SBarry Smith _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows; 773222b16d4SBarry Smith 774bd481603SBarry Smith /*MC 775bd481603SBarry Smith MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be 776bd481603SBarry Smith inserted using a local number of the rows and columns 777bd481603SBarry Smith 778bd481603SBarry Smith Synopsis: 779c1ac3661SBarry Smith PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 780bd481603SBarry Smith 781bd481603SBarry Smith Not Collective 782bd481603SBarry Smith 783bd481603SBarry Smith Input Parameters: 784784ac674SJed Brown + map - the row mapping from local numbering to global numbering 785bd481603SBarry Smith . nrows - the number of rows indicated 7861d73ed98SBarry Smith . rows - the indices of the rows 787784ac674SJed Brown . cmap - the column mapping from local to global numbering 788bd481603SBarry Smith . ncols - the number of columns in the matrix 789bd481603SBarry Smith . cols - the columns indicated 790bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 791bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 792bd481603SBarry Smith 793bd481603SBarry Smith 794bd481603SBarry Smith Level: intermediate 795bd481603SBarry Smith 796bd481603SBarry Smith Notes: 7970598bfebSBarry Smith See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details. 798bd481603SBarry Smith 7991d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 800bd481603SBarry Smith 801bd481603SBarry Smith Concepts: preallocation^Matrix 802bd481603SBarry Smith 803bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(), 804bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal() 805bd481603SBarry Smith M*/ 806784ac674SJed Brown #define MatPreallocateSetLocal(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \ 807c4f061fbSSatish Balay {\ 808c1ac3661SBarry Smith PetscInt __l;\ 809784ac674SJed Brown _4_ierr = ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\ 810784ac674SJed Brown _4_ierr = ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\ 811c4f061fbSSatish Balay for (__l=0;__l<nrows;__l++) {\ 812ef66eb69SBarry Smith _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 813c4f061fbSSatish Balay }\ 814c4f061fbSSatish Balay } 815c4f061fbSSatish Balay 816bd481603SBarry Smith /*MC 817bd481603SBarry Smith MatPreallocateSymmetricSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be 818bd481603SBarry Smith inserted using a local number of the rows and columns 819bd481603SBarry Smith 820bd481603SBarry Smith Synopsis: 821c1ac3661SBarry Smith PetscErrorCode MatPreallocateSymmetricSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 822bd481603SBarry Smith 823bd481603SBarry Smith Not Collective 824bd481603SBarry Smith 825bd481603SBarry Smith Input Parameters: 826bd481603SBarry Smith + map - the mapping between local numbering and global numbering 827bd481603SBarry Smith . nrows - the number of rows indicated 8281d73ed98SBarry Smith . rows - the indices of the rows 829bd481603SBarry Smith . ncols - the number of columns in the matrix 830bd481603SBarry Smith . cols - the columns indicated 831bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 832bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 833bd481603SBarry Smith 834bd481603SBarry Smith 835bd481603SBarry Smith Level: intermediate 836bd481603SBarry Smith 837bd481603SBarry Smith Notes: 8380598bfebSBarry Smith See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details. 839bd481603SBarry Smith 840bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 841bd481603SBarry Smith 842bd481603SBarry Smith Concepts: preallocation^Matrix 843bd481603SBarry Smith 844bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(), 845bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal() 846bd481603SBarry Smith M*/ 847d3d32019SBarry Smith #define MatPreallocateSymmetricSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\ 848d3d32019SBarry Smith {\ 849c1ac3661SBarry Smith PetscInt __l;\ 850d3d32019SBarry Smith _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\ 851d3d32019SBarry Smith _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\ 852d3d32019SBarry Smith for (__l=0;__l<nrows;__l++) {\ 853d3d32019SBarry Smith _4_ierr = MatPreallocateSymmetricSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 854d3d32019SBarry Smith }\ 855d3d32019SBarry Smith } 856d3d32019SBarry Smith 857bd481603SBarry Smith /*MC 858bd481603SBarry Smith MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be 859bd481603SBarry Smith inserted using a local number of the rows and columns 860bd481603SBarry Smith 861bd481603SBarry Smith Synopsis: 862c1ac3661SBarry Smith PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 863bd481603SBarry Smith 864bd481603SBarry Smith Not Collective 865bd481603SBarry Smith 866bd481603SBarry Smith Input Parameters: 86764075487SBarry Smith + row - the row 868bd481603SBarry Smith . ncols - the number of columns in the matrix 869a50f8bf6SHong Zhang - cols - the columns indicated 870a50f8bf6SHong Zhang 871a50f8bf6SHong Zhang Output Parameters: 872a50f8bf6SHong Zhang + dnz - the array that will be passed to the matrix preallocation routines 873bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 874bd481603SBarry Smith 875bd481603SBarry Smith 876bd481603SBarry Smith Level: intermediate 877bd481603SBarry Smith 878bd481603SBarry Smith Notes: 8790598bfebSBarry Smith See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details. 880bd481603SBarry Smith 881bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 882bd481603SBarry Smith 8831620fd73SBarry Smith This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize(). 8841620fd73SBarry Smith 885bd481603SBarry Smith Concepts: preallocation^Matrix 886bd481603SBarry Smith 887bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(), 888bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal() 889bd481603SBarry Smith M*/ 890c4f061fbSSatish Balay #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\ 891c1ac3661SBarry Smith { PetscInt __i; \ 892e32f2f54SBarry 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);\ 893e32f2f54SBarry 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);\ 8947c922b88SBarry Smith for (__i=0; __i<nc; __i++) {\ 89564075487SBarry Smith if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \ 8967cc688d7SBarry Smith else dnz[row - __rstart]++;\ 8977c922b88SBarry Smith }\ 8987c922b88SBarry Smith } 8997c922b88SBarry Smith 900bd481603SBarry Smith /*MC 901bd481603SBarry Smith MatPreallocateSymmetricSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be 902bd481603SBarry Smith inserted using a local number of the rows and columns 903bd481603SBarry Smith 904bd481603SBarry Smith Synopsis: 905c1ac3661SBarry Smith PetscErrorCode MatPreallocateSymmetricSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 906bd481603SBarry Smith 907bd481603SBarry Smith Not Collective 908bd481603SBarry Smith 909bd481603SBarry Smith Input Parameters: 910bd481603SBarry Smith + nrows - the number of rows indicated 9111d73ed98SBarry Smith . rows - the indices of the rows 912bd481603SBarry Smith . ncols - the number of columns in the matrix 913bd481603SBarry Smith . cols - the columns indicated 914bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 915bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 916bd481603SBarry Smith 917bd481603SBarry Smith 918bd481603SBarry Smith Level: intermediate 919bd481603SBarry Smith 920bd481603SBarry Smith Notes: 9210598bfebSBarry Smith See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details. 922bd481603SBarry Smith 923bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 924bd481603SBarry Smith 9251620fd73SBarry Smith This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize(). 9261620fd73SBarry Smith 927bd481603SBarry Smith Concepts: preallocation^Matrix 928bd481603SBarry Smith 929bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(), 930bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal() 931bd481603SBarry Smith M*/ 932d3d32019SBarry Smith #define MatPreallocateSymmetricSet(row,nc,cols,dnz,onz) 0;\ 933c1ac3661SBarry Smith { PetscInt __i; \ 934d3d32019SBarry Smith for (__i=0; __i<nc; __i++) {\ 935d3d32019SBarry Smith if (cols[__i] >= __end) onz[row - __rstart]++; \ 936d3d32019SBarry Smith else if (cols[__i] >= row) dnz[row - __rstart]++;\ 937d3d32019SBarry Smith }\ 938d3d32019SBarry Smith } 939d3d32019SBarry Smith 940bd481603SBarry Smith /*MC 94116371a99SBarry Smith MatPreallocateLocation - An alternative to MatPreallocationSet() that puts the nonzero locations into the matrix if it exists 94216371a99SBarry Smith 94316371a99SBarry Smith Synopsis: 94416371a99SBarry Smith PetscErrorCode MatPreallocateLocations(Mat A,PetscInt row,PetscInt ncols,PetscInt *cols,PetscInt *dnz,PetscInt *onz) 94516371a99SBarry Smith 94616371a99SBarry Smith Not Collective 94716371a99SBarry Smith 94816371a99SBarry Smith Input Parameters: 94916371a99SBarry Smith . A - matrix 95016371a99SBarry Smith . row - row where values exist (must be local to this process) 95116371a99SBarry Smith . ncols - number of columns 95216371a99SBarry Smith . cols - columns with nonzeros 95316371a99SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 95416371a99SBarry Smith - ozn - the other array passed to the matrix preallocation routines 95516371a99SBarry Smith 95616371a99SBarry Smith 95716371a99SBarry Smith Level: intermediate 95816371a99SBarry Smith 95916371a99SBarry Smith Notes: 9600598bfebSBarry Smith See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details. 96116371a99SBarry Smith 96216371a99SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 96316371a99SBarry Smith 96416371a99SBarry Smith This is a MACRO not a function because it uses a bunch of variables private to the MatPreallocation.... routines. 96516371a99SBarry Smith 96616371a99SBarry Smith Concepts: preallocation^Matrix 96716371a99SBarry Smith 96816371a99SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(), 96916371a99SBarry Smith MatPreallocateSymmetricInitialize(), MatPreallocateSymmetricSetLocal() 97016371a99SBarry Smith M*/ 97116371a99SBarry Smith #define MatPreallocateLocation(A,row,ncols,cols,dnz,onz) 0;if (A) {ierr = MatSetValues(A,1,&row,ncols,cols,PETSC_NULL,INSERT_VALUES);CHKERRQ(ierr);} else {ierr = MatPreallocateSet(row,ncols,cols,dnz,onz);CHKERRQ(ierr);} 97216371a99SBarry Smith 97316371a99SBarry Smith 97416371a99SBarry Smith /*MC 9750d2a7ff2SSatish Balay MatPreallocateFinalize - Ends the block of code that will count the number of nonzeros per 976bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 977bd481603SBarry Smith 978bd481603SBarry Smith Synopsis: 979c1ac3661SBarry Smith PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz) 980bd481603SBarry Smith 981bd481603SBarry Smith Collective on MPI_Comm 982bd481603SBarry Smith 983bd481603SBarry Smith Input Parameters: 98416371a99SBarry Smith + dnz - the array that was be passed to the matrix preallocation routines 985bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 986bd481603SBarry Smith 987bd481603SBarry Smith 988bd481603SBarry Smith Level: intermediate 989bd481603SBarry Smith 990bd481603SBarry Smith Notes: 9910598bfebSBarry Smith See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details. 992bd481603SBarry Smith 993bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 994bd481603SBarry Smith 9951620fd73SBarry Smith This is a MACRO not a function because it closes the { started in MatPreallocateInitialize(). 9961620fd73SBarry Smith 997bd481603SBarry Smith Concepts: preallocation^Matrix 998bd481603SBarry Smith 999bd481603SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(), 1000bd481603SBarry Smith MatPreallocateSymmetricInitialize(), MatPreallocateSymmetricSetLocal() 1001bd481603SBarry Smith M*/ 1002a89bc211SBarry Smith #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree2(dnz,onz);CHKERRQ(_4_ierr);} 10037c922b88SBarry Smith 1004bd481603SBarry Smith 1005bd481603SBarry Smith 10067b80b807SBarry Smith /* Routines unique to particular data structures */ 10078fe8eb27SJed Brown extern PetscErrorCode MatShellGetContext(Mat,void *); 1008ba966639SSatish Balay PetscPolymorphicFunction(MatShellGetContext,(Mat A),(A,&t),void*,t) 1009698d4c6aSKris Buschelman 10107087cfbeSBarry Smith extern PetscErrorCode MatInodeAdjustForInodes(Mat,IS*,IS*); 10117087cfbeSBarry Smith extern PetscErrorCode MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *); 1012698d4c6aSKris Buschelman 10137087cfbeSBarry Smith extern PetscErrorCode MatSeqAIJSetColumnIndices(Mat,PetscInt[]); 10147087cfbeSBarry Smith extern PetscErrorCode MatSeqBAIJSetColumnIndices(Mat,PetscInt[]); 10157087cfbeSBarry Smith extern PetscErrorCode MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 10167087cfbeSBarry Smith extern PetscErrorCode MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 10177087cfbeSBarry Smith extern PetscErrorCode MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 10187b80b807SBarry Smith 1019a96a251dSBarry Smith #define MAT_SKIP_ALLOCATION -4 1020a96a251dSBarry Smith 10217087cfbeSBarry Smith extern PetscErrorCode MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]); 1022ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz)) 10237087cfbeSBarry Smith extern PetscErrorCode MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]); 1024ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqSBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz)) 10257087cfbeSBarry Smith extern PetscErrorCode MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]); 1026ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqAIJSetPreallocation,(Mat A,const PetscInt nnz[]),(A,0,nnz)) 1027273d9f13SBarry Smith 10287087cfbeSBarry Smith extern PetscErrorCode MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 1029ba966639SSatish Balay PetscPolymorphicSubroutine(MatMPIBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[],const PetscInt onz[]),(A,bs,0,nnz,0,onz)) 10307087cfbeSBarry Smith extern PetscErrorCode MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 10317087cfbeSBarry Smith extern PetscErrorCode MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 10327087cfbeSBarry Smith extern PetscErrorCode MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []); 10337087cfbeSBarry Smith extern PetscErrorCode MatSeqBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 10347087cfbeSBarry Smith extern PetscErrorCode MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]); 10357087cfbeSBarry Smith extern PetscErrorCode MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 10367087cfbeSBarry Smith extern PetscErrorCode MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]); 10377087cfbeSBarry Smith extern PetscErrorCode MatMPIDenseSetPreallocation(Mat,PetscScalar[]); 10387087cfbeSBarry Smith extern PetscErrorCode MatSeqDenseSetPreallocation(Mat,PetscScalar[]); 10397087cfbeSBarry Smith extern PetscErrorCode MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,PetscInt*[]); 10407087cfbeSBarry Smith extern PetscErrorCode MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,PetscInt*[]); 10417087cfbeSBarry Smith extern PetscErrorCode MatAdicSetLocalFunction(Mat,void (*)(void)); 1042273d9f13SBarry Smith 10437087cfbeSBarry Smith extern PetscErrorCode MatSeqDenseSetLDA(Mat,PetscInt); 10447087cfbeSBarry Smith extern PetscErrorCode MatDenseGetLocalMatrix(Mat,Mat*); 10451b807ce4Svictorle 10467087cfbeSBarry Smith extern PetscErrorCode MatStoreValues(Mat); 10477087cfbeSBarry Smith extern PetscErrorCode MatRetrieveValues(Mat); 10482e8a6d31SBarry Smith 10497087cfbeSBarry Smith extern PetscErrorCode MatDAADSetCtx(Mat,void*); 10503a7fca6bSBarry Smith 1051b3a44c85SBarry Smith extern PetscErrorCode MatFindNonzeroRows(Mat,IS*); 10527b80b807SBarry Smith /* 10537b80b807SBarry Smith These routines are not usually accessed directly, rather solving is 105494b7f48cSBarry Smith done through the KSP and PC interfaces. 10557b80b807SBarry Smith */ 10567b80b807SBarry Smith 1057d9274352SBarry Smith /*E 1058d9274352SBarry Smith MatOrderingType - String with the name of a PETSc matrix ordering or the creation function 1059d9274352SBarry Smith with an optional dynamic library name, for example 1060d9274352SBarry Smith http://www.mcs.anl.gov/petsc/lib.a:orderingcreate() 1061d9274352SBarry Smith 1062d9274352SBarry Smith Level: beginner 1063d9274352SBarry Smith 1064e5a9bf91SBarry Smith Cannot use const because the PC objects manipulate the string 1065e5a9bf91SBarry Smith 1066d9274352SBarry Smith .seealso: MatGetOrdering() 1067d9274352SBarry Smith E*/ 10683eaad2caSSatish Balay #define MatOrderingType char* 10692692d6eeSBarry Smith #define MATORDERINGNATURAL "natural" 10702692d6eeSBarry Smith #define MATORDERINGND "nd" 10712692d6eeSBarry Smith #define MATORDERING1WD "1wd" 10722692d6eeSBarry Smith #define MATORDERINGRCM "rcm" 10732692d6eeSBarry Smith #define MATORDERINGQMD "qmd" 10742692d6eeSBarry Smith #define MATORDERINGROWLENGTH "rowlength" 10752692d6eeSBarry Smith #define MATORDERINGDSC_ND "dsc_nd" /* these three are only for DSCPACK, see its documentation for details */ 10762692d6eeSBarry Smith #define MATORDERINGDSC_MMD "dsc_mmd" 10772692d6eeSBarry Smith #define MATORDERINGDSC_MDF "dsc_mdf" 1078312e372aSBarry Smith #define MATORDERINGAMD "amd" /* only works if UMFPACK is installed with PETSc */ 1079b12f92e5SBarry Smith 10807087cfbeSBarry Smith extern PetscErrorCode MatGetOrdering(Mat,const MatOrderingType,IS*,IS*); 10817087cfbeSBarry Smith extern PetscErrorCode MatGetOrderingList(PetscFList *list); 10827087cfbeSBarry Smith extern PetscErrorCode MatOrderingRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,const MatOrderingType,IS*,IS*)); 108330de9b25SBarry Smith 108430de9b25SBarry Smith /*MC 10851890ba74SBarry Smith MatOrderingRegisterDynamic - Adds a new sparse matrix ordering to the matrix package. 108630de9b25SBarry Smith 108730de9b25SBarry Smith Synopsis: 10881890ba74SBarry Smith PetscErrorCode MatOrderingRegisterDynamic(const char *name_ordering,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatOrdering)) 108930de9b25SBarry Smith 109030de9b25SBarry Smith Not Collective 109130de9b25SBarry Smith 109230de9b25SBarry Smith Input Parameters: 10932692d6eeSBarry Smith + sname - name of ordering (for example MATORDERINGND) 109430de9b25SBarry Smith . path - location of library where creation routine is 109530de9b25SBarry Smith . name - name of function that creates the ordering type,a string 109630de9b25SBarry Smith - function - function pointer that creates the ordering 109730de9b25SBarry Smith 109830de9b25SBarry Smith Level: developer 109930de9b25SBarry Smith 110030de9b25SBarry Smith If dynamic libraries are used, then the fourth input argument (function) 110130de9b25SBarry Smith is ignored. 110230de9b25SBarry Smith 110330de9b25SBarry Smith Sample usage: 110430de9b25SBarry Smith .vb 110530de9b25SBarry Smith MatOrderingRegisterDynamic("my_order",/home/username/my_lib/lib/libO/solaris/mylib.a, 110630de9b25SBarry Smith "MyOrder",MyOrder); 110730de9b25SBarry Smith .ve 110830de9b25SBarry Smith 110930de9b25SBarry Smith Then, your partitioner can be chosen with the procedural interface via 111030de9b25SBarry Smith $ MatOrderingSetType(part,"my_order) 111130de9b25SBarry Smith or at runtime via the option 11122401956bSBarry Smith $ -pc_factor_mat_ordering_type my_order 111330de9b25SBarry Smith 1114ab901514SBarry Smith ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values. 111530de9b25SBarry Smith 111630de9b25SBarry Smith .keywords: matrix, ordering, register 111730de9b25SBarry Smith 111830de9b25SBarry Smith .seealso: MatOrderingRegisterDestroy(), MatOrderingRegisterAll() 111930de9b25SBarry Smith M*/ 1120aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 1121f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,0) 1122b12f92e5SBarry Smith #else 1123f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,d) 1124b12f92e5SBarry Smith #endif 112530de9b25SBarry Smith 11267087cfbeSBarry Smith extern PetscErrorCode MatOrderingRegisterDestroy(void); 11277087cfbeSBarry Smith extern PetscErrorCode MatOrderingRegisterAll(const char[]); 1128ace3abfcSBarry Smith extern PetscBool MatOrderingRegisterAllCalled; 1129b0a32e0cSBarry Smith extern PetscFList MatOrderingList; 1130d4fbbf0eSBarry Smith 11317087cfbeSBarry Smith extern PetscErrorCode MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS); 1132a2ce50c7SBarry Smith 1133d91e6319SBarry Smith /*S 1134d90ac83dSHong Zhang MatFactorShiftType - Numeric Shift. 1135d90ac83dSHong Zhang 1136d90ac83dSHong Zhang Level: beginner 1137d90ac83dSHong Zhang 1138d90ac83dSHong Zhang S*/ 1139d90ac83dSHong Zhang typedef enum {MAT_SHIFT_NONE,MAT_SHIFT_NONZERO,MAT_SHIFT_POSITIVE_DEFINITE,MAT_SHIFT_INBLOCKS} MatFactorShiftType; 1140d90ac83dSHong Zhang extern const char *MatFactorShiftTypes[]; 1141d90ac83dSHong Zhang 1142d90ac83dSHong Zhang /*S 11432401956bSBarry Smith MatFactorInfo - Data passed into the matrix factorization routines 1144b00f7748SHong Zhang 114561cad860SBarry Smith In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE, that is use 114661cad860SBarry Smith $ MatFactorInfo info(MAT_FACTORINFO_SIZE) 1147b00f7748SHong Zhang 114815e8a5b3SHong Zhang Notes: These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC. 1149b00f7748SHong Zhang 115061cad860SBarry Smith You can use MatFactorInfoInitialize() to set default values. 115161cad860SBarry Smith 1152b00f7748SHong Zhang Level: developer 1153b00f7748SHong Zhang 1154d7d82daaSBarry Smith .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(), 1155d7d82daaSBarry Smith MatFactorInfoInitialize() 1156b00f7748SHong Zhang 1157b00f7748SHong Zhang S*/ 1158b00f7748SHong Zhang typedef struct { 115915e8a5b3SHong Zhang PetscReal diagonal_fill; /* force diagonal to fill in if initially not filled */ 116085317021SBarry Smith PetscReal usedt; 116115e8a5b3SHong Zhang PetscReal dt; /* drop tolerance */ 1162b00f7748SHong Zhang PetscReal dtcol; /* tolerance for pivoting */ 116315e8a5b3SHong Zhang PetscReal dtcount; /* maximum nonzeros to be allowed per row */ 116467e28bfeSBarry Smith PetscReal fill; /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */ 1165348344bbSBarry Smith PetscReal levels; /* ICC/ILU(levels) */ 1166bcd9e38bSBarry Smith PetscReal pivotinblocks; /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0 1167bcd9e38bSBarry Smith factorization may be faster if do not pivot */ 116815e8a5b3SHong Zhang PetscReal zeropivot; /* pivot is called zero if less than this */ 1169f4db908eSBarry Smith PetscReal shifttype; /* type of shift added to matrix factor to prevent zero pivots */ 1170f4db908eSBarry Smith PetscReal shiftamount; /* how large the shift is */ 117115e8a5b3SHong Zhang } MatFactorInfo; 1172ffa6d0a5SLois Curfman McInnes 11737087cfbeSBarry Smith extern PetscErrorCode MatFactorInfoInitialize(MatFactorInfo*); 11747087cfbeSBarry Smith extern PetscErrorCode MatCholeskyFactor(Mat,IS,const MatFactorInfo*); 11757087cfbeSBarry Smith extern PetscErrorCode MatCholeskyFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*); 11767087cfbeSBarry Smith extern PetscErrorCode MatCholeskyFactorNumeric(Mat,Mat,const MatFactorInfo*); 11777087cfbeSBarry Smith extern PetscErrorCode MatLUFactor(Mat,IS,IS,const MatFactorInfo*); 11787087cfbeSBarry Smith extern PetscErrorCode MatILUFactor(Mat,IS,IS,const MatFactorInfo*); 11797087cfbeSBarry Smith extern PetscErrorCode MatLUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*); 11807087cfbeSBarry Smith extern PetscErrorCode MatILUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*); 11817087cfbeSBarry Smith extern PetscErrorCode MatICCFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*); 11827087cfbeSBarry Smith extern PetscErrorCode MatICCFactor(Mat,IS,const MatFactorInfo*); 11837087cfbeSBarry Smith extern PetscErrorCode MatLUFactorNumeric(Mat,Mat,const MatFactorInfo*); 11847087cfbeSBarry Smith extern PetscErrorCode MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*); 11857087cfbeSBarry Smith extern PetscErrorCode MatSolve(Mat,Vec,Vec); 11867087cfbeSBarry Smith extern PetscErrorCode MatForwardSolve(Mat,Vec,Vec); 11877087cfbeSBarry Smith extern PetscErrorCode MatBackwardSolve(Mat,Vec,Vec); 11887087cfbeSBarry Smith extern PetscErrorCode MatSolveAdd(Mat,Vec,Vec,Vec); 11897087cfbeSBarry Smith extern PetscErrorCode MatSolveTranspose(Mat,Vec,Vec); 11907087cfbeSBarry Smith extern PetscErrorCode MatSolveTransposeAdd(Mat,Vec,Vec,Vec); 11917087cfbeSBarry Smith extern PetscErrorCode MatSolves(Mat,Vecs,Vecs); 11928ed539a5SBarry Smith 11937087cfbeSBarry Smith extern PetscErrorCode MatSetUnfactored(Mat); 1194bb5a7306SBarry Smith 1195d91e6319SBarry Smith /*E 1196d91e6319SBarry Smith MatSORType - What type of (S)SOR to perform 1197bb1eb677SSatish Balay 1198d91e6319SBarry Smith Level: beginner 1199d91e6319SBarry Smith 1200d9274352SBarry Smith May be bitwise ORd together 1201d9274352SBarry Smith 1202d91e6319SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 1203d91e6319SBarry Smith 12044e7234bfSBarry Smith MatSORType may be bitwise ORd together, so do not change the numbers 12054e7234bfSBarry Smith 120641f059aeSBarry Smith .seealso: MatSOR() 1207d91e6319SBarry Smith E*/ 1208ee50ffe9SBarry Smith typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3, 1209ee50ffe9SBarry Smith SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8, 1210ee50ffe9SBarry Smith SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16, 121184cb2905SBarry Smith SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType; 12127087cfbeSBarry Smith extern PetscErrorCode MatSOR(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec); 12138ed539a5SBarry Smith 1214d4fbbf0eSBarry Smith /* 1215639f9d9dSBarry Smith These routines are for efficiently computing Jacobians via finite differences. 1216639f9d9dSBarry Smith */ 1217b12f92e5SBarry Smith 1218d9274352SBarry Smith /*E 1219d9274352SBarry Smith MatColoringType - String with the name of a PETSc matrix coloring or the creation function 1220d9274352SBarry Smith with an optional dynamic library name, for example 1221d9274352SBarry Smith http://www.mcs.anl.gov/petsc/lib.a:coloringcreate() 1222d9274352SBarry Smith 1223d9274352SBarry Smith Level: beginner 1224d9274352SBarry Smith 1225d9274352SBarry Smith .seealso: MatGetColoring() 1226d9274352SBarry Smith E*/ 1227a313700dSBarry Smith #define MatColoringType char* 12282692d6eeSBarry Smith #define MATCOLORINGNATURAL "natural" 12292692d6eeSBarry Smith #define MATCOLORINGSL "sl" 12302692d6eeSBarry Smith #define MATCOLORINGLF "lf" 12312692d6eeSBarry Smith #define MATCOLORINGID "id" 1232b12f92e5SBarry Smith 12337087cfbeSBarry Smith extern PetscErrorCode MatGetColoring(Mat,const MatColoringType,ISColoring*); 12347087cfbeSBarry Smith extern PetscErrorCode MatColoringRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,MatColoringType,ISColoring *)); 123530de9b25SBarry Smith 123630de9b25SBarry Smith /*MC 123730de9b25SBarry Smith MatColoringRegisterDynamic - Adds a new sparse matrix coloring to the 123830de9b25SBarry Smith matrix package. 123930de9b25SBarry Smith 124030de9b25SBarry Smith Synopsis: 12411890ba74SBarry Smith PetscErrorCode MatColoringRegisterDynamic(const char *name_coloring,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatColoring)) 124230de9b25SBarry Smith 124330de9b25SBarry Smith Not Collective 124430de9b25SBarry Smith 124530de9b25SBarry Smith Input Parameters: 12462692d6eeSBarry Smith + sname - name of Coloring (for example MATCOLORINGSL) 124730de9b25SBarry Smith . path - location of library where creation routine is 124830de9b25SBarry Smith . name - name of function that creates the Coloring type, a string 124930de9b25SBarry Smith - function - function pointer that creates the coloring 125030de9b25SBarry Smith 125130de9b25SBarry Smith Level: developer 125230de9b25SBarry Smith 125330de9b25SBarry Smith If dynamic libraries are used, then the fourth input argument (function) 125430de9b25SBarry Smith is ignored. 125530de9b25SBarry Smith 125630de9b25SBarry Smith Sample usage: 125730de9b25SBarry Smith .vb 125830de9b25SBarry Smith MatColoringRegisterDynamic("my_color",/home/username/my_lib/lib/libO/solaris/mylib.a, 125930de9b25SBarry Smith "MyColor",MyColor); 126030de9b25SBarry Smith .ve 126130de9b25SBarry Smith 126230de9b25SBarry Smith Then, your partitioner can be chosen with the procedural interface via 126330de9b25SBarry Smith $ MatColoringSetType(part,"my_color") 126430de9b25SBarry Smith or at runtime via the option 126530de9b25SBarry Smith $ -mat_coloring_type my_color 126630de9b25SBarry Smith 1267ab901514SBarry Smith $PETSC_ARCH occuring in pathname will be replaced with appropriate values. 126830de9b25SBarry Smith 126930de9b25SBarry Smith .keywords: matrix, Coloring, register 127030de9b25SBarry Smith 127130de9b25SBarry Smith .seealso: MatColoringRegisterDestroy(), MatColoringRegisterAll() 127230de9b25SBarry Smith M*/ 1273aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 1274f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,0) 1275b12f92e5SBarry Smith #else 1276f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,d) 1277b12f92e5SBarry Smith #endif 127830de9b25SBarry Smith 1279ace3abfcSBarry Smith extern PetscBool MatColoringRegisterAllCalled; 1280f1144a30SSatish Balay 12817087cfbeSBarry Smith extern PetscErrorCode MatColoringRegisterAll(const char[]); 12827087cfbeSBarry Smith extern PetscErrorCode MatColoringRegisterDestroy(void); 12837087cfbeSBarry Smith extern PetscErrorCode MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*); 1284639f9d9dSBarry Smith 1285d9274352SBarry Smith /*S 1286d9274352SBarry Smith MatFDColoring - Object for computing a sparse Jacobian via finite differences 1287d9274352SBarry Smith and coloring 1288639f9d9dSBarry Smith 1289d9274352SBarry Smith Level: beginner 1290d9274352SBarry Smith 1291d9274352SBarry Smith Concepts: coloring, sparse Jacobian, finite differences 1292d9274352SBarry Smith 1293d9274352SBarry Smith .seealso: MatFDColoringCreate() 1294d9274352SBarry Smith S*/ 1295e2a1c21fSSatish Balay typedef struct _p_MatFDColoring* MatFDColoring; 1296639f9d9dSBarry Smith 12977087cfbeSBarry Smith extern PetscErrorCode MatFDColoringCreate(Mat,ISColoring,MatFDColoring *); 1298fcfd50ebSBarry Smith extern PetscErrorCode MatFDColoringDestroy(MatFDColoring*); 12997087cfbeSBarry Smith extern PetscErrorCode MatFDColoringView(MatFDColoring,PetscViewer); 13007087cfbeSBarry Smith extern PetscErrorCode MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*); 13017087cfbeSBarry Smith extern PetscErrorCode MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**); 13027087cfbeSBarry Smith extern PetscErrorCode MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal); 13037087cfbeSBarry Smith extern PetscErrorCode MatFDColoringSetFromOptions(MatFDColoring); 13047087cfbeSBarry Smith extern PetscErrorCode MatFDColoringApply(Mat,MatFDColoring,Vec,MatStructure*,void *); 13057087cfbeSBarry Smith extern PetscErrorCode MatFDColoringApplyTS(Mat,MatFDColoring,PetscReal,Vec,MatStructure*,void *); 13067087cfbeSBarry Smith extern PetscErrorCode MatFDColoringSetF(MatFDColoring,Vec); 13077087cfbeSBarry Smith extern PetscErrorCode MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,PetscInt*[]); 1308639f9d9dSBarry Smith /* 13090752156aSBarry Smith These routines are for partitioning matrices: currently used only 13103eda8832SBarry Smith for adjacency matrix, MatCreateMPIAdj(). 13110752156aSBarry Smith */ 1312ca161407SBarry Smith 1313d9274352SBarry Smith /*S 1314d9274352SBarry Smith MatPartitioning - Object for managing the partitioning of a matrix or graph 1315d9274352SBarry Smith 1316d9274352SBarry Smith Level: beginner 1317d9274352SBarry Smith 1318d9274352SBarry Smith Concepts: partitioning 1319d9274352SBarry Smith 1320743a1026Svictorle .seealso: MatPartitioningCreate(), MatPartitioningType 1321d9274352SBarry Smith S*/ 132291e9ee9fSBarry Smith typedef struct _p_MatPartitioning* MatPartitioning; 1323d9274352SBarry Smith 1324d9274352SBarry Smith /*E 13255bc6e7ccSvictorle MatPartitioningType - String with the name of a PETSc matrix partitioning or the creation function 1326d9274352SBarry Smith with an optional dynamic library name, for example 1327d9274352SBarry Smith http://www.mcs.anl.gov/petsc/lib.a:partitioningcreate() 1328d9274352SBarry Smith 1329d9274352SBarry Smith Level: beginner 1330d9274352SBarry Smith 1331b547a13bSvictorle .seealso: MatPartitioningCreate(), MatPartitioning 1332d9274352SBarry Smith E*/ 1333a313700dSBarry Smith #define MatPartitioningType char* 13342692d6eeSBarry Smith #define MATPARTITIONINGCURRENT "current" 13352692d6eeSBarry Smith #define MATPARTITIONINGSQUARE "square" 13362692d6eeSBarry Smith #define MATPARTITIONINGPARMETIS "parmetis" 13372692d6eeSBarry Smith #define MATPARTITIONINGCHACO "chaco" 13382692d6eeSBarry Smith #define MATPARTITIONINGJOSTLE "jostle" 13392692d6eeSBarry Smith #define MATPARTITIONINGPARTY "party" 13402692d6eeSBarry Smith #define MATPARTITIONINGSCOTCH "scotch" 134171306c60Sjroman 1342ca161407SBarry Smith 13437087cfbeSBarry Smith extern PetscErrorCode MatPartitioningCreate(MPI_Comm,MatPartitioning*); 13447087cfbeSBarry Smith extern PetscErrorCode MatPartitioningSetType(MatPartitioning,const MatPartitioningType); 13457087cfbeSBarry Smith extern PetscErrorCode MatPartitioningSetNParts(MatPartitioning,PetscInt); 13467087cfbeSBarry Smith extern PetscErrorCode MatPartitioningSetAdjacency(MatPartitioning,Mat); 13477087cfbeSBarry Smith extern PetscErrorCode MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]); 13487087cfbeSBarry Smith extern PetscErrorCode MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []); 13497087cfbeSBarry Smith extern PetscErrorCode MatPartitioningApply(MatPartitioning,IS*); 1350fcfd50ebSBarry Smith extern PetscErrorCode MatPartitioningDestroy(MatPartitioning*); 13512aabb6bbSBarry Smith 13527087cfbeSBarry Smith extern PetscErrorCode MatPartitioningRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatPartitioning)); 135330de9b25SBarry Smith 135430de9b25SBarry Smith /*MC 135530de9b25SBarry Smith MatPartitioningRegisterDynamic - Adds a new sparse matrix partitioning to the 135630de9b25SBarry Smith matrix package. 135730de9b25SBarry Smith 135830de9b25SBarry Smith Synopsis: 13591890ba74SBarry Smith PetscErrorCode MatPartitioningRegisterDynamic(const char *name_partitioning,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatPartitioning)) 136030de9b25SBarry Smith 136130de9b25SBarry Smith Not Collective 136230de9b25SBarry Smith 136330de9b25SBarry Smith Input Parameters: 13642692d6eeSBarry Smith + sname - name of partitioning (for example MATPARTITIONINGCURRENT) or parmetis 136530de9b25SBarry Smith . path - location of library where creation routine is 136630de9b25SBarry Smith . name - name of function that creates the partitioning type, a string 136730de9b25SBarry Smith - function - function pointer that creates the partitioning type 136830de9b25SBarry Smith 136930de9b25SBarry Smith Level: developer 137030de9b25SBarry Smith 137130de9b25SBarry Smith If dynamic libraries are used, then the fourth input argument (function) 137230de9b25SBarry Smith is ignored. 137330de9b25SBarry Smith 137430de9b25SBarry Smith Sample usage: 137530de9b25SBarry Smith .vb 137630de9b25SBarry Smith MatPartitioningRegisterDynamic("my_part",/home/username/my_lib/lib/libO/solaris/mylib.a, 137730de9b25SBarry Smith "MyPartCreate",MyPartCreate); 137830de9b25SBarry Smith .ve 137930de9b25SBarry Smith 138030de9b25SBarry Smith Then, your partitioner can be chosen with the procedural interface via 138130de9b25SBarry Smith $ MatPartitioningSetType(part,"my_part") 138230de9b25SBarry Smith or at runtime via the option 138330de9b25SBarry Smith $ -mat_partitioning_type my_part 138430de9b25SBarry Smith 1385ab901514SBarry Smith $PETSC_ARCH occuring in pathname will be replaced with appropriate values. 138630de9b25SBarry Smith 138730de9b25SBarry Smith .keywords: matrix, partitioning, register 138830de9b25SBarry Smith 138930de9b25SBarry Smith .seealso: MatPartitioningRegisterDestroy(), MatPartitioningRegisterAll() 139030de9b25SBarry Smith M*/ 1391aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 1392f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,0) 13932aabb6bbSBarry Smith #else 1394f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,d) 13952aabb6bbSBarry Smith #endif 13962aabb6bbSBarry Smith 1397ace3abfcSBarry Smith extern PetscBool MatPartitioningRegisterAllCalled; 1398f1144a30SSatish Balay 13997087cfbeSBarry Smith extern PetscErrorCode MatPartitioningRegisterAll(const char[]); 14007087cfbeSBarry Smith extern PetscErrorCode MatPartitioningRegisterDestroy(void); 14012bad1931SBarry Smith 14027087cfbeSBarry Smith extern PetscErrorCode MatPartitioningView(MatPartitioning,PetscViewer); 14037087cfbeSBarry Smith extern PetscErrorCode MatPartitioningSetFromOptions(MatPartitioning); 14047087cfbeSBarry Smith extern PetscErrorCode MatPartitioningGetType(MatPartitioning,const MatPartitioningType*); 1405ca161407SBarry Smith 14067087cfbeSBarry Smith extern PetscErrorCode MatPartitioningParmetisSetCoarseSequential(MatPartitioning); 14077087cfbeSBarry Smith extern PetscErrorCode MatPartitioningParmetisGetEdgeCut(MatPartitioning, PetscInt *); 14080752156aSBarry Smith 14097087cfbeSBarry Smith extern PetscErrorCode MatPartitioningJostleSetCoarseLevel(MatPartitioning,PetscReal); 14107087cfbeSBarry Smith extern PetscErrorCode MatPartitioningJostleSetCoarseSequential(MatPartitioning); 141171306c60Sjroman 1412954b3ec6SBarry Smith typedef enum {MP_CHACO_MULTILEVEL_KL,MP_CHACO_SPECTRAL,MP_CHACO_LINEAR,MP_CHACO_RANDOM, MP_CHACO_SCATTERED} MPChacoGlobalType; 14137087cfbeSBarry Smith extern PetscErrorCode MatPartitioningChacoSetGlobal(MatPartitioning, MPChacoGlobalType); 141471306c60Sjroman typedef enum { MP_CHACO_KERNIGHAN_LIN, MP_CHACO_NONE } MPChacoLocalType; 14157087cfbeSBarry Smith extern PetscErrorCode MatPartitioningChacoSetLocal(MatPartitioning, MPChacoLocalType); 14167087cfbeSBarry Smith extern PetscErrorCode MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal); 141771306c60Sjroman typedef enum { MP_CHACO_LANCZOS, MP_CHACO_RQI_SYMMLQ } MPChacoEigenType; 14187087cfbeSBarry Smith extern PetscErrorCode MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType); 14197087cfbeSBarry Smith extern PetscErrorCode MatPartitioningChacoSetEigenTol(MatPartitioning, PetscReal); 14207087cfbeSBarry Smith extern PetscErrorCode MatPartitioningChacoSetEigenNumber(MatPartitioning, PetscInt); 142171306c60Sjroman 142271306c60Sjroman #define MP_PARTY_OPT "opt" 142371306c60Sjroman #define MP_PARTY_LIN "lin" 142471306c60Sjroman #define MP_PARTY_SCA "sca" 142571306c60Sjroman #define MP_PARTY_RAN "ran" 142671306c60Sjroman #define MP_PARTY_GBF "gbf" 142771306c60Sjroman #define MP_PARTY_GCF "gcf" 142871306c60Sjroman #define MP_PARTY_BUB "bub" 142971306c60Sjroman #define MP_PARTY_DEF "def" 14307087cfbeSBarry Smith extern PetscErrorCode MatPartitioningPartySetGlobal(MatPartitioning, const char*); 143171306c60Sjroman #define MP_PARTY_HELPFUL_SETS "hs" 143271306c60Sjroman #define MP_PARTY_KERNIGHAN_LIN "kl" 143371306c60Sjroman #define MP_PARTY_NONE "no" 14347087cfbeSBarry Smith extern PetscErrorCode MatPartitioningPartySetLocal(MatPartitioning, const char*); 14357087cfbeSBarry Smith extern PetscErrorCode MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal); 14367087cfbeSBarry Smith extern PetscErrorCode MatPartitioningPartySetBipart(MatPartitioning,PetscBool ); 14377087cfbeSBarry Smith extern PetscErrorCode MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscBool ); 143871306c60Sjroman 143971306c60Sjroman typedef enum { MP_SCOTCH_GREEDY, MP_SCOTCH_GPS, MP_SCOTCH_GR_GPS } MPScotchGlobalType; 14407087cfbeSBarry Smith extern PetscErrorCode MatPartitioningScotchSetArch(MatPartitioning,const char*); 14417087cfbeSBarry Smith extern PetscErrorCode MatPartitioningScotchSetMultilevel(MatPartitioning); 14427087cfbeSBarry Smith extern PetscErrorCode MatPartitioningScotchSetGlobal(MatPartitioning,MPScotchGlobalType); 14437087cfbeSBarry Smith extern PetscErrorCode MatPartitioningScotchSetCoarseLevel(MatPartitioning,PetscReal); 14447087cfbeSBarry Smith extern PetscErrorCode MatPartitioningScotchSetHostList(MatPartitioning,const char*); 144571306c60Sjroman typedef enum { MP_SCOTCH_KERNIGHAN_LIN, MP_SCOTCH_NONE } MPScotchLocalType; 14467087cfbeSBarry Smith extern PetscErrorCode MatPartitioningScotchSetLocal(MatPartitioning,MPScotchLocalType); 14477087cfbeSBarry Smith extern PetscErrorCode MatPartitioningScotchSetMapping(MatPartitioning); 14487087cfbeSBarry Smith extern PetscErrorCode MatPartitioningScotchSetStrategy(MatPartitioning,char*); 144971306c60Sjroman 145009573ac7SBarry Smith extern PetscErrorCode MatMeshToVertexGraph(Mat,PetscInt,Mat*); 145109573ac7SBarry Smith extern PetscErrorCode MatMeshToCellGraph(Mat,PetscInt,Mat*); 1452591294e4SBarry Smith 14530752156aSBarry Smith /* 14540a835dfdSSatish Balay If you add entries here you must also add them to finclude/petscmat.h 1455d4fbbf0eSBarry Smith */ 14561c1c02c0SLois Curfman McInnes typedef enum { MATOP_SET_VALUES=0, 14571c1c02c0SLois Curfman McInnes MATOP_GET_ROW=1, 14581c1c02c0SLois Curfman McInnes MATOP_RESTORE_ROW=2, 14591c1c02c0SLois Curfman McInnes MATOP_MULT=3, 14601c1c02c0SLois Curfman McInnes MATOP_MULT_ADD=4, 14617c922b88SBarry Smith MATOP_MULT_TRANSPOSE=5, 14627c922b88SBarry Smith MATOP_MULT_TRANSPOSE_ADD=6, 14631c1c02c0SLois Curfman McInnes MATOP_SOLVE=7, 14641c1c02c0SLois Curfman McInnes MATOP_SOLVE_ADD=8, 14657c922b88SBarry Smith MATOP_SOLVE_TRANSPOSE=9, 14667c922b88SBarry Smith MATOP_SOLVE_TRANSPOSE_ADD=10, 14671c1c02c0SLois Curfman McInnes MATOP_LUFACTOR=11, 14681c1c02c0SLois Curfman McInnes MATOP_CHOLESKYFACTOR=12, 1469a714c06dSBarry Smith MATOP_SOR=13, 14701c1c02c0SLois Curfman McInnes MATOP_TRANSPOSE=14, 14711c1c02c0SLois Curfman McInnes MATOP_GETINFO=15, 14721c1c02c0SLois Curfman McInnes MATOP_EQUAL=16, 14731c1c02c0SLois Curfman McInnes MATOP_GET_DIAGONAL=17, 14741c1c02c0SLois Curfman McInnes MATOP_DIAGONAL_SCALE=18, 14751c1c02c0SLois Curfman McInnes MATOP_NORM=19, 14761c1c02c0SLois Curfman McInnes MATOP_ASSEMBLY_BEGIN=20, 14771c1c02c0SLois Curfman McInnes MATOP_ASSEMBLY_END=21, 1478d519adbfSMatthew Knepley MATOP_SET_OPTION=22, 1479d519adbfSMatthew Knepley MATOP_ZERO_ENTRIES=23, 1480d519adbfSMatthew Knepley MATOP_ZERO_ROWS=24, 1481d519adbfSMatthew Knepley MATOP_LUFACTOR_SYMBOLIC=25, 1482d519adbfSMatthew Knepley MATOP_LUFACTOR_NUMERIC=26, 1483d519adbfSMatthew Knepley MATOP_CHOLESKY_FACTOR_SYMBOLIC=27, 1484d519adbfSMatthew Knepley MATOP_CHOLESKY_FACTOR_NUMERIC=28, 1485d519adbfSMatthew Knepley MATOP_SETUP_PREALLOCATION=29, 1486d519adbfSMatthew Knepley MATOP_ILUFACTOR_SYMBOLIC=30, 1487d519adbfSMatthew Knepley MATOP_ICCFACTOR_SYMBOLIC=31, 1488d519adbfSMatthew Knepley MATOP_GET_ARRAY=32, 1489d519adbfSMatthew Knepley MATOP_RESTORE_ARRAY=33, 1490d519adbfSMatthew Knepley MATOP_DUPLICATE=34, 1491d519adbfSMatthew Knepley MATOP_FORWARD_SOLVE=35, 1492d519adbfSMatthew Knepley MATOP_BACKWARD_SOLVE=36, 1493d519adbfSMatthew Knepley MATOP_ILUFACTOR=37, 1494d519adbfSMatthew Knepley MATOP_ICCFACTOR=38, 1495d519adbfSMatthew Knepley MATOP_AXPY=39, 1496d519adbfSMatthew Knepley MATOP_GET_SUBMATRICES=40, 1497d519adbfSMatthew Knepley MATOP_INCREASE_OVERLAP=41, 1498d519adbfSMatthew Knepley MATOP_GET_VALUES=42, 1499d519adbfSMatthew Knepley MATOP_COPY=43, 1500d519adbfSMatthew Knepley MATOP_GET_ROW_MAX=44, 1501d519adbfSMatthew Knepley MATOP_SCALE=45, 1502d519adbfSMatthew Knepley MATOP_SHIFT=46, 150335153367SBarry Smith MATOP_DIAGONAL_SET=47, 1504d519adbfSMatthew Knepley MATOP_ILUDT_FACTOR=48, 1505d519adbfSMatthew Knepley MATOP_SET_BLOCK_SIZE=49, 1506d519adbfSMatthew Knepley MATOP_GET_ROW_IJ=50, 1507d519adbfSMatthew Knepley MATOP_RESTORE_ROW_IJ=51, 1508d519adbfSMatthew Knepley MATOP_GET_COLUMN_IJ=52, 1509d519adbfSMatthew Knepley MATOP_RESTORE_COLUMN_IJ=53, 1510d519adbfSMatthew Knepley MATOP_FDCOLORING_CREATE=54, 1511d519adbfSMatthew Knepley MATOP_COLORING_PATCH=55, 1512d519adbfSMatthew Knepley MATOP_SET_UNFACTORED=56, 1513d519adbfSMatthew Knepley MATOP_PERMUTE=57, 1514d519adbfSMatthew Knepley MATOP_SET_VALUES_BLOCKED=58, 1515d519adbfSMatthew Knepley MATOP_GET_SUBMATRIX=59, 1516d519adbfSMatthew Knepley MATOP_DESTROY=60, 1517d519adbfSMatthew Knepley MATOP_VIEW=61, 1518d519adbfSMatthew Knepley MATOP_CONVERT_FROM=62, 1519d519adbfSMatthew Knepley MATOP_USE_SCALED_FORM=63, 1520d519adbfSMatthew Knepley MATOP_SCALE_SYSTEM=64, 1521d519adbfSMatthew Knepley MATOP_UNSCALE_SYSTEM=65, 1522d519adbfSMatthew Knepley MATOP_SET_LOCAL_TO_GLOBAL_MAP=66, 1523d519adbfSMatthew Knepley MATOP_SET_VALUES_LOCAL=67, 1524d519adbfSMatthew Knepley MATOP_ZERO_ROWS_LOCAL=68, 1525d519adbfSMatthew Knepley MATOP_GET_ROW_MAX_ABS=69, 1526d519adbfSMatthew Knepley MATOP_GET_ROW_MIN_ABS=70, 1527d519adbfSMatthew Knepley MATOP_CONVERT=71, 1528d519adbfSMatthew Knepley MATOP_SET_COLORING=72, 1529d519adbfSMatthew Knepley MATOP_SET_VALUES_ADIC=73, 1530d519adbfSMatthew Knepley MATOP_SET_VALUES_ADIFOR=74, 1531d519adbfSMatthew Knepley MATOP_FD_COLORING_APPLY=75, 1532d519adbfSMatthew Knepley MATOP_SET_FROM_OPTIONS=76, 1533d519adbfSMatthew Knepley MATOP_MULT_CON=77, 1534d519adbfSMatthew Knepley MATOP_MULT_TRANSPOSE_CON=78, 1535d519adbfSMatthew Knepley MATOP_PERMUTE_SPARSIFY=79, 1536d519adbfSMatthew Knepley MATOP_MULT_MULTIPLE=80, 1537d519adbfSMatthew Knepley MATOP_SOLVE_MULTIPLE=81, 1538d519adbfSMatthew Knepley MATOP_GET_INERTIA=82, 1539d519adbfSMatthew Knepley MATOP_LOAD=83, 1540d519adbfSMatthew Knepley MATOP_IS_SYMMETRIC=84, 1541d519adbfSMatthew Knepley MATOP_IS_HERMITIAN=85, 1542d519adbfSMatthew Knepley MATOP_IS_STRUCTURALLY_SYMMETRIC=86, 154341f059aeSBarry Smith MATOP_DUMMY=87, 1544d519adbfSMatthew Knepley MATOP_GET_VECS=88, 1545d519adbfSMatthew Knepley MATOP_MAT_MULT=89, 1546d519adbfSMatthew Knepley MATOP_MAT_MULT_SYMBOLIC=90, 1547d519adbfSMatthew Knepley MATOP_MAT_MULT_NUMERIC=91, 1548d519adbfSMatthew Knepley MATOP_PTAP=92, 1549d519adbfSMatthew Knepley MATOP_PTAP_SYMBOLIC=93, 1550d519adbfSMatthew Knepley MATOP_PTAP_NUMERIC=94, 1551d519adbfSMatthew Knepley MATOP_MAT_MULTTRANSPOSE=95, 15521763ac4eSSatish Balay MATOP_MAT_MULTTRANSPOSE_SYM=96, 15531763ac4eSSatish Balay MATOP_MAT_MULTTRANSPOSE_NUM=97, 1554d519adbfSMatthew Knepley MATOP_PTAP_SYMBOLIC_SEQAIJ=98, 1555d519adbfSMatthew Knepley MATOP_PTAP_NUMERIC_SEQAIJ=99, 1556d519adbfSMatthew Knepley MATOP_PTAP_SYMBOLIC_MPIAIJ=100, 1557d519adbfSMatthew Knepley MATOP_PTAP_NUMERIC_MPIAIJ=101, 1558d519adbfSMatthew Knepley MATOP_CONJUGATE=102, 1559d519adbfSMatthew Knepley MATOP_SET_SIZES=103, 1560d519adbfSMatthew Knepley MATOP_SET_VALUES_ROW=104, 1561d519adbfSMatthew Knepley MATOP_REAL_PART=105, 1562d519adbfSMatthew Knepley MATOP_IMAG_PART=106, 1563d519adbfSMatthew Knepley MATOP_GET_ROW_UTRIANGULAR=107, 1564d519adbfSMatthew Knepley MATOP_RESTORE_ROW_UTRIANGULAR=108, 1565d519adbfSMatthew Knepley MATOP_MATSOLVE=109, 1566d519adbfSMatthew Knepley MATOP_GET_REDUNDANTMATRIX=110, 1567d519adbfSMatthew Knepley MATOP_GET_ROW_MIN=111, 1568d519adbfSMatthew Knepley MATOP_GET_COLUMN_VEC=112, 1569d519adbfSMatthew Knepley MATOP_MISSING_DIAGONAL=113, 1570d519adbfSMatthew Knepley MATOP_MATGETSEQNONZEROSTRUCTURE=114, 157189c6957cSBarry Smith MATOP_CREATE=115, 157289c6957cSBarry Smith MATOP_GET_GHOSTS=116, 1573ea38565cSJed Brown MATOP_GET_LOCALSUBMATRIX=117, 1574ea38565cSJed Brown MATOP_RESTORE_LOCALSUBMATRIX=118, 1575eeffb40dSHong Zhang MATOP_MULT_DIAGONAL_BLOCK=119, 1576547795f9SHong Zhang MATOP_HERMITIANTRANSPOSE=120, 1577547795f9SHong Zhang MATOP_MULTHERMITIANTRANSPOSE=121, 1578fbdbba38SShri Abhyankar MATOP_MULTHERMITIANTRANSPOSEADD=122, 1579b9614d88SDmitry Karpeev MATOP_GETMULTIPROCBLOCK=123, 15800716a85fSBarry Smith MATOP_GETCOLUMNNORMS=125, 1581b9614d88SDmitry Karpeev MATOP_GET_SUBMATRICES_PARALLEL=128 1582fae171e0SBarry Smith } MatOperation; 15837087cfbeSBarry Smith extern PetscErrorCode MatHasOperation(Mat,MatOperation,PetscBool *); 15847087cfbeSBarry Smith extern PetscErrorCode MatShellSetOperation(Mat,MatOperation,void(*)(void)); 15857087cfbeSBarry Smith extern PetscErrorCode MatShellGetOperation(Mat,MatOperation,void(**)(void)); 15867087cfbeSBarry Smith extern PetscErrorCode MatShellSetContext(Mat,void*); 1587112a2221SBarry Smith 158890ace30eSBarry Smith /* 158990ace30eSBarry Smith Codes for matrices stored on disk. By default they are 159090ace30eSBarry Smith stored in a universal format. By changing the format with 15917973ad04SBarry Smith PetscViewerSetFormat(viewer,PETSC_VIEWER_NATIVE); the matrices will 159290ace30eSBarry Smith be stored in a way natural for the matrix, for example dense matrices 159390ace30eSBarry Smith would be stored as dense. Matrices stored this way may only be 1594f4403165SShri Abhyankar read into matrices of the same type. 159590ace30eSBarry Smith */ 159690ace30eSBarry Smith #define MATRIX_BINARY_FORMAT_DENSE -1 159790ace30eSBarry Smith 15987087cfbeSBarry Smith extern PetscErrorCode MatMPIBAIJSetHashTableFactor(Mat,PetscReal); 15997087cfbeSBarry Smith extern PetscErrorCode MatISGetLocalMat(Mat,Mat*); 16001f4e1ec7SBarry Smith 1601d9274352SBarry Smith /*S 1602d9274352SBarry Smith MatNullSpace - Object that removes a null space from a vector, i.e. 1603d9274352SBarry Smith orthogonalizes the vector to a subsapce 1604d9274352SBarry Smith 1605f7a9e4ceSBarry Smith Level: advanced 1606d9274352SBarry Smith 1607d9274352SBarry Smith Concepts: matrix; linear operator, null space 1608d9274352SBarry Smith 16096e1639daSBarry Smith Users manual sections: 16106e1639daSBarry Smith . sec_singular 16116e1639daSBarry Smith 1612d9274352SBarry Smith .seealso: MatNullSpaceCreate() 1613d9274352SBarry Smith S*/ 161474637425SBarry Smith typedef struct _p_MatNullSpace* MatNullSpace; 1615d9274352SBarry Smith 16167087cfbeSBarry Smith extern PetscErrorCode MatNullSpaceCreate(MPI_Comm,PetscBool ,PetscInt,const Vec[],MatNullSpace*); 16177087cfbeSBarry Smith extern PetscErrorCode MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(MatNullSpace,Vec,void*),void*); 1618d34fcf5fSBarry Smith extern PetscErrorCode MatNullSpaceDestroy(MatNullSpace*); 16197087cfbeSBarry Smith extern PetscErrorCode MatNullSpaceRemove(MatNullSpace,Vec,Vec*); 16207087cfbeSBarry Smith extern PetscErrorCode MatNullSpaceAttach(Mat,MatNullSpace); 16217087cfbeSBarry Smith extern PetscErrorCode MatNullSpaceTest(MatNullSpace,Mat,PetscBool *); 1622b717e993SJed Brown extern PetscErrorCode MatNullSpaceView(MatNullSpace,PetscViewer); 162374637425SBarry Smith 16247087cfbeSBarry Smith extern PetscErrorCode MatReorderingSeqSBAIJ(Mat,IS); 16257087cfbeSBarry Smith extern PetscErrorCode MatMPISBAIJSetHashTableFactor(Mat,PetscReal); 16267087cfbeSBarry Smith extern PetscErrorCode MatSeqSBAIJSetColumnIndices(Mat,PetscInt *); 16277087cfbeSBarry Smith extern PetscErrorCode MatSeqBAIJInvertBlockDiagonal(Mat); 16283f1d51d7SBarry Smith 16297087cfbeSBarry Smith extern PetscErrorCode MatCreateMAIJ(Mat,PetscInt,Mat*); 16307087cfbeSBarry Smith extern PetscErrorCode MatMAIJRedimension(Mat,PetscInt,Mat*); 16317087cfbeSBarry Smith extern PetscErrorCode MatMAIJGetAIJ(Mat,Mat*); 1632c4f061fbSSatish Balay 16337087cfbeSBarry Smith extern PetscErrorCode MatComputeExplicitOperator(Mat,Mat*); 1634b0a32e0cSBarry Smith 16357087cfbeSBarry Smith extern PetscErrorCode MatDiagonalScaleLocal(Mat,Vec); 163604f1ad80SBarry Smith 16377087cfbeSBarry Smith extern PetscErrorCode MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*); 16387087cfbeSBarry Smith extern PetscErrorCode MatMFFDSetBase(Mat,Vec,Vec); 16397087cfbeSBarry Smith extern PetscErrorCode MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*); 16407087cfbeSBarry Smith extern PetscErrorCode MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*)); 16417087cfbeSBarry Smith extern PetscErrorCode MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec)); 16427087cfbeSBarry Smith extern PetscErrorCode MatMFFDAddNullSpace(Mat,MatNullSpace); 16437087cfbeSBarry Smith extern PetscErrorCode MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt); 16447087cfbeSBarry Smith extern PetscErrorCode MatMFFDResetHHistory(Mat); 16457087cfbeSBarry Smith extern PetscErrorCode MatMFFDSetFunctionError(Mat,PetscReal); 16467087cfbeSBarry Smith extern PetscErrorCode MatMFFDSetPeriod(Mat,PetscInt); 16477087cfbeSBarry Smith extern PetscErrorCode MatMFFDGetH(Mat,PetscScalar *); 16487087cfbeSBarry Smith extern PetscErrorCode MatMFFDSetOptionsPrefix(Mat,const char[]); 16497087cfbeSBarry Smith extern PetscErrorCode MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*); 16507087cfbeSBarry Smith extern PetscErrorCode MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*); 1651e884886fSBarry Smith 16526370053bSBarry Smith /*S 16536370053bSBarry Smith MatMFFD - A data structured used to manage the computation of the h differencing parameter for matrix-free 16546370053bSBarry Smith Jacobian vector products 1655e884886fSBarry Smith 16566370053bSBarry Smith Notes: MATMFFD is a specific MatType which uses the MatMFFD data structure 16576370053bSBarry Smith 16586370053bSBarry Smith MatMFFD*() methods actually take the Mat as their first argument. Not a MatMFFD data structure 16596370053bSBarry Smith 16606370053bSBarry Smith Level: developer 16616370053bSBarry Smith 16626370053bSBarry Smith .seealso: MATMFFD, MatCreateMFFD(), MatMFFDSetFuction(), MatMFFDSetType(), MatMFFDRegister() 16636370053bSBarry Smith S*/ 1664e884886fSBarry Smith typedef struct _p_MatMFFD* MatMFFD; 1665e884886fSBarry Smith 1666e884886fSBarry Smith /*E 1667e884886fSBarry Smith MatMFFDType - algorithm used to compute the h used in computing matrix-vector products via differencing of the function 1668e884886fSBarry Smith 1669e884886fSBarry Smith Level: beginner 1670e884886fSBarry Smith 1671e884886fSBarry Smith .seealso: MatMFFDSetType(), MatMFFDRegister() 1672e884886fSBarry Smith E*/ 1673a313700dSBarry Smith #define MatMFFDType char* 1674e884886fSBarry Smith #define MATMFFD_DS "ds" 1675e884886fSBarry Smith #define MATMFFD_WP "wp" 1676e884886fSBarry Smith 16777087cfbeSBarry Smith extern PetscErrorCode MatMFFDSetType(Mat,const MatMFFDType); 16787087cfbeSBarry Smith extern PetscErrorCode MatMFFDRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatMFFD)); 1679e884886fSBarry Smith 1680e884886fSBarry Smith /*MC 1681e884886fSBarry Smith MatMFFDRegisterDynamic - Adds a method to the MatMFFD registry. 1682e884886fSBarry Smith 1683e884886fSBarry Smith Synopsis: 16841890ba74SBarry Smith PetscErrorCode MatMFFDRegisterDynamic(const char *name_solver,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatMFFD)) 1685e884886fSBarry Smith 1686e884886fSBarry Smith Not Collective 1687e884886fSBarry Smith 1688e884886fSBarry Smith Input Parameters: 1689e884886fSBarry Smith + name_solver - name of a new user-defined compute-h module 1690e884886fSBarry Smith . path - path (either absolute or relative) the library containing this solver 1691e884886fSBarry Smith . name_create - name of routine to create method context 1692e884886fSBarry Smith - routine_create - routine to create method context 1693e884886fSBarry Smith 1694e884886fSBarry Smith Level: developer 1695e884886fSBarry Smith 1696e884886fSBarry Smith Notes: 1697e884886fSBarry Smith MatMFFDRegisterDynamic() may be called multiple times to add several user-defined solvers. 1698e884886fSBarry Smith 1699e884886fSBarry Smith If dynamic libraries are used, then the fourth input argument (routine_create) 1700e884886fSBarry Smith is ignored. 1701e884886fSBarry Smith 1702e884886fSBarry Smith Sample usage: 1703e884886fSBarry Smith .vb 1704e884886fSBarry Smith MatMFFDRegisterDynamic("my_h",/home/username/my_lib/lib/libO/solaris/mylib.a, 1705e884886fSBarry Smith "MyHCreate",MyHCreate); 1706e884886fSBarry Smith .ve 1707e884886fSBarry Smith 1708e884886fSBarry Smith Then, your solver can be chosen with the procedural interface via 1709e884886fSBarry Smith $ MatMFFDSetType(mfctx,"my_h") 1710e884886fSBarry Smith or at runtime via the option 1711e884886fSBarry Smith $ -snes_mf_type my_h 1712e884886fSBarry Smith 1713e884886fSBarry Smith .keywords: MatMFFD, register 1714e884886fSBarry Smith 1715e884886fSBarry Smith .seealso: MatMFFDRegisterAll(), MatMFFDRegisterDestroy() 1716e884886fSBarry Smith M*/ 1717e884886fSBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 1718e884886fSBarry Smith #define MatMFFDRegisterDynamic(a,b,c,d) MatMFFDRegister(a,b,c,0) 1719e884886fSBarry Smith #else 1720e884886fSBarry Smith #define MatMFFDRegisterDynamic(a,b,c,d) MatMFFDRegister(a,b,c,d) 1721e884886fSBarry Smith #endif 1722e884886fSBarry Smith 17237087cfbeSBarry Smith extern PetscErrorCode MatMFFDRegisterAll(const char[]); 17247087cfbeSBarry Smith extern PetscErrorCode MatMFFDRegisterDestroy(void); 17257087cfbeSBarry Smith extern PetscErrorCode MatMFFDDSSetUmin(Mat,PetscReal); 17267087cfbeSBarry Smith extern PetscErrorCode MatMFFDWPSetComputeNormU(Mat,PetscBool ); 1727e884886fSBarry Smith 1728e884886fSBarry Smith 17297087cfbeSBarry Smith extern PetscErrorCode PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *); 17307087cfbeSBarry Smith extern PetscErrorCode PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *); 17317dbadf16SMatthew Knepley 173297969023SHong Zhang /* 173397969023SHong Zhang PETSc interface to MUMPS 173497969023SHong Zhang */ 173597969023SHong Zhang #ifdef PETSC_HAVE_MUMPS 1736f250808bSBarry Smith extern PetscErrorCode MatMumpsSetIcntl(Mat,PetscInt,PetscInt); 173797969023SHong Zhang #endif 173823a5497aSJed Brown 1739d954db57SHong Zhang /* 1740d954db57SHong Zhang PETSc interface to SUPERLU 1741d954db57SHong Zhang */ 1742d954db57SHong Zhang #ifdef PETSC_HAVE_SUPERLU 1743f250808bSBarry Smith extern PetscErrorCode MatSuperluSetILUDropTol(Mat,PetscReal); 1744d954db57SHong Zhang #endif 1745d954db57SHong Zhang 174690c53902SBarry Smith #if defined(PETSC_HAVE_CUSP) 174790c53902SBarry Smith extern PetscErrorCode MatCreateSeqAIJCUSP(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 174890c53902SBarry Smith extern PetscErrorCode MatCreateMPIAIJCUSP(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 174990c53902SBarry Smith #endif 175090c53902SBarry Smith 175154efbe56SHong Zhang /* 175254efbe56SHong Zhang PETSc interface to FFTW 175354efbe56SHong Zhang */ 175454efbe56SHong Zhang #if defined(PETSC_HAVE_FFTW) 175574a26884SAmlan Barua extern PetscErrorCode VecScatterPetscToFFTW(Mat,Vec,Vec); 175674a26884SAmlan Barua extern PetscErrorCode VecScatterFFTWToPetsc(Mat,Vec,Vec); 17574be45526SHong Zhang extern PetscErrorCode MatGetVecsFFTW(Mat,Vec*,Vec*,Vec*); 175854efbe56SHong Zhang #endif 175954efbe56SHong Zhang 1760c8883902SJed Brown extern PetscErrorCode MatCreateNest(MPI_Comm,PetscInt,const IS[],PetscInt,const IS[],const Mat[],Mat*); 1761c8883902SJed Brown extern PetscErrorCode MatNestGetSize(Mat,PetscInt*,PetscInt*); 1762c8883902SJed Brown extern PetscErrorCode MatNestGetSubMats(Mat,PetscInt*,PetscInt*,Mat***); 1763c8883902SJed Brown extern PetscErrorCode MatNestGetSubMat(Mat,PetscInt,PetscInt,Mat*); 17647087cfbeSBarry Smith extern PetscErrorCode MatNestSetVecType(Mat,const VecType); 1765c8883902SJed Brown extern PetscErrorCode MatNestSetSubMats(Mat,PetscInt,const IS[],PetscInt,const IS[],const Mat[]); 17660782ca92SJed Brown extern PetscErrorCode MatNestSetSubMat(Mat,PetscInt,PetscInt,Mat); 1767d8588912SDave May 176823a5497aSJed Brown PETSC_EXTERN_CXX_END 176923a5497aSJed Brown #endif 1770