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" 31273d9f13SBarry Smith #define MATSEQMAIJ "seqmaij" 32273d9f13SBarry Smith #define MATMPIMAIJ "mpimaij" 33209238afSKris Buschelman #define MATMAIJ "maij" 34273d9f13SBarry Smith #define MATIS "is" 35273d9f13SBarry Smith #define MATMPIROWBS "mpirowbs" 36273d9f13SBarry Smith #define MATSEQAIJ "seqaij" 37273d9f13SBarry Smith #define MATMPIAIJ "mpiaij" 38209238afSKris Buschelman #define MATAIJ "aij" 39273d9f13SBarry Smith #define MATSHELL "shell" 40209238afSKris Buschelman #define MATSEQDENSE "seqdense" 41273d9f13SBarry Smith #define MATMPIDENSE "mpidense" 42209238afSKris Buschelman #define MATDENSE "dense" 43273d9f13SBarry Smith #define MATSEQBAIJ "seqbaij" 44273d9f13SBarry Smith #define MATMPIBAIJ "mpibaij" 45209238afSKris Buschelman #define MATBAIJ "baij" 46273d9f13SBarry Smith #define MATMPIADJ "mpiadj" 47273d9f13SBarry Smith #define MATSEQSBAIJ "seqsbaij" 48273d9f13SBarry Smith #define MATMPISBAIJ "mpisbaij" 49209238afSKris Buschelman #define MATSBAIJ "sbaij" 50cebc7f6cSBarry Smith #define MATDAAD "daad" 51cebc7f6cSBarry Smith #define MATMFFD "mffd" 52c8a8475eSBarry Smith #define MATNORMAL "normal" 53ab92ecdeSBarry Smith #define MATLRC "lrc" 544b8d99adSRichard Tran Mills #define MATSEQCSRPERM "seqcsrperm" 5518def556SRichard Tran Mills #define MATMPICSRPERM "mpicsrperm" 56814655b8SBarry Smith #define MATCSRPERM "csrperm" 5781824310SBarry Smith #define MATSEQCRL "seqcrl" 58c02b24c5SBarry Smith #define MATMPICRL "mpicrl" 59c02b24c5SBarry Smith #define MATCRL "crl" 602a6744ebSBarry Smith #define MATSCATTER "scatter" 61421e10b8SBarry Smith #define MATBLOCKMAT "blockmat" 62793850ffSBarry Smith #define MATCOMPOSITE "composite" 635a7f1df3SHong Zhang #define MATSEQFFTW "seqfftw" 64557cca28SSatish Balay #define MATTRANSPOSEMAT "transpose" 6572ca8751SBarry Smith #define MATSCHURCOMPLEMENT "schurcomplement" 66d91e6319SBarry Smith 679989ab13SBarry Smith /*E 68c7393fdbSBarry Smith MatSolverPackage - String with the name of a PETSc matrix solver type. 699989ab13SBarry Smith 709989ab13SBarry Smith For example: "petsc" indicates what PETSc provides, "superlu" indicates either 719989ab13SBarry Smith SuperLU or SuperLU_Dist etc. 729989ab13SBarry Smith 739989ab13SBarry Smith 749989ab13SBarry Smith Level: beginner 759989ab13SBarry Smith 765c9eb25fSBarry Smith .seealso: MatGetFactor(), Mat, MatSetType(), MatType 779989ab13SBarry Smith E*/ 78c7393fdbSBarry Smith #define MatSolverPackage char* 79b24902e0SBarry Smith #define MAT_SOLVER_SPOOLES "spooles" 80b24902e0SBarry Smith #define MAT_SOLVER_SUPERLU "superlu" 81b24902e0SBarry Smith #define MAT_SOLVER_SUPERLU_DIST "superlu_dist" 82b24902e0SBarry Smith #define MAT_SOLVER_UMFPACK "umfpack" 83b24902e0SBarry Smith #define MAT_SOLVER_ESSL "essl" 84b24902e0SBarry Smith #define MAT_SOLVER_LUSOL "lusol" 85b24902e0SBarry Smith #define MAT_SOLVER_MUMPS "mumps" 863bf14a46SMatthew Knepley #define MAT_SOLVER_PASTIX "pastix" 87b24902e0SBarry Smith #define MAT_SOLVER_DSCPACK "dscpack" 88b24902e0SBarry Smith #define MAT_SOLVER_MATLAB "matlab" 8943244d56SBarry Smith #define MAT_SOLVER_PETSC "petsc" 90*b6806ab0SHong Zhang #define MAT_SOLVER_PLAPACK "plapack" 91b24902e0SBarry Smith 92b24902e0SBarry Smith /*E 93b24902e0SBarry Smith MatFactorType - indicates what type of factorization is requested 94b24902e0SBarry Smith 95b24902e0SBarry Smith Level: beginner 96b24902e0SBarry Smith 97b24902e0SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 98b24902e0SBarry Smith 99c7393fdbSBarry Smith .seealso: MatSolverPackage, MatGetFactor() 100b24902e0SBarry Smith E*/ 1015c9eb25fSBarry Smith typedef enum {MAT_FACTOR_NONE, MAT_FACTOR_LU, MAT_FACTOR_CHOLESKY, MAT_FACTOR_ILU, MAT_FACTOR_ICC} MatFactorType; 102c7393fdbSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetFactor(Mat,const MatSolverPackage,MatFactorType,Mat*); 103db4efbfdSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetFactorAvailable(Mat,const MatSolverPackage,MatFactorType,PetscTruth*); 10435bd34faSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFactorGetSolverPackage(Mat,const MatSolverPackage*); 105b24902e0SBarry Smith 1069989ab13SBarry Smith 107c06d978dSMatthew Knepley /* Logging support */ 108552e946dSBarry Smith #define MAT_FILE_COOKIE 1211216 /* used to indicate matrices in binary files */ 109be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_COOKIE; 110be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_FDCOLORING_COOKIE; 111be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_PARTITIONING_COOKIE; 112be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_NULLSPACE_COOKIE; 113e884886fSBarry Smith extern PetscCookie PETSCMAT_DLLEXPORT MATMFFD_COOKIE; 114c06d978dSMatthew Knepley 115ceb03754SKris Buschelman /*E 116ceb03754SKris Buschelman MatReuse - Indicates if matrices obtained from a previous call to MatGetSubMatrices() 117ceb03754SKris Buschelman or MatGetSubMatrix() are to be reused to store the new matrix values. 118ceb03754SKris Buschelman 119ceb03754SKris Buschelman Level: beginner 120ceb03754SKris Buschelman 121ceb03754SKris Buschelman Any additions/changes here MUST also be made in include/finclude/petscmat.h 122ceb03754SKris Buschelman 1230c451bc4SBarry Smith .seealso: MatGetSubMatrices(), MatGetSubMatrix(), MatDestroyMatrices(), MatConvert() 124ceb03754SKris Buschelman E*/ 125ceb03754SKris Buschelman typedef enum {MAT_INITIAL_MATRIX,MAT_REUSE_MATRIX} MatReuse; 126ceb03754SKris Buschelman 1275494a064SHong Zhang /*E 1285494a064SHong Zhang MatGetSubMatrixOption - Indicates if matrices obtained from a call to MatGetSubMatrices() 129829201f2SHong Zhang include the matrix values. Currently it is only used by MatGetSeqNonzerostructure(). 1305494a064SHong Zhang 1315494a064SHong Zhang Level: beginner 1325494a064SHong Zhang 133829201f2SHong Zhang .seealso: MatGetSeqNonzerostructure() 1345494a064SHong Zhang E*/ 1355494a064SHong Zhang typedef enum {MAT_DO_NOT_GET_VALUES,MAT_GET_VALUES} MatGetSubMatrixOption; 1365494a064SHong Zhang 137e5bd5246SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInitializePackage(const char[]); 138c06d978dSMatthew Knepley 139f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreate(MPI_Comm,Mat*); 140a69a960dSMatthew Knepley PetscPolymorphicFunction(MatCreate,(MPI_Comm comm),(comm,&A),Mat,A) 141a69a960dSMatthew Knepley PetscPolymorphicFunction(MatCreate,(),(PETSC_COMM_WORLD,&A),Mat,A) 142f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetSizes(Mat,PetscInt,PetscInt,PetscInt,PetscInt); 143a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetType(Mat,const MatType); 144be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetFromOptions(Mat); 145be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUpPreallocation(Mat); 146be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRegisterAll(const char[]); 147be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat)); 14830de9b25SBarry Smith 149f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetOptionsPrefix(Mat,const char[]); 150f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAppendOptionsPrefix(Mat,const char[]); 151f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOptionsPrefix(Mat,const char*[]); 152f69a0ea3SMatthew Knepley 15330de9b25SBarry Smith /*MC 15430de9b25SBarry Smith MatRegisterDynamic - Adds a new matrix type 15530de9b25SBarry Smith 15630de9b25SBarry Smith Synopsis: 157c1ac3661SBarry Smith PetscErrorCode MatRegisterDynamic(char *name,char *path,char *name_create,PetscErrorCode (*routine_create)(Mat)) 15830de9b25SBarry Smith 15930de9b25SBarry Smith Not Collective 16030de9b25SBarry Smith 16130de9b25SBarry Smith Input Parameters: 16230de9b25SBarry Smith + name - name of a new user-defined matrix type 16330de9b25SBarry Smith . path - path (either absolute or relative) the library containing this solver 16430de9b25SBarry Smith . name_create - name of routine to create method context 16530de9b25SBarry Smith - routine_create - routine to create method context 16630de9b25SBarry Smith 16730de9b25SBarry Smith Notes: 16830de9b25SBarry Smith MatRegisterDynamic() may be called multiple times to add several user-defined solvers. 16930de9b25SBarry Smith 17030de9b25SBarry Smith If dynamic libraries are used, then the fourth input argument (routine_create) 17130de9b25SBarry Smith is ignored. 17230de9b25SBarry Smith 17330de9b25SBarry Smith Sample usage: 17430de9b25SBarry Smith .vb 17530de9b25SBarry Smith MatRegisterDynamic("my_mat",/home/username/my_lib/lib/libO/solaris/mylib.a, 17630de9b25SBarry Smith "MyMatCreate",MyMatCreate); 17730de9b25SBarry Smith .ve 17830de9b25SBarry Smith 17930de9b25SBarry Smith Then, your solver can be chosen with the procedural interface via 18030de9b25SBarry Smith $ MatSetType(Mat,"my_mat") 18130de9b25SBarry Smith or at runtime via the option 18230de9b25SBarry Smith $ -mat_type my_mat 18330de9b25SBarry Smith 18430de9b25SBarry Smith Level: advanced 18530de9b25SBarry Smith 186ab901514SBarry Smith Notes: ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values. 18730de9b25SBarry Smith If your function is not being put into a shared library then use VecRegister() instead 18830de9b25SBarry Smith 18930de9b25SBarry Smith .keywords: Mat, register 19030de9b25SBarry Smith 19130de9b25SBarry Smith .seealso: MatRegisterAll(), MatRegisterDestroy() 19230de9b25SBarry Smith 19330de9b25SBarry Smith M*/ 194273d9f13SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 195273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,0) 196273d9f13SBarry Smith #else 197273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,d) 19830de9b25SBarry Smith #endif 19930de9b25SBarry Smith 200273d9f13SBarry Smith extern PetscTruth MatRegisterAllCalled; 201b0a32e0cSBarry Smith extern PetscFList MatList; 20228988994SBarry Smith 2033b224e63SBarry Smith /*E 2043b224e63SBarry Smith MatStructure - Indicates if the matrix has the same nonzero structure 2053b224e63SBarry Smith 2063b224e63SBarry Smith Level: beginner 2073b224e63SBarry Smith 2083b224e63SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 2093b224e63SBarry Smith 2103b224e63SBarry Smith .seealso: MatCopy(), KSPSetOperators(), PCSetOperators() 2113b224e63SBarry Smith E*/ 2123b224e63SBarry Smith typedef enum {SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER,SUBSET_NONZERO_PATTERN} MatStructure; 2133b224e63SBarry Smith 214be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqDense(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*); 215be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIDense(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*); 216be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 217ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,nz,nnz,&A),Mat,A) 218ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,&A),Mat,A) 219ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,0,nnz,&A),Mat,A) 220ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,&A),Mat,A) 221ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,A)) 222ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,m,n,0,nnz,A)) 223ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,A)) 224be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 225ba966639SSatish 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) 226ba966639SSatish 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) 227ba966639SSatish 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) 228ba966639SSatish 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) 229ba966639SSatish 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)) 230ba966639SSatish 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)) 231ba966639SSatish 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)) 232ba966639SSatish 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) 233ba966639SSatish 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) 234ba966639SSatish 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) 235ba966639SSatish 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) 236ba966639SSatish 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)) 237ba966639SSatish 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)) 238ba966639SSatish 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)) 23963ab93bfSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *); 2408d7a6e47SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJWithSplitArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],PetscInt[],PetscInt[],PetscScalar[],Mat*); 2418d7a6e47SBarry Smith 242be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIRowbs(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 243be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 244ba966639SSatish 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) 245ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A) 246ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A) 247ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A) 248ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A)) 249ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A)) 250ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A)) 251be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 252ba966639SSatish 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) 253ba966639SSatish 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) 254ba966639SSatish 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) 255ba966639SSatish 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) 256ba966639SSatish 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)) 257ba966639SSatish 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)) 258ba966639SSatish 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)) 259ba966639SSatish 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) 260ba966639SSatish 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) 261ba966639SSatish 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) 262ba966639SSatish 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) 263ba966639SSatish 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)) 264ba966639SSatish 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)) 265ba966639SSatish 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)) 266be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAdj(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscInt[],Mat*); 267be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 268ba966639SSatish 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) 269ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A) 270ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A) 271ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A) 272ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A)) 273ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A)) 274ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A)) 275be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPISBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 276ba966639SSatish 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) 277ba966639SSatish 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) 278ba966639SSatish 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) 279ba966639SSatish 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) 280ba966639SSatish 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)) 281ba966639SSatish 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)) 282ba966639SSatish 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)) 283ba966639SSatish 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) 284ba966639SSatish 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) 285ba966639SSatish 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) 286ba966639SSatish 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) 287ba966639SSatish 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)) 288ba966639SSatish 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)) 289ba966639SSatish 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)) 290be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateShell(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,void *,Mat*); 291ba966639SSatish 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) 292ba966639SSatish Balay PetscPolymorphicFunction(MatCreateShell,(PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx),(PETSC_COMM_WORLD,m,n,M,N,ctx,&A),Mat,A) 293be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateAdic(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,void (*)(void),Mat*); 294be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateNormal(Mat,Mat*); 295ba966639SSatish Balay PetscPolymorphicFunction(MatCreateNormal,(Mat mat),(mat,&A),Mat,A) 296ab92ecdeSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateLRC(Mat,Mat,Mat,Mat*); 297284134d9SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateIS(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,ISLocalToGlobalMapping,Mat*); 2981472f72bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 2991472f72bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPICRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 3002a6744ebSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateScatter(MPI_Comm,VecScatter,Mat*); 3012a6744ebSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScatterSetVecScatter(Mat,VecScatter); 3022a6744ebSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScatterGetVecScatter(Mat,VecScatter*); 3038cdf2d9bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateBlockMat(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,Mat*); 304793850ffSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompositeAddMat(Mat,Mat); 305b52f573bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompositeMerge(Mat); 306793850ffSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateComposite(MPI_Comm,PetscInt,const Mat*,Mat*); 3076d7c1e57SBarry Smith typedef enum {MAT_COMPOSITE_ADDITIVE,MAT_COMPOSITE_MULTIPLICATIVE} MatCompositeType; 3086d7c1e57SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompositeSetType(Mat,MatCompositeType); 3096d7c1e57SBarry Smith 3105a7f1df3SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqFFTW(MPI_Comm,PetscInt,const PetscInt[],Mat*); 311f20085c4SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateTranspose(Mat,Mat*); 3121472f72bSBarry Smith 313f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUp(Mat); 314be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroy(Mat); 31521c89e3eSBarry Smith 3160c451bc4SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConjugate(Mat); 31799cafbc1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRealPart(Mat); 31899cafbc1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatImaginaryPart(Mat); 31999cafbc1SBarry Smith 3208ed539a5SBarry Smith /* ------------------------------------------------------------*/ 321be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 322be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 32387d4246cSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesRow(Mat,PetscInt,const PetscScalar[]); 324f259bd47SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesRowLocal(Mat,PetscInt,const PetscScalar[]); 32584cb2905SBarry Smith 3262ef4de8bSBarry Smith /*S 3272ef4de8bSBarry Smith MatStencil - Data structure (C struct) for storing information about a single row or 3282ef4de8bSBarry Smith column of a matrix as index on an associated grid. 3292ef4de8bSBarry Smith 3302ef4de8bSBarry Smith Level: beginner 3312ef4de8bSBarry Smith 3322ef4de8bSBarry Smith Concepts: matrix; linear operator 3332ef4de8bSBarry Smith 334d7bd93baSBarry Smith .seealso: MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockStencil() 3352ef4de8bSBarry Smith S*/ 336435da068SBarry Smith typedef struct { 337c1ac3661SBarry Smith PetscInt k,j,i,c; 338435da068SBarry Smith } MatStencil; 3392ef4de8bSBarry Smith 340be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode); 341be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode); 342be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt); 343435da068SBarry Smith 344be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetColoring(Mat,ISColoring); 345be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesAdic(Mat,void*); 346be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesAdifor(Mat,PetscInt,void*); 3473a7fca6bSBarry Smith 348d91e6319SBarry Smith /*E 349d91e6319SBarry Smith MatAssemblyType - Indicates if the matrix is now to be used, or if you plan 350d91e6319SBarry Smith to continue to add values to it 351d91e6319SBarry Smith 352d91e6319SBarry Smith Level: beginner 353d91e6319SBarry Smith 354d91e6319SBarry Smith .seealso: MatAssemblyBegin(), MatAssemblyEnd() 355d91e6319SBarry Smith E*/ 3566d4a8577SBarry Smith typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType; 357be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssemblyBegin(Mat,MatAssemblyType); 358be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssemblyEnd(Mat,MatAssemblyType); 359be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssembled(Mat,PetscTruth*); 3604f9c727eSBarry Smith 3611d73ed98SBarry Smith 36230de9b25SBarry Smith 363d91e6319SBarry Smith /*E 364d91e6319SBarry Smith MatOption - Options that may be set for a matrix and its behavior or storage 365d91e6319SBarry Smith 366d91e6319SBarry Smith Level: beginner 367d91e6319SBarry Smith 3680a835dfdSSatish Balay Any additions/changes here MUST also be made in include/finclude/petscmat.h 369d91e6319SBarry Smith 370d91e6319SBarry Smith .seealso: MatSetOption() 371d91e6319SBarry Smith E*/ 372512a5fc5SBarry Smith typedef enum {MAT_ROW_ORIENTED,MAT_NEW_NONZERO_LOCATIONS, 3734e0d8c25SBarry Smith MAT_SYMMETRIC, 3744e0d8c25SBarry Smith MAT_STRUCTURALLY_SYMMETRIC, 3758047e77bSBarry Smith MAT_NEW_DIAGONALS,MAT_IGNORE_OFF_PROC_ENTRIES, 3764e0d8c25SBarry Smith MAT_NEW_NONZERO_LOCATION_ERR, 3774e0d8c25SBarry Smith MAT_NEW_NONZERO_ALLOCATION_ERR,MAT_USE_HASH_TABLE, 3784e0d8c25SBarry Smith MAT_KEEP_ZEROED_ROWS,MAT_IGNORE_ZERO_ENTRIES, 3794e0d8c25SBarry Smith MAT_USE_INODES, 3804e0d8c25SBarry Smith MAT_HERMITIAN, 3814e0d8c25SBarry Smith MAT_SYMMETRY_ETERNAL, 3824e0d8c25SBarry Smith MAT_USE_COMPRESSEDROW, 3834e0d8c25SBarry Smith MAT_IGNORE_LOWER_TRIANGULAR,MAT_ERROR_LOWER_TRIANGULAR, 38428b2fa4aSMatthew Knepley MAT_GETROW_UPPERTRIANGULAR,MAT_UNUSED_NONZERO_LOCATION_ERR, 38528b2fa4aSMatthew Knepley NUM_MAT_OPTIONS} MatOption; 386290bbb0aSBarry Smith extern const char *MatOptions[]; 3874e0d8c25SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetOption(Mat,MatOption,PetscTruth); 388a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetType(Mat,const MatType*); 389a313700dSBarry Smith PetscPolymorphicFunction(MatGetType,(Mat mat),(mat,&t),const MatType,t) 39084cb2905SBarry Smith 391be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]); 392be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 393be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 394f5edf698SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowUpperTriangular(Mat); 395f5edf698SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRowUpperTriangular(Mat); 396be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 397be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 398be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumnVector(Mat,Vec,PetscInt); 399be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetArray(Mat,PetscScalar *[]); 400ba966639SSatish Balay PetscPolymorphicFunction(MatGetArray,(Mat mat),(mat,&a),PetscScalar*,a) 401be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreArray(Mat,PetscScalar *[]); 402be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBlockSize(Mat,PetscInt *); 403ba966639SSatish Balay PetscPolymorphicFunction(MatGetBlockSize,(Mat mat),(mat,&a),PetscInt,a) 404be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetBlockSize(Mat,PetscInt); 4057b80b807SBarry Smith 4061620fd73SBarry Smith 407be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMult(Mat,Vec,Vec); 408be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultAdd(Mat,Vec,Vec,Vec); 409ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultAdd,(Mat A,Vec x,Vec y),(A,x,y,y)) 410be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTranspose(Mat,Vec,Vec); 411be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsTranspose(Mat,Mat,PetscReal,PetscTruth*); 412ba966639SSatish Balay PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B,PetscReal tol),(A,B,tol,&t),PetscTruth,t) 413e37f64cdSBarry Smith PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B),(A,B,0,&t),PetscTruth,t) 4141cbb95d3SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitianTranspose(Mat,Mat,PetscReal,PetscTruth*); 415be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeAdd(Mat,Vec,Vec,Vec); 416ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultTransposeAdd,(Mat A,Vec x,Vec y),(A,x,y,y)) 417be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultConstrained(Mat,Vec,Vec); 418be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeConstrained(Mat,Vec,Vec); 4192bebee5dSHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatSolve(Mat,Mat,Mat); 420f5edf698SHong Zhang 421d91e6319SBarry Smith /*E 422d91e6319SBarry Smith MatDuplicateOption - Indicates if a duplicated sparse matrix should have 423d91e6319SBarry Smith its numerical values copied over or just its nonzero structure. 424d91e6319SBarry Smith 425d91e6319SBarry Smith Level: beginner 426d91e6319SBarry Smith 427d91e6319SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 428d91e6319SBarry Smith 429d91e6319SBarry Smith .seealso: MatDuplicate() 430d91e6319SBarry Smith E*/ 4312e8a6d31SBarry Smith typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES} MatDuplicateOption; 4322e8a6d31SBarry Smith 433a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert(Mat,const MatType,MatReuse,Mat*); 434a313700dSBarry Smith PetscPolymorphicFunction(MatConvert,(Mat A,const MatType t),(A,t,MAT_INITIAL_MATRIX,&a),Mat,a) 435be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDuplicate(Mat,MatDuplicateOption,Mat*); 436ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A,MatDuplicateOption o),(A,o,&a),Mat,a) 437ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A),(A,MAT_COPY_VALUES,&a),Mat,a) 43894a9d846SBarry Smith 439cb5b572fSBarry Smith 440be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCopy(Mat,Mat,MatStructure); 441be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatView(Mat,PetscViewer); 442be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsSymmetric(Mat,PetscReal,PetscTruth*); 443ba966639SSatish Balay PetscPolymorphicFunction(MatIsSymmetric,(Mat A,PetscReal tol),(A,tol,&t),PetscTruth,t) 444e37f64cdSBarry Smith PetscPolymorphicFunction(MatIsSymmetric,(Mat A),(A,0,&t),PetscTruth,t) 445be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsStructurallySymmetric(Mat,PetscTruth*); 446ba966639SSatish Balay PetscPolymorphicFunction(MatIsStructurallySymmetric,(Mat A),(A,&t),PetscTruth,t) 4471cbb95d3SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitian(Mat,PetscReal,PetscTruth*); 4481cbb95d3SBarry Smith PetscPolymorphicFunction(MatIsHermitian,(Mat A),(A,0,&t),PetscTruth,t) 449be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsSymmetricKnown(Mat,PetscTruth*,PetscTruth*); 450be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitianKnown(Mat,PetscTruth*,PetscTruth*); 45164c62002SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMissingDiagonal(Mat,PetscTruth *,PetscInt *); 452a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLoad(PetscViewer,const MatType,Mat*); 453a313700dSBarry Smith PetscPolymorphicFunction(MatLoad,(PetscViewer v,const MatType t),(v,t,&a),Mat,a) 4547b80b807SBarry Smith 4558f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *); 4568f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRowIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *); 4578f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumnIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *); 4588f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreColumnIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *); 459d4fbbf0eSBarry Smith 460d91e6319SBarry Smith /*S 461d91e6319SBarry Smith MatInfo - Context of matrix information, used with MatGetInfo() 462d91e6319SBarry Smith 463d91e6319SBarry Smith In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE 464d91e6319SBarry Smith 465d91e6319SBarry Smith Level: intermediate 466d91e6319SBarry Smith 467d91e6319SBarry Smith Concepts: matrix^nonzero information 468d91e6319SBarry Smith 469d9274352SBarry Smith .seealso: MatGetInfo(), MatInfoType 470d91e6319SBarry Smith S*/ 4714e220ebcSLois Curfman McInnes typedef struct { 472b0a32e0cSBarry Smith PetscLogDouble block_size; /* block size */ 473b0a32e0cSBarry Smith PetscLogDouble nz_allocated,nz_used,nz_unneeded; /* number of nonzeros */ 474b0a32e0cSBarry Smith PetscLogDouble memory; /* memory allocated */ 475b0a32e0cSBarry Smith PetscLogDouble assemblies; /* number of matrix assemblies called */ 476b0a32e0cSBarry Smith PetscLogDouble mallocs; /* number of mallocs during MatSetValues() */ 477b0a32e0cSBarry Smith PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */ 478b0a32e0cSBarry Smith PetscLogDouble factor_mallocs; /* number of mallocs during factorization */ 4794e220ebcSLois Curfman McInnes } MatInfo; 4804e220ebcSLois Curfman McInnes 481d9274352SBarry Smith /*E 482d9274352SBarry Smith MatInfoType - Indicates if you want information about the local part of the matrix, 483d9274352SBarry Smith the entire parallel matrix or the maximum over all the local parts. 484d9274352SBarry Smith 485d9274352SBarry Smith Level: beginner 486d9274352SBarry Smith 487d9274352SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 488d9274352SBarry Smith 489d9274352SBarry Smith .seealso: MatGetInfo(), MatInfo 490d9274352SBarry Smith E*/ 4917b80b807SBarry Smith typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType; 492be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetInfo(Mat,MatInfoType,MatInfo*); 493be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatValid(Mat,PetscTruth*); 494be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetDiagonal(Mat,Vec); 495985db425SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMax(Mat,Vec,PetscInt[]); 496985db425SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMin(Mat,Vec,PetscInt[]); 497985db425SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMaxAbs(Mat,Vec,PetscInt[]); 498c87e5d42SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMinAbs(Mat,Vec,PetscInt[]); 49986697f06SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowSum(Mat,Vec); 500fc4dec0aSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatTranspose(Mat,MatReuse,Mat*); 501fc4dec0aSBarry Smith PetscPolymorphicFunction(MatTranspose,(Mat A),(A,MAT_INITIAL_MATRIX,&t),Mat,t) 502be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPermute(Mat,IS,IS,Mat *); 503ba966639SSatish Balay PetscPolymorphicFunction(MatPermute,(Mat A,IS is1,IS is2),(A,is1,is2,&t),Mat,t) 504be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPermuteSparsify(Mat,PetscInt,PetscReal,PetscReal,IS,IS,Mat *); 505be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalScale(Mat,Vec,Vec); 506be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalSet(Mat,Vec,InsertMode); 507be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatEqual(Mat,Mat,PetscTruth*); 508ba966639SSatish Balay PetscPolymorphicFunction(MatEqual,(Mat A,Mat B),(A,B,&t),PetscTruth,t) 509be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultEqual(Mat,Mat,PetscInt,PetscTruth*); 510be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultAddEqual(Mat,Mat,PetscInt,PetscTruth*); 511be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeEqual(Mat,Mat,PetscInt,PetscTruth*); 512be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscTruth*); 5137b80b807SBarry Smith 514be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNorm(Mat,NormType,PetscReal *); 515ba966639SSatish Balay PetscPolymorphicFunction(MatNorm,(Mat A,NormType t),(A,t,&n),PetscReal,n) 516be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroEntries(Mat); 517f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar); 518f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsIS(Mat,IS,PetscScalar); 519f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroColumns(Mat,PetscInt,const PetscInt [],const PetscScalar*); 520f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroColumnsIS(Mat,IS,const PetscScalar*); 5217b80b807SBarry Smith 522be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatUseScaledForm(Mat,PetscTruth); 523be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScaleSystem(Mat,Vec,Vec); 524be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatUnScaleSystem(Mat,Vec,Vec); 5255ef9f2a5SBarry Smith 526be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSize(Mat,PetscInt*,PetscInt*); 527be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalSize(Mat,PetscInt*,PetscInt*); 528be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRange(Mat,PetscInt*,PetscInt*); 5298ee2e534SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRanges(Mat,const PetscInt**); 530f4b8409dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRangeColumn(Mat,PetscInt*,PetscInt*); 53172ca8751SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRangesColumn(Mat,const PetscInt**); 5327b80b807SBarry Smith 533be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]); 534be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroyMatrices(PetscInt,Mat *[]); 535be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrix(Mat,IS,IS,PetscInt,MatReuse,Mat *); 536abf0e977SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrixRaw(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscInt,MatReuse,Mat *); 537829201f2SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSeqNonzeroStructure(Mat,Mat *[]); 538829201f2SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroySeqNonzeroStructure(Mat *[]); 5395494a064SHong Zhang 540be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge(MPI_Comm,Mat,PetscInt,MatReuse,Mat*); 541be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPI(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*); 542be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPISymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*); 543be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPINumeric(Mat,Mat); 544be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroy_MPIAIJ_SeqsToMPI(Mat); 545be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalMat(Mat,MatReuse,Mat*); 546be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*); 547be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,PetscInt*,Mat*); 548dd6ea824SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBrowsOfAoCols(Mat,Mat,MatReuse,PetscInt**,MatScalar**,Mat*); 54943eb5e2fSMatthew Knepley #if defined (PETSC_USE_CTABLE) 550cd116e26SSatish Balay #include "petscctable.h" 55143eb5e2fSMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetCommunicationStructs(Mat, Vec *, PetscTable *, VecScatter *); 55243eb5e2fSMatthew Knepley #else 55343eb5e2fSMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetCommunicationStructs(Mat, Vec *, PetscInt *[], VecScatter *); 55443eb5e2fSMatthew Knepley #endif 5558efafbd8SBarry Smith 556be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt); 5577b80b807SBarry Smith 558be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*); 559be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultSymbolic(Mat,Mat,PetscReal,Mat*); 560be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultNumeric(Mat,Mat,Mat); 56122440eb1SKris Buschelman 562be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*); 563be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAPSymbolic(Mat,Mat,PetscReal,Mat*); 564be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAPNumeric(Mat,Mat,Mat); 56522440eb1SKris Buschelman 566be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTranspose(Mat,Mat,MatReuse,PetscReal,Mat*); 567be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTransposeSymbolic(Mat,Mat,PetscReal,Mat*); 568be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTransposeNumeric(Mat,Mat,Mat); 569bc011b1eSHong Zhang 570f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAXPY(Mat,PetscScalar,Mat,MatStructure); 57104aac2b0SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAYPX(Mat,PetscScalar,Mat,MatStructure); 572be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompress(Mat); 5731c741599SHong Zhang 574f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScale(Mat,PetscScalar); 575f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShift(Mat,PetscScalar); 5767b80b807SBarry Smith 577be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping); 578be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetLocalToGlobalMappingBlock(Mat,ISLocalToGlobalMapping); 579f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar); 580f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsLocalIS(Mat,IS,PetscScalar); 581be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 582be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 583052efed2SBarry Smith 584be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStashSetInitialSize(Mat,PetscInt,PetscInt); 585be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*); 58690f02eecSBarry Smith 587be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInterpolate(Mat,Vec,Vec); 5880c451bc4SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInterpolateAdd(Mat,Vec,Vec,Vec); 589ba966639SSatish Balay PetscPolymorphicSubroutine(MatInterpolateAdd,(Mat A,Vec x,Vec y),(A,x,y,y)) 590be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestrict(Mat,Vec,Vec); 591be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetVecs(Mat,Vec*,Vec*); 59269db28dcSHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRedundantMatrix(Mat,PetscInt,MPI_Comm,PetscInt,MatReuse,Mat*); 593bd481603SBarry Smith 594bd481603SBarry Smith /*MC 5951620fd73SBarry Smith MatSetValue - Set a single entry into a matrix. 5961620fd73SBarry Smith 5971620fd73SBarry Smith Not collective 5981620fd73SBarry Smith 5991620fd73SBarry Smith Input Parameters: 6001620fd73SBarry Smith + m - the matrix 6011620fd73SBarry Smith . row - the row location of the entry 6021620fd73SBarry Smith . col - the column location of the entry 6031620fd73SBarry Smith . value - the value to insert 6041620fd73SBarry Smith - mode - either INSERT_VALUES or ADD_VALUES 6051620fd73SBarry Smith 6061620fd73SBarry Smith Notes: 6071620fd73SBarry Smith For efficiency one should use MatSetValues() and set several or many 6081620fd73SBarry Smith values simultaneously if possible. 6091620fd73SBarry Smith 6101620fd73SBarry Smith Level: beginner 6111620fd73SBarry Smith 6121620fd73SBarry Smith .seealso: MatSetValues(), MatSetValueLocal() 6131620fd73SBarry Smith M*/ 6141620fd73SBarry 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);} 6151620fd73SBarry Smith 6161620fd73SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatGetValue(Mat v,PetscInt i,PetscInt j,PetscScalar *va) {return MatGetValues(v,1,&i,1,&j,va);} 6171620fd73SBarry Smith 6181620fd73SBarry 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);} 6191620fd73SBarry Smith 6201620fd73SBarry Smith /*MC 6210d2a7ff2SSatish Balay MatPreallocateInitialize - Begins the block of code that will count the number of nonzeros per 622bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 623bd481603SBarry Smith 624bd481603SBarry Smith Synopsis: 625c1ac3661SBarry Smith PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz) 626bd481603SBarry Smith 627bd481603SBarry Smith Collective on MPI_Comm 628bd481603SBarry Smith 629bd481603SBarry Smith Input Parameters: 630bd481603SBarry Smith + comm - the communicator that will share the eventually allocated matrix 631859fcb39SBarry Smith . nrows - the number of LOCAL rows in the matrix 632859fcb39SBarry Smith - ncols - the number of LOCAL columns in the matrix 633bd481603SBarry Smith 634bd481603SBarry Smith Output Parameters: 635bd481603SBarry Smith + dnz - the array that will be passed to the matrix preallocation routines 636bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 637bd481603SBarry Smith 638bd481603SBarry Smith 639bd481603SBarry Smith Level: intermediate 640bd481603SBarry Smith 641bd481603SBarry Smith Notes: 642bd481603SBarry Smith See the chapter in the users manual on performance for more details 643bd481603SBarry Smith 6441d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 645bd481603SBarry Smith 646bd481603SBarry Smith Use MatPreallocateInitializeSymmetric() for symmetric matrices (MPISBAIJ matrices) 647bd481603SBarry Smith 6481620fd73SBarry Smith This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize(). 6491620fd73SBarry Smith 650bd481603SBarry Smith Concepts: preallocation^Matrix 651bd481603SBarry Smith 652bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(), 653bd481603SBarry Smith MatPreallocateInitializeSymmetric(), MatPreallocateSymmetricSetLocal() 654bd481603SBarry Smith M*/ 655c4f061fbSSatish Balay #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \ 6567c922b88SBarry Smith { \ 657a89bc211SBarry Smith PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__start,__end; \ 658a89bc211SBarry Smith _4_ierr = PetscMalloc2(__nrows,PetscInt,&dnz,__nrows,PetscInt,&onz);CHKERRQ(_4_ierr); \ 659a89bc211SBarry Smith _4_ierr = PetscMemzero(dnz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\ 660a89bc211SBarry Smith _4_ierr = PetscMemzero(onz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\ 661c1ac3661SBarry Smith _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ctmp;\ 662a89bc211SBarry Smith _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows; 6637c922b88SBarry Smith 664bd481603SBarry Smith /*MC 6650d2a7ff2SSatish Balay MatPreallocateSymmetricInitialize - Begins the block of code that will count the number of nonzeros per 666bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 667bd481603SBarry Smith 668bd481603SBarry Smith Synopsis: 669c1ac3661SBarry Smith PetscErrorCode MatPreallocateSymmetricInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz) 670bd481603SBarry Smith 671bd481603SBarry Smith Collective on MPI_Comm 672bd481603SBarry Smith 673bd481603SBarry Smith Input Parameters: 674bd481603SBarry Smith + comm - the communicator that will share the eventually allocated matrix 675859fcb39SBarry Smith . nrows - the number of LOCAL rows in the matrix 676859fcb39SBarry Smith - ncols - the number of LOCAL columns in the matrix 677bd481603SBarry Smith 678bd481603SBarry Smith Output Parameters: 679bd481603SBarry Smith + dnz - the array that will be passed to the matrix preallocation routines 680bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 681bd481603SBarry Smith 682bd481603SBarry Smith 683bd481603SBarry Smith Level: intermediate 684bd481603SBarry Smith 685bd481603SBarry Smith Notes: 686bd481603SBarry Smith See the chapter in the users manual on performance for more details 687bd481603SBarry Smith 6881d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 689bd481603SBarry Smith 6901620fd73SBarry Smith This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize(). 6911620fd73SBarry Smith 692bd481603SBarry Smith Concepts: preallocation^Matrix 693bd481603SBarry Smith 694bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(), 695bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal() 696bd481603SBarry Smith M*/ 697222b16d4SBarry Smith #define MatPreallocateSymmetricInitialize(comm,nrows,ncols,dnz,onz) 0; \ 698222b16d4SBarry Smith { \ 699a89bc211SBarry Smith PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__end; \ 700a89bc211SBarry Smith _4_ierr = PetscMalloc2(__nrows,PetscInt,&dnz,__nrows,PetscInt,&onz);CHKERRQ(_4_ierr); \ 701a89bc211SBarry Smith _4_ierr = PetscMemzero(dnz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\ 702a89bc211SBarry Smith _4_ierr = PetscMemzero(onz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\ 703c1ac3661SBarry Smith _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr);\ 704a89bc211SBarry Smith _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows; 705222b16d4SBarry Smith 706bd481603SBarry Smith /*MC 707bd481603SBarry Smith MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be 708bd481603SBarry Smith inserted using a local number of the rows and columns 709bd481603SBarry Smith 710bd481603SBarry Smith Synopsis: 711c1ac3661SBarry Smith PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 712bd481603SBarry Smith 713bd481603SBarry Smith Not Collective 714bd481603SBarry Smith 715bd481603SBarry Smith Input Parameters: 716bd481603SBarry Smith + map - the mapping between local numbering and global numbering 717bd481603SBarry Smith . nrows - the number of rows indicated 7181d73ed98SBarry Smith . rows - the indices of the rows 719bd481603SBarry Smith . ncols - the number of columns in the matrix 720bd481603SBarry Smith . cols - the columns indicated 721bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 722bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 723bd481603SBarry Smith 724bd481603SBarry Smith 725bd481603SBarry Smith Level: intermediate 726bd481603SBarry Smith 727bd481603SBarry Smith Notes: 728bd481603SBarry Smith See the chapter in the users manual on performance for more details 729bd481603SBarry Smith 7301d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 731bd481603SBarry Smith 732bd481603SBarry Smith Concepts: preallocation^Matrix 733bd481603SBarry Smith 734bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(), 735bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal() 736bd481603SBarry Smith M*/ 737c4f061fbSSatish Balay #define MatPreallocateSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\ 738c4f061fbSSatish Balay {\ 739c1ac3661SBarry Smith PetscInt __l;\ 740ef66eb69SBarry Smith _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\ 741ef66eb69SBarry Smith _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\ 742c4f061fbSSatish Balay for (__l=0;__l<nrows;__l++) {\ 743ef66eb69SBarry Smith _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 744c4f061fbSSatish Balay }\ 745c4f061fbSSatish Balay } 746c4f061fbSSatish Balay 747bd481603SBarry Smith /*MC 748bd481603SBarry Smith MatPreallocateSymmetricSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be 749bd481603SBarry Smith inserted using a local number of the rows and columns 750bd481603SBarry Smith 751bd481603SBarry Smith Synopsis: 752c1ac3661SBarry Smith PetscErrorCode MatPreallocateSymmetricSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 753bd481603SBarry Smith 754bd481603SBarry Smith Not Collective 755bd481603SBarry Smith 756bd481603SBarry Smith Input Parameters: 757bd481603SBarry Smith + map - the mapping between local numbering and global numbering 758bd481603SBarry Smith . nrows - the number of rows indicated 7591d73ed98SBarry Smith . rows - the indices of the rows 760bd481603SBarry Smith . ncols - the number of columns in the matrix 761bd481603SBarry Smith . cols - the columns indicated 762bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 763bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 764bd481603SBarry Smith 765bd481603SBarry Smith 766bd481603SBarry Smith Level: intermediate 767bd481603SBarry Smith 768bd481603SBarry Smith Notes: 769bd481603SBarry Smith See the chapter in the users manual on performance for more details 770bd481603SBarry Smith 771bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 772bd481603SBarry Smith 773bd481603SBarry Smith Concepts: preallocation^Matrix 774bd481603SBarry Smith 775bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(), 776bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal() 777bd481603SBarry Smith M*/ 778d3d32019SBarry Smith #define MatPreallocateSymmetricSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\ 779d3d32019SBarry Smith {\ 780c1ac3661SBarry Smith PetscInt __l;\ 781d3d32019SBarry Smith _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\ 782d3d32019SBarry Smith _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\ 783d3d32019SBarry Smith for (__l=0;__l<nrows;__l++) {\ 784d3d32019SBarry Smith _4_ierr = MatPreallocateSymmetricSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 785d3d32019SBarry Smith }\ 786d3d32019SBarry Smith } 787d3d32019SBarry Smith 788bd481603SBarry Smith /*MC 789bd481603SBarry Smith MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be 790bd481603SBarry Smith inserted using a local number of the rows and columns 791bd481603SBarry Smith 792bd481603SBarry Smith Synopsis: 793c1ac3661SBarry Smith PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 794bd481603SBarry Smith 795bd481603SBarry Smith Not Collective 796bd481603SBarry Smith 797bd481603SBarry Smith Input Parameters: 79864075487SBarry Smith + row - the row 799bd481603SBarry Smith . ncols - the number of columns in the matrix 800a50f8bf6SHong Zhang - cols - the columns indicated 801a50f8bf6SHong Zhang 802a50f8bf6SHong Zhang Output Parameters: 803a50f8bf6SHong Zhang + dnz - the array that will be passed to the matrix preallocation routines 804bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 805bd481603SBarry Smith 806bd481603SBarry Smith 807bd481603SBarry Smith Level: intermediate 808bd481603SBarry Smith 809bd481603SBarry Smith Notes: 810bd481603SBarry Smith See the chapter in the users manual on performance for more details 811bd481603SBarry Smith 812bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 813bd481603SBarry Smith 8141620fd73SBarry Smith This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize(). 8151620fd73SBarry Smith 816bd481603SBarry Smith Concepts: preallocation^Matrix 817bd481603SBarry Smith 818bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(), 819bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal() 820bd481603SBarry Smith M*/ 821c4f061fbSSatish Balay #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\ 822c1ac3661SBarry Smith { PetscInt __i; \ 8238ee2e534SBarry Smith if (row < __rstart) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Trying to set preallocation for row %D less than first local row %D",row,__rstart);\ 824a89bc211SBarry Smith if (row >= __rstart+__nrows) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Trying to set preallocation for row %D greater than last local row %D",row,__rstart+__nrows-1);\ 8257c922b88SBarry Smith for (__i=0; __i<nc; __i++) {\ 82664075487SBarry Smith if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \ 8277cc688d7SBarry Smith else dnz[row - __rstart]++;\ 8287c922b88SBarry Smith }\ 8297c922b88SBarry Smith } 8307c922b88SBarry Smith 831bd481603SBarry Smith /*MC 832bd481603SBarry Smith MatPreallocateSymmetricSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be 833bd481603SBarry Smith inserted using a local number of the rows and columns 834bd481603SBarry Smith 835bd481603SBarry Smith Synopsis: 836c1ac3661SBarry Smith PetscErrorCode MatPreallocateSymmetricSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 837bd481603SBarry Smith 838bd481603SBarry Smith Not Collective 839bd481603SBarry Smith 840bd481603SBarry Smith Input Parameters: 841bd481603SBarry Smith + nrows - the number of rows indicated 8421d73ed98SBarry Smith . rows - the indices of the rows 843bd481603SBarry Smith . ncols - the number of columns in the matrix 844bd481603SBarry Smith . cols - the columns indicated 845bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 846bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 847bd481603SBarry Smith 848bd481603SBarry Smith 849bd481603SBarry Smith Level: intermediate 850bd481603SBarry Smith 851bd481603SBarry Smith Notes: 852bd481603SBarry Smith See the chapter in the users manual on performance for more details 853bd481603SBarry Smith 854bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 855bd481603SBarry Smith 8561620fd73SBarry Smith This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize(). 8571620fd73SBarry Smith 858bd481603SBarry Smith Concepts: preallocation^Matrix 859bd481603SBarry Smith 860bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(), 861bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal() 862bd481603SBarry Smith M*/ 863d3d32019SBarry Smith #define MatPreallocateSymmetricSet(row,nc,cols,dnz,onz) 0;\ 864c1ac3661SBarry Smith { PetscInt __i; \ 865d3d32019SBarry Smith for (__i=0; __i<nc; __i++) {\ 866d3d32019SBarry Smith if (cols[__i] >= __end) onz[row - __rstart]++; \ 867d3d32019SBarry Smith else if (cols[__i] >= row) dnz[row - __rstart]++;\ 868d3d32019SBarry Smith }\ 869d3d32019SBarry Smith } 870d3d32019SBarry Smith 871bd481603SBarry Smith /*MC 87216371a99SBarry Smith MatPreallocateLocation - An alternative to MatPreallocationSet() that puts the nonzero locations into the matrix if it exists 87316371a99SBarry Smith 87416371a99SBarry Smith Synopsis: 87516371a99SBarry Smith PetscErrorCode MatPreallocateLocations(Mat A,PetscInt row,PetscInt ncols,PetscInt *cols,PetscInt *dnz,PetscInt *onz) 87616371a99SBarry Smith 87716371a99SBarry Smith Not Collective 87816371a99SBarry Smith 87916371a99SBarry Smith Input Parameters: 88016371a99SBarry Smith . A - matrix 88116371a99SBarry Smith . row - row where values exist (must be local to this process) 88216371a99SBarry Smith . ncols - number of columns 88316371a99SBarry Smith . cols - columns with nonzeros 88416371a99SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 88516371a99SBarry Smith - ozn - the other array passed to the matrix preallocation routines 88616371a99SBarry Smith 88716371a99SBarry Smith 88816371a99SBarry Smith Level: intermediate 88916371a99SBarry Smith 89016371a99SBarry Smith Notes: 89116371a99SBarry Smith See the chapter in the users manual on performance for more details 89216371a99SBarry Smith 89316371a99SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 89416371a99SBarry Smith 89516371a99SBarry Smith This is a MACRO not a function because it uses a bunch of variables private to the MatPreallocation.... routines. 89616371a99SBarry Smith 89716371a99SBarry Smith Concepts: preallocation^Matrix 89816371a99SBarry Smith 89916371a99SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(), 90016371a99SBarry Smith MatPreallocateSymmetricInitialize(), MatPreallocateSymmetricSetLocal() 90116371a99SBarry Smith M*/ 90216371a99SBarry 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);} 90316371a99SBarry Smith 90416371a99SBarry Smith 90516371a99SBarry Smith /*MC 9060d2a7ff2SSatish Balay MatPreallocateFinalize - Ends the block of code that will count the number of nonzeros per 907bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 908bd481603SBarry Smith 909bd481603SBarry Smith Synopsis: 910c1ac3661SBarry Smith PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz) 911bd481603SBarry Smith 912bd481603SBarry Smith Collective on MPI_Comm 913bd481603SBarry Smith 914bd481603SBarry Smith Input Parameters: 91516371a99SBarry Smith + dnz - the array that was be passed to the matrix preallocation routines 916bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 917bd481603SBarry Smith 918bd481603SBarry Smith 919bd481603SBarry Smith Level: intermediate 920bd481603SBarry Smith 921bd481603SBarry Smith Notes: 922bd481603SBarry Smith See the chapter in the users manual on performance for more details 923bd481603SBarry Smith 924bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 925bd481603SBarry Smith 9261620fd73SBarry Smith This is a MACRO not a function because it closes the { started in MatPreallocateInitialize(). 9271620fd73SBarry Smith 928bd481603SBarry Smith Concepts: preallocation^Matrix 929bd481603SBarry Smith 930bd481603SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(), 931bd481603SBarry Smith MatPreallocateSymmetricInitialize(), MatPreallocateSymmetricSetLocal() 932bd481603SBarry Smith M*/ 933a89bc211SBarry Smith #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree2(dnz,onz);CHKERRQ(_4_ierr);} 9347c922b88SBarry Smith 935bd481603SBarry Smith 936bd481603SBarry Smith 9377b80b807SBarry Smith /* Routines unique to particular data structures */ 938be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellGetContext(Mat,void **); 939ba966639SSatish Balay PetscPolymorphicFunction(MatShellGetContext,(Mat A),(A,&t),void*,t) 940698d4c6aSKris Buschelman 941be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInodeAdjustForInodes(Mat,IS*,IS*); 942be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *); 943698d4c6aSKris Buschelman 944be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetColumnIndices(Mat,PetscInt[]); 945be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetColumnIndices(Mat,PetscInt[]); 946be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 947c75a6043SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 948c75a6043SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 9497b80b807SBarry Smith 950a96a251dSBarry Smith #define MAT_SKIP_ALLOCATION -4 951a96a251dSBarry Smith 952be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]); 953ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz)) 954be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]); 955ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqSBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz)) 956be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]); 957ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqAIJSetPreallocation,(Mat A,const PetscInt nnz[]),(A,0,nnz)) 958273d9f13SBarry Smith 959be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 960ba966639SSatish Balay PetscPolymorphicSubroutine(MatMPIBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[],const PetscInt onz[]),(A,bs,0,nnz,0,onz)) 961be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 962be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 96353803f26SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []); 964725b52f3SLisandro Dalcin EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 965be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]); 966be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 967be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]); 968be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIDenseSetPreallocation(Mat,PetscScalar[]); 969284134d9SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation(Mat,PetscScalar[]); 970be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIRowbsSetPreallocation(Mat,PetscInt,const PetscInt[]); 971be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,PetscInt*[]); 972be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,PetscInt*[]); 973be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAdicSetLocalFunction(Mat,void (*)(void)); 974273d9f13SBarry Smith 975be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetLDA(Mat,PetscInt); 976ab92ecdeSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDenseGetLocalMatrix(Mat,Mat*); 9771b807ce4Svictorle 978be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStoreValues(Mat); 979be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRetrieveValues(Mat); 9802e8a6d31SBarry Smith 981be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDAADSetCtx(Mat,void*); 9823a7fca6bSBarry Smith 9837b80b807SBarry Smith /* 9847b80b807SBarry Smith These routines are not usually accessed directly, rather solving is 98594b7f48cSBarry Smith done through the KSP and PC interfaces. 9867b80b807SBarry Smith */ 9877b80b807SBarry Smith 988d9274352SBarry Smith /*E 989d9274352SBarry Smith MatOrderingType - String with the name of a PETSc matrix ordering or the creation function 990d9274352SBarry Smith with an optional dynamic library name, for example 991d9274352SBarry Smith http://www.mcs.anl.gov/petsc/lib.a:orderingcreate() 992d9274352SBarry Smith 993d9274352SBarry Smith Level: beginner 994d9274352SBarry Smith 995e5a9bf91SBarry Smith Cannot use const because the PC objects manipulate the string 996e5a9bf91SBarry Smith 997d9274352SBarry Smith .seealso: MatGetOrdering() 998d9274352SBarry Smith E*/ 9993eaad2caSSatish Balay #define MatOrderingType char* 1000b12f92e5SBarry Smith #define MATORDERING_NATURAL "natural" 1001b12f92e5SBarry Smith #define MATORDERING_ND "nd" 1002b12f92e5SBarry Smith #define MATORDERING_1WD "1wd" 1003b12f92e5SBarry Smith #define MATORDERING_RCM "rcm" 1004b12f92e5SBarry Smith #define MATORDERING_QMD "qmd" 1005b12f92e5SBarry Smith #define MATORDERING_ROWLENGTH "rowlength" 100662152c8bSBarry Smith #define MATORDERING_DSC_ND "dsc_nd" 100762152c8bSBarry Smith #define MATORDERING_DSC_MMD "dsc_mmd" 100862152c8bSBarry Smith #define MATORDERING_DSC_MDF "dsc_mdf" 1009c06d978dSMatthew Knepley #define MATORDERING_CONSTRAINED "constrained" 1010c06d978dSMatthew Knepley #define MATORDERING_IDENTITY "identity" 1011c06d978dSMatthew Knepley #define MATORDERING_REVERSE "reverse" 1012b12f92e5SBarry Smith 1013f2416520SSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOrdering(Mat,const MatOrderingType,IS*,IS*); 1014be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOrderingList(PetscFList *list); 1015f2416520SSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,const MatOrderingType,IS*,IS*)); 101630de9b25SBarry Smith 101730de9b25SBarry Smith /*MC 101830de9b25SBarry Smith MatOrderingRegisterDynamic - Adds a new sparse matrix ordering to the 101930de9b25SBarry Smith matrix package. 102030de9b25SBarry Smith 102130de9b25SBarry Smith Synopsis: 1022c1ac3661SBarry Smith PetscErrorCode MatOrderingRegisterDynamic(char *name_ordering,char *path,char *name_create,PetscErrorCode (*routine_create)(MatOrdering)) 102330de9b25SBarry Smith 102430de9b25SBarry Smith Not Collective 102530de9b25SBarry Smith 102630de9b25SBarry Smith Input Parameters: 102730de9b25SBarry Smith + sname - name of ordering (for example MATORDERING_ND) 102830de9b25SBarry Smith . path - location of library where creation routine is 102930de9b25SBarry Smith . name - name of function that creates the ordering type,a string 103030de9b25SBarry Smith - function - function pointer that creates the ordering 103130de9b25SBarry Smith 103230de9b25SBarry Smith Level: developer 103330de9b25SBarry Smith 103430de9b25SBarry Smith If dynamic libraries are used, then the fourth input argument (function) 103530de9b25SBarry Smith is ignored. 103630de9b25SBarry Smith 103730de9b25SBarry Smith Sample usage: 103830de9b25SBarry Smith .vb 103930de9b25SBarry Smith MatOrderingRegisterDynamic("my_order",/home/username/my_lib/lib/libO/solaris/mylib.a, 104030de9b25SBarry Smith "MyOrder",MyOrder); 104130de9b25SBarry Smith .ve 104230de9b25SBarry Smith 104330de9b25SBarry Smith Then, your partitioner can be chosen with the procedural interface via 104430de9b25SBarry Smith $ MatOrderingSetType(part,"my_order) 104530de9b25SBarry Smith or at runtime via the option 10462401956bSBarry Smith $ -pc_factor_mat_ordering_type my_order 104730de9b25SBarry Smith 1048ab901514SBarry Smith ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values. 104930de9b25SBarry Smith 105030de9b25SBarry Smith .keywords: matrix, ordering, register 105130de9b25SBarry Smith 105230de9b25SBarry Smith .seealso: MatOrderingRegisterDestroy(), MatOrderingRegisterAll() 105330de9b25SBarry Smith M*/ 1054aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 1055f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,0) 1056b12f92e5SBarry Smith #else 1057f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,d) 1058b12f92e5SBarry Smith #endif 105930de9b25SBarry Smith 1060be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegisterDestroy(void); 1061be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegisterAll(const char[]); 10622bad1931SBarry Smith extern PetscTruth MatOrderingRegisterAllCalled; 1063b0a32e0cSBarry Smith extern PetscFList MatOrderingList; 1064d4fbbf0eSBarry Smith 1065be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS); 1066a2ce50c7SBarry Smith 1067d91e6319SBarry Smith /*S 10682401956bSBarry Smith MatFactorInfo - Data passed into the matrix factorization routines 1069b00f7748SHong Zhang 107061cad860SBarry Smith In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE, that is use 107161cad860SBarry Smith $ MatFactorInfo info(MAT_FACTORINFO_SIZE) 1072b00f7748SHong Zhang 107315e8a5b3SHong Zhang Notes: These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC. 1074b00f7748SHong Zhang 107561cad860SBarry Smith You can use MatFactorInfoInitialize() to set default values. 107661cad860SBarry Smith 1077b00f7748SHong Zhang Level: developer 1078b00f7748SHong Zhang 1079d7d82daaSBarry Smith .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(), 1080d7d82daaSBarry Smith MatFactorInfoInitialize() 1081b00f7748SHong Zhang 1082b00f7748SHong Zhang S*/ 1083b00f7748SHong Zhang typedef struct { 10840a29876aSHong Zhang PetscReal shiftnz; /* scaling of identity added to matrix to prevent zero pivots */ 1085fbf22428SSatish Balay PetscReal shiftpd; /* if true, shift until positive pivots */ 108615e8a5b3SHong Zhang PetscReal diagonal_fill; /* force diagonal to fill in if initially not filled */ 108785317021SBarry Smith PetscReal usedt; 108815e8a5b3SHong Zhang PetscReal dt; /* drop tolerance */ 1089b00f7748SHong Zhang PetscReal dtcol; /* tolerance for pivoting */ 109015e8a5b3SHong Zhang PetscReal dtcount; /* maximum nonzeros to be allowed per row */ 109167e28bfeSBarry Smith PetscReal fill; /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */ 1092348344bbSBarry Smith PetscReal levels; /* ICC/ILU(levels) */ 1093bcd9e38bSBarry Smith PetscReal pivotinblocks; /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0 1094bcd9e38bSBarry Smith factorization may be faster if do not pivot */ 109562bba022SBarry Smith PetscReal shiftinblocks; /* if block in block factorization has zero pivot then shift diagonal until non-singular */ 109615e8a5b3SHong Zhang PetscReal zeropivot; /* pivot is called zero if less than this */ 109715e8a5b3SHong Zhang } MatFactorInfo; 1098ffa6d0a5SLois Curfman McInnes 1099be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFactorInfoInitialize(MatFactorInfo*); 11000481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactor(Mat,IS,const MatFactorInfo*); 11010481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*); 11020481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactorNumeric(Mat,Mat,const MatFactorInfo*); 11030481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactor(Mat,IS,IS,const MatFactorInfo*); 11040481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUFactor(Mat,IS,IS,const MatFactorInfo*); 11050481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*); 11060481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*); 11070481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatICCFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*); 11080481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatICCFactor(Mat,IS,const MatFactorInfo*); 11090481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactorNumeric(Mat,Mat,const MatFactorInfo*); 11100481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUDTFactor(Mat,IS,IS,const MatFactorInfo*,Mat *); 1111be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*); 1112be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolve(Mat,Vec,Vec); 1113be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatForwardSolve(Mat,Vec,Vec); 1114be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatBackwardSolve(Mat,Vec,Vec); 1115be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveAdd(Mat,Vec,Vec,Vec); 1116be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveTranspose(Mat,Vec,Vec); 1117be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveTransposeAdd(Mat,Vec,Vec,Vec); 1118be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolves(Mat,Vecs,Vecs); 11198ed539a5SBarry Smith 1120be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUnfactored(Mat); 1121bb5a7306SBarry Smith 1122d91e6319SBarry Smith /*E 1123d91e6319SBarry Smith MatSORType - What type of (S)SOR to perform 1124bb1eb677SSatish Balay 1125d91e6319SBarry Smith Level: beginner 1126d91e6319SBarry Smith 1127d9274352SBarry Smith May be bitwise ORd together 1128d9274352SBarry Smith 1129d91e6319SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 1130d91e6319SBarry Smith 11314e7234bfSBarry Smith MatSORType may be bitwise ORd together, so do not change the numbers 11324e7234bfSBarry Smith 1133a1d92eedSBarry Smith .seealso: MatRelax(), MatPBRelax() 1134d91e6319SBarry Smith E*/ 1135ee50ffe9SBarry Smith typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3, 1136ee50ffe9SBarry Smith SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8, 1137ee50ffe9SBarry Smith SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16, 113884cb2905SBarry Smith SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType; 1139be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRelax(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec); 1140be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPBRelax(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec); 11418ed539a5SBarry Smith 1142d4fbbf0eSBarry Smith /* 1143639f9d9dSBarry Smith These routines are for efficiently computing Jacobians via finite differences. 1144639f9d9dSBarry Smith */ 1145b12f92e5SBarry Smith 1146d9274352SBarry Smith /*E 1147d9274352SBarry Smith MatColoringType - String with the name of a PETSc matrix coloring or the creation function 1148d9274352SBarry Smith with an optional dynamic library name, for example 1149d9274352SBarry Smith http://www.mcs.anl.gov/petsc/lib.a:coloringcreate() 1150d9274352SBarry Smith 1151d9274352SBarry Smith Level: beginner 1152d9274352SBarry Smith 1153d9274352SBarry Smith .seealso: MatGetColoring() 1154d9274352SBarry Smith E*/ 1155a313700dSBarry Smith #define MatColoringType char* 1156b12f92e5SBarry Smith #define MATCOLORING_NATURAL "natural" 1157b12f92e5SBarry Smith #define MATCOLORING_SL "sl" 1158b12f92e5SBarry Smith #define MATCOLORING_LF "lf" 1159b12f92e5SBarry Smith #define MATCOLORING_ID "id" 1160b12f92e5SBarry Smith 1161ab970d9bSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColoring(Mat,const MatColoringType,ISColoring*); 11622e90c967SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,MatColoringType,ISColoring *)); 116330de9b25SBarry Smith 116430de9b25SBarry Smith /*MC 116530de9b25SBarry Smith MatColoringRegisterDynamic - Adds a new sparse matrix coloring to the 116630de9b25SBarry Smith matrix package. 116730de9b25SBarry Smith 116830de9b25SBarry Smith Synopsis: 1169c1ac3661SBarry Smith PetscErrorCode MatColoringRegisterDynamic(char *name_coloring,char *path,char *name_create,PetscErrorCode (*routine_create)(MatColoring)) 117030de9b25SBarry Smith 117130de9b25SBarry Smith Not Collective 117230de9b25SBarry Smith 117330de9b25SBarry Smith Input Parameters: 117430de9b25SBarry Smith + sname - name of Coloring (for example MATCOLORING_SL) 117530de9b25SBarry Smith . path - location of library where creation routine is 117630de9b25SBarry Smith . name - name of function that creates the Coloring type, a string 117730de9b25SBarry Smith - function - function pointer that creates the coloring 117830de9b25SBarry Smith 117930de9b25SBarry Smith Level: developer 118030de9b25SBarry Smith 118130de9b25SBarry Smith If dynamic libraries are used, then the fourth input argument (function) 118230de9b25SBarry Smith is ignored. 118330de9b25SBarry Smith 118430de9b25SBarry Smith Sample usage: 118530de9b25SBarry Smith .vb 118630de9b25SBarry Smith MatColoringRegisterDynamic("my_color",/home/username/my_lib/lib/libO/solaris/mylib.a, 118730de9b25SBarry Smith "MyColor",MyColor); 118830de9b25SBarry Smith .ve 118930de9b25SBarry Smith 119030de9b25SBarry Smith Then, your partitioner can be chosen with the procedural interface via 119130de9b25SBarry Smith $ MatColoringSetType(part,"my_color") 119230de9b25SBarry Smith or at runtime via the option 119330de9b25SBarry Smith $ -mat_coloring_type my_color 119430de9b25SBarry Smith 1195ab901514SBarry Smith $PETSC_ARCH occuring in pathname will be replaced with appropriate values. 119630de9b25SBarry Smith 119730de9b25SBarry Smith .keywords: matrix, Coloring, register 119830de9b25SBarry Smith 119930de9b25SBarry Smith .seealso: MatColoringRegisterDestroy(), MatColoringRegisterAll() 120030de9b25SBarry Smith M*/ 1201aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 1202f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,0) 1203b12f92e5SBarry Smith #else 1204f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,d) 1205b12f92e5SBarry Smith #endif 120630de9b25SBarry Smith 12072bad1931SBarry Smith extern PetscTruth MatColoringRegisterAllCalled; 1208f1144a30SSatish Balay 1209be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegisterAll(const char[]); 1210be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegisterDestroy(void); 1211be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*); 1212639f9d9dSBarry Smith 1213d9274352SBarry Smith /*S 1214d9274352SBarry Smith MatFDColoring - Object for computing a sparse Jacobian via finite differences 1215d9274352SBarry Smith and coloring 1216639f9d9dSBarry Smith 1217d9274352SBarry Smith Level: beginner 1218d9274352SBarry Smith 1219d9274352SBarry Smith Concepts: coloring, sparse Jacobian, finite differences 1220d9274352SBarry Smith 1221d9274352SBarry Smith .seealso: MatFDColoringCreate() 1222d9274352SBarry Smith S*/ 1223e2a1c21fSSatish Balay typedef struct _p_MatFDColoring* MatFDColoring; 1224639f9d9dSBarry Smith 1225be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringCreate(Mat,ISColoring,MatFDColoring *); 1226be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringDestroy(MatFDColoring); 1227be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringView(MatFDColoring,PetscViewer); 1228be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*); 12294a9d489dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**); 1230be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal); 1231be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFromOptions(MatFDColoring); 1232be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApply(Mat,MatFDColoring,Vec,MatStructure*,void *); 1233be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApplyTS(Mat,MatFDColoring,PetscReal,Vec,MatStructure*,void *); 1234be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetF(MatFDColoring,Vec); 1235be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,PetscInt*[]); 1236639f9d9dSBarry Smith /* 12370752156aSBarry Smith These routines are for partitioning matrices: currently used only 12383eda8832SBarry Smith for adjacency matrix, MatCreateMPIAdj(). 12390752156aSBarry Smith */ 1240ca161407SBarry Smith 1241d9274352SBarry Smith /*S 1242d9274352SBarry Smith MatPartitioning - Object for managing the partitioning of a matrix or graph 1243d9274352SBarry Smith 1244d9274352SBarry Smith Level: beginner 1245d9274352SBarry Smith 1246d9274352SBarry Smith Concepts: partitioning 1247d9274352SBarry Smith 1248743a1026Svictorle .seealso: MatPartitioningCreate(), MatPartitioningType 1249d9274352SBarry Smith S*/ 125091e9ee9fSBarry Smith typedef struct _p_MatPartitioning* MatPartitioning; 1251d9274352SBarry Smith 1252d9274352SBarry Smith /*E 12535bc6e7ccSvictorle MatPartitioningType - String with the name of a PETSc matrix partitioning or the creation function 1254d9274352SBarry Smith with an optional dynamic library name, for example 1255d9274352SBarry Smith http://www.mcs.anl.gov/petsc/lib.a:partitioningcreate() 1256d9274352SBarry Smith 1257d9274352SBarry Smith Level: beginner 1258d9274352SBarry Smith 1259b547a13bSvictorle .seealso: MatPartitioningCreate(), MatPartitioning 1260d9274352SBarry Smith E*/ 1261a313700dSBarry Smith #define MatPartitioningType char* 12628ba1e511SMatthew Knepley #define MAT_PARTITIONING_CURRENT "current" 1263c2ec2a73SHong Zhang #define MAT_PARTITIONING_SQUARE "square" 12648ba1e511SMatthew Knepley #define MAT_PARTITIONING_PARMETIS "parmetis" 126571306c60Sjroman #define MAT_PARTITIONING_CHACO "chaco" 126671306c60Sjroman #define MAT_PARTITIONING_JOSTLE "jostle" 126771306c60Sjroman #define MAT_PARTITIONING_PARTY "party" 126871306c60Sjroman #define MAT_PARTITIONING_SCOTCH "scotch" 126971306c60Sjroman 1270ca161407SBarry Smith 1271be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningCreate(MPI_Comm,MatPartitioning*); 1272a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetType(MatPartitioning,const MatPartitioningType); 1273be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetNParts(MatPartitioning,PetscInt); 1274be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetAdjacency(MatPartitioning,Mat); 1275be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]); 1276be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []); 1277be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningApply(MatPartitioning,IS*); 1278be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningDestroy(MatPartitioning); 12792aabb6bbSBarry Smith 1280be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatPartitioning)); 128130de9b25SBarry Smith 128230de9b25SBarry Smith /*MC 128330de9b25SBarry Smith MatPartitioningRegisterDynamic - Adds a new sparse matrix partitioning to the 128430de9b25SBarry Smith matrix package. 128530de9b25SBarry Smith 128630de9b25SBarry Smith Synopsis: 1287c1ac3661SBarry Smith PetscErrorCode MatPartitioningRegisterDynamic(char *name_partitioning,char *path,char *name_create,PetscErrorCode (*routine_create)(MatPartitioning)) 128830de9b25SBarry Smith 128930de9b25SBarry Smith Not Collective 129030de9b25SBarry Smith 129130de9b25SBarry Smith Input Parameters: 129230de9b25SBarry Smith + sname - name of partitioning (for example MAT_PARTITIONING_CURRENT) or parmetis 129330de9b25SBarry Smith . path - location of library where creation routine is 129430de9b25SBarry Smith . name - name of function that creates the partitioning type, a string 129530de9b25SBarry Smith - function - function pointer that creates the partitioning type 129630de9b25SBarry Smith 129730de9b25SBarry Smith Level: developer 129830de9b25SBarry Smith 129930de9b25SBarry Smith If dynamic libraries are used, then the fourth input argument (function) 130030de9b25SBarry Smith is ignored. 130130de9b25SBarry Smith 130230de9b25SBarry Smith Sample usage: 130330de9b25SBarry Smith .vb 130430de9b25SBarry Smith MatPartitioningRegisterDynamic("my_part",/home/username/my_lib/lib/libO/solaris/mylib.a, 130530de9b25SBarry Smith "MyPartCreate",MyPartCreate); 130630de9b25SBarry Smith .ve 130730de9b25SBarry Smith 130830de9b25SBarry Smith Then, your partitioner can be chosen with the procedural interface via 130930de9b25SBarry Smith $ MatPartitioningSetType(part,"my_part") 131030de9b25SBarry Smith or at runtime via the option 131130de9b25SBarry Smith $ -mat_partitioning_type my_part 131230de9b25SBarry Smith 1313ab901514SBarry Smith $PETSC_ARCH occuring in pathname will be replaced with appropriate values. 131430de9b25SBarry Smith 131530de9b25SBarry Smith .keywords: matrix, partitioning, register 131630de9b25SBarry Smith 131730de9b25SBarry Smith .seealso: MatPartitioningRegisterDestroy(), MatPartitioningRegisterAll() 131830de9b25SBarry Smith M*/ 1319aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 1320f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,0) 13212aabb6bbSBarry Smith #else 1322f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,d) 13232aabb6bbSBarry Smith #endif 13242aabb6bbSBarry Smith 13252bad1931SBarry Smith extern PetscTruth MatPartitioningRegisterAllCalled; 1326f1144a30SSatish Balay 1327be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegisterAll(const char[]); 1328be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegisterDestroy(void); 13292bad1931SBarry Smith 1330be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningView(MatPartitioning,PetscViewer); 1331be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetFromOptions(MatPartitioning); 1332a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningGetType(MatPartitioning,const MatPartitioningType*); 1333ca161407SBarry Smith 1334be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningParmetisSetCoarseSequential(MatPartitioning); 13350e8af3caSMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningParmetisGetEdgeCut(MatPartitioning, PetscInt *); 13360752156aSBarry Smith 1337be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningJostleSetCoarseLevel(MatPartitioning,PetscReal); 1338be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningJostleSetCoarseSequential(MatPartitioning); 133971306c60Sjroman 1340954b3ec6SBarry Smith typedef enum {MP_CHACO_MULTILEVEL_KL,MP_CHACO_SPECTRAL,MP_CHACO_LINEAR,MP_CHACO_RANDOM, MP_CHACO_SCATTERED} MPChacoGlobalType; 1341be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetGlobal(MatPartitioning, MPChacoGlobalType); 134271306c60Sjroman typedef enum { MP_CHACO_KERNIGHAN_LIN, MP_CHACO_NONE } MPChacoLocalType; 1343be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetLocal(MatPartitioning, MPChacoLocalType); 1344be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal); 134571306c60Sjroman typedef enum { MP_CHACO_LANCZOS, MP_CHACO_RQI_SYMMLQ } MPChacoEigenType; 1346be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType); 1347be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenTol(MatPartitioning, PetscReal); 1348be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenNumber(MatPartitioning, PetscInt); 134971306c60Sjroman 135071306c60Sjroman #define MP_PARTY_OPT "opt" 135171306c60Sjroman #define MP_PARTY_LIN "lin" 135271306c60Sjroman #define MP_PARTY_SCA "sca" 135371306c60Sjroman #define MP_PARTY_RAN "ran" 135471306c60Sjroman #define MP_PARTY_GBF "gbf" 135571306c60Sjroman #define MP_PARTY_GCF "gcf" 135671306c60Sjroman #define MP_PARTY_BUB "bub" 135771306c60Sjroman #define MP_PARTY_DEF "def" 1358be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetGlobal(MatPartitioning, const char*); 135971306c60Sjroman #define MP_PARTY_HELPFUL_SETS "hs" 136071306c60Sjroman #define MP_PARTY_KERNIGHAN_LIN "kl" 136171306c60Sjroman #define MP_PARTY_NONE "no" 1362be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetLocal(MatPartitioning, const char*); 1363be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal); 1364be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetBipart(MatPartitioning,PetscTruth); 1365be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscTruth); 136671306c60Sjroman 136771306c60Sjroman typedef enum { MP_SCOTCH_GREEDY, MP_SCOTCH_GPS, MP_SCOTCH_GR_GPS } MPScotchGlobalType; 1368be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetArch(MatPartitioning,const char*); 1369be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetMultilevel(MatPartitioning); 1370be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetGlobal(MatPartitioning,MPScotchGlobalType); 1371be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetCoarseLevel(MatPartitioning,PetscReal); 1372be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetHostList(MatPartitioning,const char*); 137371306c60Sjroman typedef enum { MP_SCOTCH_KERNIGHAN_LIN, MP_SCOTCH_NONE } MPScotchLocalType; 1374be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetLocal(MatPartitioning,MPScotchLocalType); 1375be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetMapping(MatPartitioning); 1376be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetStrategy(MatPartitioning,char*); 137771306c60Sjroman 13780752156aSBarry Smith /* 13790a835dfdSSatish Balay If you add entries here you must also add them to finclude/petscmat.h 1380d4fbbf0eSBarry Smith */ 13811c1c02c0SLois Curfman McInnes typedef enum { MATOP_SET_VALUES=0, 13821c1c02c0SLois Curfman McInnes MATOP_GET_ROW=1, 13831c1c02c0SLois Curfman McInnes MATOP_RESTORE_ROW=2, 13841c1c02c0SLois Curfman McInnes MATOP_MULT=3, 13851c1c02c0SLois Curfman McInnes MATOP_MULT_ADD=4, 13867c922b88SBarry Smith MATOP_MULT_TRANSPOSE=5, 13877c922b88SBarry Smith MATOP_MULT_TRANSPOSE_ADD=6, 13881c1c02c0SLois Curfman McInnes MATOP_SOLVE=7, 13891c1c02c0SLois Curfman McInnes MATOP_SOLVE_ADD=8, 13907c922b88SBarry Smith MATOP_SOLVE_TRANSPOSE=9, 13917c922b88SBarry Smith MATOP_SOLVE_TRANSPOSE_ADD=10, 13921c1c02c0SLois Curfman McInnes MATOP_LUFACTOR=11, 13931c1c02c0SLois Curfman McInnes MATOP_CHOLESKYFACTOR=12, 13941c1c02c0SLois Curfman McInnes MATOP_RELAX=13, 13951c1c02c0SLois Curfman McInnes MATOP_TRANSPOSE=14, 13961c1c02c0SLois Curfman McInnes MATOP_GETINFO=15, 13971c1c02c0SLois Curfman McInnes MATOP_EQUAL=16, 13981c1c02c0SLois Curfman McInnes MATOP_GET_DIAGONAL=17, 13991c1c02c0SLois Curfman McInnes MATOP_DIAGONAL_SCALE=18, 14001c1c02c0SLois Curfman McInnes MATOP_NORM=19, 14011c1c02c0SLois Curfman McInnes MATOP_ASSEMBLY_BEGIN=20, 14021c1c02c0SLois Curfman McInnes MATOP_ASSEMBLY_END=21, 14031c1c02c0SLois Curfman McInnes MATOP_COMPRESS=22, 14041c1c02c0SLois Curfman McInnes MATOP_SET_OPTION=23, 14051c1c02c0SLois Curfman McInnes MATOP_ZERO_ENTRIES=24, 14061c1c02c0SLois Curfman McInnes MATOP_ZERO_ROWS=25, 14071c1c02c0SLois Curfman McInnes MATOP_LUFACTOR_SYMBOLIC=26, 14081c1c02c0SLois Curfman McInnes MATOP_LUFACTOR_NUMERIC=27, 14091c1c02c0SLois Curfman McInnes MATOP_CHOLESKY_FACTOR_SYMBOLIC=28, 14101c1c02c0SLois Curfman McInnes MATOP_CHOLESKY_FACTOR_NUMERIC=29, 1411d643ce63SMatthew Knepley MATOP_SETUP_PREALLOCATION=30, 1412d643ce63SMatthew Knepley MATOP_ILUFACTOR_SYMBOLIC=31, 1413d643ce63SMatthew Knepley MATOP_ICCFACTOR_SYMBOLIC=32, 1414d643ce63SMatthew Knepley MATOP_GET_ARRAY=33, 1415d643ce63SMatthew Knepley MATOP_RESTORE_ARRAY=34, 1416e26f1c3cSBarry Smith MATOP_DUPLICATE=35, 1417d643ce63SMatthew Knepley MATOP_FORWARD_SOLVE=36, 1418d643ce63SMatthew Knepley MATOP_BACKWARD_SOLVE=37, 1419d643ce63SMatthew Knepley MATOP_ILUFACTOR=38, 1420d643ce63SMatthew Knepley MATOP_ICCFACTOR=39, 1421d643ce63SMatthew Knepley MATOP_AXPY=40, 1422d643ce63SMatthew Knepley MATOP_GET_SUBMATRICES=41, 1423d643ce63SMatthew Knepley MATOP_INCREASE_OVERLAP=42, 1424d643ce63SMatthew Knepley MATOP_GET_VALUES=43, 1425d643ce63SMatthew Knepley MATOP_COPY=44, 1426ff1cced0SMatthew Knepley MATOP_GET_ROW_MAX=45, 1427d643ce63SMatthew Knepley MATOP_SCALE=46, 1428d643ce63SMatthew Knepley MATOP_SHIFT=47, 1429d643ce63SMatthew Knepley MATOP_DIAGONAL_SHIFT=48, 1430d643ce63SMatthew Knepley MATOP_ILUDT_FACTOR=49, 1431ff1cced0SMatthew Knepley MATOP_SET_BLOCK_SIZE=50, 1432d643ce63SMatthew Knepley MATOP_GET_ROW_IJ=51, 1433d643ce63SMatthew Knepley MATOP_RESTORE_ROW_IJ=52, 1434d643ce63SMatthew Knepley MATOP_GET_COLUMN_IJ=53, 1435d643ce63SMatthew Knepley MATOP_RESTORE_COLUMN_IJ=54, 1436d643ce63SMatthew Knepley MATOP_FDCOLORING_CREATE=55, 1437d643ce63SMatthew Knepley MATOP_COLORING_PATCH=56, 1438d643ce63SMatthew Knepley MATOP_SET_UNFACTORED=57, 1439d643ce63SMatthew Knepley MATOP_PERMUTE=58, 1440d643ce63SMatthew Knepley MATOP_SET_VALUES_BLOCKED=59, 1441d643ce63SMatthew Knepley MATOP_GET_SUBMATRIX=60, 1442d643ce63SMatthew Knepley MATOP_DESTROY=61, 1443d643ce63SMatthew Knepley MATOP_VIEW=62, 1444ff1cced0SMatthew Knepley MATOP_CONVERT_FROM=63, 1445d643ce63SMatthew Knepley MATOP_USE_SCALED_FORM=64, 1446d643ce63SMatthew Knepley MATOP_SCALE_SYSTEM=65, 1447d643ce63SMatthew Knepley MATOP_UNSCALE_SYSTEM=66, 1448ba16a8cbSKris Buschelman MATOP_SET_LOCAL_TO_GLOBAL_MAP=67, 1449d643ce63SMatthew Knepley MATOP_SET_VALUES_LOCAL=68, 1450d643ce63SMatthew Knepley MATOP_ZERO_ROWS_LOCAL=69, 1451ff1cced0SMatthew Knepley MATOP_GET_ROW_MAX_ABS=70, 1452c87e5d42SMatthew Knepley MATOP_GET_ROW_MIN_ABS=71, 1453c87e5d42SMatthew Knepley MATOP_CONVERT=72, 1454c87e5d42SMatthew Knepley MATOP_SET_COLORING=73, 1455c87e5d42SMatthew Knepley MATOP_SET_VALUES_ADIC=74, 1456c87e5d42SMatthew Knepley MATOP_SET_VALUES_ADIFOR=75, 1457c87e5d42SMatthew Knepley MATOP_FD_COLORING_APPLY=76, 1458c87e5d42SMatthew Knepley MATOP_SET_FROM_OPTIONS=77, 1459c87e5d42SMatthew Knepley MATOP_MULT_CON=78, 1460c87e5d42SMatthew Knepley MATOP_MULT_TRANSPOSE_CON=79, 1461d643ce63SMatthew Knepley MATOP_PERMUTE_SPARSIFY=80, 1462d643ce63SMatthew Knepley MATOP_MULT_MULTIPLE=81, 146341acf15aSKris Buschelman MATOP_SOLVE_MULTIPLE=82, 146441acf15aSKris Buschelman MATOP_GET_INERTIA=83, 146541acf15aSKris Buschelman MATOP_LOAD=84, 146641acf15aSKris Buschelman MATOP_IS_SYMMETRIC=85, 146741acf15aSKris Buschelman MATOP_IS_HERMITIAN=86, 146841acf15aSKris Buschelman MATOP_IS_STRUCTURALLY_SYMMETRIC=87, 146941acf15aSKris Buschelman MATOP_PB_RELAX=88, 147041acf15aSKris Buschelman MATOP_GET_VECS=89, 147141acf15aSKris Buschelman MATOP_MAT_MULT=90, 147241acf15aSKris Buschelman MATOP_MAT_MULT_SYMBOLIC=91, 147341acf15aSKris Buschelman MATOP_MAT_MULT_NUMERIC=92, 147441acf15aSKris Buschelman MATOP_PTAP=93, 147541acf15aSKris Buschelman MATOP_PTAP_SYMBOLIC=94, 1476bc011b1eSHong Zhang MATOP_PTAP_NUMERIC=95, 1477bc011b1eSHong Zhang MATOP_MAT_MULTTRANSPOSE=96, 1478bc011b1eSHong Zhang MATOP_MAT_MULTTRANSPOSE_SYMBOLIC=97, 14797ba1a0bfSKris Buschelman MATOP_MAT_MULTTRANSPOSE_NUMERIC=98, 14807ba1a0bfSKris Buschelman MATOP_PTAP_SYMBOLIC_SEQAIJ=99, 14817ba1a0bfSKris Buschelman MATOP_PTAP_NUMERIC_SEQAIJ=100, 14827ba1a0bfSKris Buschelman MATOP_PTAP_SYMBOLIC_MPIAIJ=101, 148387d4246cSBarry Smith MATOP_PTAP_NUMERIC_MPIAIJ=102, 1484ff1cced0SMatthew Knepley MATOP_CONJUGATE=103, 1485ff1cced0SMatthew Knepley MATOP_SET_SIZES=104, 1486f5edf698SHong Zhang MATOP_SET_VALUES_ROW = 105, 1487ff1cced0SMatthew Knepley MATOP_REAL_PART=106, 1488ff1cced0SMatthew Knepley MATOP_IMAG_PART=107, 1489d5c63f83SSatish Balay MATOP_GET_ROW_UTRIANGULAR=108, 14902bebee5dSHong Zhang MATOP_RESTORE_ROW_UTRIANGULAR=109, 149169db28dcSHong Zhang MATOP_MATSOLVE=110, 1492829201f2SHong Zhang MATOP_GET_REDUNDANTMATRIX=111, 1493ff1cced0SMatthew Knepley MATOP_GET_ROW_MIN=112, 1494ff1cced0SMatthew Knepley MATOP_GET_COLUMN_VEC=113, 1495ff1cced0SMatthew Knepley MATOP_MISSING_DIAGONAL=114, 1496ff1cced0SMatthew Knepley MATOP_MATGETSEQNONZEROSTRUCTURE=115, 1497ff1cced0SMatthew Knepley MATOP_DESTROY_SOLVER=116 1498fae171e0SBarry Smith } MatOperation; 1499be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatHasOperation(Mat,MatOperation,PetscTruth*); 1500be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellSetOperation(Mat,MatOperation,void(*)(void)); 1501be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellGetOperation(Mat,MatOperation,void(**)(void)); 1502be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellSetContext(Mat,void*); 1503112a2221SBarry Smith 150490ace30eSBarry Smith /* 150590ace30eSBarry Smith Codes for matrices stored on disk. By default they are 150690ace30eSBarry Smith stored in a universal format. By changing the format with 15077973ad04SBarry Smith PetscViewerSetFormat(viewer,PETSC_VIEWER_NATIVE); the matrices will 150890ace30eSBarry Smith be stored in a way natural for the matrix, for example dense matrices 150990ace30eSBarry Smith would be stored as dense. Matrices stored this way may only be 151090ace30eSBarry Smith read into matrices of the same time. 151190ace30eSBarry Smith */ 151290ace30eSBarry Smith #define MATRIX_BINARY_FORMAT_DENSE -1 151390ace30eSBarry Smith 1514be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetHashTableFactor(Mat,PetscReal); 1515be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatISGetLocalMat(Mat,Mat*); 15161f4e1ec7SBarry Smith 1517d9274352SBarry Smith /*S 1518d9274352SBarry Smith MatNullSpace - Object that removes a null space from a vector, i.e. 1519d9274352SBarry Smith orthogonalizes the vector to a subsapce 1520d9274352SBarry Smith 1521f7a9e4ceSBarry Smith Level: advanced 1522d9274352SBarry Smith 1523d9274352SBarry Smith Concepts: matrix; linear operator, null space 1524d9274352SBarry Smith 15256e1639daSBarry Smith Users manual sections: 15266e1639daSBarry Smith . sec_singular 15276e1639daSBarry Smith 1528d9274352SBarry Smith .seealso: MatNullSpaceCreate() 1529d9274352SBarry Smith S*/ 153074637425SBarry Smith typedef struct _p_MatNullSpace* MatNullSpace; 1531d9274352SBarry Smith 1532be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceCreate(MPI_Comm,PetscTruth,PetscInt,const Vec[],MatNullSpace*); 1533281d1b2eSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(Vec,void*),void*); 1534be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceDestroy(MatNullSpace); 1535be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceRemove(MatNullSpace,Vec,Vec*); 1536be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceAttach(Mat,MatNullSpace); 153795902228SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceTest(MatNullSpace,Mat,PetscTruth *); 153874637425SBarry Smith 1539be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatReorderingSeqSBAIJ(Mat,IS); 1540be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPISBAIJSetHashTableFactor(Mat,PetscReal); 1541be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetColumnIndices(Mat,PetscInt *); 154246129b97SKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJInvertBlockDiagonal(Mat); 15433f1d51d7SBarry Smith 1544be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMAIJ(Mat,PetscInt,Mat*); 1545be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMAIJRedimension(Mat,PetscInt,Mat*); 1546be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMAIJGetAIJ(Mat,Mat*); 1547c4f061fbSSatish Balay 1548be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatComputeExplicitOperator(Mat,Mat*); 1549b0a32e0cSBarry Smith 1550be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalScaleLocal(Mat,Vec); 155104f1ad80SBarry Smith 1552fef1beadSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*); 15533ec795f1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetBase(Mat,Vec,Vec); 15543ec795f1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*); 1555e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*)); 1556e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec)); 1557e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDAddNullSpace(Mat,MatNullSpace); 1558e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt); 1559e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDResetHHistory(Mat); 1560e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctionError(Mat,PetscReal); 1561e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetPeriod(Mat,PetscInt); 1562e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDGetH(Mat,PetscScalar *); 1563e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFromOptions(Mat); 1564e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*); 1565e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*); 1566e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDDSSetUmin(Mat,PetscReal); 1567e884886fSBarry Smith 1568e884886fSBarry Smith 1569e884886fSBarry Smith typedef struct _p_MatMFFD* MatMFFD; 1570e884886fSBarry Smith 1571e884886fSBarry Smith /*E 1572e884886fSBarry Smith MatMFFDType - algorithm used to compute the h used in computing matrix-vector products via differencing of the function 1573e884886fSBarry Smith 1574e884886fSBarry Smith Level: beginner 1575e884886fSBarry Smith 1576e884886fSBarry Smith .seealso: MatMFFDSetType(), MatMFFDRegister() 1577e884886fSBarry Smith E*/ 1578a313700dSBarry Smith #define MatMFFDType char* 1579e884886fSBarry Smith #define MATMFFD_DS "ds" 1580e884886fSBarry Smith #define MATMFFD_WP "wp" 1581e884886fSBarry Smith 1582a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetType(Mat,const MatMFFDType); 1583e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatMFFD)); 1584e884886fSBarry Smith 1585e884886fSBarry Smith /*MC 1586e884886fSBarry Smith MatMFFDRegisterDynamic - Adds a method to the MatMFFD registry. 1587e884886fSBarry Smith 1588e884886fSBarry Smith Synopsis: 1589e884886fSBarry Smith PetscErrorCode MatMFFDRegisterDynamic(char *name_solver,char *path,char *name_create,PetscErrorCode (*routine_create)(MatMFFD)) 1590e884886fSBarry Smith 1591e884886fSBarry Smith Not Collective 1592e884886fSBarry Smith 1593e884886fSBarry Smith Input Parameters: 1594e884886fSBarry Smith + name_solver - name of a new user-defined compute-h module 1595e884886fSBarry Smith . path - path (either absolute or relative) the library containing this solver 1596e884886fSBarry Smith . name_create - name of routine to create method context 1597e884886fSBarry Smith - routine_create - routine to create method context 1598e884886fSBarry Smith 1599e884886fSBarry Smith Level: developer 1600e884886fSBarry Smith 1601e884886fSBarry Smith Notes: 1602e884886fSBarry Smith MatMFFDRegisterDynamic() may be called multiple times to add several user-defined solvers. 1603e884886fSBarry Smith 1604e884886fSBarry Smith If dynamic libraries are used, then the fourth input argument (routine_create) 1605e884886fSBarry Smith is ignored. 1606e884886fSBarry Smith 1607e884886fSBarry Smith Sample usage: 1608e884886fSBarry Smith .vb 1609e884886fSBarry Smith MatMFFDRegisterDynamic("my_h",/home/username/my_lib/lib/libO/solaris/mylib.a, 1610e884886fSBarry Smith "MyHCreate",MyHCreate); 1611e884886fSBarry Smith .ve 1612e884886fSBarry Smith 1613e884886fSBarry Smith Then, your solver can be chosen with the procedural interface via 1614e884886fSBarry Smith $ MatMFFDSetType(mfctx,"my_h") 1615e884886fSBarry Smith or at runtime via the option 1616e884886fSBarry Smith $ -snes_mf_type my_h 1617e884886fSBarry Smith 1618e884886fSBarry Smith .keywords: MatMFFD, register 1619e884886fSBarry Smith 1620e884886fSBarry Smith .seealso: MatMFFDRegisterAll(), MatMFFDRegisterDestroy() 1621e884886fSBarry Smith M*/ 1622e884886fSBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 1623e884886fSBarry Smith #define MatMFFDRegisterDynamic(a,b,c,d) MatMFFDRegister(a,b,c,0) 1624e884886fSBarry Smith #else 1625e884886fSBarry Smith #define MatMFFDRegisterDynamic(a,b,c,d) MatMFFDRegister(a,b,c,d) 1626e884886fSBarry Smith #endif 1627e884886fSBarry Smith 1628e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDRegisterAll(const char[]); 1629e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDRegisterDestroy(void); 1630e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDDefaultSetUmin(Mat,PetscReal); 1631e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDWPSetComputeNormU(Mat,PetscTruth); 1632e884886fSBarry Smith 1633e884886fSBarry Smith 1634be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *); 1635be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *); 16367dbadf16SMatthew Knepley 1637e9fa29b7SSatish Balay PETSC_EXTERN_CXX_END 16382eac72dbSBarry Smith #endif 1639