12eac72dbSBarry Smith /* 22eac72dbSBarry Smith Include file for the matrix component of PETSc 32eac72dbSBarry Smith */ 40a835dfdSSatish Balay #ifndef __PETSCMAT_H 50a835dfdSSatish Balay #define __PETSCMAT_H 60a835dfdSSatish Balay #include "petscvec.h" 7e9fa29b7SSatish Balay PETSC_EXTERN_CXX_BEGIN 82eac72dbSBarry Smith 9d9274352SBarry Smith /*S 10d9274352SBarry Smith Mat - Abstract PETSc matrix object 112eac72dbSBarry Smith 12d91e6319SBarry Smith Level: beginner 13d91e6319SBarry Smith 14d9274352SBarry Smith Concepts: matrix; linear operator 15d9274352SBarry Smith 16d9274352SBarry Smith .seealso: MatCreate(), MatType, MatSetType() 17d9274352SBarry Smith S*/ 18d9274352SBarry Smith typedef struct _p_Mat* Mat; 19d9274352SBarry Smith 20d9274352SBarry Smith /*E 21d9274352SBarry Smith MatType - String with the name of a PETSc matrix or the creation function 22d9274352SBarry Smith with an optional dynamic library name, for example 23d9274352SBarry Smith http://www.mcs.anl.gov/petsc/lib.a:mymatcreate() 24d9274352SBarry Smith 25d9274352SBarry Smith Level: beginner 26d9274352SBarry Smith 27c7393fdbSBarry Smith .seealso: MatSetType(), Mat, MatSolverPackage 28d91e6319SBarry Smith E*/ 29a313700dSBarry Smith #define MatType char* 30273d9f13SBarry Smith #define MATSAME "same" 31273d9f13SBarry Smith #define MATSEQMAIJ "seqmaij" 32273d9f13SBarry Smith #define MATMPIMAIJ "mpimaij" 33209238afSKris Buschelman #define MATMAIJ "maij" 34273d9f13SBarry Smith #define MATIS "is" 35273d9f13SBarry Smith #define MATSEQAIJ "seqaij" 36273d9f13SBarry Smith #define MATMPIAIJ "mpiaij" 37209238afSKris Buschelman #define MATAIJ "aij" 38273d9f13SBarry Smith #define MATSHELL "shell" 39209238afSKris Buschelman #define MATSEQDENSE "seqdense" 40273d9f13SBarry Smith #define MATMPIDENSE "mpidense" 41209238afSKris Buschelman #define MATDENSE "dense" 42273d9f13SBarry Smith #define MATSEQBAIJ "seqbaij" 43273d9f13SBarry Smith #define MATMPIBAIJ "mpibaij" 44209238afSKris Buschelman #define MATBAIJ "baij" 45273d9f13SBarry Smith #define MATMPIADJ "mpiadj" 46273d9f13SBarry Smith #define MATSEQSBAIJ "seqsbaij" 47273d9f13SBarry Smith #define MATMPISBAIJ "mpisbaij" 48209238afSKris Buschelman #define MATSBAIJ "sbaij" 49cebc7f6cSBarry Smith #define MATDAAD "daad" 50cebc7f6cSBarry Smith #define MATMFFD "mffd" 51c8a8475eSBarry Smith #define MATNORMAL "normal" 52ab92ecdeSBarry Smith #define MATLRC "lrc" 534b8d99adSRichard Tran Mills #define MATSEQCSRPERM "seqcsrperm" 5418def556SRichard Tran Mills #define MATMPICSRPERM "mpicsrperm" 55814655b8SBarry Smith #define MATCSRPERM "csrperm" 5681824310SBarry Smith #define MATSEQCRL "seqcrl" 57c02b24c5SBarry Smith #define MATMPICRL "mpicrl" 58c02b24c5SBarry Smith #define MATCRL "crl" 592a6744ebSBarry Smith #define MATSCATTER "scatter" 60421e10b8SBarry Smith #define MATBLOCKMAT "blockmat" 61793850ffSBarry Smith #define MATCOMPOSITE "composite" 625a7f1df3SHong Zhang #define MATSEQFFTW "seqfftw" 63557cca28SSatish Balay #define MATTRANSPOSEMAT "transpose" 6472ca8751SBarry Smith #define MATSCHURCOMPLEMENT "schurcomplement" 651d6018f0SLisandro Dalcin #define MATPYTHON "python" 66f91d8e95SBarry Smith #define MATHYPRESTRUCT "hyprestruct" 67a9e6138eSGlenn Hammond #define MATHYPRESSTRUCT "hypresstruct" 68ee1cef2cSJed Brown #define MATSUBMATRIX "submatrix" 69d91e6319SBarry Smith 70164db6ccSDmitry Karpeev 719f8a2776SDmitry Karpeev #define MATFWK "matfwk" 72164db6ccSDmitry Karpeev 73773941d6SBarry Smith 749989ab13SBarry Smith /*E 75c7393fdbSBarry Smith MatSolverPackage - String with the name of a PETSc matrix solver type. 769989ab13SBarry Smith 779989ab13SBarry Smith For example: "petsc" indicates what PETSc provides, "superlu" indicates either 789989ab13SBarry Smith SuperLU or SuperLU_Dist etc. 799989ab13SBarry Smith 809989ab13SBarry Smith 819989ab13SBarry Smith Level: beginner 829989ab13SBarry Smith 835c9eb25fSBarry Smith .seealso: MatGetFactor(), Mat, MatSetType(), MatType 849989ab13SBarry Smith E*/ 85c7393fdbSBarry Smith #define MatSolverPackage char* 86b24902e0SBarry Smith #define MAT_SOLVER_SPOOLES "spooles" 87b24902e0SBarry Smith #define MAT_SOLVER_SUPERLU "superlu" 88b24902e0SBarry Smith #define MAT_SOLVER_SUPERLU_DIST "superlu_dist" 89b24902e0SBarry Smith #define MAT_SOLVER_UMFPACK "umfpack" 90b24902e0SBarry Smith #define MAT_SOLVER_ESSL "essl" 91b24902e0SBarry Smith #define MAT_SOLVER_LUSOL "lusol" 92b24902e0SBarry Smith #define MAT_SOLVER_MUMPS "mumps" 933bf14a46SMatthew Knepley #define MAT_SOLVER_PASTIX "pastix" 94b24902e0SBarry Smith #define MAT_SOLVER_DSCPACK "dscpack" 95b24902e0SBarry Smith #define MAT_SOLVER_MATLAB "matlab" 9643244d56SBarry Smith #define MAT_SOLVER_PETSC "petsc" 97b6806ab0SHong Zhang #define MAT_SOLVER_PLAPACK "plapack" 98845006b9SBarry Smith #define MAT_SOLVER_BAS "bas" 99773941d6SBarry Smith 100b24902e0SBarry Smith /*E 101b24902e0SBarry Smith MatFactorType - indicates what type of factorization is requested 102b24902e0SBarry Smith 103b24902e0SBarry Smith Level: beginner 104b24902e0SBarry Smith 105b24902e0SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 106b24902e0SBarry Smith 107c7393fdbSBarry Smith .seealso: MatSolverPackage, MatGetFactor() 108b24902e0SBarry Smith E*/ 109599ef60dSHong Zhang typedef enum {MAT_FACTOR_NONE, MAT_FACTOR_LU, MAT_FACTOR_CHOLESKY, MAT_FACTOR_ILU, MAT_FACTOR_ICC,MAT_FACTOR_ILUDT} MatFactorType; 110c7393fdbSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetFactor(Mat,const MatSolverPackage,MatFactorType,Mat*); 111db4efbfdSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetFactorAvailable(Mat,const MatSolverPackage,MatFactorType,PetscTruth*); 11235bd34faSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFactorGetSolverPackage(Mat,const MatSolverPackage*); 113b24902e0SBarry Smith 1149989ab13SBarry Smith 115c06d978dSMatthew Knepley /* Logging support */ 116*0700a824SBarry Smith #define MAT_FILE_CLASSID 1211216 /* used to indicate matrices in binary files */ 117*0700a824SBarry Smith extern PetscClassId PETSCMAT_DLLEXPORT MAT_CLASSID; 118*0700a824SBarry Smith extern PetscClassId PETSCMAT_DLLEXPORT MAT_FDCOLORING_CLASSID; 119*0700a824SBarry Smith extern PetscClassId PETSCMAT_DLLEXPORT MAT_PARTITIONING_CLASSID; 120*0700a824SBarry Smith extern PetscClassId PETSCMAT_DLLEXPORT MAT_NULLSPACE_CLASSID; 121*0700a824SBarry Smith extern PetscClassId PETSCMAT_DLLEXPORT MATMFFD_CLASSID; 122c06d978dSMatthew Knepley 123ceb03754SKris Buschelman /*E 124ceb03754SKris Buschelman MatReuse - Indicates if matrices obtained from a previous call to MatGetSubMatrices() 125ceb03754SKris Buschelman or MatGetSubMatrix() are to be reused to store the new matrix values. 126ceb03754SKris Buschelman 127ceb03754SKris Buschelman Level: beginner 128ceb03754SKris Buschelman 129ceb03754SKris Buschelman Any additions/changes here MUST also be made in include/finclude/petscmat.h 130ceb03754SKris Buschelman 1310c451bc4SBarry Smith .seealso: MatGetSubMatrices(), MatGetSubMatrix(), MatDestroyMatrices(), MatConvert() 132ceb03754SKris Buschelman E*/ 133dfe085dbSJed Brown typedef enum {MAT_INITIAL_MATRIX,MAT_REUSE_MATRIX,MAT_IGNORE_MATRIX} MatReuse; 134ceb03754SKris Buschelman 1355494a064SHong Zhang /*E 1365494a064SHong Zhang MatGetSubMatrixOption - Indicates if matrices obtained from a call to MatGetSubMatrices() 137829201f2SHong Zhang include the matrix values. Currently it is only used by MatGetSeqNonzerostructure(). 1385494a064SHong Zhang 1395494a064SHong Zhang Level: beginner 1405494a064SHong Zhang 141829201f2SHong Zhang .seealso: MatGetSeqNonzerostructure() 1425494a064SHong Zhang E*/ 1435494a064SHong Zhang typedef enum {MAT_DO_NOT_GET_VALUES,MAT_GET_VALUES} MatGetSubMatrixOption; 1445494a064SHong Zhang 145e5bd5246SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInitializePackage(const char[]); 146c06d978dSMatthew Knepley 147f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreate(MPI_Comm,Mat*); 148a69a960dSMatthew Knepley PetscPolymorphicFunction(MatCreate,(MPI_Comm comm),(comm,&A),Mat,A) 149a69a960dSMatthew Knepley PetscPolymorphicFunction(MatCreate,(),(PETSC_COMM_WORLD,&A),Mat,A) 150f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetSizes(Mat,PetscInt,PetscInt,PetscInt,PetscInt); 151a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetType(Mat,const MatType); 152be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetFromOptions(Mat); 153be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUpPreallocation(Mat); 154be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRegisterAll(const char[]); 155be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat)); 15630de9b25SBarry Smith 157f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetOptionsPrefix(Mat,const char[]); 158f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAppendOptionsPrefix(Mat,const char[]); 159f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOptionsPrefix(Mat,const char*[]); 160f69a0ea3SMatthew Knepley 16130de9b25SBarry Smith /*MC 16230de9b25SBarry Smith MatRegisterDynamic - Adds a new matrix type 16330de9b25SBarry Smith 16430de9b25SBarry Smith Synopsis: 165c1ac3661SBarry Smith PetscErrorCode MatRegisterDynamic(char *name,char *path,char *name_create,PetscErrorCode (*routine_create)(Mat)) 16630de9b25SBarry Smith 16730de9b25SBarry Smith Not Collective 16830de9b25SBarry Smith 16930de9b25SBarry Smith Input Parameters: 17030de9b25SBarry Smith + name - name of a new user-defined matrix type 17130de9b25SBarry Smith . path - path (either absolute or relative) the library containing this solver 17230de9b25SBarry Smith . name_create - name of routine to create method context 17330de9b25SBarry Smith - routine_create - routine to create method context 17430de9b25SBarry Smith 17530de9b25SBarry Smith Notes: 17630de9b25SBarry Smith MatRegisterDynamic() may be called multiple times to add several user-defined solvers. 17730de9b25SBarry Smith 17830de9b25SBarry Smith If dynamic libraries are used, then the fourth input argument (routine_create) 17930de9b25SBarry Smith is ignored. 18030de9b25SBarry Smith 18130de9b25SBarry Smith Sample usage: 18230de9b25SBarry Smith .vb 18330de9b25SBarry Smith MatRegisterDynamic("my_mat",/home/username/my_lib/lib/libO/solaris/mylib.a, 18430de9b25SBarry Smith "MyMatCreate",MyMatCreate); 18530de9b25SBarry Smith .ve 18630de9b25SBarry Smith 18730de9b25SBarry Smith Then, your solver can be chosen with the procedural interface via 18830de9b25SBarry Smith $ MatSetType(Mat,"my_mat") 18930de9b25SBarry Smith or at runtime via the option 19030de9b25SBarry Smith $ -mat_type my_mat 19130de9b25SBarry Smith 19230de9b25SBarry Smith Level: advanced 19330de9b25SBarry Smith 194ab901514SBarry Smith Notes: ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values. 19530de9b25SBarry Smith If your function is not being put into a shared library then use VecRegister() instead 19630de9b25SBarry Smith 19730de9b25SBarry Smith .keywords: Mat, register 19830de9b25SBarry Smith 19930de9b25SBarry Smith .seealso: MatRegisterAll(), MatRegisterDestroy() 20030de9b25SBarry Smith 20130de9b25SBarry Smith M*/ 202273d9f13SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 203273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,0) 204273d9f13SBarry Smith #else 205273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,d) 20630de9b25SBarry Smith #endif 20730de9b25SBarry Smith 208273d9f13SBarry Smith extern PetscTruth MatRegisterAllCalled; 209b0a32e0cSBarry Smith extern PetscFList MatList; 210b022a5c1SBarry Smith extern PetscFList MatColoringList; 211b022a5c1SBarry Smith extern PetscFList MatPartitioningList; 21228988994SBarry Smith 2133b224e63SBarry Smith /*E 2143b224e63SBarry Smith MatStructure - Indicates if the matrix has the same nonzero structure 2153b224e63SBarry Smith 2163b224e63SBarry Smith Level: beginner 2173b224e63SBarry Smith 2183b224e63SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 2193b224e63SBarry Smith 2203b224e63SBarry Smith .seealso: MatCopy(), KSPSetOperators(), PCSetOperators() 2213b224e63SBarry Smith E*/ 2223b224e63SBarry Smith typedef enum {SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER,SUBSET_NONZERO_PATTERN} MatStructure; 2233b224e63SBarry Smith 224be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqDense(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*); 225be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIDense(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*); 226be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 227ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,nz,nnz,&A),Mat,A) 228ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,&A),Mat,A) 229ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,0,nnz,&A),Mat,A) 230ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,&A),Mat,A) 231ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,A)) 232ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,m,n,0,nnz,A)) 233ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,A)) 234be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 235ba966639SSatish 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) 236ba966639SSatish 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) 237ba966639SSatish 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) 238ba966639SSatish 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) 239ba966639SSatish 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)) 240ba966639SSatish 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)) 241ba966639SSatish 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)) 242ba966639SSatish 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) 243ba966639SSatish 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) 244ba966639SSatish 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) 245ba966639SSatish 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) 246ba966639SSatish 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)) 247ba966639SSatish 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)) 248ba966639SSatish 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)) 24963ab93bfSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *); 2508d7a6e47SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJWithSplitArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],PetscInt[],PetscInt[],PetscScalar[],Mat*); 2518d7a6e47SBarry Smith 252be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 253ba966639SSatish 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) 254ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A) 255ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A) 256ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A) 257ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A)) 258ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A)) 259ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A)) 260be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 261ba966639SSatish 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) 262ba966639SSatish 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) 263ba966639SSatish 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) 264ba966639SSatish 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) 265ba966639SSatish 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)) 266ba966639SSatish 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)) 267ba966639SSatish 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)) 268ba966639SSatish 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) 269ba966639SSatish 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) 270ba966639SSatish 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) 271ba966639SSatish 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) 272ba966639SSatish 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)) 273ba966639SSatish 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)) 274ba966639SSatish 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)) 275be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAdj(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscInt[],Mat*); 276be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 277ba966639SSatish 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) 278ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A) 279ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A) 280ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A) 281ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A)) 282ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A)) 283ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A)) 284be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPISBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 285ba966639SSatish 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) 286ba966639SSatish 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) 287ba966639SSatish 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) 288ba966639SSatish 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) 289ba966639SSatish 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)) 290ba966639SSatish 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)) 291ba966639SSatish 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)) 292ba966639SSatish 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) 293ba966639SSatish 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) 294ba966639SSatish 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) 295ba966639SSatish 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) 296ba966639SSatish 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)) 297ba966639SSatish 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)) 298ba966639SSatish 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)) 299be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateShell(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,void *,Mat*); 300ba966639SSatish 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) 301ba966639SSatish Balay PetscPolymorphicFunction(MatCreateShell,(PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx),(PETSC_COMM_WORLD,m,n,M,N,ctx,&A),Mat,A) 302be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateAdic(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,void (*)(void),Mat*); 303be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateNormal(Mat,Mat*); 304ba966639SSatish Balay PetscPolymorphicFunction(MatCreateNormal,(Mat mat),(mat,&A),Mat,A) 305ab92ecdeSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateLRC(Mat,Mat,Mat,Mat*); 306284134d9SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateIS(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,ISLocalToGlobalMapping,Mat*); 3071472f72bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 3081472f72bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPICRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 3092a6744ebSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateScatter(MPI_Comm,VecScatter,Mat*); 3102a6744ebSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScatterSetVecScatter(Mat,VecScatter); 3112a6744ebSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScatterGetVecScatter(Mat,VecScatter*); 3128cdf2d9bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateBlockMat(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,Mat*); 313793850ffSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompositeAddMat(Mat,Mat); 314b52f573bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompositeMerge(Mat); 315793850ffSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateComposite(MPI_Comm,PetscInt,const Mat*,Mat*); 3166d7c1e57SBarry Smith typedef enum {MAT_COMPOSITE_ADDITIVE,MAT_COMPOSITE_MULTIPLICATIVE} MatCompositeType; 3176d7c1e57SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompositeSetType(Mat,MatCompositeType); 3186d7c1e57SBarry Smith 319164db6ccSDmitry Karpeev EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFwkAIJSetPreallocation(Mat A,PetscInt nz,PetscInt *nnz); 320164db6ccSDmitry Karpeev EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFwkGetDefaltBlockType(Mat A, const MatType *type); 321164db6ccSDmitry Karpeev EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFwkSetDefaltBlockType(Mat A, const MatType type); 322164db6ccSDmitry Karpeev EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFwkSetPreallocation_AIJ(Mat A,PetscInt nz,PetscInt *nnz); 323164db6ccSDmitry Karpeev EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFwkAddBlock(Mat A, PetscInt rowblock, PetscInt colblock, const MatType type, Mat *block); 324ad4c773dSDmitry Karpeev EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFwkSetBlock(Mat A, PetscInt rowblock, PetscInt colblock, Mat block); 325164db6ccSDmitry Karpeev EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFwkGetBlock(Mat A, PetscInt rowblock, PetscInt colblock, Mat *block); 326164db6ccSDmitry Karpeev EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFwkSetScatter(Mat A, Mat S, PetscInt blockcount, PetscInt *blocksizes); 327164db6ccSDmitry Karpeev EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFwkSetGather(Mat A, Mat G, PetscInt blockcount, PetscInt *blocksizes); 328a6a5cd3fSDmitry Karpeev 3295a7f1df3SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqFFTW(MPI_Comm,PetscInt,const PetscInt[],Mat*); 330f20085c4SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateTranspose(Mat,Mat*); 331ee1cef2cSJed Brown EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSubMatrix(Mat,IS,IS,Mat*); 332ee1cef2cSJed Brown EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSubMatrixUpdate(Mat,Mat,IS,IS); 3331472f72bSBarry Smith 3341d6018f0SLisandro Dalcin EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreatePython(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const char[],Mat*); 3351d6018f0SLisandro Dalcin EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPythonSetType(Mat,const char[]); 3361d6018f0SLisandro Dalcin 3371d6018f0SLisandro Dalcin 338f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUp(Mat); 339be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroy(Mat); 34021c89e3eSBarry Smith 3410c451bc4SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConjugate(Mat); 34299cafbc1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRealPart(Mat); 34399cafbc1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatImaginaryPart(Mat); 344bfc7d00aSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetDiagonalBlock(Mat,PetscTruth*,MatReuse,Mat*); 345fe97e370SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetTrace(Mat,PetscScalar*); 34699cafbc1SBarry Smith 3478ed539a5SBarry Smith /* ------------------------------------------------------------*/ 348be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 349be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 35087d4246cSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesRow(Mat,PetscInt,const PetscScalar[]); 351f259bd47SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesRowLocal(Mat,PetscInt,const PetscScalar[]); 35284cb2905SBarry Smith 3532ef4de8bSBarry Smith /*S 3542ef4de8bSBarry Smith MatStencil - Data structure (C struct) for storing information about a single row or 3552ef4de8bSBarry Smith column of a matrix as index on an associated grid. 3562ef4de8bSBarry Smith 3572ef4de8bSBarry Smith Level: beginner 3582ef4de8bSBarry Smith 3592ef4de8bSBarry Smith Concepts: matrix; linear operator 3602ef4de8bSBarry Smith 361d7bd93baSBarry Smith .seealso: MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockStencil() 3622ef4de8bSBarry Smith S*/ 363435da068SBarry Smith typedef struct { 364c1ac3661SBarry Smith PetscInt k,j,i,c; 365435da068SBarry Smith } MatStencil; 3662ef4de8bSBarry Smith 367be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode); 368be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode); 369be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt); 370435da068SBarry Smith 371be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetColoring(Mat,ISColoring); 372be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesAdic(Mat,void*); 373be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesAdifor(Mat,PetscInt,void*); 3743a7fca6bSBarry Smith 375d91e6319SBarry Smith /*E 376d91e6319SBarry Smith MatAssemblyType - Indicates if the matrix is now to be used, or if you plan 377d91e6319SBarry Smith to continue to add values to it 378d91e6319SBarry Smith 379d91e6319SBarry Smith Level: beginner 380d91e6319SBarry Smith 381d91e6319SBarry Smith .seealso: MatAssemblyBegin(), MatAssemblyEnd() 382d91e6319SBarry Smith E*/ 3836d4a8577SBarry Smith typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType; 384be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssemblyBegin(Mat,MatAssemblyType); 385be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssemblyEnd(Mat,MatAssemblyType); 386be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssembled(Mat,PetscTruth*); 3874f9c727eSBarry Smith 3881d73ed98SBarry Smith 38930de9b25SBarry Smith 390d91e6319SBarry Smith /*E 391d91e6319SBarry Smith MatOption - Options that may be set for a matrix and its behavior or storage 392d91e6319SBarry Smith 393d91e6319SBarry Smith Level: beginner 394d91e6319SBarry Smith 3950a835dfdSSatish Balay Any additions/changes here MUST also be made in include/finclude/petscmat.h 396d91e6319SBarry Smith 397d91e6319SBarry Smith .seealso: MatSetOption() 398d91e6319SBarry Smith E*/ 399512a5fc5SBarry Smith typedef enum {MAT_ROW_ORIENTED,MAT_NEW_NONZERO_LOCATIONS, 4004e0d8c25SBarry Smith MAT_SYMMETRIC, 4014e0d8c25SBarry Smith MAT_STRUCTURALLY_SYMMETRIC, 4028047e77bSBarry Smith MAT_NEW_DIAGONALS,MAT_IGNORE_OFF_PROC_ENTRIES, 4034e0d8c25SBarry Smith MAT_NEW_NONZERO_LOCATION_ERR, 4044e0d8c25SBarry Smith MAT_NEW_NONZERO_ALLOCATION_ERR,MAT_USE_HASH_TABLE, 405a9817697SBarry Smith MAT_KEEP_NONZERO_PATTERN,MAT_IGNORE_ZERO_ENTRIES, 4064e0d8c25SBarry Smith MAT_USE_INODES, 4074e0d8c25SBarry Smith MAT_HERMITIAN, 4084e0d8c25SBarry Smith MAT_SYMMETRY_ETERNAL, 4094e0d8c25SBarry Smith MAT_USE_COMPRESSEDROW, 4104e0d8c25SBarry Smith MAT_IGNORE_LOWER_TRIANGULAR,MAT_ERROR_LOWER_TRIANGULAR, 41128b2fa4aSMatthew Knepley MAT_GETROW_UPPERTRIANGULAR,MAT_UNUSED_NONZERO_LOCATION_ERR, 41228b2fa4aSMatthew Knepley NUM_MAT_OPTIONS} MatOption; 413290bbb0aSBarry Smith extern const char *MatOptions[]; 4144e0d8c25SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetOption(Mat,MatOption,PetscTruth); 415a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetType(Mat,const MatType*); 416a313700dSBarry Smith PetscPolymorphicFunction(MatGetType,(Mat mat),(mat,&t),const MatType,t) 41784cb2905SBarry Smith 418be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]); 419be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 420be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 421f5edf698SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowUpperTriangular(Mat); 422f5edf698SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRowUpperTriangular(Mat); 423be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 424be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 425be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumnVector(Mat,Vec,PetscInt); 426be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetArray(Mat,PetscScalar *[]); 427ba966639SSatish Balay PetscPolymorphicFunction(MatGetArray,(Mat mat),(mat,&a),PetscScalar*,a) 428be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreArray(Mat,PetscScalar *[]); 429be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBlockSize(Mat,PetscInt *); 430ba966639SSatish Balay PetscPolymorphicFunction(MatGetBlockSize,(Mat mat),(mat,&a),PetscInt,a) 431be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetBlockSize(Mat,PetscInt); 4327b80b807SBarry Smith 4331620fd73SBarry Smith 434be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMult(Mat,Vec,Vec); 43589c6957cSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultDiagonalBlock(Mat,Vec,Vec); 436be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultAdd(Mat,Vec,Vec,Vec); 437ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultAdd,(Mat A,Vec x,Vec y),(A,x,y,y)) 438be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTranspose(Mat,Vec,Vec); 439547795f9SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultHermitianTranspose(Mat,Vec,Vec); 440be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsTranspose(Mat,Mat,PetscReal,PetscTruth*); 441ba966639SSatish Balay PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B,PetscReal tol),(A,B,tol,&t),PetscTruth,t) 442e37f64cdSBarry Smith PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B),(A,B,0,&t),PetscTruth,t) 4431cbb95d3SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitianTranspose(Mat,Mat,PetscReal,PetscTruth*); 444be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeAdd(Mat,Vec,Vec,Vec); 4455cb05578SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultHermitianTransposeAdd(Mat,Vec,Vec,Vec); 446ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultTransposeAdd,(Mat A,Vec x,Vec y),(A,x,y,y)) 447be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultConstrained(Mat,Vec,Vec); 448be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeConstrained(Mat,Vec,Vec); 4492bebee5dSHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatSolve(Mat,Mat,Mat); 450f5edf698SHong Zhang 451d91e6319SBarry Smith /*E 452d91e6319SBarry Smith MatDuplicateOption - Indicates if a duplicated sparse matrix should have 453d91e6319SBarry Smith its numerical values copied over or just its nonzero structure. 454d91e6319SBarry Smith 455d91e6319SBarry Smith Level: beginner 456d91e6319SBarry Smith 457d91e6319SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 458d91e6319SBarry Smith 45970dcbbb9SBarry Smith $ MAT_SHARE_NONZERO_PATTERN - the i and j arrays in the new matrix will be shared with the original matrix 46070dcbbb9SBarry Smith $ this also triggers the MAT_DO_NOT_COPY_VALUES option. This is used when you 46170dcbbb9SBarry Smith $ have several matrices with the same nonzero pattern. 46270dcbbb9SBarry Smith 463d91e6319SBarry Smith .seealso: MatDuplicate() 464d91e6319SBarry Smith E*/ 46570dcbbb9SBarry Smith typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES,MAT_SHARE_NONZERO_PATTERN} MatDuplicateOption; 4662e8a6d31SBarry Smith 467a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert(Mat,const MatType,MatReuse,Mat*); 468a313700dSBarry Smith PetscPolymorphicFunction(MatConvert,(Mat A,const MatType t),(A,t,MAT_INITIAL_MATRIX,&a),Mat,a) 469be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDuplicate(Mat,MatDuplicateOption,Mat*); 470ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A,MatDuplicateOption o),(A,o,&a),Mat,a) 471ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A),(A,MAT_COPY_VALUES,&a),Mat,a) 47294a9d846SBarry Smith 473cb5b572fSBarry Smith 474be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCopy(Mat,Mat,MatStructure); 475be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatView(Mat,PetscViewer); 476be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsSymmetric(Mat,PetscReal,PetscTruth*); 477ba966639SSatish Balay PetscPolymorphicFunction(MatIsSymmetric,(Mat A,PetscReal tol),(A,tol,&t),PetscTruth,t) 478e37f64cdSBarry Smith PetscPolymorphicFunction(MatIsSymmetric,(Mat A),(A,0,&t),PetscTruth,t) 479be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsStructurallySymmetric(Mat,PetscTruth*); 480ba966639SSatish Balay PetscPolymorphicFunction(MatIsStructurallySymmetric,(Mat A),(A,&t),PetscTruth,t) 4811cbb95d3SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitian(Mat,PetscReal,PetscTruth*); 4821cbb95d3SBarry Smith PetscPolymorphicFunction(MatIsHermitian,(Mat A),(A,0,&t),PetscTruth,t) 483be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsSymmetricKnown(Mat,PetscTruth*,PetscTruth*); 484be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitianKnown(Mat,PetscTruth*,PetscTruth*); 48564c62002SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMissingDiagonal(Mat,PetscTruth *,PetscInt *); 486a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLoad(PetscViewer,const MatType,Mat*); 487a313700dSBarry Smith PetscPolymorphicFunction(MatLoad,(PetscViewer v,const MatType t),(v,t,&a),Mat,a) 4887b80b807SBarry Smith 4898f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *); 4908f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRowIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *); 4918f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumnIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *); 4928f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreColumnIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *); 493d4fbbf0eSBarry Smith 494d91e6319SBarry Smith /*S 495d91e6319SBarry Smith MatInfo - Context of matrix information, used with MatGetInfo() 496d91e6319SBarry Smith 497d91e6319SBarry Smith In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE 498d91e6319SBarry Smith 499d91e6319SBarry Smith Level: intermediate 500d91e6319SBarry Smith 501d91e6319SBarry Smith Concepts: matrix^nonzero information 502d91e6319SBarry Smith 503d9274352SBarry Smith .seealso: MatGetInfo(), MatInfoType 504d91e6319SBarry Smith S*/ 5054e220ebcSLois Curfman McInnes typedef struct { 506b0a32e0cSBarry Smith PetscLogDouble block_size; /* block size */ 507b0a32e0cSBarry Smith PetscLogDouble nz_allocated,nz_used,nz_unneeded; /* number of nonzeros */ 508b0a32e0cSBarry Smith PetscLogDouble memory; /* memory allocated */ 509b0a32e0cSBarry Smith PetscLogDouble assemblies; /* number of matrix assemblies called */ 510b0a32e0cSBarry Smith PetscLogDouble mallocs; /* number of mallocs during MatSetValues() */ 511b0a32e0cSBarry Smith PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */ 512b0a32e0cSBarry Smith PetscLogDouble factor_mallocs; /* number of mallocs during factorization */ 5134e220ebcSLois Curfman McInnes } MatInfo; 5144e220ebcSLois Curfman McInnes 515d9274352SBarry Smith /*E 516d9274352SBarry Smith MatInfoType - Indicates if you want information about the local part of the matrix, 517d9274352SBarry Smith the entire parallel matrix or the maximum over all the local parts. 518d9274352SBarry Smith 519d9274352SBarry Smith Level: beginner 520d9274352SBarry Smith 521d9274352SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 522d9274352SBarry Smith 523d9274352SBarry Smith .seealso: MatGetInfo(), MatInfo 524d9274352SBarry Smith E*/ 5257b80b807SBarry Smith typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType; 526be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetInfo(Mat,MatInfoType,MatInfo*); 527be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatValid(Mat,PetscTruth*); 528be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetDiagonal(Mat,Vec); 529985db425SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMax(Mat,Vec,PetscInt[]); 530985db425SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMin(Mat,Vec,PetscInt[]); 531985db425SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMaxAbs(Mat,Vec,PetscInt[]); 532c87e5d42SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMinAbs(Mat,Vec,PetscInt[]); 53386697f06SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowSum(Mat,Vec); 534fc4dec0aSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatTranspose(Mat,MatReuse,Mat*); 535fc4dec0aSBarry Smith PetscPolymorphicFunction(MatTranspose,(Mat A),(A,MAT_INITIAL_MATRIX,&t),Mat,t) 536eeffb40dSHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatHermitianTranspose(Mat,MatReuse,Mat*); 537be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPermute(Mat,IS,IS,Mat *); 538ba966639SSatish Balay PetscPolymorphicFunction(MatPermute,(Mat A,IS is1,IS is2),(A,is1,is2,&t),Mat,t) 539be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPermuteSparsify(Mat,PetscInt,PetscReal,PetscReal,IS,IS,Mat *); 540be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalScale(Mat,Vec,Vec); 541be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalSet(Mat,Vec,InsertMode); 542be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatEqual(Mat,Mat,PetscTruth*); 543ba966639SSatish Balay PetscPolymorphicFunction(MatEqual,(Mat A,Mat B),(A,B,&t),PetscTruth,t) 544be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultEqual(Mat,Mat,PetscInt,PetscTruth*); 545be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultAddEqual(Mat,Mat,PetscInt,PetscTruth*); 546be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeEqual(Mat,Mat,PetscInt,PetscTruth*); 547be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscTruth*); 5487b80b807SBarry Smith 549be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNorm(Mat,NormType,PetscReal *); 550ba966639SSatish Balay PetscPolymorphicFunction(MatNorm,(Mat A,NormType t),(A,t,&n),PetscReal,n) 551242f1d38SBarry Smith EXTERN PetscErrorCode MatGetColumnNorms(Mat,NormType,PetscReal *); 552be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroEntries(Mat); 553f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar); 554f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsIS(Mat,IS,PetscScalar); 555f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroColumns(Mat,PetscInt,const PetscInt [],const PetscScalar*); 556f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroColumnsIS(Mat,IS,const PetscScalar*); 5577b80b807SBarry Smith 558be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatUseScaledForm(Mat,PetscTruth); 559be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScaleSystem(Mat,Vec,Vec); 560be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatUnScaleSystem(Mat,Vec,Vec); 5615ef9f2a5SBarry Smith 562be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSize(Mat,PetscInt*,PetscInt*); 563be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalSize(Mat,PetscInt*,PetscInt*); 564be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRange(Mat,PetscInt*,PetscInt*); 5658ee2e534SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRanges(Mat,const PetscInt**); 566f4b8409dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRangeColumn(Mat,PetscInt*,PetscInt*); 56772ca8751SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRangesColumn(Mat,const PetscInt**); 5687b80b807SBarry Smith 569be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]); 570be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroyMatrices(PetscInt,Mat *[]); 5714aa3045dSJed Brown EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrix(Mat,IS,IS,MatReuse,Mat *); 572f6d58c54SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSeqNonzeroStructure(Mat,Mat*); 573f6d58c54SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroySeqNonzeroStructure(Mat*); 5745494a064SHong Zhang 575be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge(MPI_Comm,Mat,PetscInt,MatReuse,Mat*); 576be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPI(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*); 577be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPISymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*); 578be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPINumeric(Mat,Mat); 579be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroy_MPIAIJ_SeqsToMPI(Mat); 580be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalMat(Mat,MatReuse,Mat*); 581be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*); 582be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,PetscInt*,Mat*); 5831d79065fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBrowsOfAoCols(Mat,Mat,MatReuse,PetscInt**,PetscInt**,MatScalar**,Mat*); 58443eb5e2fSMatthew Knepley #if defined (PETSC_USE_CTABLE) 585cd116e26SSatish Balay #include "petscctable.h" 58643eb5e2fSMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetCommunicationStructs(Mat, Vec *, PetscTable *, VecScatter *); 58743eb5e2fSMatthew Knepley #else 58843eb5e2fSMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetCommunicationStructs(Mat, Vec *, PetscInt *[], VecScatter *); 58943eb5e2fSMatthew Knepley #endif 5908c7482ecSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetGhosts(Mat, PetscInt *,const PetscInt *[]); 5918efafbd8SBarry Smith 592be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt); 5937b80b807SBarry Smith 594be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*); 595be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultSymbolic(Mat,Mat,PetscReal,Mat*); 596be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultNumeric(Mat,Mat,Mat); 59722440eb1SKris Buschelman 598be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*); 599be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAPSymbolic(Mat,Mat,PetscReal,Mat*); 600be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAPNumeric(Mat,Mat,Mat); 60122440eb1SKris Buschelman 602be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTranspose(Mat,Mat,MatReuse,PetscReal,Mat*); 603be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTransposeSymbolic(Mat,Mat,PetscReal,Mat*); 604be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTransposeNumeric(Mat,Mat,Mat); 605bc011b1eSHong Zhang 606f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAXPY(Mat,PetscScalar,Mat,MatStructure); 60704aac2b0SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAYPX(Mat,PetscScalar,Mat,MatStructure); 6081c741599SHong Zhang 609f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScale(Mat,PetscScalar); 610f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShift(Mat,PetscScalar); 6117b80b807SBarry Smith 612be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping); 613be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetLocalToGlobalMappingBlock(Mat,ISLocalToGlobalMapping); 614f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar); 615f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsLocalIS(Mat,IS,PetscScalar); 616be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 617be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 618052efed2SBarry Smith 619be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStashSetInitialSize(Mat,PetscInt,PetscInt); 620be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*); 62190f02eecSBarry Smith 622be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInterpolate(Mat,Vec,Vec); 6230c451bc4SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInterpolateAdd(Mat,Vec,Vec,Vec); 624ba966639SSatish Balay PetscPolymorphicSubroutine(MatInterpolateAdd,(Mat A,Vec x,Vec y),(A,x,y,y)) 625be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestrict(Mat,Vec,Vec); 626be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetVecs(Mat,Vec*,Vec*); 62769db28dcSHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRedundantMatrix(Mat,PetscInt,MPI_Comm,PetscInt,MatReuse,Mat*); 628bd481603SBarry Smith 629bd481603SBarry Smith /*MC 6301620fd73SBarry Smith MatSetValue - Set a single entry into a matrix. 6311620fd73SBarry Smith 6321620fd73SBarry Smith Not collective 6331620fd73SBarry Smith 6341620fd73SBarry Smith Input Parameters: 6351620fd73SBarry Smith + m - the matrix 6361620fd73SBarry Smith . row - the row location of the entry 6371620fd73SBarry Smith . col - the column location of the entry 6381620fd73SBarry Smith . value - the value to insert 6391620fd73SBarry Smith - mode - either INSERT_VALUES or ADD_VALUES 6401620fd73SBarry Smith 6411620fd73SBarry Smith Notes: 6421620fd73SBarry Smith For efficiency one should use MatSetValues() and set several or many 6431620fd73SBarry Smith values simultaneously if possible. 6441620fd73SBarry Smith 6451620fd73SBarry Smith Level: beginner 6461620fd73SBarry Smith 6471620fd73SBarry Smith .seealso: MatSetValues(), MatSetValueLocal() 6481620fd73SBarry Smith M*/ 6491620fd73SBarry 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);} 6501620fd73SBarry Smith 6511620fd73SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatGetValue(Mat v,PetscInt i,PetscInt j,PetscScalar *va) {return MatGetValues(v,1,&i,1,&j,va);} 6521620fd73SBarry Smith 6531620fd73SBarry 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);} 6541620fd73SBarry Smith 6551620fd73SBarry Smith /*MC 6560d2a7ff2SSatish Balay MatPreallocateInitialize - Begins the block of code that will count the number of nonzeros per 657bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 658bd481603SBarry Smith 659bd481603SBarry Smith Synopsis: 660c1ac3661SBarry Smith PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz) 661bd481603SBarry Smith 662bd481603SBarry Smith Collective on MPI_Comm 663bd481603SBarry Smith 664bd481603SBarry Smith Input Parameters: 665bd481603SBarry Smith + comm - the communicator that will share the eventually allocated matrix 666859fcb39SBarry Smith . nrows - the number of LOCAL rows in the matrix 667859fcb39SBarry Smith - ncols - the number of LOCAL columns in the matrix 668bd481603SBarry Smith 669bd481603SBarry Smith Output Parameters: 670bd481603SBarry Smith + dnz - the array that will be passed to the matrix preallocation routines 671bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 672bd481603SBarry Smith 673bd481603SBarry Smith 674bd481603SBarry Smith Level: intermediate 675bd481603SBarry Smith 676bd481603SBarry Smith Notes: 677bd481603SBarry Smith See the chapter in the users manual on performance for more details 678bd481603SBarry Smith 6791d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 680bd481603SBarry Smith 681bd481603SBarry Smith Use MatPreallocateInitializeSymmetric() for symmetric matrices (MPISBAIJ matrices) 682bd481603SBarry Smith 6831620fd73SBarry Smith This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize(). 6841620fd73SBarry Smith 685bd481603SBarry Smith Concepts: preallocation^Matrix 686bd481603SBarry Smith 687bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(), 688bd481603SBarry Smith MatPreallocateInitializeSymmetric(), MatPreallocateSymmetricSetLocal() 689bd481603SBarry Smith M*/ 690c4f061fbSSatish Balay #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \ 6917c922b88SBarry Smith { \ 692a89bc211SBarry Smith PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__start,__end; \ 693a89bc211SBarry Smith _4_ierr = PetscMalloc2(__nrows,PetscInt,&dnz,__nrows,PetscInt,&onz);CHKERRQ(_4_ierr); \ 694a89bc211SBarry Smith _4_ierr = PetscMemzero(dnz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\ 695a89bc211SBarry Smith _4_ierr = PetscMemzero(onz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\ 696c1ac3661SBarry Smith _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ctmp;\ 697a89bc211SBarry Smith _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows; 6987c922b88SBarry Smith 699bd481603SBarry Smith /*MC 7000d2a7ff2SSatish Balay MatPreallocateSymmetricInitialize - Begins the block of code that will count the number of nonzeros per 701bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 702bd481603SBarry Smith 703bd481603SBarry Smith Synopsis: 704c1ac3661SBarry Smith PetscErrorCode MatPreallocateSymmetricInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz) 705bd481603SBarry Smith 706bd481603SBarry Smith Collective on MPI_Comm 707bd481603SBarry Smith 708bd481603SBarry Smith Input Parameters: 709bd481603SBarry Smith + comm - the communicator that will share the eventually allocated matrix 710859fcb39SBarry Smith . nrows - the number of LOCAL rows in the matrix 711859fcb39SBarry Smith - ncols - the number of LOCAL columns in the matrix 712bd481603SBarry Smith 713bd481603SBarry Smith Output Parameters: 714bd481603SBarry Smith + dnz - the array that will be passed to the matrix preallocation routines 715bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 716bd481603SBarry Smith 717bd481603SBarry Smith 718bd481603SBarry Smith Level: intermediate 719bd481603SBarry Smith 720bd481603SBarry Smith Notes: 721bd481603SBarry Smith See the chapter in the users manual on performance for more details 722bd481603SBarry Smith 7231d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 724bd481603SBarry Smith 7251620fd73SBarry Smith This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize(). 7261620fd73SBarry Smith 727bd481603SBarry Smith Concepts: preallocation^Matrix 728bd481603SBarry Smith 729bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(), 730bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal() 731bd481603SBarry Smith M*/ 732222b16d4SBarry Smith #define MatPreallocateSymmetricInitialize(comm,nrows,ncols,dnz,onz) 0; \ 733222b16d4SBarry Smith { \ 734a89bc211SBarry Smith PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__end; \ 735a89bc211SBarry Smith _4_ierr = PetscMalloc2(__nrows,PetscInt,&dnz,__nrows,PetscInt,&onz);CHKERRQ(_4_ierr); \ 736a89bc211SBarry Smith _4_ierr = PetscMemzero(dnz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\ 737a89bc211SBarry Smith _4_ierr = PetscMemzero(onz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\ 738c1ac3661SBarry Smith _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr);\ 739a89bc211SBarry Smith _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows; 740222b16d4SBarry Smith 741bd481603SBarry Smith /*MC 742bd481603SBarry Smith MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be 743bd481603SBarry Smith inserted using a local number of the rows and columns 744bd481603SBarry Smith 745bd481603SBarry Smith Synopsis: 746c1ac3661SBarry Smith PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 747bd481603SBarry Smith 748bd481603SBarry Smith Not Collective 749bd481603SBarry Smith 750bd481603SBarry Smith Input Parameters: 751bd481603SBarry Smith + map - the mapping between local numbering and global numbering 752bd481603SBarry Smith . nrows - the number of rows indicated 7531d73ed98SBarry Smith . rows - the indices of the rows 754bd481603SBarry Smith . ncols - the number of columns in the matrix 755bd481603SBarry Smith . cols - the columns indicated 756bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 757bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 758bd481603SBarry Smith 759bd481603SBarry Smith 760bd481603SBarry Smith Level: intermediate 761bd481603SBarry Smith 762bd481603SBarry Smith Notes: 763bd481603SBarry Smith See the chapter in the users manual on performance for more details 764bd481603SBarry Smith 7651d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 766bd481603SBarry Smith 767bd481603SBarry Smith Concepts: preallocation^Matrix 768bd481603SBarry Smith 769bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(), 770bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal() 771bd481603SBarry Smith M*/ 772c4f061fbSSatish Balay #define MatPreallocateSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\ 773c4f061fbSSatish Balay {\ 774c1ac3661SBarry Smith PetscInt __l;\ 775ef66eb69SBarry Smith _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\ 776ef66eb69SBarry Smith _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\ 777c4f061fbSSatish Balay for (__l=0;__l<nrows;__l++) {\ 778ef66eb69SBarry Smith _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 779c4f061fbSSatish Balay }\ 780c4f061fbSSatish Balay } 781c4f061fbSSatish Balay 782bd481603SBarry Smith /*MC 783bd481603SBarry Smith MatPreallocateSymmetricSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be 784bd481603SBarry Smith inserted using a local number of the rows and columns 785bd481603SBarry Smith 786bd481603SBarry Smith Synopsis: 787c1ac3661SBarry Smith PetscErrorCode MatPreallocateSymmetricSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 788bd481603SBarry Smith 789bd481603SBarry Smith Not Collective 790bd481603SBarry Smith 791bd481603SBarry Smith Input Parameters: 792bd481603SBarry Smith + map - the mapping between local numbering and global numbering 793bd481603SBarry Smith . nrows - the number of rows indicated 7941d73ed98SBarry Smith . rows - the indices of the rows 795bd481603SBarry Smith . ncols - the number of columns in the matrix 796bd481603SBarry Smith . cols - the columns indicated 797bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 798bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 799bd481603SBarry Smith 800bd481603SBarry Smith 801bd481603SBarry Smith Level: intermediate 802bd481603SBarry Smith 803bd481603SBarry Smith Notes: 804bd481603SBarry Smith See the chapter in the users manual on performance for more details 805bd481603SBarry Smith 806bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 807bd481603SBarry Smith 808bd481603SBarry Smith Concepts: preallocation^Matrix 809bd481603SBarry Smith 810bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(), 811bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal() 812bd481603SBarry Smith M*/ 813d3d32019SBarry Smith #define MatPreallocateSymmetricSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\ 814d3d32019SBarry Smith {\ 815c1ac3661SBarry Smith PetscInt __l;\ 816d3d32019SBarry Smith _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\ 817d3d32019SBarry Smith _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\ 818d3d32019SBarry Smith for (__l=0;__l<nrows;__l++) {\ 819d3d32019SBarry Smith _4_ierr = MatPreallocateSymmetricSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 820d3d32019SBarry Smith }\ 821d3d32019SBarry Smith } 822d3d32019SBarry Smith 823bd481603SBarry Smith /*MC 824bd481603SBarry Smith MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be 825bd481603SBarry Smith inserted using a local number of the rows and columns 826bd481603SBarry Smith 827bd481603SBarry Smith Synopsis: 828c1ac3661SBarry Smith PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 829bd481603SBarry Smith 830bd481603SBarry Smith Not Collective 831bd481603SBarry Smith 832bd481603SBarry Smith Input Parameters: 83364075487SBarry Smith + row - the row 834bd481603SBarry Smith . ncols - the number of columns in the matrix 835a50f8bf6SHong Zhang - cols - the columns indicated 836a50f8bf6SHong Zhang 837a50f8bf6SHong Zhang Output Parameters: 838a50f8bf6SHong Zhang + dnz - the array that will be passed to the matrix preallocation routines 839bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 840bd481603SBarry Smith 841bd481603SBarry Smith 842bd481603SBarry Smith Level: intermediate 843bd481603SBarry Smith 844bd481603SBarry Smith Notes: 845bd481603SBarry Smith See the chapter in the users manual on performance for more details 846bd481603SBarry Smith 847bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 848bd481603SBarry Smith 8491620fd73SBarry Smith This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize(). 8501620fd73SBarry Smith 851bd481603SBarry Smith Concepts: preallocation^Matrix 852bd481603SBarry Smith 853bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(), 854bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal() 855bd481603SBarry Smith M*/ 856c4f061fbSSatish Balay #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\ 857c1ac3661SBarry Smith { PetscInt __i; \ 8588ee2e534SBarry Smith if (row < __rstart) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Trying to set preallocation for row %D less than first local row %D",row,__rstart);\ 859a89bc211SBarry Smith if (row >= __rstart+__nrows) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Trying to set preallocation for row %D greater than last local row %D",row,__rstart+__nrows-1);\ 8607c922b88SBarry Smith for (__i=0; __i<nc; __i++) {\ 86164075487SBarry Smith if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \ 8627cc688d7SBarry Smith else dnz[row - __rstart]++;\ 8637c922b88SBarry Smith }\ 8647c922b88SBarry Smith } 8657c922b88SBarry Smith 866bd481603SBarry Smith /*MC 867bd481603SBarry Smith MatPreallocateSymmetricSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be 868bd481603SBarry Smith inserted using a local number of the rows and columns 869bd481603SBarry Smith 870bd481603SBarry Smith Synopsis: 871c1ac3661SBarry Smith PetscErrorCode MatPreallocateSymmetricSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 872bd481603SBarry Smith 873bd481603SBarry Smith Not Collective 874bd481603SBarry Smith 875bd481603SBarry Smith Input Parameters: 876bd481603SBarry Smith + nrows - the number of rows indicated 8771d73ed98SBarry Smith . rows - the indices of the rows 878bd481603SBarry Smith . ncols - the number of columns in the matrix 879bd481603SBarry Smith . cols - the columns indicated 880bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 881bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 882bd481603SBarry Smith 883bd481603SBarry Smith 884bd481603SBarry Smith Level: intermediate 885bd481603SBarry Smith 886bd481603SBarry Smith Notes: 887bd481603SBarry Smith See the chapter in the users manual on performance for more details 888bd481603SBarry Smith 889bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 890bd481603SBarry Smith 8911620fd73SBarry Smith This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize(). 8921620fd73SBarry Smith 893bd481603SBarry Smith Concepts: preallocation^Matrix 894bd481603SBarry Smith 895bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(), 896bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal() 897bd481603SBarry Smith M*/ 898d3d32019SBarry Smith #define MatPreallocateSymmetricSet(row,nc,cols,dnz,onz) 0;\ 899c1ac3661SBarry Smith { PetscInt __i; \ 900d3d32019SBarry Smith for (__i=0; __i<nc; __i++) {\ 901d3d32019SBarry Smith if (cols[__i] >= __end) onz[row - __rstart]++; \ 902d3d32019SBarry Smith else if (cols[__i] >= row) dnz[row - __rstart]++;\ 903d3d32019SBarry Smith }\ 904d3d32019SBarry Smith } 905d3d32019SBarry Smith 906bd481603SBarry Smith /*MC 90716371a99SBarry Smith MatPreallocateLocation - An alternative to MatPreallocationSet() that puts the nonzero locations into the matrix if it exists 90816371a99SBarry Smith 90916371a99SBarry Smith Synopsis: 91016371a99SBarry Smith PetscErrorCode MatPreallocateLocations(Mat A,PetscInt row,PetscInt ncols,PetscInt *cols,PetscInt *dnz,PetscInt *onz) 91116371a99SBarry Smith 91216371a99SBarry Smith Not Collective 91316371a99SBarry Smith 91416371a99SBarry Smith Input Parameters: 91516371a99SBarry Smith . A - matrix 91616371a99SBarry Smith . row - row where values exist (must be local to this process) 91716371a99SBarry Smith . ncols - number of columns 91816371a99SBarry Smith . cols - columns with nonzeros 91916371a99SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 92016371a99SBarry Smith - ozn - the other array passed to the matrix preallocation routines 92116371a99SBarry Smith 92216371a99SBarry Smith 92316371a99SBarry Smith Level: intermediate 92416371a99SBarry Smith 92516371a99SBarry Smith Notes: 92616371a99SBarry Smith See the chapter in the users manual on performance for more details 92716371a99SBarry Smith 92816371a99SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 92916371a99SBarry Smith 93016371a99SBarry Smith This is a MACRO not a function because it uses a bunch of variables private to the MatPreallocation.... routines. 93116371a99SBarry Smith 93216371a99SBarry Smith Concepts: preallocation^Matrix 93316371a99SBarry Smith 93416371a99SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(), 93516371a99SBarry Smith MatPreallocateSymmetricInitialize(), MatPreallocateSymmetricSetLocal() 93616371a99SBarry Smith M*/ 93716371a99SBarry Smith #define MatPreallocateLocation(A,row,ncols,cols,dnz,onz) 0;if (A) {ierr = MatSetValues(A,1,&row,ncols,cols,PETSC_NULL,INSERT_VALUES);CHKERRQ(ierr);} else {ierr = MatPreallocateSet(row,ncols,cols,dnz,onz);CHKERRQ(ierr);} 93816371a99SBarry Smith 93916371a99SBarry Smith 94016371a99SBarry Smith /*MC 9410d2a7ff2SSatish Balay MatPreallocateFinalize - Ends the block of code that will count the number of nonzeros per 942bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 943bd481603SBarry Smith 944bd481603SBarry Smith Synopsis: 945c1ac3661SBarry Smith PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz) 946bd481603SBarry Smith 947bd481603SBarry Smith Collective on MPI_Comm 948bd481603SBarry Smith 949bd481603SBarry Smith Input Parameters: 95016371a99SBarry Smith + dnz - the array that was be passed to the matrix preallocation routines 951bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 952bd481603SBarry Smith 953bd481603SBarry Smith 954bd481603SBarry Smith Level: intermediate 955bd481603SBarry Smith 956bd481603SBarry Smith Notes: 957bd481603SBarry Smith See the chapter in the users manual on performance for more details 958bd481603SBarry Smith 959bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 960bd481603SBarry Smith 9611620fd73SBarry Smith This is a MACRO not a function because it closes the { started in MatPreallocateInitialize(). 9621620fd73SBarry Smith 963bd481603SBarry Smith Concepts: preallocation^Matrix 964bd481603SBarry Smith 965bd481603SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(), 966bd481603SBarry Smith MatPreallocateSymmetricInitialize(), MatPreallocateSymmetricSetLocal() 967bd481603SBarry Smith M*/ 968a89bc211SBarry Smith #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree2(dnz,onz);CHKERRQ(_4_ierr);} 9697c922b88SBarry Smith 970bd481603SBarry Smith 971bd481603SBarry Smith 9727b80b807SBarry Smith /* Routines unique to particular data structures */ 973be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellGetContext(Mat,void **); 974ba966639SSatish Balay PetscPolymorphicFunction(MatShellGetContext,(Mat A),(A,&t),void*,t) 975698d4c6aSKris Buschelman 976be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInodeAdjustForInodes(Mat,IS*,IS*); 977be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *); 978698d4c6aSKris Buschelman 979be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetColumnIndices(Mat,PetscInt[]); 980be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetColumnIndices(Mat,PetscInt[]); 981be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 982c75a6043SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 983c75a6043SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 9847b80b807SBarry Smith 985a96a251dSBarry Smith #define MAT_SKIP_ALLOCATION -4 986a96a251dSBarry Smith 987be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]); 988ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz)) 989be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]); 990ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqSBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz)) 991be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]); 992ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqAIJSetPreallocation,(Mat A,const PetscInt nnz[]),(A,0,nnz)) 993273d9f13SBarry Smith 994be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 995ba966639SSatish Balay PetscPolymorphicSubroutine(MatMPIBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[],const PetscInt onz[]),(A,bs,0,nnz,0,onz)) 996be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 997be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 99853803f26SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []); 999725b52f3SLisandro Dalcin EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 1000be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]); 1001be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 1002be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]); 1003be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIDenseSetPreallocation(Mat,PetscScalar[]); 1004284134d9SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation(Mat,PetscScalar[]); 1005be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,PetscInt*[]); 1006be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,PetscInt*[]); 1007be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAdicSetLocalFunction(Mat,void (*)(void)); 1008273d9f13SBarry Smith 1009be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetLDA(Mat,PetscInt); 1010ab92ecdeSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDenseGetLocalMatrix(Mat,Mat*); 10111b807ce4Svictorle 1012be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStoreValues(Mat); 1013be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRetrieveValues(Mat); 10142e8a6d31SBarry Smith 1015be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDAADSetCtx(Mat,void*); 10163a7fca6bSBarry Smith 10177b80b807SBarry Smith /* 10187b80b807SBarry Smith These routines are not usually accessed directly, rather solving is 101994b7f48cSBarry Smith done through the KSP and PC interfaces. 10207b80b807SBarry Smith */ 10217b80b807SBarry Smith 1022d9274352SBarry Smith /*E 1023d9274352SBarry Smith MatOrderingType - String with the name of a PETSc matrix ordering or the creation function 1024d9274352SBarry Smith with an optional dynamic library name, for example 1025d9274352SBarry Smith http://www.mcs.anl.gov/petsc/lib.a:orderingcreate() 1026d9274352SBarry Smith 1027d9274352SBarry Smith Level: beginner 1028d9274352SBarry Smith 1029e5a9bf91SBarry Smith Cannot use const because the PC objects manipulate the string 1030e5a9bf91SBarry Smith 1031d9274352SBarry Smith .seealso: MatGetOrdering() 1032d9274352SBarry Smith E*/ 10333eaad2caSSatish Balay #define MatOrderingType char* 1034b12f92e5SBarry Smith #define MATORDERING_NATURAL "natural" 1035b12f92e5SBarry Smith #define MATORDERING_ND "nd" 1036b12f92e5SBarry Smith #define MATORDERING_1WD "1wd" 1037b12f92e5SBarry Smith #define MATORDERING_RCM "rcm" 1038b12f92e5SBarry Smith #define MATORDERING_QMD "qmd" 1039b12f92e5SBarry Smith #define MATORDERING_ROWLENGTH "rowlength" 1040e106eecfSBarry Smith #define MATORDERING_DSC_ND "dsc_nd" /* these three are only for DSCPACK, see its documentation for details */ 104162152c8bSBarry Smith #define MATORDERING_DSC_MMD "dsc_mmd" 104262152c8bSBarry Smith #define MATORDERING_DSC_MDF "dsc_mdf" 1043c06d978dSMatthew Knepley #define MATORDERING_CONSTRAINED "constrained" 1044c06d978dSMatthew Knepley #define MATORDERING_IDENTITY "identity" 1045c06d978dSMatthew Knepley #define MATORDERING_REVERSE "reverse" 10468fa24674SBarry Smith #define MATORDERING_FLOW "flow" 104799671722SJed Brown #define MATORDERING_AMD "amd" 1048b12f92e5SBarry Smith 1049f2416520SSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOrdering(Mat,const MatOrderingType,IS*,IS*); 1050be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOrderingList(PetscFList *list); 1051f2416520SSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,const MatOrderingType,IS*,IS*)); 105230de9b25SBarry Smith 105330de9b25SBarry Smith /*MC 105430de9b25SBarry Smith MatOrderingRegisterDynamic - Adds a new sparse matrix ordering to the 105530de9b25SBarry Smith matrix package. 105630de9b25SBarry Smith 105730de9b25SBarry Smith Synopsis: 1058c1ac3661SBarry Smith PetscErrorCode MatOrderingRegisterDynamic(char *name_ordering,char *path,char *name_create,PetscErrorCode (*routine_create)(MatOrdering)) 105930de9b25SBarry Smith 106030de9b25SBarry Smith Not Collective 106130de9b25SBarry Smith 106230de9b25SBarry Smith Input Parameters: 106330de9b25SBarry Smith + sname - name of ordering (for example MATORDERING_ND) 106430de9b25SBarry Smith . path - location of library where creation routine is 106530de9b25SBarry Smith . name - name of function that creates the ordering type,a string 106630de9b25SBarry Smith - function - function pointer that creates the ordering 106730de9b25SBarry Smith 106830de9b25SBarry Smith Level: developer 106930de9b25SBarry Smith 107030de9b25SBarry Smith If dynamic libraries are used, then the fourth input argument (function) 107130de9b25SBarry Smith is ignored. 107230de9b25SBarry Smith 107330de9b25SBarry Smith Sample usage: 107430de9b25SBarry Smith .vb 107530de9b25SBarry Smith MatOrderingRegisterDynamic("my_order",/home/username/my_lib/lib/libO/solaris/mylib.a, 107630de9b25SBarry Smith "MyOrder",MyOrder); 107730de9b25SBarry Smith .ve 107830de9b25SBarry Smith 107930de9b25SBarry Smith Then, your partitioner can be chosen with the procedural interface via 108030de9b25SBarry Smith $ MatOrderingSetType(part,"my_order) 108130de9b25SBarry Smith or at runtime via the option 10822401956bSBarry Smith $ -pc_factor_mat_ordering_type my_order 108330de9b25SBarry Smith 1084ab901514SBarry Smith ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values. 108530de9b25SBarry Smith 108630de9b25SBarry Smith .keywords: matrix, ordering, register 108730de9b25SBarry Smith 108830de9b25SBarry Smith .seealso: MatOrderingRegisterDestroy(), MatOrderingRegisterAll() 108930de9b25SBarry Smith M*/ 1090aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 1091f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,0) 1092b12f92e5SBarry Smith #else 1093f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,d) 1094b12f92e5SBarry Smith #endif 109530de9b25SBarry Smith 1096be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegisterDestroy(void); 1097be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegisterAll(const char[]); 10982bad1931SBarry Smith extern PetscTruth MatOrderingRegisterAllCalled; 1099b0a32e0cSBarry Smith extern PetscFList MatOrderingList; 1100d4fbbf0eSBarry Smith 1101be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS); 1102a2ce50c7SBarry Smith 1103d91e6319SBarry Smith /*S 1104d90ac83dSHong Zhang MatFactorShiftType - Numeric Shift. 1105d90ac83dSHong Zhang 1106d90ac83dSHong Zhang Level: beginner 1107d90ac83dSHong Zhang 1108d90ac83dSHong Zhang S*/ 1109d90ac83dSHong Zhang typedef enum {MAT_SHIFT_NONE,MAT_SHIFT_NONZERO,MAT_SHIFT_POSITIVE_DEFINITE,MAT_SHIFT_INBLOCKS} MatFactorShiftType; 1110d90ac83dSHong Zhang extern const char *MatFactorShiftTypes[]; 1111d90ac83dSHong Zhang 1112d90ac83dSHong Zhang /*S 11132401956bSBarry Smith MatFactorInfo - Data passed into the matrix factorization routines 1114b00f7748SHong Zhang 111561cad860SBarry Smith In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE, that is use 111661cad860SBarry Smith $ MatFactorInfo info(MAT_FACTORINFO_SIZE) 1117b00f7748SHong Zhang 111815e8a5b3SHong Zhang Notes: These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC. 1119b00f7748SHong Zhang 112061cad860SBarry Smith You can use MatFactorInfoInitialize() to set default values. 112161cad860SBarry Smith 1122b00f7748SHong Zhang Level: developer 1123b00f7748SHong Zhang 1124d7d82daaSBarry Smith .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(), 1125d7d82daaSBarry Smith MatFactorInfoInitialize() 1126b00f7748SHong Zhang 1127b00f7748SHong Zhang S*/ 1128b00f7748SHong Zhang typedef struct { 112915e8a5b3SHong Zhang PetscReal diagonal_fill; /* force diagonal to fill in if initially not filled */ 113085317021SBarry Smith PetscReal usedt; 113115e8a5b3SHong Zhang PetscReal dt; /* drop tolerance */ 1132b00f7748SHong Zhang PetscReal dtcol; /* tolerance for pivoting */ 113315e8a5b3SHong Zhang PetscReal dtcount; /* maximum nonzeros to be allowed per row */ 113467e28bfeSBarry Smith PetscReal fill; /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */ 1135348344bbSBarry Smith PetscReal levels; /* ICC/ILU(levels) */ 1136bcd9e38bSBarry Smith PetscReal pivotinblocks; /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0 1137bcd9e38bSBarry Smith factorization may be faster if do not pivot */ 113815e8a5b3SHong Zhang PetscReal zeropivot; /* pivot is called zero if less than this */ 1139f4db908eSBarry Smith PetscReal shifttype; /* type of shift added to matrix factor to prevent zero pivots */ 1140f4db908eSBarry Smith PetscReal shiftamount; /* how large the shift is */ 114115e8a5b3SHong Zhang } MatFactorInfo; 1142ffa6d0a5SLois Curfman McInnes 1143be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFactorInfoInitialize(MatFactorInfo*); 11440481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactor(Mat,IS,const MatFactorInfo*); 11450481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*); 11460481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactorNumeric(Mat,Mat,const MatFactorInfo*); 11470481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactor(Mat,IS,IS,const MatFactorInfo*); 11480481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUFactor(Mat,IS,IS,const MatFactorInfo*); 11490481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*); 11500481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*); 11510481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatICCFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*); 11520481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatICCFactor(Mat,IS,const MatFactorInfo*); 11530481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactorNumeric(Mat,Mat,const MatFactorInfo*); 1154be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*); 1155be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolve(Mat,Vec,Vec); 1156be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatForwardSolve(Mat,Vec,Vec); 1157be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatBackwardSolve(Mat,Vec,Vec); 1158be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveAdd(Mat,Vec,Vec,Vec); 1159be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveTranspose(Mat,Vec,Vec); 1160be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveTransposeAdd(Mat,Vec,Vec,Vec); 1161be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolves(Mat,Vecs,Vecs); 11628ed539a5SBarry Smith 1163be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUnfactored(Mat); 1164bb5a7306SBarry Smith 1165d91e6319SBarry Smith /*E 1166d91e6319SBarry Smith MatSORType - What type of (S)SOR to perform 1167bb1eb677SSatish Balay 1168d91e6319SBarry Smith Level: beginner 1169d91e6319SBarry Smith 1170d9274352SBarry Smith May be bitwise ORd together 1171d9274352SBarry Smith 1172d91e6319SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 1173d91e6319SBarry Smith 11744e7234bfSBarry Smith MatSORType may be bitwise ORd together, so do not change the numbers 11754e7234bfSBarry Smith 117641f059aeSBarry Smith .seealso: MatSOR() 1177d91e6319SBarry Smith E*/ 1178ee50ffe9SBarry Smith typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3, 1179ee50ffe9SBarry Smith SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8, 1180ee50ffe9SBarry Smith SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16, 118184cb2905SBarry Smith SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType; 118241f059aeSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSOR(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec); 11838ed539a5SBarry Smith 1184d4fbbf0eSBarry Smith /* 1185639f9d9dSBarry Smith These routines are for efficiently computing Jacobians via finite differences. 1186639f9d9dSBarry Smith */ 1187b12f92e5SBarry Smith 1188d9274352SBarry Smith /*E 1189d9274352SBarry Smith MatColoringType - String with the name of a PETSc matrix coloring or the creation function 1190d9274352SBarry Smith with an optional dynamic library name, for example 1191d9274352SBarry Smith http://www.mcs.anl.gov/petsc/lib.a:coloringcreate() 1192d9274352SBarry Smith 1193d9274352SBarry Smith Level: beginner 1194d9274352SBarry Smith 1195d9274352SBarry Smith .seealso: MatGetColoring() 1196d9274352SBarry Smith E*/ 1197a313700dSBarry Smith #define MatColoringType char* 1198b12f92e5SBarry Smith #define MATCOLORING_NATURAL "natural" 1199b12f92e5SBarry Smith #define MATCOLORING_SL "sl" 1200b12f92e5SBarry Smith #define MATCOLORING_LF "lf" 1201b12f92e5SBarry Smith #define MATCOLORING_ID "id" 1202b12f92e5SBarry Smith 1203ab970d9bSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColoring(Mat,const MatColoringType,ISColoring*); 12042e90c967SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,MatColoringType,ISColoring *)); 120530de9b25SBarry Smith 120630de9b25SBarry Smith /*MC 120730de9b25SBarry Smith MatColoringRegisterDynamic - Adds a new sparse matrix coloring to the 120830de9b25SBarry Smith matrix package. 120930de9b25SBarry Smith 121030de9b25SBarry Smith Synopsis: 1211c1ac3661SBarry Smith PetscErrorCode MatColoringRegisterDynamic(char *name_coloring,char *path,char *name_create,PetscErrorCode (*routine_create)(MatColoring)) 121230de9b25SBarry Smith 121330de9b25SBarry Smith Not Collective 121430de9b25SBarry Smith 121530de9b25SBarry Smith Input Parameters: 121630de9b25SBarry Smith + sname - name of Coloring (for example MATCOLORING_SL) 121730de9b25SBarry Smith . path - location of library where creation routine is 121830de9b25SBarry Smith . name - name of function that creates the Coloring type, a string 121930de9b25SBarry Smith - function - function pointer that creates the coloring 122030de9b25SBarry Smith 122130de9b25SBarry Smith Level: developer 122230de9b25SBarry Smith 122330de9b25SBarry Smith If dynamic libraries are used, then the fourth input argument (function) 122430de9b25SBarry Smith is ignored. 122530de9b25SBarry Smith 122630de9b25SBarry Smith Sample usage: 122730de9b25SBarry Smith .vb 122830de9b25SBarry Smith MatColoringRegisterDynamic("my_color",/home/username/my_lib/lib/libO/solaris/mylib.a, 122930de9b25SBarry Smith "MyColor",MyColor); 123030de9b25SBarry Smith .ve 123130de9b25SBarry Smith 123230de9b25SBarry Smith Then, your partitioner can be chosen with the procedural interface via 123330de9b25SBarry Smith $ MatColoringSetType(part,"my_color") 123430de9b25SBarry Smith or at runtime via the option 123530de9b25SBarry Smith $ -mat_coloring_type my_color 123630de9b25SBarry Smith 1237ab901514SBarry Smith $PETSC_ARCH occuring in pathname will be replaced with appropriate values. 123830de9b25SBarry Smith 123930de9b25SBarry Smith .keywords: matrix, Coloring, register 124030de9b25SBarry Smith 124130de9b25SBarry Smith .seealso: MatColoringRegisterDestroy(), MatColoringRegisterAll() 124230de9b25SBarry Smith M*/ 1243aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 1244f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,0) 1245b12f92e5SBarry Smith #else 1246f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,d) 1247b12f92e5SBarry Smith #endif 124830de9b25SBarry Smith 12492bad1931SBarry Smith extern PetscTruth MatColoringRegisterAllCalled; 1250f1144a30SSatish Balay 1251be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegisterAll(const char[]); 1252be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegisterDestroy(void); 1253be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*); 1254639f9d9dSBarry Smith 1255d9274352SBarry Smith /*S 1256d9274352SBarry Smith MatFDColoring - Object for computing a sparse Jacobian via finite differences 1257d9274352SBarry Smith and coloring 1258639f9d9dSBarry Smith 1259d9274352SBarry Smith Level: beginner 1260d9274352SBarry Smith 1261d9274352SBarry Smith Concepts: coloring, sparse Jacobian, finite differences 1262d9274352SBarry Smith 1263d9274352SBarry Smith .seealso: MatFDColoringCreate() 1264d9274352SBarry Smith S*/ 1265e2a1c21fSSatish Balay typedef struct _p_MatFDColoring* MatFDColoring; 1266639f9d9dSBarry Smith 1267be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringCreate(Mat,ISColoring,MatFDColoring *); 1268be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringDestroy(MatFDColoring); 1269be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringView(MatFDColoring,PetscViewer); 1270be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*); 12714a9d489dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**); 1272be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal); 1273be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFromOptions(MatFDColoring); 1274be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApply(Mat,MatFDColoring,Vec,MatStructure*,void *); 1275be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApplyTS(Mat,MatFDColoring,PetscReal,Vec,MatStructure*,void *); 1276be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetF(MatFDColoring,Vec); 1277be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,PetscInt*[]); 1278639f9d9dSBarry Smith /* 12790752156aSBarry Smith These routines are for partitioning matrices: currently used only 12803eda8832SBarry Smith for adjacency matrix, MatCreateMPIAdj(). 12810752156aSBarry Smith */ 1282ca161407SBarry Smith 1283d9274352SBarry Smith /*S 1284d9274352SBarry Smith MatPartitioning - Object for managing the partitioning of a matrix or graph 1285d9274352SBarry Smith 1286d9274352SBarry Smith Level: beginner 1287d9274352SBarry Smith 1288d9274352SBarry Smith Concepts: partitioning 1289d9274352SBarry Smith 1290743a1026Svictorle .seealso: MatPartitioningCreate(), MatPartitioningType 1291d9274352SBarry Smith S*/ 129291e9ee9fSBarry Smith typedef struct _p_MatPartitioning* MatPartitioning; 1293d9274352SBarry Smith 1294d9274352SBarry Smith /*E 12955bc6e7ccSvictorle MatPartitioningType - String with the name of a PETSc matrix partitioning or the creation function 1296d9274352SBarry Smith with an optional dynamic library name, for example 1297d9274352SBarry Smith http://www.mcs.anl.gov/petsc/lib.a:partitioningcreate() 1298d9274352SBarry Smith 1299d9274352SBarry Smith Level: beginner 1300d9274352SBarry Smith 1301b547a13bSvictorle .seealso: MatPartitioningCreate(), MatPartitioning 1302d9274352SBarry Smith E*/ 1303a313700dSBarry Smith #define MatPartitioningType char* 13048ba1e511SMatthew Knepley #define MAT_PARTITIONING_CURRENT "current" 1305c2ec2a73SHong Zhang #define MAT_PARTITIONING_SQUARE "square" 13068ba1e511SMatthew Knepley #define MAT_PARTITIONING_PARMETIS "parmetis" 130771306c60Sjroman #define MAT_PARTITIONING_CHACO "chaco" 130871306c60Sjroman #define MAT_PARTITIONING_JOSTLE "jostle" 130971306c60Sjroman #define MAT_PARTITIONING_PARTY "party" 131071306c60Sjroman #define MAT_PARTITIONING_SCOTCH "scotch" 131171306c60Sjroman 1312ca161407SBarry Smith 1313be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningCreate(MPI_Comm,MatPartitioning*); 1314a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetType(MatPartitioning,const MatPartitioningType); 1315be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetNParts(MatPartitioning,PetscInt); 1316be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetAdjacency(MatPartitioning,Mat); 1317be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]); 1318be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []); 1319be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningApply(MatPartitioning,IS*); 1320be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningDestroy(MatPartitioning); 13212aabb6bbSBarry Smith 1322be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatPartitioning)); 132330de9b25SBarry Smith 132430de9b25SBarry Smith /*MC 132530de9b25SBarry Smith MatPartitioningRegisterDynamic - Adds a new sparse matrix partitioning to the 132630de9b25SBarry Smith matrix package. 132730de9b25SBarry Smith 132830de9b25SBarry Smith Synopsis: 1329c1ac3661SBarry Smith PetscErrorCode MatPartitioningRegisterDynamic(char *name_partitioning,char *path,char *name_create,PetscErrorCode (*routine_create)(MatPartitioning)) 133030de9b25SBarry Smith 133130de9b25SBarry Smith Not Collective 133230de9b25SBarry Smith 133330de9b25SBarry Smith Input Parameters: 133430de9b25SBarry Smith + sname - name of partitioning (for example MAT_PARTITIONING_CURRENT) or parmetis 133530de9b25SBarry Smith . path - location of library where creation routine is 133630de9b25SBarry Smith . name - name of function that creates the partitioning type, a string 133730de9b25SBarry Smith - function - function pointer that creates the partitioning type 133830de9b25SBarry Smith 133930de9b25SBarry Smith Level: developer 134030de9b25SBarry Smith 134130de9b25SBarry Smith If dynamic libraries are used, then the fourth input argument (function) 134230de9b25SBarry Smith is ignored. 134330de9b25SBarry Smith 134430de9b25SBarry Smith Sample usage: 134530de9b25SBarry Smith .vb 134630de9b25SBarry Smith MatPartitioningRegisterDynamic("my_part",/home/username/my_lib/lib/libO/solaris/mylib.a, 134730de9b25SBarry Smith "MyPartCreate",MyPartCreate); 134830de9b25SBarry Smith .ve 134930de9b25SBarry Smith 135030de9b25SBarry Smith Then, your partitioner can be chosen with the procedural interface via 135130de9b25SBarry Smith $ MatPartitioningSetType(part,"my_part") 135230de9b25SBarry Smith or at runtime via the option 135330de9b25SBarry Smith $ -mat_partitioning_type my_part 135430de9b25SBarry Smith 1355ab901514SBarry Smith $PETSC_ARCH occuring in pathname will be replaced with appropriate values. 135630de9b25SBarry Smith 135730de9b25SBarry Smith .keywords: matrix, partitioning, register 135830de9b25SBarry Smith 135930de9b25SBarry Smith .seealso: MatPartitioningRegisterDestroy(), MatPartitioningRegisterAll() 136030de9b25SBarry Smith M*/ 1361aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 1362f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,0) 13632aabb6bbSBarry Smith #else 1364f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,d) 13652aabb6bbSBarry Smith #endif 13662aabb6bbSBarry Smith 13672bad1931SBarry Smith extern PetscTruth MatPartitioningRegisterAllCalled; 1368f1144a30SSatish Balay 1369be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegisterAll(const char[]); 1370be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegisterDestroy(void); 13712bad1931SBarry Smith 1372be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningView(MatPartitioning,PetscViewer); 1373be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetFromOptions(MatPartitioning); 1374a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningGetType(MatPartitioning,const MatPartitioningType*); 1375ca161407SBarry Smith 1376be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningParmetisSetCoarseSequential(MatPartitioning); 13770e8af3caSMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningParmetisGetEdgeCut(MatPartitioning, PetscInt *); 13780752156aSBarry Smith 1379be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningJostleSetCoarseLevel(MatPartitioning,PetscReal); 1380be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningJostleSetCoarseSequential(MatPartitioning); 138171306c60Sjroman 1382954b3ec6SBarry Smith typedef enum {MP_CHACO_MULTILEVEL_KL,MP_CHACO_SPECTRAL,MP_CHACO_LINEAR,MP_CHACO_RANDOM, MP_CHACO_SCATTERED} MPChacoGlobalType; 1383be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetGlobal(MatPartitioning, MPChacoGlobalType); 138471306c60Sjroman typedef enum { MP_CHACO_KERNIGHAN_LIN, MP_CHACO_NONE } MPChacoLocalType; 1385be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetLocal(MatPartitioning, MPChacoLocalType); 1386be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal); 138771306c60Sjroman typedef enum { MP_CHACO_LANCZOS, MP_CHACO_RQI_SYMMLQ } MPChacoEigenType; 1388be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType); 1389be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenTol(MatPartitioning, PetscReal); 1390be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenNumber(MatPartitioning, PetscInt); 139171306c60Sjroman 139271306c60Sjroman #define MP_PARTY_OPT "opt" 139371306c60Sjroman #define MP_PARTY_LIN "lin" 139471306c60Sjroman #define MP_PARTY_SCA "sca" 139571306c60Sjroman #define MP_PARTY_RAN "ran" 139671306c60Sjroman #define MP_PARTY_GBF "gbf" 139771306c60Sjroman #define MP_PARTY_GCF "gcf" 139871306c60Sjroman #define MP_PARTY_BUB "bub" 139971306c60Sjroman #define MP_PARTY_DEF "def" 1400be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetGlobal(MatPartitioning, const char*); 140171306c60Sjroman #define MP_PARTY_HELPFUL_SETS "hs" 140271306c60Sjroman #define MP_PARTY_KERNIGHAN_LIN "kl" 140371306c60Sjroman #define MP_PARTY_NONE "no" 1404be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetLocal(MatPartitioning, const char*); 1405be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal); 1406be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetBipart(MatPartitioning,PetscTruth); 1407be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscTruth); 140871306c60Sjroman 140971306c60Sjroman typedef enum { MP_SCOTCH_GREEDY, MP_SCOTCH_GPS, MP_SCOTCH_GR_GPS } MPScotchGlobalType; 1410be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetArch(MatPartitioning,const char*); 1411be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetMultilevel(MatPartitioning); 1412be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetGlobal(MatPartitioning,MPScotchGlobalType); 1413be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetCoarseLevel(MatPartitioning,PetscReal); 1414be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetHostList(MatPartitioning,const char*); 141571306c60Sjroman typedef enum { MP_SCOTCH_KERNIGHAN_LIN, MP_SCOTCH_NONE } MPScotchLocalType; 1416be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetLocal(MatPartitioning,MPScotchLocalType); 1417be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetMapping(MatPartitioning); 1418be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetStrategy(MatPartitioning,char*); 141971306c60Sjroman 1420591294e4SBarry Smith EXTERN PetscErrorCode MatMeshToVertexGraph(Mat,PetscInt,Mat*); 142114196391SBarry Smith EXTERN PetscErrorCode MatMeshToCellGraph(Mat,PetscInt,Mat*); 1422591294e4SBarry Smith 14230752156aSBarry Smith /* 14240a835dfdSSatish Balay If you add entries here you must also add them to finclude/petscmat.h 1425d4fbbf0eSBarry Smith */ 14261c1c02c0SLois Curfman McInnes typedef enum { MATOP_SET_VALUES=0, 14271c1c02c0SLois Curfman McInnes MATOP_GET_ROW=1, 14281c1c02c0SLois Curfman McInnes MATOP_RESTORE_ROW=2, 14291c1c02c0SLois Curfman McInnes MATOP_MULT=3, 14301c1c02c0SLois Curfman McInnes MATOP_MULT_ADD=4, 14317c922b88SBarry Smith MATOP_MULT_TRANSPOSE=5, 14327c922b88SBarry Smith MATOP_MULT_TRANSPOSE_ADD=6, 14331c1c02c0SLois Curfman McInnes MATOP_SOLVE=7, 14341c1c02c0SLois Curfman McInnes MATOP_SOLVE_ADD=8, 14357c922b88SBarry Smith MATOP_SOLVE_TRANSPOSE=9, 14367c922b88SBarry Smith MATOP_SOLVE_TRANSPOSE_ADD=10, 14371c1c02c0SLois Curfman McInnes MATOP_LUFACTOR=11, 14381c1c02c0SLois Curfman McInnes MATOP_CHOLESKYFACTOR=12, 1439a714c06dSBarry Smith MATOP_SOR=13, 14401c1c02c0SLois Curfman McInnes MATOP_TRANSPOSE=14, 14411c1c02c0SLois Curfman McInnes MATOP_GETINFO=15, 14421c1c02c0SLois Curfman McInnes MATOP_EQUAL=16, 14431c1c02c0SLois Curfman McInnes MATOP_GET_DIAGONAL=17, 14441c1c02c0SLois Curfman McInnes MATOP_DIAGONAL_SCALE=18, 14451c1c02c0SLois Curfman McInnes MATOP_NORM=19, 14461c1c02c0SLois Curfman McInnes MATOP_ASSEMBLY_BEGIN=20, 14471c1c02c0SLois Curfman McInnes MATOP_ASSEMBLY_END=21, 1448d519adbfSMatthew Knepley MATOP_SET_OPTION=22, 1449d519adbfSMatthew Knepley MATOP_ZERO_ENTRIES=23, 1450d519adbfSMatthew Knepley MATOP_ZERO_ROWS=24, 1451d519adbfSMatthew Knepley MATOP_LUFACTOR_SYMBOLIC=25, 1452d519adbfSMatthew Knepley MATOP_LUFACTOR_NUMERIC=26, 1453d519adbfSMatthew Knepley MATOP_CHOLESKY_FACTOR_SYMBOLIC=27, 1454d519adbfSMatthew Knepley MATOP_CHOLESKY_FACTOR_NUMERIC=28, 1455d519adbfSMatthew Knepley MATOP_SETUP_PREALLOCATION=29, 1456d519adbfSMatthew Knepley MATOP_ILUFACTOR_SYMBOLIC=30, 1457d519adbfSMatthew Knepley MATOP_ICCFACTOR_SYMBOLIC=31, 1458d519adbfSMatthew Knepley MATOP_GET_ARRAY=32, 1459d519adbfSMatthew Knepley MATOP_RESTORE_ARRAY=33, 1460d519adbfSMatthew Knepley MATOP_DUPLICATE=34, 1461d519adbfSMatthew Knepley MATOP_FORWARD_SOLVE=35, 1462d519adbfSMatthew Knepley MATOP_BACKWARD_SOLVE=36, 1463d519adbfSMatthew Knepley MATOP_ILUFACTOR=37, 1464d519adbfSMatthew Knepley MATOP_ICCFACTOR=38, 1465d519adbfSMatthew Knepley MATOP_AXPY=39, 1466d519adbfSMatthew Knepley MATOP_GET_SUBMATRICES=40, 1467d519adbfSMatthew Knepley MATOP_INCREASE_OVERLAP=41, 1468d519adbfSMatthew Knepley MATOP_GET_VALUES=42, 1469d519adbfSMatthew Knepley MATOP_COPY=43, 1470d519adbfSMatthew Knepley MATOP_GET_ROW_MAX=44, 1471d519adbfSMatthew Knepley MATOP_SCALE=45, 1472d519adbfSMatthew Knepley MATOP_SHIFT=46, 147335153367SBarry Smith MATOP_DIAGONAL_SET=47, 1474d519adbfSMatthew Knepley MATOP_ILUDT_FACTOR=48, 1475d519adbfSMatthew Knepley MATOP_SET_BLOCK_SIZE=49, 1476d519adbfSMatthew Knepley MATOP_GET_ROW_IJ=50, 1477d519adbfSMatthew Knepley MATOP_RESTORE_ROW_IJ=51, 1478d519adbfSMatthew Knepley MATOP_GET_COLUMN_IJ=52, 1479d519adbfSMatthew Knepley MATOP_RESTORE_COLUMN_IJ=53, 1480d519adbfSMatthew Knepley MATOP_FDCOLORING_CREATE=54, 1481d519adbfSMatthew Knepley MATOP_COLORING_PATCH=55, 1482d519adbfSMatthew Knepley MATOP_SET_UNFACTORED=56, 1483d519adbfSMatthew Knepley MATOP_PERMUTE=57, 1484d519adbfSMatthew Knepley MATOP_SET_VALUES_BLOCKED=58, 1485d519adbfSMatthew Knepley MATOP_GET_SUBMATRIX=59, 1486d519adbfSMatthew Knepley MATOP_DESTROY=60, 1487d519adbfSMatthew Knepley MATOP_VIEW=61, 1488d519adbfSMatthew Knepley MATOP_CONVERT_FROM=62, 1489d519adbfSMatthew Knepley MATOP_USE_SCALED_FORM=63, 1490d519adbfSMatthew Knepley MATOP_SCALE_SYSTEM=64, 1491d519adbfSMatthew Knepley MATOP_UNSCALE_SYSTEM=65, 1492d519adbfSMatthew Knepley MATOP_SET_LOCAL_TO_GLOBAL_MAP=66, 1493d519adbfSMatthew Knepley MATOP_SET_VALUES_LOCAL=67, 1494d519adbfSMatthew Knepley MATOP_ZERO_ROWS_LOCAL=68, 1495d519adbfSMatthew Knepley MATOP_GET_ROW_MAX_ABS=69, 1496d519adbfSMatthew Knepley MATOP_GET_ROW_MIN_ABS=70, 1497d519adbfSMatthew Knepley MATOP_CONVERT=71, 1498d519adbfSMatthew Knepley MATOP_SET_COLORING=72, 1499d519adbfSMatthew Knepley MATOP_SET_VALUES_ADIC=73, 1500d519adbfSMatthew Knepley MATOP_SET_VALUES_ADIFOR=74, 1501d519adbfSMatthew Knepley MATOP_FD_COLORING_APPLY=75, 1502d519adbfSMatthew Knepley MATOP_SET_FROM_OPTIONS=76, 1503d519adbfSMatthew Knepley MATOP_MULT_CON=77, 1504d519adbfSMatthew Knepley MATOP_MULT_TRANSPOSE_CON=78, 1505d519adbfSMatthew Knepley MATOP_PERMUTE_SPARSIFY=79, 1506d519adbfSMatthew Knepley MATOP_MULT_MULTIPLE=80, 1507d519adbfSMatthew Knepley MATOP_SOLVE_MULTIPLE=81, 1508d519adbfSMatthew Knepley MATOP_GET_INERTIA=82, 1509d519adbfSMatthew Knepley MATOP_LOAD=83, 1510d519adbfSMatthew Knepley MATOP_IS_SYMMETRIC=84, 1511d519adbfSMatthew Knepley MATOP_IS_HERMITIAN=85, 1512d519adbfSMatthew Knepley MATOP_IS_STRUCTURALLY_SYMMETRIC=86, 151341f059aeSBarry Smith MATOP_DUMMY=87, 1514d519adbfSMatthew Knepley MATOP_GET_VECS=88, 1515d519adbfSMatthew Knepley MATOP_MAT_MULT=89, 1516d519adbfSMatthew Knepley MATOP_MAT_MULT_SYMBOLIC=90, 1517d519adbfSMatthew Knepley MATOP_MAT_MULT_NUMERIC=91, 1518d519adbfSMatthew Knepley MATOP_PTAP=92, 1519d519adbfSMatthew Knepley MATOP_PTAP_SYMBOLIC=93, 1520d519adbfSMatthew Knepley MATOP_PTAP_NUMERIC=94, 1521d519adbfSMatthew Knepley MATOP_MAT_MULTTRANSPOSE=95, 15221763ac4eSSatish Balay MATOP_MAT_MULTTRANSPOSE_SYM=96, 15231763ac4eSSatish Balay MATOP_MAT_MULTTRANSPOSE_NUM=97, 1524d519adbfSMatthew Knepley MATOP_PTAP_SYMBOLIC_SEQAIJ=98, 1525d519adbfSMatthew Knepley MATOP_PTAP_NUMERIC_SEQAIJ=99, 1526d519adbfSMatthew Knepley MATOP_PTAP_SYMBOLIC_MPIAIJ=100, 1527d519adbfSMatthew Knepley MATOP_PTAP_NUMERIC_MPIAIJ=101, 1528d519adbfSMatthew Knepley MATOP_CONJUGATE=102, 1529d519adbfSMatthew Knepley MATOP_SET_SIZES=103, 1530d519adbfSMatthew Knepley MATOP_SET_VALUES_ROW=104, 1531d519adbfSMatthew Knepley MATOP_REAL_PART=105, 1532d519adbfSMatthew Knepley MATOP_IMAG_PART=106, 1533d519adbfSMatthew Knepley MATOP_GET_ROW_UTRIANGULAR=107, 1534d519adbfSMatthew Knepley MATOP_RESTORE_ROW_UTRIANGULAR=108, 1535d519adbfSMatthew Knepley MATOP_MATSOLVE=109, 1536d519adbfSMatthew Knepley MATOP_GET_REDUNDANTMATRIX=110, 1537d519adbfSMatthew Knepley MATOP_GET_ROW_MIN=111, 1538d519adbfSMatthew Knepley MATOP_GET_COLUMN_VEC=112, 1539d519adbfSMatthew Knepley MATOP_MISSING_DIAGONAL=113, 1540d519adbfSMatthew Knepley MATOP_MATGETSEQNONZEROSTRUCTURE=114, 154189c6957cSBarry Smith MATOP_CREATE=115, 154289c6957cSBarry Smith MATOP_GET_GHOSTS=116, 1543eeffb40dSHong Zhang MATOP_MULT_DIAGONAL_BLOCK=119, 1544547795f9SHong Zhang MATOP_HERMITIANTRANSPOSE=120, 1545547795f9SHong Zhang MATOP_MULTHERMITIANTRANSPOSE=121, 1546547795f9SHong Zhang MATOP_MULTHERMITIANTRANSPOSEADD=122 1547fae171e0SBarry Smith } MatOperation; 1548be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatHasOperation(Mat,MatOperation,PetscTruth*); 1549be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellSetOperation(Mat,MatOperation,void(*)(void)); 1550be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellGetOperation(Mat,MatOperation,void(**)(void)); 1551be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellSetContext(Mat,void*); 1552112a2221SBarry Smith 155390ace30eSBarry Smith /* 155490ace30eSBarry Smith Codes for matrices stored on disk. By default they are 155590ace30eSBarry Smith stored in a universal format. By changing the format with 15567973ad04SBarry Smith PetscViewerSetFormat(viewer,PETSC_VIEWER_NATIVE); the matrices will 155790ace30eSBarry Smith be stored in a way natural for the matrix, for example dense matrices 155890ace30eSBarry Smith would be stored as dense. Matrices stored this way may only be 155990ace30eSBarry Smith read into matrices of the same time. 156090ace30eSBarry Smith */ 156190ace30eSBarry Smith #define MATRIX_BINARY_FORMAT_DENSE -1 156290ace30eSBarry Smith 1563be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetHashTableFactor(Mat,PetscReal); 1564be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatISGetLocalMat(Mat,Mat*); 15651f4e1ec7SBarry Smith 1566d9274352SBarry Smith /*S 1567d9274352SBarry Smith MatNullSpace - Object that removes a null space from a vector, i.e. 1568d9274352SBarry Smith orthogonalizes the vector to a subsapce 1569d9274352SBarry Smith 1570f7a9e4ceSBarry Smith Level: advanced 1571d9274352SBarry Smith 1572d9274352SBarry Smith Concepts: matrix; linear operator, null space 1573d9274352SBarry Smith 15746e1639daSBarry Smith Users manual sections: 15756e1639daSBarry Smith . sec_singular 15766e1639daSBarry Smith 1577d9274352SBarry Smith .seealso: MatNullSpaceCreate() 1578d9274352SBarry Smith S*/ 157974637425SBarry Smith typedef struct _p_MatNullSpace* MatNullSpace; 1580d9274352SBarry Smith 1581be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceCreate(MPI_Comm,PetscTruth,PetscInt,const Vec[],MatNullSpace*); 1582b22b330cSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(MatNullSpace,Vec,void*),void*); 1583be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceDestroy(MatNullSpace); 1584be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceRemove(MatNullSpace,Vec,Vec*); 1585be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceAttach(Mat,MatNullSpace); 158695902228SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceTest(MatNullSpace,Mat,PetscTruth *); 158774637425SBarry Smith 1588be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatReorderingSeqSBAIJ(Mat,IS); 1589be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPISBAIJSetHashTableFactor(Mat,PetscReal); 1590be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetColumnIndices(Mat,PetscInt *); 159146129b97SKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJInvertBlockDiagonal(Mat); 15923f1d51d7SBarry Smith 1593be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMAIJ(Mat,PetscInt,Mat*); 1594be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMAIJRedimension(Mat,PetscInt,Mat*); 1595be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMAIJGetAIJ(Mat,Mat*); 1596c4f061fbSSatish Balay 1597be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatComputeExplicitOperator(Mat,Mat*); 1598b0a32e0cSBarry Smith 1599be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalScaleLocal(Mat,Vec); 160004f1ad80SBarry Smith 1601fef1beadSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*); 16023ec795f1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetBase(Mat,Vec,Vec); 16033ec795f1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*); 1604e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*)); 1605e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec)); 1606e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDAddNullSpace(Mat,MatNullSpace); 1607e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt); 1608e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDResetHHistory(Mat); 1609e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctionError(Mat,PetscReal); 1610e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetPeriod(Mat,PetscInt); 1611e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDGetH(Mat,PetscScalar *); 16126aa9148fSLisandro Dalcin EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetOptionsPrefix(Mat,const char[]); 1613e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFromOptions(Mat); 1614e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*); 1615e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*); 1616e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDDSSetUmin(Mat,PetscReal); 1617e884886fSBarry Smith 1618e884886fSBarry Smith 1619e884886fSBarry Smith typedef struct _p_MatMFFD* MatMFFD; 1620e884886fSBarry Smith 1621e884886fSBarry Smith /*E 1622e884886fSBarry Smith MatMFFDType - algorithm used to compute the h used in computing matrix-vector products via differencing of the function 1623e884886fSBarry Smith 1624e884886fSBarry Smith Level: beginner 1625e884886fSBarry Smith 1626e884886fSBarry Smith .seealso: MatMFFDSetType(), MatMFFDRegister() 1627e884886fSBarry Smith E*/ 1628a313700dSBarry Smith #define MatMFFDType char* 1629e884886fSBarry Smith #define MATMFFD_DS "ds" 1630e884886fSBarry Smith #define MATMFFD_WP "wp" 1631e884886fSBarry Smith 1632a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetType(Mat,const MatMFFDType); 1633e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatMFFD)); 1634e884886fSBarry Smith 1635e884886fSBarry Smith /*MC 1636e884886fSBarry Smith MatMFFDRegisterDynamic - Adds a method to the MatMFFD registry. 1637e884886fSBarry Smith 1638e884886fSBarry Smith Synopsis: 1639e884886fSBarry Smith PetscErrorCode MatMFFDRegisterDynamic(char *name_solver,char *path,char *name_create,PetscErrorCode (*routine_create)(MatMFFD)) 1640e884886fSBarry Smith 1641e884886fSBarry Smith Not Collective 1642e884886fSBarry Smith 1643e884886fSBarry Smith Input Parameters: 1644e884886fSBarry Smith + name_solver - name of a new user-defined compute-h module 1645e884886fSBarry Smith . path - path (either absolute or relative) the library containing this solver 1646e884886fSBarry Smith . name_create - name of routine to create method context 1647e884886fSBarry Smith - routine_create - routine to create method context 1648e884886fSBarry Smith 1649e884886fSBarry Smith Level: developer 1650e884886fSBarry Smith 1651e884886fSBarry Smith Notes: 1652e884886fSBarry Smith MatMFFDRegisterDynamic() may be called multiple times to add several user-defined solvers. 1653e884886fSBarry Smith 1654e884886fSBarry Smith If dynamic libraries are used, then the fourth input argument (routine_create) 1655e884886fSBarry Smith is ignored. 1656e884886fSBarry Smith 1657e884886fSBarry Smith Sample usage: 1658e884886fSBarry Smith .vb 1659e884886fSBarry Smith MatMFFDRegisterDynamic("my_h",/home/username/my_lib/lib/libO/solaris/mylib.a, 1660e884886fSBarry Smith "MyHCreate",MyHCreate); 1661e884886fSBarry Smith .ve 1662e884886fSBarry Smith 1663e884886fSBarry Smith Then, your solver can be chosen with the procedural interface via 1664e884886fSBarry Smith $ MatMFFDSetType(mfctx,"my_h") 1665e884886fSBarry Smith or at runtime via the option 1666e884886fSBarry Smith $ -snes_mf_type my_h 1667e884886fSBarry Smith 1668e884886fSBarry Smith .keywords: MatMFFD, register 1669e884886fSBarry Smith 1670e884886fSBarry Smith .seealso: MatMFFDRegisterAll(), MatMFFDRegisterDestroy() 1671e884886fSBarry Smith M*/ 1672e884886fSBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 1673e884886fSBarry Smith #define MatMFFDRegisterDynamic(a,b,c,d) MatMFFDRegister(a,b,c,0) 1674e884886fSBarry Smith #else 1675e884886fSBarry Smith #define MatMFFDRegisterDynamic(a,b,c,d) MatMFFDRegister(a,b,c,d) 1676e884886fSBarry Smith #endif 1677e884886fSBarry Smith 1678e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDRegisterAll(const char[]); 1679e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDRegisterDestroy(void); 1680e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDDefaultSetUmin(Mat,PetscReal); 1681e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDWPSetComputeNormU(Mat,PetscTruth); 1682e884886fSBarry Smith 1683e884886fSBarry Smith 1684be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *); 1685be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *); 16867dbadf16SMatthew Knepley 168797969023SHong Zhang /* 168897969023SHong Zhang PETSc interface to MUMPS 168997969023SHong Zhang */ 169097969023SHong Zhang #ifdef PETSC_HAVE_MUMPS 169197969023SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMumpsSetIcntl(Mat,PetscInt,PetscInt); 169297969023SHong Zhang #endif 169323a5497aSJed Brown 169423a5497aSJed Brown PETSC_EXTERN_CXX_END 169523a5497aSJed Brown #endif 1696