12eac72dbSBarry Smith /* 22eac72dbSBarry Smith Include file for the matrix component of PETSc 32eac72dbSBarry Smith */ 40a835dfdSSatish Balay #ifndef __PETSCMAT_H 50a835dfdSSatish Balay #define __PETSCMAT_H 60a835dfdSSatish Balay #include "petscvec.h" 7e9fa29b7SSatish Balay PETSC_EXTERN_CXX_BEGIN 82eac72dbSBarry Smith 9d9274352SBarry Smith /*S 10d9274352SBarry Smith Mat - Abstract PETSc matrix object 112eac72dbSBarry Smith 12d91e6319SBarry Smith Level: beginner 13d91e6319SBarry Smith 14d9274352SBarry Smith Concepts: matrix; linear operator 15d9274352SBarry Smith 16d9274352SBarry Smith .seealso: MatCreate(), MatType, MatSetType() 17d9274352SBarry Smith S*/ 18d9274352SBarry Smith typedef struct _p_Mat* Mat; 19d9274352SBarry Smith 20d9274352SBarry Smith /*E 21d9274352SBarry Smith MatType - String with the name of a PETSc matrix or the creation function 22d9274352SBarry Smith with an optional dynamic library name, for example 23d9274352SBarry Smith http://www.mcs.anl.gov/petsc/lib.a:mymatcreate() 24d9274352SBarry Smith 25d9274352SBarry Smith Level: beginner 26d9274352SBarry Smith 27d9274352SBarry Smith .seealso: MatSetType(), Mat 28d91e6319SBarry Smith E*/ 29e5a9bf91SBarry Smith #define MatType const char* 30273d9f13SBarry Smith #define MATSAME "same" 31273d9f13SBarry Smith #define MATSEQMAIJ "seqmaij" 32273d9f13SBarry Smith #define MATMPIMAIJ "mpimaij" 33209238afSKris Buschelman #define MATMAIJ "maij" 34273d9f13SBarry Smith #define MATIS "is" 35273d9f13SBarry Smith #define MATMPIROWBS "mpirowbs" 36273d9f13SBarry Smith #define MATSEQAIJ "seqaij" 37273d9f13SBarry Smith #define MATMPIAIJ "mpiaij" 38209238afSKris Buschelman #define MATAIJ "aij" 39273d9f13SBarry Smith #define MATSHELL "shell" 40209238afSKris Buschelman #define MATSEQDENSE "seqdense" 41273d9f13SBarry Smith #define MATMPIDENSE "mpidense" 42209238afSKris Buschelman #define MATDENSE "dense" 43273d9f13SBarry Smith #define MATSEQBAIJ "seqbaij" 44273d9f13SBarry Smith #define MATMPIBAIJ "mpibaij" 45209238afSKris Buschelman #define MATBAIJ "baij" 46273d9f13SBarry Smith #define MATMPIADJ "mpiadj" 47273d9f13SBarry Smith #define MATSEQSBAIJ "seqsbaij" 48273d9f13SBarry Smith #define MATMPISBAIJ "mpisbaij" 49209238afSKris Buschelman #define MATSBAIJ "sbaij" 50cebc7f6cSBarry Smith #define MATDAAD "daad" 51cebc7f6cSBarry Smith #define MATMFFD "mffd" 52c8a8475eSBarry Smith #define MATNORMAL "normal" 53ab92ecdeSBarry Smith #define MATLRC "lrc" 54b3a1e11cSKris Buschelman #define MATSEQAIJSPOOLES "seqaijspooles" 55d10c748bSKris Buschelman #define MATMPIAIJSPOOLES "mpiaijspooles" 569abb65ffSKris Buschelman #define MATSEQSBAIJSPOOLES "seqsbaijspooles" 5722191285SKris Buschelman #define MATMPISBAIJSPOOLES "mpisbaijspooles" 58bb4d4166SHong Zhang #define MATAIJSPOOLES "aijspooles" 59bb4d4166SHong Zhang #define MATSBAIJSPOOLES "sbaijspooles" 6014b396bbSKris Buschelman #define MATSUPERLU "superlu" 611d96aa28SKris Buschelman #define MATSUPERLU_DIST "superlu_dist" 621677d0b8SKris Buschelman #define MATUMFPACK "umfpack" 63e8aa55a4SKris Buschelman #define MATESSL "essl" 644eb8e494SKris Buschelman #define MATLUSOL "lusol" 65397b6df1SKris Buschelman #define MATAIJMUMPS "aijmumps" 66397b6df1SKris Buschelman #define MATSBAIJMUMPS "sbaijmumps" 678da957c5SKris Buschelman #define MATDSCPACK "dscpack" 687065e2aaSKris Buschelman #define MATMATLAB "matlab" 694b8d99adSRichard Tran Mills #define MATSEQCSRPERM "seqcsrperm" 7018def556SRichard Tran Mills #define MATMPICSRPERM "mpicsrperm" 71814655b8SBarry Smith #define MATCSRPERM "csrperm" 7281824310SBarry Smith #define MATSEQCRL "seqcrl" 73c02b24c5SBarry Smith #define MATMPICRL "mpicrl" 74c02b24c5SBarry Smith #define MATCRL "crl" 752a6744ebSBarry Smith #define MATSCATTER "scatter" 76421e10b8SBarry Smith #define MATBLOCKMAT "blockmat" 77793850ffSBarry Smith #define MATCOMPOSITE "composite" 785a7f1df3SHong Zhang #define MATSEQFFTW "seqfftw" 79d91e6319SBarry Smith 809989ab13SBarry Smith /*E 819989ab13SBarry Smith MatSolverType - String with the name of a PETSc matrix solver type. 829989ab13SBarry Smith 839989ab13SBarry Smith For example: "petsc" indicates what PETSc provides, "superlu" indicates either 849989ab13SBarry Smith SuperLU or SuperLU_Dist etc. 859989ab13SBarry Smith 869989ab13SBarry Smith 879989ab13SBarry Smith Level: beginner 889989ab13SBarry Smith 899989ab13SBarry Smith .seealso: MatSetSolverType(), Mat, MatSetType(), MatType 909989ab13SBarry Smith E*/ 919989ab13SBarry Smith #define MatSolverType char* 929989ab13SBarry Smith 93c06d978dSMatthew Knepley /* Logging support */ 94552e946dSBarry Smith #define MAT_FILE_COOKIE 1211216 /* used to indicate matrices in binary files */ 95be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_COOKIE; 96be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_FDCOLORING_COOKIE; 97be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_PARTITIONING_COOKIE; 98be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_NULLSPACE_COOKIE; 99e884886fSBarry Smith extern PetscCookie PETSCMAT_DLLEXPORT MATMFFD_COOKIE; 100c06d978dSMatthew Knepley 101ceb03754SKris Buschelman /*E 102ceb03754SKris Buschelman MatReuse - Indicates if matrices obtained from a previous call to MatGetSubMatrices() 103ceb03754SKris Buschelman or MatGetSubMatrix() are to be reused to store the new matrix values. 104ceb03754SKris Buschelman 105ceb03754SKris Buschelman Level: beginner 106ceb03754SKris Buschelman 107ceb03754SKris Buschelman Any additions/changes here MUST also be made in include/finclude/petscmat.h 108ceb03754SKris Buschelman 1090c451bc4SBarry Smith .seealso: MatGetSubMatrices(), MatGetSubMatrix(), MatDestroyMatrices(), MatConvert() 110ceb03754SKris Buschelman E*/ 111ceb03754SKris Buschelman typedef enum {MAT_INITIAL_MATRIX,MAT_REUSE_MATRIX} MatReuse; 112ceb03754SKris Buschelman 1135494a064SHong Zhang /*E 1145494a064SHong Zhang MatGetSubMatrixOption - Indicates if matrices obtained from a call to MatGetSubMatrices() 115829201f2SHong Zhang include the matrix values. Currently it is only used by MatGetSeqNonzerostructure(). 1165494a064SHong Zhang 1175494a064SHong Zhang Level: beginner 1185494a064SHong Zhang 119829201f2SHong Zhang .seealso: MatGetSeqNonzerostructure() 1205494a064SHong Zhang E*/ 1215494a064SHong Zhang typedef enum {MAT_DO_NOT_GET_VALUES,MAT_GET_VALUES} MatGetSubMatrixOption; 1225494a064SHong Zhang 123e5bd5246SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInitializePackage(const char[]); 124c06d978dSMatthew Knepley 125f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreate(MPI_Comm,Mat*); 126a69a960dSMatthew Knepley PetscPolymorphicFunction(MatCreate,(MPI_Comm comm),(comm,&A),Mat,A) 127a69a960dSMatthew Knepley PetscPolymorphicFunction(MatCreate,(),(PETSC_COMM_WORLD,&A),Mat,A) 128f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetSizes(Mat,PetscInt,PetscInt,PetscInt,PetscInt); 129f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetType(Mat,MatType); 1309989ab13SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetSolverType(Mat,MatSolverType); 1319989ab13SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSolverType(Mat,MatSolverType*); 132be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetFromOptions(Mat); 133be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUpPreallocation(Mat); 134be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRegisterAll(const char[]); 135be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat)); 13630de9b25SBarry Smith 137f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetOptionsPrefix(Mat,const char[]); 138f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAppendOptionsPrefix(Mat,const char[]); 139f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOptionsPrefix(Mat,const char*[]); 140f69a0ea3SMatthew Knepley 14130de9b25SBarry Smith /*MC 14230de9b25SBarry Smith MatRegisterDynamic - Adds a new matrix type 14330de9b25SBarry Smith 14430de9b25SBarry Smith Synopsis: 145c1ac3661SBarry Smith PetscErrorCode MatRegisterDynamic(char *name,char *path,char *name_create,PetscErrorCode (*routine_create)(Mat)) 14630de9b25SBarry Smith 14730de9b25SBarry Smith Not Collective 14830de9b25SBarry Smith 14930de9b25SBarry Smith Input Parameters: 15030de9b25SBarry Smith + name - name of a new user-defined matrix type 15130de9b25SBarry Smith . path - path (either absolute or relative) the library containing this solver 15230de9b25SBarry Smith . name_create - name of routine to create method context 15330de9b25SBarry Smith - routine_create - routine to create method context 15430de9b25SBarry Smith 15530de9b25SBarry Smith Notes: 15630de9b25SBarry Smith MatRegisterDynamic() may be called multiple times to add several user-defined solvers. 15730de9b25SBarry Smith 15830de9b25SBarry Smith If dynamic libraries are used, then the fourth input argument (routine_create) 15930de9b25SBarry Smith is ignored. 16030de9b25SBarry Smith 16130de9b25SBarry Smith Sample usage: 16230de9b25SBarry Smith .vb 16330de9b25SBarry Smith MatRegisterDynamic("my_mat",/home/username/my_lib/lib/libO/solaris/mylib.a, 16430de9b25SBarry Smith "MyMatCreate",MyMatCreate); 16530de9b25SBarry Smith .ve 16630de9b25SBarry Smith 16730de9b25SBarry Smith Then, your solver can be chosen with the procedural interface via 16830de9b25SBarry Smith $ MatSetType(Mat,"my_mat") 16930de9b25SBarry Smith or at runtime via the option 17030de9b25SBarry Smith $ -mat_type my_mat 17130de9b25SBarry Smith 17230de9b25SBarry Smith Level: advanced 17330de9b25SBarry Smith 174ab901514SBarry Smith Notes: ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values. 17530de9b25SBarry Smith If your function is not being put into a shared library then use VecRegister() instead 17630de9b25SBarry Smith 17730de9b25SBarry Smith .keywords: Mat, register 17830de9b25SBarry Smith 17930de9b25SBarry Smith .seealso: MatRegisterAll(), MatRegisterDestroy() 18030de9b25SBarry Smith 18130de9b25SBarry Smith M*/ 182273d9f13SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 183273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,0) 184273d9f13SBarry Smith #else 185273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,d) 18630de9b25SBarry Smith #endif 18730de9b25SBarry Smith 1889c666560SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRegisterAll(const char[]); 189273d9f13SBarry Smith extern PetscTruth MatRegisterAllCalled; 190b0a32e0cSBarry Smith extern PetscFList MatList; 19128988994SBarry Smith 192be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqDense(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*); 193be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIDense(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*); 194be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 195ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,nz,nnz,&A),Mat,A) 196ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,&A),Mat,A) 197ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,0,nnz,&A),Mat,A) 198ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,&A),Mat,A) 199ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,A)) 200ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,m,n,0,nnz,A)) 201ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,A)) 202be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 203ba966639SSatish 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) 204ba966639SSatish 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) 205ba966639SSatish 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) 206ba966639SSatish 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) 207ba966639SSatish 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)) 208ba966639SSatish 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)) 209ba966639SSatish 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)) 210ba966639SSatish 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) 211ba966639SSatish 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) 212ba966639SSatish 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) 213ba966639SSatish 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) 214ba966639SSatish 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)) 215ba966639SSatish 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)) 216ba966639SSatish 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)) 21763ab93bfSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *); 2188d7a6e47SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJWithSplitArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],PetscInt[],PetscInt[],PetscScalar[],Mat*); 2198d7a6e47SBarry Smith 220be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIRowbs(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 221be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 222ba966639SSatish 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) 223ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A) 224ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A) 225ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A) 226ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A)) 227ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A)) 228ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A)) 229be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 230ba966639SSatish 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) 231ba966639SSatish 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) 232ba966639SSatish 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) 233ba966639SSatish 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) 234ba966639SSatish 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)) 235ba966639SSatish 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)) 236ba966639SSatish 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)) 237ba966639SSatish 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) 238ba966639SSatish 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) 239ba966639SSatish 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) 240ba966639SSatish 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) 241ba966639SSatish 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)) 242ba966639SSatish 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)) 243ba966639SSatish 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)) 244be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAdj(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscInt[],Mat*); 245be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 246ba966639SSatish 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) 247ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A) 248ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A) 249ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A) 250ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A)) 251ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A)) 252ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A)) 253be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPISBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 254ba966639SSatish 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) 255ba966639SSatish 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) 256ba966639SSatish 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) 257ba966639SSatish 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) 258ba966639SSatish 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)) 259ba966639SSatish 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)) 260ba966639SSatish 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)) 261ba966639SSatish 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) 262ba966639SSatish 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) 263ba966639SSatish 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) 264ba966639SSatish 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) 265ba966639SSatish 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)) 266ba966639SSatish 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)) 267ba966639SSatish 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)) 268be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateShell(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,void *,Mat*); 269ba966639SSatish 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) 270ba966639SSatish Balay PetscPolymorphicFunction(MatCreateShell,(PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx),(PETSC_COMM_WORLD,m,n,M,N,ctx,&A),Mat,A) 271be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateAdic(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,void (*)(void),Mat*); 272be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateNormal(Mat,Mat*); 273ba966639SSatish Balay PetscPolymorphicFunction(MatCreateNormal,(Mat mat),(mat,&A),Mat,A) 274ab92ecdeSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateLRC(Mat,Mat,Mat,Mat*); 275284134d9SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateIS(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,ISLocalToGlobalMapping,Mat*); 2761472f72bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 2771472f72bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPICRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 2782a6744ebSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateScatter(MPI_Comm,VecScatter,Mat*); 2792a6744ebSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScatterSetVecScatter(Mat,VecScatter); 2802a6744ebSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScatterGetVecScatter(Mat,VecScatter*); 2818cdf2d9bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateBlockMat(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,Mat*); 282793850ffSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompositeAddMat(Mat,Mat); 283b52f573bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompositeMerge(Mat); 284793850ffSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateComposite(MPI_Comm,PetscInt,const Mat*,Mat*); 2855a7f1df3SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqFFTW(MPI_Comm,PetscInt,const PetscInt[],Mat*); 2861472f72bSBarry Smith 287f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUp(Mat); 288be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroy(Mat); 28921c89e3eSBarry Smith 2900c451bc4SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConjugate(Mat); 29199cafbc1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRealPart(Mat); 29299cafbc1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatImaginaryPart(Mat); 29399cafbc1SBarry Smith 2948ed539a5SBarry Smith /* ------------------------------------------------------------*/ 295be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 296be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 29787d4246cSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesRow(Mat,PetscInt,const PetscScalar[]); 298f259bd47SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesRowLocal(Mat,PetscInt,const PetscScalar[]); 29984cb2905SBarry Smith 3002ef4de8bSBarry Smith /*S 3012ef4de8bSBarry Smith MatStencil - Data structure (C struct) for storing information about a single row or 3022ef4de8bSBarry Smith column of a matrix as index on an associated grid. 3032ef4de8bSBarry Smith 3042ef4de8bSBarry Smith Level: beginner 3052ef4de8bSBarry Smith 3062ef4de8bSBarry Smith Concepts: matrix; linear operator 3072ef4de8bSBarry Smith 308d7bd93baSBarry Smith .seealso: MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockStencil() 3092ef4de8bSBarry Smith S*/ 310435da068SBarry Smith typedef struct { 311c1ac3661SBarry Smith PetscInt k,j,i,c; 312435da068SBarry Smith } MatStencil; 3132ef4de8bSBarry Smith 314be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode); 315be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode); 316be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt); 317435da068SBarry Smith 318be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetColoring(Mat,ISColoring); 319be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesAdic(Mat,void*); 320be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesAdifor(Mat,PetscInt,void*); 3213a7fca6bSBarry Smith 322d91e6319SBarry Smith /*E 323d91e6319SBarry Smith MatAssemblyType - Indicates if the matrix is now to be used, or if you plan 324d91e6319SBarry Smith to continue to add values to it 325d91e6319SBarry Smith 326d91e6319SBarry Smith Level: beginner 327d91e6319SBarry Smith 328d91e6319SBarry Smith .seealso: MatAssemblyBegin(), MatAssemblyEnd() 329d91e6319SBarry Smith E*/ 3306d4a8577SBarry Smith typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType; 331be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssemblyBegin(Mat,MatAssemblyType); 332be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssemblyEnd(Mat,MatAssemblyType); 333be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssembled(Mat,PetscTruth*); 3344f9c727eSBarry Smith 3351d73ed98SBarry Smith 33630de9b25SBarry Smith 337d91e6319SBarry Smith /*E 338d91e6319SBarry Smith MatOption - Options that may be set for a matrix and its behavior or storage 339d91e6319SBarry Smith 340d91e6319SBarry Smith Level: beginner 341d91e6319SBarry Smith 3420a835dfdSSatish Balay Any additions/changes here MUST also be made in include/finclude/petscmat.h 343d91e6319SBarry Smith 344d91e6319SBarry Smith .seealso: MatSetOption() 345d91e6319SBarry Smith E*/ 346512a5fc5SBarry Smith typedef enum {MAT_ROW_ORIENTED,MAT_NEW_NONZERO_LOCATIONS, 3474e0d8c25SBarry Smith MAT_SYMMETRIC, 3484e0d8c25SBarry Smith MAT_STRUCTURALLY_SYMMETRIC, 3498047e77bSBarry Smith MAT_NEW_DIAGONALS,MAT_IGNORE_OFF_PROC_ENTRIES, 3504e0d8c25SBarry Smith MAT_NEW_NONZERO_LOCATION_ERR, 3514e0d8c25SBarry Smith MAT_NEW_NONZERO_ALLOCATION_ERR,MAT_USE_HASH_TABLE, 3524e0d8c25SBarry Smith MAT_KEEP_ZEROED_ROWS,MAT_IGNORE_ZERO_ENTRIES, 3534e0d8c25SBarry Smith MAT_USE_INODES, 3544e0d8c25SBarry Smith MAT_HERMITIAN, 3554e0d8c25SBarry Smith MAT_SYMMETRY_ETERNAL, 3564e0d8c25SBarry Smith MAT_USE_COMPRESSEDROW, 3574e0d8c25SBarry Smith MAT_IGNORE_LOWER_TRIANGULAR,MAT_ERROR_LOWER_TRIANGULAR, 3584e0d8c25SBarry Smith MAT_GETROW_UPPERTRIANGULAR} MatOption; 359290bbb0aSBarry Smith extern const char *MatOptions[]; 3604e0d8c25SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetOption(Mat,MatOption,PetscTruth); 361be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetType(Mat,MatType*); 362ba966639SSatish Balay PetscPolymorphicFunction(MatGetType,(Mat mat),(mat,&t),MatType,t) 36384cb2905SBarry Smith 364be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]); 365be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 366be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 367f5edf698SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowUpperTriangular(Mat); 368f5edf698SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRowUpperTriangular(Mat); 369be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 370be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 371be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumnVector(Mat,Vec,PetscInt); 372be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetArray(Mat,PetscScalar *[]); 373ba966639SSatish Balay PetscPolymorphicFunction(MatGetArray,(Mat mat),(mat,&a),PetscScalar*,a) 374be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreArray(Mat,PetscScalar *[]); 375be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBlockSize(Mat,PetscInt *); 376ba966639SSatish Balay PetscPolymorphicFunction(MatGetBlockSize,(Mat mat),(mat,&a),PetscInt,a) 377be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetBlockSize(Mat,PetscInt); 3787b80b807SBarry Smith 3791620fd73SBarry Smith 380be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMult(Mat,Vec,Vec); 381be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultAdd(Mat,Vec,Vec,Vec); 382ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultAdd,(Mat A,Vec x,Vec y),(A,x,y,y)) 383be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTranspose(Mat,Vec,Vec); 384be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsTranspose(Mat,Mat,PetscReal,PetscTruth*); 385ba966639SSatish Balay PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B,PetscReal tol),(A,B,tol,&t),PetscTruth,t) 386e37f64cdSBarry Smith PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B),(A,B,0,&t),PetscTruth,t) 3871cbb95d3SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitianTranspose(Mat,Mat,PetscReal,PetscTruth*); 388be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeAdd(Mat,Vec,Vec,Vec); 389ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultTransposeAdd,(Mat A,Vec x,Vec y),(A,x,y,y)) 390be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultConstrained(Mat,Vec,Vec); 391be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeConstrained(Mat,Vec,Vec); 3922bebee5dSHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatSolve(Mat,Mat,Mat); 393f5edf698SHong Zhang 394d91e6319SBarry Smith /*E 395d91e6319SBarry Smith MatDuplicateOption - Indicates if a duplicated sparse matrix should have 396d91e6319SBarry Smith its numerical values copied over or just its nonzero structure. 397d91e6319SBarry Smith 398d91e6319SBarry Smith Level: beginner 399d91e6319SBarry Smith 400d91e6319SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 401d91e6319SBarry Smith 402d91e6319SBarry Smith .seealso: MatDuplicate() 403d91e6319SBarry Smith E*/ 4042e8a6d31SBarry Smith typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES} MatDuplicateOption; 4052e8a6d31SBarry Smith 406f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert(Mat,MatType,MatReuse,Mat*); 407ba966639SSatish Balay PetscPolymorphicFunction(MatConvert,(Mat A,MatType t),(A,t,MAT_INITIAL_MATRIX,&a),Mat,a) 408be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDuplicate(Mat,MatDuplicateOption,Mat*); 409ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A,MatDuplicateOption o),(A,o,&a),Mat,a) 410ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A),(A,MAT_COPY_VALUES,&a),Mat,a) 41194a9d846SBarry Smith 412d91e6319SBarry Smith /*E 413d91e6319SBarry Smith MatStructure - Indicates if the matrix has the same nonzero structure 414d91e6319SBarry Smith 415d91e6319SBarry Smith Level: beginner 416d91e6319SBarry Smith 417d91e6319SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 418d91e6319SBarry Smith 41994b7f48cSBarry Smith .seealso: MatCopy(), KSPSetOperators(), PCSetOperators() 420d91e6319SBarry Smith E*/ 421c537a176SHong Zhang typedef enum {SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER,SUBSET_NONZERO_PATTERN} MatStructure; 422cb5b572fSBarry Smith 423be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCopy(Mat,Mat,MatStructure); 424be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatView(Mat,PetscViewer); 425be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsSymmetric(Mat,PetscReal,PetscTruth*); 426ba966639SSatish Balay PetscPolymorphicFunction(MatIsSymmetric,(Mat A,PetscReal tol),(A,tol,&t),PetscTruth,t) 427e37f64cdSBarry Smith PetscPolymorphicFunction(MatIsSymmetric,(Mat A),(A,0,&t),PetscTruth,t) 428be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsStructurallySymmetric(Mat,PetscTruth*); 429ba966639SSatish Balay PetscPolymorphicFunction(MatIsStructurallySymmetric,(Mat A),(A,&t),PetscTruth,t) 4301cbb95d3SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitian(Mat,PetscReal,PetscTruth*); 4311cbb95d3SBarry Smith PetscPolymorphicFunction(MatIsHermitian,(Mat A),(A,0,&t),PetscTruth,t) 432be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsSymmetricKnown(Mat,PetscTruth*,PetscTruth*); 433be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitianKnown(Mat,PetscTruth*,PetscTruth*); 43464c62002SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMissingDiagonal(Mat,PetscTruth *,PetscInt *); 435f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLoad(PetscViewer,MatType,Mat*); 436ba966639SSatish Balay PetscPolymorphicFunction(MatLoad,(PetscViewer v,MatType t),(v,t,&a),Mat,a) 4377b80b807SBarry Smith 4388f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *); 4398f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRowIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *); 4408f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumnIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *); 4418f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreColumnIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *); 442d4fbbf0eSBarry Smith 443d91e6319SBarry Smith /*S 444d91e6319SBarry Smith MatInfo - Context of matrix information, used with MatGetInfo() 445d91e6319SBarry Smith 446d91e6319SBarry Smith In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE 447d91e6319SBarry Smith 448d91e6319SBarry Smith Level: intermediate 449d91e6319SBarry Smith 450d91e6319SBarry Smith Concepts: matrix^nonzero information 451d91e6319SBarry Smith 452d9274352SBarry Smith .seealso: MatGetInfo(), MatInfoType 453d91e6319SBarry Smith S*/ 4544e220ebcSLois Curfman McInnes typedef struct { 455b0a32e0cSBarry Smith PetscLogDouble rows_global,columns_global; /* number of global rows and columns */ 456b0a32e0cSBarry Smith PetscLogDouble rows_local,columns_local; /* number of local rows and columns */ 457b0a32e0cSBarry Smith PetscLogDouble block_size; /* block size */ 458b0a32e0cSBarry Smith PetscLogDouble nz_allocated,nz_used,nz_unneeded; /* number of nonzeros */ 459b0a32e0cSBarry Smith PetscLogDouble memory; /* memory allocated */ 460b0a32e0cSBarry Smith PetscLogDouble assemblies; /* number of matrix assemblies called */ 461b0a32e0cSBarry Smith PetscLogDouble mallocs; /* number of mallocs during MatSetValues() */ 462b0a32e0cSBarry Smith PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */ 463b0a32e0cSBarry Smith PetscLogDouble factor_mallocs; /* number of mallocs during factorization */ 4644e220ebcSLois Curfman McInnes } MatInfo; 4654e220ebcSLois Curfman McInnes 466d9274352SBarry Smith /*E 467d9274352SBarry Smith MatInfoType - Indicates if you want information about the local part of the matrix, 468d9274352SBarry Smith the entire parallel matrix or the maximum over all the local parts. 469d9274352SBarry Smith 470d9274352SBarry Smith Level: beginner 471d9274352SBarry Smith 472d9274352SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 473d9274352SBarry Smith 474d9274352SBarry Smith .seealso: MatGetInfo(), MatInfo 475d9274352SBarry Smith E*/ 4767b80b807SBarry Smith typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType; 477be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetInfo(Mat,MatInfoType,MatInfo*); 478be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatValid(Mat,PetscTruth*); 479be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetDiagonal(Mat,Vec); 480985db425SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMax(Mat,Vec,PetscInt[]); 481985db425SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMin(Mat,Vec,PetscInt[]); 482985db425SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMaxAbs(Mat,Vec,PetscInt[]); 48386697f06SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowSum(Mat,Vec); 484fc4dec0aSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatTranspose(Mat,MatReuse,Mat*); 485fc4dec0aSBarry Smith PetscPolymorphicFunction(MatTranspose,(Mat A),(A,MAT_INITIAL_MATRIX,&t),Mat,t) 486be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPermute(Mat,IS,IS,Mat *); 487ba966639SSatish Balay PetscPolymorphicFunction(MatPermute,(Mat A,IS is1,IS is2),(A,is1,is2,&t),Mat,t) 488be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPermuteSparsify(Mat,PetscInt,PetscReal,PetscReal,IS,IS,Mat *); 489be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalScale(Mat,Vec,Vec); 490be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalSet(Mat,Vec,InsertMode); 491be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatEqual(Mat,Mat,PetscTruth*); 492ba966639SSatish Balay PetscPolymorphicFunction(MatEqual,(Mat A,Mat B),(A,B,&t),PetscTruth,t) 493be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultEqual(Mat,Mat,PetscInt,PetscTruth*); 494be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultAddEqual(Mat,Mat,PetscInt,PetscTruth*); 495be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeEqual(Mat,Mat,PetscInt,PetscTruth*); 496be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscTruth*); 4977b80b807SBarry Smith 498be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNorm(Mat,NormType,PetscReal *); 499ba966639SSatish Balay PetscPolymorphicFunction(MatNorm,(Mat A,NormType t),(A,t,&n),PetscReal,n) 500be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroEntries(Mat); 501f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar); 502f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsIS(Mat,IS,PetscScalar); 503f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroColumns(Mat,PetscInt,const PetscInt [],const PetscScalar*); 504f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroColumnsIS(Mat,IS,const PetscScalar*); 5057b80b807SBarry Smith 506be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatUseScaledForm(Mat,PetscTruth); 507be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScaleSystem(Mat,Vec,Vec); 508be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatUnScaleSystem(Mat,Vec,Vec); 5095ef9f2a5SBarry Smith 510be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSize(Mat,PetscInt*,PetscInt*); 511be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalSize(Mat,PetscInt*,PetscInt*); 512be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRange(Mat,PetscInt*,PetscInt*); 5138ee2e534SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRanges(Mat,const PetscInt**); 514f4b8409dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRangeColumn(Mat,PetscInt*,PetscInt*); 5157b80b807SBarry Smith 516be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]); 517be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroyMatrices(PetscInt,Mat *[]); 518be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrix(Mat,IS,IS,PetscInt,MatReuse,Mat *); 519abf0e977SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrixRaw(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscInt,MatReuse,Mat *); 520829201f2SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSeqNonzeroStructure(Mat,Mat *[]); 521829201f2SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroySeqNonzeroStructure(Mat *[]); 5225494a064SHong Zhang 523be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge(MPI_Comm,Mat,PetscInt,MatReuse,Mat*); 524be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPI(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*); 525be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPISymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*); 526be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPINumeric(Mat,Mat); 527be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroy_MPIAIJ_SeqsToMPI(Mat); 528be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalMat(Mat,MatReuse,Mat*); 529be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*); 530be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,PetscInt*,Mat*); 531dd6ea824SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBrowsOfAoCols(Mat,Mat,MatReuse,PetscInt**,MatScalar**,Mat*); 53243eb5e2fSMatthew Knepley #if defined (PETSC_USE_CTABLE) 533cd116e26SSatish Balay #include "petscctable.h" 53443eb5e2fSMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetCommunicationStructs(Mat, Vec *, PetscTable *, VecScatter *); 53543eb5e2fSMatthew Knepley #else 53643eb5e2fSMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetCommunicationStructs(Mat, Vec *, PetscInt *[], VecScatter *); 53743eb5e2fSMatthew Knepley #endif 5388efafbd8SBarry Smith 539be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt); 5407b80b807SBarry Smith 541be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*); 542be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultSymbolic(Mat,Mat,PetscReal,Mat*); 543be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultNumeric(Mat,Mat,Mat); 54422440eb1SKris Buschelman 545be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*); 546be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAPSymbolic(Mat,Mat,PetscReal,Mat*); 547be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAPNumeric(Mat,Mat,Mat); 54822440eb1SKris Buschelman 549be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTranspose(Mat,Mat,MatReuse,PetscReal,Mat*); 550be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTransposeSymbolic(Mat,Mat,PetscReal,Mat*); 551be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTransposeNumeric(Mat,Mat,Mat); 552bc011b1eSHong Zhang 553f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAXPY(Mat,PetscScalar,Mat,MatStructure); 55404aac2b0SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAYPX(Mat,PetscScalar,Mat,MatStructure); 555be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompress(Mat); 5561c741599SHong Zhang 557f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScale(Mat,PetscScalar); 558f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShift(Mat,PetscScalar); 5597b80b807SBarry Smith 560be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping); 561be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetLocalToGlobalMappingBlock(Mat,ISLocalToGlobalMapping); 562f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar); 563f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsLocalIS(Mat,IS,PetscScalar); 564be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 565be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 566052efed2SBarry Smith 567be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStashSetInitialSize(Mat,PetscInt,PetscInt); 568be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*); 56990f02eecSBarry Smith 570be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInterpolate(Mat,Vec,Vec); 5710c451bc4SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInterpolateAdd(Mat,Vec,Vec,Vec); 572ba966639SSatish Balay PetscPolymorphicSubroutine(MatInterpolateAdd,(Mat A,Vec x,Vec y),(A,x,y,y)) 573be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestrict(Mat,Vec,Vec); 574be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetVecs(Mat,Vec*,Vec*); 57569db28dcSHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRedundantMatrix(Mat,PetscInt,MPI_Comm,PetscInt,MatReuse,Mat*); 576bd481603SBarry Smith 577bd481603SBarry Smith /*MC 5781620fd73SBarry Smith MatSetValue - Set a single entry into a matrix. 5791620fd73SBarry Smith 5801620fd73SBarry Smith Not collective 5811620fd73SBarry Smith 5821620fd73SBarry Smith Input Parameters: 5831620fd73SBarry Smith + m - the matrix 5841620fd73SBarry Smith . row - the row location of the entry 5851620fd73SBarry Smith . col - the column location of the entry 5861620fd73SBarry Smith . value - the value to insert 5871620fd73SBarry Smith - mode - either INSERT_VALUES or ADD_VALUES 5881620fd73SBarry Smith 5891620fd73SBarry Smith Notes: 5901620fd73SBarry Smith For efficiency one should use MatSetValues() and set several or many 5911620fd73SBarry Smith values simultaneously if possible. 5921620fd73SBarry Smith 5931620fd73SBarry Smith Level: beginner 5941620fd73SBarry Smith 5951620fd73SBarry Smith .seealso: MatSetValues(), MatSetValueLocal() 5961620fd73SBarry Smith M*/ 5971620fd73SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatSetValue(Mat v,PetscInt i,PetscInt j,PetscScalar va,InsertMode mode) {return MatSetValues(v,1,&i,1,&j,&va,mode);} 5981620fd73SBarry Smith 5991620fd73SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatGetValue(Mat v,PetscInt i,PetscInt j,PetscScalar *va) {return MatGetValues(v,1,&i,1,&j,va);} 6001620fd73SBarry Smith 6011620fd73SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatSetValueLocal(Mat v,PetscInt i,PetscInt j,PetscScalar va,InsertMode mode) {return MatSetValuesLocal(v,1,&i,1,&j,&va,mode);} 6021620fd73SBarry Smith 6031620fd73SBarry Smith /*MC 6040d2a7ff2SSatish Balay MatPreallocateInitialize - Begins the block of code that will count the number of nonzeros per 605bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 606bd481603SBarry Smith 607bd481603SBarry Smith Synopsis: 608c1ac3661SBarry Smith PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz) 609bd481603SBarry Smith 610bd481603SBarry Smith Collective on MPI_Comm 611bd481603SBarry Smith 612bd481603SBarry Smith Input Parameters: 613bd481603SBarry Smith + comm - the communicator that will share the eventually allocated matrix 614859fcb39SBarry Smith . nrows - the number of LOCAL rows in the matrix 615859fcb39SBarry Smith - ncols - the number of LOCAL columns in the matrix 616bd481603SBarry Smith 617bd481603SBarry Smith Output Parameters: 618bd481603SBarry Smith + dnz - the array that will be passed to the matrix preallocation routines 619bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 620bd481603SBarry Smith 621bd481603SBarry Smith 622bd481603SBarry Smith Level: intermediate 623bd481603SBarry Smith 624bd481603SBarry Smith Notes: 625bd481603SBarry Smith See the chapter in the users manual on performance for more details 626bd481603SBarry Smith 6271d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 628bd481603SBarry Smith 629bd481603SBarry Smith Use MatPreallocateInitializeSymmetric() for symmetric matrices (MPISBAIJ matrices) 630bd481603SBarry Smith 6311620fd73SBarry Smith This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize(). 6321620fd73SBarry Smith 633bd481603SBarry Smith Concepts: preallocation^Matrix 634bd481603SBarry Smith 635bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(), 636bd481603SBarry Smith MatPreallocateInitializeSymmetric(), MatPreallocateSymmetricSetLocal() 637bd481603SBarry Smith M*/ 638c4f061fbSSatish Balay #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \ 6397c922b88SBarry Smith { \ 640c1ac3661SBarry Smith PetscErrorCode _4_ierr; PetscInt __tmp = (nrows),__ctmp = (ncols),__rstart,__start,__end; \ 641c1ac3661SBarry Smith _4_ierr = PetscMalloc(2*__tmp*sizeof(PetscInt),&dnz);CHKERRQ(_4_ierr);onz = dnz + __tmp;\ 642c1ac3661SBarry Smith _4_ierr = PetscMemzero(dnz,2*__tmp*sizeof(PetscInt));CHKERRQ(_4_ierr);\ 643c1ac3661SBarry Smith _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ctmp;\ 644c1ac3661SBarry Smith _4_ierr = MPI_Scan(&__tmp,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __tmp; 6457c922b88SBarry Smith 646bd481603SBarry Smith /*MC 6470d2a7ff2SSatish Balay MatPreallocateSymmetricInitialize - Begins the block of code that will count the number of nonzeros per 648bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 649bd481603SBarry Smith 650bd481603SBarry Smith Synopsis: 651c1ac3661SBarry Smith PetscErrorCode MatPreallocateSymmetricInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz) 652bd481603SBarry Smith 653bd481603SBarry Smith Collective on MPI_Comm 654bd481603SBarry Smith 655bd481603SBarry Smith Input Parameters: 656bd481603SBarry Smith + comm - the communicator that will share the eventually allocated matrix 657859fcb39SBarry Smith . nrows - the number of LOCAL rows in the matrix 658859fcb39SBarry Smith - ncols - the number of LOCAL columns in the matrix 659bd481603SBarry Smith 660bd481603SBarry Smith Output Parameters: 661bd481603SBarry Smith + dnz - the array that will be passed to the matrix preallocation routines 662bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 663bd481603SBarry Smith 664bd481603SBarry Smith 665bd481603SBarry Smith Level: intermediate 666bd481603SBarry Smith 667bd481603SBarry Smith Notes: 668bd481603SBarry Smith See the chapter in the users manual on performance for more details 669bd481603SBarry Smith 6701d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 671bd481603SBarry Smith 6721620fd73SBarry Smith This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize(). 6731620fd73SBarry Smith 674bd481603SBarry Smith Concepts: preallocation^Matrix 675bd481603SBarry Smith 676bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(), 677bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal() 678bd481603SBarry Smith M*/ 679222b16d4SBarry Smith #define MatPreallocateSymmetricInitialize(comm,nrows,ncols,dnz,onz) 0; \ 680222b16d4SBarry Smith { \ 681c1ac3661SBarry Smith PetscErrorCode _4_ierr; PetscInt __tmp = (nrows),__ctmp = (ncols),__rstart,__end; \ 682c1ac3661SBarry Smith _4_ierr = PetscMalloc(2*__tmp*sizeof(PetscInt),&dnz);CHKERRQ(_4_ierr);onz = dnz + __tmp;\ 683c1ac3661SBarry Smith _4_ierr = PetscMemzero(dnz,2*__tmp*sizeof(PetscInt));CHKERRQ(_4_ierr);\ 684c1ac3661SBarry Smith _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr);\ 685c1ac3661SBarry Smith _4_ierr = MPI_Scan(&__tmp,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __tmp; 686222b16d4SBarry Smith 687bd481603SBarry Smith /*MC 688bd481603SBarry Smith MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be 689bd481603SBarry Smith inserted using a local number of the rows and columns 690bd481603SBarry Smith 691bd481603SBarry Smith Synopsis: 692c1ac3661SBarry Smith PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 693bd481603SBarry Smith 694bd481603SBarry Smith Not Collective 695bd481603SBarry Smith 696bd481603SBarry Smith Input Parameters: 697bd481603SBarry Smith + map - the mapping between local numbering and global numbering 698bd481603SBarry Smith . nrows - the number of rows indicated 6991d73ed98SBarry Smith . rows - the indices of the rows 700bd481603SBarry Smith . ncols - the number of columns in the matrix 701bd481603SBarry Smith . cols - the columns indicated 702bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 703bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 704bd481603SBarry Smith 705bd481603SBarry Smith 706bd481603SBarry Smith Level: intermediate 707bd481603SBarry Smith 708bd481603SBarry Smith Notes: 709bd481603SBarry Smith See the chapter in the users manual on performance for more details 710bd481603SBarry Smith 7111d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 712bd481603SBarry Smith 713bd481603SBarry Smith Concepts: preallocation^Matrix 714bd481603SBarry Smith 715bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(), 716bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal() 717bd481603SBarry Smith M*/ 718c4f061fbSSatish Balay #define MatPreallocateSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\ 719c4f061fbSSatish Balay {\ 720c1ac3661SBarry Smith PetscInt __l;\ 721ef66eb69SBarry Smith _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\ 722ef66eb69SBarry Smith _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\ 723c4f061fbSSatish Balay for (__l=0;__l<nrows;__l++) {\ 724ef66eb69SBarry Smith _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 725c4f061fbSSatish Balay }\ 726c4f061fbSSatish Balay } 727c4f061fbSSatish Balay 728bd481603SBarry Smith /*MC 729bd481603SBarry Smith MatPreallocateSymmetricSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be 730bd481603SBarry Smith inserted using a local number of the rows and columns 731bd481603SBarry Smith 732bd481603SBarry Smith Synopsis: 733c1ac3661SBarry Smith PetscErrorCode MatPreallocateSymmetricSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 734bd481603SBarry Smith 735bd481603SBarry Smith Not Collective 736bd481603SBarry Smith 737bd481603SBarry Smith Input Parameters: 738bd481603SBarry Smith + map - the mapping between local numbering and global numbering 739bd481603SBarry Smith . nrows - the number of rows indicated 7401d73ed98SBarry Smith . rows - the indices of the rows 741bd481603SBarry Smith . ncols - the number of columns in the matrix 742bd481603SBarry Smith . cols - the columns indicated 743bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 744bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 745bd481603SBarry Smith 746bd481603SBarry Smith 747bd481603SBarry Smith Level: intermediate 748bd481603SBarry Smith 749bd481603SBarry Smith Notes: 750bd481603SBarry Smith See the chapter in the users manual on performance for more details 751bd481603SBarry Smith 752bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 753bd481603SBarry Smith 754bd481603SBarry Smith Concepts: preallocation^Matrix 755bd481603SBarry Smith 756bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(), 757bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal() 758bd481603SBarry Smith M*/ 759d3d32019SBarry Smith #define MatPreallocateSymmetricSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\ 760d3d32019SBarry Smith {\ 761c1ac3661SBarry Smith PetscInt __l;\ 762d3d32019SBarry Smith _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\ 763d3d32019SBarry Smith _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\ 764d3d32019SBarry Smith for (__l=0;__l<nrows;__l++) {\ 765d3d32019SBarry Smith _4_ierr = MatPreallocateSymmetricSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 766d3d32019SBarry Smith }\ 767d3d32019SBarry Smith } 768d3d32019SBarry Smith 769bd481603SBarry Smith /*MC 770bd481603SBarry Smith MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be 771bd481603SBarry Smith inserted using a local number of the rows and columns 772bd481603SBarry Smith 773bd481603SBarry Smith Synopsis: 774c1ac3661SBarry Smith PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 775bd481603SBarry Smith 776bd481603SBarry Smith Not Collective 777bd481603SBarry Smith 778bd481603SBarry Smith Input Parameters: 77964075487SBarry Smith + row - the row 780bd481603SBarry Smith . ncols - the number of columns in the matrix 781a50f8bf6SHong Zhang - cols - the columns indicated 782a50f8bf6SHong Zhang 783a50f8bf6SHong Zhang Output Parameters: 784a50f8bf6SHong Zhang + dnz - the array that will be passed to the matrix preallocation routines 785bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 786bd481603SBarry Smith 787bd481603SBarry Smith 788bd481603SBarry Smith Level: intermediate 789bd481603SBarry Smith 790bd481603SBarry Smith Notes: 791bd481603SBarry Smith See the chapter in the users manual on performance for more details 792bd481603SBarry Smith 793bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 794bd481603SBarry Smith 7951620fd73SBarry Smith This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize(). 7961620fd73SBarry Smith 797bd481603SBarry Smith Concepts: preallocation^Matrix 798bd481603SBarry Smith 799bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(), 800bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal() 801bd481603SBarry Smith M*/ 802c4f061fbSSatish Balay #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\ 803c1ac3661SBarry Smith { PetscInt __i; \ 8048ee2e534SBarry Smith if (row < __rstart) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Trying to set preallocation for row %D less than first local row %D",row,__rstart);\ 8058ee2e534SBarry 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);\ 8067c922b88SBarry Smith for (__i=0; __i<nc; __i++) {\ 80764075487SBarry Smith if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \ 8087cc688d7SBarry Smith else dnz[row - __rstart]++;\ 8097c922b88SBarry Smith }\ 8107c922b88SBarry Smith } 8117c922b88SBarry Smith 812bd481603SBarry Smith /*MC 813bd481603SBarry Smith MatPreallocateSymmetricSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be 814bd481603SBarry Smith inserted using a local number of the rows and columns 815bd481603SBarry Smith 816bd481603SBarry Smith Synopsis: 817c1ac3661SBarry Smith PetscErrorCode MatPreallocateSymmetricSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 818bd481603SBarry Smith 819bd481603SBarry Smith Not Collective 820bd481603SBarry Smith 821bd481603SBarry Smith Input Parameters: 822bd481603SBarry Smith + nrows - the number of rows indicated 8231d73ed98SBarry Smith . rows - the indices of the rows 824bd481603SBarry Smith . ncols - the number of columns in the matrix 825bd481603SBarry Smith . cols - the columns indicated 826bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 827bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 828bd481603SBarry Smith 829bd481603SBarry Smith 830bd481603SBarry Smith Level: intermediate 831bd481603SBarry Smith 832bd481603SBarry Smith Notes: 833bd481603SBarry Smith See the chapter in the users manual on performance for more details 834bd481603SBarry Smith 835bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 836bd481603SBarry Smith 8371620fd73SBarry Smith This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize(). 8381620fd73SBarry Smith 839bd481603SBarry Smith Concepts: preallocation^Matrix 840bd481603SBarry Smith 841bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(), 842bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal() 843bd481603SBarry Smith M*/ 844d3d32019SBarry Smith #define MatPreallocateSymmetricSet(row,nc,cols,dnz,onz) 0;\ 845c1ac3661SBarry Smith { PetscInt __i; \ 846d3d32019SBarry Smith for (__i=0; __i<nc; __i++) {\ 847d3d32019SBarry Smith if (cols[__i] >= __end) onz[row - __rstart]++; \ 848d3d32019SBarry Smith else if (cols[__i] >= row) dnz[row - __rstart]++;\ 849d3d32019SBarry Smith }\ 850d3d32019SBarry Smith } 851d3d32019SBarry Smith 852bd481603SBarry Smith /*MC 8530d2a7ff2SSatish Balay MatPreallocateFinalize - Ends the block of code that will count the number of nonzeros per 854bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 855bd481603SBarry Smith 856bd481603SBarry Smith Synopsis: 857c1ac3661SBarry Smith PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz) 858bd481603SBarry Smith 859bd481603SBarry Smith Collective on MPI_Comm 860bd481603SBarry Smith 861bd481603SBarry Smith Input Parameters: 862bd481603SBarry Smith + dnz - the array that will be passed to the matrix preallocation routines 863bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 864bd481603SBarry Smith 865bd481603SBarry Smith 866bd481603SBarry Smith Level: intermediate 867bd481603SBarry Smith 868bd481603SBarry Smith Notes: 869bd481603SBarry Smith See the chapter in the users manual on performance for more details 870bd481603SBarry Smith 871bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 872bd481603SBarry Smith 8731620fd73SBarry Smith This is a MACRO not a function because it closes the { started in MatPreallocateInitialize(). 8741620fd73SBarry Smith 875bd481603SBarry Smith Concepts: preallocation^Matrix 876bd481603SBarry Smith 877bd481603SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(), 878bd481603SBarry Smith MatPreallocateSymmetricInitialize(), MatPreallocateSymmetricSetLocal() 879bd481603SBarry Smith M*/ 880ef66eb69SBarry Smith #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree(dnz);CHKERRQ(_4_ierr);} 8817c922b88SBarry Smith 882bd481603SBarry Smith 883bd481603SBarry Smith 8847b80b807SBarry Smith /* Routines unique to particular data structures */ 885be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellGetContext(Mat,void **); 886ba966639SSatish Balay PetscPolymorphicFunction(MatShellGetContext,(Mat A),(A,&t),void*,t) 887698d4c6aSKris Buschelman 888be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInodeAdjustForInodes(Mat,IS*,IS*); 889be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *); 890698d4c6aSKris Buschelman 891be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetColumnIndices(Mat,PetscInt[]); 892be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetColumnIndices(Mat,PetscInt[]); 893be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 894c75a6043SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 895c75a6043SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 8967b80b807SBarry Smith 897a96a251dSBarry Smith #define MAT_SKIP_ALLOCATION -4 898a96a251dSBarry Smith 899be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]); 900ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz)) 901be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]); 902ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqSBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz)) 903be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]); 904ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqAIJSetPreallocation,(Mat A,const PetscInt nnz[]),(A,0,nnz)) 905be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDensePreallocation(Mat,PetscScalar[]); 906273d9f13SBarry Smith 907be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 908ba966639SSatish Balay PetscPolymorphicSubroutine(MatMPIBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[],const PetscInt onz[]),(A,bs,0,nnz,0,onz)) 909be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 910be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 91153803f26SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []); 912*725b52f3SLisandro Dalcin EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 913be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]); 914be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 915be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]); 916be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIDenseSetPreallocation(Mat,PetscScalar[]); 917284134d9SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation(Mat,PetscScalar[]); 918be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIRowbsSetPreallocation(Mat,PetscInt,const PetscInt[]); 919be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,PetscInt*[]); 920be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,PetscInt*[]); 921be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAdicSetLocalFunction(Mat,void (*)(void)); 922273d9f13SBarry Smith 923be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetLDA(Mat,PetscInt); 924ab92ecdeSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDenseGetLocalMatrix(Mat,Mat*); 9251b807ce4Svictorle 926be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStoreValues(Mat); 927be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRetrieveValues(Mat); 9282e8a6d31SBarry Smith 929be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDAADSetCtx(Mat,void*); 9303a7fca6bSBarry Smith 9317b80b807SBarry Smith /* 9327b80b807SBarry Smith These routines are not usually accessed directly, rather solving is 93394b7f48cSBarry Smith done through the KSP and PC interfaces. 9347b80b807SBarry Smith */ 9357b80b807SBarry Smith 936d9274352SBarry Smith /*E 937d9274352SBarry Smith MatOrderingType - String with the name of a PETSc matrix ordering or the creation function 938d9274352SBarry Smith with an optional dynamic library name, for example 939d9274352SBarry Smith http://www.mcs.anl.gov/petsc/lib.a:orderingcreate() 940d9274352SBarry Smith 941d9274352SBarry Smith Level: beginner 942d9274352SBarry Smith 943e5a9bf91SBarry Smith Cannot use const because the PC objects manipulate the string 944e5a9bf91SBarry Smith 945d9274352SBarry Smith .seealso: MatGetOrdering() 946d9274352SBarry Smith E*/ 94749773a63SBarry Smith #define MatOrderingType char* 948b12f92e5SBarry Smith #define MATORDERING_NATURAL "natural" 949b12f92e5SBarry Smith #define MATORDERING_ND "nd" 950b12f92e5SBarry Smith #define MATORDERING_1WD "1wd" 951b12f92e5SBarry Smith #define MATORDERING_RCM "rcm" 952b12f92e5SBarry Smith #define MATORDERING_QMD "qmd" 953b12f92e5SBarry Smith #define MATORDERING_ROWLENGTH "rowlength" 95462152c8bSBarry Smith #define MATORDERING_DSC_ND "dsc_nd" 95562152c8bSBarry Smith #define MATORDERING_DSC_MMD "dsc_mmd" 95662152c8bSBarry Smith #define MATORDERING_DSC_MDF "dsc_mdf" 957c06d978dSMatthew Knepley #define MATORDERING_CONSTRAINED "constrained" 958c06d978dSMatthew Knepley #define MATORDERING_IDENTITY "identity" 959c06d978dSMatthew Knepley #define MATORDERING_REVERSE "reverse" 960b12f92e5SBarry Smith 961be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOrdering(Mat,const MatOrderingType,IS*,IS*); 962be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOrderingList(PetscFList *list); 963be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,const MatOrderingType,IS*,IS*)); 96430de9b25SBarry Smith 96530de9b25SBarry Smith /*MC 96630de9b25SBarry Smith MatOrderingRegisterDynamic - Adds a new sparse matrix ordering to the 96730de9b25SBarry Smith matrix package. 96830de9b25SBarry Smith 96930de9b25SBarry Smith Synopsis: 970c1ac3661SBarry Smith PetscErrorCode MatOrderingRegisterDynamic(char *name_ordering,char *path,char *name_create,PetscErrorCode (*routine_create)(MatOrdering)) 97130de9b25SBarry Smith 97230de9b25SBarry Smith Not Collective 97330de9b25SBarry Smith 97430de9b25SBarry Smith Input Parameters: 97530de9b25SBarry Smith + sname - name of ordering (for example MATORDERING_ND) 97630de9b25SBarry Smith . path - location of library where creation routine is 97730de9b25SBarry Smith . name - name of function that creates the ordering type,a string 97830de9b25SBarry Smith - function - function pointer that creates the ordering 97930de9b25SBarry Smith 98030de9b25SBarry Smith Level: developer 98130de9b25SBarry Smith 98230de9b25SBarry Smith If dynamic libraries are used, then the fourth input argument (function) 98330de9b25SBarry Smith is ignored. 98430de9b25SBarry Smith 98530de9b25SBarry Smith Sample usage: 98630de9b25SBarry Smith .vb 98730de9b25SBarry Smith MatOrderingRegisterDynamic("my_order",/home/username/my_lib/lib/libO/solaris/mylib.a, 98830de9b25SBarry Smith "MyOrder",MyOrder); 98930de9b25SBarry Smith .ve 99030de9b25SBarry Smith 99130de9b25SBarry Smith Then, your partitioner can be chosen with the procedural interface via 99230de9b25SBarry Smith $ MatOrderingSetType(part,"my_order) 99330de9b25SBarry Smith or at runtime via the option 9942401956bSBarry Smith $ -pc_factor_mat_ordering_type my_order 99530de9b25SBarry Smith 996ab901514SBarry Smith ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values. 99730de9b25SBarry Smith 99830de9b25SBarry Smith .keywords: matrix, ordering, register 99930de9b25SBarry Smith 100030de9b25SBarry Smith .seealso: MatOrderingRegisterDestroy(), MatOrderingRegisterAll() 100130de9b25SBarry Smith M*/ 1002aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 1003f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,0) 1004b12f92e5SBarry Smith #else 1005f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,d) 1006b12f92e5SBarry Smith #endif 100730de9b25SBarry Smith 1008be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegisterDestroy(void); 1009be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegisterAll(const char[]); 10102bad1931SBarry Smith extern PetscTruth MatOrderingRegisterAllCalled; 1011b0a32e0cSBarry Smith extern PetscFList MatOrderingList; 1012d4fbbf0eSBarry Smith 1013be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS); 1014a2ce50c7SBarry Smith 1015d91e6319SBarry Smith /*S 10162401956bSBarry Smith MatFactorInfo - Data passed into the matrix factorization routines 1017b00f7748SHong Zhang 101861cad860SBarry Smith In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE, that is use 101961cad860SBarry Smith $ MatFactorInfo info(MAT_FACTORINFO_SIZE) 1020b00f7748SHong Zhang 102115e8a5b3SHong Zhang Notes: These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC. 1022b00f7748SHong Zhang 102361cad860SBarry Smith You can use MatFactorInfoInitialize() to set default values. 102461cad860SBarry Smith 1025b00f7748SHong Zhang Level: developer 1026b00f7748SHong Zhang 1027d7d82daaSBarry Smith .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(), 1028d7d82daaSBarry Smith MatFactorInfoInitialize() 1029b00f7748SHong Zhang 1030b00f7748SHong Zhang S*/ 1031b00f7748SHong Zhang typedef struct { 10320a29876aSHong Zhang PetscReal shiftnz; /* scaling of identity added to matrix to prevent zero pivots */ 1033fbf22428SSatish Balay PetscReal shiftpd; /* if true, shift until positive pivots */ 10342cea7109SBarry Smith PetscReal shift_fraction; /* record shift fraction taken */ 103515e8a5b3SHong Zhang PetscReal diagonal_fill; /* force diagonal to fill in if initially not filled */ 103615e8a5b3SHong Zhang PetscReal dt; /* drop tolerance */ 1037b00f7748SHong Zhang PetscReal dtcol; /* tolerance for pivoting */ 103815e8a5b3SHong Zhang PetscReal dtcount; /* maximum nonzeros to be allowed per row */ 103967e28bfeSBarry Smith PetscReal fill; /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */ 1040348344bbSBarry Smith PetscReal levels; /* ICC/ILU(levels) */ 1041bcd9e38bSBarry Smith PetscReal pivotinblocks; /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0 1042bcd9e38bSBarry Smith factorization may be faster if do not pivot */ 104362bba022SBarry Smith PetscReal shiftinblocks; /* if block in block factorization has zero pivot then shift diagonal until non-singular */ 104415e8a5b3SHong Zhang PetscReal zeropivot; /* pivot is called zero if less than this */ 104515e8a5b3SHong Zhang } MatFactorInfo; 1046ffa6d0a5SLois Curfman McInnes 1047be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFactorInfoInitialize(MatFactorInfo*); 1048be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactor(Mat,IS,MatFactorInfo*); 1049be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactorSymbolic(Mat,IS,MatFactorInfo*,Mat*); 1050be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactorNumeric(Mat,MatFactorInfo*,Mat*); 1051be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactor(Mat,IS,IS,MatFactorInfo*); 1052be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUFactor(Mat,IS,IS,MatFactorInfo*); 1053be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactorSymbolic(Mat,IS,IS,MatFactorInfo*,Mat*); 1054be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUFactorSymbolic(Mat,IS,IS,MatFactorInfo*,Mat*); 1055be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatICCFactorSymbolic(Mat,IS,MatFactorInfo*,Mat*); 1056be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatICCFactor(Mat,IS,MatFactorInfo*); 1057be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactorNumeric(Mat,MatFactorInfo*,Mat*); 1058be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUDTFactor(Mat,IS,IS,MatFactorInfo*,Mat *); 1059be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*); 1060be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolve(Mat,Vec,Vec); 1061be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatForwardSolve(Mat,Vec,Vec); 1062be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatBackwardSolve(Mat,Vec,Vec); 1063be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveAdd(Mat,Vec,Vec,Vec); 1064be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveTranspose(Mat,Vec,Vec); 1065be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveTransposeAdd(Mat,Vec,Vec,Vec); 1066be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolves(Mat,Vecs,Vecs); 10678ed539a5SBarry Smith 1068be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUnfactored(Mat); 1069bb5a7306SBarry Smith 1070d91e6319SBarry Smith /*E 1071d91e6319SBarry Smith MatSORType - What type of (S)SOR to perform 1072bb1eb677SSatish Balay 1073d91e6319SBarry Smith Level: beginner 1074d91e6319SBarry Smith 1075d9274352SBarry Smith May be bitwise ORd together 1076d9274352SBarry Smith 1077d91e6319SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 1078d91e6319SBarry Smith 10794e7234bfSBarry Smith MatSORType may be bitwise ORd together, so do not change the numbers 10804e7234bfSBarry Smith 1081a1d92eedSBarry Smith .seealso: MatRelax(), MatPBRelax() 1082d91e6319SBarry Smith E*/ 1083ee50ffe9SBarry Smith typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3, 1084ee50ffe9SBarry Smith SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8, 1085ee50ffe9SBarry Smith SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16, 108684cb2905SBarry Smith SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType; 1087be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRelax(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec); 1088be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPBRelax(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec); 10898ed539a5SBarry Smith 1090d4fbbf0eSBarry Smith /* 1091639f9d9dSBarry Smith These routines are for efficiently computing Jacobians via finite differences. 1092639f9d9dSBarry Smith */ 1093b12f92e5SBarry Smith 1094d9274352SBarry Smith /*E 1095d9274352SBarry Smith MatColoringType - String with the name of a PETSc matrix coloring or the creation function 1096d9274352SBarry Smith with an optional dynamic library name, for example 1097d9274352SBarry Smith http://www.mcs.anl.gov/petsc/lib.a:coloringcreate() 1098d9274352SBarry Smith 1099d9274352SBarry Smith Level: beginner 1100d9274352SBarry Smith 1101d9274352SBarry Smith .seealso: MatGetColoring() 1102d9274352SBarry Smith E*/ 1103e5a9bf91SBarry Smith #define MatColoringType const char* 1104b12f92e5SBarry Smith #define MATCOLORING_NATURAL "natural" 1105b12f92e5SBarry Smith #define MATCOLORING_SL "sl" 1106b12f92e5SBarry Smith #define MATCOLORING_LF "lf" 1107b12f92e5SBarry Smith #define MATCOLORING_ID "id" 1108b12f92e5SBarry Smith 11092e90c967SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColoring(Mat,MatColoringType,ISColoring*); 11102e90c967SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,MatColoringType,ISColoring *)); 111130de9b25SBarry Smith 111230de9b25SBarry Smith /*MC 111330de9b25SBarry Smith MatColoringRegisterDynamic - Adds a new sparse matrix coloring to the 111430de9b25SBarry Smith matrix package. 111530de9b25SBarry Smith 111630de9b25SBarry Smith Synopsis: 1117c1ac3661SBarry Smith PetscErrorCode MatColoringRegisterDynamic(char *name_coloring,char *path,char *name_create,PetscErrorCode (*routine_create)(MatColoring)) 111830de9b25SBarry Smith 111930de9b25SBarry Smith Not Collective 112030de9b25SBarry Smith 112130de9b25SBarry Smith Input Parameters: 112230de9b25SBarry Smith + sname - name of Coloring (for example MATCOLORING_SL) 112330de9b25SBarry Smith . path - location of library where creation routine is 112430de9b25SBarry Smith . name - name of function that creates the Coloring type, a string 112530de9b25SBarry Smith - function - function pointer that creates the coloring 112630de9b25SBarry Smith 112730de9b25SBarry Smith Level: developer 112830de9b25SBarry Smith 112930de9b25SBarry Smith If dynamic libraries are used, then the fourth input argument (function) 113030de9b25SBarry Smith is ignored. 113130de9b25SBarry Smith 113230de9b25SBarry Smith Sample usage: 113330de9b25SBarry Smith .vb 113430de9b25SBarry Smith MatColoringRegisterDynamic("my_color",/home/username/my_lib/lib/libO/solaris/mylib.a, 113530de9b25SBarry Smith "MyColor",MyColor); 113630de9b25SBarry Smith .ve 113730de9b25SBarry Smith 113830de9b25SBarry Smith Then, your partitioner can be chosen with the procedural interface via 113930de9b25SBarry Smith $ MatColoringSetType(part,"my_color") 114030de9b25SBarry Smith or at runtime via the option 114130de9b25SBarry Smith $ -mat_coloring_type my_color 114230de9b25SBarry Smith 1143ab901514SBarry Smith $PETSC_ARCH occuring in pathname will be replaced with appropriate values. 114430de9b25SBarry Smith 114530de9b25SBarry Smith .keywords: matrix, Coloring, register 114630de9b25SBarry Smith 114730de9b25SBarry Smith .seealso: MatColoringRegisterDestroy(), MatColoringRegisterAll() 114830de9b25SBarry Smith M*/ 1149aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 1150f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,0) 1151b12f92e5SBarry Smith #else 1152f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,d) 1153b12f92e5SBarry Smith #endif 115430de9b25SBarry Smith 11552bad1931SBarry Smith extern PetscTruth MatColoringRegisterAllCalled; 1156f1144a30SSatish Balay 1157be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegisterAll(const char[]); 1158be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegisterDestroy(void); 1159be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*); 1160639f9d9dSBarry Smith 1161d9274352SBarry Smith /*S 1162d9274352SBarry Smith MatFDColoring - Object for computing a sparse Jacobian via finite differences 1163d9274352SBarry Smith and coloring 1164639f9d9dSBarry Smith 1165d9274352SBarry Smith Level: beginner 1166d9274352SBarry Smith 1167d9274352SBarry Smith Concepts: coloring, sparse Jacobian, finite differences 1168d9274352SBarry Smith 1169d9274352SBarry Smith .seealso: MatFDColoringCreate() 1170d9274352SBarry Smith S*/ 1171e2a1c21fSSatish Balay typedef struct _p_MatFDColoring* MatFDColoring; 1172639f9d9dSBarry Smith 1173be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringCreate(Mat,ISColoring,MatFDColoring *); 1174be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringDestroy(MatFDColoring); 1175be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringView(MatFDColoring,PetscViewer); 1176be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*); 11774a9d489dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**); 1178be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal); 1179be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFrequency(MatFDColoring,PetscInt); 1180be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetFrequency(MatFDColoring,PetscInt*); 1181be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFromOptions(MatFDColoring); 1182be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApply(Mat,MatFDColoring,Vec,MatStructure*,void *); 1183be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApplyTS(Mat,MatFDColoring,PetscReal,Vec,MatStructure*,void *); 1184be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetRecompute(MatFDColoring); 1185be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetF(MatFDColoring,Vec); 1186be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,PetscInt*[]); 1187639f9d9dSBarry Smith /* 11880752156aSBarry Smith These routines are for partitioning matrices: currently used only 11893eda8832SBarry Smith for adjacency matrix, MatCreateMPIAdj(). 11900752156aSBarry Smith */ 1191ca161407SBarry Smith 1192d9274352SBarry Smith /*S 1193d9274352SBarry Smith MatPartitioning - Object for managing the partitioning of a matrix or graph 1194d9274352SBarry Smith 1195d9274352SBarry Smith Level: beginner 1196d9274352SBarry Smith 1197d9274352SBarry Smith Concepts: partitioning 1198d9274352SBarry Smith 1199743a1026Svictorle .seealso: MatPartitioningCreate(), MatPartitioningType 1200d9274352SBarry Smith S*/ 120191e9ee9fSBarry Smith typedef struct _p_MatPartitioning* MatPartitioning; 1202d9274352SBarry Smith 1203d9274352SBarry Smith /*E 12045bc6e7ccSvictorle MatPartitioningType - String with the name of a PETSc matrix partitioning or the creation function 1205d9274352SBarry Smith with an optional dynamic library name, for example 1206d9274352SBarry Smith http://www.mcs.anl.gov/petsc/lib.a:partitioningcreate() 1207d9274352SBarry Smith 1208d9274352SBarry Smith Level: beginner 1209d9274352SBarry Smith 1210b547a13bSvictorle .seealso: MatPartitioningCreate(), MatPartitioning 1211d9274352SBarry Smith E*/ 1212e5a9bf91SBarry Smith #define MatPartitioningType const char* 12138ba1e511SMatthew Knepley #define MAT_PARTITIONING_CURRENT "current" 1214c2ec2a73SHong Zhang #define MAT_PARTITIONING_SQUARE "square" 12158ba1e511SMatthew Knepley #define MAT_PARTITIONING_PARMETIS "parmetis" 121671306c60Sjroman #define MAT_PARTITIONING_CHACO "chaco" 121771306c60Sjroman #define MAT_PARTITIONING_JOSTLE "jostle" 121871306c60Sjroman #define MAT_PARTITIONING_PARTY "party" 121971306c60Sjroman #define MAT_PARTITIONING_SCOTCH "scotch" 122071306c60Sjroman 1221ca161407SBarry Smith 1222be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningCreate(MPI_Comm,MatPartitioning*); 12232e90c967SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetType(MatPartitioning,MatPartitioningType); 1224be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetNParts(MatPartitioning,PetscInt); 1225be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetAdjacency(MatPartitioning,Mat); 1226be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]); 1227be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []); 1228be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningApply(MatPartitioning,IS*); 1229be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningDestroy(MatPartitioning); 12302aabb6bbSBarry Smith 1231be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatPartitioning)); 123230de9b25SBarry Smith 123330de9b25SBarry Smith /*MC 123430de9b25SBarry Smith MatPartitioningRegisterDynamic - Adds a new sparse matrix partitioning to the 123530de9b25SBarry Smith matrix package. 123630de9b25SBarry Smith 123730de9b25SBarry Smith Synopsis: 1238c1ac3661SBarry Smith PetscErrorCode MatPartitioningRegisterDynamic(char *name_partitioning,char *path,char *name_create,PetscErrorCode (*routine_create)(MatPartitioning)) 123930de9b25SBarry Smith 124030de9b25SBarry Smith Not Collective 124130de9b25SBarry Smith 124230de9b25SBarry Smith Input Parameters: 124330de9b25SBarry Smith + sname - name of partitioning (for example MAT_PARTITIONING_CURRENT) or parmetis 124430de9b25SBarry Smith . path - location of library where creation routine is 124530de9b25SBarry Smith . name - name of function that creates the partitioning type, a string 124630de9b25SBarry Smith - function - function pointer that creates the partitioning type 124730de9b25SBarry Smith 124830de9b25SBarry Smith Level: developer 124930de9b25SBarry Smith 125030de9b25SBarry Smith If dynamic libraries are used, then the fourth input argument (function) 125130de9b25SBarry Smith is ignored. 125230de9b25SBarry Smith 125330de9b25SBarry Smith Sample usage: 125430de9b25SBarry Smith .vb 125530de9b25SBarry Smith MatPartitioningRegisterDynamic("my_part",/home/username/my_lib/lib/libO/solaris/mylib.a, 125630de9b25SBarry Smith "MyPartCreate",MyPartCreate); 125730de9b25SBarry Smith .ve 125830de9b25SBarry Smith 125930de9b25SBarry Smith Then, your partitioner can be chosen with the procedural interface via 126030de9b25SBarry Smith $ MatPartitioningSetType(part,"my_part") 126130de9b25SBarry Smith or at runtime via the option 126230de9b25SBarry Smith $ -mat_partitioning_type my_part 126330de9b25SBarry Smith 1264ab901514SBarry Smith $PETSC_ARCH occuring in pathname will be replaced with appropriate values. 126530de9b25SBarry Smith 126630de9b25SBarry Smith .keywords: matrix, partitioning, register 126730de9b25SBarry Smith 126830de9b25SBarry Smith .seealso: MatPartitioningRegisterDestroy(), MatPartitioningRegisterAll() 126930de9b25SBarry Smith M*/ 1270aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 1271f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,0) 12722aabb6bbSBarry Smith #else 1273f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,d) 12742aabb6bbSBarry Smith #endif 12752aabb6bbSBarry Smith 12762bad1931SBarry Smith extern PetscTruth MatPartitioningRegisterAllCalled; 1277f1144a30SSatish Balay 1278be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegisterAll(const char[]); 1279be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegisterDestroy(void); 12802bad1931SBarry Smith 1281be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningView(MatPartitioning,PetscViewer); 1282be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetFromOptions(MatPartitioning); 1283be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningGetType(MatPartitioning,MatPartitioningType*); 1284ca161407SBarry Smith 1285be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningParmetisSetCoarseSequential(MatPartitioning); 12860752156aSBarry Smith 1287be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningJostleSetCoarseLevel(MatPartitioning,PetscReal); 1288be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningJostleSetCoarseSequential(MatPartitioning); 128971306c60Sjroman 1290954b3ec6SBarry Smith typedef enum {MP_CHACO_MULTILEVEL_KL,MP_CHACO_SPECTRAL,MP_CHACO_LINEAR,MP_CHACO_RANDOM, MP_CHACO_SCATTERED} MPChacoGlobalType; 1291be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetGlobal(MatPartitioning, MPChacoGlobalType); 129271306c60Sjroman typedef enum { MP_CHACO_KERNIGHAN_LIN, MP_CHACO_NONE } MPChacoLocalType; 1293be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetLocal(MatPartitioning, MPChacoLocalType); 1294be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal); 129571306c60Sjroman typedef enum { MP_CHACO_LANCZOS, MP_CHACO_RQI_SYMMLQ } MPChacoEigenType; 1296be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType); 1297be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenTol(MatPartitioning, PetscReal); 1298be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenNumber(MatPartitioning, PetscInt); 129971306c60Sjroman 130071306c60Sjroman #define MP_PARTY_OPT "opt" 130171306c60Sjroman #define MP_PARTY_LIN "lin" 130271306c60Sjroman #define MP_PARTY_SCA "sca" 130371306c60Sjroman #define MP_PARTY_RAN "ran" 130471306c60Sjroman #define MP_PARTY_GBF "gbf" 130571306c60Sjroman #define MP_PARTY_GCF "gcf" 130671306c60Sjroman #define MP_PARTY_BUB "bub" 130771306c60Sjroman #define MP_PARTY_DEF "def" 1308be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetGlobal(MatPartitioning, const char*); 130971306c60Sjroman #define MP_PARTY_HELPFUL_SETS "hs" 131071306c60Sjroman #define MP_PARTY_KERNIGHAN_LIN "kl" 131171306c60Sjroman #define MP_PARTY_NONE "no" 1312be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetLocal(MatPartitioning, const char*); 1313be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal); 1314be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetBipart(MatPartitioning,PetscTruth); 1315be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscTruth); 131671306c60Sjroman 131771306c60Sjroman typedef enum { MP_SCOTCH_GREEDY, MP_SCOTCH_GPS, MP_SCOTCH_GR_GPS } MPScotchGlobalType; 1318be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetArch(MatPartitioning,const char*); 1319be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetMultilevel(MatPartitioning); 1320be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetGlobal(MatPartitioning,MPScotchGlobalType); 1321be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetCoarseLevel(MatPartitioning,PetscReal); 1322be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetHostList(MatPartitioning,const char*); 132371306c60Sjroman typedef enum { MP_SCOTCH_KERNIGHAN_LIN, MP_SCOTCH_NONE } MPScotchLocalType; 1324be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetLocal(MatPartitioning,MPScotchLocalType); 1325be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetMapping(MatPartitioning); 1326be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetStrategy(MatPartitioning,char*); 132771306c60Sjroman 13280752156aSBarry Smith /* 13290a835dfdSSatish Balay If you add entries here you must also add them to finclude/petscmat.h 1330d4fbbf0eSBarry Smith */ 13311c1c02c0SLois Curfman McInnes typedef enum { MATOP_SET_VALUES=0, 13321c1c02c0SLois Curfman McInnes MATOP_GET_ROW=1, 13331c1c02c0SLois Curfman McInnes MATOP_RESTORE_ROW=2, 13341c1c02c0SLois Curfman McInnes MATOP_MULT=3, 13351c1c02c0SLois Curfman McInnes MATOP_MULT_ADD=4, 13367c922b88SBarry Smith MATOP_MULT_TRANSPOSE=5, 13377c922b88SBarry Smith MATOP_MULT_TRANSPOSE_ADD=6, 13381c1c02c0SLois Curfman McInnes MATOP_SOLVE=7, 13391c1c02c0SLois Curfman McInnes MATOP_SOLVE_ADD=8, 13407c922b88SBarry Smith MATOP_SOLVE_TRANSPOSE=9, 13417c922b88SBarry Smith MATOP_SOLVE_TRANSPOSE_ADD=10, 13421c1c02c0SLois Curfman McInnes MATOP_LUFACTOR=11, 13431c1c02c0SLois Curfman McInnes MATOP_CHOLESKYFACTOR=12, 13441c1c02c0SLois Curfman McInnes MATOP_RELAX=13, 13451c1c02c0SLois Curfman McInnes MATOP_TRANSPOSE=14, 13461c1c02c0SLois Curfman McInnes MATOP_GETINFO=15, 13471c1c02c0SLois Curfman McInnes MATOP_EQUAL=16, 13481c1c02c0SLois Curfman McInnes MATOP_GET_DIAGONAL=17, 13491c1c02c0SLois Curfman McInnes MATOP_DIAGONAL_SCALE=18, 13501c1c02c0SLois Curfman McInnes MATOP_NORM=19, 13511c1c02c0SLois Curfman McInnes MATOP_ASSEMBLY_BEGIN=20, 13521c1c02c0SLois Curfman McInnes MATOP_ASSEMBLY_END=21, 13531c1c02c0SLois Curfman McInnes MATOP_COMPRESS=22, 13541c1c02c0SLois Curfman McInnes MATOP_SET_OPTION=23, 13551c1c02c0SLois Curfman McInnes MATOP_ZERO_ENTRIES=24, 13561c1c02c0SLois Curfman McInnes MATOP_ZERO_ROWS=25, 13571c1c02c0SLois Curfman McInnes MATOP_LUFACTOR_SYMBOLIC=26, 13581c1c02c0SLois Curfman McInnes MATOP_LUFACTOR_NUMERIC=27, 13591c1c02c0SLois Curfman McInnes MATOP_CHOLESKY_FACTOR_SYMBOLIC=28, 13601c1c02c0SLois Curfman McInnes MATOP_CHOLESKY_FACTOR_NUMERIC=29, 1361d643ce63SMatthew Knepley MATOP_SETUP_PREALLOCATION=30, 1362d643ce63SMatthew Knepley MATOP_ILUFACTOR_SYMBOLIC=31, 1363d643ce63SMatthew Knepley MATOP_ICCFACTOR_SYMBOLIC=32, 1364d643ce63SMatthew Knepley MATOP_GET_ARRAY=33, 1365d643ce63SMatthew Knepley MATOP_RESTORE_ARRAY=34, 1366e26f1c3cSBarry Smith MATOP_DUPLICATE=35, 1367d643ce63SMatthew Knepley MATOP_FORWARD_SOLVE=36, 1368d643ce63SMatthew Knepley MATOP_BACKWARD_SOLVE=37, 1369d643ce63SMatthew Knepley MATOP_ILUFACTOR=38, 1370d643ce63SMatthew Knepley MATOP_ICCFACTOR=39, 1371d643ce63SMatthew Knepley MATOP_AXPY=40, 1372d643ce63SMatthew Knepley MATOP_GET_SUBMATRICES=41, 1373d643ce63SMatthew Knepley MATOP_INCREASE_OVERLAP=42, 1374d643ce63SMatthew Knepley MATOP_GET_VALUES=43, 1375d643ce63SMatthew Knepley MATOP_COPY=44, 1376d643ce63SMatthew Knepley MATOP_PRINT_HELP=45, 1377d643ce63SMatthew Knepley MATOP_SCALE=46, 1378d643ce63SMatthew Knepley MATOP_SHIFT=47, 1379d643ce63SMatthew Knepley MATOP_DIAGONAL_SHIFT=48, 1380d643ce63SMatthew Knepley MATOP_ILUDT_FACTOR=49, 1381d643ce63SMatthew Knepley MATOP_GET_BLOCK_SIZE=50, 1382d643ce63SMatthew Knepley MATOP_GET_ROW_IJ=51, 1383d643ce63SMatthew Knepley MATOP_RESTORE_ROW_IJ=52, 1384d643ce63SMatthew Knepley MATOP_GET_COLUMN_IJ=53, 1385d643ce63SMatthew Knepley MATOP_RESTORE_COLUMN_IJ=54, 1386d643ce63SMatthew Knepley MATOP_FDCOLORING_CREATE=55, 1387d643ce63SMatthew Knepley MATOP_COLORING_PATCH=56, 1388d643ce63SMatthew Knepley MATOP_SET_UNFACTORED=57, 1389d643ce63SMatthew Knepley MATOP_PERMUTE=58, 1390d643ce63SMatthew Knepley MATOP_SET_VALUES_BLOCKED=59, 1391d643ce63SMatthew Knepley MATOP_GET_SUBMATRIX=60, 1392d643ce63SMatthew Knepley MATOP_DESTROY=61, 1393d643ce63SMatthew Knepley MATOP_VIEW=62, 1394d643ce63SMatthew Knepley MATOP_GET_MAPS=63, 1395d643ce63SMatthew Knepley MATOP_USE_SCALED_FORM=64, 1396d643ce63SMatthew Knepley MATOP_SCALE_SYSTEM=65, 1397d643ce63SMatthew Knepley MATOP_UNSCALE_SYSTEM=66, 1398ba16a8cbSKris Buschelman MATOP_SET_LOCAL_TO_GLOBAL_MAP=67, 1399d643ce63SMatthew Knepley MATOP_SET_VALUES_LOCAL=68, 1400d643ce63SMatthew Knepley MATOP_ZERO_ROWS_LOCAL=69, 1401d643ce63SMatthew Knepley MATOP_GET_ROW_MAX=70, 1402d643ce63SMatthew Knepley MATOP_CONVERT=71, 1403d643ce63SMatthew Knepley MATOP_SET_COLORING=72, 1404d643ce63SMatthew Knepley MATOP_SET_VALUES_ADIC=73, 1405d643ce63SMatthew Knepley MATOP_SET_VALUES_ADIFOR=74, 1406d643ce63SMatthew Knepley MATOP_FD_COLORING_APPLY=75, 1407d643ce63SMatthew Knepley MATOP_SET_FROM_OPTIONS=76, 1408ba16a8cbSKris Buschelman MATOP_MULT_CON=77, 1409ba16a8cbSKris Buschelman MATOP_MULT_TRANSPOSE_CON=78, 1410ba16a8cbSKris Buschelman MATOP_ILU_FACTOR_SYMBOLIC_CON=79, 1411d643ce63SMatthew Knepley MATOP_PERMUTE_SPARSIFY=80, 1412d643ce63SMatthew Knepley MATOP_MULT_MULTIPLE=81, 141341acf15aSKris Buschelman MATOP_SOLVE_MULTIPLE=82, 141441acf15aSKris Buschelman MATOP_GET_INERTIA=83, 141541acf15aSKris Buschelman MATOP_LOAD=84, 141641acf15aSKris Buschelman MATOP_IS_SYMMETRIC=85, 141741acf15aSKris Buschelman MATOP_IS_HERMITIAN=86, 141841acf15aSKris Buschelman MATOP_IS_STRUCTURALLY_SYMMETRIC=87, 141941acf15aSKris Buschelman MATOP_PB_RELAX=88, 142041acf15aSKris Buschelman MATOP_GET_VECS=89, 142141acf15aSKris Buschelman MATOP_MAT_MULT=90, 142241acf15aSKris Buschelman MATOP_MAT_MULT_SYMBOLIC=91, 142341acf15aSKris Buschelman MATOP_MAT_MULT_NUMERIC=92, 142441acf15aSKris Buschelman MATOP_PTAP=93, 142541acf15aSKris Buschelman MATOP_PTAP_SYMBOLIC=94, 1426bc011b1eSHong Zhang MATOP_PTAP_NUMERIC=95, 1427bc011b1eSHong Zhang MATOP_MAT_MULTTRANSPOSE=96, 1428bc011b1eSHong Zhang MATOP_MAT_MULTTRANSPOSE_SYMBOLIC=97, 14297ba1a0bfSKris Buschelman MATOP_MAT_MULTTRANSPOSE_NUMERIC=98, 14307ba1a0bfSKris Buschelman MATOP_PTAP_SYMBOLIC_SEQAIJ=99, 14317ba1a0bfSKris Buschelman MATOP_PTAP_NUMERIC_SEQAIJ=100, 14327ba1a0bfSKris Buschelman MATOP_PTAP_SYMBOLIC_MPIAIJ=101, 143387d4246cSBarry Smith MATOP_PTAP_NUMERIC_MPIAIJ=102, 1434f5edf698SHong Zhang MATOP_SET_VALUES_ROW = 105, 1435d5c63f83SSatish Balay MATOP_GET_ROW_UTRIANGULAR=108, 14362bebee5dSHong Zhang MATOP_RESTORE_ROW_UTRIANGULAR=109, 143769db28dcSHong Zhang MATOP_MATSOLVE=110, 1438829201f2SHong Zhang MATOP_GET_REDUNDANTMATRIX=111, 1439829201f2SHong Zhang MATOP_MATGETSEQNONZEROSTRUCTURE=115 1440fae171e0SBarry Smith } MatOperation; 1441be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatHasOperation(Mat,MatOperation,PetscTruth*); 1442be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellSetOperation(Mat,MatOperation,void(*)(void)); 1443be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellGetOperation(Mat,MatOperation,void(**)(void)); 1444be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellSetContext(Mat,void*); 1445112a2221SBarry Smith 144690ace30eSBarry Smith /* 144790ace30eSBarry Smith Codes for matrices stored on disk. By default they are 144890ace30eSBarry Smith stored in a universal format. By changing the format with 1449fb9695e5SSatish Balay PetscViewerSetFormat(viewer,PETSC_VIEWER_BINARY_NATIVE); the matrices will 145090ace30eSBarry Smith be stored in a way natural for the matrix, for example dense matrices 145190ace30eSBarry Smith would be stored as dense. Matrices stored this way may only be 145290ace30eSBarry Smith read into matrices of the same time. 145390ace30eSBarry Smith */ 145490ace30eSBarry Smith #define MATRIX_BINARY_FORMAT_DENSE -1 145590ace30eSBarry Smith 1456be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetHashTableFactor(Mat,PetscReal); 1457be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatISGetLocalMat(Mat,Mat*); 14581f4e1ec7SBarry Smith 1459d9274352SBarry Smith /*S 1460d9274352SBarry Smith MatNullSpace - Object that removes a null space from a vector, i.e. 1461d9274352SBarry Smith orthogonalizes the vector to a subsapce 1462d9274352SBarry Smith 1463f7a9e4ceSBarry Smith Level: advanced 1464d9274352SBarry Smith 1465d9274352SBarry Smith Concepts: matrix; linear operator, null space 1466d9274352SBarry Smith 14676e1639daSBarry Smith Users manual sections: 14686e1639daSBarry Smith . sec_singular 14696e1639daSBarry Smith 1470d9274352SBarry Smith .seealso: MatNullSpaceCreate() 1471d9274352SBarry Smith S*/ 147274637425SBarry Smith typedef struct _p_MatNullSpace* MatNullSpace; 1473d9274352SBarry Smith 1474be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceCreate(MPI_Comm,PetscTruth,PetscInt,const Vec[],MatNullSpace*); 1475281d1b2eSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(Vec,void*),void*); 1476be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceDestroy(MatNullSpace); 1477be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceRemove(MatNullSpace,Vec,Vec*); 1478be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceAttach(Mat,MatNullSpace); 1479be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceTest(MatNullSpace,Mat); 148074637425SBarry Smith 1481be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatReorderingSeqSBAIJ(Mat,IS); 1482be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPISBAIJSetHashTableFactor(Mat,PetscReal); 1483be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetColumnIndices(Mat,PetscInt *); 148446129b97SKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJInvertBlockDiagonal(Mat); 14853f1d51d7SBarry Smith 1486be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMAIJ(Mat,PetscInt,Mat*); 1487be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMAIJRedimension(Mat,PetscInt,Mat*); 1488be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMAIJGetAIJ(Mat,Mat*); 1489c4f061fbSSatish Balay 1490be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatComputeExplicitOperator(Mat,Mat*); 1491b0a32e0cSBarry Smith 1492be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalScaleLocal(Mat,Vec); 149304f1ad80SBarry Smith 1494fef1beadSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*); 14953ec795f1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetBase(Mat,Vec,Vec); 14963ec795f1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*); 1497e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*)); 1498e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec)); 1499e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDAddNullSpace(Mat,MatNullSpace); 1500e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt); 1501e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDResetHHistory(Mat); 1502e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctionError(Mat,PetscReal); 1503e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetPeriod(Mat,PetscInt); 1504e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDGetH(Mat,PetscScalar *); 1505e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFromOptions(Mat); 1506e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*); 1507e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*); 1508e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDDSSetUmin(Mat,PetscReal); 1509e884886fSBarry Smith 1510e884886fSBarry Smith 1511e884886fSBarry Smith typedef struct _p_MatMFFD* MatMFFD; 1512e884886fSBarry Smith 1513e884886fSBarry Smith /*E 1514e884886fSBarry Smith MatMFFDType - algorithm used to compute the h used in computing matrix-vector products via differencing of the function 1515e884886fSBarry Smith 1516e884886fSBarry Smith Level: beginner 1517e884886fSBarry Smith 1518e884886fSBarry Smith .seealso: MatMFFDSetType(), MatMFFDRegister() 1519e884886fSBarry Smith E*/ 1520e884886fSBarry Smith #define MatMFFDType const char* 1521e884886fSBarry Smith #define MATMFFD_DS "ds" 1522e884886fSBarry Smith #define MATMFFD_WP "wp" 1523e884886fSBarry Smith 1524e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetType(Mat,MatMFFDType); 1525e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatMFFD)); 1526e884886fSBarry Smith 1527e884886fSBarry Smith /*MC 1528e884886fSBarry Smith MatMFFDRegisterDynamic - Adds a method to the MatMFFD registry. 1529e884886fSBarry Smith 1530e884886fSBarry Smith Synopsis: 1531e884886fSBarry Smith PetscErrorCode MatMFFDRegisterDynamic(char *name_solver,char *path,char *name_create,PetscErrorCode (*routine_create)(MatMFFD)) 1532e884886fSBarry Smith 1533e884886fSBarry Smith Not Collective 1534e884886fSBarry Smith 1535e884886fSBarry Smith Input Parameters: 1536e884886fSBarry Smith + name_solver - name of a new user-defined compute-h module 1537e884886fSBarry Smith . path - path (either absolute or relative) the library containing this solver 1538e884886fSBarry Smith . name_create - name of routine to create method context 1539e884886fSBarry Smith - routine_create - routine to create method context 1540e884886fSBarry Smith 1541e884886fSBarry Smith Level: developer 1542e884886fSBarry Smith 1543e884886fSBarry Smith Notes: 1544e884886fSBarry Smith MatMFFDRegisterDynamic() may be called multiple times to add several user-defined solvers. 1545e884886fSBarry Smith 1546e884886fSBarry Smith If dynamic libraries are used, then the fourth input argument (routine_create) 1547e884886fSBarry Smith is ignored. 1548e884886fSBarry Smith 1549e884886fSBarry Smith Sample usage: 1550e884886fSBarry Smith .vb 1551e884886fSBarry Smith MatMFFDRegisterDynamic("my_h",/home/username/my_lib/lib/libO/solaris/mylib.a, 1552e884886fSBarry Smith "MyHCreate",MyHCreate); 1553e884886fSBarry Smith .ve 1554e884886fSBarry Smith 1555e884886fSBarry Smith Then, your solver can be chosen with the procedural interface via 1556e884886fSBarry Smith $ MatMFFDSetType(mfctx,"my_h") 1557e884886fSBarry Smith or at runtime via the option 1558e884886fSBarry Smith $ -snes_mf_type my_h 1559e884886fSBarry Smith 1560e884886fSBarry Smith .keywords: MatMFFD, register 1561e884886fSBarry Smith 1562e884886fSBarry Smith .seealso: MatMFFDRegisterAll(), MatMFFDRegisterDestroy() 1563e884886fSBarry Smith M*/ 1564e884886fSBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 1565e884886fSBarry Smith #define MatMFFDRegisterDynamic(a,b,c,d) MatMFFDRegister(a,b,c,0) 1566e884886fSBarry Smith #else 1567e884886fSBarry Smith #define MatMFFDRegisterDynamic(a,b,c,d) MatMFFDRegister(a,b,c,d) 1568e884886fSBarry Smith #endif 1569e884886fSBarry Smith 1570e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDRegisterAll(const char[]); 1571e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDRegisterDestroy(void); 1572e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDDefaultSetUmin(Mat,PetscReal); 1573e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDWPSetComputeNormU(Mat,PetscTruth); 1574e884886fSBarry Smith 1575e884886fSBarry Smith 1576be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *); 1577be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *); 15787dbadf16SMatthew Knepley 1579e9fa29b7SSatish Balay PETSC_EXTERN_CXX_END 15802eac72dbSBarry Smith #endif 1581