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*/ 29*e5a9bf91SBarry 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" 40273d9f13SBarry Smith #define MATSEQBDIAG "seqbdiag" 41273d9f13SBarry Smith #define MATMPIBDIAG "mpibdiag" 42209238afSKris Buschelman #define MATBDIAG "bdiag" 43209238afSKris Buschelman #define MATSEQDENSE "seqdense" 44273d9f13SBarry Smith #define MATMPIDENSE "mpidense" 45209238afSKris Buschelman #define MATDENSE "dense" 46273d9f13SBarry Smith #define MATSEQBAIJ "seqbaij" 47273d9f13SBarry Smith #define MATMPIBAIJ "mpibaij" 48209238afSKris Buschelman #define MATBAIJ "baij" 49273d9f13SBarry Smith #define MATMPIADJ "mpiadj" 50273d9f13SBarry Smith #define MATSEQSBAIJ "seqsbaij" 51273d9f13SBarry Smith #define MATMPISBAIJ "mpisbaij" 52209238afSKris Buschelman #define MATSBAIJ "sbaij" 53cebc7f6cSBarry Smith #define MATDAAD "daad" 54cebc7f6cSBarry Smith #define MATMFFD "mffd" 55c8a8475eSBarry Smith #define MATNORMAL "normal" 56ab92ecdeSBarry Smith #define MATLRC "lrc" 57b3a1e11cSKris Buschelman #define MATSEQAIJSPOOLES "seqaijspooles" 58d10c748bSKris Buschelman #define MATMPIAIJSPOOLES "mpiaijspooles" 599abb65ffSKris Buschelman #define MATSEQSBAIJSPOOLES "seqsbaijspooles" 6022191285SKris Buschelman #define MATMPISBAIJSPOOLES "mpisbaijspooles" 61bb4d4166SHong Zhang #define MATAIJSPOOLES "aijspooles" 62bb4d4166SHong Zhang #define MATSBAIJSPOOLES "sbaijspooles" 6314b396bbSKris Buschelman #define MATSUPERLU "superlu" 641d96aa28SKris Buschelman #define MATSUPERLU_DIST "superlu_dist" 651677d0b8SKris Buschelman #define MATUMFPACK "umfpack" 66e8aa55a4SKris Buschelman #define MATESSL "essl" 674eb8e494SKris Buschelman #define MATLUSOL "lusol" 68397b6df1SKris Buschelman #define MATAIJMUMPS "aijmumps" 69397b6df1SKris Buschelman #define MATSBAIJMUMPS "sbaijmumps" 708da957c5SKris Buschelman #define MATDSCPACK "dscpack" 717065e2aaSKris Buschelman #define MATMATLAB "matlab" 724b8d99adSRichard Tran Mills #define MATSEQCSRPERM "seqcsrperm" 7318def556SRichard Tran Mills #define MATMPICSRPERM "mpicsrperm" 74814655b8SBarry Smith #define MATCSRPERM "csrperm" 7581824310SBarry Smith #define MATSEQCRL "seqcrl" 76c02b24c5SBarry Smith #define MATMPICRL "mpicrl" 77c02b24c5SBarry Smith #define MATCRL "crl" 78c0aa2d19SHong Zhang #define MATPLAPACK "plapack" 792a6744ebSBarry Smith #define MATSCATTER "scatter" 80421e10b8SBarry Smith #define MATBLOCKMAT "blockmat" 81793850ffSBarry Smith #define MATCOMPOSITE "composite" 825a7f1df3SHong Zhang #define MATSEQFFTW "seqfftw" 83d91e6319SBarry Smith 84c06d978dSMatthew Knepley /* Logging support */ 85552e946dSBarry Smith #define MAT_FILE_COOKIE 1211216 /* used to indicate matrices in binary files */ 86be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_COOKIE; 87be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_FDCOLORING_COOKIE; 88be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_PARTITIONING_COOKIE; 89be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_NULLSPACE_COOKIE; 90c06d978dSMatthew Knepley 91ceb03754SKris Buschelman /*E 92ceb03754SKris Buschelman MatReuse - Indicates if matrices obtained from a previous call to MatGetSubMatrices() 93ceb03754SKris Buschelman or MatGetSubMatrix() are to be reused to store the new matrix values. 94ceb03754SKris Buschelman 95ceb03754SKris Buschelman Level: beginner 96ceb03754SKris Buschelman 97ceb03754SKris Buschelman Any additions/changes here MUST also be made in include/finclude/petscmat.h 98ceb03754SKris Buschelman 990c451bc4SBarry Smith .seealso: MatGetSubMatrices(), MatGetSubMatrix(), MatDestroyMatrices(), MatConvert() 100ceb03754SKris Buschelman E*/ 101ceb03754SKris Buschelman typedef enum {MAT_INITIAL_MATRIX,MAT_REUSE_MATRIX} MatReuse; 102ceb03754SKris Buschelman 103be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInitializePackage(char *); 104c06d978dSMatthew Knepley 105f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreate(MPI_Comm,Mat*); 106a69a960dSMatthew Knepley PetscPolymorphicFunction(MatCreate,(MPI_Comm comm),(comm,&A),Mat,A) 107a69a960dSMatthew Knepley PetscPolymorphicFunction(MatCreate,(),(PETSC_COMM_WORLD,&A),Mat,A) 108f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetSizes(Mat,PetscInt,PetscInt,PetscInt,PetscInt); 109f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetType(Mat,MatType); 110be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetFromOptions(Mat); 111be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUpPreallocation(Mat); 112be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRegisterAll(const char[]); 113be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat)); 11430de9b25SBarry Smith 115f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetOptionsPrefix(Mat,const char[]); 116f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAppendOptionsPrefix(Mat,const char[]); 117f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOptionsPrefix(Mat,const char*[]); 118f69a0ea3SMatthew Knepley 11930de9b25SBarry Smith /*MC 12030de9b25SBarry Smith MatRegisterDynamic - Adds a new matrix type 12130de9b25SBarry Smith 12230de9b25SBarry Smith Synopsis: 123c1ac3661SBarry Smith PetscErrorCode MatRegisterDynamic(char *name,char *path,char *name_create,PetscErrorCode (*routine_create)(Mat)) 12430de9b25SBarry Smith 12530de9b25SBarry Smith Not Collective 12630de9b25SBarry Smith 12730de9b25SBarry Smith Input Parameters: 12830de9b25SBarry Smith + name - name of a new user-defined matrix type 12930de9b25SBarry Smith . path - path (either absolute or relative) the library containing this solver 13030de9b25SBarry Smith . name_create - name of routine to create method context 13130de9b25SBarry Smith - routine_create - routine to create method context 13230de9b25SBarry Smith 13330de9b25SBarry Smith Notes: 13430de9b25SBarry Smith MatRegisterDynamic() may be called multiple times to add several user-defined solvers. 13530de9b25SBarry Smith 13630de9b25SBarry Smith If dynamic libraries are used, then the fourth input argument (routine_create) 13730de9b25SBarry Smith is ignored. 13830de9b25SBarry Smith 13930de9b25SBarry Smith Sample usage: 14030de9b25SBarry Smith .vb 14130de9b25SBarry Smith MatRegisterDynamic("my_mat",/home/username/my_lib/lib/libO/solaris/mylib.a, 14230de9b25SBarry Smith "MyMatCreate",MyMatCreate); 14330de9b25SBarry Smith .ve 14430de9b25SBarry Smith 14530de9b25SBarry Smith Then, your solver can be chosen with the procedural interface via 14630de9b25SBarry Smith $ MatSetType(Mat,"my_mat") 14730de9b25SBarry Smith or at runtime via the option 14830de9b25SBarry Smith $ -mat_type my_mat 14930de9b25SBarry Smith 15030de9b25SBarry Smith Level: advanced 15130de9b25SBarry Smith 152ab901514SBarry Smith Notes: ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values. 15330de9b25SBarry Smith If your function is not being put into a shared library then use VecRegister() instead 15430de9b25SBarry Smith 15530de9b25SBarry Smith .keywords: Mat, register 15630de9b25SBarry Smith 15730de9b25SBarry Smith .seealso: MatRegisterAll(), MatRegisterDestroy() 15830de9b25SBarry Smith 15930de9b25SBarry Smith M*/ 160273d9f13SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 161273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,0) 162273d9f13SBarry Smith #else 163273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,d) 16430de9b25SBarry Smith #endif 16530de9b25SBarry Smith 166273d9f13SBarry Smith extern PetscTruth MatRegisterAllCalled; 167b0a32e0cSBarry Smith extern PetscFList MatList; 16828988994SBarry Smith 169be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqDense(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*); 170be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIDense(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*); 171be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 172ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,nz,nnz,&A),Mat,A) 173ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,&A),Mat,A) 174ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,0,nnz,&A),Mat,A) 175ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,&A),Mat,A) 176ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,A)) 177ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,m,n,0,nnz,A)) 178ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,A)) 179be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 180ba966639SSatish 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) 181ba966639SSatish 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) 182ba966639SSatish 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) 183ba966639SSatish 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) 184ba966639SSatish 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)) 185ba966639SSatish 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)) 186ba966639SSatish 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)) 187ba966639SSatish 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) 188ba966639SSatish 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) 189ba966639SSatish 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) 190ba966639SSatish 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) 191ba966639SSatish 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)) 192ba966639SSatish 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)) 193ba966639SSatish 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)) 19463ab93bfSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *); 1958d7a6e47SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJWithSplitArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],PetscInt[],PetscInt[],PetscScalar[],Mat*); 1968d7a6e47SBarry Smith 197be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIRowbs(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 198be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqBDiag(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscScalar*[],Mat*); 199be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIBDiag(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscScalar*[],Mat*); 200be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 201ba966639SSatish 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) 202ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A) 203ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A) 204ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A) 205ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A)) 206ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A)) 207ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A)) 208be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 209ba966639SSatish 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) 210ba966639SSatish 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) 211ba966639SSatish 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) 212ba966639SSatish 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) 213ba966639SSatish 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)) 214ba966639SSatish 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)) 215ba966639SSatish 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)) 216ba966639SSatish 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) 217ba966639SSatish 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) 218ba966639SSatish 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) 219ba966639SSatish 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) 220ba966639SSatish 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)) 221ba966639SSatish 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)) 222ba966639SSatish 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)) 223be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAdj(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscInt[],Mat*); 224be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 225ba966639SSatish 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) 226ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A) 227ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A) 228ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A) 229ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A)) 230ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A)) 231ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A)) 232be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPISBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 233ba966639SSatish 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) 234ba966639SSatish 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) 235ba966639SSatish 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) 236ba966639SSatish 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) 237ba966639SSatish 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)) 238ba966639SSatish 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)) 239ba966639SSatish 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)) 240ba966639SSatish 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) 241ba966639SSatish 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) 242ba966639SSatish 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) 243ba966639SSatish 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) 244ba966639SSatish 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)) 245ba966639SSatish 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)) 246ba966639SSatish 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)) 247be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateShell(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,void *,Mat*); 248ba966639SSatish 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) 249ba966639SSatish Balay PetscPolymorphicFunction(MatCreateShell,(PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx),(PETSC_COMM_WORLD,m,n,M,N,ctx,&A),Mat,A) 250be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateAdic(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,void (*)(void),Mat*); 251be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateNormal(Mat,Mat*); 252ba966639SSatish Balay PetscPolymorphicFunction(MatCreateNormal,(Mat mat),(mat,&A),Mat,A) 253ab92ecdeSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateLRC(Mat,Mat,Mat,Mat*); 254284134d9SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateIS(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,ISLocalToGlobalMapping,Mat*); 2551472f72bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 2561472f72bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPICRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 2572a6744ebSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateScatter(MPI_Comm,VecScatter,Mat*); 2582a6744ebSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScatterSetVecScatter(Mat,VecScatter); 2592a6744ebSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScatterGetVecScatter(Mat,VecScatter*); 2608cdf2d9bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateBlockMat(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,Mat*); 261793850ffSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompositeAddMat(Mat,Mat); 262793850ffSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateComposite(MPI_Comm,PetscInt,const Mat*,Mat*); 2635a7f1df3SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqFFTW(MPI_Comm,PetscInt,const PetscInt[],Mat*); 2641472f72bSBarry Smith 265f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUp(Mat); 266be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroy(Mat); 26721c89e3eSBarry Smith 2680c451bc4SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConjugate(Mat); 26999cafbc1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRealPart(Mat); 27099cafbc1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatImaginaryPart(Mat); 27199cafbc1SBarry Smith 2728ed539a5SBarry Smith /* ------------------------------------------------------------*/ 273be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 274be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 27587d4246cSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesRow(Mat,PetscInt,const PetscScalar[]); 276f259bd47SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesRowLocal(Mat,PetscInt,const PetscScalar[]); 27784cb2905SBarry Smith 2782ef4de8bSBarry Smith /*S 2792ef4de8bSBarry Smith MatStencil - Data structure (C struct) for storing information about a single row or 2802ef4de8bSBarry Smith column of a matrix as index on an associated grid. 2812ef4de8bSBarry Smith 2822ef4de8bSBarry Smith Level: beginner 2832ef4de8bSBarry Smith 2842ef4de8bSBarry Smith Concepts: matrix; linear operator 2852ef4de8bSBarry Smith 286d7bd93baSBarry Smith .seealso: MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockStencil() 2872ef4de8bSBarry Smith S*/ 288435da068SBarry Smith typedef struct { 289c1ac3661SBarry Smith PetscInt k,j,i,c; 290435da068SBarry Smith } MatStencil; 2912ef4de8bSBarry Smith 292be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode); 293be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode); 294be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt); 295435da068SBarry Smith 296be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetColoring(Mat,ISColoring); 297be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesAdic(Mat,void*); 298be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesAdifor(Mat,PetscInt,void*); 2993a7fca6bSBarry Smith 300d91e6319SBarry Smith /*E 301d91e6319SBarry Smith MatAssemblyType - Indicates if the matrix is now to be used, or if you plan 302d91e6319SBarry Smith to continue to add values to it 303d91e6319SBarry Smith 304d91e6319SBarry Smith Level: beginner 305d91e6319SBarry Smith 306d91e6319SBarry Smith .seealso: MatAssemblyBegin(), MatAssemblyEnd() 307d91e6319SBarry Smith E*/ 3086d4a8577SBarry Smith typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType; 309be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssemblyBegin(Mat,MatAssemblyType); 310be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssemblyEnd(Mat,MatAssemblyType); 311be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssembled(Mat,PetscTruth*); 3124f9c727eSBarry Smith 313be1d678aSKris Buschelman extern PETSCMAT_DLLEXPORT PetscInt MatSetValue_Row; 314be1d678aSKris Buschelman extern PETSCMAT_DLLEXPORT PetscInt MatSetValue_Column; 315be1d678aSKris Buschelman extern PETSCMAT_DLLEXPORT PetscScalar MatSetValue_Value; 3161d73ed98SBarry Smith 31730de9b25SBarry Smith /*MC 31830de9b25SBarry Smith MatSetValue - Set a single entry into a matrix. 31930de9b25SBarry Smith 32030de9b25SBarry Smith Synopsis: 321c1ac3661SBarry Smith PetscErrorCode MatSetValue(Mat m,PetscInt row,PetscInt col,PetscScalar value,InsertMode mode); 32230de9b25SBarry Smith 32330de9b25SBarry Smith Not collective 32430de9b25SBarry Smith 32530de9b25SBarry Smith Input Parameters: 32630de9b25SBarry Smith + m - the matrix 32730de9b25SBarry Smith . row - the row location of the entry 32830de9b25SBarry Smith . col - the column location of the entry 32930de9b25SBarry Smith . value - the value to insert 33030de9b25SBarry Smith - mode - either INSERT_VALUES or ADD_VALUES 33130de9b25SBarry Smith 33230de9b25SBarry Smith Notes: 33330de9b25SBarry Smith For efficiency one should use MatSetValues() and set several or many 33430de9b25SBarry Smith values simultaneously if possible. 33530de9b25SBarry Smith 33630de9b25SBarry Smith Level: beginner 33730de9b25SBarry Smith 33830de9b25SBarry Smith .seealso: MatSetValues(), MatSetValueLocal() 33930de9b25SBarry Smith M*/ 340b951964fSBarry Smith #define MatSetValue(v,i,j,va,mode) \ 341f1144a30SSatish Balay ((MatSetValue_Row = i,MatSetValue_Column = j,MatSetValue_Value = va,0) || \ 3421d73ed98SBarry Smith MatSetValues(v,1,&MatSetValue_Row,1,&MatSetValue_Column,&MatSetValue_Value,mode)) 34330de9b25SBarry Smith 344ea06a074SBarry Smith #define MatGetValue(v,i,j,va) \ 345f1144a30SSatish Balay ((MatSetValue_Row = i,MatSetValue_Column = j,0) || \ 3461d73ed98SBarry Smith MatGetValues(v,1,&MatSetValue_Row,1,&MatSetValue_Column,&va)) 34730de9b25SBarry Smith 348d91e6319SBarry Smith #define MatSetValueLocal(v,i,j,va,mode) \ 349f1144a30SSatish Balay ((MatSetValue_Row = i,MatSetValue_Column = j,MatSetValue_Value = va,0) || \ 3501d73ed98SBarry Smith MatSetValuesLocal(v,1,&MatSetValue_Row,1,&MatSetValue_Column,&MatSetValue_Value,mode)) 35130de9b25SBarry Smith 352d91e6319SBarry Smith /*E 353d91e6319SBarry Smith MatOption - Options that may be set for a matrix and its behavior or storage 354d91e6319SBarry Smith 355d91e6319SBarry Smith Level: beginner 356d91e6319SBarry Smith 3570a835dfdSSatish Balay Any additions/changes here MUST also be made in include/finclude/petscmat.h 358d91e6319SBarry Smith 359d91e6319SBarry Smith .seealso: MatSetOption() 360d91e6319SBarry Smith E*/ 361290bbb0aSBarry Smith typedef enum {MAT_ROW_ORIENTED,MAT_COLUMN_ORIENTED,MAT_ROWS_SORTED, 362290bbb0aSBarry Smith MAT_COLUMNS_SORTED,MAT_NO_NEW_NONZERO_LOCATIONS, 363290bbb0aSBarry Smith MAT_YES_NEW_NONZERO_LOCATIONS,MAT_SYMMETRIC, 364290bbb0aSBarry Smith MAT_STRUCTURALLY_SYMMETRIC,MAT_NO_NEW_DIAGONALS, 365290bbb0aSBarry Smith MAT_YES_NEW_DIAGONALS,MAT_INODE_LIMIT_1,MAT_INODE_LIMIT_2, 366290bbb0aSBarry Smith MAT_INODE_LIMIT_3,MAT_INODE_LIMIT_4,MAT_INODE_LIMIT_5, 367290bbb0aSBarry Smith MAT_IGNORE_OFF_PROC_ENTRIES,MAT_ROWS_UNSORTED, 368290bbb0aSBarry Smith MAT_COLUMNS_UNSORTED,MAT_NEW_NONZERO_LOCATION_ERR, 369290bbb0aSBarry Smith MAT_NEW_NONZERO_ALLOCATION_ERR,MAT_USE_HASH_TABLE, 370290bbb0aSBarry Smith MAT_KEEP_ZEROED_ROWS,MAT_IGNORE_ZERO_ENTRIES,MAT_USE_INODES, 371290bbb0aSBarry Smith MAT_DO_NOT_USE_INODES,MAT_NOT_SYMMETRIC,MAT_HERMITIAN, 372290bbb0aSBarry Smith MAT_NOT_STRUCTURALLY_SYMMETRIC,MAT_NOT_HERMITIAN, 373290bbb0aSBarry Smith MAT_SYMMETRY_ETERNAL,MAT_NOT_SYMMETRY_ETERNAL, 374290bbb0aSBarry Smith MAT_USE_COMPRESSEDROW,MAT_DO_NOT_USE_COMPRESSEDROW, 375290bbb0aSBarry Smith MAT_IGNORE_LOWER_TRIANGULAR,MAT_ERROR_LOWER_TRIANGULAR,MAT_GETROW_UPPERTRIANGULAR} MatOption; 376290bbb0aSBarry Smith extern const char *MatOptions[]; 377be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetOption(Mat,MatOption); 378be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetType(Mat,MatType*); 379ba966639SSatish Balay PetscPolymorphicFunction(MatGetType,(Mat mat),(mat,&t),MatType,t) 38084cb2905SBarry Smith 381be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]); 382be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 383be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 384f5edf698SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowUpperTriangular(Mat); 385f5edf698SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRowUpperTriangular(Mat); 386be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 387be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 388be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumnVector(Mat,Vec,PetscInt); 389be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetArray(Mat,PetscScalar *[]); 390ba966639SSatish Balay PetscPolymorphicFunction(MatGetArray,(Mat mat),(mat,&a),PetscScalar*,a) 391be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreArray(Mat,PetscScalar *[]); 392be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBlockSize(Mat,PetscInt *); 393ba966639SSatish Balay PetscPolymorphicFunction(MatGetBlockSize,(Mat mat),(mat,&a),PetscInt,a) 394be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetBlockSize(Mat,PetscInt); 3957b80b807SBarry Smith 396be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMult(Mat,Vec,Vec); 397be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultAdd(Mat,Vec,Vec,Vec); 398ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultAdd,(Mat A,Vec x,Vec y),(A,x,y,y)) 399be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTranspose(Mat,Vec,Vec); 400be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsTranspose(Mat,Mat,PetscReal,PetscTruth*); 401ba966639SSatish Balay PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B,PetscReal tol),(A,B,tol,&t),PetscTruth,t) 402e37f64cdSBarry Smith PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B),(A,B,0,&t),PetscTruth,t) 403be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeAdd(Mat,Vec,Vec,Vec); 404ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultTransposeAdd,(Mat A,Vec x,Vec y),(A,x,y,y)) 405be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultConstrained(Mat,Vec,Vec); 406be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeConstrained(Mat,Vec,Vec); 4072bebee5dSHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatSolve(Mat,Mat,Mat); 408f5edf698SHong Zhang 409d91e6319SBarry Smith /*E 410d91e6319SBarry Smith MatDuplicateOption - Indicates if a duplicated sparse matrix should have 411d91e6319SBarry Smith its numerical values copied over or just its nonzero structure. 412d91e6319SBarry Smith 413d91e6319SBarry Smith Level: beginner 414d91e6319SBarry Smith 415d91e6319SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 416d91e6319SBarry Smith 417d91e6319SBarry Smith .seealso: MatDuplicate() 418d91e6319SBarry Smith E*/ 4192e8a6d31SBarry Smith typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES} MatDuplicateOption; 4202e8a6d31SBarry Smith 421f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert(Mat,MatType,MatReuse,Mat*); 422ba966639SSatish Balay PetscPolymorphicFunction(MatConvert,(Mat A,MatType t),(A,t,MAT_INITIAL_MATRIX,&a),Mat,a) 423be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDuplicate(Mat,MatDuplicateOption,Mat*); 424ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A,MatDuplicateOption o),(A,o,&a),Mat,a) 425ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A),(A,MAT_COPY_VALUES,&a),Mat,a) 42694a9d846SBarry Smith 427d91e6319SBarry Smith /*E 428d91e6319SBarry Smith MatStructure - Indicates if the matrix has the same nonzero structure 429d91e6319SBarry Smith 430d91e6319SBarry Smith Level: beginner 431d91e6319SBarry Smith 432d91e6319SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 433d91e6319SBarry Smith 43494b7f48cSBarry Smith .seealso: MatCopy(), KSPSetOperators(), PCSetOperators() 435d91e6319SBarry Smith E*/ 436c537a176SHong Zhang typedef enum {SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER,SUBSET_NONZERO_PATTERN} MatStructure; 437cb5b572fSBarry Smith 438be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCopy(Mat,Mat,MatStructure); 439be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatView(Mat,PetscViewer); 440be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsSymmetric(Mat,PetscReal,PetscTruth*); 441ba966639SSatish Balay PetscPolymorphicFunction(MatIsSymmetric,(Mat A,PetscReal tol),(A,tol,&t),PetscTruth,t) 442e37f64cdSBarry Smith PetscPolymorphicFunction(MatIsSymmetric,(Mat A),(A,0,&t),PetscTruth,t) 443be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsStructurallySymmetric(Mat,PetscTruth*); 444ba966639SSatish Balay PetscPolymorphicFunction(MatIsStructurallySymmetric,(Mat A),(A,&t),PetscTruth,t) 445be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitian(Mat,PetscTruth*); 446ba966639SSatish Balay PetscPolymorphicFunction(MatIsHermitian,(Mat A),(A,&t),PetscTruth,t) 447be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsSymmetricKnown(Mat,PetscTruth*,PetscTruth*); 448be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitianKnown(Mat,PetscTruth*,PetscTruth*); 449f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLoad(PetscViewer,MatType,Mat*); 450ba966639SSatish Balay PetscPolymorphicFunction(MatLoad,(PetscViewer v,MatType t),(v,t,&a),Mat,a) 4517b80b807SBarry Smith 452be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowIJ(Mat,PetscInt,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *); 453be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRowIJ(Mat,PetscInt,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *); 454be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumnIJ(Mat,PetscInt,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *); 455be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreColumnIJ(Mat,PetscInt,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *); 456d4fbbf0eSBarry Smith 457d91e6319SBarry Smith /*S 458d91e6319SBarry Smith MatInfo - Context of matrix information, used with MatGetInfo() 459d91e6319SBarry Smith 460d91e6319SBarry Smith In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE 461d91e6319SBarry Smith 462d91e6319SBarry Smith Level: intermediate 463d91e6319SBarry Smith 464d91e6319SBarry Smith Concepts: matrix^nonzero information 465d91e6319SBarry Smith 466d9274352SBarry Smith .seealso: MatGetInfo(), MatInfoType 467d91e6319SBarry Smith S*/ 4684e220ebcSLois Curfman McInnes typedef struct { 469b0a32e0cSBarry Smith PetscLogDouble rows_global,columns_global; /* number of global rows and columns */ 470b0a32e0cSBarry Smith PetscLogDouble rows_local,columns_local; /* number of local rows and columns */ 471b0a32e0cSBarry Smith PetscLogDouble block_size; /* block size */ 472b0a32e0cSBarry Smith PetscLogDouble nz_allocated,nz_used,nz_unneeded; /* number of nonzeros */ 473b0a32e0cSBarry Smith PetscLogDouble memory; /* memory allocated */ 474b0a32e0cSBarry Smith PetscLogDouble assemblies; /* number of matrix assemblies called */ 475b0a32e0cSBarry Smith PetscLogDouble mallocs; /* number of mallocs during MatSetValues() */ 476b0a32e0cSBarry Smith PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */ 477b0a32e0cSBarry Smith PetscLogDouble factor_mallocs; /* number of mallocs during factorization */ 4784e220ebcSLois Curfman McInnes } MatInfo; 4794e220ebcSLois Curfman McInnes 480d9274352SBarry Smith /*E 481d9274352SBarry Smith MatInfoType - Indicates if you want information about the local part of the matrix, 482d9274352SBarry Smith the entire parallel matrix or the maximum over all the local parts. 483d9274352SBarry Smith 484d9274352SBarry Smith Level: beginner 485d9274352SBarry Smith 486d9274352SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 487d9274352SBarry Smith 488d9274352SBarry Smith .seealso: MatGetInfo(), MatInfo 489d9274352SBarry Smith E*/ 4907b80b807SBarry Smith typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType; 491be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetInfo(Mat,MatInfoType,MatInfo*); 492be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatValid(Mat,PetscTruth*); 493be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetDiagonal(Mat,Vec); 494be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMax(Mat,Vec); 495be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatTranspose(Mat,Mat*); 496ba966639SSatish Balay PetscPolymorphicFunction(MatTranspose,(Mat A),(A,&t),Mat,t) 497be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPermute(Mat,IS,IS,Mat *); 498ba966639SSatish Balay PetscPolymorphicFunction(MatPermute,(Mat A,IS is1,IS is2),(A,is1,is2,&t),Mat,t) 499be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPermuteSparsify(Mat,PetscInt,PetscReal,PetscReal,IS,IS,Mat *); 500be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalScale(Mat,Vec,Vec); 501be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalSet(Mat,Vec,InsertMode); 502be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatEqual(Mat,Mat,PetscTruth*); 503ba966639SSatish Balay PetscPolymorphicFunction(MatEqual,(Mat A,Mat B),(A,B,&t),PetscTruth,t) 504be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultEqual(Mat,Mat,PetscInt,PetscTruth*); 505be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultAddEqual(Mat,Mat,PetscInt,PetscTruth*); 506be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeEqual(Mat,Mat,PetscInt,PetscTruth*); 507be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscTruth*); 5087b80b807SBarry Smith 509be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNorm(Mat,NormType,PetscReal *); 510ba966639SSatish Balay PetscPolymorphicFunction(MatNorm,(Mat A,NormType t),(A,t,&n),PetscReal,n) 511be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroEntries(Mat); 512f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar); 513f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsIS(Mat,IS,PetscScalar); 514f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroColumns(Mat,PetscInt,const PetscInt [],const PetscScalar*); 515f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroColumnsIS(Mat,IS,const PetscScalar*); 5167b80b807SBarry Smith 517be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatUseScaledForm(Mat,PetscTruth); 518be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScaleSystem(Mat,Vec,Vec); 519be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatUnScaleSystem(Mat,Vec,Vec); 5205ef9f2a5SBarry Smith 521be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSize(Mat,PetscInt*,PetscInt*); 522be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalSize(Mat,PetscInt*,PetscInt*); 523be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRange(Mat,PetscInt*,PetscInt*); 5248ee2e534SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRanges(Mat,const PetscInt**); 5257b80b807SBarry Smith 526be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]); 527be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroyMatrices(PetscInt,Mat *[]); 528be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrix(Mat,IS,IS,PetscInt,MatReuse,Mat *); 529abf0e977SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrixRaw(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscInt,MatReuse,Mat *); 530be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge(MPI_Comm,Mat,PetscInt,MatReuse,Mat*); 531be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPI(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*); 532be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPISymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*); 533be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPINumeric(Mat,Mat); 534be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroy_MPIAIJ_SeqsToMPI(Mat); 535be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalMat(Mat,MatReuse,Mat*); 536be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*); 537be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,PetscInt*,Mat*); 538be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBrowsOfAoCols(Mat,Mat,MatReuse,PetscInt**,PetscScalar**,Mat*); 53943eb5e2fSMatthew Knepley #if defined (PETSC_USE_CTABLE) 540cd116e26SSatish Balay #include "petscctable.h" 54143eb5e2fSMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetCommunicationStructs(Mat, Vec *, PetscTable *, VecScatter *); 54243eb5e2fSMatthew Knepley #else 54343eb5e2fSMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetCommunicationStructs(Mat, Vec *, PetscInt *[], VecScatter *); 54443eb5e2fSMatthew Knepley #endif 5458efafbd8SBarry Smith 546be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt); 5477b80b807SBarry Smith 548be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*); 549be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultSymbolic(Mat,Mat,PetscReal,Mat*); 550be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultNumeric(Mat,Mat,Mat); 55122440eb1SKris Buschelman 552be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*); 553be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAPSymbolic(Mat,Mat,PetscReal,Mat*); 554be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAPNumeric(Mat,Mat,Mat); 55522440eb1SKris Buschelman 556be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTranspose(Mat,Mat,MatReuse,PetscReal,Mat*); 557be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTransposeSymbolic(Mat,Mat,PetscReal,Mat*); 558be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTransposeNumeric(Mat,Mat,Mat); 559bc011b1eSHong Zhang 560f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAXPY(Mat,PetscScalar,Mat,MatStructure); 56104aac2b0SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAYPX(Mat,PetscScalar,Mat,MatStructure); 562be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompress(Mat); 5631c741599SHong Zhang 564f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScale(Mat,PetscScalar); 565f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShift(Mat,PetscScalar); 5667b80b807SBarry Smith 567be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping); 568be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetLocalToGlobalMappingBlock(Mat,ISLocalToGlobalMapping); 569f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar); 570f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsLocalIS(Mat,IS,PetscScalar); 571be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 572be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 573052efed2SBarry Smith 574be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStashSetInitialSize(Mat,PetscInt,PetscInt); 575be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*); 57690f02eecSBarry Smith 577be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInterpolate(Mat,Vec,Vec); 5780c451bc4SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInterpolateAdd(Mat,Vec,Vec,Vec); 579ba966639SSatish Balay PetscPolymorphicSubroutine(MatInterpolateAdd,(Mat A,Vec x,Vec y),(A,x,y,y)) 580be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestrict(Mat,Vec,Vec); 581be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetVecs(Mat,Vec*,Vec*); 58269db28dcSHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRedundantMatrix(Mat,PetscInt,MPI_Comm,PetscInt,MatReuse,Mat*); 583bd481603SBarry Smith 584bd481603SBarry Smith /*MC 5850d2a7ff2SSatish Balay MatPreallocateInitialize - Begins the block of code that will count the number of nonzeros per 586bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 587bd481603SBarry Smith 588bd481603SBarry Smith Synopsis: 589c1ac3661SBarry Smith PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz) 590bd481603SBarry Smith 591bd481603SBarry Smith Collective on MPI_Comm 592bd481603SBarry Smith 593bd481603SBarry Smith Input Parameters: 594bd481603SBarry Smith + comm - the communicator that will share the eventually allocated matrix 595859fcb39SBarry Smith . nrows - the number of LOCAL rows in the matrix 596859fcb39SBarry Smith - ncols - the number of LOCAL columns in the matrix 597bd481603SBarry Smith 598bd481603SBarry Smith Output Parameters: 599bd481603SBarry Smith + dnz - the array that will be passed to the matrix preallocation routines 600bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 601bd481603SBarry Smith 602bd481603SBarry Smith 603bd481603SBarry Smith Level: intermediate 604bd481603SBarry Smith 605bd481603SBarry Smith Notes: 606bd481603SBarry Smith See the chapter in the users manual on performance for more details 607bd481603SBarry Smith 6081d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 609bd481603SBarry Smith 610bd481603SBarry Smith Use MatPreallocateInitializeSymmetric() for symmetric matrices (MPISBAIJ matrices) 611bd481603SBarry Smith 612bd481603SBarry Smith Concepts: preallocation^Matrix 613bd481603SBarry Smith 614bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(), 615bd481603SBarry Smith MatPreallocateInitializeSymmetric(), MatPreallocateSymmetricSetLocal() 616bd481603SBarry Smith M*/ 617c4f061fbSSatish Balay #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \ 6187c922b88SBarry Smith { \ 619c1ac3661SBarry Smith PetscErrorCode _4_ierr; PetscInt __tmp = (nrows),__ctmp = (ncols),__rstart,__start,__end; \ 620c1ac3661SBarry Smith _4_ierr = PetscMalloc(2*__tmp*sizeof(PetscInt),&dnz);CHKERRQ(_4_ierr);onz = dnz + __tmp;\ 621c1ac3661SBarry Smith _4_ierr = PetscMemzero(dnz,2*__tmp*sizeof(PetscInt));CHKERRQ(_4_ierr);\ 622c1ac3661SBarry Smith _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ctmp;\ 623c1ac3661SBarry Smith _4_ierr = MPI_Scan(&__tmp,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __tmp; 6247c922b88SBarry Smith 625bd481603SBarry Smith /*MC 6260d2a7ff2SSatish Balay MatPreallocateSymmetricInitialize - Begins the block of code that will count the number of nonzeros per 627bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 628bd481603SBarry Smith 629bd481603SBarry Smith Synopsis: 630c1ac3661SBarry Smith PetscErrorCode MatPreallocateSymmetricInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz) 631bd481603SBarry Smith 632bd481603SBarry Smith Collective on MPI_Comm 633bd481603SBarry Smith 634bd481603SBarry Smith Input Parameters: 635bd481603SBarry Smith + comm - the communicator that will share the eventually allocated matrix 636859fcb39SBarry Smith . nrows - the number of LOCAL rows in the matrix 637859fcb39SBarry Smith - ncols - the number of LOCAL columns in the matrix 638bd481603SBarry Smith 639bd481603SBarry Smith Output Parameters: 640bd481603SBarry Smith + dnz - the array that will be passed to the matrix preallocation routines 641bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 642bd481603SBarry Smith 643bd481603SBarry Smith 644bd481603SBarry Smith Level: intermediate 645bd481603SBarry Smith 646bd481603SBarry Smith Notes: 647bd481603SBarry Smith See the chapter in the users manual on performance for more details 648bd481603SBarry Smith 6491d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 650bd481603SBarry Smith 651bd481603SBarry Smith Concepts: preallocation^Matrix 652bd481603SBarry Smith 653bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(), 654bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal() 655bd481603SBarry Smith M*/ 656222b16d4SBarry Smith #define MatPreallocateSymmetricInitialize(comm,nrows,ncols,dnz,onz) 0; \ 657222b16d4SBarry Smith { \ 658c1ac3661SBarry Smith PetscErrorCode _4_ierr; PetscInt __tmp = (nrows),__ctmp = (ncols),__rstart,__end; \ 659c1ac3661SBarry Smith _4_ierr = PetscMalloc(2*__tmp*sizeof(PetscInt),&dnz);CHKERRQ(_4_ierr);onz = dnz + __tmp;\ 660c1ac3661SBarry Smith _4_ierr = PetscMemzero(dnz,2*__tmp*sizeof(PetscInt));CHKERRQ(_4_ierr);\ 661c1ac3661SBarry Smith _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr);\ 662c1ac3661SBarry Smith _4_ierr = MPI_Scan(&__tmp,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __tmp; 663222b16d4SBarry Smith 664bd481603SBarry Smith /*MC 665bd481603SBarry Smith MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be 666bd481603SBarry Smith inserted using a local number of the rows and columns 667bd481603SBarry Smith 668bd481603SBarry Smith Synopsis: 669c1ac3661SBarry Smith PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 670bd481603SBarry Smith 671bd481603SBarry Smith Not Collective 672bd481603SBarry Smith 673bd481603SBarry Smith Input Parameters: 674bd481603SBarry Smith + map - the mapping between local numbering and global numbering 675bd481603SBarry Smith . nrows - the number of rows indicated 6761d73ed98SBarry Smith . rows - the indices of the rows 677bd481603SBarry Smith . ncols - the number of columns in the matrix 678bd481603SBarry Smith . cols - the columns indicated 679bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 680bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 681bd481603SBarry Smith 682bd481603SBarry Smith 683bd481603SBarry Smith Level: intermediate 684bd481603SBarry Smith 685bd481603SBarry Smith Notes: 686bd481603SBarry Smith See the chapter in the users manual on performance for more details 687bd481603SBarry Smith 6881d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 689bd481603SBarry Smith 690bd481603SBarry Smith Concepts: preallocation^Matrix 691bd481603SBarry Smith 692bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(), 693bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal() 694bd481603SBarry Smith M*/ 695c4f061fbSSatish Balay #define MatPreallocateSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\ 696c4f061fbSSatish Balay {\ 697c1ac3661SBarry Smith PetscInt __l;\ 698ef66eb69SBarry Smith _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\ 699ef66eb69SBarry Smith _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\ 700c4f061fbSSatish Balay for (__l=0;__l<nrows;__l++) {\ 701ef66eb69SBarry Smith _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 702c4f061fbSSatish Balay }\ 703c4f061fbSSatish Balay } 704c4f061fbSSatish Balay 705bd481603SBarry Smith /*MC 706bd481603SBarry Smith MatPreallocateSymmetricSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be 707bd481603SBarry Smith inserted using a local number of the rows and columns 708bd481603SBarry Smith 709bd481603SBarry Smith Synopsis: 710c1ac3661SBarry Smith PetscErrorCode MatPreallocateSymmetricSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 711bd481603SBarry Smith 712bd481603SBarry Smith Not Collective 713bd481603SBarry Smith 714bd481603SBarry Smith Input Parameters: 715bd481603SBarry Smith + map - the mapping between local numbering and global numbering 716bd481603SBarry Smith . nrows - the number of rows indicated 7171d73ed98SBarry Smith . rows - the indices of the rows 718bd481603SBarry Smith . ncols - the number of columns in the matrix 719bd481603SBarry Smith . cols - the columns indicated 720bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 721bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 722bd481603SBarry Smith 723bd481603SBarry Smith 724bd481603SBarry Smith Level: intermediate 725bd481603SBarry Smith 726bd481603SBarry Smith Notes: 727bd481603SBarry Smith See the chapter in the users manual on performance for more details 728bd481603SBarry Smith 729bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 730bd481603SBarry Smith 731bd481603SBarry Smith Concepts: preallocation^Matrix 732bd481603SBarry Smith 733bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(), 734bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal() 735bd481603SBarry Smith M*/ 736d3d32019SBarry Smith #define MatPreallocateSymmetricSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\ 737d3d32019SBarry Smith {\ 738c1ac3661SBarry Smith PetscInt __l;\ 739d3d32019SBarry Smith _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\ 740d3d32019SBarry Smith _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\ 741d3d32019SBarry Smith for (__l=0;__l<nrows;__l++) {\ 742d3d32019SBarry Smith _4_ierr = MatPreallocateSymmetricSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 743d3d32019SBarry Smith }\ 744d3d32019SBarry Smith } 745d3d32019SBarry Smith 746bd481603SBarry Smith /*MC 747bd481603SBarry Smith MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be 748bd481603SBarry Smith inserted using a local number of the rows and columns 749bd481603SBarry Smith 750bd481603SBarry Smith Synopsis: 751c1ac3661SBarry Smith PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 752bd481603SBarry Smith 753bd481603SBarry Smith Not Collective 754bd481603SBarry Smith 755bd481603SBarry Smith Input Parameters: 75664075487SBarry Smith + row - the row 757bd481603SBarry Smith . ncols - the number of columns in the matrix 758a50f8bf6SHong Zhang - cols - the columns indicated 759a50f8bf6SHong Zhang 760a50f8bf6SHong Zhang Output Parameters: 761a50f8bf6SHong Zhang + dnz - the array that will be passed to the matrix preallocation routines 762bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 763bd481603SBarry Smith 764bd481603SBarry Smith 765bd481603SBarry Smith Level: intermediate 766bd481603SBarry Smith 767bd481603SBarry Smith Notes: 768bd481603SBarry Smith See the chapter in the users manual on performance for more details 769bd481603SBarry Smith 770bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 771bd481603SBarry Smith 772bd481603SBarry Smith Concepts: preallocation^Matrix 773bd481603SBarry Smith 774bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(), 775bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal() 776bd481603SBarry Smith M*/ 777c4f061fbSSatish Balay #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\ 778c1ac3661SBarry Smith { PetscInt __i; \ 7798ee2e534SBarry Smith if (row < __rstart) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Trying to set preallocation for row %D less than first local row %D",row,__rstart);\ 7808ee2e534SBarry 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);\ 7817c922b88SBarry Smith for (__i=0; __i<nc; __i++) {\ 78264075487SBarry Smith if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \ 7837cc688d7SBarry Smith else dnz[row - __rstart]++;\ 7847c922b88SBarry Smith }\ 7857c922b88SBarry Smith } 7867c922b88SBarry Smith 787bd481603SBarry Smith /*MC 788bd481603SBarry Smith MatPreallocateSymmetricSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be 789bd481603SBarry Smith inserted using a local number of the rows and columns 790bd481603SBarry Smith 791bd481603SBarry Smith Synopsis: 792c1ac3661SBarry Smith PetscErrorCode MatPreallocateSymmetricSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 793bd481603SBarry Smith 794bd481603SBarry Smith Not Collective 795bd481603SBarry Smith 796bd481603SBarry Smith Input Parameters: 797bd481603SBarry Smith + nrows - the number of rows indicated 7981d73ed98SBarry Smith . rows - the indices of the rows 799bd481603SBarry Smith . ncols - the number of columns in the matrix 800bd481603SBarry Smith . cols - the columns indicated 801bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 802bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 803bd481603SBarry Smith 804bd481603SBarry Smith 805bd481603SBarry Smith Level: intermediate 806bd481603SBarry Smith 807bd481603SBarry Smith Notes: 808bd481603SBarry Smith See the chapter in the users manual on performance for more details 809bd481603SBarry Smith 810bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 811bd481603SBarry Smith 812bd481603SBarry Smith Concepts: preallocation^Matrix 813bd481603SBarry Smith 814bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(), 815bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal() 816bd481603SBarry Smith M*/ 817d3d32019SBarry Smith #define MatPreallocateSymmetricSet(row,nc,cols,dnz,onz) 0;\ 818c1ac3661SBarry Smith { PetscInt __i; \ 819d3d32019SBarry Smith for (__i=0; __i<nc; __i++) {\ 820d3d32019SBarry Smith if (cols[__i] >= __end) onz[row - __rstart]++; \ 821d3d32019SBarry Smith else if (cols[__i] >= row) dnz[row - __rstart]++;\ 822d3d32019SBarry Smith }\ 823d3d32019SBarry Smith } 824d3d32019SBarry Smith 825bd481603SBarry Smith /*MC 8260d2a7ff2SSatish Balay MatPreallocateFinalize - Ends the block of code that will count the number of nonzeros per 827bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 828bd481603SBarry Smith 829bd481603SBarry Smith Synopsis: 830c1ac3661SBarry Smith PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz) 831bd481603SBarry Smith 832bd481603SBarry Smith Collective on MPI_Comm 833bd481603SBarry Smith 834bd481603SBarry Smith Input Parameters: 835bd481603SBarry Smith + dnz - the array that will be passed to the matrix preallocation routines 836bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 837bd481603SBarry Smith 838bd481603SBarry Smith 839bd481603SBarry Smith Level: intermediate 840bd481603SBarry Smith 841bd481603SBarry Smith Notes: 842bd481603SBarry Smith See the chapter in the users manual on performance for more details 843bd481603SBarry Smith 844bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 845bd481603SBarry Smith 846bd481603SBarry Smith Concepts: preallocation^Matrix 847bd481603SBarry Smith 848bd481603SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(), 849bd481603SBarry Smith MatPreallocateSymmetricInitialize(), MatPreallocateSymmetricSetLocal() 850bd481603SBarry Smith M*/ 851ef66eb69SBarry Smith #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree(dnz);CHKERRQ(_4_ierr);} 8527c922b88SBarry Smith 853bd481603SBarry Smith 854bd481603SBarry Smith 8557b80b807SBarry Smith /* Routines unique to particular data structures */ 856be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellGetContext(Mat,void **); 857ba966639SSatish Balay PetscPolymorphicFunction(MatShellGetContext,(Mat A),(A,&t),void*,t) 858be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatBDiagGetData(Mat,PetscInt*,PetscInt*,PetscInt*[],PetscInt*[],PetscScalar***); 859698d4c6aSKris Buschelman 860be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInodeAdjustForInodes(Mat,IS*,IS*); 861be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *); 862698d4c6aSKris Buschelman 863be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetColumnIndices(Mat,PetscInt[]); 864be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetColumnIndices(Mat,PetscInt[]); 865be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 866c75a6043SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 867c75a6043SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 8687b80b807SBarry Smith 869a96a251dSBarry Smith #define MAT_SKIP_ALLOCATION -4 870a96a251dSBarry Smith 871be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]); 872ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz)) 873be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]); 874ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqSBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz)) 875be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]); 876ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqAIJSetPreallocation,(Mat A,const PetscInt nnz[]),(A,0,nnz)) 877be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDensePreallocation(Mat,PetscScalar[]); 878be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBDiagSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscScalar*[]); 879273d9f13SBarry Smith 880be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 881ba966639SSatish Balay PetscPolymorphicSubroutine(MatMPIBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[],const PetscInt onz[]),(A,bs,0,nnz,0,onz)) 882be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 883be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 88453803f26SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []); 885be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]); 886be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 887be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIDensePreallocation(Mat,PetscScalar[]); 888be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBDiagSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscScalar*[]); 889be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]); 890be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIDenseSetPreallocation(Mat,PetscScalar[]); 891284134d9SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation(Mat,PetscScalar[]); 892be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIRowbsSetPreallocation(Mat,PetscInt,const PetscInt[]); 893be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,PetscInt*[]); 894be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,PetscInt*[]); 895be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAdicSetLocalFunction(Mat,void (*)(void)); 896273d9f13SBarry Smith 897be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetLDA(Mat,PetscInt); 898ab92ecdeSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDenseGetLocalMatrix(Mat,Mat*); 8991b807ce4Svictorle 900be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStoreValues(Mat); 901be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRetrieveValues(Mat); 9022e8a6d31SBarry Smith 903be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDAADSetCtx(Mat,void*); 9043a7fca6bSBarry Smith 9057b80b807SBarry Smith /* 9067b80b807SBarry Smith These routines are not usually accessed directly, rather solving is 90794b7f48cSBarry Smith done through the KSP and PC interfaces. 9087b80b807SBarry Smith */ 9097b80b807SBarry Smith 910d9274352SBarry Smith /*E 911d9274352SBarry Smith MatOrderingType - String with the name of a PETSc matrix ordering or the creation function 912d9274352SBarry Smith with an optional dynamic library name, for example 913d9274352SBarry Smith http://www.mcs.anl.gov/petsc/lib.a:orderingcreate() 914d9274352SBarry Smith 915d9274352SBarry Smith Level: beginner 916d9274352SBarry Smith 917*e5a9bf91SBarry Smith Cannot use const because the PC objects manipulate the string 918*e5a9bf91SBarry Smith 919d9274352SBarry Smith .seealso: MatGetOrdering() 920d9274352SBarry Smith E*/ 92149773a63SBarry Smith #define MatOrderingType char* 922b12f92e5SBarry Smith #define MATORDERING_NATURAL "natural" 923b12f92e5SBarry Smith #define MATORDERING_ND "nd" 924b12f92e5SBarry Smith #define MATORDERING_1WD "1wd" 925b12f92e5SBarry Smith #define MATORDERING_RCM "rcm" 926b12f92e5SBarry Smith #define MATORDERING_QMD "qmd" 927b12f92e5SBarry Smith #define MATORDERING_ROWLENGTH "rowlength" 92862152c8bSBarry Smith #define MATORDERING_DSC_ND "dsc_nd" 92962152c8bSBarry Smith #define MATORDERING_DSC_MMD "dsc_mmd" 93062152c8bSBarry Smith #define MATORDERING_DSC_MDF "dsc_mdf" 931c06d978dSMatthew Knepley #define MATORDERING_CONSTRAINED "constrained" 932c06d978dSMatthew Knepley #define MATORDERING_IDENTITY "identity" 933c06d978dSMatthew Knepley #define MATORDERING_REVERSE "reverse" 934b12f92e5SBarry Smith 935be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOrdering(Mat,const MatOrderingType,IS*,IS*); 936be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOrderingList(PetscFList *list); 937be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,const MatOrderingType,IS*,IS*)); 93830de9b25SBarry Smith 93930de9b25SBarry Smith /*MC 94030de9b25SBarry Smith MatOrderingRegisterDynamic - Adds a new sparse matrix ordering to the 94130de9b25SBarry Smith matrix package. 94230de9b25SBarry Smith 94330de9b25SBarry Smith Synopsis: 944c1ac3661SBarry Smith PetscErrorCode MatOrderingRegisterDynamic(char *name_ordering,char *path,char *name_create,PetscErrorCode (*routine_create)(MatOrdering)) 94530de9b25SBarry Smith 94630de9b25SBarry Smith Not Collective 94730de9b25SBarry Smith 94830de9b25SBarry Smith Input Parameters: 94930de9b25SBarry Smith + sname - name of ordering (for example MATORDERING_ND) 95030de9b25SBarry Smith . path - location of library where creation routine is 95130de9b25SBarry Smith . name - name of function that creates the ordering type,a string 95230de9b25SBarry Smith - function - function pointer that creates the ordering 95330de9b25SBarry Smith 95430de9b25SBarry Smith Level: developer 95530de9b25SBarry Smith 95630de9b25SBarry Smith If dynamic libraries are used, then the fourth input argument (function) 95730de9b25SBarry Smith is ignored. 95830de9b25SBarry Smith 95930de9b25SBarry Smith Sample usage: 96030de9b25SBarry Smith .vb 96130de9b25SBarry Smith MatOrderingRegisterDynamic("my_order",/home/username/my_lib/lib/libO/solaris/mylib.a, 96230de9b25SBarry Smith "MyOrder",MyOrder); 96330de9b25SBarry Smith .ve 96430de9b25SBarry Smith 96530de9b25SBarry Smith Then, your partitioner can be chosen with the procedural interface via 96630de9b25SBarry Smith $ MatOrderingSetType(part,"my_order) 96730de9b25SBarry Smith or at runtime via the option 9682401956bSBarry Smith $ -pc_factor_mat_ordering_type my_order 96930de9b25SBarry Smith 970ab901514SBarry Smith ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values. 97130de9b25SBarry Smith 97230de9b25SBarry Smith .keywords: matrix, ordering, register 97330de9b25SBarry Smith 97430de9b25SBarry Smith .seealso: MatOrderingRegisterDestroy(), MatOrderingRegisterAll() 97530de9b25SBarry Smith M*/ 976aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 977f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,0) 978b12f92e5SBarry Smith #else 979f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,d) 980b12f92e5SBarry Smith #endif 98130de9b25SBarry Smith 982be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegisterDestroy(void); 983be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegisterAll(const char[]); 9842bad1931SBarry Smith extern PetscTruth MatOrderingRegisterAllCalled; 985b0a32e0cSBarry Smith extern PetscFList MatOrderingList; 986d4fbbf0eSBarry Smith 987be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS); 988a2ce50c7SBarry Smith 989d91e6319SBarry Smith /*S 9902401956bSBarry Smith MatFactorInfo - Data passed into the matrix factorization routines 991b00f7748SHong Zhang 99215e8a5b3SHong Zhang In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE 993b00f7748SHong Zhang 99415e8a5b3SHong Zhang Notes: These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC. 995b00f7748SHong Zhang 996b00f7748SHong Zhang Level: developer 997b00f7748SHong Zhang 998d7d82daaSBarry Smith .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(), 999d7d82daaSBarry Smith MatFactorInfoInitialize() 1000b00f7748SHong Zhang 1001b00f7748SHong Zhang S*/ 1002b00f7748SHong Zhang typedef struct { 10030a29876aSHong Zhang PetscReal shiftnz; /* scaling of identity added to matrix to prevent zero pivots */ 1004fbf22428SSatish Balay PetscReal shiftpd; /* if true, shift until positive pivots */ 10052cea7109SBarry Smith PetscReal shift_fraction; /* record shift fraction taken */ 100615e8a5b3SHong Zhang PetscReal diagonal_fill; /* force diagonal to fill in if initially not filled */ 100715e8a5b3SHong Zhang PetscReal dt; /* drop tolerance */ 1008b00f7748SHong Zhang PetscReal dtcol; /* tolerance for pivoting */ 100915e8a5b3SHong Zhang PetscReal dtcount; /* maximum nonzeros to be allowed per row */ 1010f6275e2eSBarry Smith PetscReal fill; /* expected fill; nonzeros in factored matrix/nonzeros in original matrix*/ 1011348344bbSBarry Smith PetscReal levels; /* ICC/ILU(levels) */ 1012bcd9e38bSBarry Smith PetscReal pivotinblocks; /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0 1013bcd9e38bSBarry Smith factorization may be faster if do not pivot */ 101415e8a5b3SHong Zhang PetscReal zeropivot; /* pivot is called zero if less than this */ 101515e8a5b3SHong Zhang } MatFactorInfo; 1016ffa6d0a5SLois Curfman McInnes 1017be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFactorInfoInitialize(MatFactorInfo*); 1018be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactor(Mat,IS,MatFactorInfo*); 1019be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactorSymbolic(Mat,IS,MatFactorInfo*,Mat*); 1020be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactorNumeric(Mat,MatFactorInfo*,Mat*); 1021be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactor(Mat,IS,IS,MatFactorInfo*); 1022be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUFactor(Mat,IS,IS,MatFactorInfo*); 1023be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactorSymbolic(Mat,IS,IS,MatFactorInfo*,Mat*); 1024be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUFactorSymbolic(Mat,IS,IS,MatFactorInfo*,Mat*); 1025be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatICCFactorSymbolic(Mat,IS,MatFactorInfo*,Mat*); 1026be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatICCFactor(Mat,IS,MatFactorInfo*); 1027be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactorNumeric(Mat,MatFactorInfo*,Mat*); 1028be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUDTFactor(Mat,IS,IS,MatFactorInfo*,Mat *); 1029be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*); 1030be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolve(Mat,Vec,Vec); 1031be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatForwardSolve(Mat,Vec,Vec); 1032be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatBackwardSolve(Mat,Vec,Vec); 1033be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveAdd(Mat,Vec,Vec,Vec); 1034be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveTranspose(Mat,Vec,Vec); 1035be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveTransposeAdd(Mat,Vec,Vec,Vec); 1036be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolves(Mat,Vecs,Vecs); 10378ed539a5SBarry Smith 1038be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUnfactored(Mat); 1039bb5a7306SBarry Smith 1040d91e6319SBarry Smith /*E 1041d91e6319SBarry Smith MatSORType - What type of (S)SOR to perform 1042bb1eb677SSatish Balay 1043d91e6319SBarry Smith Level: beginner 1044d91e6319SBarry Smith 1045d9274352SBarry Smith May be bitwise ORd together 1046d9274352SBarry Smith 1047d91e6319SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 1048d91e6319SBarry Smith 10494e7234bfSBarry Smith MatSORType may be bitwise ORd together, so do not change the numbers 10504e7234bfSBarry Smith 1051a1d92eedSBarry Smith .seealso: MatRelax(), MatPBRelax() 1052d91e6319SBarry Smith E*/ 1053ee50ffe9SBarry Smith typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3, 1054ee50ffe9SBarry Smith SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8, 1055ee50ffe9SBarry Smith SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16, 105684cb2905SBarry Smith SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType; 1057be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRelax(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec); 1058be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPBRelax(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec); 10598ed539a5SBarry Smith 1060d4fbbf0eSBarry Smith /* 1061639f9d9dSBarry Smith These routines are for efficiently computing Jacobians via finite differences. 1062639f9d9dSBarry Smith */ 1063b12f92e5SBarry Smith 1064d9274352SBarry Smith /*E 1065d9274352SBarry Smith MatColoringType - String with the name of a PETSc matrix coloring or the creation function 1066d9274352SBarry Smith with an optional dynamic library name, for example 1067d9274352SBarry Smith http://www.mcs.anl.gov/petsc/lib.a:coloringcreate() 1068d9274352SBarry Smith 1069d9274352SBarry Smith Level: beginner 1070d9274352SBarry Smith 1071d9274352SBarry Smith .seealso: MatGetColoring() 1072d9274352SBarry Smith E*/ 1073*e5a9bf91SBarry Smith #define MatColoringType const char* 1074b12f92e5SBarry Smith #define MATCOLORING_NATURAL "natural" 1075b12f92e5SBarry Smith #define MATCOLORING_SL "sl" 1076b12f92e5SBarry Smith #define MATCOLORING_LF "lf" 1077b12f92e5SBarry Smith #define MATCOLORING_ID "id" 1078b12f92e5SBarry Smith 1079be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColoring(Mat,const MatColoringType,ISColoring*); 1080be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,const MatColoringType,ISColoring *)); 108130de9b25SBarry Smith 108230de9b25SBarry Smith /*MC 108330de9b25SBarry Smith MatColoringRegisterDynamic - Adds a new sparse matrix coloring to the 108430de9b25SBarry Smith matrix package. 108530de9b25SBarry Smith 108630de9b25SBarry Smith Synopsis: 1087c1ac3661SBarry Smith PetscErrorCode MatColoringRegisterDynamic(char *name_coloring,char *path,char *name_create,PetscErrorCode (*routine_create)(MatColoring)) 108830de9b25SBarry Smith 108930de9b25SBarry Smith Not Collective 109030de9b25SBarry Smith 109130de9b25SBarry Smith Input Parameters: 109230de9b25SBarry Smith + sname - name of Coloring (for example MATCOLORING_SL) 109330de9b25SBarry Smith . path - location of library where creation routine is 109430de9b25SBarry Smith . name - name of function that creates the Coloring type, a string 109530de9b25SBarry Smith - function - function pointer that creates the coloring 109630de9b25SBarry Smith 109730de9b25SBarry Smith Level: developer 109830de9b25SBarry Smith 109930de9b25SBarry Smith If dynamic libraries are used, then the fourth input argument (function) 110030de9b25SBarry Smith is ignored. 110130de9b25SBarry Smith 110230de9b25SBarry Smith Sample usage: 110330de9b25SBarry Smith .vb 110430de9b25SBarry Smith MatColoringRegisterDynamic("my_color",/home/username/my_lib/lib/libO/solaris/mylib.a, 110530de9b25SBarry Smith "MyColor",MyColor); 110630de9b25SBarry Smith .ve 110730de9b25SBarry Smith 110830de9b25SBarry Smith Then, your partitioner can be chosen with the procedural interface via 110930de9b25SBarry Smith $ MatColoringSetType(part,"my_color") 111030de9b25SBarry Smith or at runtime via the option 111130de9b25SBarry Smith $ -mat_coloring_type my_color 111230de9b25SBarry Smith 1113ab901514SBarry Smith $PETSC_ARCH occuring in pathname will be replaced with appropriate values. 111430de9b25SBarry Smith 111530de9b25SBarry Smith .keywords: matrix, Coloring, register 111630de9b25SBarry Smith 111730de9b25SBarry Smith .seealso: MatColoringRegisterDestroy(), MatColoringRegisterAll() 111830de9b25SBarry Smith M*/ 1119aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 1120f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,0) 1121b12f92e5SBarry Smith #else 1122f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,d) 1123b12f92e5SBarry Smith #endif 112430de9b25SBarry Smith 11252bad1931SBarry Smith extern PetscTruth MatColoringRegisterAllCalled; 1126f1144a30SSatish Balay 1127be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegisterAll(const char[]); 1128be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegisterDestroy(void); 1129be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*); 1130639f9d9dSBarry Smith 1131d9274352SBarry Smith /*S 1132d9274352SBarry Smith MatFDColoring - Object for computing a sparse Jacobian via finite differences 1133d9274352SBarry Smith and coloring 1134639f9d9dSBarry Smith 1135d9274352SBarry Smith Level: beginner 1136d9274352SBarry Smith 1137d9274352SBarry Smith Concepts: coloring, sparse Jacobian, finite differences 1138d9274352SBarry Smith 1139d9274352SBarry Smith .seealso: MatFDColoringCreate() 1140d9274352SBarry Smith S*/ 1141e2a1c21fSSatish Balay typedef struct _p_MatFDColoring* MatFDColoring; 1142639f9d9dSBarry Smith 1143be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringCreate(Mat,ISColoring,MatFDColoring *); 1144be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringDestroy(MatFDColoring); 1145be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringView(MatFDColoring,PetscViewer); 1146be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*); 11474a9d489dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**); 1148be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal); 1149be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFrequency(MatFDColoring,PetscInt); 1150be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetFrequency(MatFDColoring,PetscInt*); 1151be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFromOptions(MatFDColoring); 1152be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApply(Mat,MatFDColoring,Vec,MatStructure*,void *); 1153be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApplyTS(Mat,MatFDColoring,PetscReal,Vec,MatStructure*,void *); 1154be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetRecompute(MatFDColoring); 1155be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetF(MatFDColoring,Vec); 1156be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,PetscInt*[]); 1157639f9d9dSBarry Smith /* 11580752156aSBarry Smith These routines are for partitioning matrices: currently used only 11593eda8832SBarry Smith for adjacency matrix, MatCreateMPIAdj(). 11600752156aSBarry Smith */ 1161ca161407SBarry Smith 1162d9274352SBarry Smith /*S 1163d9274352SBarry Smith MatPartitioning - Object for managing the partitioning of a matrix or graph 1164d9274352SBarry Smith 1165d9274352SBarry Smith Level: beginner 1166d9274352SBarry Smith 1167d9274352SBarry Smith Concepts: partitioning 1168d9274352SBarry Smith 1169743a1026Svictorle .seealso: MatPartitioningCreate(), MatPartitioningType 1170d9274352SBarry Smith S*/ 117191e9ee9fSBarry Smith typedef struct _p_MatPartitioning* MatPartitioning; 1172d9274352SBarry Smith 1173d9274352SBarry Smith /*E 11745bc6e7ccSvictorle MatPartitioningType - String with the name of a PETSc matrix partitioning or the creation function 1175d9274352SBarry Smith with an optional dynamic library name, for example 1176d9274352SBarry Smith http://www.mcs.anl.gov/petsc/lib.a:partitioningcreate() 1177d9274352SBarry Smith 1178d9274352SBarry Smith Level: beginner 1179d9274352SBarry Smith 1180b547a13bSvictorle .seealso: MatPartitioningCreate(), MatPartitioning 1181d9274352SBarry Smith E*/ 1182*e5a9bf91SBarry Smith #define MatPartitioningType const char* 11838ba1e511SMatthew Knepley #define MAT_PARTITIONING_CURRENT "current" 1184c2ec2a73SHong Zhang #define MAT_PARTITIONING_SQUARE "square" 11858ba1e511SMatthew Knepley #define MAT_PARTITIONING_PARMETIS "parmetis" 118671306c60Sjroman #define MAT_PARTITIONING_CHACO "chaco" 118771306c60Sjroman #define MAT_PARTITIONING_JOSTLE "jostle" 118871306c60Sjroman #define MAT_PARTITIONING_PARTY "party" 118971306c60Sjroman #define MAT_PARTITIONING_SCOTCH "scotch" 119071306c60Sjroman 1191ca161407SBarry Smith 1192be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningCreate(MPI_Comm,MatPartitioning*); 1193be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetType(MatPartitioning,const MatPartitioningType); 1194be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetNParts(MatPartitioning,PetscInt); 1195be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetAdjacency(MatPartitioning,Mat); 1196be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]); 1197be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []); 1198be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningApply(MatPartitioning,IS*); 1199be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningDestroy(MatPartitioning); 12002aabb6bbSBarry Smith 1201be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatPartitioning)); 120230de9b25SBarry Smith 120330de9b25SBarry Smith /*MC 120430de9b25SBarry Smith MatPartitioningRegisterDynamic - Adds a new sparse matrix partitioning to the 120530de9b25SBarry Smith matrix package. 120630de9b25SBarry Smith 120730de9b25SBarry Smith Synopsis: 1208c1ac3661SBarry Smith PetscErrorCode MatPartitioningRegisterDynamic(char *name_partitioning,char *path,char *name_create,PetscErrorCode (*routine_create)(MatPartitioning)) 120930de9b25SBarry Smith 121030de9b25SBarry Smith Not Collective 121130de9b25SBarry Smith 121230de9b25SBarry Smith Input Parameters: 121330de9b25SBarry Smith + sname - name of partitioning (for example MAT_PARTITIONING_CURRENT) or parmetis 121430de9b25SBarry Smith . path - location of library where creation routine is 121530de9b25SBarry Smith . name - name of function that creates the partitioning type, a string 121630de9b25SBarry Smith - function - function pointer that creates the partitioning type 121730de9b25SBarry Smith 121830de9b25SBarry Smith Level: developer 121930de9b25SBarry Smith 122030de9b25SBarry Smith If dynamic libraries are used, then the fourth input argument (function) 122130de9b25SBarry Smith is ignored. 122230de9b25SBarry Smith 122330de9b25SBarry Smith Sample usage: 122430de9b25SBarry Smith .vb 122530de9b25SBarry Smith MatPartitioningRegisterDynamic("my_part",/home/username/my_lib/lib/libO/solaris/mylib.a, 122630de9b25SBarry Smith "MyPartCreate",MyPartCreate); 122730de9b25SBarry Smith .ve 122830de9b25SBarry Smith 122930de9b25SBarry Smith Then, your partitioner can be chosen with the procedural interface via 123030de9b25SBarry Smith $ MatPartitioningSetType(part,"my_part") 123130de9b25SBarry Smith or at runtime via the option 123230de9b25SBarry Smith $ -mat_partitioning_type my_part 123330de9b25SBarry Smith 1234ab901514SBarry Smith $PETSC_ARCH occuring in pathname will be replaced with appropriate values. 123530de9b25SBarry Smith 123630de9b25SBarry Smith .keywords: matrix, partitioning, register 123730de9b25SBarry Smith 123830de9b25SBarry Smith .seealso: MatPartitioningRegisterDestroy(), MatPartitioningRegisterAll() 123930de9b25SBarry Smith M*/ 1240aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 1241f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,0) 12422aabb6bbSBarry Smith #else 1243f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,d) 12442aabb6bbSBarry Smith #endif 12452aabb6bbSBarry Smith 12462bad1931SBarry Smith extern PetscTruth MatPartitioningRegisterAllCalled; 1247f1144a30SSatish Balay 1248be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegisterAll(const char[]); 1249be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegisterDestroy(void); 12502bad1931SBarry Smith 1251be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningView(MatPartitioning,PetscViewer); 1252be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetFromOptions(MatPartitioning); 1253be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningGetType(MatPartitioning,MatPartitioningType*); 1254ca161407SBarry Smith 1255be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningParmetisSetCoarseSequential(MatPartitioning); 12560752156aSBarry Smith 1257be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningJostleSetCoarseLevel(MatPartitioning,PetscReal); 1258be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningJostleSetCoarseSequential(MatPartitioning); 125971306c60Sjroman 1260954b3ec6SBarry Smith typedef enum {MP_CHACO_MULTILEVEL_KL,MP_CHACO_SPECTRAL,MP_CHACO_LINEAR,MP_CHACO_RANDOM, MP_CHACO_SCATTERED} MPChacoGlobalType; 1261be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetGlobal(MatPartitioning, MPChacoGlobalType); 126271306c60Sjroman typedef enum { MP_CHACO_KERNIGHAN_LIN, MP_CHACO_NONE } MPChacoLocalType; 1263be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetLocal(MatPartitioning, MPChacoLocalType); 1264be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal); 126571306c60Sjroman typedef enum { MP_CHACO_LANCZOS, MP_CHACO_RQI_SYMMLQ } MPChacoEigenType; 1266be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType); 1267be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenTol(MatPartitioning, PetscReal); 1268be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenNumber(MatPartitioning, PetscInt); 126971306c60Sjroman 127071306c60Sjroman #define MP_PARTY_OPT "opt" 127171306c60Sjroman #define MP_PARTY_LIN "lin" 127271306c60Sjroman #define MP_PARTY_SCA "sca" 127371306c60Sjroman #define MP_PARTY_RAN "ran" 127471306c60Sjroman #define MP_PARTY_GBF "gbf" 127571306c60Sjroman #define MP_PARTY_GCF "gcf" 127671306c60Sjroman #define MP_PARTY_BUB "bub" 127771306c60Sjroman #define MP_PARTY_DEF "def" 1278be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetGlobal(MatPartitioning, const char*); 127971306c60Sjroman #define MP_PARTY_HELPFUL_SETS "hs" 128071306c60Sjroman #define MP_PARTY_KERNIGHAN_LIN "kl" 128171306c60Sjroman #define MP_PARTY_NONE "no" 1282be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetLocal(MatPartitioning, const char*); 1283be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal); 1284be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetBipart(MatPartitioning,PetscTruth); 1285be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscTruth); 128671306c60Sjroman 128771306c60Sjroman typedef enum { MP_SCOTCH_GREEDY, MP_SCOTCH_GPS, MP_SCOTCH_GR_GPS } MPScotchGlobalType; 1288be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetArch(MatPartitioning,const char*); 1289be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetMultilevel(MatPartitioning); 1290be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetGlobal(MatPartitioning,MPScotchGlobalType); 1291be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetCoarseLevel(MatPartitioning,PetscReal); 1292be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetHostList(MatPartitioning,const char*); 129371306c60Sjroman typedef enum { MP_SCOTCH_KERNIGHAN_LIN, MP_SCOTCH_NONE } MPScotchLocalType; 1294be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetLocal(MatPartitioning,MPScotchLocalType); 1295be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetMapping(MatPartitioning); 1296be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetStrategy(MatPartitioning,char*); 129771306c60Sjroman 12980752156aSBarry Smith /* 12990a835dfdSSatish Balay If you add entries here you must also add them to finclude/petscmat.h 1300d4fbbf0eSBarry Smith */ 13011c1c02c0SLois Curfman McInnes typedef enum { MATOP_SET_VALUES=0, 13021c1c02c0SLois Curfman McInnes MATOP_GET_ROW=1, 13031c1c02c0SLois Curfman McInnes MATOP_RESTORE_ROW=2, 13041c1c02c0SLois Curfman McInnes MATOP_MULT=3, 13051c1c02c0SLois Curfman McInnes MATOP_MULT_ADD=4, 13067c922b88SBarry Smith MATOP_MULT_TRANSPOSE=5, 13077c922b88SBarry Smith MATOP_MULT_TRANSPOSE_ADD=6, 13081c1c02c0SLois Curfman McInnes MATOP_SOLVE=7, 13091c1c02c0SLois Curfman McInnes MATOP_SOLVE_ADD=8, 13107c922b88SBarry Smith MATOP_SOLVE_TRANSPOSE=9, 13117c922b88SBarry Smith MATOP_SOLVE_TRANSPOSE_ADD=10, 13121c1c02c0SLois Curfman McInnes MATOP_LUFACTOR=11, 13131c1c02c0SLois Curfman McInnes MATOP_CHOLESKYFACTOR=12, 13141c1c02c0SLois Curfman McInnes MATOP_RELAX=13, 13151c1c02c0SLois Curfman McInnes MATOP_TRANSPOSE=14, 13161c1c02c0SLois Curfman McInnes MATOP_GETINFO=15, 13171c1c02c0SLois Curfman McInnes MATOP_EQUAL=16, 13181c1c02c0SLois Curfman McInnes MATOP_GET_DIAGONAL=17, 13191c1c02c0SLois Curfman McInnes MATOP_DIAGONAL_SCALE=18, 13201c1c02c0SLois Curfman McInnes MATOP_NORM=19, 13211c1c02c0SLois Curfman McInnes MATOP_ASSEMBLY_BEGIN=20, 13221c1c02c0SLois Curfman McInnes MATOP_ASSEMBLY_END=21, 13231c1c02c0SLois Curfman McInnes MATOP_COMPRESS=22, 13241c1c02c0SLois Curfman McInnes MATOP_SET_OPTION=23, 13251c1c02c0SLois Curfman McInnes MATOP_ZERO_ENTRIES=24, 13261c1c02c0SLois Curfman McInnes MATOP_ZERO_ROWS=25, 13271c1c02c0SLois Curfman McInnes MATOP_LUFACTOR_SYMBOLIC=26, 13281c1c02c0SLois Curfman McInnes MATOP_LUFACTOR_NUMERIC=27, 13291c1c02c0SLois Curfman McInnes MATOP_CHOLESKY_FACTOR_SYMBOLIC=28, 13301c1c02c0SLois Curfman McInnes MATOP_CHOLESKY_FACTOR_NUMERIC=29, 1331d643ce63SMatthew Knepley MATOP_SETUP_PREALLOCATION=30, 1332d643ce63SMatthew Knepley MATOP_ILUFACTOR_SYMBOLIC=31, 1333d643ce63SMatthew Knepley MATOP_ICCFACTOR_SYMBOLIC=32, 1334d643ce63SMatthew Knepley MATOP_GET_ARRAY=33, 1335d643ce63SMatthew Knepley MATOP_RESTORE_ARRAY=34, 1336e26f1c3cSBarry Smith MATOP_DUPLICATE=35, 1337d643ce63SMatthew Knepley MATOP_FORWARD_SOLVE=36, 1338d643ce63SMatthew Knepley MATOP_BACKWARD_SOLVE=37, 1339d643ce63SMatthew Knepley MATOP_ILUFACTOR=38, 1340d643ce63SMatthew Knepley MATOP_ICCFACTOR=39, 1341d643ce63SMatthew Knepley MATOP_AXPY=40, 1342d643ce63SMatthew Knepley MATOP_GET_SUBMATRICES=41, 1343d643ce63SMatthew Knepley MATOP_INCREASE_OVERLAP=42, 1344d643ce63SMatthew Knepley MATOP_GET_VALUES=43, 1345d643ce63SMatthew Knepley MATOP_COPY=44, 1346d643ce63SMatthew Knepley MATOP_PRINT_HELP=45, 1347d643ce63SMatthew Knepley MATOP_SCALE=46, 1348d643ce63SMatthew Knepley MATOP_SHIFT=47, 1349d643ce63SMatthew Knepley MATOP_DIAGONAL_SHIFT=48, 1350d643ce63SMatthew Knepley MATOP_ILUDT_FACTOR=49, 1351d643ce63SMatthew Knepley MATOP_GET_BLOCK_SIZE=50, 1352d643ce63SMatthew Knepley MATOP_GET_ROW_IJ=51, 1353d643ce63SMatthew Knepley MATOP_RESTORE_ROW_IJ=52, 1354d643ce63SMatthew Knepley MATOP_GET_COLUMN_IJ=53, 1355d643ce63SMatthew Knepley MATOP_RESTORE_COLUMN_IJ=54, 1356d643ce63SMatthew Knepley MATOP_FDCOLORING_CREATE=55, 1357d643ce63SMatthew Knepley MATOP_COLORING_PATCH=56, 1358d643ce63SMatthew Knepley MATOP_SET_UNFACTORED=57, 1359d643ce63SMatthew Knepley MATOP_PERMUTE=58, 1360d643ce63SMatthew Knepley MATOP_SET_VALUES_BLOCKED=59, 1361d643ce63SMatthew Knepley MATOP_GET_SUBMATRIX=60, 1362d643ce63SMatthew Knepley MATOP_DESTROY=61, 1363d643ce63SMatthew Knepley MATOP_VIEW=62, 1364d643ce63SMatthew Knepley MATOP_GET_MAPS=63, 1365d643ce63SMatthew Knepley MATOP_USE_SCALED_FORM=64, 1366d643ce63SMatthew Knepley MATOP_SCALE_SYSTEM=65, 1367d643ce63SMatthew Knepley MATOP_UNSCALE_SYSTEM=66, 1368ba16a8cbSKris Buschelman MATOP_SET_LOCAL_TO_GLOBAL_MAP=67, 1369d643ce63SMatthew Knepley MATOP_SET_VALUES_LOCAL=68, 1370d643ce63SMatthew Knepley MATOP_ZERO_ROWS_LOCAL=69, 1371d643ce63SMatthew Knepley MATOP_GET_ROW_MAX=70, 1372d643ce63SMatthew Knepley MATOP_CONVERT=71, 1373d643ce63SMatthew Knepley MATOP_SET_COLORING=72, 1374d643ce63SMatthew Knepley MATOP_SET_VALUES_ADIC=73, 1375d643ce63SMatthew Knepley MATOP_SET_VALUES_ADIFOR=74, 1376d643ce63SMatthew Knepley MATOP_FD_COLORING_APPLY=75, 1377d643ce63SMatthew Knepley MATOP_SET_FROM_OPTIONS=76, 1378ba16a8cbSKris Buschelman MATOP_MULT_CON=77, 1379ba16a8cbSKris Buschelman MATOP_MULT_TRANSPOSE_CON=78, 1380ba16a8cbSKris Buschelman MATOP_ILU_FACTOR_SYMBOLIC_CON=79, 1381d643ce63SMatthew Knepley MATOP_PERMUTE_SPARSIFY=80, 1382d643ce63SMatthew Knepley MATOP_MULT_MULTIPLE=81, 138341acf15aSKris Buschelman MATOP_SOLVE_MULTIPLE=82, 138441acf15aSKris Buschelman MATOP_GET_INERTIA=83, 138541acf15aSKris Buschelman MATOP_LOAD=84, 138641acf15aSKris Buschelman MATOP_IS_SYMMETRIC=85, 138741acf15aSKris Buschelman MATOP_IS_HERMITIAN=86, 138841acf15aSKris Buschelman MATOP_IS_STRUCTURALLY_SYMMETRIC=87, 138941acf15aSKris Buschelman MATOP_PB_RELAX=88, 139041acf15aSKris Buschelman MATOP_GET_VECS=89, 139141acf15aSKris Buschelman MATOP_MAT_MULT=90, 139241acf15aSKris Buschelman MATOP_MAT_MULT_SYMBOLIC=91, 139341acf15aSKris Buschelman MATOP_MAT_MULT_NUMERIC=92, 139441acf15aSKris Buschelman MATOP_PTAP=93, 139541acf15aSKris Buschelman MATOP_PTAP_SYMBOLIC=94, 1396bc011b1eSHong Zhang MATOP_PTAP_NUMERIC=95, 1397bc011b1eSHong Zhang MATOP_MAT_MULTTRANSPOSE=96, 1398bc011b1eSHong Zhang MATOP_MAT_MULTTRANSPOSE_SYMBOLIC=97, 13997ba1a0bfSKris Buschelman MATOP_MAT_MULTTRANSPOSE_NUMERIC=98, 14007ba1a0bfSKris Buschelman MATOP_PTAP_SYMBOLIC_SEQAIJ=99, 14017ba1a0bfSKris Buschelman MATOP_PTAP_NUMERIC_SEQAIJ=100, 14027ba1a0bfSKris Buschelman MATOP_PTAP_SYMBOLIC_MPIAIJ=101, 140387d4246cSBarry Smith MATOP_PTAP_NUMERIC_MPIAIJ=102, 1404f5edf698SHong Zhang MATOP_SET_VALUES_ROW = 105, 1405d5c63f83SSatish Balay MATOP_GET_ROW_UTRIANGULAR=108, 14062bebee5dSHong Zhang MATOP_RESTORE_ROW_UTRIANGULAR=109, 140769db28dcSHong Zhang MATOP_MATSOLVE=110, 140869db28dcSHong Zhang MATOP_GET_REDUNDANTMATRIX=111 1409fae171e0SBarry Smith } MatOperation; 1410be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatHasOperation(Mat,MatOperation,PetscTruth*); 1411be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellSetOperation(Mat,MatOperation,void(*)(void)); 1412be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellGetOperation(Mat,MatOperation,void(**)(void)); 1413be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellSetContext(Mat,void*); 1414112a2221SBarry Smith 141590ace30eSBarry Smith /* 141690ace30eSBarry Smith Codes for matrices stored on disk. By default they are 141790ace30eSBarry Smith stored in a universal format. By changing the format with 1418fb9695e5SSatish Balay PetscViewerSetFormat(viewer,PETSC_VIEWER_BINARY_NATIVE); the matrices will 141990ace30eSBarry Smith be stored in a way natural for the matrix, for example dense matrices 142090ace30eSBarry Smith would be stored as dense. Matrices stored this way may only be 142190ace30eSBarry Smith read into matrices of the same time. 142290ace30eSBarry Smith */ 142390ace30eSBarry Smith #define MATRIX_BINARY_FORMAT_DENSE -1 142490ace30eSBarry Smith 1425be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetHashTableFactor(Mat,PetscReal); 1426be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatISGetLocalMat(Mat,Mat*); 14271f4e1ec7SBarry Smith 1428d9274352SBarry Smith /*S 1429d9274352SBarry Smith MatNullSpace - Object that removes a null space from a vector, i.e. 1430d9274352SBarry Smith orthogonalizes the vector to a subsapce 1431d9274352SBarry Smith 1432f7a9e4ceSBarry Smith Level: advanced 1433d9274352SBarry Smith 1434d9274352SBarry Smith Concepts: matrix; linear operator, null space 1435d9274352SBarry Smith 14366e1639daSBarry Smith Users manual sections: 14376e1639daSBarry Smith . sec_singular 14386e1639daSBarry Smith 1439d9274352SBarry Smith .seealso: MatNullSpaceCreate() 1440d9274352SBarry Smith S*/ 144174637425SBarry Smith typedef struct _p_MatNullSpace* MatNullSpace; 1442d9274352SBarry Smith 1443be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceCreate(MPI_Comm,PetscTruth,PetscInt,const Vec[],MatNullSpace*); 1444281d1b2eSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(Vec,void*),void*); 1445be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceDestroy(MatNullSpace); 1446be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceRemove(MatNullSpace,Vec,Vec*); 1447be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceAttach(Mat,MatNullSpace); 1448be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceTest(MatNullSpace,Mat); 144974637425SBarry Smith 1450be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatReorderingSeqSBAIJ(Mat,IS); 1451be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPISBAIJSetHashTableFactor(Mat,PetscReal); 1452be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetColumnIndices(Mat,PetscInt *); 145346129b97SKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJInvertBlockDiagonal(Mat); 14543f1d51d7SBarry Smith 1455be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMAIJ(Mat,PetscInt,Mat*); 1456be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMAIJRedimension(Mat,PetscInt,Mat*); 1457be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMAIJGetAIJ(Mat,Mat*); 1458c4f061fbSSatish Balay 1459be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatComputeExplicitOperator(Mat,Mat*); 1460b0a32e0cSBarry Smith 1461be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalScaleLocal(Mat,Vec); 146204f1ad80SBarry Smith 1463be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *); 1464be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *); 14657dbadf16SMatthew Knepley 1466e9fa29b7SSatish Balay PETSC_EXTERN_CXX_END 14672eac72dbSBarry Smith #endif 1468