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" 37*51d315f7SKerry 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*/ 2413b224e63SBarry Smith typedef enum {SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER,SUBSET_NONZERO_PATTERN} 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*); 36899cafbc1SBarry Smith 3698ed539a5SBarry Smith /* ------------------------------------------------------------*/ 3707087cfbeSBarry Smith extern PetscErrorCode MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 3717087cfbeSBarry Smith extern PetscErrorCode MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 3727087cfbeSBarry Smith extern PetscErrorCode MatSetValuesRow(Mat,PetscInt,const PetscScalar[]); 3737087cfbeSBarry Smith extern PetscErrorCode MatSetValuesRowLocal(Mat,PetscInt,const PetscScalar[]); 37484cb2905SBarry Smith 3752ef4de8bSBarry Smith /*S 3762ef4de8bSBarry Smith MatStencil - Data structure (C struct) for storing information about a single row or 3772ef4de8bSBarry Smith column of a matrix as index on an associated grid. 3782ef4de8bSBarry Smith 3792ef4de8bSBarry Smith Level: beginner 3802ef4de8bSBarry Smith 3812ef4de8bSBarry Smith Concepts: matrix; linear operator 3822ef4de8bSBarry Smith 383d7bd93baSBarry Smith .seealso: MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockStencil() 3842ef4de8bSBarry Smith S*/ 385435da068SBarry Smith typedef struct { 386c1ac3661SBarry Smith PetscInt k,j,i,c; 387435da068SBarry Smith } MatStencil; 3882ef4de8bSBarry Smith 3897087cfbeSBarry Smith extern PetscErrorCode MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode); 3907087cfbeSBarry Smith extern PetscErrorCode MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode); 3917087cfbeSBarry Smith extern PetscErrorCode MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt); 392435da068SBarry Smith 3937087cfbeSBarry Smith extern PetscErrorCode MatSetColoring(Mat,ISColoring); 3947087cfbeSBarry Smith extern PetscErrorCode MatSetValuesAdic(Mat,void*); 3957087cfbeSBarry Smith extern PetscErrorCode MatSetValuesAdifor(Mat,PetscInt,void*); 3963a7fca6bSBarry Smith 397d91e6319SBarry Smith /*E 398d91e6319SBarry Smith MatAssemblyType - Indicates if the matrix is now to be used, or if you plan 399d91e6319SBarry Smith to continue to add values to it 400d91e6319SBarry Smith 401d91e6319SBarry Smith Level: beginner 402d91e6319SBarry Smith 403d91e6319SBarry Smith .seealso: MatAssemblyBegin(), MatAssemblyEnd() 404d91e6319SBarry Smith E*/ 4056d4a8577SBarry Smith typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType; 4067087cfbeSBarry Smith extern PetscErrorCode MatAssemblyBegin(Mat,MatAssemblyType); 4077087cfbeSBarry Smith extern PetscErrorCode MatAssemblyEnd(Mat,MatAssemblyType); 4087087cfbeSBarry Smith extern PetscErrorCode MatAssembled(Mat,PetscBool *); 4094f9c727eSBarry Smith 4101d73ed98SBarry Smith 41130de9b25SBarry Smith 412d91e6319SBarry Smith /*E 413d91e6319SBarry Smith MatOption - Options that may be set for a matrix and its behavior or storage 414d91e6319SBarry Smith 415d91e6319SBarry Smith Level: beginner 416d91e6319SBarry Smith 4170a835dfdSSatish Balay Any additions/changes here MUST also be made in include/finclude/petscmat.h 418d91e6319SBarry Smith 419d91e6319SBarry Smith .seealso: MatSetOption() 420d91e6319SBarry Smith E*/ 421512a5fc5SBarry Smith typedef enum {MAT_ROW_ORIENTED,MAT_NEW_NONZERO_LOCATIONS, 4224e0d8c25SBarry Smith MAT_SYMMETRIC, 4234e0d8c25SBarry Smith MAT_STRUCTURALLY_SYMMETRIC, 4248047e77bSBarry Smith MAT_NEW_DIAGONALS,MAT_IGNORE_OFF_PROC_ENTRIES, 4254e0d8c25SBarry Smith MAT_NEW_NONZERO_LOCATION_ERR, 4264e0d8c25SBarry Smith MAT_NEW_NONZERO_ALLOCATION_ERR,MAT_USE_HASH_TABLE, 427a9817697SBarry Smith MAT_KEEP_NONZERO_PATTERN,MAT_IGNORE_ZERO_ENTRIES, 4284e0d8c25SBarry Smith MAT_USE_INODES, 4294e0d8c25SBarry Smith MAT_HERMITIAN, 4304e0d8c25SBarry Smith MAT_SYMMETRY_ETERNAL, 431cd6b891eSBarry Smith MAT_CHECK_COMPRESSED_ROW, 4324e0d8c25SBarry Smith MAT_IGNORE_LOWER_TRIANGULAR,MAT_ERROR_LOWER_TRIANGULAR, 43328b2fa4aSMatthew Knepley MAT_GETROW_UPPERTRIANGULAR,MAT_UNUSED_NONZERO_LOCATION_ERR, 4344cb17eb5SBarry Smith MAT_SPD,MAT_NO_OFF_PROC_ENTRIES,MAT_NO_OFF_PROC_ZERO_ROWS, 43528b2fa4aSMatthew Knepley NUM_MAT_OPTIONS} MatOption; 436290bbb0aSBarry Smith extern const char *MatOptions[]; 4377087cfbeSBarry Smith extern PetscErrorCode MatSetOption(Mat,MatOption,PetscBool ); 4387087cfbeSBarry Smith extern PetscErrorCode MatGetType(Mat,const MatType*); 439a313700dSBarry Smith PetscPolymorphicFunction(MatGetType,(Mat mat),(mat,&t),const MatType,t) 44084cb2905SBarry Smith 4417087cfbeSBarry Smith extern PetscErrorCode MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]); 4427087cfbeSBarry Smith extern PetscErrorCode MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 4437087cfbeSBarry Smith extern PetscErrorCode MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 4447087cfbeSBarry Smith extern PetscErrorCode MatGetRowUpperTriangular(Mat); 4457087cfbeSBarry Smith extern PetscErrorCode MatRestoreRowUpperTriangular(Mat); 4467087cfbeSBarry Smith extern PetscErrorCode MatGetColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 4477087cfbeSBarry Smith extern PetscErrorCode MatRestoreColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 4487087cfbeSBarry Smith extern PetscErrorCode MatGetColumnVector(Mat,Vec,PetscInt); 4497087cfbeSBarry Smith extern PetscErrorCode MatGetArray(Mat,PetscScalar *[]); 450ba966639SSatish Balay PetscPolymorphicFunction(MatGetArray,(Mat mat),(mat,&a),PetscScalar*,a) 4517087cfbeSBarry Smith extern PetscErrorCode MatRestoreArray(Mat,PetscScalar *[]); 4527087cfbeSBarry Smith extern PetscErrorCode MatGetBlockSize(Mat,PetscInt *); 453ba966639SSatish Balay PetscPolymorphicFunction(MatGetBlockSize,(Mat mat),(mat,&a),PetscInt,a) 4547087cfbeSBarry Smith extern PetscErrorCode MatSetBlockSize(Mat,PetscInt); 4557b80b807SBarry Smith 4561620fd73SBarry Smith 4577087cfbeSBarry Smith extern PetscErrorCode MatMult(Mat,Vec,Vec); 4587087cfbeSBarry Smith extern PetscErrorCode MatMultDiagonalBlock(Mat,Vec,Vec); 4597087cfbeSBarry Smith extern PetscErrorCode MatMultAdd(Mat,Vec,Vec,Vec); 460ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultAdd,(Mat A,Vec x,Vec y),(A,x,y,y)) 4617087cfbeSBarry Smith extern PetscErrorCode MatMultTranspose(Mat,Vec,Vec); 4627087cfbeSBarry Smith extern PetscErrorCode MatMultHermitianTranspose(Mat,Vec,Vec); 4637087cfbeSBarry Smith extern PetscErrorCode MatIsTranspose(Mat,Mat,PetscReal,PetscBool *); 464ace3abfcSBarry Smith PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B,PetscReal tol),(A,B,tol,&t),PetscBool ,t) 465ace3abfcSBarry Smith PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B),(A,B,0,&t),PetscBool ,t) 4667087cfbeSBarry Smith extern PetscErrorCode MatIsHermitianTranspose(Mat,Mat,PetscReal,PetscBool *); 4677087cfbeSBarry Smith extern PetscErrorCode MatMultTransposeAdd(Mat,Vec,Vec,Vec); 4687087cfbeSBarry Smith extern PetscErrorCode MatMultHermitianTransposeAdd(Mat,Vec,Vec,Vec); 469ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultTransposeAdd,(Mat A,Vec x,Vec y),(A,x,y,y)) 4707087cfbeSBarry Smith extern PetscErrorCode MatMultConstrained(Mat,Vec,Vec); 4717087cfbeSBarry Smith extern PetscErrorCode MatMultTransposeConstrained(Mat,Vec,Vec); 4727087cfbeSBarry Smith extern PetscErrorCode MatMatSolve(Mat,Mat,Mat); 473f5edf698SHong Zhang 474d91e6319SBarry Smith /*E 475d91e6319SBarry Smith MatDuplicateOption - Indicates if a duplicated sparse matrix should have 476d91e6319SBarry Smith its numerical values copied over or just its nonzero structure. 477d91e6319SBarry Smith 478d91e6319SBarry Smith Level: beginner 479d91e6319SBarry Smith 480d91e6319SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 481d91e6319SBarry Smith 48270dcbbb9SBarry Smith $ MAT_SHARE_NONZERO_PATTERN - the i and j arrays in the new matrix will be shared with the original matrix 48370dcbbb9SBarry Smith $ this also triggers the MAT_DO_NOT_COPY_VALUES option. This is used when you 48470dcbbb9SBarry Smith $ have several matrices with the same nonzero pattern. 48570dcbbb9SBarry Smith 486d91e6319SBarry Smith .seealso: MatDuplicate() 487d91e6319SBarry Smith E*/ 48870dcbbb9SBarry Smith typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES,MAT_SHARE_NONZERO_PATTERN} MatDuplicateOption; 4892e8a6d31SBarry Smith 4907087cfbeSBarry Smith extern PetscErrorCode MatConvert(Mat,const MatType,MatReuse,Mat*); 491a313700dSBarry Smith PetscPolymorphicFunction(MatConvert,(Mat A,const MatType t),(A,t,MAT_INITIAL_MATRIX,&a),Mat,a) 4927087cfbeSBarry Smith extern PetscErrorCode MatDuplicate(Mat,MatDuplicateOption,Mat*); 493ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A,MatDuplicateOption o),(A,o,&a),Mat,a) 494ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A),(A,MAT_COPY_VALUES,&a),Mat,a) 49594a9d846SBarry Smith 496cb5b572fSBarry Smith 4977087cfbeSBarry Smith extern PetscErrorCode MatCopy(Mat,Mat,MatStructure); 4987087cfbeSBarry Smith extern PetscErrorCode MatView(Mat,PetscViewer); 4997087cfbeSBarry Smith extern PetscErrorCode MatIsSymmetric(Mat,PetscReal,PetscBool *); 500ace3abfcSBarry Smith PetscPolymorphicFunction(MatIsSymmetric,(Mat A,PetscReal tol),(A,tol,&t),PetscBool ,t) 501ace3abfcSBarry Smith PetscPolymorphicFunction(MatIsSymmetric,(Mat A),(A,0,&t),PetscBool ,t) 5027087cfbeSBarry Smith extern PetscErrorCode MatIsStructurallySymmetric(Mat,PetscBool *); 503ace3abfcSBarry Smith PetscPolymorphicFunction(MatIsStructurallySymmetric,(Mat A),(A,&t),PetscBool ,t) 5047087cfbeSBarry Smith extern PetscErrorCode MatIsHermitian(Mat,PetscReal,PetscBool *); 505ace3abfcSBarry Smith PetscPolymorphicFunction(MatIsHermitian,(Mat A),(A,0,&t),PetscBool ,t) 5067087cfbeSBarry Smith extern PetscErrorCode MatIsSymmetricKnown(Mat,PetscBool *,PetscBool *); 5077087cfbeSBarry Smith extern PetscErrorCode MatIsHermitianKnown(Mat,PetscBool *,PetscBool *); 5087087cfbeSBarry Smith extern PetscErrorCode MatMissingDiagonal(Mat,PetscBool *,PetscInt *); 5097087cfbeSBarry Smith extern PetscErrorCode MatLoad(Mat, PetscViewer); 5107b80b807SBarry Smith 5117087cfbeSBarry Smith extern PetscErrorCode MatGetRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,PetscInt *[],PetscInt *[],PetscBool *); 5127087cfbeSBarry Smith extern PetscErrorCode MatRestoreRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,PetscInt *[],PetscInt *[],PetscBool *); 5137087cfbeSBarry Smith extern PetscErrorCode MatGetColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,PetscInt *[],PetscInt *[],PetscBool *); 5147087cfbeSBarry Smith extern PetscErrorCode MatRestoreColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,PetscInt *[],PetscInt *[],PetscBool *); 515d4fbbf0eSBarry Smith 516d91e6319SBarry Smith /*S 517d91e6319SBarry Smith MatInfo - Context of matrix information, used with MatGetInfo() 518d91e6319SBarry Smith 519d91e6319SBarry Smith In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE 520d91e6319SBarry Smith 521d91e6319SBarry Smith Level: intermediate 522d91e6319SBarry Smith 523d91e6319SBarry Smith Concepts: matrix^nonzero information 524d91e6319SBarry Smith 525d9274352SBarry Smith .seealso: MatGetInfo(), MatInfoType 526d91e6319SBarry Smith S*/ 5274e220ebcSLois Curfman McInnes typedef struct { 528b0a32e0cSBarry Smith PetscLogDouble block_size; /* block size */ 529b0a32e0cSBarry Smith PetscLogDouble nz_allocated,nz_used,nz_unneeded; /* number of nonzeros */ 530b0a32e0cSBarry Smith PetscLogDouble memory; /* memory allocated */ 531b0a32e0cSBarry Smith PetscLogDouble assemblies; /* number of matrix assemblies called */ 532b0a32e0cSBarry Smith PetscLogDouble mallocs; /* number of mallocs during MatSetValues() */ 533b0a32e0cSBarry Smith PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */ 534b0a32e0cSBarry Smith PetscLogDouble factor_mallocs; /* number of mallocs during factorization */ 5354e220ebcSLois Curfman McInnes } MatInfo; 5364e220ebcSLois Curfman McInnes 537d9274352SBarry Smith /*E 538d9274352SBarry Smith MatInfoType - Indicates if you want information about the local part of the matrix, 539d9274352SBarry Smith the entire parallel matrix or the maximum over all the local parts. 540d9274352SBarry Smith 541d9274352SBarry Smith Level: beginner 542d9274352SBarry Smith 543d9274352SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 544d9274352SBarry Smith 545d9274352SBarry Smith .seealso: MatGetInfo(), MatInfo 546d9274352SBarry Smith E*/ 5477b80b807SBarry Smith typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType; 5487087cfbeSBarry Smith extern PetscErrorCode MatGetInfo(Mat,MatInfoType,MatInfo*); 5497087cfbeSBarry Smith extern PetscErrorCode MatGetDiagonal(Mat,Vec); 5507087cfbeSBarry Smith extern PetscErrorCode MatGetRowMax(Mat,Vec,PetscInt[]); 5517087cfbeSBarry Smith extern PetscErrorCode MatGetRowMin(Mat,Vec,PetscInt[]); 5527087cfbeSBarry Smith extern PetscErrorCode MatGetRowMaxAbs(Mat,Vec,PetscInt[]); 5537087cfbeSBarry Smith extern PetscErrorCode MatGetRowMinAbs(Mat,Vec,PetscInt[]); 5547087cfbeSBarry Smith extern PetscErrorCode MatGetRowSum(Mat,Vec); 5557087cfbeSBarry Smith extern PetscErrorCode MatTranspose(Mat,MatReuse,Mat*); 556fc4dec0aSBarry Smith PetscPolymorphicFunction(MatTranspose,(Mat A),(A,MAT_INITIAL_MATRIX,&t),Mat,t) 5577087cfbeSBarry Smith extern PetscErrorCode MatHermitianTranspose(Mat,MatReuse,Mat*); 5587087cfbeSBarry Smith extern PetscErrorCode MatPermute(Mat,IS,IS,Mat *); 559ba966639SSatish Balay PetscPolymorphicFunction(MatPermute,(Mat A,IS is1,IS is2),(A,is1,is2,&t),Mat,t) 5607087cfbeSBarry Smith extern PetscErrorCode MatDiagonalScale(Mat,Vec,Vec); 5617087cfbeSBarry Smith extern PetscErrorCode MatDiagonalSet(Mat,Vec,InsertMode); 5627087cfbeSBarry Smith extern PetscErrorCode MatEqual(Mat,Mat,PetscBool *); 563ace3abfcSBarry Smith PetscPolymorphicFunction(MatEqual,(Mat A,Mat B),(A,B,&t),PetscBool ,t) 5647087cfbeSBarry Smith extern PetscErrorCode MatMultEqual(Mat,Mat,PetscInt,PetscBool *); 5657087cfbeSBarry Smith extern PetscErrorCode MatMultAddEqual(Mat,Mat,PetscInt,PetscBool *); 5667087cfbeSBarry Smith extern PetscErrorCode MatMultTransposeEqual(Mat,Mat,PetscInt,PetscBool *); 5677087cfbeSBarry Smith extern PetscErrorCode MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscBool *); 5687b80b807SBarry Smith 5697087cfbeSBarry Smith extern PetscErrorCode MatNorm(Mat,NormType,PetscReal *); 570ba966639SSatish Balay PetscPolymorphicFunction(MatNorm,(Mat A,NormType t),(A,t,&n),PetscReal,n) 57109573ac7SBarry Smith extern PetscErrorCode MatGetColumnNorms(Mat,NormType,PetscReal *); 5727087cfbeSBarry Smith extern PetscErrorCode MatZeroEntries(Mat); 5737087cfbeSBarry Smith extern PetscErrorCode MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 5747087cfbeSBarry Smith extern PetscErrorCode MatZeroRowsIS(Mat,IS,PetscScalar,Vec,Vec); 5757087cfbeSBarry Smith extern PetscErrorCode MatZeroRowsStencil(Mat,PetscInt,const MatStencil [],PetscScalar,Vec,Vec); 5767087cfbeSBarry Smith extern PetscErrorCode MatZeroRowsColumns(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 5777087cfbeSBarry Smith extern PetscErrorCode MatZeroRowsColumnsIS(Mat,IS,PetscScalar,Vec,Vec); 5787b80b807SBarry Smith 5797087cfbeSBarry Smith extern PetscErrorCode MatUseScaledForm(Mat,PetscBool ); 5807087cfbeSBarry Smith extern PetscErrorCode MatScaleSystem(Mat,Vec,Vec); 5817087cfbeSBarry Smith extern PetscErrorCode MatUnScaleSystem(Mat,Vec,Vec); 5825ef9f2a5SBarry Smith 5837087cfbeSBarry Smith extern PetscErrorCode MatGetSize(Mat,PetscInt*,PetscInt*); 5847087cfbeSBarry Smith extern PetscErrorCode MatGetLocalSize(Mat,PetscInt*,PetscInt*); 5857087cfbeSBarry Smith extern PetscErrorCode MatGetOwnershipRange(Mat,PetscInt*,PetscInt*); 5867087cfbeSBarry Smith extern PetscErrorCode MatGetOwnershipRanges(Mat,const PetscInt**); 5877087cfbeSBarry Smith extern PetscErrorCode MatGetOwnershipRangeColumn(Mat,PetscInt*,PetscInt*); 5887087cfbeSBarry Smith extern PetscErrorCode MatGetOwnershipRangesColumn(Mat,const PetscInt**); 5897b80b807SBarry Smith 5907087cfbeSBarry Smith extern PetscErrorCode MatGetSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]); 5917087cfbeSBarry Smith extern PetscErrorCode MatGetSubMatricesParallel(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]); 5927087cfbeSBarry Smith extern PetscErrorCode MatDestroyMatrices(PetscInt,Mat *[]); 5937087cfbeSBarry Smith extern PetscErrorCode MatGetSubMatrix(Mat,IS,IS,MatReuse,Mat *); 5947087cfbeSBarry Smith extern PetscErrorCode MatGetLocalSubMatrix(Mat,IS,IS,Mat*); 5957087cfbeSBarry Smith extern PetscErrorCode MatRestoreLocalSubMatrix(Mat,IS,IS,Mat*); 5967087cfbeSBarry Smith extern PetscErrorCode MatGetSeqNonzeroStructure(Mat,Mat*); 5977087cfbeSBarry Smith extern PetscErrorCode MatDestroySeqNonzeroStructure(Mat*); 5985494a064SHong Zhang 5997087cfbeSBarry Smith extern PetscErrorCode MatMerge(MPI_Comm,Mat,PetscInt,MatReuse,Mat*); 6007087cfbeSBarry Smith extern PetscErrorCode MatMerge_SeqsToMPI(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*); 6017087cfbeSBarry Smith extern PetscErrorCode MatMerge_SeqsToMPISymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*); 6027087cfbeSBarry Smith extern PetscErrorCode MatMerge_SeqsToMPINumeric(Mat,Mat); 6034a2b5492SBarry Smith extern PetscErrorCode MatMPIAIJGetLocalMat(Mat,MatReuse,Mat*); 6044a2b5492SBarry Smith extern PetscErrorCode MatMPIAIJGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*); 6057087cfbeSBarry Smith extern PetscErrorCode MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,PetscInt*,Mat*); 6067087cfbeSBarry Smith extern PetscErrorCode MatGetBrowsOfAoCols(Mat,Mat,MatReuse,PetscInt**,PetscInt**,MatScalar**,Mat*); 60743eb5e2fSMatthew Knepley #if defined (PETSC_USE_CTABLE) 608cd116e26SSatish Balay #include "petscctable.h" 6097087cfbeSBarry Smith extern PetscErrorCode MatGetCommunicationStructs(Mat, Vec *, PetscTable *, VecScatter *); 61043eb5e2fSMatthew Knepley #else 6117087cfbeSBarry Smith extern PetscErrorCode MatGetCommunicationStructs(Mat, Vec *, PetscInt *[], VecScatter *); 61243eb5e2fSMatthew Knepley #endif 6137087cfbeSBarry Smith extern PetscErrorCode MatGetGhosts(Mat, PetscInt *,const PetscInt *[]); 6148efafbd8SBarry Smith 6157087cfbeSBarry Smith extern PetscErrorCode MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt); 6167b80b807SBarry Smith 6177087cfbeSBarry Smith extern PetscErrorCode MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*); 6187087cfbeSBarry Smith extern PetscErrorCode MatMatMultSymbolic(Mat,Mat,PetscReal,Mat*); 6197087cfbeSBarry Smith extern PetscErrorCode MatMatMultNumeric(Mat,Mat,Mat); 62022440eb1SKris Buschelman 6217087cfbeSBarry Smith extern PetscErrorCode MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*); 6227087cfbeSBarry Smith extern PetscErrorCode MatPtAPSymbolic(Mat,Mat,PetscReal,Mat*); 6237087cfbeSBarry Smith extern PetscErrorCode MatPtAPNumeric(Mat,Mat,Mat); 62422440eb1SKris Buschelman 6257087cfbeSBarry Smith extern PetscErrorCode MatMatMultTranspose(Mat,Mat,MatReuse,PetscReal,Mat*); 6267087cfbeSBarry Smith extern PetscErrorCode MatMatMultTransposeSymbolic(Mat,Mat,PetscReal,Mat*); 6277087cfbeSBarry Smith extern PetscErrorCode MatMatMultTransposeNumeric(Mat,Mat,Mat); 628bc011b1eSHong Zhang 6297087cfbeSBarry Smith extern PetscErrorCode MatAXPY(Mat,PetscScalar,Mat,MatStructure); 6307087cfbeSBarry Smith extern PetscErrorCode MatAYPX(Mat,PetscScalar,Mat,MatStructure); 6311c741599SHong Zhang 6327087cfbeSBarry Smith extern PetscErrorCode MatScale(Mat,PetscScalar); 6337087cfbeSBarry Smith extern PetscErrorCode MatShift(Mat,PetscScalar); 6347b80b807SBarry Smith 6357087cfbeSBarry Smith extern PetscErrorCode MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping); 6367087cfbeSBarry Smith extern PetscErrorCode MatSetLocalToGlobalMappingBlock(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping); 6377087cfbeSBarry Smith extern PetscErrorCode MatGetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*); 6387087cfbeSBarry Smith extern PetscErrorCode MatGetLocalToGlobalMappingBlock(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*); 6397087cfbeSBarry Smith extern PetscErrorCode MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 6407087cfbeSBarry Smith extern PetscErrorCode MatZeroRowsLocalIS(Mat,IS,PetscScalar,Vec,Vec); 6417087cfbeSBarry Smith extern PetscErrorCode MatZeroRowsColumnsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 6427087cfbeSBarry Smith extern PetscErrorCode MatZeroRowsColumnsLocalIS(Mat,IS,PetscScalar,Vec,Vec); 6437087cfbeSBarry Smith extern PetscErrorCode MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 6447087cfbeSBarry Smith extern PetscErrorCode MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 645052efed2SBarry Smith 6467087cfbeSBarry Smith extern PetscErrorCode MatStashSetInitialSize(Mat,PetscInt,PetscInt); 6477087cfbeSBarry Smith extern PetscErrorCode MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*); 64890f02eecSBarry Smith 6497087cfbeSBarry Smith extern PetscErrorCode MatInterpolate(Mat,Vec,Vec); 6507087cfbeSBarry Smith extern PetscErrorCode MatInterpolateAdd(Mat,Vec,Vec,Vec); 651ba966639SSatish Balay PetscPolymorphicSubroutine(MatInterpolateAdd,(Mat A,Vec x,Vec y),(A,x,y,y)) 6527087cfbeSBarry Smith extern PetscErrorCode MatRestrict(Mat,Vec,Vec); 6537087cfbeSBarry Smith extern PetscErrorCode MatGetVecs(Mat,Vec*,Vec*); 6547087cfbeSBarry Smith extern PetscErrorCode MatGetRedundantMatrix(Mat,PetscInt,MPI_Comm,PetscInt,MatReuse,Mat*); 6557087cfbeSBarry Smith extern PetscErrorCode MatGetMultiProcBlock(Mat,MPI_Comm,Mat*); 656fc9bc008SSatish Balay extern PetscErrorCode MatFindZeroDiagonals(Mat,IS*); 657bd481603SBarry Smith 658bd481603SBarry Smith /*MC 6591620fd73SBarry Smith MatSetValue - Set a single entry into a matrix. 6601620fd73SBarry Smith 6611620fd73SBarry Smith Not collective 6621620fd73SBarry Smith 6631620fd73SBarry Smith Input Parameters: 6641620fd73SBarry Smith + m - the matrix 6651620fd73SBarry Smith . row - the row location of the entry 6661620fd73SBarry Smith . col - the column location of the entry 6671620fd73SBarry Smith . value - the value to insert 6681620fd73SBarry Smith - mode - either INSERT_VALUES or ADD_VALUES 6691620fd73SBarry Smith 6701620fd73SBarry Smith Notes: 6711620fd73SBarry Smith For efficiency one should use MatSetValues() and set several or many 6721620fd73SBarry Smith values simultaneously if possible. 6731620fd73SBarry Smith 6741620fd73SBarry Smith Level: beginner 6751620fd73SBarry Smith 6761620fd73SBarry Smith .seealso: MatSetValues(), MatSetValueLocal() 6771620fd73SBarry Smith M*/ 6781620fd73SBarry 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);} 6791620fd73SBarry Smith 6801620fd73SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatGetValue(Mat v,PetscInt i,PetscInt j,PetscScalar *va) {return MatGetValues(v,1,&i,1,&j,va);} 6811620fd73SBarry Smith 6821620fd73SBarry 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);} 6831620fd73SBarry Smith 6841620fd73SBarry Smith /*MC 6850d2a7ff2SSatish Balay MatPreallocateInitialize - Begins the block of code that will count the number of nonzeros per 686bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 687bd481603SBarry Smith 688bd481603SBarry Smith Synopsis: 689c1ac3661SBarry Smith PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz) 690bd481603SBarry Smith 691bd481603SBarry Smith Collective on MPI_Comm 692bd481603SBarry Smith 693bd481603SBarry Smith Input Parameters: 694bd481603SBarry Smith + comm - the communicator that will share the eventually allocated matrix 695859fcb39SBarry Smith . nrows - the number of LOCAL rows in the matrix 696859fcb39SBarry Smith - ncols - the number of LOCAL columns in the matrix 697bd481603SBarry Smith 698bd481603SBarry Smith Output Parameters: 699bd481603SBarry Smith + dnz - the array that will be passed to the matrix preallocation routines 700bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 701bd481603SBarry Smith 702bd481603SBarry Smith 703bd481603SBarry Smith Level: intermediate 704bd481603SBarry Smith 705bd481603SBarry Smith Notes: 7060598bfebSBarry Smith See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details. 707bd481603SBarry Smith 7081d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 709bd481603SBarry Smith 710bd481603SBarry Smith Use MatPreallocateInitializeSymmetric() for symmetric matrices (MPISBAIJ matrices) 711bd481603SBarry Smith 7121620fd73SBarry Smith This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize(). 7131620fd73SBarry Smith 714bd481603SBarry Smith Concepts: preallocation^Matrix 715bd481603SBarry Smith 716bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(), 717bd481603SBarry Smith MatPreallocateInitializeSymmetric(), MatPreallocateSymmetricSetLocal() 718bd481603SBarry Smith M*/ 719c4f061fbSSatish Balay #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \ 7207c922b88SBarry Smith { \ 721a89bc211SBarry Smith PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__start,__end; \ 722a89bc211SBarry Smith _4_ierr = PetscMalloc2(__nrows,PetscInt,&dnz,__nrows,PetscInt,&onz);CHKERRQ(_4_ierr); \ 723a89bc211SBarry Smith _4_ierr = PetscMemzero(dnz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\ 724a89bc211SBarry Smith _4_ierr = PetscMemzero(onz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\ 725c1ac3661SBarry Smith _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ctmp;\ 726a89bc211SBarry Smith _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows; 7277c922b88SBarry Smith 728bd481603SBarry Smith /*MC 7290d2a7ff2SSatish Balay MatPreallocateSymmetricInitialize - Begins the block of code that will count the number of nonzeros per 730bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 731bd481603SBarry Smith 732bd481603SBarry Smith Synopsis: 733c1ac3661SBarry Smith PetscErrorCode MatPreallocateSymmetricInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz) 734bd481603SBarry Smith 735bd481603SBarry Smith Collective on MPI_Comm 736bd481603SBarry Smith 737bd481603SBarry Smith Input Parameters: 738bd481603SBarry Smith + comm - the communicator that will share the eventually allocated matrix 739859fcb39SBarry Smith . nrows - the number of LOCAL rows in the matrix 740859fcb39SBarry Smith - ncols - the number of LOCAL columns in the matrix 741bd481603SBarry Smith 742bd481603SBarry Smith Output Parameters: 743bd481603SBarry Smith + dnz - the array that will be passed to the matrix preallocation routines 744bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 745bd481603SBarry Smith 746bd481603SBarry Smith 747bd481603SBarry Smith Level: intermediate 748bd481603SBarry Smith 749bd481603SBarry Smith Notes: 7500598bfebSBarry Smith See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details. 751bd481603SBarry Smith 7521d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 753bd481603SBarry Smith 7541620fd73SBarry Smith This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize(). 7551620fd73SBarry Smith 756bd481603SBarry Smith Concepts: preallocation^Matrix 757bd481603SBarry Smith 758bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(), 759bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal() 760bd481603SBarry Smith M*/ 761222b16d4SBarry Smith #define MatPreallocateSymmetricInitialize(comm,nrows,ncols,dnz,onz) 0; \ 762222b16d4SBarry Smith { \ 763a89bc211SBarry Smith PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__end; \ 764a89bc211SBarry Smith _4_ierr = PetscMalloc2(__nrows,PetscInt,&dnz,__nrows,PetscInt,&onz);CHKERRQ(_4_ierr); \ 765a89bc211SBarry Smith _4_ierr = PetscMemzero(dnz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\ 766a89bc211SBarry Smith _4_ierr = PetscMemzero(onz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\ 767c1ac3661SBarry Smith _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr);\ 768a89bc211SBarry Smith _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows; 769222b16d4SBarry Smith 770bd481603SBarry Smith /*MC 771bd481603SBarry Smith MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be 772bd481603SBarry Smith inserted using a local number of the rows and columns 773bd481603SBarry Smith 774bd481603SBarry Smith Synopsis: 775c1ac3661SBarry Smith PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 776bd481603SBarry Smith 777bd481603SBarry Smith Not Collective 778bd481603SBarry Smith 779bd481603SBarry Smith Input Parameters: 780784ac674SJed Brown + map - the row mapping from local numbering to global numbering 781bd481603SBarry Smith . nrows - the number of rows indicated 7821d73ed98SBarry Smith . rows - the indices of the rows 783784ac674SJed Brown . cmap - the column mapping from local to global numbering 784bd481603SBarry Smith . ncols - the number of columns in the matrix 785bd481603SBarry Smith . cols - the columns indicated 786bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 787bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 788bd481603SBarry Smith 789bd481603SBarry Smith 790bd481603SBarry Smith Level: intermediate 791bd481603SBarry Smith 792bd481603SBarry Smith Notes: 7930598bfebSBarry Smith See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details. 794bd481603SBarry Smith 7951d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 796bd481603SBarry Smith 797bd481603SBarry Smith Concepts: preallocation^Matrix 798bd481603SBarry Smith 799bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(), 800bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal() 801bd481603SBarry Smith M*/ 802784ac674SJed Brown #define MatPreallocateSetLocal(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \ 803c4f061fbSSatish Balay {\ 804c1ac3661SBarry Smith PetscInt __l;\ 805784ac674SJed Brown _4_ierr = ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\ 806784ac674SJed Brown _4_ierr = ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\ 807c4f061fbSSatish Balay for (__l=0;__l<nrows;__l++) {\ 808ef66eb69SBarry Smith _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 809c4f061fbSSatish Balay }\ 810c4f061fbSSatish Balay } 811c4f061fbSSatish Balay 812bd481603SBarry Smith /*MC 813bd481603SBarry Smith MatPreallocateSymmetricSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be 814bd481603SBarry Smith inserted using a local number of the rows and columns 815bd481603SBarry Smith 816bd481603SBarry Smith Synopsis: 817c1ac3661SBarry Smith PetscErrorCode MatPreallocateSymmetricSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 818bd481603SBarry Smith 819bd481603SBarry Smith Not Collective 820bd481603SBarry Smith 821bd481603SBarry Smith Input Parameters: 822bd481603SBarry Smith + map - the mapping between local numbering and global numbering 823bd481603SBarry Smith . nrows - the number of rows indicated 8241d73ed98SBarry Smith . rows - the indices of the rows 825bd481603SBarry Smith . ncols - the number of columns in the matrix 826bd481603SBarry Smith . cols - the columns indicated 827bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 828bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 829bd481603SBarry Smith 830bd481603SBarry Smith 831bd481603SBarry Smith Level: intermediate 832bd481603SBarry Smith 833bd481603SBarry Smith Notes: 8340598bfebSBarry Smith See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details. 835bd481603SBarry Smith 836bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 837bd481603SBarry Smith 838bd481603SBarry Smith Concepts: preallocation^Matrix 839bd481603SBarry Smith 840bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(), 841bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal() 842bd481603SBarry Smith M*/ 843d3d32019SBarry Smith #define MatPreallocateSymmetricSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\ 844d3d32019SBarry Smith {\ 845c1ac3661SBarry Smith PetscInt __l;\ 846d3d32019SBarry Smith _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\ 847d3d32019SBarry Smith _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\ 848d3d32019SBarry Smith for (__l=0;__l<nrows;__l++) {\ 849d3d32019SBarry Smith _4_ierr = MatPreallocateSymmetricSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 850d3d32019SBarry Smith }\ 851d3d32019SBarry Smith } 852d3d32019SBarry Smith 853bd481603SBarry Smith /*MC 854bd481603SBarry Smith MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be 855bd481603SBarry Smith inserted using a local number of the rows and columns 856bd481603SBarry Smith 857bd481603SBarry Smith Synopsis: 858c1ac3661SBarry Smith PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 859bd481603SBarry Smith 860bd481603SBarry Smith Not Collective 861bd481603SBarry Smith 862bd481603SBarry Smith Input Parameters: 86364075487SBarry Smith + row - the row 864bd481603SBarry Smith . ncols - the number of columns in the matrix 865a50f8bf6SHong Zhang - cols - the columns indicated 866a50f8bf6SHong Zhang 867a50f8bf6SHong Zhang Output Parameters: 868a50f8bf6SHong Zhang + dnz - the array that will be passed to the matrix preallocation routines 869bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 870bd481603SBarry Smith 871bd481603SBarry Smith 872bd481603SBarry Smith Level: intermediate 873bd481603SBarry Smith 874bd481603SBarry Smith Notes: 8750598bfebSBarry Smith See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details. 876bd481603SBarry Smith 877bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 878bd481603SBarry Smith 8791620fd73SBarry Smith This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize(). 8801620fd73SBarry Smith 881bd481603SBarry Smith Concepts: preallocation^Matrix 882bd481603SBarry Smith 883bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(), 884bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal() 885bd481603SBarry Smith M*/ 886c4f061fbSSatish Balay #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\ 887c1ac3661SBarry Smith { PetscInt __i; \ 888e32f2f54SBarry Smith if (row < __rstart) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Trying to set preallocation for row %D less than first local row %D",row,__rstart);\ 889e32f2f54SBarry Smith if (row >= __rstart+__nrows) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Trying to set preallocation for row %D greater than last local row %D",row,__rstart+__nrows-1);\ 8907c922b88SBarry Smith for (__i=0; __i<nc; __i++) {\ 89164075487SBarry Smith if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \ 8927cc688d7SBarry Smith else dnz[row - __rstart]++;\ 8937c922b88SBarry Smith }\ 8947c922b88SBarry Smith } 8957c922b88SBarry Smith 896bd481603SBarry Smith /*MC 897bd481603SBarry Smith MatPreallocateSymmetricSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be 898bd481603SBarry Smith inserted using a local number of the rows and columns 899bd481603SBarry Smith 900bd481603SBarry Smith Synopsis: 901c1ac3661SBarry Smith PetscErrorCode MatPreallocateSymmetricSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 902bd481603SBarry Smith 903bd481603SBarry Smith Not Collective 904bd481603SBarry Smith 905bd481603SBarry Smith Input Parameters: 906bd481603SBarry Smith + nrows - the number of rows indicated 9071d73ed98SBarry Smith . rows - the indices of the rows 908bd481603SBarry Smith . ncols - the number of columns in the matrix 909bd481603SBarry Smith . cols - the columns indicated 910bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 911bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 912bd481603SBarry Smith 913bd481603SBarry Smith 914bd481603SBarry Smith Level: intermediate 915bd481603SBarry Smith 916bd481603SBarry Smith Notes: 9170598bfebSBarry Smith See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details. 918bd481603SBarry Smith 919bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 920bd481603SBarry Smith 9211620fd73SBarry Smith This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize(). 9221620fd73SBarry Smith 923bd481603SBarry Smith Concepts: preallocation^Matrix 924bd481603SBarry Smith 925bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(), 926bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal() 927bd481603SBarry Smith M*/ 928d3d32019SBarry Smith #define MatPreallocateSymmetricSet(row,nc,cols,dnz,onz) 0;\ 929c1ac3661SBarry Smith { PetscInt __i; \ 930d3d32019SBarry Smith for (__i=0; __i<nc; __i++) {\ 931d3d32019SBarry Smith if (cols[__i] >= __end) onz[row - __rstart]++; \ 932d3d32019SBarry Smith else if (cols[__i] >= row) dnz[row - __rstart]++;\ 933d3d32019SBarry Smith }\ 934d3d32019SBarry Smith } 935d3d32019SBarry Smith 936bd481603SBarry Smith /*MC 93716371a99SBarry Smith MatPreallocateLocation - An alternative to MatPreallocationSet() that puts the nonzero locations into the matrix if it exists 93816371a99SBarry Smith 93916371a99SBarry Smith Synopsis: 94016371a99SBarry Smith PetscErrorCode MatPreallocateLocations(Mat A,PetscInt row,PetscInt ncols,PetscInt *cols,PetscInt *dnz,PetscInt *onz) 94116371a99SBarry Smith 94216371a99SBarry Smith Not Collective 94316371a99SBarry Smith 94416371a99SBarry Smith Input Parameters: 94516371a99SBarry Smith . A - matrix 94616371a99SBarry Smith . row - row where values exist (must be local to this process) 94716371a99SBarry Smith . ncols - number of columns 94816371a99SBarry Smith . cols - columns with nonzeros 94916371a99SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 95016371a99SBarry Smith - ozn - the other array passed to the matrix preallocation routines 95116371a99SBarry Smith 95216371a99SBarry Smith 95316371a99SBarry Smith Level: intermediate 95416371a99SBarry Smith 95516371a99SBarry Smith Notes: 9560598bfebSBarry Smith See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details. 95716371a99SBarry Smith 95816371a99SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 95916371a99SBarry Smith 96016371a99SBarry Smith This is a MACRO not a function because it uses a bunch of variables private to the MatPreallocation.... routines. 96116371a99SBarry Smith 96216371a99SBarry Smith Concepts: preallocation^Matrix 96316371a99SBarry Smith 96416371a99SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(), 96516371a99SBarry Smith MatPreallocateSymmetricInitialize(), MatPreallocateSymmetricSetLocal() 96616371a99SBarry Smith M*/ 96716371a99SBarry 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);} 96816371a99SBarry Smith 96916371a99SBarry Smith 97016371a99SBarry Smith /*MC 9710d2a7ff2SSatish Balay MatPreallocateFinalize - Ends the block of code that will count the number of nonzeros per 972bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 973bd481603SBarry Smith 974bd481603SBarry Smith Synopsis: 975c1ac3661SBarry Smith PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz) 976bd481603SBarry Smith 977bd481603SBarry Smith Collective on MPI_Comm 978bd481603SBarry Smith 979bd481603SBarry Smith Input Parameters: 98016371a99SBarry Smith + dnz - the array that was be passed to the matrix preallocation routines 981bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 982bd481603SBarry Smith 983bd481603SBarry Smith 984bd481603SBarry Smith Level: intermediate 985bd481603SBarry Smith 986bd481603SBarry Smith Notes: 9870598bfebSBarry Smith See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details. 988bd481603SBarry Smith 989bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 990bd481603SBarry Smith 9911620fd73SBarry Smith This is a MACRO not a function because it closes the { started in MatPreallocateInitialize(). 9921620fd73SBarry Smith 993bd481603SBarry Smith Concepts: preallocation^Matrix 994bd481603SBarry Smith 995bd481603SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(), 996bd481603SBarry Smith MatPreallocateSymmetricInitialize(), MatPreallocateSymmetricSetLocal() 997bd481603SBarry Smith M*/ 998a89bc211SBarry Smith #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree2(dnz,onz);CHKERRQ(_4_ierr);} 9997c922b88SBarry Smith 1000bd481603SBarry Smith 1001bd481603SBarry Smith 10027b80b807SBarry Smith /* Routines unique to particular data structures */ 10037087cfbeSBarry Smith extern PetscErrorCode MatShellGetContext(Mat,void **); 1004ba966639SSatish Balay PetscPolymorphicFunction(MatShellGetContext,(Mat A),(A,&t),void*,t) 1005698d4c6aSKris Buschelman 10067087cfbeSBarry Smith extern PetscErrorCode MatInodeAdjustForInodes(Mat,IS*,IS*); 10077087cfbeSBarry Smith extern PetscErrorCode MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *); 1008698d4c6aSKris Buschelman 10097087cfbeSBarry Smith extern PetscErrorCode MatSeqAIJSetColumnIndices(Mat,PetscInt[]); 10107087cfbeSBarry Smith extern PetscErrorCode MatSeqBAIJSetColumnIndices(Mat,PetscInt[]); 10117087cfbeSBarry Smith extern PetscErrorCode MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 10127087cfbeSBarry Smith extern PetscErrorCode MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 10137087cfbeSBarry Smith extern PetscErrorCode MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 10147b80b807SBarry Smith 1015a96a251dSBarry Smith #define MAT_SKIP_ALLOCATION -4 1016a96a251dSBarry Smith 10177087cfbeSBarry Smith extern PetscErrorCode MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]); 1018ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz)) 10197087cfbeSBarry Smith extern PetscErrorCode MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]); 1020ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqSBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz)) 10217087cfbeSBarry Smith extern PetscErrorCode MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]); 1022ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqAIJSetPreallocation,(Mat A,const PetscInt nnz[]),(A,0,nnz)) 1023273d9f13SBarry Smith 10247087cfbeSBarry Smith extern PetscErrorCode MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 1025ba966639SSatish Balay PetscPolymorphicSubroutine(MatMPIBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[],const PetscInt onz[]),(A,bs,0,nnz,0,onz)) 10267087cfbeSBarry Smith extern PetscErrorCode MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 10277087cfbeSBarry Smith extern PetscErrorCode MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 10287087cfbeSBarry Smith extern PetscErrorCode MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []); 10297087cfbeSBarry Smith extern PetscErrorCode MatSeqBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 10307087cfbeSBarry Smith extern PetscErrorCode MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]); 10317087cfbeSBarry Smith extern PetscErrorCode MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 10327087cfbeSBarry Smith extern PetscErrorCode MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]); 10337087cfbeSBarry Smith extern PetscErrorCode MatMPIDenseSetPreallocation(Mat,PetscScalar[]); 10347087cfbeSBarry Smith extern PetscErrorCode MatSeqDenseSetPreallocation(Mat,PetscScalar[]); 10357087cfbeSBarry Smith extern PetscErrorCode MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,PetscInt*[]); 10367087cfbeSBarry Smith extern PetscErrorCode MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,PetscInt*[]); 10377087cfbeSBarry Smith extern PetscErrorCode MatAdicSetLocalFunction(Mat,void (*)(void)); 1038273d9f13SBarry Smith 10397087cfbeSBarry Smith extern PetscErrorCode MatSeqDenseSetLDA(Mat,PetscInt); 10407087cfbeSBarry Smith extern PetscErrorCode MatDenseGetLocalMatrix(Mat,Mat*); 10411b807ce4Svictorle 10427087cfbeSBarry Smith extern PetscErrorCode MatStoreValues(Mat); 10437087cfbeSBarry Smith extern PetscErrorCode MatRetrieveValues(Mat); 10442e8a6d31SBarry Smith 10457087cfbeSBarry Smith extern PetscErrorCode MatDAADSetCtx(Mat,void*); 10463a7fca6bSBarry Smith 1047b3a44c85SBarry Smith extern PetscErrorCode MatFindNonzeroRows(Mat,IS*); 10487b80b807SBarry Smith /* 10497b80b807SBarry Smith These routines are not usually accessed directly, rather solving is 105094b7f48cSBarry Smith done through the KSP and PC interfaces. 10517b80b807SBarry Smith */ 10527b80b807SBarry Smith 1053d9274352SBarry Smith /*E 1054d9274352SBarry Smith MatOrderingType - String with the name of a PETSc matrix ordering or the creation function 1055d9274352SBarry Smith with an optional dynamic library name, for example 1056d9274352SBarry Smith http://www.mcs.anl.gov/petsc/lib.a:orderingcreate() 1057d9274352SBarry Smith 1058d9274352SBarry Smith Level: beginner 1059d9274352SBarry Smith 1060e5a9bf91SBarry Smith Cannot use const because the PC objects manipulate the string 1061e5a9bf91SBarry Smith 1062d9274352SBarry Smith .seealso: MatGetOrdering() 1063d9274352SBarry Smith E*/ 10643eaad2caSSatish Balay #define MatOrderingType char* 10652692d6eeSBarry Smith #define MATORDERINGNATURAL "natural" 10662692d6eeSBarry Smith #define MATORDERINGND "nd" 10672692d6eeSBarry Smith #define MATORDERING1WD "1wd" 10682692d6eeSBarry Smith #define MATORDERINGRCM "rcm" 10692692d6eeSBarry Smith #define MATORDERINGQMD "qmd" 10702692d6eeSBarry Smith #define MATORDERINGROWLENGTH "rowlength" 10712692d6eeSBarry Smith #define MATORDERINGDSC_ND "dsc_nd" /* these three are only for DSCPACK, see its documentation for details */ 10722692d6eeSBarry Smith #define MATORDERINGDSC_MMD "dsc_mmd" 10732692d6eeSBarry Smith #define MATORDERINGDSC_MDF "dsc_mdf" 1074312e372aSBarry Smith #define MATORDERINGAMD "amd" /* only works if UMFPACK is installed with PETSc */ 1075b12f92e5SBarry Smith 10767087cfbeSBarry Smith extern PetscErrorCode MatGetOrdering(Mat,const MatOrderingType,IS*,IS*); 10777087cfbeSBarry Smith extern PetscErrorCode MatGetOrderingList(PetscFList *list); 10787087cfbeSBarry Smith extern PetscErrorCode MatOrderingRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,const MatOrderingType,IS*,IS*)); 107930de9b25SBarry Smith 108030de9b25SBarry Smith /*MC 10811890ba74SBarry Smith MatOrderingRegisterDynamic - Adds a new sparse matrix ordering to the matrix package. 108230de9b25SBarry Smith 108330de9b25SBarry Smith Synopsis: 10841890ba74SBarry Smith PetscErrorCode MatOrderingRegisterDynamic(const char *name_ordering,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatOrdering)) 108530de9b25SBarry Smith 108630de9b25SBarry Smith Not Collective 108730de9b25SBarry Smith 108830de9b25SBarry Smith Input Parameters: 10892692d6eeSBarry Smith + sname - name of ordering (for example MATORDERINGND) 109030de9b25SBarry Smith . path - location of library where creation routine is 109130de9b25SBarry Smith . name - name of function that creates the ordering type,a string 109230de9b25SBarry Smith - function - function pointer that creates the ordering 109330de9b25SBarry Smith 109430de9b25SBarry Smith Level: developer 109530de9b25SBarry Smith 109630de9b25SBarry Smith If dynamic libraries are used, then the fourth input argument (function) 109730de9b25SBarry Smith is ignored. 109830de9b25SBarry Smith 109930de9b25SBarry Smith Sample usage: 110030de9b25SBarry Smith .vb 110130de9b25SBarry Smith MatOrderingRegisterDynamic("my_order",/home/username/my_lib/lib/libO/solaris/mylib.a, 110230de9b25SBarry Smith "MyOrder",MyOrder); 110330de9b25SBarry Smith .ve 110430de9b25SBarry Smith 110530de9b25SBarry Smith Then, your partitioner can be chosen with the procedural interface via 110630de9b25SBarry Smith $ MatOrderingSetType(part,"my_order) 110730de9b25SBarry Smith or at runtime via the option 11082401956bSBarry Smith $ -pc_factor_mat_ordering_type my_order 110930de9b25SBarry Smith 1110ab901514SBarry Smith ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values. 111130de9b25SBarry Smith 111230de9b25SBarry Smith .keywords: matrix, ordering, register 111330de9b25SBarry Smith 111430de9b25SBarry Smith .seealso: MatOrderingRegisterDestroy(), MatOrderingRegisterAll() 111530de9b25SBarry Smith M*/ 1116aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 1117f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,0) 1118b12f92e5SBarry Smith #else 1119f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,d) 1120b12f92e5SBarry Smith #endif 112130de9b25SBarry Smith 11227087cfbeSBarry Smith extern PetscErrorCode MatOrderingRegisterDestroy(void); 11237087cfbeSBarry Smith extern PetscErrorCode MatOrderingRegisterAll(const char[]); 1124ace3abfcSBarry Smith extern PetscBool MatOrderingRegisterAllCalled; 1125b0a32e0cSBarry Smith extern PetscFList MatOrderingList; 1126d4fbbf0eSBarry Smith 11277087cfbeSBarry Smith extern PetscErrorCode MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS); 1128a2ce50c7SBarry Smith 1129d91e6319SBarry Smith /*S 1130d90ac83dSHong Zhang MatFactorShiftType - Numeric Shift. 1131d90ac83dSHong Zhang 1132d90ac83dSHong Zhang Level: beginner 1133d90ac83dSHong Zhang 1134d90ac83dSHong Zhang S*/ 1135d90ac83dSHong Zhang typedef enum {MAT_SHIFT_NONE,MAT_SHIFT_NONZERO,MAT_SHIFT_POSITIVE_DEFINITE,MAT_SHIFT_INBLOCKS} MatFactorShiftType; 1136d90ac83dSHong Zhang extern const char *MatFactorShiftTypes[]; 1137d90ac83dSHong Zhang 1138d90ac83dSHong Zhang /*S 11392401956bSBarry Smith MatFactorInfo - Data passed into the matrix factorization routines 1140b00f7748SHong Zhang 114161cad860SBarry Smith In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE, that is use 114261cad860SBarry Smith $ MatFactorInfo info(MAT_FACTORINFO_SIZE) 1143b00f7748SHong Zhang 114415e8a5b3SHong Zhang Notes: These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC. 1145b00f7748SHong Zhang 114661cad860SBarry Smith You can use MatFactorInfoInitialize() to set default values. 114761cad860SBarry Smith 1148b00f7748SHong Zhang Level: developer 1149b00f7748SHong Zhang 1150d7d82daaSBarry Smith .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(), 1151d7d82daaSBarry Smith MatFactorInfoInitialize() 1152b00f7748SHong Zhang 1153b00f7748SHong Zhang S*/ 1154b00f7748SHong Zhang typedef struct { 115515e8a5b3SHong Zhang PetscReal diagonal_fill; /* force diagonal to fill in if initially not filled */ 115685317021SBarry Smith PetscReal usedt; 115715e8a5b3SHong Zhang PetscReal dt; /* drop tolerance */ 1158b00f7748SHong Zhang PetscReal dtcol; /* tolerance for pivoting */ 115915e8a5b3SHong Zhang PetscReal dtcount; /* maximum nonzeros to be allowed per row */ 116067e28bfeSBarry Smith PetscReal fill; /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */ 1161348344bbSBarry Smith PetscReal levels; /* ICC/ILU(levels) */ 1162bcd9e38bSBarry Smith PetscReal pivotinblocks; /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0 1163bcd9e38bSBarry Smith factorization may be faster if do not pivot */ 116415e8a5b3SHong Zhang PetscReal zeropivot; /* pivot is called zero if less than this */ 1165f4db908eSBarry Smith PetscReal shifttype; /* type of shift added to matrix factor to prevent zero pivots */ 1166f4db908eSBarry Smith PetscReal shiftamount; /* how large the shift is */ 116715e8a5b3SHong Zhang } MatFactorInfo; 1168ffa6d0a5SLois Curfman McInnes 11697087cfbeSBarry Smith extern PetscErrorCode MatFactorInfoInitialize(MatFactorInfo*); 11707087cfbeSBarry Smith extern PetscErrorCode MatCholeskyFactor(Mat,IS,const MatFactorInfo*); 11717087cfbeSBarry Smith extern PetscErrorCode MatCholeskyFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*); 11727087cfbeSBarry Smith extern PetscErrorCode MatCholeskyFactorNumeric(Mat,Mat,const MatFactorInfo*); 11737087cfbeSBarry Smith extern PetscErrorCode MatLUFactor(Mat,IS,IS,const MatFactorInfo*); 11747087cfbeSBarry Smith extern PetscErrorCode MatILUFactor(Mat,IS,IS,const MatFactorInfo*); 11757087cfbeSBarry Smith extern PetscErrorCode MatLUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*); 11767087cfbeSBarry Smith extern PetscErrorCode MatILUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*); 11777087cfbeSBarry Smith extern PetscErrorCode MatICCFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*); 11787087cfbeSBarry Smith extern PetscErrorCode MatICCFactor(Mat,IS,const MatFactorInfo*); 11797087cfbeSBarry Smith extern PetscErrorCode MatLUFactorNumeric(Mat,Mat,const MatFactorInfo*); 11807087cfbeSBarry Smith extern PetscErrorCode MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*); 11817087cfbeSBarry Smith extern PetscErrorCode MatSolve(Mat,Vec,Vec); 11827087cfbeSBarry Smith extern PetscErrorCode MatForwardSolve(Mat,Vec,Vec); 11837087cfbeSBarry Smith extern PetscErrorCode MatBackwardSolve(Mat,Vec,Vec); 11847087cfbeSBarry Smith extern PetscErrorCode MatSolveAdd(Mat,Vec,Vec,Vec); 11857087cfbeSBarry Smith extern PetscErrorCode MatSolveTranspose(Mat,Vec,Vec); 11867087cfbeSBarry Smith extern PetscErrorCode MatSolveTransposeAdd(Mat,Vec,Vec,Vec); 11877087cfbeSBarry Smith extern PetscErrorCode MatSolves(Mat,Vecs,Vecs); 11888ed539a5SBarry Smith 11897087cfbeSBarry Smith extern PetscErrorCode MatSetUnfactored(Mat); 1190bb5a7306SBarry Smith 1191d91e6319SBarry Smith /*E 1192d91e6319SBarry Smith MatSORType - What type of (S)SOR to perform 1193bb1eb677SSatish Balay 1194d91e6319SBarry Smith Level: beginner 1195d91e6319SBarry Smith 1196d9274352SBarry Smith May be bitwise ORd together 1197d9274352SBarry Smith 1198d91e6319SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 1199d91e6319SBarry Smith 12004e7234bfSBarry Smith MatSORType may be bitwise ORd together, so do not change the numbers 12014e7234bfSBarry Smith 120241f059aeSBarry Smith .seealso: MatSOR() 1203d91e6319SBarry Smith E*/ 1204ee50ffe9SBarry Smith typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3, 1205ee50ffe9SBarry Smith SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8, 1206ee50ffe9SBarry Smith SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16, 120784cb2905SBarry Smith SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType; 12087087cfbeSBarry Smith extern PetscErrorCode MatSOR(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec); 12098ed539a5SBarry Smith 1210d4fbbf0eSBarry Smith /* 1211639f9d9dSBarry Smith These routines are for efficiently computing Jacobians via finite differences. 1212639f9d9dSBarry Smith */ 1213b12f92e5SBarry Smith 1214d9274352SBarry Smith /*E 1215d9274352SBarry Smith MatColoringType - String with the name of a PETSc matrix coloring or the creation function 1216d9274352SBarry Smith with an optional dynamic library name, for example 1217d9274352SBarry Smith http://www.mcs.anl.gov/petsc/lib.a:coloringcreate() 1218d9274352SBarry Smith 1219d9274352SBarry Smith Level: beginner 1220d9274352SBarry Smith 1221d9274352SBarry Smith .seealso: MatGetColoring() 1222d9274352SBarry Smith E*/ 1223a313700dSBarry Smith #define MatColoringType char* 12242692d6eeSBarry Smith #define MATCOLORINGNATURAL "natural" 12252692d6eeSBarry Smith #define MATCOLORINGSL "sl" 12262692d6eeSBarry Smith #define MATCOLORINGLF "lf" 12272692d6eeSBarry Smith #define MATCOLORINGID "id" 1228b12f92e5SBarry Smith 12297087cfbeSBarry Smith extern PetscErrorCode MatGetColoring(Mat,const MatColoringType,ISColoring*); 12307087cfbeSBarry Smith extern PetscErrorCode MatColoringRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,MatColoringType,ISColoring *)); 123130de9b25SBarry Smith 123230de9b25SBarry Smith /*MC 123330de9b25SBarry Smith MatColoringRegisterDynamic - Adds a new sparse matrix coloring to the 123430de9b25SBarry Smith matrix package. 123530de9b25SBarry Smith 123630de9b25SBarry Smith Synopsis: 12371890ba74SBarry Smith PetscErrorCode MatColoringRegisterDynamic(const char *name_coloring,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatColoring)) 123830de9b25SBarry Smith 123930de9b25SBarry Smith Not Collective 124030de9b25SBarry Smith 124130de9b25SBarry Smith Input Parameters: 12422692d6eeSBarry Smith + sname - name of Coloring (for example MATCOLORINGSL) 124330de9b25SBarry Smith . path - location of library where creation routine is 124430de9b25SBarry Smith . name - name of function that creates the Coloring type, a string 124530de9b25SBarry Smith - function - function pointer that creates the coloring 124630de9b25SBarry Smith 124730de9b25SBarry Smith Level: developer 124830de9b25SBarry Smith 124930de9b25SBarry Smith If dynamic libraries are used, then the fourth input argument (function) 125030de9b25SBarry Smith is ignored. 125130de9b25SBarry Smith 125230de9b25SBarry Smith Sample usage: 125330de9b25SBarry Smith .vb 125430de9b25SBarry Smith MatColoringRegisterDynamic("my_color",/home/username/my_lib/lib/libO/solaris/mylib.a, 125530de9b25SBarry Smith "MyColor",MyColor); 125630de9b25SBarry Smith .ve 125730de9b25SBarry Smith 125830de9b25SBarry Smith Then, your partitioner can be chosen with the procedural interface via 125930de9b25SBarry Smith $ MatColoringSetType(part,"my_color") 126030de9b25SBarry Smith or at runtime via the option 126130de9b25SBarry Smith $ -mat_coloring_type my_color 126230de9b25SBarry Smith 1263ab901514SBarry Smith $PETSC_ARCH occuring in pathname will be replaced with appropriate values. 126430de9b25SBarry Smith 126530de9b25SBarry Smith .keywords: matrix, Coloring, register 126630de9b25SBarry Smith 126730de9b25SBarry Smith .seealso: MatColoringRegisterDestroy(), MatColoringRegisterAll() 126830de9b25SBarry Smith M*/ 1269aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 1270f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,0) 1271b12f92e5SBarry Smith #else 1272f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,d) 1273b12f92e5SBarry Smith #endif 127430de9b25SBarry Smith 1275ace3abfcSBarry Smith extern PetscBool MatColoringRegisterAllCalled; 1276f1144a30SSatish Balay 12777087cfbeSBarry Smith extern PetscErrorCode MatColoringRegisterAll(const char[]); 12787087cfbeSBarry Smith extern PetscErrorCode MatColoringRegisterDestroy(void); 12797087cfbeSBarry Smith extern PetscErrorCode MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*); 1280639f9d9dSBarry Smith 1281d9274352SBarry Smith /*S 1282d9274352SBarry Smith MatFDColoring - Object for computing a sparse Jacobian via finite differences 1283d9274352SBarry Smith and coloring 1284639f9d9dSBarry Smith 1285d9274352SBarry Smith Level: beginner 1286d9274352SBarry Smith 1287d9274352SBarry Smith Concepts: coloring, sparse Jacobian, finite differences 1288d9274352SBarry Smith 1289d9274352SBarry Smith .seealso: MatFDColoringCreate() 1290d9274352SBarry Smith S*/ 1291e2a1c21fSSatish Balay typedef struct _p_MatFDColoring* MatFDColoring; 1292639f9d9dSBarry Smith 12937087cfbeSBarry Smith extern PetscErrorCode MatFDColoringCreate(Mat,ISColoring,MatFDColoring *); 1294fcfd50ebSBarry Smith extern PetscErrorCode MatFDColoringDestroy(MatFDColoring*); 12957087cfbeSBarry Smith extern PetscErrorCode MatFDColoringView(MatFDColoring,PetscViewer); 12967087cfbeSBarry Smith extern PetscErrorCode MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*); 12977087cfbeSBarry Smith extern PetscErrorCode MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**); 12987087cfbeSBarry Smith extern PetscErrorCode MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal); 12997087cfbeSBarry Smith extern PetscErrorCode MatFDColoringSetFromOptions(MatFDColoring); 13007087cfbeSBarry Smith extern PetscErrorCode MatFDColoringApply(Mat,MatFDColoring,Vec,MatStructure*,void *); 13017087cfbeSBarry Smith extern PetscErrorCode MatFDColoringApplyTS(Mat,MatFDColoring,PetscReal,Vec,MatStructure*,void *); 13027087cfbeSBarry Smith extern PetscErrorCode MatFDColoringSetF(MatFDColoring,Vec); 13037087cfbeSBarry Smith extern PetscErrorCode MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,PetscInt*[]); 1304639f9d9dSBarry Smith /* 13050752156aSBarry Smith These routines are for partitioning matrices: currently used only 13063eda8832SBarry Smith for adjacency matrix, MatCreateMPIAdj(). 13070752156aSBarry Smith */ 1308ca161407SBarry Smith 1309d9274352SBarry Smith /*S 1310d9274352SBarry Smith MatPartitioning - Object for managing the partitioning of a matrix or graph 1311d9274352SBarry Smith 1312d9274352SBarry Smith Level: beginner 1313d9274352SBarry Smith 1314d9274352SBarry Smith Concepts: partitioning 1315d9274352SBarry Smith 1316743a1026Svictorle .seealso: MatPartitioningCreate(), MatPartitioningType 1317d9274352SBarry Smith S*/ 131891e9ee9fSBarry Smith typedef struct _p_MatPartitioning* MatPartitioning; 1319d9274352SBarry Smith 1320d9274352SBarry Smith /*E 13215bc6e7ccSvictorle MatPartitioningType - String with the name of a PETSc matrix partitioning or the creation function 1322d9274352SBarry Smith with an optional dynamic library name, for example 1323d9274352SBarry Smith http://www.mcs.anl.gov/petsc/lib.a:partitioningcreate() 1324d9274352SBarry Smith 1325d9274352SBarry Smith Level: beginner 1326d9274352SBarry Smith 1327b547a13bSvictorle .seealso: MatPartitioningCreate(), MatPartitioning 1328d9274352SBarry Smith E*/ 1329a313700dSBarry Smith #define MatPartitioningType char* 13302692d6eeSBarry Smith #define MATPARTITIONINGCURRENT "current" 13312692d6eeSBarry Smith #define MATPARTITIONINGSQUARE "square" 13322692d6eeSBarry Smith #define MATPARTITIONINGPARMETIS "parmetis" 13332692d6eeSBarry Smith #define MATPARTITIONINGCHACO "chaco" 13342692d6eeSBarry Smith #define MATPARTITIONINGJOSTLE "jostle" 13352692d6eeSBarry Smith #define MATPARTITIONINGPARTY "party" 13362692d6eeSBarry Smith #define MATPARTITIONINGSCOTCH "scotch" 133771306c60Sjroman 1338ca161407SBarry Smith 13397087cfbeSBarry Smith extern PetscErrorCode MatPartitioningCreate(MPI_Comm,MatPartitioning*); 13407087cfbeSBarry Smith extern PetscErrorCode MatPartitioningSetType(MatPartitioning,const MatPartitioningType); 13417087cfbeSBarry Smith extern PetscErrorCode MatPartitioningSetNParts(MatPartitioning,PetscInt); 13427087cfbeSBarry Smith extern PetscErrorCode MatPartitioningSetAdjacency(MatPartitioning,Mat); 13437087cfbeSBarry Smith extern PetscErrorCode MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]); 13447087cfbeSBarry Smith extern PetscErrorCode MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []); 13457087cfbeSBarry Smith extern PetscErrorCode MatPartitioningApply(MatPartitioning,IS*); 1346fcfd50ebSBarry Smith extern PetscErrorCode MatPartitioningDestroy(MatPartitioning*); 13472aabb6bbSBarry Smith 13487087cfbeSBarry Smith extern PetscErrorCode MatPartitioningRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatPartitioning)); 134930de9b25SBarry Smith 135030de9b25SBarry Smith /*MC 135130de9b25SBarry Smith MatPartitioningRegisterDynamic - Adds a new sparse matrix partitioning to the 135230de9b25SBarry Smith matrix package. 135330de9b25SBarry Smith 135430de9b25SBarry Smith Synopsis: 13551890ba74SBarry Smith PetscErrorCode MatPartitioningRegisterDynamic(const char *name_partitioning,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatPartitioning)) 135630de9b25SBarry Smith 135730de9b25SBarry Smith Not Collective 135830de9b25SBarry Smith 135930de9b25SBarry Smith Input Parameters: 13602692d6eeSBarry Smith + sname - name of partitioning (for example MATPARTITIONINGCURRENT) or parmetis 136130de9b25SBarry Smith . path - location of library where creation routine is 136230de9b25SBarry Smith . name - name of function that creates the partitioning type, a string 136330de9b25SBarry Smith - function - function pointer that creates the partitioning type 136430de9b25SBarry Smith 136530de9b25SBarry Smith Level: developer 136630de9b25SBarry Smith 136730de9b25SBarry Smith If dynamic libraries are used, then the fourth input argument (function) 136830de9b25SBarry Smith is ignored. 136930de9b25SBarry Smith 137030de9b25SBarry Smith Sample usage: 137130de9b25SBarry Smith .vb 137230de9b25SBarry Smith MatPartitioningRegisterDynamic("my_part",/home/username/my_lib/lib/libO/solaris/mylib.a, 137330de9b25SBarry Smith "MyPartCreate",MyPartCreate); 137430de9b25SBarry Smith .ve 137530de9b25SBarry Smith 137630de9b25SBarry Smith Then, your partitioner can be chosen with the procedural interface via 137730de9b25SBarry Smith $ MatPartitioningSetType(part,"my_part") 137830de9b25SBarry Smith or at runtime via the option 137930de9b25SBarry Smith $ -mat_partitioning_type my_part 138030de9b25SBarry Smith 1381ab901514SBarry Smith $PETSC_ARCH occuring in pathname will be replaced with appropriate values. 138230de9b25SBarry Smith 138330de9b25SBarry Smith .keywords: matrix, partitioning, register 138430de9b25SBarry Smith 138530de9b25SBarry Smith .seealso: MatPartitioningRegisterDestroy(), MatPartitioningRegisterAll() 138630de9b25SBarry Smith M*/ 1387aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 1388f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,0) 13892aabb6bbSBarry Smith #else 1390f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,d) 13912aabb6bbSBarry Smith #endif 13922aabb6bbSBarry Smith 1393ace3abfcSBarry Smith extern PetscBool MatPartitioningRegisterAllCalled; 1394f1144a30SSatish Balay 13957087cfbeSBarry Smith extern PetscErrorCode MatPartitioningRegisterAll(const char[]); 13967087cfbeSBarry Smith extern PetscErrorCode MatPartitioningRegisterDestroy(void); 13972bad1931SBarry Smith 13987087cfbeSBarry Smith extern PetscErrorCode MatPartitioningView(MatPartitioning,PetscViewer); 13997087cfbeSBarry Smith extern PetscErrorCode MatPartitioningSetFromOptions(MatPartitioning); 14007087cfbeSBarry Smith extern PetscErrorCode MatPartitioningGetType(MatPartitioning,const MatPartitioningType*); 1401ca161407SBarry Smith 14027087cfbeSBarry Smith extern PetscErrorCode MatPartitioningParmetisSetCoarseSequential(MatPartitioning); 14037087cfbeSBarry Smith extern PetscErrorCode MatPartitioningParmetisGetEdgeCut(MatPartitioning, PetscInt *); 14040752156aSBarry Smith 14057087cfbeSBarry Smith extern PetscErrorCode MatPartitioningJostleSetCoarseLevel(MatPartitioning,PetscReal); 14067087cfbeSBarry Smith extern PetscErrorCode MatPartitioningJostleSetCoarseSequential(MatPartitioning); 140771306c60Sjroman 1408954b3ec6SBarry Smith typedef enum {MP_CHACO_MULTILEVEL_KL,MP_CHACO_SPECTRAL,MP_CHACO_LINEAR,MP_CHACO_RANDOM, MP_CHACO_SCATTERED} MPChacoGlobalType; 14097087cfbeSBarry Smith extern PetscErrorCode MatPartitioningChacoSetGlobal(MatPartitioning, MPChacoGlobalType); 141071306c60Sjroman typedef enum { MP_CHACO_KERNIGHAN_LIN, MP_CHACO_NONE } MPChacoLocalType; 14117087cfbeSBarry Smith extern PetscErrorCode MatPartitioningChacoSetLocal(MatPartitioning, MPChacoLocalType); 14127087cfbeSBarry Smith extern PetscErrorCode MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal); 141371306c60Sjroman typedef enum { MP_CHACO_LANCZOS, MP_CHACO_RQI_SYMMLQ } MPChacoEigenType; 14147087cfbeSBarry Smith extern PetscErrorCode MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType); 14157087cfbeSBarry Smith extern PetscErrorCode MatPartitioningChacoSetEigenTol(MatPartitioning, PetscReal); 14167087cfbeSBarry Smith extern PetscErrorCode MatPartitioningChacoSetEigenNumber(MatPartitioning, PetscInt); 141771306c60Sjroman 141871306c60Sjroman #define MP_PARTY_OPT "opt" 141971306c60Sjroman #define MP_PARTY_LIN "lin" 142071306c60Sjroman #define MP_PARTY_SCA "sca" 142171306c60Sjroman #define MP_PARTY_RAN "ran" 142271306c60Sjroman #define MP_PARTY_GBF "gbf" 142371306c60Sjroman #define MP_PARTY_GCF "gcf" 142471306c60Sjroman #define MP_PARTY_BUB "bub" 142571306c60Sjroman #define MP_PARTY_DEF "def" 14267087cfbeSBarry Smith extern PetscErrorCode MatPartitioningPartySetGlobal(MatPartitioning, const char*); 142771306c60Sjroman #define MP_PARTY_HELPFUL_SETS "hs" 142871306c60Sjroman #define MP_PARTY_KERNIGHAN_LIN "kl" 142971306c60Sjroman #define MP_PARTY_NONE "no" 14307087cfbeSBarry Smith extern PetscErrorCode MatPartitioningPartySetLocal(MatPartitioning, const char*); 14317087cfbeSBarry Smith extern PetscErrorCode MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal); 14327087cfbeSBarry Smith extern PetscErrorCode MatPartitioningPartySetBipart(MatPartitioning,PetscBool ); 14337087cfbeSBarry Smith extern PetscErrorCode MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscBool ); 143471306c60Sjroman 143571306c60Sjroman typedef enum { MP_SCOTCH_GREEDY, MP_SCOTCH_GPS, MP_SCOTCH_GR_GPS } MPScotchGlobalType; 14367087cfbeSBarry Smith extern PetscErrorCode MatPartitioningScotchSetArch(MatPartitioning,const char*); 14377087cfbeSBarry Smith extern PetscErrorCode MatPartitioningScotchSetMultilevel(MatPartitioning); 14387087cfbeSBarry Smith extern PetscErrorCode MatPartitioningScotchSetGlobal(MatPartitioning,MPScotchGlobalType); 14397087cfbeSBarry Smith extern PetscErrorCode MatPartitioningScotchSetCoarseLevel(MatPartitioning,PetscReal); 14407087cfbeSBarry Smith extern PetscErrorCode MatPartitioningScotchSetHostList(MatPartitioning,const char*); 144171306c60Sjroman typedef enum { MP_SCOTCH_KERNIGHAN_LIN, MP_SCOTCH_NONE } MPScotchLocalType; 14427087cfbeSBarry Smith extern PetscErrorCode MatPartitioningScotchSetLocal(MatPartitioning,MPScotchLocalType); 14437087cfbeSBarry Smith extern PetscErrorCode MatPartitioningScotchSetMapping(MatPartitioning); 14447087cfbeSBarry Smith extern PetscErrorCode MatPartitioningScotchSetStrategy(MatPartitioning,char*); 144571306c60Sjroman 144609573ac7SBarry Smith extern PetscErrorCode MatMeshToVertexGraph(Mat,PetscInt,Mat*); 144709573ac7SBarry Smith extern PetscErrorCode MatMeshToCellGraph(Mat,PetscInt,Mat*); 1448591294e4SBarry Smith 14490752156aSBarry Smith /* 14500a835dfdSSatish Balay If you add entries here you must also add them to finclude/petscmat.h 1451d4fbbf0eSBarry Smith */ 14521c1c02c0SLois Curfman McInnes typedef enum { MATOP_SET_VALUES=0, 14531c1c02c0SLois Curfman McInnes MATOP_GET_ROW=1, 14541c1c02c0SLois Curfman McInnes MATOP_RESTORE_ROW=2, 14551c1c02c0SLois Curfman McInnes MATOP_MULT=3, 14561c1c02c0SLois Curfman McInnes MATOP_MULT_ADD=4, 14577c922b88SBarry Smith MATOP_MULT_TRANSPOSE=5, 14587c922b88SBarry Smith MATOP_MULT_TRANSPOSE_ADD=6, 14591c1c02c0SLois Curfman McInnes MATOP_SOLVE=7, 14601c1c02c0SLois Curfman McInnes MATOP_SOLVE_ADD=8, 14617c922b88SBarry Smith MATOP_SOLVE_TRANSPOSE=9, 14627c922b88SBarry Smith MATOP_SOLVE_TRANSPOSE_ADD=10, 14631c1c02c0SLois Curfman McInnes MATOP_LUFACTOR=11, 14641c1c02c0SLois Curfman McInnes MATOP_CHOLESKYFACTOR=12, 1465a714c06dSBarry Smith MATOP_SOR=13, 14661c1c02c0SLois Curfman McInnes MATOP_TRANSPOSE=14, 14671c1c02c0SLois Curfman McInnes MATOP_GETINFO=15, 14681c1c02c0SLois Curfman McInnes MATOP_EQUAL=16, 14691c1c02c0SLois Curfman McInnes MATOP_GET_DIAGONAL=17, 14701c1c02c0SLois Curfman McInnes MATOP_DIAGONAL_SCALE=18, 14711c1c02c0SLois Curfman McInnes MATOP_NORM=19, 14721c1c02c0SLois Curfman McInnes MATOP_ASSEMBLY_BEGIN=20, 14731c1c02c0SLois Curfman McInnes MATOP_ASSEMBLY_END=21, 1474d519adbfSMatthew Knepley MATOP_SET_OPTION=22, 1475d519adbfSMatthew Knepley MATOP_ZERO_ENTRIES=23, 1476d519adbfSMatthew Knepley MATOP_ZERO_ROWS=24, 1477d519adbfSMatthew Knepley MATOP_LUFACTOR_SYMBOLIC=25, 1478d519adbfSMatthew Knepley MATOP_LUFACTOR_NUMERIC=26, 1479d519adbfSMatthew Knepley MATOP_CHOLESKY_FACTOR_SYMBOLIC=27, 1480d519adbfSMatthew Knepley MATOP_CHOLESKY_FACTOR_NUMERIC=28, 1481d519adbfSMatthew Knepley MATOP_SETUP_PREALLOCATION=29, 1482d519adbfSMatthew Knepley MATOP_ILUFACTOR_SYMBOLIC=30, 1483d519adbfSMatthew Knepley MATOP_ICCFACTOR_SYMBOLIC=31, 1484d519adbfSMatthew Knepley MATOP_GET_ARRAY=32, 1485d519adbfSMatthew Knepley MATOP_RESTORE_ARRAY=33, 1486d519adbfSMatthew Knepley MATOP_DUPLICATE=34, 1487d519adbfSMatthew Knepley MATOP_FORWARD_SOLVE=35, 1488d519adbfSMatthew Knepley MATOP_BACKWARD_SOLVE=36, 1489d519adbfSMatthew Knepley MATOP_ILUFACTOR=37, 1490d519adbfSMatthew Knepley MATOP_ICCFACTOR=38, 1491d519adbfSMatthew Knepley MATOP_AXPY=39, 1492d519adbfSMatthew Knepley MATOP_GET_SUBMATRICES=40, 1493d519adbfSMatthew Knepley MATOP_INCREASE_OVERLAP=41, 1494d519adbfSMatthew Knepley MATOP_GET_VALUES=42, 1495d519adbfSMatthew Knepley MATOP_COPY=43, 1496d519adbfSMatthew Knepley MATOP_GET_ROW_MAX=44, 1497d519adbfSMatthew Knepley MATOP_SCALE=45, 1498d519adbfSMatthew Knepley MATOP_SHIFT=46, 149935153367SBarry Smith MATOP_DIAGONAL_SET=47, 1500d519adbfSMatthew Knepley MATOP_ILUDT_FACTOR=48, 1501d519adbfSMatthew Knepley MATOP_SET_BLOCK_SIZE=49, 1502d519adbfSMatthew Knepley MATOP_GET_ROW_IJ=50, 1503d519adbfSMatthew Knepley MATOP_RESTORE_ROW_IJ=51, 1504d519adbfSMatthew Knepley MATOP_GET_COLUMN_IJ=52, 1505d519adbfSMatthew Knepley MATOP_RESTORE_COLUMN_IJ=53, 1506d519adbfSMatthew Knepley MATOP_FDCOLORING_CREATE=54, 1507d519adbfSMatthew Knepley MATOP_COLORING_PATCH=55, 1508d519adbfSMatthew Knepley MATOP_SET_UNFACTORED=56, 1509d519adbfSMatthew Knepley MATOP_PERMUTE=57, 1510d519adbfSMatthew Knepley MATOP_SET_VALUES_BLOCKED=58, 1511d519adbfSMatthew Knepley MATOP_GET_SUBMATRIX=59, 1512d519adbfSMatthew Knepley MATOP_DESTROY=60, 1513d519adbfSMatthew Knepley MATOP_VIEW=61, 1514d519adbfSMatthew Knepley MATOP_CONVERT_FROM=62, 1515d519adbfSMatthew Knepley MATOP_USE_SCALED_FORM=63, 1516d519adbfSMatthew Knepley MATOP_SCALE_SYSTEM=64, 1517d519adbfSMatthew Knepley MATOP_UNSCALE_SYSTEM=65, 1518d519adbfSMatthew Knepley MATOP_SET_LOCAL_TO_GLOBAL_MAP=66, 1519d519adbfSMatthew Knepley MATOP_SET_VALUES_LOCAL=67, 1520d519adbfSMatthew Knepley MATOP_ZERO_ROWS_LOCAL=68, 1521d519adbfSMatthew Knepley MATOP_GET_ROW_MAX_ABS=69, 1522d519adbfSMatthew Knepley MATOP_GET_ROW_MIN_ABS=70, 1523d519adbfSMatthew Knepley MATOP_CONVERT=71, 1524d519adbfSMatthew Knepley MATOP_SET_COLORING=72, 1525d519adbfSMatthew Knepley MATOP_SET_VALUES_ADIC=73, 1526d519adbfSMatthew Knepley MATOP_SET_VALUES_ADIFOR=74, 1527d519adbfSMatthew Knepley MATOP_FD_COLORING_APPLY=75, 1528d519adbfSMatthew Knepley MATOP_SET_FROM_OPTIONS=76, 1529d519adbfSMatthew Knepley MATOP_MULT_CON=77, 1530d519adbfSMatthew Knepley MATOP_MULT_TRANSPOSE_CON=78, 1531d519adbfSMatthew Knepley MATOP_PERMUTE_SPARSIFY=79, 1532d519adbfSMatthew Knepley MATOP_MULT_MULTIPLE=80, 1533d519adbfSMatthew Knepley MATOP_SOLVE_MULTIPLE=81, 1534d519adbfSMatthew Knepley MATOP_GET_INERTIA=82, 1535d519adbfSMatthew Knepley MATOP_LOAD=83, 1536d519adbfSMatthew Knepley MATOP_IS_SYMMETRIC=84, 1537d519adbfSMatthew Knepley MATOP_IS_HERMITIAN=85, 1538d519adbfSMatthew Knepley MATOP_IS_STRUCTURALLY_SYMMETRIC=86, 153941f059aeSBarry Smith MATOP_DUMMY=87, 1540d519adbfSMatthew Knepley MATOP_GET_VECS=88, 1541d519adbfSMatthew Knepley MATOP_MAT_MULT=89, 1542d519adbfSMatthew Knepley MATOP_MAT_MULT_SYMBOLIC=90, 1543d519adbfSMatthew Knepley MATOP_MAT_MULT_NUMERIC=91, 1544d519adbfSMatthew Knepley MATOP_PTAP=92, 1545d519adbfSMatthew Knepley MATOP_PTAP_SYMBOLIC=93, 1546d519adbfSMatthew Knepley MATOP_PTAP_NUMERIC=94, 1547d519adbfSMatthew Knepley MATOP_MAT_MULTTRANSPOSE=95, 15481763ac4eSSatish Balay MATOP_MAT_MULTTRANSPOSE_SYM=96, 15491763ac4eSSatish Balay MATOP_MAT_MULTTRANSPOSE_NUM=97, 1550d519adbfSMatthew Knepley MATOP_PTAP_SYMBOLIC_SEQAIJ=98, 1551d519adbfSMatthew Knepley MATOP_PTAP_NUMERIC_SEQAIJ=99, 1552d519adbfSMatthew Knepley MATOP_PTAP_SYMBOLIC_MPIAIJ=100, 1553d519adbfSMatthew Knepley MATOP_PTAP_NUMERIC_MPIAIJ=101, 1554d519adbfSMatthew Knepley MATOP_CONJUGATE=102, 1555d519adbfSMatthew Knepley MATOP_SET_SIZES=103, 1556d519adbfSMatthew Knepley MATOP_SET_VALUES_ROW=104, 1557d519adbfSMatthew Knepley MATOP_REAL_PART=105, 1558d519adbfSMatthew Knepley MATOP_IMAG_PART=106, 1559d519adbfSMatthew Knepley MATOP_GET_ROW_UTRIANGULAR=107, 1560d519adbfSMatthew Knepley MATOP_RESTORE_ROW_UTRIANGULAR=108, 1561d519adbfSMatthew Knepley MATOP_MATSOLVE=109, 1562d519adbfSMatthew Knepley MATOP_GET_REDUNDANTMATRIX=110, 1563d519adbfSMatthew Knepley MATOP_GET_ROW_MIN=111, 1564d519adbfSMatthew Knepley MATOP_GET_COLUMN_VEC=112, 1565d519adbfSMatthew Knepley MATOP_MISSING_DIAGONAL=113, 1566d519adbfSMatthew Knepley MATOP_MATGETSEQNONZEROSTRUCTURE=114, 156789c6957cSBarry Smith MATOP_CREATE=115, 156889c6957cSBarry Smith MATOP_GET_GHOSTS=116, 1569ea38565cSJed Brown MATOP_GET_LOCALSUBMATRIX=117, 1570ea38565cSJed Brown MATOP_RESTORE_LOCALSUBMATRIX=118, 1571eeffb40dSHong Zhang MATOP_MULT_DIAGONAL_BLOCK=119, 1572547795f9SHong Zhang MATOP_HERMITIANTRANSPOSE=120, 1573547795f9SHong Zhang MATOP_MULTHERMITIANTRANSPOSE=121, 1574fbdbba38SShri Abhyankar MATOP_MULTHERMITIANTRANSPOSEADD=122, 1575b9614d88SDmitry Karpeev MATOP_GETMULTIPROCBLOCK=123, 15760716a85fSBarry Smith MATOP_GETCOLUMNNORMS=125, 1577b9614d88SDmitry Karpeev MATOP_GET_SUBMATRICES_PARALLEL=128 1578fae171e0SBarry Smith } MatOperation; 15797087cfbeSBarry Smith extern PetscErrorCode MatHasOperation(Mat,MatOperation,PetscBool *); 15807087cfbeSBarry Smith extern PetscErrorCode MatShellSetOperation(Mat,MatOperation,void(*)(void)); 15817087cfbeSBarry Smith extern PetscErrorCode MatShellGetOperation(Mat,MatOperation,void(**)(void)); 15827087cfbeSBarry Smith extern PetscErrorCode MatShellSetContext(Mat,void*); 1583112a2221SBarry Smith 158490ace30eSBarry Smith /* 158590ace30eSBarry Smith Codes for matrices stored on disk. By default they are 158690ace30eSBarry Smith stored in a universal format. By changing the format with 15877973ad04SBarry Smith PetscViewerSetFormat(viewer,PETSC_VIEWER_NATIVE); the matrices will 158890ace30eSBarry Smith be stored in a way natural for the matrix, for example dense matrices 158990ace30eSBarry Smith would be stored as dense. Matrices stored this way may only be 1590f4403165SShri Abhyankar read into matrices of the same type. 159190ace30eSBarry Smith */ 159290ace30eSBarry Smith #define MATRIX_BINARY_FORMAT_DENSE -1 159390ace30eSBarry Smith 15947087cfbeSBarry Smith extern PetscErrorCode MatMPIBAIJSetHashTableFactor(Mat,PetscReal); 15957087cfbeSBarry Smith extern PetscErrorCode MatISGetLocalMat(Mat,Mat*); 15961f4e1ec7SBarry Smith 1597d9274352SBarry Smith /*S 1598d9274352SBarry Smith MatNullSpace - Object that removes a null space from a vector, i.e. 1599d9274352SBarry Smith orthogonalizes the vector to a subsapce 1600d9274352SBarry Smith 1601f7a9e4ceSBarry Smith Level: advanced 1602d9274352SBarry Smith 1603d9274352SBarry Smith Concepts: matrix; linear operator, null space 1604d9274352SBarry Smith 16056e1639daSBarry Smith Users manual sections: 16066e1639daSBarry Smith . sec_singular 16076e1639daSBarry Smith 1608d9274352SBarry Smith .seealso: MatNullSpaceCreate() 1609d9274352SBarry Smith S*/ 161074637425SBarry Smith typedef struct _p_MatNullSpace* MatNullSpace; 1611d9274352SBarry Smith 16127087cfbeSBarry Smith extern PetscErrorCode MatNullSpaceCreate(MPI_Comm,PetscBool ,PetscInt,const Vec[],MatNullSpace*); 16137087cfbeSBarry Smith extern PetscErrorCode MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(MatNullSpace,Vec,void*),void*); 1614d34fcf5fSBarry Smith extern PetscErrorCode MatNullSpaceDestroy(MatNullSpace*); 16157087cfbeSBarry Smith extern PetscErrorCode MatNullSpaceRemove(MatNullSpace,Vec,Vec*); 16167087cfbeSBarry Smith extern PetscErrorCode MatNullSpaceAttach(Mat,MatNullSpace); 16177087cfbeSBarry Smith extern PetscErrorCode MatNullSpaceTest(MatNullSpace,Mat,PetscBool *); 1618b717e993SJed Brown extern PetscErrorCode MatNullSpaceView(MatNullSpace,PetscViewer); 161974637425SBarry Smith 16207087cfbeSBarry Smith extern PetscErrorCode MatReorderingSeqSBAIJ(Mat,IS); 16217087cfbeSBarry Smith extern PetscErrorCode MatMPISBAIJSetHashTableFactor(Mat,PetscReal); 16227087cfbeSBarry Smith extern PetscErrorCode MatSeqSBAIJSetColumnIndices(Mat,PetscInt *); 16237087cfbeSBarry Smith extern PetscErrorCode MatSeqBAIJInvertBlockDiagonal(Mat); 16243f1d51d7SBarry Smith 16257087cfbeSBarry Smith extern PetscErrorCode MatCreateMAIJ(Mat,PetscInt,Mat*); 16267087cfbeSBarry Smith extern PetscErrorCode MatMAIJRedimension(Mat,PetscInt,Mat*); 16277087cfbeSBarry Smith extern PetscErrorCode MatMAIJGetAIJ(Mat,Mat*); 1628c4f061fbSSatish Balay 16297087cfbeSBarry Smith extern PetscErrorCode MatComputeExplicitOperator(Mat,Mat*); 1630b0a32e0cSBarry Smith 16317087cfbeSBarry Smith extern PetscErrorCode MatDiagonalScaleLocal(Mat,Vec); 163204f1ad80SBarry Smith 16337087cfbeSBarry Smith extern PetscErrorCode MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*); 16347087cfbeSBarry Smith extern PetscErrorCode MatMFFDSetBase(Mat,Vec,Vec); 16357087cfbeSBarry Smith extern PetscErrorCode MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*); 16367087cfbeSBarry Smith extern PetscErrorCode MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*)); 16377087cfbeSBarry Smith extern PetscErrorCode MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec)); 16387087cfbeSBarry Smith extern PetscErrorCode MatMFFDAddNullSpace(Mat,MatNullSpace); 16397087cfbeSBarry Smith extern PetscErrorCode MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt); 16407087cfbeSBarry Smith extern PetscErrorCode MatMFFDResetHHistory(Mat); 16417087cfbeSBarry Smith extern PetscErrorCode MatMFFDSetFunctionError(Mat,PetscReal); 16427087cfbeSBarry Smith extern PetscErrorCode MatMFFDSetPeriod(Mat,PetscInt); 16437087cfbeSBarry Smith extern PetscErrorCode MatMFFDGetH(Mat,PetscScalar *); 16447087cfbeSBarry Smith extern PetscErrorCode MatMFFDSetOptionsPrefix(Mat,const char[]); 16457087cfbeSBarry Smith extern PetscErrorCode MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*); 16467087cfbeSBarry Smith extern PetscErrorCode MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*); 1647e884886fSBarry Smith 16486370053bSBarry Smith /*S 16496370053bSBarry Smith MatMFFD - A data structured used to manage the computation of the h differencing parameter for matrix-free 16506370053bSBarry Smith Jacobian vector products 1651e884886fSBarry Smith 16526370053bSBarry Smith Notes: MATMFFD is a specific MatType which uses the MatMFFD data structure 16536370053bSBarry Smith 16546370053bSBarry Smith MatMFFD*() methods actually take the Mat as their first argument. Not a MatMFFD data structure 16556370053bSBarry Smith 16566370053bSBarry Smith Level: developer 16576370053bSBarry Smith 16586370053bSBarry Smith .seealso: MATMFFD, MatCreateMFFD(), MatMFFDSetFuction(), MatMFFDSetType(), MatMFFDRegister() 16596370053bSBarry Smith S*/ 1660e884886fSBarry Smith typedef struct _p_MatMFFD* MatMFFD; 1661e884886fSBarry Smith 1662e884886fSBarry Smith /*E 1663e884886fSBarry Smith MatMFFDType - algorithm used to compute the h used in computing matrix-vector products via differencing of the function 1664e884886fSBarry Smith 1665e884886fSBarry Smith Level: beginner 1666e884886fSBarry Smith 1667e884886fSBarry Smith .seealso: MatMFFDSetType(), MatMFFDRegister() 1668e884886fSBarry Smith E*/ 1669a313700dSBarry Smith #define MatMFFDType char* 1670e884886fSBarry Smith #define MATMFFD_DS "ds" 1671e884886fSBarry Smith #define MATMFFD_WP "wp" 1672e884886fSBarry Smith 16737087cfbeSBarry Smith extern PetscErrorCode MatMFFDSetType(Mat,const MatMFFDType); 16747087cfbeSBarry Smith extern PetscErrorCode MatMFFDRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatMFFD)); 1675e884886fSBarry Smith 1676e884886fSBarry Smith /*MC 1677e884886fSBarry Smith MatMFFDRegisterDynamic - Adds a method to the MatMFFD registry. 1678e884886fSBarry Smith 1679e884886fSBarry Smith Synopsis: 16801890ba74SBarry Smith PetscErrorCode MatMFFDRegisterDynamic(const char *name_solver,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatMFFD)) 1681e884886fSBarry Smith 1682e884886fSBarry Smith Not Collective 1683e884886fSBarry Smith 1684e884886fSBarry Smith Input Parameters: 1685e884886fSBarry Smith + name_solver - name of a new user-defined compute-h module 1686e884886fSBarry Smith . path - path (either absolute or relative) the library containing this solver 1687e884886fSBarry Smith . name_create - name of routine to create method context 1688e884886fSBarry Smith - routine_create - routine to create method context 1689e884886fSBarry Smith 1690e884886fSBarry Smith Level: developer 1691e884886fSBarry Smith 1692e884886fSBarry Smith Notes: 1693e884886fSBarry Smith MatMFFDRegisterDynamic() may be called multiple times to add several user-defined solvers. 1694e884886fSBarry Smith 1695e884886fSBarry Smith If dynamic libraries are used, then the fourth input argument (routine_create) 1696e884886fSBarry Smith is ignored. 1697e884886fSBarry Smith 1698e884886fSBarry Smith Sample usage: 1699e884886fSBarry Smith .vb 1700e884886fSBarry Smith MatMFFDRegisterDynamic("my_h",/home/username/my_lib/lib/libO/solaris/mylib.a, 1701e884886fSBarry Smith "MyHCreate",MyHCreate); 1702e884886fSBarry Smith .ve 1703e884886fSBarry Smith 1704e884886fSBarry Smith Then, your solver can be chosen with the procedural interface via 1705e884886fSBarry Smith $ MatMFFDSetType(mfctx,"my_h") 1706e884886fSBarry Smith or at runtime via the option 1707e884886fSBarry Smith $ -snes_mf_type my_h 1708e884886fSBarry Smith 1709e884886fSBarry Smith .keywords: MatMFFD, register 1710e884886fSBarry Smith 1711e884886fSBarry Smith .seealso: MatMFFDRegisterAll(), MatMFFDRegisterDestroy() 1712e884886fSBarry Smith M*/ 1713e884886fSBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 1714e884886fSBarry Smith #define MatMFFDRegisterDynamic(a,b,c,d) MatMFFDRegister(a,b,c,0) 1715e884886fSBarry Smith #else 1716e884886fSBarry Smith #define MatMFFDRegisterDynamic(a,b,c,d) MatMFFDRegister(a,b,c,d) 1717e884886fSBarry Smith #endif 1718e884886fSBarry Smith 17197087cfbeSBarry Smith extern PetscErrorCode MatMFFDRegisterAll(const char[]); 17207087cfbeSBarry Smith extern PetscErrorCode MatMFFDRegisterDestroy(void); 17217087cfbeSBarry Smith extern PetscErrorCode MatMFFDDSSetUmin(Mat,PetscReal); 17227087cfbeSBarry Smith extern PetscErrorCode MatMFFDWPSetComputeNormU(Mat,PetscBool ); 1723e884886fSBarry Smith 1724e884886fSBarry Smith 17257087cfbeSBarry Smith extern PetscErrorCode PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *); 17267087cfbeSBarry Smith extern PetscErrorCode PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *); 17277dbadf16SMatthew Knepley 172897969023SHong Zhang /* 172997969023SHong Zhang PETSc interface to MUMPS 173097969023SHong Zhang */ 173197969023SHong Zhang #ifdef PETSC_HAVE_MUMPS 1732f250808bSBarry Smith extern PetscErrorCode MatMumpsSetIcntl(Mat,PetscInt,PetscInt); 173397969023SHong Zhang #endif 173423a5497aSJed Brown 1735d954db57SHong Zhang /* 1736d954db57SHong Zhang PETSc interface to SUPERLU 1737d954db57SHong Zhang */ 1738d954db57SHong Zhang #ifdef PETSC_HAVE_SUPERLU 1739f250808bSBarry Smith extern PetscErrorCode MatSuperluSetILUDropTol(Mat,PetscReal); 1740d954db57SHong Zhang #endif 1741d954db57SHong Zhang 1742c8883902SJed Brown extern PetscErrorCode MatCreateNest(MPI_Comm,PetscInt,const IS[],PetscInt,const IS[],const Mat[],Mat*); 1743c8883902SJed Brown extern PetscErrorCode MatNestGetSize(Mat,PetscInt*,PetscInt*); 1744c8883902SJed Brown extern PetscErrorCode MatNestGetSubMats(Mat,PetscInt*,PetscInt*,Mat***); 1745c8883902SJed Brown extern PetscErrorCode MatNestGetSubMat(Mat,PetscInt,PetscInt,Mat*); 17467087cfbeSBarry Smith extern PetscErrorCode MatNestSetVecType(Mat,const VecType); 1747c8883902SJed Brown extern PetscErrorCode MatNestSetSubMats(Mat,PetscInt,const IS[],PetscInt,const IS[],const Mat[]); 1748d8588912SDave May 174923a5497aSJed Brown PETSC_EXTERN_CXX_END 175023a5497aSJed Brown #endif 1751