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 27d9274352SBarry Smith .seealso: MatSetType(), Mat 28d91e6319SBarry Smith E*/ 29e5a9bf91SBarry Smith #define MatType const 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" 54b3a1e11cSKris Buschelman #define MATSEQAIJSPOOLES "seqaijspooles" 55d10c748bSKris Buschelman #define MATMPIAIJSPOOLES "mpiaijspooles" 569abb65ffSKris Buschelman #define MATSEQSBAIJSPOOLES "seqsbaijspooles" 5722191285SKris Buschelman #define MATMPISBAIJSPOOLES "mpisbaijspooles" 58bb4d4166SHong Zhang #define MATAIJSPOOLES "aijspooles" 59bb4d4166SHong Zhang #define MATSBAIJSPOOLES "sbaijspooles" 6014b396bbSKris Buschelman #define MATSUPERLU "superlu" 611d96aa28SKris Buschelman #define MATSUPERLU_DIST "superlu_dist" 621677d0b8SKris Buschelman #define MATUMFPACK "umfpack" 63e8aa55a4SKris Buschelman #define MATESSL "essl" 644eb8e494SKris Buschelman #define MATLUSOL "lusol" 65397b6df1SKris Buschelman #define MATAIJMUMPS "aijmumps" 66397b6df1SKris Buschelman #define MATSBAIJMUMPS "sbaijmumps" 678da957c5SKris Buschelman #define MATDSCPACK "dscpack" 687065e2aaSKris Buschelman #define MATMATLAB "matlab" 694b8d99adSRichard Tran Mills #define MATSEQCSRPERM "seqcsrperm" 7018def556SRichard Tran Mills #define MATMPICSRPERM "mpicsrperm" 71814655b8SBarry Smith #define MATCSRPERM "csrperm" 7281824310SBarry Smith #define MATSEQCRL "seqcrl" 73c02b24c5SBarry Smith #define MATMPICRL "mpicrl" 74c02b24c5SBarry Smith #define MATCRL "crl" 75c0aa2d19SHong Zhang #define MATPLAPACK "plapack" 762a6744ebSBarry Smith #define MATSCATTER "scatter" 77421e10b8SBarry Smith #define MATBLOCKMAT "blockmat" 78793850ffSBarry Smith #define MATCOMPOSITE "composite" 795a7f1df3SHong Zhang #define MATSEQFFTW "seqfftw" 80d91e6319SBarry Smith 81c06d978dSMatthew Knepley /* Logging support */ 82552e946dSBarry Smith #define MAT_FILE_COOKIE 1211216 /* used to indicate matrices in binary files */ 83be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_COOKIE; 84be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_FDCOLORING_COOKIE; 85be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_PARTITIONING_COOKIE; 86be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_NULLSPACE_COOKIE; 87e884886fSBarry Smith extern PetscCookie PETSCMAT_DLLEXPORT MATMFFD_COOKIE; 88c06d978dSMatthew Knepley 89ceb03754SKris Buschelman /*E 90ceb03754SKris Buschelman MatReuse - Indicates if matrices obtained from a previous call to MatGetSubMatrices() 91ceb03754SKris Buschelman or MatGetSubMatrix() are to be reused to store the new matrix values. 92ceb03754SKris Buschelman 93ceb03754SKris Buschelman Level: beginner 94ceb03754SKris Buschelman 95ceb03754SKris Buschelman Any additions/changes here MUST also be made in include/finclude/petscmat.h 96ceb03754SKris Buschelman 970c451bc4SBarry Smith .seealso: MatGetSubMatrices(), MatGetSubMatrix(), MatDestroyMatrices(), MatConvert() 98ceb03754SKris Buschelman E*/ 99ceb03754SKris Buschelman typedef enum {MAT_INITIAL_MATRIX,MAT_REUSE_MATRIX} MatReuse; 100ceb03754SKris Buschelman 1015494a064SHong Zhang /*E 1025494a064SHong Zhang MatGetSubMatrixOption - Indicates if matrices obtained from a call to MatGetSubMatrices() 103829201f2SHong Zhang include the matrix values. Currently it is only used by MatGetSeqNonzerostructure(). 1045494a064SHong Zhang 1055494a064SHong Zhang Level: beginner 1065494a064SHong Zhang 107829201f2SHong Zhang .seealso: MatGetSeqNonzerostructure() 1085494a064SHong Zhang E*/ 1095494a064SHong Zhang typedef enum {MAT_DO_NOT_GET_VALUES,MAT_GET_VALUES} MatGetSubMatrixOption; 1105494a064SHong Zhang 111e5bd5246SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInitializePackage(const char[]); 112c06d978dSMatthew Knepley 113f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreate(MPI_Comm,Mat*); 114a69a960dSMatthew Knepley PetscPolymorphicFunction(MatCreate,(MPI_Comm comm),(comm,&A),Mat,A) 115a69a960dSMatthew Knepley PetscPolymorphicFunction(MatCreate,(),(PETSC_COMM_WORLD,&A),Mat,A) 116f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetSizes(Mat,PetscInt,PetscInt,PetscInt,PetscInt); 117f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetType(Mat,MatType); 118be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetFromOptions(Mat); 119be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUpPreallocation(Mat); 120be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRegisterAll(const char[]); 121be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat)); 12230de9b25SBarry Smith 123f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetOptionsPrefix(Mat,const char[]); 124f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAppendOptionsPrefix(Mat,const char[]); 125f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOptionsPrefix(Mat,const char*[]); 126f69a0ea3SMatthew Knepley 12730de9b25SBarry Smith /*MC 12830de9b25SBarry Smith MatRegisterDynamic - Adds a new matrix type 12930de9b25SBarry Smith 13030de9b25SBarry Smith Synopsis: 131c1ac3661SBarry Smith PetscErrorCode MatRegisterDynamic(char *name,char *path,char *name_create,PetscErrorCode (*routine_create)(Mat)) 13230de9b25SBarry Smith 13330de9b25SBarry Smith Not Collective 13430de9b25SBarry Smith 13530de9b25SBarry Smith Input Parameters: 13630de9b25SBarry Smith + name - name of a new user-defined matrix type 13730de9b25SBarry Smith . path - path (either absolute or relative) the library containing this solver 13830de9b25SBarry Smith . name_create - name of routine to create method context 13930de9b25SBarry Smith - routine_create - routine to create method context 14030de9b25SBarry Smith 14130de9b25SBarry Smith Notes: 14230de9b25SBarry Smith MatRegisterDynamic() may be called multiple times to add several user-defined solvers. 14330de9b25SBarry Smith 14430de9b25SBarry Smith If dynamic libraries are used, then the fourth input argument (routine_create) 14530de9b25SBarry Smith is ignored. 14630de9b25SBarry Smith 14730de9b25SBarry Smith Sample usage: 14830de9b25SBarry Smith .vb 14930de9b25SBarry Smith MatRegisterDynamic("my_mat",/home/username/my_lib/lib/libO/solaris/mylib.a, 15030de9b25SBarry Smith "MyMatCreate",MyMatCreate); 15130de9b25SBarry Smith .ve 15230de9b25SBarry Smith 15330de9b25SBarry Smith Then, your solver can be chosen with the procedural interface via 15430de9b25SBarry Smith $ MatSetType(Mat,"my_mat") 15530de9b25SBarry Smith or at runtime via the option 15630de9b25SBarry Smith $ -mat_type my_mat 15730de9b25SBarry Smith 15830de9b25SBarry Smith Level: advanced 15930de9b25SBarry Smith 160ab901514SBarry Smith Notes: ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values. 16130de9b25SBarry Smith If your function is not being put into a shared library then use VecRegister() instead 16230de9b25SBarry Smith 16330de9b25SBarry Smith .keywords: Mat, register 16430de9b25SBarry Smith 16530de9b25SBarry Smith .seealso: MatRegisterAll(), MatRegisterDestroy() 16630de9b25SBarry Smith 16730de9b25SBarry Smith M*/ 168273d9f13SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 169273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,0) 170273d9f13SBarry Smith #else 171273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,d) 17230de9b25SBarry Smith #endif 17330de9b25SBarry Smith 1749c666560SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRegisterAll(const char[]); 175273d9f13SBarry Smith extern PetscTruth MatRegisterAllCalled; 176b0a32e0cSBarry Smith extern PetscFList MatList; 17728988994SBarry Smith 178be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqDense(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*); 179be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIDense(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*); 180be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 181ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,nz,nnz,&A),Mat,A) 182ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,&A),Mat,A) 183ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,0,nnz,&A),Mat,A) 184ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,&A),Mat,A) 185ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,A)) 186ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,m,n,0,nnz,A)) 187ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,A)) 188be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 189ba966639SSatish 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) 190ba966639SSatish 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) 191ba966639SSatish 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) 192ba966639SSatish 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) 193ba966639SSatish 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)) 194ba966639SSatish 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)) 195ba966639SSatish 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)) 196ba966639SSatish 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) 197ba966639SSatish 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) 198ba966639SSatish 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) 199ba966639SSatish 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) 200ba966639SSatish 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)) 201ba966639SSatish 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)) 202ba966639SSatish 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)) 20363ab93bfSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *); 2048d7a6e47SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJWithSplitArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],PetscInt[],PetscInt[],PetscScalar[],Mat*); 2058d7a6e47SBarry Smith 206be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIRowbs(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 207be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 208ba966639SSatish 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) 209ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A) 210ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A) 211ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A) 212ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A)) 213ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A)) 214ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A)) 215be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 216ba966639SSatish 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) 217ba966639SSatish 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) 218ba966639SSatish 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) 219ba966639SSatish 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) 220ba966639SSatish 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)) 221ba966639SSatish 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)) 222ba966639SSatish 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)) 223ba966639SSatish 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) 224ba966639SSatish 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) 225ba966639SSatish 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) 226ba966639SSatish 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) 227ba966639SSatish 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)) 228ba966639SSatish 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)) 229ba966639SSatish 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)) 230be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAdj(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscInt[],Mat*); 231be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 232ba966639SSatish 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) 233ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A) 234ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A) 235ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A) 236ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A)) 237ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A)) 238ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A)) 239be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPISBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 240ba966639SSatish 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) 241ba966639SSatish 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) 242ba966639SSatish 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) 243ba966639SSatish 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) 244ba966639SSatish 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)) 245ba966639SSatish 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)) 246ba966639SSatish 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)) 247ba966639SSatish 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) 248ba966639SSatish 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) 249ba966639SSatish 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) 250ba966639SSatish 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) 251ba966639SSatish 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)) 252ba966639SSatish 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)) 253ba966639SSatish 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)) 254be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateShell(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,void *,Mat*); 255ba966639SSatish 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) 256ba966639SSatish Balay PetscPolymorphicFunction(MatCreateShell,(PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx),(PETSC_COMM_WORLD,m,n,M,N,ctx,&A),Mat,A) 257be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateAdic(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,void (*)(void),Mat*); 258be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateNormal(Mat,Mat*); 259ba966639SSatish Balay PetscPolymorphicFunction(MatCreateNormal,(Mat mat),(mat,&A),Mat,A) 260ab92ecdeSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateLRC(Mat,Mat,Mat,Mat*); 261284134d9SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateIS(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,ISLocalToGlobalMapping,Mat*); 2621472f72bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 2631472f72bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPICRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 2642a6744ebSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateScatter(MPI_Comm,VecScatter,Mat*); 2652a6744ebSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScatterSetVecScatter(Mat,VecScatter); 2662a6744ebSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScatterGetVecScatter(Mat,VecScatter*); 2678cdf2d9bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateBlockMat(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,Mat*); 268793850ffSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompositeAddMat(Mat,Mat); 269b52f573bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompositeMerge(Mat); 270793850ffSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateComposite(MPI_Comm,PetscInt,const Mat*,Mat*); 2715a7f1df3SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqFFTW(MPI_Comm,PetscInt,const PetscInt[],Mat*); 2721472f72bSBarry Smith 273f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUp(Mat); 274be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroy(Mat); 27521c89e3eSBarry Smith 2760c451bc4SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConjugate(Mat); 27799cafbc1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRealPart(Mat); 27899cafbc1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatImaginaryPart(Mat); 27999cafbc1SBarry Smith 2808ed539a5SBarry Smith /* ------------------------------------------------------------*/ 281be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 282be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 28387d4246cSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesRow(Mat,PetscInt,const PetscScalar[]); 284f259bd47SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesRowLocal(Mat,PetscInt,const PetscScalar[]); 28584cb2905SBarry Smith 2862ef4de8bSBarry Smith /*S 2872ef4de8bSBarry Smith MatStencil - Data structure (C struct) for storing information about a single row or 2882ef4de8bSBarry Smith column of a matrix as index on an associated grid. 2892ef4de8bSBarry Smith 2902ef4de8bSBarry Smith Level: beginner 2912ef4de8bSBarry Smith 2922ef4de8bSBarry Smith Concepts: matrix; linear operator 2932ef4de8bSBarry Smith 294d7bd93baSBarry Smith .seealso: MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockStencil() 2952ef4de8bSBarry Smith S*/ 296435da068SBarry Smith typedef struct { 297c1ac3661SBarry Smith PetscInt k,j,i,c; 298435da068SBarry Smith } MatStencil; 2992ef4de8bSBarry Smith 300be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode); 301be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode); 302be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt); 303435da068SBarry Smith 304be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetColoring(Mat,ISColoring); 305be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesAdic(Mat,void*); 306be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesAdifor(Mat,PetscInt,void*); 3073a7fca6bSBarry Smith 308d91e6319SBarry Smith /*E 309d91e6319SBarry Smith MatAssemblyType - Indicates if the matrix is now to be used, or if you plan 310d91e6319SBarry Smith to continue to add values to it 311d91e6319SBarry Smith 312d91e6319SBarry Smith Level: beginner 313d91e6319SBarry Smith 314d91e6319SBarry Smith .seealso: MatAssemblyBegin(), MatAssemblyEnd() 315d91e6319SBarry Smith E*/ 3166d4a8577SBarry Smith typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType; 317be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssemblyBegin(Mat,MatAssemblyType); 318be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssemblyEnd(Mat,MatAssemblyType); 319be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssembled(Mat,PetscTruth*); 3204f9c727eSBarry Smith 3211d73ed98SBarry Smith 32230de9b25SBarry Smith 323d91e6319SBarry Smith /*E 324d91e6319SBarry Smith MatOption - Options that may be set for a matrix and its behavior or storage 325d91e6319SBarry Smith 326d91e6319SBarry Smith Level: beginner 327d91e6319SBarry Smith 3280a835dfdSSatish Balay Any additions/changes here MUST also be made in include/finclude/petscmat.h 329d91e6319SBarry Smith 330d91e6319SBarry Smith .seealso: MatSetOption() 331d91e6319SBarry Smith E*/ 332512a5fc5SBarry Smith typedef enum {MAT_ROW_ORIENTED,MAT_NEW_NONZERO_LOCATIONS, 3334e0d8c25SBarry Smith MAT_SYMMETRIC, 3344e0d8c25SBarry Smith MAT_STRUCTURALLY_SYMMETRIC, 3358047e77bSBarry Smith MAT_NEW_DIAGONALS,MAT_IGNORE_OFF_PROC_ENTRIES, 3364e0d8c25SBarry Smith MAT_NEW_NONZERO_LOCATION_ERR, 3374e0d8c25SBarry Smith MAT_NEW_NONZERO_ALLOCATION_ERR,MAT_USE_HASH_TABLE, 3384e0d8c25SBarry Smith MAT_KEEP_ZEROED_ROWS,MAT_IGNORE_ZERO_ENTRIES, 3394e0d8c25SBarry Smith MAT_USE_INODES, 3404e0d8c25SBarry Smith MAT_HERMITIAN, 3414e0d8c25SBarry Smith MAT_SYMMETRY_ETERNAL, 3424e0d8c25SBarry Smith MAT_USE_COMPRESSEDROW, 3434e0d8c25SBarry Smith MAT_IGNORE_LOWER_TRIANGULAR,MAT_ERROR_LOWER_TRIANGULAR, 3444e0d8c25SBarry Smith MAT_GETROW_UPPERTRIANGULAR} MatOption; 345290bbb0aSBarry Smith extern const char *MatOptions[]; 3464e0d8c25SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetOption(Mat,MatOption,PetscTruth); 347be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetType(Mat,MatType*); 348ba966639SSatish Balay PetscPolymorphicFunction(MatGetType,(Mat mat),(mat,&t),MatType,t) 34984cb2905SBarry Smith 350be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]); 351be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 352be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 353f5edf698SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowUpperTriangular(Mat); 354f5edf698SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRowUpperTriangular(Mat); 355be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 356be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 357be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumnVector(Mat,Vec,PetscInt); 358be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetArray(Mat,PetscScalar *[]); 359ba966639SSatish Balay PetscPolymorphicFunction(MatGetArray,(Mat mat),(mat,&a),PetscScalar*,a) 360be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreArray(Mat,PetscScalar *[]); 361be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBlockSize(Mat,PetscInt *); 362ba966639SSatish Balay PetscPolymorphicFunction(MatGetBlockSize,(Mat mat),(mat,&a),PetscInt,a) 363be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetBlockSize(Mat,PetscInt); 3647b80b807SBarry Smith 3651620fd73SBarry Smith 366be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMult(Mat,Vec,Vec); 367be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultAdd(Mat,Vec,Vec,Vec); 368ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultAdd,(Mat A,Vec x,Vec y),(A,x,y,y)) 369be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTranspose(Mat,Vec,Vec); 370be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsTranspose(Mat,Mat,PetscReal,PetscTruth*); 371ba966639SSatish Balay PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B,PetscReal tol),(A,B,tol,&t),PetscTruth,t) 372e37f64cdSBarry Smith PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B),(A,B,0,&t),PetscTruth,t) 3731cbb95d3SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitianTranspose(Mat,Mat,PetscReal,PetscTruth*); 374be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeAdd(Mat,Vec,Vec,Vec); 375ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultTransposeAdd,(Mat A,Vec x,Vec y),(A,x,y,y)) 376be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultConstrained(Mat,Vec,Vec); 377be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeConstrained(Mat,Vec,Vec); 3782bebee5dSHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatSolve(Mat,Mat,Mat); 379f5edf698SHong Zhang 380d91e6319SBarry Smith /*E 381d91e6319SBarry Smith MatDuplicateOption - Indicates if a duplicated sparse matrix should have 382d91e6319SBarry Smith its numerical values copied over or just its nonzero structure. 383d91e6319SBarry Smith 384d91e6319SBarry Smith Level: beginner 385d91e6319SBarry Smith 386d91e6319SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 387d91e6319SBarry Smith 388d91e6319SBarry Smith .seealso: MatDuplicate() 389d91e6319SBarry Smith E*/ 3902e8a6d31SBarry Smith typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES} MatDuplicateOption; 3912e8a6d31SBarry Smith 392f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert(Mat,MatType,MatReuse,Mat*); 393ba966639SSatish Balay PetscPolymorphicFunction(MatConvert,(Mat A,MatType t),(A,t,MAT_INITIAL_MATRIX,&a),Mat,a) 394be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDuplicate(Mat,MatDuplicateOption,Mat*); 395ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A,MatDuplicateOption o),(A,o,&a),Mat,a) 396ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A),(A,MAT_COPY_VALUES,&a),Mat,a) 39794a9d846SBarry Smith 398d91e6319SBarry Smith /*E 399d91e6319SBarry Smith MatStructure - Indicates if the matrix has the same nonzero structure 400d91e6319SBarry Smith 401d91e6319SBarry Smith Level: beginner 402d91e6319SBarry Smith 403d91e6319SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 404d91e6319SBarry Smith 40594b7f48cSBarry Smith .seealso: MatCopy(), KSPSetOperators(), PCSetOperators() 406d91e6319SBarry Smith E*/ 407c537a176SHong Zhang typedef enum {SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER,SUBSET_NONZERO_PATTERN} MatStructure; 408cb5b572fSBarry Smith 409be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCopy(Mat,Mat,MatStructure); 410be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatView(Mat,PetscViewer); 411be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsSymmetric(Mat,PetscReal,PetscTruth*); 412ba966639SSatish Balay PetscPolymorphicFunction(MatIsSymmetric,(Mat A,PetscReal tol),(A,tol,&t),PetscTruth,t) 413e37f64cdSBarry Smith PetscPolymorphicFunction(MatIsSymmetric,(Mat A),(A,0,&t),PetscTruth,t) 414be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsStructurallySymmetric(Mat,PetscTruth*); 415ba966639SSatish Balay PetscPolymorphicFunction(MatIsStructurallySymmetric,(Mat A),(A,&t),PetscTruth,t) 4161cbb95d3SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitian(Mat,PetscReal,PetscTruth*); 4171cbb95d3SBarry Smith PetscPolymorphicFunction(MatIsHermitian,(Mat A),(A,0,&t),PetscTruth,t) 418be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsSymmetricKnown(Mat,PetscTruth*,PetscTruth*); 419be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitianKnown(Mat,PetscTruth*,PetscTruth*); 42064c62002SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMissingDiagonal(Mat,PetscTruth *,PetscInt *); 421f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLoad(PetscViewer,MatType,Mat*); 422ba966639SSatish Balay PetscPolymorphicFunction(MatLoad,(PetscViewer v,MatType t),(v,t,&a),Mat,a) 4237b80b807SBarry Smith 4248f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *); 4258f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRowIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *); 4268f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumnIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *); 4278f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreColumnIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *); 428d4fbbf0eSBarry Smith 429d91e6319SBarry Smith /*S 430d91e6319SBarry Smith MatInfo - Context of matrix information, used with MatGetInfo() 431d91e6319SBarry Smith 432d91e6319SBarry Smith In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE 433d91e6319SBarry Smith 434d91e6319SBarry Smith Level: intermediate 435d91e6319SBarry Smith 436d91e6319SBarry Smith Concepts: matrix^nonzero information 437d91e6319SBarry Smith 438d9274352SBarry Smith .seealso: MatGetInfo(), MatInfoType 439d91e6319SBarry Smith S*/ 4404e220ebcSLois Curfman McInnes typedef struct { 441b0a32e0cSBarry Smith PetscLogDouble rows_global,columns_global; /* number of global rows and columns */ 442b0a32e0cSBarry Smith PetscLogDouble rows_local,columns_local; /* number of local rows and columns */ 443b0a32e0cSBarry Smith PetscLogDouble block_size; /* block size */ 444b0a32e0cSBarry Smith PetscLogDouble nz_allocated,nz_used,nz_unneeded; /* number of nonzeros */ 445b0a32e0cSBarry Smith PetscLogDouble memory; /* memory allocated */ 446b0a32e0cSBarry Smith PetscLogDouble assemblies; /* number of matrix assemblies called */ 447b0a32e0cSBarry Smith PetscLogDouble mallocs; /* number of mallocs during MatSetValues() */ 448b0a32e0cSBarry Smith PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */ 449b0a32e0cSBarry Smith PetscLogDouble factor_mallocs; /* number of mallocs during factorization */ 4504e220ebcSLois Curfman McInnes } MatInfo; 4514e220ebcSLois Curfman McInnes 452d9274352SBarry Smith /*E 453d9274352SBarry Smith MatInfoType - Indicates if you want information about the local part of the matrix, 454d9274352SBarry Smith the entire parallel matrix or the maximum over all the local parts. 455d9274352SBarry Smith 456d9274352SBarry Smith Level: beginner 457d9274352SBarry Smith 458d9274352SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 459d9274352SBarry Smith 460d9274352SBarry Smith .seealso: MatGetInfo(), MatInfo 461d9274352SBarry Smith E*/ 4627b80b807SBarry Smith typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType; 463be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetInfo(Mat,MatInfoType,MatInfo*); 464be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatValid(Mat,PetscTruth*); 465be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetDiagonal(Mat,Vec); 466985db425SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMax(Mat,Vec,PetscInt[]); 467985db425SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMin(Mat,Vec,PetscInt[]); 468985db425SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMaxAbs(Mat,Vec,PetscInt[]); 46986697f06SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowSum(Mat,Vec); 470*fc4dec0aSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatTranspose(Mat,MatReuse,Mat*); 471*fc4dec0aSBarry Smith PetscPolymorphicFunction(MatTranspose,(Mat A),(A,MAT_INITIAL_MATRIX,&t),Mat,t) 472be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPermute(Mat,IS,IS,Mat *); 473ba966639SSatish Balay PetscPolymorphicFunction(MatPermute,(Mat A,IS is1,IS is2),(A,is1,is2,&t),Mat,t) 474be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPermuteSparsify(Mat,PetscInt,PetscReal,PetscReal,IS,IS,Mat *); 475be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalScale(Mat,Vec,Vec); 476be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalSet(Mat,Vec,InsertMode); 477be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatEqual(Mat,Mat,PetscTruth*); 478ba966639SSatish Balay PetscPolymorphicFunction(MatEqual,(Mat A,Mat B),(A,B,&t),PetscTruth,t) 479be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultEqual(Mat,Mat,PetscInt,PetscTruth*); 480be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultAddEqual(Mat,Mat,PetscInt,PetscTruth*); 481be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeEqual(Mat,Mat,PetscInt,PetscTruth*); 482be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscTruth*); 4837b80b807SBarry Smith 484be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNorm(Mat,NormType,PetscReal *); 485ba966639SSatish Balay PetscPolymorphicFunction(MatNorm,(Mat A,NormType t),(A,t,&n),PetscReal,n) 486be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroEntries(Mat); 487f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar); 488f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsIS(Mat,IS,PetscScalar); 489f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroColumns(Mat,PetscInt,const PetscInt [],const PetscScalar*); 490f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroColumnsIS(Mat,IS,const PetscScalar*); 4917b80b807SBarry Smith 492be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatUseScaledForm(Mat,PetscTruth); 493be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScaleSystem(Mat,Vec,Vec); 494be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatUnScaleSystem(Mat,Vec,Vec); 4955ef9f2a5SBarry Smith 496be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSize(Mat,PetscInt*,PetscInt*); 497be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalSize(Mat,PetscInt*,PetscInt*); 498be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRange(Mat,PetscInt*,PetscInt*); 4998ee2e534SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRanges(Mat,const PetscInt**); 500f4b8409dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRangeColumn(Mat,PetscInt*,PetscInt*); 5017b80b807SBarry Smith 502be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]); 503be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroyMatrices(PetscInt,Mat *[]); 504be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrix(Mat,IS,IS,PetscInt,MatReuse,Mat *); 505abf0e977SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrixRaw(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscInt,MatReuse,Mat *); 506829201f2SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSeqNonzeroStructure(Mat,Mat *[]); 507829201f2SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroySeqNonzeroStructure(Mat *[]); 5085494a064SHong Zhang 509be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge(MPI_Comm,Mat,PetscInt,MatReuse,Mat*); 510be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPI(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*); 511be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPISymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*); 512be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPINumeric(Mat,Mat); 513be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroy_MPIAIJ_SeqsToMPI(Mat); 514be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalMat(Mat,MatReuse,Mat*); 515be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*); 516be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,PetscInt*,Mat*); 517dd6ea824SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBrowsOfAoCols(Mat,Mat,MatReuse,PetscInt**,MatScalar**,Mat*); 51843eb5e2fSMatthew Knepley #if defined (PETSC_USE_CTABLE) 519cd116e26SSatish Balay #include "petscctable.h" 52043eb5e2fSMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetCommunicationStructs(Mat, Vec *, PetscTable *, VecScatter *); 52143eb5e2fSMatthew Knepley #else 52243eb5e2fSMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetCommunicationStructs(Mat, Vec *, PetscInt *[], VecScatter *); 52343eb5e2fSMatthew Knepley #endif 5248efafbd8SBarry Smith 525be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt); 5267b80b807SBarry Smith 527be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*); 528be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultSymbolic(Mat,Mat,PetscReal,Mat*); 529be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultNumeric(Mat,Mat,Mat); 53022440eb1SKris Buschelman 531be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*); 532be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAPSymbolic(Mat,Mat,PetscReal,Mat*); 533be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAPNumeric(Mat,Mat,Mat); 53422440eb1SKris Buschelman 535be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTranspose(Mat,Mat,MatReuse,PetscReal,Mat*); 536be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTransposeSymbolic(Mat,Mat,PetscReal,Mat*); 537be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTransposeNumeric(Mat,Mat,Mat); 538bc011b1eSHong Zhang 539f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAXPY(Mat,PetscScalar,Mat,MatStructure); 54004aac2b0SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAYPX(Mat,PetscScalar,Mat,MatStructure); 541be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompress(Mat); 5421c741599SHong Zhang 543f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScale(Mat,PetscScalar); 544f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShift(Mat,PetscScalar); 5457b80b807SBarry Smith 546be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping); 547be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetLocalToGlobalMappingBlock(Mat,ISLocalToGlobalMapping); 548f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar); 549f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsLocalIS(Mat,IS,PetscScalar); 550be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 551be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 552052efed2SBarry Smith 553be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStashSetInitialSize(Mat,PetscInt,PetscInt); 554be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*); 55590f02eecSBarry Smith 556be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInterpolate(Mat,Vec,Vec); 5570c451bc4SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInterpolateAdd(Mat,Vec,Vec,Vec); 558ba966639SSatish Balay PetscPolymorphicSubroutine(MatInterpolateAdd,(Mat A,Vec x,Vec y),(A,x,y,y)) 559be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestrict(Mat,Vec,Vec); 560be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetVecs(Mat,Vec*,Vec*); 56169db28dcSHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRedundantMatrix(Mat,PetscInt,MPI_Comm,PetscInt,MatReuse,Mat*); 562bd481603SBarry Smith 563bd481603SBarry Smith /*MC 5641620fd73SBarry Smith MatSetValue - Set a single entry into a matrix. 5651620fd73SBarry Smith 5661620fd73SBarry Smith Not collective 5671620fd73SBarry Smith 5681620fd73SBarry Smith Input Parameters: 5691620fd73SBarry Smith + m - the matrix 5701620fd73SBarry Smith . row - the row location of the entry 5711620fd73SBarry Smith . col - the column location of the entry 5721620fd73SBarry Smith . value - the value to insert 5731620fd73SBarry Smith - mode - either INSERT_VALUES or ADD_VALUES 5741620fd73SBarry Smith 5751620fd73SBarry Smith Notes: 5761620fd73SBarry Smith For efficiency one should use MatSetValues() and set several or many 5771620fd73SBarry Smith values simultaneously if possible. 5781620fd73SBarry Smith 5791620fd73SBarry Smith Level: beginner 5801620fd73SBarry Smith 5811620fd73SBarry Smith .seealso: MatSetValues(), MatSetValueLocal() 5821620fd73SBarry Smith M*/ 5831620fd73SBarry 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);} 5841620fd73SBarry Smith 5851620fd73SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatGetValue(Mat v,PetscInt i,PetscInt j,PetscScalar *va) {return MatGetValues(v,1,&i,1,&j,va);} 5861620fd73SBarry Smith 5871620fd73SBarry 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);} 5881620fd73SBarry Smith 5891620fd73SBarry Smith /*MC 5900d2a7ff2SSatish Balay MatPreallocateInitialize - Begins the block of code that will count the number of nonzeros per 591bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 592bd481603SBarry Smith 593bd481603SBarry Smith Synopsis: 594c1ac3661SBarry Smith PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz) 595bd481603SBarry Smith 596bd481603SBarry Smith Collective on MPI_Comm 597bd481603SBarry Smith 598bd481603SBarry Smith Input Parameters: 599bd481603SBarry Smith + comm - the communicator that will share the eventually allocated matrix 600859fcb39SBarry Smith . nrows - the number of LOCAL rows in the matrix 601859fcb39SBarry Smith - ncols - the number of LOCAL columns in the matrix 602bd481603SBarry Smith 603bd481603SBarry Smith Output Parameters: 604bd481603SBarry Smith + dnz - the array that will be passed to the matrix preallocation routines 605bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 606bd481603SBarry Smith 607bd481603SBarry Smith 608bd481603SBarry Smith Level: intermediate 609bd481603SBarry Smith 610bd481603SBarry Smith Notes: 611bd481603SBarry Smith See the chapter in the users manual on performance for more details 612bd481603SBarry Smith 6131d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 614bd481603SBarry Smith 615bd481603SBarry Smith Use MatPreallocateInitializeSymmetric() for symmetric matrices (MPISBAIJ matrices) 616bd481603SBarry Smith 6171620fd73SBarry Smith This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize(). 6181620fd73SBarry Smith 619bd481603SBarry Smith Concepts: preallocation^Matrix 620bd481603SBarry Smith 621bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(), 622bd481603SBarry Smith MatPreallocateInitializeSymmetric(), MatPreallocateSymmetricSetLocal() 623bd481603SBarry Smith M*/ 624c4f061fbSSatish Balay #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \ 6257c922b88SBarry Smith { \ 626c1ac3661SBarry Smith PetscErrorCode _4_ierr; PetscInt __tmp = (nrows),__ctmp = (ncols),__rstart,__start,__end; \ 627c1ac3661SBarry Smith _4_ierr = PetscMalloc(2*__tmp*sizeof(PetscInt),&dnz);CHKERRQ(_4_ierr);onz = dnz + __tmp;\ 628c1ac3661SBarry Smith _4_ierr = PetscMemzero(dnz,2*__tmp*sizeof(PetscInt));CHKERRQ(_4_ierr);\ 629c1ac3661SBarry Smith _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ctmp;\ 630c1ac3661SBarry Smith _4_ierr = MPI_Scan(&__tmp,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __tmp; 6317c922b88SBarry Smith 632bd481603SBarry Smith /*MC 6330d2a7ff2SSatish Balay MatPreallocateSymmetricInitialize - Begins the block of code that will count the number of nonzeros per 634bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 635bd481603SBarry Smith 636bd481603SBarry Smith Synopsis: 637c1ac3661SBarry Smith PetscErrorCode MatPreallocateSymmetricInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz) 638bd481603SBarry Smith 639bd481603SBarry Smith Collective on MPI_Comm 640bd481603SBarry Smith 641bd481603SBarry Smith Input Parameters: 642bd481603SBarry Smith + comm - the communicator that will share the eventually allocated matrix 643859fcb39SBarry Smith . nrows - the number of LOCAL rows in the matrix 644859fcb39SBarry Smith - ncols - the number of LOCAL columns in the matrix 645bd481603SBarry Smith 646bd481603SBarry Smith Output Parameters: 647bd481603SBarry Smith + dnz - the array that will be passed to the matrix preallocation routines 648bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 649bd481603SBarry Smith 650bd481603SBarry Smith 651bd481603SBarry Smith Level: intermediate 652bd481603SBarry Smith 653bd481603SBarry Smith Notes: 654bd481603SBarry Smith See the chapter in the users manual on performance for more details 655bd481603SBarry Smith 6561d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 657bd481603SBarry Smith 6581620fd73SBarry Smith This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize(). 6591620fd73SBarry Smith 660bd481603SBarry Smith Concepts: preallocation^Matrix 661bd481603SBarry Smith 662bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(), 663bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal() 664bd481603SBarry Smith M*/ 665222b16d4SBarry Smith #define MatPreallocateSymmetricInitialize(comm,nrows,ncols,dnz,onz) 0; \ 666222b16d4SBarry Smith { \ 667c1ac3661SBarry Smith PetscErrorCode _4_ierr; PetscInt __tmp = (nrows),__ctmp = (ncols),__rstart,__end; \ 668c1ac3661SBarry Smith _4_ierr = PetscMalloc(2*__tmp*sizeof(PetscInt),&dnz);CHKERRQ(_4_ierr);onz = dnz + __tmp;\ 669c1ac3661SBarry Smith _4_ierr = PetscMemzero(dnz,2*__tmp*sizeof(PetscInt));CHKERRQ(_4_ierr);\ 670c1ac3661SBarry Smith _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr);\ 671c1ac3661SBarry Smith _4_ierr = MPI_Scan(&__tmp,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __tmp; 672222b16d4SBarry Smith 673bd481603SBarry Smith /*MC 674bd481603SBarry Smith MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be 675bd481603SBarry Smith inserted using a local number of the rows and columns 676bd481603SBarry Smith 677bd481603SBarry Smith Synopsis: 678c1ac3661SBarry Smith PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 679bd481603SBarry Smith 680bd481603SBarry Smith Not Collective 681bd481603SBarry Smith 682bd481603SBarry Smith Input Parameters: 683bd481603SBarry Smith + map - the mapping between local numbering and global numbering 684bd481603SBarry Smith . nrows - the number of rows indicated 6851d73ed98SBarry Smith . rows - the indices of the rows 686bd481603SBarry Smith . ncols - the number of columns in the matrix 687bd481603SBarry Smith . cols - the columns indicated 688bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 689bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 690bd481603SBarry Smith 691bd481603SBarry Smith 692bd481603SBarry Smith Level: intermediate 693bd481603SBarry Smith 694bd481603SBarry Smith Notes: 695bd481603SBarry Smith See the chapter in the users manual on performance for more details 696bd481603SBarry Smith 6971d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 698bd481603SBarry Smith 699bd481603SBarry Smith Concepts: preallocation^Matrix 700bd481603SBarry Smith 701bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(), 702bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal() 703bd481603SBarry Smith M*/ 704c4f061fbSSatish Balay #define MatPreallocateSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\ 705c4f061fbSSatish Balay {\ 706c1ac3661SBarry Smith PetscInt __l;\ 707ef66eb69SBarry Smith _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\ 708ef66eb69SBarry Smith _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\ 709c4f061fbSSatish Balay for (__l=0;__l<nrows;__l++) {\ 710ef66eb69SBarry Smith _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 711c4f061fbSSatish Balay }\ 712c4f061fbSSatish Balay } 713c4f061fbSSatish Balay 714bd481603SBarry Smith /*MC 715bd481603SBarry Smith MatPreallocateSymmetricSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be 716bd481603SBarry Smith inserted using a local number of the rows and columns 717bd481603SBarry Smith 718bd481603SBarry Smith Synopsis: 719c1ac3661SBarry Smith PetscErrorCode MatPreallocateSymmetricSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 720bd481603SBarry Smith 721bd481603SBarry Smith Not Collective 722bd481603SBarry Smith 723bd481603SBarry Smith Input Parameters: 724bd481603SBarry Smith + map - the mapping between local numbering and global numbering 725bd481603SBarry Smith . nrows - the number of rows indicated 7261d73ed98SBarry Smith . rows - the indices of the rows 727bd481603SBarry Smith . ncols - the number of columns in the matrix 728bd481603SBarry Smith . cols - the columns indicated 729bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 730bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 731bd481603SBarry Smith 732bd481603SBarry Smith 733bd481603SBarry Smith Level: intermediate 734bd481603SBarry Smith 735bd481603SBarry Smith Notes: 736bd481603SBarry Smith See the chapter in the users manual on performance for more details 737bd481603SBarry Smith 738bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 739bd481603SBarry Smith 740bd481603SBarry Smith Concepts: preallocation^Matrix 741bd481603SBarry Smith 742bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(), 743bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal() 744bd481603SBarry Smith M*/ 745d3d32019SBarry Smith #define MatPreallocateSymmetricSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\ 746d3d32019SBarry Smith {\ 747c1ac3661SBarry Smith PetscInt __l;\ 748d3d32019SBarry Smith _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\ 749d3d32019SBarry Smith _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\ 750d3d32019SBarry Smith for (__l=0;__l<nrows;__l++) {\ 751d3d32019SBarry Smith _4_ierr = MatPreallocateSymmetricSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 752d3d32019SBarry Smith }\ 753d3d32019SBarry Smith } 754d3d32019SBarry Smith 755bd481603SBarry Smith /*MC 756bd481603SBarry Smith MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be 757bd481603SBarry Smith inserted using a local number of the rows and columns 758bd481603SBarry Smith 759bd481603SBarry Smith Synopsis: 760c1ac3661SBarry Smith PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 761bd481603SBarry Smith 762bd481603SBarry Smith Not Collective 763bd481603SBarry Smith 764bd481603SBarry Smith Input Parameters: 76564075487SBarry Smith + row - the row 766bd481603SBarry Smith . ncols - the number of columns in the matrix 767a50f8bf6SHong Zhang - cols - the columns indicated 768a50f8bf6SHong Zhang 769a50f8bf6SHong Zhang Output Parameters: 770a50f8bf6SHong Zhang + dnz - the array that will be passed to the matrix preallocation routines 771bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 772bd481603SBarry Smith 773bd481603SBarry Smith 774bd481603SBarry Smith Level: intermediate 775bd481603SBarry Smith 776bd481603SBarry Smith Notes: 777bd481603SBarry Smith See the chapter in the users manual on performance for more details 778bd481603SBarry Smith 779bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 780bd481603SBarry Smith 7811620fd73SBarry Smith This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize(). 7821620fd73SBarry Smith 783bd481603SBarry Smith Concepts: preallocation^Matrix 784bd481603SBarry Smith 785bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(), 786bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal() 787bd481603SBarry Smith M*/ 788c4f061fbSSatish Balay #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\ 789c1ac3661SBarry Smith { PetscInt __i; \ 7908ee2e534SBarry Smith if (row < __rstart) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Trying to set preallocation for row %D less than first local row %D",row,__rstart);\ 7918ee2e534SBarry Smith if (row >= __rstart+__tmp) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Trying to set preallocation for row %D greater than last local row %D",row,__rstart+__tmp-1);\ 7927c922b88SBarry Smith for (__i=0; __i<nc; __i++) {\ 79364075487SBarry Smith if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \ 7947cc688d7SBarry Smith else dnz[row - __rstart]++;\ 7957c922b88SBarry Smith }\ 7967c922b88SBarry Smith } 7977c922b88SBarry Smith 798bd481603SBarry Smith /*MC 799bd481603SBarry Smith MatPreallocateSymmetricSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be 800bd481603SBarry Smith inserted using a local number of the rows and columns 801bd481603SBarry Smith 802bd481603SBarry Smith Synopsis: 803c1ac3661SBarry Smith PetscErrorCode MatPreallocateSymmetricSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 804bd481603SBarry Smith 805bd481603SBarry Smith Not Collective 806bd481603SBarry Smith 807bd481603SBarry Smith Input Parameters: 808bd481603SBarry Smith + nrows - the number of rows indicated 8091d73ed98SBarry Smith . rows - the indices of the rows 810bd481603SBarry Smith . ncols - the number of columns in the matrix 811bd481603SBarry Smith . cols - the columns indicated 812bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 813bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 814bd481603SBarry Smith 815bd481603SBarry Smith 816bd481603SBarry Smith Level: intermediate 817bd481603SBarry Smith 818bd481603SBarry Smith Notes: 819bd481603SBarry Smith See the chapter in the users manual on performance for more details 820bd481603SBarry Smith 821bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 822bd481603SBarry Smith 8231620fd73SBarry Smith This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize(). 8241620fd73SBarry Smith 825bd481603SBarry Smith Concepts: preallocation^Matrix 826bd481603SBarry Smith 827bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(), 828bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal() 829bd481603SBarry Smith M*/ 830d3d32019SBarry Smith #define MatPreallocateSymmetricSet(row,nc,cols,dnz,onz) 0;\ 831c1ac3661SBarry Smith { PetscInt __i; \ 832d3d32019SBarry Smith for (__i=0; __i<nc; __i++) {\ 833d3d32019SBarry Smith if (cols[__i] >= __end) onz[row - __rstart]++; \ 834d3d32019SBarry Smith else if (cols[__i] >= row) dnz[row - __rstart]++;\ 835d3d32019SBarry Smith }\ 836d3d32019SBarry Smith } 837d3d32019SBarry Smith 838bd481603SBarry Smith /*MC 8390d2a7ff2SSatish Balay MatPreallocateFinalize - Ends the block of code that will count the number of nonzeros per 840bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 841bd481603SBarry Smith 842bd481603SBarry Smith Synopsis: 843c1ac3661SBarry Smith PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz) 844bd481603SBarry Smith 845bd481603SBarry Smith Collective on MPI_Comm 846bd481603SBarry Smith 847bd481603SBarry Smith Input Parameters: 848bd481603SBarry Smith + dnz - the array that will be passed to the matrix preallocation routines 849bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 850bd481603SBarry Smith 851bd481603SBarry Smith 852bd481603SBarry Smith Level: intermediate 853bd481603SBarry Smith 854bd481603SBarry Smith Notes: 855bd481603SBarry Smith See the chapter in the users manual on performance for more details 856bd481603SBarry Smith 857bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 858bd481603SBarry Smith 8591620fd73SBarry Smith This is a MACRO not a function because it closes the { started in MatPreallocateInitialize(). 8601620fd73SBarry Smith 861bd481603SBarry Smith Concepts: preallocation^Matrix 862bd481603SBarry Smith 863bd481603SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(), 864bd481603SBarry Smith MatPreallocateSymmetricInitialize(), MatPreallocateSymmetricSetLocal() 865bd481603SBarry Smith M*/ 866ef66eb69SBarry Smith #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree(dnz);CHKERRQ(_4_ierr);} 8677c922b88SBarry Smith 868bd481603SBarry Smith 869bd481603SBarry Smith 8707b80b807SBarry Smith /* Routines unique to particular data structures */ 871be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellGetContext(Mat,void **); 872ba966639SSatish Balay PetscPolymorphicFunction(MatShellGetContext,(Mat A),(A,&t),void*,t) 873698d4c6aSKris Buschelman 874be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInodeAdjustForInodes(Mat,IS*,IS*); 875be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *); 876698d4c6aSKris Buschelman 877be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetColumnIndices(Mat,PetscInt[]); 878be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetColumnIndices(Mat,PetscInt[]); 879be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 880c75a6043SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 881c75a6043SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 8827b80b807SBarry Smith 883a96a251dSBarry Smith #define MAT_SKIP_ALLOCATION -4 884a96a251dSBarry Smith 885be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]); 886ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz)) 887be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]); 888ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqSBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz)) 889be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]); 890ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqAIJSetPreallocation,(Mat A,const PetscInt nnz[]),(A,0,nnz)) 891be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDensePreallocation(Mat,PetscScalar[]); 892273d9f13SBarry Smith 893be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 894ba966639SSatish Balay PetscPolymorphicSubroutine(MatMPIBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[],const PetscInt onz[]),(A,bs,0,nnz,0,onz)) 895be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 896be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 89753803f26SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []); 898be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]); 899be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 900be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]); 901be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIDenseSetPreallocation(Mat,PetscScalar[]); 902284134d9SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation(Mat,PetscScalar[]); 903be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIRowbsSetPreallocation(Mat,PetscInt,const PetscInt[]); 904be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,PetscInt*[]); 905be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,PetscInt*[]); 906be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAdicSetLocalFunction(Mat,void (*)(void)); 907273d9f13SBarry Smith 908be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetLDA(Mat,PetscInt); 909ab92ecdeSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDenseGetLocalMatrix(Mat,Mat*); 9101b807ce4Svictorle 911be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStoreValues(Mat); 912be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRetrieveValues(Mat); 9132e8a6d31SBarry Smith 914be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDAADSetCtx(Mat,void*); 9153a7fca6bSBarry Smith 9167b80b807SBarry Smith /* 9177b80b807SBarry Smith These routines are not usually accessed directly, rather solving is 91894b7f48cSBarry Smith done through the KSP and PC interfaces. 9197b80b807SBarry Smith */ 9207b80b807SBarry Smith 921d9274352SBarry Smith /*E 922d9274352SBarry Smith MatOrderingType - String with the name of a PETSc matrix ordering or the creation function 923d9274352SBarry Smith with an optional dynamic library name, for example 924d9274352SBarry Smith http://www.mcs.anl.gov/petsc/lib.a:orderingcreate() 925d9274352SBarry Smith 926d9274352SBarry Smith Level: beginner 927d9274352SBarry Smith 928e5a9bf91SBarry Smith Cannot use const because the PC objects manipulate the string 929e5a9bf91SBarry Smith 930d9274352SBarry Smith .seealso: MatGetOrdering() 931d9274352SBarry Smith E*/ 93249773a63SBarry Smith #define MatOrderingType char* 933b12f92e5SBarry Smith #define MATORDERING_NATURAL "natural" 934b12f92e5SBarry Smith #define MATORDERING_ND "nd" 935b12f92e5SBarry Smith #define MATORDERING_1WD "1wd" 936b12f92e5SBarry Smith #define MATORDERING_RCM "rcm" 937b12f92e5SBarry Smith #define MATORDERING_QMD "qmd" 938b12f92e5SBarry Smith #define MATORDERING_ROWLENGTH "rowlength" 93962152c8bSBarry Smith #define MATORDERING_DSC_ND "dsc_nd" 94062152c8bSBarry Smith #define MATORDERING_DSC_MMD "dsc_mmd" 94162152c8bSBarry Smith #define MATORDERING_DSC_MDF "dsc_mdf" 942c06d978dSMatthew Knepley #define MATORDERING_CONSTRAINED "constrained" 943c06d978dSMatthew Knepley #define MATORDERING_IDENTITY "identity" 944c06d978dSMatthew Knepley #define MATORDERING_REVERSE "reverse" 945b12f92e5SBarry Smith 946be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOrdering(Mat,const MatOrderingType,IS*,IS*); 947be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOrderingList(PetscFList *list); 948be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,const MatOrderingType,IS*,IS*)); 94930de9b25SBarry Smith 95030de9b25SBarry Smith /*MC 95130de9b25SBarry Smith MatOrderingRegisterDynamic - Adds a new sparse matrix ordering to the 95230de9b25SBarry Smith matrix package. 95330de9b25SBarry Smith 95430de9b25SBarry Smith Synopsis: 955c1ac3661SBarry Smith PetscErrorCode MatOrderingRegisterDynamic(char *name_ordering,char *path,char *name_create,PetscErrorCode (*routine_create)(MatOrdering)) 95630de9b25SBarry Smith 95730de9b25SBarry Smith Not Collective 95830de9b25SBarry Smith 95930de9b25SBarry Smith Input Parameters: 96030de9b25SBarry Smith + sname - name of ordering (for example MATORDERING_ND) 96130de9b25SBarry Smith . path - location of library where creation routine is 96230de9b25SBarry Smith . name - name of function that creates the ordering type,a string 96330de9b25SBarry Smith - function - function pointer that creates the ordering 96430de9b25SBarry Smith 96530de9b25SBarry Smith Level: developer 96630de9b25SBarry Smith 96730de9b25SBarry Smith If dynamic libraries are used, then the fourth input argument (function) 96830de9b25SBarry Smith is ignored. 96930de9b25SBarry Smith 97030de9b25SBarry Smith Sample usage: 97130de9b25SBarry Smith .vb 97230de9b25SBarry Smith MatOrderingRegisterDynamic("my_order",/home/username/my_lib/lib/libO/solaris/mylib.a, 97330de9b25SBarry Smith "MyOrder",MyOrder); 97430de9b25SBarry Smith .ve 97530de9b25SBarry Smith 97630de9b25SBarry Smith Then, your partitioner can be chosen with the procedural interface via 97730de9b25SBarry Smith $ MatOrderingSetType(part,"my_order) 97830de9b25SBarry Smith or at runtime via the option 9792401956bSBarry Smith $ -pc_factor_mat_ordering_type my_order 98030de9b25SBarry Smith 981ab901514SBarry Smith ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values. 98230de9b25SBarry Smith 98330de9b25SBarry Smith .keywords: matrix, ordering, register 98430de9b25SBarry Smith 98530de9b25SBarry Smith .seealso: MatOrderingRegisterDestroy(), MatOrderingRegisterAll() 98630de9b25SBarry Smith M*/ 987aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 988f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,0) 989b12f92e5SBarry Smith #else 990f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,d) 991b12f92e5SBarry Smith #endif 99230de9b25SBarry Smith 993be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegisterDestroy(void); 994be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegisterAll(const char[]); 9952bad1931SBarry Smith extern PetscTruth MatOrderingRegisterAllCalled; 996b0a32e0cSBarry Smith extern PetscFList MatOrderingList; 997d4fbbf0eSBarry Smith 998be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS); 999a2ce50c7SBarry Smith 1000d91e6319SBarry Smith /*S 10012401956bSBarry Smith MatFactorInfo - Data passed into the matrix factorization routines 1002b00f7748SHong Zhang 100361cad860SBarry Smith In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE, that is use 100461cad860SBarry Smith $ MatFactorInfo info(MAT_FACTORINFO_SIZE) 1005b00f7748SHong Zhang 100615e8a5b3SHong Zhang Notes: These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC. 1007b00f7748SHong Zhang 100861cad860SBarry Smith You can use MatFactorInfoInitialize() to set default values. 100961cad860SBarry Smith 1010b00f7748SHong Zhang Level: developer 1011b00f7748SHong Zhang 1012d7d82daaSBarry Smith .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(), 1013d7d82daaSBarry Smith MatFactorInfoInitialize() 1014b00f7748SHong Zhang 1015b00f7748SHong Zhang S*/ 1016b00f7748SHong Zhang typedef struct { 10170a29876aSHong Zhang PetscReal shiftnz; /* scaling of identity added to matrix to prevent zero pivots */ 1018fbf22428SSatish Balay PetscReal shiftpd; /* if true, shift until positive pivots */ 10192cea7109SBarry Smith PetscReal shift_fraction; /* record shift fraction taken */ 102015e8a5b3SHong Zhang PetscReal diagonal_fill; /* force diagonal to fill in if initially not filled */ 102115e8a5b3SHong Zhang PetscReal dt; /* drop tolerance */ 1022b00f7748SHong Zhang PetscReal dtcol; /* tolerance for pivoting */ 102315e8a5b3SHong Zhang PetscReal dtcount; /* maximum nonzeros to be allowed per row */ 102467e28bfeSBarry Smith PetscReal fill; /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */ 1025348344bbSBarry Smith PetscReal levels; /* ICC/ILU(levels) */ 1026bcd9e38bSBarry Smith PetscReal pivotinblocks; /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0 1027bcd9e38bSBarry Smith factorization may be faster if do not pivot */ 102862bba022SBarry Smith PetscReal shiftinblocks; /* if block in block factorization has zero pivot then shift diagonal until non-singular */ 102915e8a5b3SHong Zhang PetscReal zeropivot; /* pivot is called zero if less than this */ 103015e8a5b3SHong Zhang } MatFactorInfo; 1031ffa6d0a5SLois Curfman McInnes 1032be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFactorInfoInitialize(MatFactorInfo*); 1033be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactor(Mat,IS,MatFactorInfo*); 1034be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactorSymbolic(Mat,IS,MatFactorInfo*,Mat*); 1035be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactorNumeric(Mat,MatFactorInfo*,Mat*); 1036be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactor(Mat,IS,IS,MatFactorInfo*); 1037be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUFactor(Mat,IS,IS,MatFactorInfo*); 1038be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactorSymbolic(Mat,IS,IS,MatFactorInfo*,Mat*); 1039be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUFactorSymbolic(Mat,IS,IS,MatFactorInfo*,Mat*); 1040be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatICCFactorSymbolic(Mat,IS,MatFactorInfo*,Mat*); 1041be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatICCFactor(Mat,IS,MatFactorInfo*); 1042be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactorNumeric(Mat,MatFactorInfo*,Mat*); 1043be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUDTFactor(Mat,IS,IS,MatFactorInfo*,Mat *); 1044be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*); 1045be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolve(Mat,Vec,Vec); 1046be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatForwardSolve(Mat,Vec,Vec); 1047be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatBackwardSolve(Mat,Vec,Vec); 1048be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveAdd(Mat,Vec,Vec,Vec); 1049be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveTranspose(Mat,Vec,Vec); 1050be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveTransposeAdd(Mat,Vec,Vec,Vec); 1051be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolves(Mat,Vecs,Vecs); 10528ed539a5SBarry Smith 1053be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUnfactored(Mat); 1054bb5a7306SBarry Smith 1055d91e6319SBarry Smith /*E 1056d91e6319SBarry Smith MatSORType - What type of (S)SOR to perform 1057bb1eb677SSatish Balay 1058d91e6319SBarry Smith Level: beginner 1059d91e6319SBarry Smith 1060d9274352SBarry Smith May be bitwise ORd together 1061d9274352SBarry Smith 1062d91e6319SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 1063d91e6319SBarry Smith 10644e7234bfSBarry Smith MatSORType may be bitwise ORd together, so do not change the numbers 10654e7234bfSBarry Smith 1066a1d92eedSBarry Smith .seealso: MatRelax(), MatPBRelax() 1067d91e6319SBarry Smith E*/ 1068ee50ffe9SBarry Smith typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3, 1069ee50ffe9SBarry Smith SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8, 1070ee50ffe9SBarry Smith SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16, 107184cb2905SBarry Smith SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType; 1072be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRelax(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec); 1073be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPBRelax(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec); 10748ed539a5SBarry Smith 1075d4fbbf0eSBarry Smith /* 1076639f9d9dSBarry Smith These routines are for efficiently computing Jacobians via finite differences. 1077639f9d9dSBarry Smith */ 1078b12f92e5SBarry Smith 1079d9274352SBarry Smith /*E 1080d9274352SBarry Smith MatColoringType - String with the name of a PETSc matrix coloring or the creation function 1081d9274352SBarry Smith with an optional dynamic library name, for example 1082d9274352SBarry Smith http://www.mcs.anl.gov/petsc/lib.a:coloringcreate() 1083d9274352SBarry Smith 1084d9274352SBarry Smith Level: beginner 1085d9274352SBarry Smith 1086d9274352SBarry Smith .seealso: MatGetColoring() 1087d9274352SBarry Smith E*/ 1088e5a9bf91SBarry Smith #define MatColoringType const char* 1089b12f92e5SBarry Smith #define MATCOLORING_NATURAL "natural" 1090b12f92e5SBarry Smith #define MATCOLORING_SL "sl" 1091b12f92e5SBarry Smith #define MATCOLORING_LF "lf" 1092b12f92e5SBarry Smith #define MATCOLORING_ID "id" 1093b12f92e5SBarry Smith 10942e90c967SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColoring(Mat,MatColoringType,ISColoring*); 10952e90c967SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,MatColoringType,ISColoring *)); 109630de9b25SBarry Smith 109730de9b25SBarry Smith /*MC 109830de9b25SBarry Smith MatColoringRegisterDynamic - Adds a new sparse matrix coloring to the 109930de9b25SBarry Smith matrix package. 110030de9b25SBarry Smith 110130de9b25SBarry Smith Synopsis: 1102c1ac3661SBarry Smith PetscErrorCode MatColoringRegisterDynamic(char *name_coloring,char *path,char *name_create,PetscErrorCode (*routine_create)(MatColoring)) 110330de9b25SBarry Smith 110430de9b25SBarry Smith Not Collective 110530de9b25SBarry Smith 110630de9b25SBarry Smith Input Parameters: 110730de9b25SBarry Smith + sname - name of Coloring (for example MATCOLORING_SL) 110830de9b25SBarry Smith . path - location of library where creation routine is 110930de9b25SBarry Smith . name - name of function that creates the Coloring type, a string 111030de9b25SBarry Smith - function - function pointer that creates the coloring 111130de9b25SBarry Smith 111230de9b25SBarry Smith Level: developer 111330de9b25SBarry Smith 111430de9b25SBarry Smith If dynamic libraries are used, then the fourth input argument (function) 111530de9b25SBarry Smith is ignored. 111630de9b25SBarry Smith 111730de9b25SBarry Smith Sample usage: 111830de9b25SBarry Smith .vb 111930de9b25SBarry Smith MatColoringRegisterDynamic("my_color",/home/username/my_lib/lib/libO/solaris/mylib.a, 112030de9b25SBarry Smith "MyColor",MyColor); 112130de9b25SBarry Smith .ve 112230de9b25SBarry Smith 112330de9b25SBarry Smith Then, your partitioner can be chosen with the procedural interface via 112430de9b25SBarry Smith $ MatColoringSetType(part,"my_color") 112530de9b25SBarry Smith or at runtime via the option 112630de9b25SBarry Smith $ -mat_coloring_type my_color 112730de9b25SBarry Smith 1128ab901514SBarry Smith $PETSC_ARCH occuring in pathname will be replaced with appropriate values. 112930de9b25SBarry Smith 113030de9b25SBarry Smith .keywords: matrix, Coloring, register 113130de9b25SBarry Smith 113230de9b25SBarry Smith .seealso: MatColoringRegisterDestroy(), MatColoringRegisterAll() 113330de9b25SBarry Smith M*/ 1134aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 1135f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,0) 1136b12f92e5SBarry Smith #else 1137f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,d) 1138b12f92e5SBarry Smith #endif 113930de9b25SBarry Smith 11402bad1931SBarry Smith extern PetscTruth MatColoringRegisterAllCalled; 1141f1144a30SSatish Balay 1142be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegisterAll(const char[]); 1143be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegisterDestroy(void); 1144be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*); 1145639f9d9dSBarry Smith 1146d9274352SBarry Smith /*S 1147d9274352SBarry Smith MatFDColoring - Object for computing a sparse Jacobian via finite differences 1148d9274352SBarry Smith and coloring 1149639f9d9dSBarry Smith 1150d9274352SBarry Smith Level: beginner 1151d9274352SBarry Smith 1152d9274352SBarry Smith Concepts: coloring, sparse Jacobian, finite differences 1153d9274352SBarry Smith 1154d9274352SBarry Smith .seealso: MatFDColoringCreate() 1155d9274352SBarry Smith S*/ 1156e2a1c21fSSatish Balay typedef struct _p_MatFDColoring* MatFDColoring; 1157639f9d9dSBarry Smith 1158be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringCreate(Mat,ISColoring,MatFDColoring *); 1159be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringDestroy(MatFDColoring); 1160be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringView(MatFDColoring,PetscViewer); 1161be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*); 11624a9d489dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**); 1163be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal); 1164be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFrequency(MatFDColoring,PetscInt); 1165be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetFrequency(MatFDColoring,PetscInt*); 1166be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFromOptions(MatFDColoring); 1167be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApply(Mat,MatFDColoring,Vec,MatStructure*,void *); 1168be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApplyTS(Mat,MatFDColoring,PetscReal,Vec,MatStructure*,void *); 1169be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetRecompute(MatFDColoring); 1170be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetF(MatFDColoring,Vec); 1171be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,PetscInt*[]); 1172639f9d9dSBarry Smith /* 11730752156aSBarry Smith These routines are for partitioning matrices: currently used only 11743eda8832SBarry Smith for adjacency matrix, MatCreateMPIAdj(). 11750752156aSBarry Smith */ 1176ca161407SBarry Smith 1177d9274352SBarry Smith /*S 1178d9274352SBarry Smith MatPartitioning - Object for managing the partitioning of a matrix or graph 1179d9274352SBarry Smith 1180d9274352SBarry Smith Level: beginner 1181d9274352SBarry Smith 1182d9274352SBarry Smith Concepts: partitioning 1183d9274352SBarry Smith 1184743a1026Svictorle .seealso: MatPartitioningCreate(), MatPartitioningType 1185d9274352SBarry Smith S*/ 118691e9ee9fSBarry Smith typedef struct _p_MatPartitioning* MatPartitioning; 1187d9274352SBarry Smith 1188d9274352SBarry Smith /*E 11895bc6e7ccSvictorle MatPartitioningType - String with the name of a PETSc matrix partitioning or the creation function 1190d9274352SBarry Smith with an optional dynamic library name, for example 1191d9274352SBarry Smith http://www.mcs.anl.gov/petsc/lib.a:partitioningcreate() 1192d9274352SBarry Smith 1193d9274352SBarry Smith Level: beginner 1194d9274352SBarry Smith 1195b547a13bSvictorle .seealso: MatPartitioningCreate(), MatPartitioning 1196d9274352SBarry Smith E*/ 1197e5a9bf91SBarry Smith #define MatPartitioningType const char* 11988ba1e511SMatthew Knepley #define MAT_PARTITIONING_CURRENT "current" 1199c2ec2a73SHong Zhang #define MAT_PARTITIONING_SQUARE "square" 12008ba1e511SMatthew Knepley #define MAT_PARTITIONING_PARMETIS "parmetis" 120171306c60Sjroman #define MAT_PARTITIONING_CHACO "chaco" 120271306c60Sjroman #define MAT_PARTITIONING_JOSTLE "jostle" 120371306c60Sjroman #define MAT_PARTITIONING_PARTY "party" 120471306c60Sjroman #define MAT_PARTITIONING_SCOTCH "scotch" 120571306c60Sjroman 1206ca161407SBarry Smith 1207be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningCreate(MPI_Comm,MatPartitioning*); 12082e90c967SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetType(MatPartitioning,MatPartitioningType); 1209be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetNParts(MatPartitioning,PetscInt); 1210be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetAdjacency(MatPartitioning,Mat); 1211be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]); 1212be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []); 1213be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningApply(MatPartitioning,IS*); 1214be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningDestroy(MatPartitioning); 12152aabb6bbSBarry Smith 1216be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatPartitioning)); 121730de9b25SBarry Smith 121830de9b25SBarry Smith /*MC 121930de9b25SBarry Smith MatPartitioningRegisterDynamic - Adds a new sparse matrix partitioning to the 122030de9b25SBarry Smith matrix package. 122130de9b25SBarry Smith 122230de9b25SBarry Smith Synopsis: 1223c1ac3661SBarry Smith PetscErrorCode MatPartitioningRegisterDynamic(char *name_partitioning,char *path,char *name_create,PetscErrorCode (*routine_create)(MatPartitioning)) 122430de9b25SBarry Smith 122530de9b25SBarry Smith Not Collective 122630de9b25SBarry Smith 122730de9b25SBarry Smith Input Parameters: 122830de9b25SBarry Smith + sname - name of partitioning (for example MAT_PARTITIONING_CURRENT) or parmetis 122930de9b25SBarry Smith . path - location of library where creation routine is 123030de9b25SBarry Smith . name - name of function that creates the partitioning type, a string 123130de9b25SBarry Smith - function - function pointer that creates the partitioning type 123230de9b25SBarry Smith 123330de9b25SBarry Smith Level: developer 123430de9b25SBarry Smith 123530de9b25SBarry Smith If dynamic libraries are used, then the fourth input argument (function) 123630de9b25SBarry Smith is ignored. 123730de9b25SBarry Smith 123830de9b25SBarry Smith Sample usage: 123930de9b25SBarry Smith .vb 124030de9b25SBarry Smith MatPartitioningRegisterDynamic("my_part",/home/username/my_lib/lib/libO/solaris/mylib.a, 124130de9b25SBarry Smith "MyPartCreate",MyPartCreate); 124230de9b25SBarry Smith .ve 124330de9b25SBarry Smith 124430de9b25SBarry Smith Then, your partitioner can be chosen with the procedural interface via 124530de9b25SBarry Smith $ MatPartitioningSetType(part,"my_part") 124630de9b25SBarry Smith or at runtime via the option 124730de9b25SBarry Smith $ -mat_partitioning_type my_part 124830de9b25SBarry Smith 1249ab901514SBarry Smith $PETSC_ARCH occuring in pathname will be replaced with appropriate values. 125030de9b25SBarry Smith 125130de9b25SBarry Smith .keywords: matrix, partitioning, register 125230de9b25SBarry Smith 125330de9b25SBarry Smith .seealso: MatPartitioningRegisterDestroy(), MatPartitioningRegisterAll() 125430de9b25SBarry Smith M*/ 1255aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 1256f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,0) 12572aabb6bbSBarry Smith #else 1258f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,d) 12592aabb6bbSBarry Smith #endif 12602aabb6bbSBarry Smith 12612bad1931SBarry Smith extern PetscTruth MatPartitioningRegisterAllCalled; 1262f1144a30SSatish Balay 1263be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegisterAll(const char[]); 1264be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegisterDestroy(void); 12652bad1931SBarry Smith 1266be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningView(MatPartitioning,PetscViewer); 1267be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetFromOptions(MatPartitioning); 1268be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningGetType(MatPartitioning,MatPartitioningType*); 1269ca161407SBarry Smith 1270be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningParmetisSetCoarseSequential(MatPartitioning); 12710752156aSBarry Smith 1272be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningJostleSetCoarseLevel(MatPartitioning,PetscReal); 1273be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningJostleSetCoarseSequential(MatPartitioning); 127471306c60Sjroman 1275954b3ec6SBarry Smith typedef enum {MP_CHACO_MULTILEVEL_KL,MP_CHACO_SPECTRAL,MP_CHACO_LINEAR,MP_CHACO_RANDOM, MP_CHACO_SCATTERED} MPChacoGlobalType; 1276be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetGlobal(MatPartitioning, MPChacoGlobalType); 127771306c60Sjroman typedef enum { MP_CHACO_KERNIGHAN_LIN, MP_CHACO_NONE } MPChacoLocalType; 1278be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetLocal(MatPartitioning, MPChacoLocalType); 1279be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal); 128071306c60Sjroman typedef enum { MP_CHACO_LANCZOS, MP_CHACO_RQI_SYMMLQ } MPChacoEigenType; 1281be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType); 1282be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenTol(MatPartitioning, PetscReal); 1283be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenNumber(MatPartitioning, PetscInt); 128471306c60Sjroman 128571306c60Sjroman #define MP_PARTY_OPT "opt" 128671306c60Sjroman #define MP_PARTY_LIN "lin" 128771306c60Sjroman #define MP_PARTY_SCA "sca" 128871306c60Sjroman #define MP_PARTY_RAN "ran" 128971306c60Sjroman #define MP_PARTY_GBF "gbf" 129071306c60Sjroman #define MP_PARTY_GCF "gcf" 129171306c60Sjroman #define MP_PARTY_BUB "bub" 129271306c60Sjroman #define MP_PARTY_DEF "def" 1293be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetGlobal(MatPartitioning, const char*); 129471306c60Sjroman #define MP_PARTY_HELPFUL_SETS "hs" 129571306c60Sjroman #define MP_PARTY_KERNIGHAN_LIN "kl" 129671306c60Sjroman #define MP_PARTY_NONE "no" 1297be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetLocal(MatPartitioning, const char*); 1298be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal); 1299be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetBipart(MatPartitioning,PetscTruth); 1300be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscTruth); 130171306c60Sjroman 130271306c60Sjroman typedef enum { MP_SCOTCH_GREEDY, MP_SCOTCH_GPS, MP_SCOTCH_GR_GPS } MPScotchGlobalType; 1303be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetArch(MatPartitioning,const char*); 1304be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetMultilevel(MatPartitioning); 1305be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetGlobal(MatPartitioning,MPScotchGlobalType); 1306be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetCoarseLevel(MatPartitioning,PetscReal); 1307be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetHostList(MatPartitioning,const char*); 130871306c60Sjroman typedef enum { MP_SCOTCH_KERNIGHAN_LIN, MP_SCOTCH_NONE } MPScotchLocalType; 1309be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetLocal(MatPartitioning,MPScotchLocalType); 1310be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetMapping(MatPartitioning); 1311be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetStrategy(MatPartitioning,char*); 131271306c60Sjroman 13130752156aSBarry Smith /* 13140a835dfdSSatish Balay If you add entries here you must also add them to finclude/petscmat.h 1315d4fbbf0eSBarry Smith */ 13161c1c02c0SLois Curfman McInnes typedef enum { MATOP_SET_VALUES=0, 13171c1c02c0SLois Curfman McInnes MATOP_GET_ROW=1, 13181c1c02c0SLois Curfman McInnes MATOP_RESTORE_ROW=2, 13191c1c02c0SLois Curfman McInnes MATOP_MULT=3, 13201c1c02c0SLois Curfman McInnes MATOP_MULT_ADD=4, 13217c922b88SBarry Smith MATOP_MULT_TRANSPOSE=5, 13227c922b88SBarry Smith MATOP_MULT_TRANSPOSE_ADD=6, 13231c1c02c0SLois Curfman McInnes MATOP_SOLVE=7, 13241c1c02c0SLois Curfman McInnes MATOP_SOLVE_ADD=8, 13257c922b88SBarry Smith MATOP_SOLVE_TRANSPOSE=9, 13267c922b88SBarry Smith MATOP_SOLVE_TRANSPOSE_ADD=10, 13271c1c02c0SLois Curfman McInnes MATOP_LUFACTOR=11, 13281c1c02c0SLois Curfman McInnes MATOP_CHOLESKYFACTOR=12, 13291c1c02c0SLois Curfman McInnes MATOP_RELAX=13, 13301c1c02c0SLois Curfman McInnes MATOP_TRANSPOSE=14, 13311c1c02c0SLois Curfman McInnes MATOP_GETINFO=15, 13321c1c02c0SLois Curfman McInnes MATOP_EQUAL=16, 13331c1c02c0SLois Curfman McInnes MATOP_GET_DIAGONAL=17, 13341c1c02c0SLois Curfman McInnes MATOP_DIAGONAL_SCALE=18, 13351c1c02c0SLois Curfman McInnes MATOP_NORM=19, 13361c1c02c0SLois Curfman McInnes MATOP_ASSEMBLY_BEGIN=20, 13371c1c02c0SLois Curfman McInnes MATOP_ASSEMBLY_END=21, 13381c1c02c0SLois Curfman McInnes MATOP_COMPRESS=22, 13391c1c02c0SLois Curfman McInnes MATOP_SET_OPTION=23, 13401c1c02c0SLois Curfman McInnes MATOP_ZERO_ENTRIES=24, 13411c1c02c0SLois Curfman McInnes MATOP_ZERO_ROWS=25, 13421c1c02c0SLois Curfman McInnes MATOP_LUFACTOR_SYMBOLIC=26, 13431c1c02c0SLois Curfman McInnes MATOP_LUFACTOR_NUMERIC=27, 13441c1c02c0SLois Curfman McInnes MATOP_CHOLESKY_FACTOR_SYMBOLIC=28, 13451c1c02c0SLois Curfman McInnes MATOP_CHOLESKY_FACTOR_NUMERIC=29, 1346d643ce63SMatthew Knepley MATOP_SETUP_PREALLOCATION=30, 1347d643ce63SMatthew Knepley MATOP_ILUFACTOR_SYMBOLIC=31, 1348d643ce63SMatthew Knepley MATOP_ICCFACTOR_SYMBOLIC=32, 1349d643ce63SMatthew Knepley MATOP_GET_ARRAY=33, 1350d643ce63SMatthew Knepley MATOP_RESTORE_ARRAY=34, 1351e26f1c3cSBarry Smith MATOP_DUPLICATE=35, 1352d643ce63SMatthew Knepley MATOP_FORWARD_SOLVE=36, 1353d643ce63SMatthew Knepley MATOP_BACKWARD_SOLVE=37, 1354d643ce63SMatthew Knepley MATOP_ILUFACTOR=38, 1355d643ce63SMatthew Knepley MATOP_ICCFACTOR=39, 1356d643ce63SMatthew Knepley MATOP_AXPY=40, 1357d643ce63SMatthew Knepley MATOP_GET_SUBMATRICES=41, 1358d643ce63SMatthew Knepley MATOP_INCREASE_OVERLAP=42, 1359d643ce63SMatthew Knepley MATOP_GET_VALUES=43, 1360d643ce63SMatthew Knepley MATOP_COPY=44, 1361d643ce63SMatthew Knepley MATOP_PRINT_HELP=45, 1362d643ce63SMatthew Knepley MATOP_SCALE=46, 1363d643ce63SMatthew Knepley MATOP_SHIFT=47, 1364d643ce63SMatthew Knepley MATOP_DIAGONAL_SHIFT=48, 1365d643ce63SMatthew Knepley MATOP_ILUDT_FACTOR=49, 1366d643ce63SMatthew Knepley MATOP_GET_BLOCK_SIZE=50, 1367d643ce63SMatthew Knepley MATOP_GET_ROW_IJ=51, 1368d643ce63SMatthew Knepley MATOP_RESTORE_ROW_IJ=52, 1369d643ce63SMatthew Knepley MATOP_GET_COLUMN_IJ=53, 1370d643ce63SMatthew Knepley MATOP_RESTORE_COLUMN_IJ=54, 1371d643ce63SMatthew Knepley MATOP_FDCOLORING_CREATE=55, 1372d643ce63SMatthew Knepley MATOP_COLORING_PATCH=56, 1373d643ce63SMatthew Knepley MATOP_SET_UNFACTORED=57, 1374d643ce63SMatthew Knepley MATOP_PERMUTE=58, 1375d643ce63SMatthew Knepley MATOP_SET_VALUES_BLOCKED=59, 1376d643ce63SMatthew Knepley MATOP_GET_SUBMATRIX=60, 1377d643ce63SMatthew Knepley MATOP_DESTROY=61, 1378d643ce63SMatthew Knepley MATOP_VIEW=62, 1379d643ce63SMatthew Knepley MATOP_GET_MAPS=63, 1380d643ce63SMatthew Knepley MATOP_USE_SCALED_FORM=64, 1381d643ce63SMatthew Knepley MATOP_SCALE_SYSTEM=65, 1382d643ce63SMatthew Knepley MATOP_UNSCALE_SYSTEM=66, 1383ba16a8cbSKris Buschelman MATOP_SET_LOCAL_TO_GLOBAL_MAP=67, 1384d643ce63SMatthew Knepley MATOP_SET_VALUES_LOCAL=68, 1385d643ce63SMatthew Knepley MATOP_ZERO_ROWS_LOCAL=69, 1386d643ce63SMatthew Knepley MATOP_GET_ROW_MAX=70, 1387d643ce63SMatthew Knepley MATOP_CONVERT=71, 1388d643ce63SMatthew Knepley MATOP_SET_COLORING=72, 1389d643ce63SMatthew Knepley MATOP_SET_VALUES_ADIC=73, 1390d643ce63SMatthew Knepley MATOP_SET_VALUES_ADIFOR=74, 1391d643ce63SMatthew Knepley MATOP_FD_COLORING_APPLY=75, 1392d643ce63SMatthew Knepley MATOP_SET_FROM_OPTIONS=76, 1393ba16a8cbSKris Buschelman MATOP_MULT_CON=77, 1394ba16a8cbSKris Buschelman MATOP_MULT_TRANSPOSE_CON=78, 1395ba16a8cbSKris Buschelman MATOP_ILU_FACTOR_SYMBOLIC_CON=79, 1396d643ce63SMatthew Knepley MATOP_PERMUTE_SPARSIFY=80, 1397d643ce63SMatthew Knepley MATOP_MULT_MULTIPLE=81, 139841acf15aSKris Buschelman MATOP_SOLVE_MULTIPLE=82, 139941acf15aSKris Buschelman MATOP_GET_INERTIA=83, 140041acf15aSKris Buschelman MATOP_LOAD=84, 140141acf15aSKris Buschelman MATOP_IS_SYMMETRIC=85, 140241acf15aSKris Buschelman MATOP_IS_HERMITIAN=86, 140341acf15aSKris Buschelman MATOP_IS_STRUCTURALLY_SYMMETRIC=87, 140441acf15aSKris Buschelman MATOP_PB_RELAX=88, 140541acf15aSKris Buschelman MATOP_GET_VECS=89, 140641acf15aSKris Buschelman MATOP_MAT_MULT=90, 140741acf15aSKris Buschelman MATOP_MAT_MULT_SYMBOLIC=91, 140841acf15aSKris Buschelman MATOP_MAT_MULT_NUMERIC=92, 140941acf15aSKris Buschelman MATOP_PTAP=93, 141041acf15aSKris Buschelman MATOP_PTAP_SYMBOLIC=94, 1411bc011b1eSHong Zhang MATOP_PTAP_NUMERIC=95, 1412bc011b1eSHong Zhang MATOP_MAT_MULTTRANSPOSE=96, 1413bc011b1eSHong Zhang MATOP_MAT_MULTTRANSPOSE_SYMBOLIC=97, 14147ba1a0bfSKris Buschelman MATOP_MAT_MULTTRANSPOSE_NUMERIC=98, 14157ba1a0bfSKris Buschelman MATOP_PTAP_SYMBOLIC_SEQAIJ=99, 14167ba1a0bfSKris Buschelman MATOP_PTAP_NUMERIC_SEQAIJ=100, 14177ba1a0bfSKris Buschelman MATOP_PTAP_SYMBOLIC_MPIAIJ=101, 141887d4246cSBarry Smith MATOP_PTAP_NUMERIC_MPIAIJ=102, 1419f5edf698SHong Zhang MATOP_SET_VALUES_ROW = 105, 1420d5c63f83SSatish Balay MATOP_GET_ROW_UTRIANGULAR=108, 14212bebee5dSHong Zhang MATOP_RESTORE_ROW_UTRIANGULAR=109, 142269db28dcSHong Zhang MATOP_MATSOLVE=110, 1423829201f2SHong Zhang MATOP_GET_REDUNDANTMATRIX=111, 1424829201f2SHong Zhang MATOP_MATGETSEQNONZEROSTRUCTURE=115 1425fae171e0SBarry Smith } MatOperation; 1426be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatHasOperation(Mat,MatOperation,PetscTruth*); 1427be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellSetOperation(Mat,MatOperation,void(*)(void)); 1428be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellGetOperation(Mat,MatOperation,void(**)(void)); 1429be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellSetContext(Mat,void*); 1430112a2221SBarry Smith 143190ace30eSBarry Smith /* 143290ace30eSBarry Smith Codes for matrices stored on disk. By default they are 143390ace30eSBarry Smith stored in a universal format. By changing the format with 1434fb9695e5SSatish Balay PetscViewerSetFormat(viewer,PETSC_VIEWER_BINARY_NATIVE); the matrices will 143590ace30eSBarry Smith be stored in a way natural for the matrix, for example dense matrices 143690ace30eSBarry Smith would be stored as dense. Matrices stored this way may only be 143790ace30eSBarry Smith read into matrices of the same time. 143890ace30eSBarry Smith */ 143990ace30eSBarry Smith #define MATRIX_BINARY_FORMAT_DENSE -1 144090ace30eSBarry Smith 1441be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetHashTableFactor(Mat,PetscReal); 1442be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatISGetLocalMat(Mat,Mat*); 14431f4e1ec7SBarry Smith 1444d9274352SBarry Smith /*S 1445d9274352SBarry Smith MatNullSpace - Object that removes a null space from a vector, i.e. 1446d9274352SBarry Smith orthogonalizes the vector to a subsapce 1447d9274352SBarry Smith 1448f7a9e4ceSBarry Smith Level: advanced 1449d9274352SBarry Smith 1450d9274352SBarry Smith Concepts: matrix; linear operator, null space 1451d9274352SBarry Smith 14526e1639daSBarry Smith Users manual sections: 14536e1639daSBarry Smith . sec_singular 14546e1639daSBarry Smith 1455d9274352SBarry Smith .seealso: MatNullSpaceCreate() 1456d9274352SBarry Smith S*/ 145774637425SBarry Smith typedef struct _p_MatNullSpace* MatNullSpace; 1458d9274352SBarry Smith 1459be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceCreate(MPI_Comm,PetscTruth,PetscInt,const Vec[],MatNullSpace*); 1460281d1b2eSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(Vec,void*),void*); 1461be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceDestroy(MatNullSpace); 1462be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceRemove(MatNullSpace,Vec,Vec*); 1463be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceAttach(Mat,MatNullSpace); 1464be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceTest(MatNullSpace,Mat); 146574637425SBarry Smith 1466be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatReorderingSeqSBAIJ(Mat,IS); 1467be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPISBAIJSetHashTableFactor(Mat,PetscReal); 1468be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetColumnIndices(Mat,PetscInt *); 146946129b97SKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJInvertBlockDiagonal(Mat); 14703f1d51d7SBarry Smith 1471be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMAIJ(Mat,PetscInt,Mat*); 1472be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMAIJRedimension(Mat,PetscInt,Mat*); 1473be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMAIJGetAIJ(Mat,Mat*); 1474c4f061fbSSatish Balay 1475be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatComputeExplicitOperator(Mat,Mat*); 1476b0a32e0cSBarry Smith 1477be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalScaleLocal(Mat,Vec); 147804f1ad80SBarry Smith 1479fef1beadSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*); 14803ec795f1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetBase(Mat,Vec,Vec); 14813ec795f1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*); 1482e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*)); 1483e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec)); 1484e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDAddNullSpace(Mat,MatNullSpace); 1485e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt); 1486e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDResetHHistory(Mat); 1487e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctionError(Mat,PetscReal); 1488e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetPeriod(Mat,PetscInt); 1489e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDGetH(Mat,PetscScalar *); 1490e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFromOptions(Mat); 1491e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*); 1492e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*); 1493e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDDSSetUmin(Mat,PetscReal); 1494e884886fSBarry Smith 1495e884886fSBarry Smith 1496e884886fSBarry Smith typedef struct _p_MatMFFD* MatMFFD; 1497e884886fSBarry Smith 1498e884886fSBarry Smith /*E 1499e884886fSBarry Smith MatMFFDType - algorithm used to compute the h used in computing matrix-vector products via differencing of the function 1500e884886fSBarry Smith 1501e884886fSBarry Smith Level: beginner 1502e884886fSBarry Smith 1503e884886fSBarry Smith .seealso: MatMFFDSetType(), MatMFFDRegister() 1504e884886fSBarry Smith E*/ 1505e884886fSBarry Smith #define MatMFFDType const char* 1506e884886fSBarry Smith #define MATMFFD_DS "ds" 1507e884886fSBarry Smith #define MATMFFD_WP "wp" 1508e884886fSBarry Smith 1509e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetType(Mat,MatMFFDType); 1510e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatMFFD)); 1511e884886fSBarry Smith 1512e884886fSBarry Smith /*MC 1513e884886fSBarry Smith MatMFFDRegisterDynamic - Adds a method to the MatMFFD registry. 1514e884886fSBarry Smith 1515e884886fSBarry Smith Synopsis: 1516e884886fSBarry Smith PetscErrorCode MatMFFDRegisterDynamic(char *name_solver,char *path,char *name_create,PetscErrorCode (*routine_create)(MatMFFD)) 1517e884886fSBarry Smith 1518e884886fSBarry Smith Not Collective 1519e884886fSBarry Smith 1520e884886fSBarry Smith Input Parameters: 1521e884886fSBarry Smith + name_solver - name of a new user-defined compute-h module 1522e884886fSBarry Smith . path - path (either absolute or relative) the library containing this solver 1523e884886fSBarry Smith . name_create - name of routine to create method context 1524e884886fSBarry Smith - routine_create - routine to create method context 1525e884886fSBarry Smith 1526e884886fSBarry Smith Level: developer 1527e884886fSBarry Smith 1528e884886fSBarry Smith Notes: 1529e884886fSBarry Smith MatMFFDRegisterDynamic() may be called multiple times to add several user-defined solvers. 1530e884886fSBarry Smith 1531e884886fSBarry Smith If dynamic libraries are used, then the fourth input argument (routine_create) 1532e884886fSBarry Smith is ignored. 1533e884886fSBarry Smith 1534e884886fSBarry Smith Sample usage: 1535e884886fSBarry Smith .vb 1536e884886fSBarry Smith MatMFFDRegisterDynamic("my_h",/home/username/my_lib/lib/libO/solaris/mylib.a, 1537e884886fSBarry Smith "MyHCreate",MyHCreate); 1538e884886fSBarry Smith .ve 1539e884886fSBarry Smith 1540e884886fSBarry Smith Then, your solver can be chosen with the procedural interface via 1541e884886fSBarry Smith $ MatMFFDSetType(mfctx,"my_h") 1542e884886fSBarry Smith or at runtime via the option 1543e884886fSBarry Smith $ -snes_mf_type my_h 1544e884886fSBarry Smith 1545e884886fSBarry Smith .keywords: MatMFFD, register 1546e884886fSBarry Smith 1547e884886fSBarry Smith .seealso: MatMFFDRegisterAll(), MatMFFDRegisterDestroy() 1548e884886fSBarry Smith M*/ 1549e884886fSBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 1550e884886fSBarry Smith #define MatMFFDRegisterDynamic(a,b,c,d) MatMFFDRegister(a,b,c,0) 1551e884886fSBarry Smith #else 1552e884886fSBarry Smith #define MatMFFDRegisterDynamic(a,b,c,d) MatMFFDRegister(a,b,c,d) 1553e884886fSBarry Smith #endif 1554e884886fSBarry Smith 1555e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDRegisterAll(const char[]); 1556e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDRegisterDestroy(void); 1557e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDDefaultSetUmin(Mat,PetscReal); 1558e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDWPSetComputeNormU(Mat,PetscTruth); 1559e884886fSBarry Smith 1560e884886fSBarry Smith 1561be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *); 1562be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *); 15637dbadf16SMatthew Knepley 1564e9fa29b7SSatish Balay PETSC_EXTERN_CXX_END 15652eac72dbSBarry Smith #endif 1566