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 27*5c9eb25fSBarry Smith .seealso: MatSetType(), Mat, MatSolverType 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" 544b8d99adSRichard Tran Mills #define MATSEQCSRPERM "seqcsrperm" 5518def556SRichard Tran Mills #define MATMPICSRPERM "mpicsrperm" 56814655b8SBarry Smith #define MATCSRPERM "csrperm" 5781824310SBarry Smith #define MATSEQCRL "seqcrl" 58c02b24c5SBarry Smith #define MATMPICRL "mpicrl" 59c02b24c5SBarry Smith #define MATCRL "crl" 602a6744ebSBarry Smith #define MATSCATTER "scatter" 61421e10b8SBarry Smith #define MATBLOCKMAT "blockmat" 62793850ffSBarry Smith #define MATCOMPOSITE "composite" 635a7f1df3SHong Zhang #define MATSEQFFTW "seqfftw" 64d91e6319SBarry Smith 659989ab13SBarry Smith /*E 669989ab13SBarry Smith MatSolverType - String with the name of a PETSc matrix solver type. 679989ab13SBarry Smith 689989ab13SBarry Smith For example: "petsc" indicates what PETSc provides, "superlu" indicates either 699989ab13SBarry Smith SuperLU or SuperLU_Dist etc. 709989ab13SBarry Smith 719989ab13SBarry Smith 729989ab13SBarry Smith Level: beginner 739989ab13SBarry Smith 74*5c9eb25fSBarry Smith .seealso: MatGetFactor(), Mat, MatSetType(), MatType 759989ab13SBarry Smith E*/ 769989ab13SBarry Smith #define MatSolverType char* 77b24902e0SBarry Smith #define MAT_SOLVER_SPOOLES "spooles" 78b24902e0SBarry Smith #define MAT_SOLVER_SUPERLU "superlu" 79b24902e0SBarry Smith #define MAT_SOLVER_SUPERLU_DIST "superlu_dist" 80b24902e0SBarry Smith #define MAT_SOLVER_UMFPACK "umfpack" 81b24902e0SBarry Smith #define MAT_SOLVER_ESSL "essl" 82b24902e0SBarry Smith #define MAT_SOLVER_LUSOL "lusol" 83b24902e0SBarry Smith #define MAT_SOLVER_MUMPS "mumps" 84b24902e0SBarry Smith #define MAT_SOLVER_DSCPACK "dscpack" 85b24902e0SBarry Smith #define MAT_SOLVER_MATLAB "matlab" 86b24902e0SBarry Smith 87b24902e0SBarry Smith /*E 88b24902e0SBarry Smith MatFactorType - indicates what type of factorization is requested 89b24902e0SBarry Smith 90b24902e0SBarry Smith Level: beginner 91b24902e0SBarry Smith 92b24902e0SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 93b24902e0SBarry Smith 94b24902e0SBarry Smith .seealso: MatSolverType, MatGetFactor() 95b24902e0SBarry Smith E*/ 96*5c9eb25fSBarry Smith typedef enum {MAT_FACTOR_NONE, MAT_FACTOR_LU, MAT_FACTOR_CHOLESKY, MAT_FACTOR_ILU, MAT_FACTOR_ICC} MatFactorType; 97b24902e0SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetFactor(Mat,MatSolverType,MatFactorType,Mat*); 98b24902e0SBarry Smith 999989ab13SBarry Smith 100c06d978dSMatthew Knepley /* Logging support */ 101552e946dSBarry Smith #define MAT_FILE_COOKIE 1211216 /* used to indicate matrices in binary files */ 102be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_COOKIE; 103be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_FDCOLORING_COOKIE; 104be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_PARTITIONING_COOKIE; 105be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_NULLSPACE_COOKIE; 106e884886fSBarry Smith extern PetscCookie PETSCMAT_DLLEXPORT MATMFFD_COOKIE; 107c06d978dSMatthew Knepley 108ceb03754SKris Buschelman /*E 109ceb03754SKris Buschelman MatReuse - Indicates if matrices obtained from a previous call to MatGetSubMatrices() 110ceb03754SKris Buschelman or MatGetSubMatrix() are to be reused to store the new matrix values. 111ceb03754SKris Buschelman 112ceb03754SKris Buschelman Level: beginner 113ceb03754SKris Buschelman 114ceb03754SKris Buschelman Any additions/changes here MUST also be made in include/finclude/petscmat.h 115ceb03754SKris Buschelman 1160c451bc4SBarry Smith .seealso: MatGetSubMatrices(), MatGetSubMatrix(), MatDestroyMatrices(), MatConvert() 117ceb03754SKris Buschelman E*/ 118ceb03754SKris Buschelman typedef enum {MAT_INITIAL_MATRIX,MAT_REUSE_MATRIX} MatReuse; 119ceb03754SKris Buschelman 1205494a064SHong Zhang /*E 1215494a064SHong Zhang MatGetSubMatrixOption - Indicates if matrices obtained from a call to MatGetSubMatrices() 122829201f2SHong Zhang include the matrix values. Currently it is only used by MatGetSeqNonzerostructure(). 1235494a064SHong Zhang 1245494a064SHong Zhang Level: beginner 1255494a064SHong Zhang 126829201f2SHong Zhang .seealso: MatGetSeqNonzerostructure() 1275494a064SHong Zhang E*/ 1285494a064SHong Zhang typedef enum {MAT_DO_NOT_GET_VALUES,MAT_GET_VALUES} MatGetSubMatrixOption; 1295494a064SHong Zhang 130e5bd5246SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInitializePackage(const char[]); 131c06d978dSMatthew Knepley 132f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreate(MPI_Comm,Mat*); 133a69a960dSMatthew Knepley PetscPolymorphicFunction(MatCreate,(MPI_Comm comm),(comm,&A),Mat,A) 134a69a960dSMatthew Knepley PetscPolymorphicFunction(MatCreate,(),(PETSC_COMM_WORLD,&A),Mat,A) 135f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetSizes(Mat,PetscInt,PetscInt,PetscInt,PetscInt); 136f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetType(Mat,MatType); 137be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetFromOptions(Mat); 138be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUpPreallocation(Mat); 139be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRegisterAll(const char[]); 140be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat)); 14130de9b25SBarry Smith 142f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetOptionsPrefix(Mat,const char[]); 143f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAppendOptionsPrefix(Mat,const char[]); 144f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOptionsPrefix(Mat,const char*[]); 145f69a0ea3SMatthew Knepley 14630de9b25SBarry Smith /*MC 14730de9b25SBarry Smith MatRegisterDynamic - Adds a new matrix type 14830de9b25SBarry Smith 14930de9b25SBarry Smith Synopsis: 150c1ac3661SBarry Smith PetscErrorCode MatRegisterDynamic(char *name,char *path,char *name_create,PetscErrorCode (*routine_create)(Mat)) 15130de9b25SBarry Smith 15230de9b25SBarry Smith Not Collective 15330de9b25SBarry Smith 15430de9b25SBarry Smith Input Parameters: 15530de9b25SBarry Smith + name - name of a new user-defined matrix type 15630de9b25SBarry Smith . path - path (either absolute or relative) the library containing this solver 15730de9b25SBarry Smith . name_create - name of routine to create method context 15830de9b25SBarry Smith - routine_create - routine to create method context 15930de9b25SBarry Smith 16030de9b25SBarry Smith Notes: 16130de9b25SBarry Smith MatRegisterDynamic() may be called multiple times to add several user-defined solvers. 16230de9b25SBarry Smith 16330de9b25SBarry Smith If dynamic libraries are used, then the fourth input argument (routine_create) 16430de9b25SBarry Smith is ignored. 16530de9b25SBarry Smith 16630de9b25SBarry Smith Sample usage: 16730de9b25SBarry Smith .vb 16830de9b25SBarry Smith MatRegisterDynamic("my_mat",/home/username/my_lib/lib/libO/solaris/mylib.a, 16930de9b25SBarry Smith "MyMatCreate",MyMatCreate); 17030de9b25SBarry Smith .ve 17130de9b25SBarry Smith 17230de9b25SBarry Smith Then, your solver can be chosen with the procedural interface via 17330de9b25SBarry Smith $ MatSetType(Mat,"my_mat") 17430de9b25SBarry Smith or at runtime via the option 17530de9b25SBarry Smith $ -mat_type my_mat 17630de9b25SBarry Smith 17730de9b25SBarry Smith Level: advanced 17830de9b25SBarry Smith 179ab901514SBarry Smith Notes: ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values. 18030de9b25SBarry Smith If your function is not being put into a shared library then use VecRegister() instead 18130de9b25SBarry Smith 18230de9b25SBarry Smith .keywords: Mat, register 18330de9b25SBarry Smith 18430de9b25SBarry Smith .seealso: MatRegisterAll(), MatRegisterDestroy() 18530de9b25SBarry Smith 18630de9b25SBarry Smith M*/ 187273d9f13SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 188273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,0) 189273d9f13SBarry Smith #else 190273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,d) 19130de9b25SBarry Smith #endif 19230de9b25SBarry Smith 1939c666560SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRegisterAll(const char[]); 194273d9f13SBarry Smith extern PetscTruth MatRegisterAllCalled; 195b0a32e0cSBarry Smith extern PetscFList MatList; 19628988994SBarry Smith 197be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqDense(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*); 198be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIDense(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*); 199be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 200ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,nz,nnz,&A),Mat,A) 201ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,&A),Mat,A) 202ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,0,nnz,&A),Mat,A) 203ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,&A),Mat,A) 204ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,A)) 205ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,m,n,0,nnz,A)) 206ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,A)) 207be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 208ba966639SSatish 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) 209ba966639SSatish 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) 210ba966639SSatish 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) 211ba966639SSatish 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) 212ba966639SSatish 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)) 213ba966639SSatish 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)) 214ba966639SSatish 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)) 215ba966639SSatish 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) 216ba966639SSatish 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) 217ba966639SSatish 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) 218ba966639SSatish 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) 219ba966639SSatish 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)) 220ba966639SSatish 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)) 221ba966639SSatish 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)) 22263ab93bfSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *); 2238d7a6e47SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJWithSplitArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],PetscInt[],PetscInt[],PetscScalar[],Mat*); 2248d7a6e47SBarry Smith 225be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIRowbs(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 226be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 227ba966639SSatish 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) 228ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A) 229ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A) 230ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A) 231ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A)) 232ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A)) 233ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A)) 234be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 235ba966639SSatish 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) 236ba966639SSatish 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) 237ba966639SSatish 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) 238ba966639SSatish 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) 239ba966639SSatish 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)) 240ba966639SSatish 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)) 241ba966639SSatish 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)) 242ba966639SSatish 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) 243ba966639SSatish 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) 244ba966639SSatish 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) 245ba966639SSatish 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) 246ba966639SSatish 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)) 247ba966639SSatish 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)) 248ba966639SSatish 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)) 249be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAdj(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscInt[],Mat*); 250be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 251ba966639SSatish 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) 252ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A) 253ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A) 254ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A) 255ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A)) 256ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A)) 257ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A)) 258be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPISBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 259ba966639SSatish 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) 260ba966639SSatish 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) 261ba966639SSatish 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) 262ba966639SSatish 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) 263ba966639SSatish 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)) 264ba966639SSatish 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)) 265ba966639SSatish 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)) 266ba966639SSatish 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) 267ba966639SSatish 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) 268ba966639SSatish 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) 269ba966639SSatish 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) 270ba966639SSatish 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)) 271ba966639SSatish 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)) 272ba966639SSatish 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)) 273be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateShell(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,void *,Mat*); 274ba966639SSatish 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) 275ba966639SSatish Balay PetscPolymorphicFunction(MatCreateShell,(PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx),(PETSC_COMM_WORLD,m,n,M,N,ctx,&A),Mat,A) 276be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateAdic(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,void (*)(void),Mat*); 277be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateNormal(Mat,Mat*); 278ba966639SSatish Balay PetscPolymorphicFunction(MatCreateNormal,(Mat mat),(mat,&A),Mat,A) 279ab92ecdeSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateLRC(Mat,Mat,Mat,Mat*); 280284134d9SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateIS(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,ISLocalToGlobalMapping,Mat*); 2811472f72bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 2821472f72bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPICRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 2832a6744ebSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateScatter(MPI_Comm,VecScatter,Mat*); 2842a6744ebSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScatterSetVecScatter(Mat,VecScatter); 2852a6744ebSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScatterGetVecScatter(Mat,VecScatter*); 2868cdf2d9bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateBlockMat(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,Mat*); 287793850ffSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompositeAddMat(Mat,Mat); 288b52f573bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompositeMerge(Mat); 289793850ffSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateComposite(MPI_Comm,PetscInt,const Mat*,Mat*); 2905a7f1df3SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqFFTW(MPI_Comm,PetscInt,const PetscInt[],Mat*); 2911472f72bSBarry Smith 292f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUp(Mat); 293be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroy(Mat); 29421c89e3eSBarry Smith 2950c451bc4SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConjugate(Mat); 29699cafbc1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRealPart(Mat); 29799cafbc1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatImaginaryPart(Mat); 29899cafbc1SBarry Smith 2998ed539a5SBarry Smith /* ------------------------------------------------------------*/ 300be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 301be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 30287d4246cSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesRow(Mat,PetscInt,const PetscScalar[]); 303f259bd47SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesRowLocal(Mat,PetscInt,const PetscScalar[]); 30484cb2905SBarry Smith 3052ef4de8bSBarry Smith /*S 3062ef4de8bSBarry Smith MatStencil - Data structure (C struct) for storing information about a single row or 3072ef4de8bSBarry Smith column of a matrix as index on an associated grid. 3082ef4de8bSBarry Smith 3092ef4de8bSBarry Smith Level: beginner 3102ef4de8bSBarry Smith 3112ef4de8bSBarry Smith Concepts: matrix; linear operator 3122ef4de8bSBarry Smith 313d7bd93baSBarry Smith .seealso: MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockStencil() 3142ef4de8bSBarry Smith S*/ 315435da068SBarry Smith typedef struct { 316c1ac3661SBarry Smith PetscInt k,j,i,c; 317435da068SBarry Smith } MatStencil; 3182ef4de8bSBarry Smith 319be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode); 320be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode); 321be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt); 322435da068SBarry Smith 323be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetColoring(Mat,ISColoring); 324be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesAdic(Mat,void*); 325be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesAdifor(Mat,PetscInt,void*); 3263a7fca6bSBarry Smith 327d91e6319SBarry Smith /*E 328d91e6319SBarry Smith MatAssemblyType - Indicates if the matrix is now to be used, or if you plan 329d91e6319SBarry Smith to continue to add values to it 330d91e6319SBarry Smith 331d91e6319SBarry Smith Level: beginner 332d91e6319SBarry Smith 333d91e6319SBarry Smith .seealso: MatAssemblyBegin(), MatAssemblyEnd() 334d91e6319SBarry Smith E*/ 3356d4a8577SBarry Smith typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType; 336be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssemblyBegin(Mat,MatAssemblyType); 337be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssemblyEnd(Mat,MatAssemblyType); 338be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssembled(Mat,PetscTruth*); 3394f9c727eSBarry Smith 3401d73ed98SBarry Smith 34130de9b25SBarry Smith 342d91e6319SBarry Smith /*E 343d91e6319SBarry Smith MatOption - Options that may be set for a matrix and its behavior or storage 344d91e6319SBarry Smith 345d91e6319SBarry Smith Level: beginner 346d91e6319SBarry Smith 3470a835dfdSSatish Balay Any additions/changes here MUST also be made in include/finclude/petscmat.h 348d91e6319SBarry Smith 349d91e6319SBarry Smith .seealso: MatSetOption() 350d91e6319SBarry Smith E*/ 351512a5fc5SBarry Smith typedef enum {MAT_ROW_ORIENTED,MAT_NEW_NONZERO_LOCATIONS, 3524e0d8c25SBarry Smith MAT_SYMMETRIC, 3534e0d8c25SBarry Smith MAT_STRUCTURALLY_SYMMETRIC, 3548047e77bSBarry Smith MAT_NEW_DIAGONALS,MAT_IGNORE_OFF_PROC_ENTRIES, 3554e0d8c25SBarry Smith MAT_NEW_NONZERO_LOCATION_ERR, 3564e0d8c25SBarry Smith MAT_NEW_NONZERO_ALLOCATION_ERR,MAT_USE_HASH_TABLE, 3574e0d8c25SBarry Smith MAT_KEEP_ZEROED_ROWS,MAT_IGNORE_ZERO_ENTRIES, 3584e0d8c25SBarry Smith MAT_USE_INODES, 3594e0d8c25SBarry Smith MAT_HERMITIAN, 3604e0d8c25SBarry Smith MAT_SYMMETRY_ETERNAL, 3614e0d8c25SBarry Smith MAT_USE_COMPRESSEDROW, 3624e0d8c25SBarry Smith MAT_IGNORE_LOWER_TRIANGULAR,MAT_ERROR_LOWER_TRIANGULAR, 3634e0d8c25SBarry Smith MAT_GETROW_UPPERTRIANGULAR} MatOption; 364290bbb0aSBarry Smith extern const char *MatOptions[]; 3654e0d8c25SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetOption(Mat,MatOption,PetscTruth); 366be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetType(Mat,MatType*); 367ba966639SSatish Balay PetscPolymorphicFunction(MatGetType,(Mat mat),(mat,&t),MatType,t) 36884cb2905SBarry Smith 369be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]); 370be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 371be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 372f5edf698SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowUpperTriangular(Mat); 373f5edf698SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRowUpperTriangular(Mat); 374be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 375be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 376be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumnVector(Mat,Vec,PetscInt); 377be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetArray(Mat,PetscScalar *[]); 378ba966639SSatish Balay PetscPolymorphicFunction(MatGetArray,(Mat mat),(mat,&a),PetscScalar*,a) 379be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreArray(Mat,PetscScalar *[]); 380be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBlockSize(Mat,PetscInt *); 381ba966639SSatish Balay PetscPolymorphicFunction(MatGetBlockSize,(Mat mat),(mat,&a),PetscInt,a) 382be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetBlockSize(Mat,PetscInt); 3837b80b807SBarry Smith 3841620fd73SBarry Smith 385be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMult(Mat,Vec,Vec); 386be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultAdd(Mat,Vec,Vec,Vec); 387ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultAdd,(Mat A,Vec x,Vec y),(A,x,y,y)) 388be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTranspose(Mat,Vec,Vec); 389be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsTranspose(Mat,Mat,PetscReal,PetscTruth*); 390ba966639SSatish Balay PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B,PetscReal tol),(A,B,tol,&t),PetscTruth,t) 391e37f64cdSBarry Smith PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B),(A,B,0,&t),PetscTruth,t) 3921cbb95d3SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitianTranspose(Mat,Mat,PetscReal,PetscTruth*); 393be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeAdd(Mat,Vec,Vec,Vec); 394ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultTransposeAdd,(Mat A,Vec x,Vec y),(A,x,y,y)) 395be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultConstrained(Mat,Vec,Vec); 396be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeConstrained(Mat,Vec,Vec); 3972bebee5dSHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatSolve(Mat,Mat,Mat); 398f5edf698SHong Zhang 399d91e6319SBarry Smith /*E 400d91e6319SBarry Smith MatDuplicateOption - Indicates if a duplicated sparse matrix should have 401d91e6319SBarry Smith its numerical values copied over or just its nonzero structure. 402d91e6319SBarry Smith 403d91e6319SBarry Smith Level: beginner 404d91e6319SBarry Smith 405d91e6319SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 406d91e6319SBarry Smith 407d91e6319SBarry Smith .seealso: MatDuplicate() 408d91e6319SBarry Smith E*/ 4092e8a6d31SBarry Smith typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES} MatDuplicateOption; 4102e8a6d31SBarry Smith 411f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert(Mat,MatType,MatReuse,Mat*); 412ba966639SSatish Balay PetscPolymorphicFunction(MatConvert,(Mat A,MatType t),(A,t,MAT_INITIAL_MATRIX,&a),Mat,a) 413be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDuplicate(Mat,MatDuplicateOption,Mat*); 414ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A,MatDuplicateOption o),(A,o,&a),Mat,a) 415ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A),(A,MAT_COPY_VALUES,&a),Mat,a) 41694a9d846SBarry Smith 417d91e6319SBarry Smith /*E 418d91e6319SBarry Smith MatStructure - Indicates if the matrix has the same nonzero structure 419d91e6319SBarry Smith 420d91e6319SBarry Smith Level: beginner 421d91e6319SBarry Smith 422d91e6319SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 423d91e6319SBarry Smith 42494b7f48cSBarry Smith .seealso: MatCopy(), KSPSetOperators(), PCSetOperators() 425d91e6319SBarry Smith E*/ 426c537a176SHong Zhang typedef enum {SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER,SUBSET_NONZERO_PATTERN} MatStructure; 427cb5b572fSBarry Smith 428be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCopy(Mat,Mat,MatStructure); 429be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatView(Mat,PetscViewer); 430be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsSymmetric(Mat,PetscReal,PetscTruth*); 431ba966639SSatish Balay PetscPolymorphicFunction(MatIsSymmetric,(Mat A,PetscReal tol),(A,tol,&t),PetscTruth,t) 432e37f64cdSBarry Smith PetscPolymorphicFunction(MatIsSymmetric,(Mat A),(A,0,&t),PetscTruth,t) 433be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsStructurallySymmetric(Mat,PetscTruth*); 434ba966639SSatish Balay PetscPolymorphicFunction(MatIsStructurallySymmetric,(Mat A),(A,&t),PetscTruth,t) 4351cbb95d3SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitian(Mat,PetscReal,PetscTruth*); 4361cbb95d3SBarry Smith PetscPolymorphicFunction(MatIsHermitian,(Mat A),(A,0,&t),PetscTruth,t) 437be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsSymmetricKnown(Mat,PetscTruth*,PetscTruth*); 438be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitianKnown(Mat,PetscTruth*,PetscTruth*); 43964c62002SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMissingDiagonal(Mat,PetscTruth *,PetscInt *); 440f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLoad(PetscViewer,MatType,Mat*); 441ba966639SSatish Balay PetscPolymorphicFunction(MatLoad,(PetscViewer v,MatType t),(v,t,&a),Mat,a) 4427b80b807SBarry Smith 4438f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *); 4448f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRowIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *); 4458f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumnIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *); 4468f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreColumnIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *); 447d4fbbf0eSBarry Smith 448d91e6319SBarry Smith /*S 449d91e6319SBarry Smith MatInfo - Context of matrix information, used with MatGetInfo() 450d91e6319SBarry Smith 451d91e6319SBarry Smith In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE 452d91e6319SBarry Smith 453d91e6319SBarry Smith Level: intermediate 454d91e6319SBarry Smith 455d91e6319SBarry Smith Concepts: matrix^nonzero information 456d91e6319SBarry Smith 457d9274352SBarry Smith .seealso: MatGetInfo(), MatInfoType 458d91e6319SBarry Smith S*/ 4594e220ebcSLois Curfman McInnes typedef struct { 460b0a32e0cSBarry Smith PetscLogDouble rows_global,columns_global; /* number of global rows and columns */ 461b0a32e0cSBarry Smith PetscLogDouble rows_local,columns_local; /* number of local rows and columns */ 462b0a32e0cSBarry Smith PetscLogDouble block_size; /* block size */ 463b0a32e0cSBarry Smith PetscLogDouble nz_allocated,nz_used,nz_unneeded; /* number of nonzeros */ 464b0a32e0cSBarry Smith PetscLogDouble memory; /* memory allocated */ 465b0a32e0cSBarry Smith PetscLogDouble assemblies; /* number of matrix assemblies called */ 466b0a32e0cSBarry Smith PetscLogDouble mallocs; /* number of mallocs during MatSetValues() */ 467b0a32e0cSBarry Smith PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */ 468b0a32e0cSBarry Smith PetscLogDouble factor_mallocs; /* number of mallocs during factorization */ 4694e220ebcSLois Curfman McInnes } MatInfo; 4704e220ebcSLois Curfman McInnes 471d9274352SBarry Smith /*E 472d9274352SBarry Smith MatInfoType - Indicates if you want information about the local part of the matrix, 473d9274352SBarry Smith the entire parallel matrix or the maximum over all the local parts. 474d9274352SBarry Smith 475d9274352SBarry Smith Level: beginner 476d9274352SBarry Smith 477d9274352SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 478d9274352SBarry Smith 479d9274352SBarry Smith .seealso: MatGetInfo(), MatInfo 480d9274352SBarry Smith E*/ 4817b80b807SBarry Smith typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType; 482be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetInfo(Mat,MatInfoType,MatInfo*); 483be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatValid(Mat,PetscTruth*); 484be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetDiagonal(Mat,Vec); 485985db425SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMax(Mat,Vec,PetscInt[]); 486985db425SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMin(Mat,Vec,PetscInt[]); 487985db425SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMaxAbs(Mat,Vec,PetscInt[]); 48886697f06SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowSum(Mat,Vec); 489fc4dec0aSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatTranspose(Mat,MatReuse,Mat*); 490fc4dec0aSBarry Smith PetscPolymorphicFunction(MatTranspose,(Mat A),(A,MAT_INITIAL_MATRIX,&t),Mat,t) 491be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPermute(Mat,IS,IS,Mat *); 492ba966639SSatish Balay PetscPolymorphicFunction(MatPermute,(Mat A,IS is1,IS is2),(A,is1,is2,&t),Mat,t) 493be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPermuteSparsify(Mat,PetscInt,PetscReal,PetscReal,IS,IS,Mat *); 494be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalScale(Mat,Vec,Vec); 495be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalSet(Mat,Vec,InsertMode); 496be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatEqual(Mat,Mat,PetscTruth*); 497ba966639SSatish Balay PetscPolymorphicFunction(MatEqual,(Mat A,Mat B),(A,B,&t),PetscTruth,t) 498be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultEqual(Mat,Mat,PetscInt,PetscTruth*); 499be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultAddEqual(Mat,Mat,PetscInt,PetscTruth*); 500be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeEqual(Mat,Mat,PetscInt,PetscTruth*); 501be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscTruth*); 5027b80b807SBarry Smith 503be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNorm(Mat,NormType,PetscReal *); 504ba966639SSatish Balay PetscPolymorphicFunction(MatNorm,(Mat A,NormType t),(A,t,&n),PetscReal,n) 505be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroEntries(Mat); 506f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar); 507f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsIS(Mat,IS,PetscScalar); 508f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroColumns(Mat,PetscInt,const PetscInt [],const PetscScalar*); 509f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroColumnsIS(Mat,IS,const PetscScalar*); 5107b80b807SBarry Smith 511be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatUseScaledForm(Mat,PetscTruth); 512be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScaleSystem(Mat,Vec,Vec); 513be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatUnScaleSystem(Mat,Vec,Vec); 5145ef9f2a5SBarry Smith 515be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSize(Mat,PetscInt*,PetscInt*); 516be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalSize(Mat,PetscInt*,PetscInt*); 517be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRange(Mat,PetscInt*,PetscInt*); 5188ee2e534SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRanges(Mat,const PetscInt**); 519f4b8409dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRangeColumn(Mat,PetscInt*,PetscInt*); 5207b80b807SBarry Smith 521be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]); 522be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroyMatrices(PetscInt,Mat *[]); 523be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrix(Mat,IS,IS,PetscInt,MatReuse,Mat *); 524abf0e977SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrixRaw(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscInt,MatReuse,Mat *); 525829201f2SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSeqNonzeroStructure(Mat,Mat *[]); 526829201f2SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroySeqNonzeroStructure(Mat *[]); 5275494a064SHong Zhang 528be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge(MPI_Comm,Mat,PetscInt,MatReuse,Mat*); 529be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPI(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*); 530be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPISymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*); 531be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPINumeric(Mat,Mat); 532be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroy_MPIAIJ_SeqsToMPI(Mat); 533be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalMat(Mat,MatReuse,Mat*); 534be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*); 535be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,PetscInt*,Mat*); 536dd6ea824SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBrowsOfAoCols(Mat,Mat,MatReuse,PetscInt**,MatScalar**,Mat*); 53743eb5e2fSMatthew Knepley #if defined (PETSC_USE_CTABLE) 538cd116e26SSatish Balay #include "petscctable.h" 53943eb5e2fSMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetCommunicationStructs(Mat, Vec *, PetscTable *, VecScatter *); 54043eb5e2fSMatthew Knepley #else 54143eb5e2fSMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetCommunicationStructs(Mat, Vec *, PetscInt *[], VecScatter *); 54243eb5e2fSMatthew Knepley #endif 5438efafbd8SBarry Smith 544be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt); 5457b80b807SBarry Smith 546be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*); 547be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultSymbolic(Mat,Mat,PetscReal,Mat*); 548be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultNumeric(Mat,Mat,Mat); 54922440eb1SKris Buschelman 550be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*); 551be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAPSymbolic(Mat,Mat,PetscReal,Mat*); 552be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAPNumeric(Mat,Mat,Mat); 55322440eb1SKris Buschelman 554be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTranspose(Mat,Mat,MatReuse,PetscReal,Mat*); 555be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTransposeSymbolic(Mat,Mat,PetscReal,Mat*); 556be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTransposeNumeric(Mat,Mat,Mat); 557bc011b1eSHong Zhang 558f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAXPY(Mat,PetscScalar,Mat,MatStructure); 55904aac2b0SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAYPX(Mat,PetscScalar,Mat,MatStructure); 560be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompress(Mat); 5611c741599SHong Zhang 562f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScale(Mat,PetscScalar); 563f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShift(Mat,PetscScalar); 5647b80b807SBarry Smith 565be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping); 566be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetLocalToGlobalMappingBlock(Mat,ISLocalToGlobalMapping); 567f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar); 568f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsLocalIS(Mat,IS,PetscScalar); 569be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 570be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 571052efed2SBarry Smith 572be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStashSetInitialSize(Mat,PetscInt,PetscInt); 573be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*); 57490f02eecSBarry Smith 575be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInterpolate(Mat,Vec,Vec); 5760c451bc4SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInterpolateAdd(Mat,Vec,Vec,Vec); 577ba966639SSatish Balay PetscPolymorphicSubroutine(MatInterpolateAdd,(Mat A,Vec x,Vec y),(A,x,y,y)) 578be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestrict(Mat,Vec,Vec); 579be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetVecs(Mat,Vec*,Vec*); 58069db28dcSHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRedundantMatrix(Mat,PetscInt,MPI_Comm,PetscInt,MatReuse,Mat*); 581bd481603SBarry Smith 582bd481603SBarry Smith /*MC 5831620fd73SBarry Smith MatSetValue - Set a single entry into a matrix. 5841620fd73SBarry Smith 5851620fd73SBarry Smith Not collective 5861620fd73SBarry Smith 5871620fd73SBarry Smith Input Parameters: 5881620fd73SBarry Smith + m - the matrix 5891620fd73SBarry Smith . row - the row location of the entry 5901620fd73SBarry Smith . col - the column location of the entry 5911620fd73SBarry Smith . value - the value to insert 5921620fd73SBarry Smith - mode - either INSERT_VALUES or ADD_VALUES 5931620fd73SBarry Smith 5941620fd73SBarry Smith Notes: 5951620fd73SBarry Smith For efficiency one should use MatSetValues() and set several or many 5961620fd73SBarry Smith values simultaneously if possible. 5971620fd73SBarry Smith 5981620fd73SBarry Smith Level: beginner 5991620fd73SBarry Smith 6001620fd73SBarry Smith .seealso: MatSetValues(), MatSetValueLocal() 6011620fd73SBarry Smith M*/ 6021620fd73SBarry 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);} 6031620fd73SBarry Smith 6041620fd73SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatGetValue(Mat v,PetscInt i,PetscInt j,PetscScalar *va) {return MatGetValues(v,1,&i,1,&j,va);} 6051620fd73SBarry Smith 6061620fd73SBarry 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);} 6071620fd73SBarry Smith 6081620fd73SBarry Smith /*MC 6090d2a7ff2SSatish Balay MatPreallocateInitialize - Begins the block of code that will count the number of nonzeros per 610bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 611bd481603SBarry Smith 612bd481603SBarry Smith Synopsis: 613c1ac3661SBarry Smith PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz) 614bd481603SBarry Smith 615bd481603SBarry Smith Collective on MPI_Comm 616bd481603SBarry Smith 617bd481603SBarry Smith Input Parameters: 618bd481603SBarry Smith + comm - the communicator that will share the eventually allocated matrix 619859fcb39SBarry Smith . nrows - the number of LOCAL rows in the matrix 620859fcb39SBarry Smith - ncols - the number of LOCAL columns in the matrix 621bd481603SBarry Smith 622bd481603SBarry Smith Output Parameters: 623bd481603SBarry Smith + dnz - the array that will be passed to the matrix preallocation routines 624bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 625bd481603SBarry Smith 626bd481603SBarry Smith 627bd481603SBarry Smith Level: intermediate 628bd481603SBarry Smith 629bd481603SBarry Smith Notes: 630bd481603SBarry Smith See the chapter in the users manual on performance for more details 631bd481603SBarry Smith 6321d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 633bd481603SBarry Smith 634bd481603SBarry Smith Use MatPreallocateInitializeSymmetric() for symmetric matrices (MPISBAIJ matrices) 635bd481603SBarry Smith 6361620fd73SBarry Smith This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize(). 6371620fd73SBarry Smith 638bd481603SBarry Smith Concepts: preallocation^Matrix 639bd481603SBarry Smith 640bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(), 641bd481603SBarry Smith MatPreallocateInitializeSymmetric(), MatPreallocateSymmetricSetLocal() 642bd481603SBarry Smith M*/ 643c4f061fbSSatish Balay #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \ 6447c922b88SBarry Smith { \ 645c1ac3661SBarry Smith PetscErrorCode _4_ierr; PetscInt __tmp = (nrows),__ctmp = (ncols),__rstart,__start,__end; \ 646c1ac3661SBarry Smith _4_ierr = PetscMalloc(2*__tmp*sizeof(PetscInt),&dnz);CHKERRQ(_4_ierr);onz = dnz + __tmp;\ 647c1ac3661SBarry Smith _4_ierr = PetscMemzero(dnz,2*__tmp*sizeof(PetscInt));CHKERRQ(_4_ierr);\ 648c1ac3661SBarry Smith _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ctmp;\ 649c1ac3661SBarry Smith _4_ierr = MPI_Scan(&__tmp,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __tmp; 6507c922b88SBarry Smith 651bd481603SBarry Smith /*MC 6520d2a7ff2SSatish Balay MatPreallocateSymmetricInitialize - Begins the block of code that will count the number of nonzeros per 653bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 654bd481603SBarry Smith 655bd481603SBarry Smith Synopsis: 656c1ac3661SBarry Smith PetscErrorCode MatPreallocateSymmetricInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz) 657bd481603SBarry Smith 658bd481603SBarry Smith Collective on MPI_Comm 659bd481603SBarry Smith 660bd481603SBarry Smith Input Parameters: 661bd481603SBarry Smith + comm - the communicator that will share the eventually allocated matrix 662859fcb39SBarry Smith . nrows - the number of LOCAL rows in the matrix 663859fcb39SBarry Smith - ncols - the number of LOCAL columns in the matrix 664bd481603SBarry Smith 665bd481603SBarry Smith Output Parameters: 666bd481603SBarry Smith + dnz - the array that will be passed to the matrix preallocation routines 667bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 668bd481603SBarry Smith 669bd481603SBarry Smith 670bd481603SBarry Smith Level: intermediate 671bd481603SBarry Smith 672bd481603SBarry Smith Notes: 673bd481603SBarry Smith See the chapter in the users manual on performance for more details 674bd481603SBarry Smith 6751d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 676bd481603SBarry Smith 6771620fd73SBarry Smith This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize(). 6781620fd73SBarry Smith 679bd481603SBarry Smith Concepts: preallocation^Matrix 680bd481603SBarry Smith 681bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(), 682bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal() 683bd481603SBarry Smith M*/ 684222b16d4SBarry Smith #define MatPreallocateSymmetricInitialize(comm,nrows,ncols,dnz,onz) 0; \ 685222b16d4SBarry Smith { \ 686c1ac3661SBarry Smith PetscErrorCode _4_ierr; PetscInt __tmp = (nrows),__ctmp = (ncols),__rstart,__end; \ 687c1ac3661SBarry Smith _4_ierr = PetscMalloc(2*__tmp*sizeof(PetscInt),&dnz);CHKERRQ(_4_ierr);onz = dnz + __tmp;\ 688c1ac3661SBarry Smith _4_ierr = PetscMemzero(dnz,2*__tmp*sizeof(PetscInt));CHKERRQ(_4_ierr);\ 689c1ac3661SBarry Smith _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr);\ 690c1ac3661SBarry Smith _4_ierr = MPI_Scan(&__tmp,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __tmp; 691222b16d4SBarry Smith 692bd481603SBarry Smith /*MC 693bd481603SBarry Smith MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be 694bd481603SBarry Smith inserted using a local number of the rows and columns 695bd481603SBarry Smith 696bd481603SBarry Smith Synopsis: 697c1ac3661SBarry Smith PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 698bd481603SBarry Smith 699bd481603SBarry Smith Not Collective 700bd481603SBarry Smith 701bd481603SBarry Smith Input Parameters: 702bd481603SBarry Smith + map - the mapping between local numbering and global numbering 703bd481603SBarry Smith . nrows - the number of rows indicated 7041d73ed98SBarry Smith . rows - the indices of the rows 705bd481603SBarry Smith . ncols - the number of columns in the matrix 706bd481603SBarry Smith . cols - the columns indicated 707bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 708bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 709bd481603SBarry Smith 710bd481603SBarry Smith 711bd481603SBarry Smith Level: intermediate 712bd481603SBarry Smith 713bd481603SBarry Smith Notes: 714bd481603SBarry Smith See the chapter in the users manual on performance for more details 715bd481603SBarry Smith 7161d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 717bd481603SBarry Smith 718bd481603SBarry Smith Concepts: preallocation^Matrix 719bd481603SBarry Smith 720bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(), 721bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal() 722bd481603SBarry Smith M*/ 723c4f061fbSSatish Balay #define MatPreallocateSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\ 724c4f061fbSSatish Balay {\ 725c1ac3661SBarry Smith PetscInt __l;\ 726ef66eb69SBarry Smith _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\ 727ef66eb69SBarry Smith _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\ 728c4f061fbSSatish Balay for (__l=0;__l<nrows;__l++) {\ 729ef66eb69SBarry Smith _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 730c4f061fbSSatish Balay }\ 731c4f061fbSSatish Balay } 732c4f061fbSSatish Balay 733bd481603SBarry Smith /*MC 734bd481603SBarry Smith MatPreallocateSymmetricSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be 735bd481603SBarry Smith inserted using a local number of the rows and columns 736bd481603SBarry Smith 737bd481603SBarry Smith Synopsis: 738c1ac3661SBarry Smith PetscErrorCode MatPreallocateSymmetricSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 739bd481603SBarry Smith 740bd481603SBarry Smith Not Collective 741bd481603SBarry Smith 742bd481603SBarry Smith Input Parameters: 743bd481603SBarry Smith + map - the mapping between local numbering and global numbering 744bd481603SBarry Smith . nrows - the number of rows indicated 7451d73ed98SBarry Smith . rows - the indices of the rows 746bd481603SBarry Smith . ncols - the number of columns in the matrix 747bd481603SBarry Smith . cols - the columns indicated 748bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 749bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 750bd481603SBarry Smith 751bd481603SBarry Smith 752bd481603SBarry Smith Level: intermediate 753bd481603SBarry Smith 754bd481603SBarry Smith Notes: 755bd481603SBarry Smith See the chapter in the users manual on performance for more details 756bd481603SBarry Smith 757bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 758bd481603SBarry Smith 759bd481603SBarry Smith Concepts: preallocation^Matrix 760bd481603SBarry Smith 761bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(), 762bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal() 763bd481603SBarry Smith M*/ 764d3d32019SBarry Smith #define MatPreallocateSymmetricSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\ 765d3d32019SBarry Smith {\ 766c1ac3661SBarry Smith PetscInt __l;\ 767d3d32019SBarry Smith _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\ 768d3d32019SBarry Smith _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\ 769d3d32019SBarry Smith for (__l=0;__l<nrows;__l++) {\ 770d3d32019SBarry Smith _4_ierr = MatPreallocateSymmetricSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 771d3d32019SBarry Smith }\ 772d3d32019SBarry Smith } 773d3d32019SBarry Smith 774bd481603SBarry Smith /*MC 775bd481603SBarry Smith MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be 776bd481603SBarry Smith inserted using a local number of the rows and columns 777bd481603SBarry Smith 778bd481603SBarry Smith Synopsis: 779c1ac3661SBarry Smith PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 780bd481603SBarry Smith 781bd481603SBarry Smith Not Collective 782bd481603SBarry Smith 783bd481603SBarry Smith Input Parameters: 78464075487SBarry Smith + row - the row 785bd481603SBarry Smith . ncols - the number of columns in the matrix 786a50f8bf6SHong Zhang - cols - the columns indicated 787a50f8bf6SHong Zhang 788a50f8bf6SHong Zhang Output Parameters: 789a50f8bf6SHong Zhang + dnz - the array that will be passed to the matrix preallocation routines 790bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 791bd481603SBarry Smith 792bd481603SBarry Smith 793bd481603SBarry Smith Level: intermediate 794bd481603SBarry Smith 795bd481603SBarry Smith Notes: 796bd481603SBarry Smith See the chapter in the users manual on performance for more details 797bd481603SBarry Smith 798bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 799bd481603SBarry Smith 8001620fd73SBarry Smith This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize(). 8011620fd73SBarry Smith 802bd481603SBarry Smith Concepts: preallocation^Matrix 803bd481603SBarry Smith 804bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(), 805bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal() 806bd481603SBarry Smith M*/ 807c4f061fbSSatish Balay #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\ 808c1ac3661SBarry Smith { PetscInt __i; \ 8098ee2e534SBarry Smith if (row < __rstart) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Trying to set preallocation for row %D less than first local row %D",row,__rstart);\ 8108ee2e534SBarry 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);\ 8117c922b88SBarry Smith for (__i=0; __i<nc; __i++) {\ 81264075487SBarry Smith if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \ 8137cc688d7SBarry Smith else dnz[row - __rstart]++;\ 8147c922b88SBarry Smith }\ 8157c922b88SBarry Smith } 8167c922b88SBarry Smith 817bd481603SBarry Smith /*MC 818bd481603SBarry Smith MatPreallocateSymmetricSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be 819bd481603SBarry Smith inserted using a local number of the rows and columns 820bd481603SBarry Smith 821bd481603SBarry Smith Synopsis: 822c1ac3661SBarry Smith PetscErrorCode MatPreallocateSymmetricSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 823bd481603SBarry Smith 824bd481603SBarry Smith Not Collective 825bd481603SBarry Smith 826bd481603SBarry Smith Input Parameters: 827bd481603SBarry Smith + nrows - the number of rows indicated 8281d73ed98SBarry Smith . rows - the indices of the rows 829bd481603SBarry Smith . ncols - the number of columns in the matrix 830bd481603SBarry Smith . cols - the columns indicated 831bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 832bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 833bd481603SBarry Smith 834bd481603SBarry Smith 835bd481603SBarry Smith Level: intermediate 836bd481603SBarry Smith 837bd481603SBarry Smith Notes: 838bd481603SBarry Smith See the chapter in the users manual on performance for more details 839bd481603SBarry Smith 840bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 841bd481603SBarry Smith 8421620fd73SBarry Smith This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize(). 8431620fd73SBarry Smith 844bd481603SBarry Smith Concepts: preallocation^Matrix 845bd481603SBarry Smith 846bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(), 847bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal() 848bd481603SBarry Smith M*/ 849d3d32019SBarry Smith #define MatPreallocateSymmetricSet(row,nc,cols,dnz,onz) 0;\ 850c1ac3661SBarry Smith { PetscInt __i; \ 851d3d32019SBarry Smith for (__i=0; __i<nc; __i++) {\ 852d3d32019SBarry Smith if (cols[__i] >= __end) onz[row - __rstart]++; \ 853d3d32019SBarry Smith else if (cols[__i] >= row) dnz[row - __rstart]++;\ 854d3d32019SBarry Smith }\ 855d3d32019SBarry Smith } 856d3d32019SBarry Smith 857bd481603SBarry Smith /*MC 8580d2a7ff2SSatish Balay MatPreallocateFinalize - Ends the block of code that will count the number of nonzeros per 859bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 860bd481603SBarry Smith 861bd481603SBarry Smith Synopsis: 862c1ac3661SBarry Smith PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz) 863bd481603SBarry Smith 864bd481603SBarry Smith Collective on MPI_Comm 865bd481603SBarry Smith 866bd481603SBarry Smith Input Parameters: 867bd481603SBarry Smith + dnz - the array that will be passed to the matrix preallocation routines 868bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 869bd481603SBarry Smith 870bd481603SBarry Smith 871bd481603SBarry Smith Level: intermediate 872bd481603SBarry Smith 873bd481603SBarry Smith Notes: 874bd481603SBarry Smith See the chapter in the users manual on performance for more details 875bd481603SBarry Smith 876bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 877bd481603SBarry Smith 8781620fd73SBarry Smith This is a MACRO not a function because it closes the { started in MatPreallocateInitialize(). 8791620fd73SBarry Smith 880bd481603SBarry Smith Concepts: preallocation^Matrix 881bd481603SBarry Smith 882bd481603SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(), 883bd481603SBarry Smith MatPreallocateSymmetricInitialize(), MatPreallocateSymmetricSetLocal() 884bd481603SBarry Smith M*/ 885ef66eb69SBarry Smith #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree(dnz);CHKERRQ(_4_ierr);} 8867c922b88SBarry Smith 887bd481603SBarry Smith 888bd481603SBarry Smith 8897b80b807SBarry Smith /* Routines unique to particular data structures */ 890be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellGetContext(Mat,void **); 891ba966639SSatish Balay PetscPolymorphicFunction(MatShellGetContext,(Mat A),(A,&t),void*,t) 892698d4c6aSKris Buschelman 893be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInodeAdjustForInodes(Mat,IS*,IS*); 894be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *); 895698d4c6aSKris Buschelman 896be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetColumnIndices(Mat,PetscInt[]); 897be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetColumnIndices(Mat,PetscInt[]); 898be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 899c75a6043SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 900c75a6043SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 9017b80b807SBarry Smith 902a96a251dSBarry Smith #define MAT_SKIP_ALLOCATION -4 903a96a251dSBarry Smith 904be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]); 905ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz)) 906be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]); 907ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqSBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz)) 908be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]); 909ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqAIJSetPreallocation,(Mat A,const PetscInt nnz[]),(A,0,nnz)) 910be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDensePreallocation(Mat,PetscScalar[]); 911273d9f13SBarry Smith 912be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 913ba966639SSatish Balay PetscPolymorphicSubroutine(MatMPIBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[],const PetscInt onz[]),(A,bs,0,nnz,0,onz)) 914be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 915be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 91653803f26SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []); 917be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]); 918be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 919be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]); 920be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIDenseSetPreallocation(Mat,PetscScalar[]); 921284134d9SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation(Mat,PetscScalar[]); 922be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIRowbsSetPreallocation(Mat,PetscInt,const PetscInt[]); 923be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,PetscInt*[]); 924be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,PetscInt*[]); 925be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAdicSetLocalFunction(Mat,void (*)(void)); 926273d9f13SBarry Smith 927be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetLDA(Mat,PetscInt); 928ab92ecdeSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDenseGetLocalMatrix(Mat,Mat*); 9291b807ce4Svictorle 930be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStoreValues(Mat); 931be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRetrieveValues(Mat); 9322e8a6d31SBarry Smith 933be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDAADSetCtx(Mat,void*); 9343a7fca6bSBarry Smith 9357b80b807SBarry Smith /* 9367b80b807SBarry Smith These routines are not usually accessed directly, rather solving is 93794b7f48cSBarry Smith done through the KSP and PC interfaces. 9387b80b807SBarry Smith */ 9397b80b807SBarry Smith 940d9274352SBarry Smith /*E 941d9274352SBarry Smith MatOrderingType - String with the name of a PETSc matrix ordering or the creation function 942d9274352SBarry Smith with an optional dynamic library name, for example 943d9274352SBarry Smith http://www.mcs.anl.gov/petsc/lib.a:orderingcreate() 944d9274352SBarry Smith 945d9274352SBarry Smith Level: beginner 946d9274352SBarry Smith 947e5a9bf91SBarry Smith Cannot use const because the PC objects manipulate the string 948e5a9bf91SBarry Smith 949d9274352SBarry Smith .seealso: MatGetOrdering() 950d9274352SBarry Smith E*/ 95149773a63SBarry Smith #define MatOrderingType char* 952b12f92e5SBarry Smith #define MATORDERING_NATURAL "natural" 953b12f92e5SBarry Smith #define MATORDERING_ND "nd" 954b12f92e5SBarry Smith #define MATORDERING_1WD "1wd" 955b12f92e5SBarry Smith #define MATORDERING_RCM "rcm" 956b12f92e5SBarry Smith #define MATORDERING_QMD "qmd" 957b12f92e5SBarry Smith #define MATORDERING_ROWLENGTH "rowlength" 95862152c8bSBarry Smith #define MATORDERING_DSC_ND "dsc_nd" 95962152c8bSBarry Smith #define MATORDERING_DSC_MMD "dsc_mmd" 96062152c8bSBarry Smith #define MATORDERING_DSC_MDF "dsc_mdf" 961c06d978dSMatthew Knepley #define MATORDERING_CONSTRAINED "constrained" 962c06d978dSMatthew Knepley #define MATORDERING_IDENTITY "identity" 963c06d978dSMatthew Knepley #define MATORDERING_REVERSE "reverse" 964b12f92e5SBarry Smith 965be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOrdering(Mat,const MatOrderingType,IS*,IS*); 966be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOrderingList(PetscFList *list); 967be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,const MatOrderingType,IS*,IS*)); 96830de9b25SBarry Smith 96930de9b25SBarry Smith /*MC 97030de9b25SBarry Smith MatOrderingRegisterDynamic - Adds a new sparse matrix ordering to the 97130de9b25SBarry Smith matrix package. 97230de9b25SBarry Smith 97330de9b25SBarry Smith Synopsis: 974c1ac3661SBarry Smith PetscErrorCode MatOrderingRegisterDynamic(char *name_ordering,char *path,char *name_create,PetscErrorCode (*routine_create)(MatOrdering)) 97530de9b25SBarry Smith 97630de9b25SBarry Smith Not Collective 97730de9b25SBarry Smith 97830de9b25SBarry Smith Input Parameters: 97930de9b25SBarry Smith + sname - name of ordering (for example MATORDERING_ND) 98030de9b25SBarry Smith . path - location of library where creation routine is 98130de9b25SBarry Smith . name - name of function that creates the ordering type,a string 98230de9b25SBarry Smith - function - function pointer that creates the ordering 98330de9b25SBarry Smith 98430de9b25SBarry Smith Level: developer 98530de9b25SBarry Smith 98630de9b25SBarry Smith If dynamic libraries are used, then the fourth input argument (function) 98730de9b25SBarry Smith is ignored. 98830de9b25SBarry Smith 98930de9b25SBarry Smith Sample usage: 99030de9b25SBarry Smith .vb 99130de9b25SBarry Smith MatOrderingRegisterDynamic("my_order",/home/username/my_lib/lib/libO/solaris/mylib.a, 99230de9b25SBarry Smith "MyOrder",MyOrder); 99330de9b25SBarry Smith .ve 99430de9b25SBarry Smith 99530de9b25SBarry Smith Then, your partitioner can be chosen with the procedural interface via 99630de9b25SBarry Smith $ MatOrderingSetType(part,"my_order) 99730de9b25SBarry Smith or at runtime via the option 9982401956bSBarry Smith $ -pc_factor_mat_ordering_type my_order 99930de9b25SBarry Smith 1000ab901514SBarry Smith ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values. 100130de9b25SBarry Smith 100230de9b25SBarry Smith .keywords: matrix, ordering, register 100330de9b25SBarry Smith 100430de9b25SBarry Smith .seealso: MatOrderingRegisterDestroy(), MatOrderingRegisterAll() 100530de9b25SBarry Smith M*/ 1006aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 1007f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,0) 1008b12f92e5SBarry Smith #else 1009f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,d) 1010b12f92e5SBarry Smith #endif 101130de9b25SBarry Smith 1012be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegisterDestroy(void); 1013be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegisterAll(const char[]); 10142bad1931SBarry Smith extern PetscTruth MatOrderingRegisterAllCalled; 1015b0a32e0cSBarry Smith extern PetscFList MatOrderingList; 1016d4fbbf0eSBarry Smith 1017be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS); 1018a2ce50c7SBarry Smith 1019d91e6319SBarry Smith /*S 10202401956bSBarry Smith MatFactorInfo - Data passed into the matrix factorization routines 1021b00f7748SHong Zhang 102261cad860SBarry Smith In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE, that is use 102361cad860SBarry Smith $ MatFactorInfo info(MAT_FACTORINFO_SIZE) 1024b00f7748SHong Zhang 102515e8a5b3SHong Zhang Notes: These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC. 1026b00f7748SHong Zhang 102761cad860SBarry Smith You can use MatFactorInfoInitialize() to set default values. 102861cad860SBarry Smith 1029b00f7748SHong Zhang Level: developer 1030b00f7748SHong Zhang 1031d7d82daaSBarry Smith .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(), 1032d7d82daaSBarry Smith MatFactorInfoInitialize() 1033b00f7748SHong Zhang 1034b00f7748SHong Zhang S*/ 1035b00f7748SHong Zhang typedef struct { 10360a29876aSHong Zhang PetscReal shiftnz; /* scaling of identity added to matrix to prevent zero pivots */ 1037fbf22428SSatish Balay PetscReal shiftpd; /* if true, shift until positive pivots */ 10382cea7109SBarry Smith PetscReal shift_fraction; /* record shift fraction taken */ 103915e8a5b3SHong Zhang PetscReal diagonal_fill; /* force diagonal to fill in if initially not filled */ 104015e8a5b3SHong Zhang PetscReal dt; /* drop tolerance */ 1041b00f7748SHong Zhang PetscReal dtcol; /* tolerance for pivoting */ 104215e8a5b3SHong Zhang PetscReal dtcount; /* maximum nonzeros to be allowed per row */ 104367e28bfeSBarry Smith PetscReal fill; /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */ 1044348344bbSBarry Smith PetscReal levels; /* ICC/ILU(levels) */ 1045bcd9e38bSBarry Smith PetscReal pivotinblocks; /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0 1046bcd9e38bSBarry Smith factorization may be faster if do not pivot */ 104762bba022SBarry Smith PetscReal shiftinblocks; /* if block in block factorization has zero pivot then shift diagonal until non-singular */ 104815e8a5b3SHong Zhang PetscReal zeropivot; /* pivot is called zero if less than this */ 104915e8a5b3SHong Zhang } MatFactorInfo; 1050ffa6d0a5SLois Curfman McInnes 1051be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFactorInfoInitialize(MatFactorInfo*); 1052be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactor(Mat,IS,MatFactorInfo*); 1053be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactorSymbolic(Mat,IS,MatFactorInfo*,Mat*); 1054be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactorNumeric(Mat,MatFactorInfo*,Mat*); 1055be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactor(Mat,IS,IS,MatFactorInfo*); 1056be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUFactor(Mat,IS,IS,MatFactorInfo*); 1057be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactorSymbolic(Mat,IS,IS,MatFactorInfo*,Mat*); 1058be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUFactorSymbolic(Mat,IS,IS,MatFactorInfo*,Mat*); 1059be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatICCFactorSymbolic(Mat,IS,MatFactorInfo*,Mat*); 1060be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatICCFactor(Mat,IS,MatFactorInfo*); 1061be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactorNumeric(Mat,MatFactorInfo*,Mat*); 1062be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUDTFactor(Mat,IS,IS,MatFactorInfo*,Mat *); 1063be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*); 1064be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolve(Mat,Vec,Vec); 1065be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatForwardSolve(Mat,Vec,Vec); 1066be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatBackwardSolve(Mat,Vec,Vec); 1067be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveAdd(Mat,Vec,Vec,Vec); 1068be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveTranspose(Mat,Vec,Vec); 1069be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveTransposeAdd(Mat,Vec,Vec,Vec); 1070be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolves(Mat,Vecs,Vecs); 10718ed539a5SBarry Smith 1072be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUnfactored(Mat); 1073bb5a7306SBarry Smith 1074d91e6319SBarry Smith /*E 1075d91e6319SBarry Smith MatSORType - What type of (S)SOR to perform 1076bb1eb677SSatish Balay 1077d91e6319SBarry Smith Level: beginner 1078d91e6319SBarry Smith 1079d9274352SBarry Smith May be bitwise ORd together 1080d9274352SBarry Smith 1081d91e6319SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 1082d91e6319SBarry Smith 10834e7234bfSBarry Smith MatSORType may be bitwise ORd together, so do not change the numbers 10844e7234bfSBarry Smith 1085a1d92eedSBarry Smith .seealso: MatRelax(), MatPBRelax() 1086d91e6319SBarry Smith E*/ 1087ee50ffe9SBarry Smith typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3, 1088ee50ffe9SBarry Smith SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8, 1089ee50ffe9SBarry Smith SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16, 109084cb2905SBarry Smith SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType; 1091be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRelax(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec); 1092be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPBRelax(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec); 10938ed539a5SBarry Smith 1094d4fbbf0eSBarry Smith /* 1095639f9d9dSBarry Smith These routines are for efficiently computing Jacobians via finite differences. 1096639f9d9dSBarry Smith */ 1097b12f92e5SBarry Smith 1098d9274352SBarry Smith /*E 1099d9274352SBarry Smith MatColoringType - String with the name of a PETSc matrix coloring or the creation function 1100d9274352SBarry Smith with an optional dynamic library name, for example 1101d9274352SBarry Smith http://www.mcs.anl.gov/petsc/lib.a:coloringcreate() 1102d9274352SBarry Smith 1103d9274352SBarry Smith Level: beginner 1104d9274352SBarry Smith 1105d9274352SBarry Smith .seealso: MatGetColoring() 1106d9274352SBarry Smith E*/ 1107e5a9bf91SBarry Smith #define MatColoringType const char* 1108b12f92e5SBarry Smith #define MATCOLORING_NATURAL "natural" 1109b12f92e5SBarry Smith #define MATCOLORING_SL "sl" 1110b12f92e5SBarry Smith #define MATCOLORING_LF "lf" 1111b12f92e5SBarry Smith #define MATCOLORING_ID "id" 1112b12f92e5SBarry Smith 11132e90c967SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColoring(Mat,MatColoringType,ISColoring*); 11142e90c967SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,MatColoringType,ISColoring *)); 111530de9b25SBarry Smith 111630de9b25SBarry Smith /*MC 111730de9b25SBarry Smith MatColoringRegisterDynamic - Adds a new sparse matrix coloring to the 111830de9b25SBarry Smith matrix package. 111930de9b25SBarry Smith 112030de9b25SBarry Smith Synopsis: 1121c1ac3661SBarry Smith PetscErrorCode MatColoringRegisterDynamic(char *name_coloring,char *path,char *name_create,PetscErrorCode (*routine_create)(MatColoring)) 112230de9b25SBarry Smith 112330de9b25SBarry Smith Not Collective 112430de9b25SBarry Smith 112530de9b25SBarry Smith Input Parameters: 112630de9b25SBarry Smith + sname - name of Coloring (for example MATCOLORING_SL) 112730de9b25SBarry Smith . path - location of library where creation routine is 112830de9b25SBarry Smith . name - name of function that creates the Coloring type, a string 112930de9b25SBarry Smith - function - function pointer that creates the coloring 113030de9b25SBarry Smith 113130de9b25SBarry Smith Level: developer 113230de9b25SBarry Smith 113330de9b25SBarry Smith If dynamic libraries are used, then the fourth input argument (function) 113430de9b25SBarry Smith is ignored. 113530de9b25SBarry Smith 113630de9b25SBarry Smith Sample usage: 113730de9b25SBarry Smith .vb 113830de9b25SBarry Smith MatColoringRegisterDynamic("my_color",/home/username/my_lib/lib/libO/solaris/mylib.a, 113930de9b25SBarry Smith "MyColor",MyColor); 114030de9b25SBarry Smith .ve 114130de9b25SBarry Smith 114230de9b25SBarry Smith Then, your partitioner can be chosen with the procedural interface via 114330de9b25SBarry Smith $ MatColoringSetType(part,"my_color") 114430de9b25SBarry Smith or at runtime via the option 114530de9b25SBarry Smith $ -mat_coloring_type my_color 114630de9b25SBarry Smith 1147ab901514SBarry Smith $PETSC_ARCH occuring in pathname will be replaced with appropriate values. 114830de9b25SBarry Smith 114930de9b25SBarry Smith .keywords: matrix, Coloring, register 115030de9b25SBarry Smith 115130de9b25SBarry Smith .seealso: MatColoringRegisterDestroy(), MatColoringRegisterAll() 115230de9b25SBarry Smith M*/ 1153aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 1154f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,0) 1155b12f92e5SBarry Smith #else 1156f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,d) 1157b12f92e5SBarry Smith #endif 115830de9b25SBarry Smith 11592bad1931SBarry Smith extern PetscTruth MatColoringRegisterAllCalled; 1160f1144a30SSatish Balay 1161be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegisterAll(const char[]); 1162be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegisterDestroy(void); 1163be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*); 1164639f9d9dSBarry Smith 1165d9274352SBarry Smith /*S 1166d9274352SBarry Smith MatFDColoring - Object for computing a sparse Jacobian via finite differences 1167d9274352SBarry Smith and coloring 1168639f9d9dSBarry Smith 1169d9274352SBarry Smith Level: beginner 1170d9274352SBarry Smith 1171d9274352SBarry Smith Concepts: coloring, sparse Jacobian, finite differences 1172d9274352SBarry Smith 1173d9274352SBarry Smith .seealso: MatFDColoringCreate() 1174d9274352SBarry Smith S*/ 1175e2a1c21fSSatish Balay typedef struct _p_MatFDColoring* MatFDColoring; 1176639f9d9dSBarry Smith 1177be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringCreate(Mat,ISColoring,MatFDColoring *); 1178be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringDestroy(MatFDColoring); 1179be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringView(MatFDColoring,PetscViewer); 1180be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*); 11814a9d489dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**); 1182be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal); 1183be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFrequency(MatFDColoring,PetscInt); 1184be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetFrequency(MatFDColoring,PetscInt*); 1185be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFromOptions(MatFDColoring); 1186be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApply(Mat,MatFDColoring,Vec,MatStructure*,void *); 1187be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApplyTS(Mat,MatFDColoring,PetscReal,Vec,MatStructure*,void *); 1188be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetRecompute(MatFDColoring); 1189be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetF(MatFDColoring,Vec); 1190be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,PetscInt*[]); 1191639f9d9dSBarry Smith /* 11920752156aSBarry Smith These routines are for partitioning matrices: currently used only 11933eda8832SBarry Smith for adjacency matrix, MatCreateMPIAdj(). 11940752156aSBarry Smith */ 1195ca161407SBarry Smith 1196d9274352SBarry Smith /*S 1197d9274352SBarry Smith MatPartitioning - Object for managing the partitioning of a matrix or graph 1198d9274352SBarry Smith 1199d9274352SBarry Smith Level: beginner 1200d9274352SBarry Smith 1201d9274352SBarry Smith Concepts: partitioning 1202d9274352SBarry Smith 1203743a1026Svictorle .seealso: MatPartitioningCreate(), MatPartitioningType 1204d9274352SBarry Smith S*/ 120591e9ee9fSBarry Smith typedef struct _p_MatPartitioning* MatPartitioning; 1206d9274352SBarry Smith 1207d9274352SBarry Smith /*E 12085bc6e7ccSvictorle MatPartitioningType - String with the name of a PETSc matrix partitioning or the creation function 1209d9274352SBarry Smith with an optional dynamic library name, for example 1210d9274352SBarry Smith http://www.mcs.anl.gov/petsc/lib.a:partitioningcreate() 1211d9274352SBarry Smith 1212d9274352SBarry Smith Level: beginner 1213d9274352SBarry Smith 1214b547a13bSvictorle .seealso: MatPartitioningCreate(), MatPartitioning 1215d9274352SBarry Smith E*/ 1216e5a9bf91SBarry Smith #define MatPartitioningType const char* 12178ba1e511SMatthew Knepley #define MAT_PARTITIONING_CURRENT "current" 1218c2ec2a73SHong Zhang #define MAT_PARTITIONING_SQUARE "square" 12198ba1e511SMatthew Knepley #define MAT_PARTITIONING_PARMETIS "parmetis" 122071306c60Sjroman #define MAT_PARTITIONING_CHACO "chaco" 122171306c60Sjroman #define MAT_PARTITIONING_JOSTLE "jostle" 122271306c60Sjroman #define MAT_PARTITIONING_PARTY "party" 122371306c60Sjroman #define MAT_PARTITIONING_SCOTCH "scotch" 122471306c60Sjroman 1225ca161407SBarry Smith 1226be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningCreate(MPI_Comm,MatPartitioning*); 12272e90c967SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetType(MatPartitioning,MatPartitioningType); 1228be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetNParts(MatPartitioning,PetscInt); 1229be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetAdjacency(MatPartitioning,Mat); 1230be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]); 1231be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []); 1232be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningApply(MatPartitioning,IS*); 1233be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningDestroy(MatPartitioning); 12342aabb6bbSBarry Smith 1235be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatPartitioning)); 123630de9b25SBarry Smith 123730de9b25SBarry Smith /*MC 123830de9b25SBarry Smith MatPartitioningRegisterDynamic - Adds a new sparse matrix partitioning to the 123930de9b25SBarry Smith matrix package. 124030de9b25SBarry Smith 124130de9b25SBarry Smith Synopsis: 1242c1ac3661SBarry Smith PetscErrorCode MatPartitioningRegisterDynamic(char *name_partitioning,char *path,char *name_create,PetscErrorCode (*routine_create)(MatPartitioning)) 124330de9b25SBarry Smith 124430de9b25SBarry Smith Not Collective 124530de9b25SBarry Smith 124630de9b25SBarry Smith Input Parameters: 124730de9b25SBarry Smith + sname - name of partitioning (for example MAT_PARTITIONING_CURRENT) or parmetis 124830de9b25SBarry Smith . path - location of library where creation routine is 124930de9b25SBarry Smith . name - name of function that creates the partitioning type, a string 125030de9b25SBarry Smith - function - function pointer that creates the partitioning type 125130de9b25SBarry Smith 125230de9b25SBarry Smith Level: developer 125330de9b25SBarry Smith 125430de9b25SBarry Smith If dynamic libraries are used, then the fourth input argument (function) 125530de9b25SBarry Smith is ignored. 125630de9b25SBarry Smith 125730de9b25SBarry Smith Sample usage: 125830de9b25SBarry Smith .vb 125930de9b25SBarry Smith MatPartitioningRegisterDynamic("my_part",/home/username/my_lib/lib/libO/solaris/mylib.a, 126030de9b25SBarry Smith "MyPartCreate",MyPartCreate); 126130de9b25SBarry Smith .ve 126230de9b25SBarry Smith 126330de9b25SBarry Smith Then, your partitioner can be chosen with the procedural interface via 126430de9b25SBarry Smith $ MatPartitioningSetType(part,"my_part") 126530de9b25SBarry Smith or at runtime via the option 126630de9b25SBarry Smith $ -mat_partitioning_type my_part 126730de9b25SBarry Smith 1268ab901514SBarry Smith $PETSC_ARCH occuring in pathname will be replaced with appropriate values. 126930de9b25SBarry Smith 127030de9b25SBarry Smith .keywords: matrix, partitioning, register 127130de9b25SBarry Smith 127230de9b25SBarry Smith .seealso: MatPartitioningRegisterDestroy(), MatPartitioningRegisterAll() 127330de9b25SBarry Smith M*/ 1274aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 1275f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,0) 12762aabb6bbSBarry Smith #else 1277f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,d) 12782aabb6bbSBarry Smith #endif 12792aabb6bbSBarry Smith 12802bad1931SBarry Smith extern PetscTruth MatPartitioningRegisterAllCalled; 1281f1144a30SSatish Balay 1282be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegisterAll(const char[]); 1283be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegisterDestroy(void); 12842bad1931SBarry Smith 1285be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningView(MatPartitioning,PetscViewer); 1286be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetFromOptions(MatPartitioning); 1287be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningGetType(MatPartitioning,MatPartitioningType*); 1288ca161407SBarry Smith 1289be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningParmetisSetCoarseSequential(MatPartitioning); 12900752156aSBarry Smith 1291be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningJostleSetCoarseLevel(MatPartitioning,PetscReal); 1292be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningJostleSetCoarseSequential(MatPartitioning); 129371306c60Sjroman 1294954b3ec6SBarry Smith typedef enum {MP_CHACO_MULTILEVEL_KL,MP_CHACO_SPECTRAL,MP_CHACO_LINEAR,MP_CHACO_RANDOM, MP_CHACO_SCATTERED} MPChacoGlobalType; 1295be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetGlobal(MatPartitioning, MPChacoGlobalType); 129671306c60Sjroman typedef enum { MP_CHACO_KERNIGHAN_LIN, MP_CHACO_NONE } MPChacoLocalType; 1297be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetLocal(MatPartitioning, MPChacoLocalType); 1298be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal); 129971306c60Sjroman typedef enum { MP_CHACO_LANCZOS, MP_CHACO_RQI_SYMMLQ } MPChacoEigenType; 1300be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType); 1301be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenTol(MatPartitioning, PetscReal); 1302be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenNumber(MatPartitioning, PetscInt); 130371306c60Sjroman 130471306c60Sjroman #define MP_PARTY_OPT "opt" 130571306c60Sjroman #define MP_PARTY_LIN "lin" 130671306c60Sjroman #define MP_PARTY_SCA "sca" 130771306c60Sjroman #define MP_PARTY_RAN "ran" 130871306c60Sjroman #define MP_PARTY_GBF "gbf" 130971306c60Sjroman #define MP_PARTY_GCF "gcf" 131071306c60Sjroman #define MP_PARTY_BUB "bub" 131171306c60Sjroman #define MP_PARTY_DEF "def" 1312be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetGlobal(MatPartitioning, const char*); 131371306c60Sjroman #define MP_PARTY_HELPFUL_SETS "hs" 131471306c60Sjroman #define MP_PARTY_KERNIGHAN_LIN "kl" 131571306c60Sjroman #define MP_PARTY_NONE "no" 1316be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetLocal(MatPartitioning, const char*); 1317be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal); 1318be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetBipart(MatPartitioning,PetscTruth); 1319be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscTruth); 132071306c60Sjroman 132171306c60Sjroman typedef enum { MP_SCOTCH_GREEDY, MP_SCOTCH_GPS, MP_SCOTCH_GR_GPS } MPScotchGlobalType; 1322be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetArch(MatPartitioning,const char*); 1323be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetMultilevel(MatPartitioning); 1324be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetGlobal(MatPartitioning,MPScotchGlobalType); 1325be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetCoarseLevel(MatPartitioning,PetscReal); 1326be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetHostList(MatPartitioning,const char*); 132771306c60Sjroman typedef enum { MP_SCOTCH_KERNIGHAN_LIN, MP_SCOTCH_NONE } MPScotchLocalType; 1328be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetLocal(MatPartitioning,MPScotchLocalType); 1329be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetMapping(MatPartitioning); 1330be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetStrategy(MatPartitioning,char*); 133171306c60Sjroman 13320752156aSBarry Smith /* 13330a835dfdSSatish Balay If you add entries here you must also add them to finclude/petscmat.h 1334d4fbbf0eSBarry Smith */ 13351c1c02c0SLois Curfman McInnes typedef enum { MATOP_SET_VALUES=0, 13361c1c02c0SLois Curfman McInnes MATOP_GET_ROW=1, 13371c1c02c0SLois Curfman McInnes MATOP_RESTORE_ROW=2, 13381c1c02c0SLois Curfman McInnes MATOP_MULT=3, 13391c1c02c0SLois Curfman McInnes MATOP_MULT_ADD=4, 13407c922b88SBarry Smith MATOP_MULT_TRANSPOSE=5, 13417c922b88SBarry Smith MATOP_MULT_TRANSPOSE_ADD=6, 13421c1c02c0SLois Curfman McInnes MATOP_SOLVE=7, 13431c1c02c0SLois Curfman McInnes MATOP_SOLVE_ADD=8, 13447c922b88SBarry Smith MATOP_SOLVE_TRANSPOSE=9, 13457c922b88SBarry Smith MATOP_SOLVE_TRANSPOSE_ADD=10, 13461c1c02c0SLois Curfman McInnes MATOP_LUFACTOR=11, 13471c1c02c0SLois Curfman McInnes MATOP_CHOLESKYFACTOR=12, 13481c1c02c0SLois Curfman McInnes MATOP_RELAX=13, 13491c1c02c0SLois Curfman McInnes MATOP_TRANSPOSE=14, 13501c1c02c0SLois Curfman McInnes MATOP_GETINFO=15, 13511c1c02c0SLois Curfman McInnes MATOP_EQUAL=16, 13521c1c02c0SLois Curfman McInnes MATOP_GET_DIAGONAL=17, 13531c1c02c0SLois Curfman McInnes MATOP_DIAGONAL_SCALE=18, 13541c1c02c0SLois Curfman McInnes MATOP_NORM=19, 13551c1c02c0SLois Curfman McInnes MATOP_ASSEMBLY_BEGIN=20, 13561c1c02c0SLois Curfman McInnes MATOP_ASSEMBLY_END=21, 13571c1c02c0SLois Curfman McInnes MATOP_COMPRESS=22, 13581c1c02c0SLois Curfman McInnes MATOP_SET_OPTION=23, 13591c1c02c0SLois Curfman McInnes MATOP_ZERO_ENTRIES=24, 13601c1c02c0SLois Curfman McInnes MATOP_ZERO_ROWS=25, 13611c1c02c0SLois Curfman McInnes MATOP_LUFACTOR_SYMBOLIC=26, 13621c1c02c0SLois Curfman McInnes MATOP_LUFACTOR_NUMERIC=27, 13631c1c02c0SLois Curfman McInnes MATOP_CHOLESKY_FACTOR_SYMBOLIC=28, 13641c1c02c0SLois Curfman McInnes MATOP_CHOLESKY_FACTOR_NUMERIC=29, 1365d643ce63SMatthew Knepley MATOP_SETUP_PREALLOCATION=30, 1366d643ce63SMatthew Knepley MATOP_ILUFACTOR_SYMBOLIC=31, 1367d643ce63SMatthew Knepley MATOP_ICCFACTOR_SYMBOLIC=32, 1368d643ce63SMatthew Knepley MATOP_GET_ARRAY=33, 1369d643ce63SMatthew Knepley MATOP_RESTORE_ARRAY=34, 1370e26f1c3cSBarry Smith MATOP_DUPLICATE=35, 1371d643ce63SMatthew Knepley MATOP_FORWARD_SOLVE=36, 1372d643ce63SMatthew Knepley MATOP_BACKWARD_SOLVE=37, 1373d643ce63SMatthew Knepley MATOP_ILUFACTOR=38, 1374d643ce63SMatthew Knepley MATOP_ICCFACTOR=39, 1375d643ce63SMatthew Knepley MATOP_AXPY=40, 1376d643ce63SMatthew Knepley MATOP_GET_SUBMATRICES=41, 1377d643ce63SMatthew Knepley MATOP_INCREASE_OVERLAP=42, 1378d643ce63SMatthew Knepley MATOP_GET_VALUES=43, 1379d643ce63SMatthew Knepley MATOP_COPY=44, 1380d643ce63SMatthew Knepley MATOP_PRINT_HELP=45, 1381d643ce63SMatthew Knepley MATOP_SCALE=46, 1382d643ce63SMatthew Knepley MATOP_SHIFT=47, 1383d643ce63SMatthew Knepley MATOP_DIAGONAL_SHIFT=48, 1384d643ce63SMatthew Knepley MATOP_ILUDT_FACTOR=49, 1385d643ce63SMatthew Knepley MATOP_GET_BLOCK_SIZE=50, 1386d643ce63SMatthew Knepley MATOP_GET_ROW_IJ=51, 1387d643ce63SMatthew Knepley MATOP_RESTORE_ROW_IJ=52, 1388d643ce63SMatthew Knepley MATOP_GET_COLUMN_IJ=53, 1389d643ce63SMatthew Knepley MATOP_RESTORE_COLUMN_IJ=54, 1390d643ce63SMatthew Knepley MATOP_FDCOLORING_CREATE=55, 1391d643ce63SMatthew Knepley MATOP_COLORING_PATCH=56, 1392d643ce63SMatthew Knepley MATOP_SET_UNFACTORED=57, 1393d643ce63SMatthew Knepley MATOP_PERMUTE=58, 1394d643ce63SMatthew Knepley MATOP_SET_VALUES_BLOCKED=59, 1395d643ce63SMatthew Knepley MATOP_GET_SUBMATRIX=60, 1396d643ce63SMatthew Knepley MATOP_DESTROY=61, 1397d643ce63SMatthew Knepley MATOP_VIEW=62, 1398d643ce63SMatthew Knepley MATOP_GET_MAPS=63, 1399d643ce63SMatthew Knepley MATOP_USE_SCALED_FORM=64, 1400d643ce63SMatthew Knepley MATOP_SCALE_SYSTEM=65, 1401d643ce63SMatthew Knepley MATOP_UNSCALE_SYSTEM=66, 1402ba16a8cbSKris Buschelman MATOP_SET_LOCAL_TO_GLOBAL_MAP=67, 1403d643ce63SMatthew Knepley MATOP_SET_VALUES_LOCAL=68, 1404d643ce63SMatthew Knepley MATOP_ZERO_ROWS_LOCAL=69, 1405d643ce63SMatthew Knepley MATOP_GET_ROW_MAX=70, 1406d643ce63SMatthew Knepley MATOP_CONVERT=71, 1407d643ce63SMatthew Knepley MATOP_SET_COLORING=72, 1408d643ce63SMatthew Knepley MATOP_SET_VALUES_ADIC=73, 1409d643ce63SMatthew Knepley MATOP_SET_VALUES_ADIFOR=74, 1410d643ce63SMatthew Knepley MATOP_FD_COLORING_APPLY=75, 1411d643ce63SMatthew Knepley MATOP_SET_FROM_OPTIONS=76, 1412ba16a8cbSKris Buschelman MATOP_MULT_CON=77, 1413ba16a8cbSKris Buschelman MATOP_MULT_TRANSPOSE_CON=78, 1414ba16a8cbSKris Buschelman MATOP_ILU_FACTOR_SYMBOLIC_CON=79, 1415d643ce63SMatthew Knepley MATOP_PERMUTE_SPARSIFY=80, 1416d643ce63SMatthew Knepley MATOP_MULT_MULTIPLE=81, 141741acf15aSKris Buschelman MATOP_SOLVE_MULTIPLE=82, 141841acf15aSKris Buschelman MATOP_GET_INERTIA=83, 141941acf15aSKris Buschelman MATOP_LOAD=84, 142041acf15aSKris Buschelman MATOP_IS_SYMMETRIC=85, 142141acf15aSKris Buschelman MATOP_IS_HERMITIAN=86, 142241acf15aSKris Buschelman MATOP_IS_STRUCTURALLY_SYMMETRIC=87, 142341acf15aSKris Buschelman MATOP_PB_RELAX=88, 142441acf15aSKris Buschelman MATOP_GET_VECS=89, 142541acf15aSKris Buschelman MATOP_MAT_MULT=90, 142641acf15aSKris Buschelman MATOP_MAT_MULT_SYMBOLIC=91, 142741acf15aSKris Buschelman MATOP_MAT_MULT_NUMERIC=92, 142841acf15aSKris Buschelman MATOP_PTAP=93, 142941acf15aSKris Buschelman MATOP_PTAP_SYMBOLIC=94, 1430bc011b1eSHong Zhang MATOP_PTAP_NUMERIC=95, 1431bc011b1eSHong Zhang MATOP_MAT_MULTTRANSPOSE=96, 1432bc011b1eSHong Zhang MATOP_MAT_MULTTRANSPOSE_SYMBOLIC=97, 14337ba1a0bfSKris Buschelman MATOP_MAT_MULTTRANSPOSE_NUMERIC=98, 14347ba1a0bfSKris Buschelman MATOP_PTAP_SYMBOLIC_SEQAIJ=99, 14357ba1a0bfSKris Buschelman MATOP_PTAP_NUMERIC_SEQAIJ=100, 14367ba1a0bfSKris Buschelman MATOP_PTAP_SYMBOLIC_MPIAIJ=101, 143787d4246cSBarry Smith MATOP_PTAP_NUMERIC_MPIAIJ=102, 1438f5edf698SHong Zhang MATOP_SET_VALUES_ROW = 105, 1439d5c63f83SSatish Balay MATOP_GET_ROW_UTRIANGULAR=108, 14402bebee5dSHong Zhang MATOP_RESTORE_ROW_UTRIANGULAR=109, 144169db28dcSHong Zhang MATOP_MATSOLVE=110, 1442829201f2SHong Zhang MATOP_GET_REDUNDANTMATRIX=111, 1443829201f2SHong Zhang MATOP_MATGETSEQNONZEROSTRUCTURE=115 1444fae171e0SBarry Smith } MatOperation; 1445be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatHasOperation(Mat,MatOperation,PetscTruth*); 1446be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellSetOperation(Mat,MatOperation,void(*)(void)); 1447be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellGetOperation(Mat,MatOperation,void(**)(void)); 1448be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellSetContext(Mat,void*); 1449112a2221SBarry Smith 145090ace30eSBarry Smith /* 145190ace30eSBarry Smith Codes for matrices stored on disk. By default they are 145290ace30eSBarry Smith stored in a universal format. By changing the format with 1453fb9695e5SSatish Balay PetscViewerSetFormat(viewer,PETSC_VIEWER_BINARY_NATIVE); the matrices will 145490ace30eSBarry Smith be stored in a way natural for the matrix, for example dense matrices 145590ace30eSBarry Smith would be stored as dense. Matrices stored this way may only be 145690ace30eSBarry Smith read into matrices of the same time. 145790ace30eSBarry Smith */ 145890ace30eSBarry Smith #define MATRIX_BINARY_FORMAT_DENSE -1 145990ace30eSBarry Smith 1460be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetHashTableFactor(Mat,PetscReal); 1461be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatISGetLocalMat(Mat,Mat*); 14621f4e1ec7SBarry Smith 1463d9274352SBarry Smith /*S 1464d9274352SBarry Smith MatNullSpace - Object that removes a null space from a vector, i.e. 1465d9274352SBarry Smith orthogonalizes the vector to a subsapce 1466d9274352SBarry Smith 1467f7a9e4ceSBarry Smith Level: advanced 1468d9274352SBarry Smith 1469d9274352SBarry Smith Concepts: matrix; linear operator, null space 1470d9274352SBarry Smith 14716e1639daSBarry Smith Users manual sections: 14726e1639daSBarry Smith . sec_singular 14736e1639daSBarry Smith 1474d9274352SBarry Smith .seealso: MatNullSpaceCreate() 1475d9274352SBarry Smith S*/ 147674637425SBarry Smith typedef struct _p_MatNullSpace* MatNullSpace; 1477d9274352SBarry Smith 1478be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceCreate(MPI_Comm,PetscTruth,PetscInt,const Vec[],MatNullSpace*); 1479281d1b2eSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(Vec,void*),void*); 1480be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceDestroy(MatNullSpace); 1481be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceRemove(MatNullSpace,Vec,Vec*); 1482be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceAttach(Mat,MatNullSpace); 1483be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceTest(MatNullSpace,Mat); 148474637425SBarry Smith 1485be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatReorderingSeqSBAIJ(Mat,IS); 1486be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPISBAIJSetHashTableFactor(Mat,PetscReal); 1487be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetColumnIndices(Mat,PetscInt *); 148846129b97SKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJInvertBlockDiagonal(Mat); 14893f1d51d7SBarry Smith 1490be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMAIJ(Mat,PetscInt,Mat*); 1491be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMAIJRedimension(Mat,PetscInt,Mat*); 1492be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMAIJGetAIJ(Mat,Mat*); 1493c4f061fbSSatish Balay 1494be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatComputeExplicitOperator(Mat,Mat*); 1495b0a32e0cSBarry Smith 1496be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalScaleLocal(Mat,Vec); 149704f1ad80SBarry Smith 1498fef1beadSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*); 14993ec795f1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetBase(Mat,Vec,Vec); 15003ec795f1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*); 1501e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*)); 1502e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec)); 1503e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDAddNullSpace(Mat,MatNullSpace); 1504e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt); 1505e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDResetHHistory(Mat); 1506e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctionError(Mat,PetscReal); 1507e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetPeriod(Mat,PetscInt); 1508e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDGetH(Mat,PetscScalar *); 1509e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFromOptions(Mat); 1510e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*); 1511e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*); 1512e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDDSSetUmin(Mat,PetscReal); 1513e884886fSBarry Smith 1514e884886fSBarry Smith 1515e884886fSBarry Smith typedef struct _p_MatMFFD* MatMFFD; 1516e884886fSBarry Smith 1517e884886fSBarry Smith /*E 1518e884886fSBarry Smith MatMFFDType - algorithm used to compute the h used in computing matrix-vector products via differencing of the function 1519e884886fSBarry Smith 1520e884886fSBarry Smith Level: beginner 1521e884886fSBarry Smith 1522e884886fSBarry Smith .seealso: MatMFFDSetType(), MatMFFDRegister() 1523e884886fSBarry Smith E*/ 1524e884886fSBarry Smith #define MatMFFDType const char* 1525e884886fSBarry Smith #define MATMFFD_DS "ds" 1526e884886fSBarry Smith #define MATMFFD_WP "wp" 1527e884886fSBarry Smith 1528e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetType(Mat,MatMFFDType); 1529e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatMFFD)); 1530e884886fSBarry Smith 1531e884886fSBarry Smith /*MC 1532e884886fSBarry Smith MatMFFDRegisterDynamic - Adds a method to the MatMFFD registry. 1533e884886fSBarry Smith 1534e884886fSBarry Smith Synopsis: 1535e884886fSBarry Smith PetscErrorCode MatMFFDRegisterDynamic(char *name_solver,char *path,char *name_create,PetscErrorCode (*routine_create)(MatMFFD)) 1536e884886fSBarry Smith 1537e884886fSBarry Smith Not Collective 1538e884886fSBarry Smith 1539e884886fSBarry Smith Input Parameters: 1540e884886fSBarry Smith + name_solver - name of a new user-defined compute-h module 1541e884886fSBarry Smith . path - path (either absolute or relative) the library containing this solver 1542e884886fSBarry Smith . name_create - name of routine to create method context 1543e884886fSBarry Smith - routine_create - routine to create method context 1544e884886fSBarry Smith 1545e884886fSBarry Smith Level: developer 1546e884886fSBarry Smith 1547e884886fSBarry Smith Notes: 1548e884886fSBarry Smith MatMFFDRegisterDynamic() may be called multiple times to add several user-defined solvers. 1549e884886fSBarry Smith 1550e884886fSBarry Smith If dynamic libraries are used, then the fourth input argument (routine_create) 1551e884886fSBarry Smith is ignored. 1552e884886fSBarry Smith 1553e884886fSBarry Smith Sample usage: 1554e884886fSBarry Smith .vb 1555e884886fSBarry Smith MatMFFDRegisterDynamic("my_h",/home/username/my_lib/lib/libO/solaris/mylib.a, 1556e884886fSBarry Smith "MyHCreate",MyHCreate); 1557e884886fSBarry Smith .ve 1558e884886fSBarry Smith 1559e884886fSBarry Smith Then, your solver can be chosen with the procedural interface via 1560e884886fSBarry Smith $ MatMFFDSetType(mfctx,"my_h") 1561e884886fSBarry Smith or at runtime via the option 1562e884886fSBarry Smith $ -snes_mf_type my_h 1563e884886fSBarry Smith 1564e884886fSBarry Smith .keywords: MatMFFD, register 1565e884886fSBarry Smith 1566e884886fSBarry Smith .seealso: MatMFFDRegisterAll(), MatMFFDRegisterDestroy() 1567e884886fSBarry Smith M*/ 1568e884886fSBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 1569e884886fSBarry Smith #define MatMFFDRegisterDynamic(a,b,c,d) MatMFFDRegister(a,b,c,0) 1570e884886fSBarry Smith #else 1571e884886fSBarry Smith #define MatMFFDRegisterDynamic(a,b,c,d) MatMFFDRegister(a,b,c,d) 1572e884886fSBarry Smith #endif 1573e884886fSBarry Smith 1574e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDRegisterAll(const char[]); 1575e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDRegisterDestroy(void); 1576e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDDefaultSetUmin(Mat,PetscReal); 1577e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDWPSetComputeNormU(Mat,PetscTruth); 1578e884886fSBarry Smith 1579e884886fSBarry Smith 1580be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *); 1581be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *); 15827dbadf16SMatthew Knepley 1583e9fa29b7SSatish Balay PETSC_EXTERN_CXX_END 15842eac72dbSBarry Smith #endif 1585