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 71dbc6227fSDmitry Karpeev #define MATDD "matdd" 72164db6ccSDmitry Karpeev 73ba2f8784SDmitry Karpeev #define MATIM "matim" 74ba2f8784SDmitry Karpeev 75773941d6SBarry Smith 769989ab13SBarry Smith /*E 77c7393fdbSBarry Smith MatSolverPackage - String with the name of a PETSc matrix solver type. 789989ab13SBarry Smith 799989ab13SBarry Smith For example: "petsc" indicates what PETSc provides, "superlu" indicates either 809989ab13SBarry Smith SuperLU or SuperLU_Dist etc. 819989ab13SBarry Smith 829989ab13SBarry Smith 839989ab13SBarry Smith Level: beginner 849989ab13SBarry Smith 855c9eb25fSBarry Smith .seealso: MatGetFactor(), Mat, MatSetType(), MatType 869989ab13SBarry Smith E*/ 87c7393fdbSBarry Smith #define MatSolverPackage char* 88b24902e0SBarry Smith #define MAT_SOLVER_SPOOLES "spooles" 89b24902e0SBarry Smith #define MAT_SOLVER_SUPERLU "superlu" 90b24902e0SBarry Smith #define MAT_SOLVER_SUPERLU_DIST "superlu_dist" 91b24902e0SBarry Smith #define MAT_SOLVER_UMFPACK "umfpack" 92b24902e0SBarry Smith #define MAT_SOLVER_ESSL "essl" 93b24902e0SBarry Smith #define MAT_SOLVER_LUSOL "lusol" 94b24902e0SBarry Smith #define MAT_SOLVER_MUMPS "mumps" 953bf14a46SMatthew Knepley #define MAT_SOLVER_PASTIX "pastix" 96b24902e0SBarry Smith #define MAT_SOLVER_DSCPACK "dscpack" 97b24902e0SBarry Smith #define MAT_SOLVER_MATLAB "matlab" 9843244d56SBarry Smith #define MAT_SOLVER_PETSC "petsc" 99b6806ab0SHong Zhang #define MAT_SOLVER_PLAPACK "plapack" 100845006b9SBarry Smith #define MAT_SOLVER_BAS "bas" 101773941d6SBarry Smith 102b24902e0SBarry Smith /*E 103b24902e0SBarry Smith MatFactorType - indicates what type of factorization is requested 104b24902e0SBarry Smith 105b24902e0SBarry Smith Level: beginner 106b24902e0SBarry Smith 107b24902e0SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 108b24902e0SBarry Smith 109c7393fdbSBarry Smith .seealso: MatSolverPackage, MatGetFactor() 110b24902e0SBarry Smith E*/ 111599ef60dSHong Zhang typedef enum {MAT_FACTOR_NONE, MAT_FACTOR_LU, MAT_FACTOR_CHOLESKY, MAT_FACTOR_ILU, MAT_FACTOR_ICC,MAT_FACTOR_ILUDT} MatFactorType; 112e92e720dSBarry Smith extern const char *MatFactorTypes[]; 113e92e720dSBarry Smith 114c7393fdbSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetFactor(Mat,const MatSolverPackage,MatFactorType,Mat*); 115db4efbfdSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetFactorAvailable(Mat,const MatSolverPackage,MatFactorType,PetscTruth*); 11635bd34faSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFactorGetSolverPackage(Mat,const MatSolverPackage*); 1176335c336SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetFactorType(Mat,MatFactorType*); 1189989ab13SBarry Smith 119c06d978dSMatthew Knepley /* Logging support */ 1200700a824SBarry Smith #define MAT_FILE_CLASSID 1211216 /* used to indicate matrices in binary files */ 1210700a824SBarry Smith extern PetscClassId PETSCMAT_DLLEXPORT MAT_CLASSID; 1220700a824SBarry Smith extern PetscClassId PETSCMAT_DLLEXPORT MAT_FDCOLORING_CLASSID; 1230700a824SBarry Smith extern PetscClassId PETSCMAT_DLLEXPORT MAT_PARTITIONING_CLASSID; 1240700a824SBarry Smith extern PetscClassId PETSCMAT_DLLEXPORT MAT_NULLSPACE_CLASSID; 1250700a824SBarry Smith extern PetscClassId PETSCMAT_DLLEXPORT MATMFFD_CLASSID; 126c06d978dSMatthew Knepley 127ceb03754SKris Buschelman /*E 128ceb03754SKris Buschelman MatReuse - Indicates if matrices obtained from a previous call to MatGetSubMatrices() 129ceb03754SKris Buschelman or MatGetSubMatrix() are to be reused to store the new matrix values. 130ceb03754SKris Buschelman 131ceb03754SKris Buschelman Level: beginner 132ceb03754SKris Buschelman 133ceb03754SKris Buschelman Any additions/changes here MUST also be made in include/finclude/petscmat.h 134ceb03754SKris Buschelman 1350c451bc4SBarry Smith .seealso: MatGetSubMatrices(), MatGetSubMatrix(), MatDestroyMatrices(), MatConvert() 136ceb03754SKris Buschelman E*/ 137dfe085dbSJed Brown typedef enum {MAT_INITIAL_MATRIX,MAT_REUSE_MATRIX,MAT_IGNORE_MATRIX} MatReuse; 138ceb03754SKris Buschelman 1395494a064SHong Zhang /*E 1405494a064SHong Zhang MatGetSubMatrixOption - Indicates if matrices obtained from a call to MatGetSubMatrices() 141829201f2SHong Zhang include the matrix values. Currently it is only used by MatGetSeqNonzerostructure(). 1425494a064SHong Zhang 1435494a064SHong Zhang Level: beginner 1445494a064SHong Zhang 145829201f2SHong Zhang .seealso: MatGetSeqNonzerostructure() 1465494a064SHong Zhang E*/ 1475494a064SHong Zhang typedef enum {MAT_DO_NOT_GET_VALUES,MAT_GET_VALUES} MatGetSubMatrixOption; 1485494a064SHong Zhang 149e5bd5246SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInitializePackage(const char[]); 150c06d978dSMatthew Knepley 151f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreate(MPI_Comm,Mat*); 152a69a960dSMatthew Knepley PetscPolymorphicFunction(MatCreate,(MPI_Comm comm),(comm,&A),Mat,A) 153a69a960dSMatthew Knepley PetscPolymorphicFunction(MatCreate,(),(PETSC_COMM_WORLD,&A),Mat,A) 154f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetSizes(Mat,PetscInt,PetscInt,PetscInt,PetscInt); 155a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetType(Mat,const MatType); 156be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetFromOptions(Mat); 157be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUpPreallocation(Mat); 158be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRegisterAll(const char[]); 159be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat)); 16030de9b25SBarry Smith 161f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetOptionsPrefix(Mat,const char[]); 162f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAppendOptionsPrefix(Mat,const char[]); 163f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOptionsPrefix(Mat,const char*[]); 164f69a0ea3SMatthew Knepley 16530de9b25SBarry Smith /*MC 16630de9b25SBarry Smith MatRegisterDynamic - Adds a new matrix type 16730de9b25SBarry Smith 16830de9b25SBarry Smith Synopsis: 1691890ba74SBarry Smith PetscErrorCode MatRegisterDynamic(const char *name,const char *path,const char *name_create,PetscErrorCode (*routine_create)(Mat)) 17030de9b25SBarry Smith 17130de9b25SBarry Smith Not Collective 17230de9b25SBarry Smith 17330de9b25SBarry Smith Input Parameters: 17430de9b25SBarry Smith + name - name of a new user-defined matrix type 17530de9b25SBarry Smith . path - path (either absolute or relative) the library containing this solver 17630de9b25SBarry Smith . name_create - name of routine to create method context 17730de9b25SBarry Smith - routine_create - routine to create method context 17830de9b25SBarry Smith 17930de9b25SBarry Smith Notes: 18030de9b25SBarry Smith MatRegisterDynamic() may be called multiple times to add several user-defined solvers. 18130de9b25SBarry Smith 18230de9b25SBarry Smith If dynamic libraries are used, then the fourth input argument (routine_create) 18330de9b25SBarry Smith is ignored. 18430de9b25SBarry Smith 18530de9b25SBarry Smith Sample usage: 18630de9b25SBarry Smith .vb 18730de9b25SBarry Smith MatRegisterDynamic("my_mat",/home/username/my_lib/lib/libO/solaris/mylib.a, 18830de9b25SBarry Smith "MyMatCreate",MyMatCreate); 18930de9b25SBarry Smith .ve 19030de9b25SBarry Smith 19130de9b25SBarry Smith Then, your solver can be chosen with the procedural interface via 19230de9b25SBarry Smith $ MatSetType(Mat,"my_mat") 19330de9b25SBarry Smith or at runtime via the option 19430de9b25SBarry Smith $ -mat_type my_mat 19530de9b25SBarry Smith 19630de9b25SBarry Smith Level: advanced 19730de9b25SBarry Smith 198ab901514SBarry Smith Notes: ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values. 19930de9b25SBarry Smith If your function is not being put into a shared library then use VecRegister() instead 20030de9b25SBarry Smith 20130de9b25SBarry Smith .keywords: Mat, register 20230de9b25SBarry Smith 20330de9b25SBarry Smith .seealso: MatRegisterAll(), MatRegisterDestroy() 20430de9b25SBarry Smith 20530de9b25SBarry Smith M*/ 206273d9f13SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 207273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,0) 208273d9f13SBarry Smith #else 209273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,d) 21030de9b25SBarry Smith #endif 21130de9b25SBarry Smith 212273d9f13SBarry Smith extern PetscTruth MatRegisterAllCalled; 213b0a32e0cSBarry Smith extern PetscFList MatList; 214b022a5c1SBarry Smith extern PetscFList MatColoringList; 215b022a5c1SBarry Smith extern PetscFList MatPartitioningList; 21628988994SBarry Smith 2173b224e63SBarry Smith /*E 2183b224e63SBarry Smith MatStructure - Indicates if the matrix has the same nonzero structure 2193b224e63SBarry Smith 2203b224e63SBarry Smith Level: beginner 2213b224e63SBarry Smith 2223b224e63SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 2233b224e63SBarry Smith 2243b224e63SBarry Smith .seealso: MatCopy(), KSPSetOperators(), PCSetOperators() 2253b224e63SBarry Smith E*/ 2263b224e63SBarry Smith typedef enum {SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER,SUBSET_NONZERO_PATTERN} MatStructure; 2273b224e63SBarry Smith 228be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqDense(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*); 229be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIDense(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*); 230be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 231ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,nz,nnz,&A),Mat,A) 232ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,&A),Mat,A) 233ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,0,nnz,&A),Mat,A) 234ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,&A),Mat,A) 235ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,A)) 236ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,m,n,0,nnz,A)) 237ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,A)) 238be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 239ba966639SSatish 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) 240ba966639SSatish 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) 241ba966639SSatish 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) 242ba966639SSatish 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) 243ba966639SSatish 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)) 244ba966639SSatish 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)) 245ba966639SSatish 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)) 246ba966639SSatish 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) 247ba966639SSatish 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) 248ba966639SSatish 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) 249ba966639SSatish 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) 250ba966639SSatish 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)) 251ba966639SSatish 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)) 252ba966639SSatish 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)) 25363ab93bfSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *); 2548d7a6e47SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJWithSplitArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],PetscInt[],PetscInt[],PetscScalar[],Mat*); 2558d7a6e47SBarry Smith 256be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 257ba966639SSatish 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) 258ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A) 259ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A) 260ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A) 261ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A)) 262ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A)) 263ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A)) 264be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 265ba966639SSatish 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) 266ba966639SSatish 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) 267ba966639SSatish 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) 268ba966639SSatish 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) 269ba966639SSatish 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)) 270ba966639SSatish 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)) 271ba966639SSatish 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)) 272ba966639SSatish 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) 273ba966639SSatish 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) 274ba966639SSatish 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) 275ba966639SSatish 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) 276ba966639SSatish 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)) 277ba966639SSatish 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)) 278ba966639SSatish 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)) 279be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAdj(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscInt[],Mat*); 280be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 281ba966639SSatish 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) 282ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A) 283ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A) 284ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A) 285ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A)) 286ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A)) 287ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A)) 288be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPISBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 289ba966639SSatish 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) 290ba966639SSatish 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) 291ba966639SSatish 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) 292ba966639SSatish 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) 293ba966639SSatish 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)) 294ba966639SSatish 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)) 295ba966639SSatish 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)) 296ba966639SSatish 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) 297ba966639SSatish 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) 298ba966639SSatish 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) 299ba966639SSatish 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) 300ba966639SSatish 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)) 301ba966639SSatish 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)) 302ba966639SSatish 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)) 303be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateShell(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,void *,Mat*); 304ba966639SSatish 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) 305ba966639SSatish Balay PetscPolymorphicFunction(MatCreateShell,(PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx),(PETSC_COMM_WORLD,m,n,M,N,ctx,&A),Mat,A) 306be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateAdic(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,void (*)(void),Mat*); 307be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateNormal(Mat,Mat*); 308ba966639SSatish Balay PetscPolymorphicFunction(MatCreateNormal,(Mat mat),(mat,&A),Mat,A) 309ab92ecdeSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateLRC(Mat,Mat,Mat,Mat*); 310284134d9SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateIS(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,ISLocalToGlobalMapping,Mat*); 3111472f72bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 3121472f72bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPICRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 3132a6744ebSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateScatter(MPI_Comm,VecScatter,Mat*); 3142a6744ebSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScatterSetVecScatter(Mat,VecScatter); 3152a6744ebSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScatterGetVecScatter(Mat,VecScatter*); 3168cdf2d9bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateBlockMat(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,Mat*); 317793850ffSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompositeAddMat(Mat,Mat); 318b52f573bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompositeMerge(Mat); 319793850ffSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateComposite(MPI_Comm,PetscInt,const Mat*,Mat*); 3206d7c1e57SBarry Smith typedef enum {MAT_COMPOSITE_ADDITIVE,MAT_COMPOSITE_MULTIPLICATIVE} MatCompositeType; 3216d7c1e57SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompositeSetType(Mat,MatCompositeType); 3226d7c1e57SBarry Smith 323dbc6227fSDmitry Karpeev typedef enum {MATDD_BLOCK_COMM_DEFAULT = 0, MATDD_BLOCK_COMM_SELF = -1, MATDD_BLOCK_COMM_DETERMINE = -2} MatDDBlockCommType; 324dbc6227fSDmitry Karpeev EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDDAIJSetPreallocation(Mat A,PetscInt nz,PetscInt *nnz); 325dbc6227fSDmitry Karpeev EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDDGetDefaltBlockType(Mat A, const MatType *type); 326dbc6227fSDmitry Karpeev EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDDSetDefaltBlockType(Mat A, const MatType type); 327dbc6227fSDmitry Karpeev EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDDSetPreallocation_AIJ(Mat A,PetscInt nz,PetscInt *nnz); 328dbc6227fSDmitry Karpeev EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDDAddBlock(Mat A, PetscInt rowblock, PetscInt colblock, const MatType blockmattype, MatDDBlockCommType blockcommtype, Mat *block); 329dbc6227fSDmitry Karpeev EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDDSetBlock(Mat A, PetscInt rowblock, PetscInt colblock, Mat block); 330dbc6227fSDmitry Karpeev EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDDGetBlock(Mat A, PetscInt rowblock, PetscInt colblock, Mat *block); 331dbc6227fSDmitry Karpeev EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDDSetScatters(Mat A, PetscInt blockcount, Mat scatters[]); 332dbc6227fSDmitry Karpeev EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDDSetGathers(Mat A, PetscInt blockcount, Mat gathers[]); 33362ee7b53SDmitry Karpeev 334ba2f8784SDmitry Karpeev #if defined PETSC_HAVE_MATIM 335ba2f8784SDmitry Karpeev EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIMSetIS(Mat A, IS in, IS out); 336ba2f8784SDmitry Karpeev EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIMGetIS(Mat A, IS *in, IS *out); 337ba2f8784SDmitry Karpeev #endif 338ba2f8784SDmitry Karpeev 339ba2f8784SDmitry Karpeev 340a6a5cd3fSDmitry Karpeev 3415a7f1df3SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqFFTW(MPI_Comm,PetscInt,const PetscInt[],Mat*); 342f20085c4SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateTranspose(Mat,Mat*); 343ee1cef2cSJed Brown EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSubMatrix(Mat,IS,IS,Mat*); 344ee1cef2cSJed Brown EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSubMatrixUpdate(Mat,Mat,IS,IS); 3451472f72bSBarry Smith 3461d6018f0SLisandro Dalcin EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreatePython(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const char[],Mat*); 3471d6018f0SLisandro Dalcin EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPythonSetType(Mat,const char[]); 3481d6018f0SLisandro Dalcin 3491d6018f0SLisandro Dalcin 350f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUp(Mat); 351be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroy(Mat); 35221c89e3eSBarry Smith 3530c451bc4SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConjugate(Mat); 35499cafbc1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRealPart(Mat); 35599cafbc1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatImaginaryPart(Mat); 356bfc7d00aSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetDiagonalBlock(Mat,PetscTruth*,MatReuse,Mat*); 357fe97e370SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetTrace(Mat,PetscScalar*); 35899cafbc1SBarry Smith 3598ed539a5SBarry Smith /* ------------------------------------------------------------*/ 360be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 361be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 36287d4246cSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesRow(Mat,PetscInt,const PetscScalar[]); 363f259bd47SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesRowLocal(Mat,PetscInt,const PetscScalar[]); 36484cb2905SBarry Smith 3652ef4de8bSBarry Smith /*S 3662ef4de8bSBarry Smith MatStencil - Data structure (C struct) for storing information about a single row or 3672ef4de8bSBarry Smith column of a matrix as index on an associated grid. 3682ef4de8bSBarry Smith 3692ef4de8bSBarry Smith Level: beginner 3702ef4de8bSBarry Smith 3712ef4de8bSBarry Smith Concepts: matrix; linear operator 3722ef4de8bSBarry Smith 373d7bd93baSBarry Smith .seealso: MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockStencil() 3742ef4de8bSBarry Smith S*/ 375435da068SBarry Smith typedef struct { 376c1ac3661SBarry Smith PetscInt k,j,i,c; 377435da068SBarry Smith } MatStencil; 3782ef4de8bSBarry Smith 379be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode); 380be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode); 381be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt); 382435da068SBarry Smith 383be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetColoring(Mat,ISColoring); 384be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesAdic(Mat,void*); 385be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesAdifor(Mat,PetscInt,void*); 3863a7fca6bSBarry Smith 387d91e6319SBarry Smith /*E 388d91e6319SBarry Smith MatAssemblyType - Indicates if the matrix is now to be used, or if you plan 389d91e6319SBarry Smith to continue to add values to it 390d91e6319SBarry Smith 391d91e6319SBarry Smith Level: beginner 392d91e6319SBarry Smith 393d91e6319SBarry Smith .seealso: MatAssemblyBegin(), MatAssemblyEnd() 394d91e6319SBarry Smith E*/ 3956d4a8577SBarry Smith typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType; 396be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssemblyBegin(Mat,MatAssemblyType); 397be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssemblyEnd(Mat,MatAssemblyType); 398be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssembled(Mat,PetscTruth*); 3994f9c727eSBarry Smith 4001d73ed98SBarry Smith 40130de9b25SBarry Smith 402d91e6319SBarry Smith /*E 403d91e6319SBarry Smith MatOption - Options that may be set for a matrix and its behavior or storage 404d91e6319SBarry Smith 405d91e6319SBarry Smith Level: beginner 406d91e6319SBarry Smith 4070a835dfdSSatish Balay Any additions/changes here MUST also be made in include/finclude/petscmat.h 408d91e6319SBarry Smith 409d91e6319SBarry Smith .seealso: MatSetOption() 410d91e6319SBarry Smith E*/ 411512a5fc5SBarry Smith typedef enum {MAT_ROW_ORIENTED,MAT_NEW_NONZERO_LOCATIONS, 4124e0d8c25SBarry Smith MAT_SYMMETRIC, 4134e0d8c25SBarry Smith MAT_STRUCTURALLY_SYMMETRIC, 4148047e77bSBarry Smith MAT_NEW_DIAGONALS,MAT_IGNORE_OFF_PROC_ENTRIES, 4154e0d8c25SBarry Smith MAT_NEW_NONZERO_LOCATION_ERR, 4164e0d8c25SBarry Smith MAT_NEW_NONZERO_ALLOCATION_ERR,MAT_USE_HASH_TABLE, 417a9817697SBarry Smith MAT_KEEP_NONZERO_PATTERN,MAT_IGNORE_ZERO_ENTRIES, 4184e0d8c25SBarry Smith MAT_USE_INODES, 4194e0d8c25SBarry Smith MAT_HERMITIAN, 4204e0d8c25SBarry Smith MAT_SYMMETRY_ETERNAL, 4214e0d8c25SBarry Smith MAT_USE_COMPRESSEDROW, 4224e0d8c25SBarry Smith MAT_IGNORE_LOWER_TRIANGULAR,MAT_ERROR_LOWER_TRIANGULAR, 42328b2fa4aSMatthew Knepley MAT_GETROW_UPPERTRIANGULAR,MAT_UNUSED_NONZERO_LOCATION_ERR, 42428b2fa4aSMatthew Knepley NUM_MAT_OPTIONS} MatOption; 425290bbb0aSBarry Smith extern const char *MatOptions[]; 4264e0d8c25SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetOption(Mat,MatOption,PetscTruth); 427a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetType(Mat,const MatType*); 428a313700dSBarry Smith PetscPolymorphicFunction(MatGetType,(Mat mat),(mat,&t),const MatType,t) 42984cb2905SBarry Smith 430be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]); 431be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 432be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 433f5edf698SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowUpperTriangular(Mat); 434f5edf698SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRowUpperTriangular(Mat); 435be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 436be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 437be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumnVector(Mat,Vec,PetscInt); 438be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetArray(Mat,PetscScalar *[]); 439ba966639SSatish Balay PetscPolymorphicFunction(MatGetArray,(Mat mat),(mat,&a),PetscScalar*,a) 440be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreArray(Mat,PetscScalar *[]); 441be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBlockSize(Mat,PetscInt *); 442ba966639SSatish Balay PetscPolymorphicFunction(MatGetBlockSize,(Mat mat),(mat,&a),PetscInt,a) 443be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetBlockSize(Mat,PetscInt); 4447b80b807SBarry Smith 4451620fd73SBarry Smith 446be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMult(Mat,Vec,Vec); 44789c6957cSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultDiagonalBlock(Mat,Vec,Vec); 448be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultAdd(Mat,Vec,Vec,Vec); 449ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultAdd,(Mat A,Vec x,Vec y),(A,x,y,y)) 450be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTranspose(Mat,Vec,Vec); 451547795f9SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultHermitianTranspose(Mat,Vec,Vec); 452be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsTranspose(Mat,Mat,PetscReal,PetscTruth*); 453ba966639SSatish Balay PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B,PetscReal tol),(A,B,tol,&t),PetscTruth,t) 454e37f64cdSBarry Smith PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B),(A,B,0,&t),PetscTruth,t) 4551cbb95d3SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitianTranspose(Mat,Mat,PetscReal,PetscTruth*); 456be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeAdd(Mat,Vec,Vec,Vec); 4575cb05578SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultHermitianTransposeAdd(Mat,Vec,Vec,Vec); 458ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultTransposeAdd,(Mat A,Vec x,Vec y),(A,x,y,y)) 459be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultConstrained(Mat,Vec,Vec); 460be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeConstrained(Mat,Vec,Vec); 4612bebee5dSHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatSolve(Mat,Mat,Mat); 462f5edf698SHong Zhang 463d91e6319SBarry Smith /*E 464d91e6319SBarry Smith MatDuplicateOption - Indicates if a duplicated sparse matrix should have 465d91e6319SBarry Smith its numerical values copied over or just its nonzero structure. 466d91e6319SBarry Smith 467d91e6319SBarry Smith Level: beginner 468d91e6319SBarry Smith 469d91e6319SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 470d91e6319SBarry Smith 47170dcbbb9SBarry Smith $ MAT_SHARE_NONZERO_PATTERN - the i and j arrays in the new matrix will be shared with the original matrix 47270dcbbb9SBarry Smith $ this also triggers the MAT_DO_NOT_COPY_VALUES option. This is used when you 47370dcbbb9SBarry Smith $ have several matrices with the same nonzero pattern. 47470dcbbb9SBarry Smith 475d91e6319SBarry Smith .seealso: MatDuplicate() 476d91e6319SBarry Smith E*/ 47770dcbbb9SBarry Smith typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES,MAT_SHARE_NONZERO_PATTERN} MatDuplicateOption; 4782e8a6d31SBarry Smith 479a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert(Mat,const MatType,MatReuse,Mat*); 480a313700dSBarry Smith PetscPolymorphicFunction(MatConvert,(Mat A,const MatType t),(A,t,MAT_INITIAL_MATRIX,&a),Mat,a) 481be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDuplicate(Mat,MatDuplicateOption,Mat*); 482ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A,MatDuplicateOption o),(A,o,&a),Mat,a) 483ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A),(A,MAT_COPY_VALUES,&a),Mat,a) 48494a9d846SBarry Smith 485cb5b572fSBarry Smith 486be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCopy(Mat,Mat,MatStructure); 487be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatView(Mat,PetscViewer); 488be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsSymmetric(Mat,PetscReal,PetscTruth*); 489ba966639SSatish Balay PetscPolymorphicFunction(MatIsSymmetric,(Mat A,PetscReal tol),(A,tol,&t),PetscTruth,t) 490e37f64cdSBarry Smith PetscPolymorphicFunction(MatIsSymmetric,(Mat A),(A,0,&t),PetscTruth,t) 491be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsStructurallySymmetric(Mat,PetscTruth*); 492ba966639SSatish Balay PetscPolymorphicFunction(MatIsStructurallySymmetric,(Mat A),(A,&t),PetscTruth,t) 4931cbb95d3SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitian(Mat,PetscReal,PetscTruth*); 4941cbb95d3SBarry Smith PetscPolymorphicFunction(MatIsHermitian,(Mat A),(A,0,&t),PetscTruth,t) 495be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsSymmetricKnown(Mat,PetscTruth*,PetscTruth*); 496be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitianKnown(Mat,PetscTruth*,PetscTruth*); 49764c62002SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMissingDiagonal(Mat,PetscTruth *,PetscInt *); 498a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLoad(PetscViewer,const MatType,Mat*); 499a313700dSBarry Smith PetscPolymorphicFunction(MatLoad,(PetscViewer v,const MatType t),(v,t,&a),Mat,a) 5007b80b807SBarry Smith 5018f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *); 5028f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRowIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *); 5038f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumnIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *); 5048f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreColumnIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *); 505d4fbbf0eSBarry Smith 506d91e6319SBarry Smith /*S 507d91e6319SBarry Smith MatInfo - Context of matrix information, used with MatGetInfo() 508d91e6319SBarry Smith 509d91e6319SBarry Smith In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE 510d91e6319SBarry Smith 511d91e6319SBarry Smith Level: intermediate 512d91e6319SBarry Smith 513d91e6319SBarry Smith Concepts: matrix^nonzero information 514d91e6319SBarry Smith 515d9274352SBarry Smith .seealso: MatGetInfo(), MatInfoType 516d91e6319SBarry Smith S*/ 5174e220ebcSLois Curfman McInnes typedef struct { 518b0a32e0cSBarry Smith PetscLogDouble block_size; /* block size */ 519b0a32e0cSBarry Smith PetscLogDouble nz_allocated,nz_used,nz_unneeded; /* number of nonzeros */ 520b0a32e0cSBarry Smith PetscLogDouble memory; /* memory allocated */ 521b0a32e0cSBarry Smith PetscLogDouble assemblies; /* number of matrix assemblies called */ 522b0a32e0cSBarry Smith PetscLogDouble mallocs; /* number of mallocs during MatSetValues() */ 523b0a32e0cSBarry Smith PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */ 524b0a32e0cSBarry Smith PetscLogDouble factor_mallocs; /* number of mallocs during factorization */ 5254e220ebcSLois Curfman McInnes } MatInfo; 5264e220ebcSLois Curfman McInnes 527d9274352SBarry Smith /*E 528d9274352SBarry Smith MatInfoType - Indicates if you want information about the local part of the matrix, 529d9274352SBarry Smith the entire parallel matrix or the maximum over all the local parts. 530d9274352SBarry Smith 531d9274352SBarry Smith Level: beginner 532d9274352SBarry Smith 533d9274352SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 534d9274352SBarry Smith 535d9274352SBarry Smith .seealso: MatGetInfo(), MatInfo 536d9274352SBarry Smith E*/ 5377b80b807SBarry Smith typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType; 538be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetInfo(Mat,MatInfoType,MatInfo*); 539be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatValid(Mat,PetscTruth*); 540be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetDiagonal(Mat,Vec); 541985db425SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMax(Mat,Vec,PetscInt[]); 542985db425SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMin(Mat,Vec,PetscInt[]); 543985db425SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMaxAbs(Mat,Vec,PetscInt[]); 544c87e5d42SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMinAbs(Mat,Vec,PetscInt[]); 54586697f06SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowSum(Mat,Vec); 546fc4dec0aSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatTranspose(Mat,MatReuse,Mat*); 547fc4dec0aSBarry Smith PetscPolymorphicFunction(MatTranspose,(Mat A),(A,MAT_INITIAL_MATRIX,&t),Mat,t) 548eeffb40dSHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatHermitianTranspose(Mat,MatReuse,Mat*); 549be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPermute(Mat,IS,IS,Mat *); 550ba966639SSatish Balay PetscPolymorphicFunction(MatPermute,(Mat A,IS is1,IS is2),(A,is1,is2,&t),Mat,t) 551be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPermuteSparsify(Mat,PetscInt,PetscReal,PetscReal,IS,IS,Mat *); 552be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalScale(Mat,Vec,Vec); 553be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalSet(Mat,Vec,InsertMode); 554be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatEqual(Mat,Mat,PetscTruth*); 555ba966639SSatish Balay PetscPolymorphicFunction(MatEqual,(Mat A,Mat B),(A,B,&t),PetscTruth,t) 556be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultEqual(Mat,Mat,PetscInt,PetscTruth*); 557be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultAddEqual(Mat,Mat,PetscInt,PetscTruth*); 558be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeEqual(Mat,Mat,PetscInt,PetscTruth*); 559be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscTruth*); 5607b80b807SBarry Smith 561be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNorm(Mat,NormType,PetscReal *); 562ba966639SSatish Balay PetscPolymorphicFunction(MatNorm,(Mat A,NormType t),(A,t,&n),PetscReal,n) 563242f1d38SBarry Smith EXTERN PetscErrorCode MatGetColumnNorms(Mat,NormType,PetscReal *); 564be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroEntries(Mat); 565f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar); 566f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsIS(Mat,IS,PetscScalar); 567f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroColumns(Mat,PetscInt,const PetscInt [],const PetscScalar*); 568f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroColumnsIS(Mat,IS,const PetscScalar*); 5697b80b807SBarry Smith 570be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatUseScaledForm(Mat,PetscTruth); 571be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScaleSystem(Mat,Vec,Vec); 572be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatUnScaleSystem(Mat,Vec,Vec); 5735ef9f2a5SBarry Smith 574be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSize(Mat,PetscInt*,PetscInt*); 575be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalSize(Mat,PetscInt*,PetscInt*); 576be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRange(Mat,PetscInt*,PetscInt*); 5778ee2e534SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRanges(Mat,const PetscInt**); 578f4b8409dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRangeColumn(Mat,PetscInt*,PetscInt*); 57972ca8751SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRangesColumn(Mat,const PetscInt**); 5807b80b807SBarry Smith 581be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]); 582be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroyMatrices(PetscInt,Mat *[]); 5834aa3045dSJed Brown EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrix(Mat,IS,IS,MatReuse,Mat *); 584f6d58c54SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSeqNonzeroStructure(Mat,Mat*); 585f6d58c54SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroySeqNonzeroStructure(Mat*); 5865494a064SHong Zhang 587be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge(MPI_Comm,Mat,PetscInt,MatReuse,Mat*); 588be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPI(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*); 589be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPISymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*); 590be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPINumeric(Mat,Mat); 591be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroy_MPIAIJ_SeqsToMPI(Mat); 592be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalMat(Mat,MatReuse,Mat*); 593be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*); 594be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,PetscInt*,Mat*); 5951d79065fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBrowsOfAoCols(Mat,Mat,MatReuse,PetscInt**,PetscInt**,MatScalar**,Mat*); 59643eb5e2fSMatthew Knepley #if defined (PETSC_USE_CTABLE) 597cd116e26SSatish Balay #include "petscctable.h" 59843eb5e2fSMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetCommunicationStructs(Mat, Vec *, PetscTable *, VecScatter *); 59943eb5e2fSMatthew Knepley #else 60043eb5e2fSMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetCommunicationStructs(Mat, Vec *, PetscInt *[], VecScatter *); 60143eb5e2fSMatthew Knepley #endif 6028c7482ecSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetGhosts(Mat, PetscInt *,const PetscInt *[]); 6038efafbd8SBarry Smith 604be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt); 6057b80b807SBarry Smith 606be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*); 607be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultSymbolic(Mat,Mat,PetscReal,Mat*); 608be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultNumeric(Mat,Mat,Mat); 60922440eb1SKris Buschelman 610be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*); 611be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAPSymbolic(Mat,Mat,PetscReal,Mat*); 612be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAPNumeric(Mat,Mat,Mat); 61322440eb1SKris Buschelman 614be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTranspose(Mat,Mat,MatReuse,PetscReal,Mat*); 615be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTransposeSymbolic(Mat,Mat,PetscReal,Mat*); 616be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTransposeNumeric(Mat,Mat,Mat); 617bc011b1eSHong Zhang 618f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAXPY(Mat,PetscScalar,Mat,MatStructure); 61904aac2b0SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAYPX(Mat,PetscScalar,Mat,MatStructure); 6201c741599SHong Zhang 621f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScale(Mat,PetscScalar); 622f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShift(Mat,PetscScalar); 6237b80b807SBarry Smith 624be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping); 625be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetLocalToGlobalMappingBlock(Mat,ISLocalToGlobalMapping); 626f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar); 627f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsLocalIS(Mat,IS,PetscScalar); 628be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 629be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 630052efed2SBarry Smith 631be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStashSetInitialSize(Mat,PetscInt,PetscInt); 632be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*); 63390f02eecSBarry Smith 634be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInterpolate(Mat,Vec,Vec); 6350c451bc4SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInterpolateAdd(Mat,Vec,Vec,Vec); 636ba966639SSatish Balay PetscPolymorphicSubroutine(MatInterpolateAdd,(Mat A,Vec x,Vec y),(A,x,y,y)) 637be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestrict(Mat,Vec,Vec); 638be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetVecs(Mat,Vec*,Vec*); 63969db28dcSHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRedundantMatrix(Mat,PetscInt,MPI_Comm,PetscInt,MatReuse,Mat*); 640bd481603SBarry Smith 641bd481603SBarry Smith /*MC 6421620fd73SBarry Smith MatSetValue - Set a single entry into a matrix. 6431620fd73SBarry Smith 6441620fd73SBarry Smith Not collective 6451620fd73SBarry Smith 6461620fd73SBarry Smith Input Parameters: 6471620fd73SBarry Smith + m - the matrix 6481620fd73SBarry Smith . row - the row location of the entry 6491620fd73SBarry Smith . col - the column location of the entry 6501620fd73SBarry Smith . value - the value to insert 6511620fd73SBarry Smith - mode - either INSERT_VALUES or ADD_VALUES 6521620fd73SBarry Smith 6531620fd73SBarry Smith Notes: 6541620fd73SBarry Smith For efficiency one should use MatSetValues() and set several or many 6551620fd73SBarry Smith values simultaneously if possible. 6561620fd73SBarry Smith 6571620fd73SBarry Smith Level: beginner 6581620fd73SBarry Smith 6591620fd73SBarry Smith .seealso: MatSetValues(), MatSetValueLocal() 6601620fd73SBarry Smith M*/ 6611620fd73SBarry 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);} 6621620fd73SBarry Smith 6631620fd73SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatGetValue(Mat v,PetscInt i,PetscInt j,PetscScalar *va) {return MatGetValues(v,1,&i,1,&j,va);} 6641620fd73SBarry Smith 6651620fd73SBarry 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);} 6661620fd73SBarry Smith 6671620fd73SBarry Smith /*MC 6680d2a7ff2SSatish Balay MatPreallocateInitialize - Begins the block of code that will count the number of nonzeros per 669bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 670bd481603SBarry Smith 671bd481603SBarry Smith Synopsis: 672c1ac3661SBarry Smith PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz) 673bd481603SBarry Smith 674bd481603SBarry Smith Collective on MPI_Comm 675bd481603SBarry Smith 676bd481603SBarry Smith Input Parameters: 677bd481603SBarry Smith + comm - the communicator that will share the eventually allocated matrix 678859fcb39SBarry Smith . nrows - the number of LOCAL rows in the matrix 679859fcb39SBarry Smith - ncols - the number of LOCAL columns in the matrix 680bd481603SBarry Smith 681bd481603SBarry Smith Output Parameters: 682bd481603SBarry Smith + dnz - the array that will be passed to the matrix preallocation routines 683bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 684bd481603SBarry Smith 685bd481603SBarry Smith 686bd481603SBarry Smith Level: intermediate 687bd481603SBarry Smith 688bd481603SBarry Smith Notes: 689bd481603SBarry Smith See the chapter in the users manual on performance for more details 690bd481603SBarry Smith 6911d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 692bd481603SBarry Smith 693bd481603SBarry Smith Use MatPreallocateInitializeSymmetric() for symmetric matrices (MPISBAIJ matrices) 694bd481603SBarry Smith 6951620fd73SBarry Smith This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize(). 6961620fd73SBarry Smith 697bd481603SBarry Smith Concepts: preallocation^Matrix 698bd481603SBarry Smith 699bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(), 700bd481603SBarry Smith MatPreallocateInitializeSymmetric(), MatPreallocateSymmetricSetLocal() 701bd481603SBarry Smith M*/ 702c4f061fbSSatish Balay #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \ 7037c922b88SBarry Smith { \ 704a89bc211SBarry Smith PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__start,__end; \ 705a89bc211SBarry Smith _4_ierr = PetscMalloc2(__nrows,PetscInt,&dnz,__nrows,PetscInt,&onz);CHKERRQ(_4_ierr); \ 706a89bc211SBarry Smith _4_ierr = PetscMemzero(dnz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\ 707a89bc211SBarry Smith _4_ierr = PetscMemzero(onz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\ 708c1ac3661SBarry Smith _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ctmp;\ 709a89bc211SBarry Smith _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows; 7107c922b88SBarry Smith 711bd481603SBarry Smith /*MC 7120d2a7ff2SSatish Balay MatPreallocateSymmetricInitialize - Begins the block of code that will count the number of nonzeros per 713bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 714bd481603SBarry Smith 715bd481603SBarry Smith Synopsis: 716c1ac3661SBarry Smith PetscErrorCode MatPreallocateSymmetricInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz) 717bd481603SBarry Smith 718bd481603SBarry Smith Collective on MPI_Comm 719bd481603SBarry Smith 720bd481603SBarry Smith Input Parameters: 721bd481603SBarry Smith + comm - the communicator that will share the eventually allocated matrix 722859fcb39SBarry Smith . nrows - the number of LOCAL rows in the matrix 723859fcb39SBarry Smith - ncols - the number of LOCAL columns in the matrix 724bd481603SBarry Smith 725bd481603SBarry Smith Output Parameters: 726bd481603SBarry Smith + dnz - the array that will be passed to the matrix preallocation routines 727bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 728bd481603SBarry Smith 729bd481603SBarry Smith 730bd481603SBarry Smith Level: intermediate 731bd481603SBarry Smith 732bd481603SBarry Smith Notes: 733bd481603SBarry Smith See the chapter in the users manual on performance for more details 734bd481603SBarry Smith 7351d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 736bd481603SBarry Smith 7371620fd73SBarry Smith This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize(). 7381620fd73SBarry Smith 739bd481603SBarry Smith Concepts: preallocation^Matrix 740bd481603SBarry Smith 741bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(), 742bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal() 743bd481603SBarry Smith M*/ 744222b16d4SBarry Smith #define MatPreallocateSymmetricInitialize(comm,nrows,ncols,dnz,onz) 0; \ 745222b16d4SBarry Smith { \ 746a89bc211SBarry Smith PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__end; \ 747a89bc211SBarry Smith _4_ierr = PetscMalloc2(__nrows,PetscInt,&dnz,__nrows,PetscInt,&onz);CHKERRQ(_4_ierr); \ 748a89bc211SBarry Smith _4_ierr = PetscMemzero(dnz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\ 749a89bc211SBarry Smith _4_ierr = PetscMemzero(onz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\ 750c1ac3661SBarry Smith _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr);\ 751a89bc211SBarry Smith _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows; 752222b16d4SBarry Smith 753bd481603SBarry Smith /*MC 754bd481603SBarry Smith MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be 755bd481603SBarry Smith inserted using a local number of the rows and columns 756bd481603SBarry Smith 757bd481603SBarry Smith Synopsis: 758c1ac3661SBarry Smith PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 759bd481603SBarry Smith 760bd481603SBarry Smith Not Collective 761bd481603SBarry Smith 762bd481603SBarry Smith Input Parameters: 763bd481603SBarry Smith + map - the mapping between local numbering and global numbering 764bd481603SBarry Smith . nrows - the number of rows indicated 7651d73ed98SBarry Smith . rows - the indices of the rows 766bd481603SBarry Smith . ncols - the number of columns in the matrix 767bd481603SBarry Smith . cols - the columns indicated 768bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 769bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 770bd481603SBarry Smith 771bd481603SBarry Smith 772bd481603SBarry Smith Level: intermediate 773bd481603SBarry Smith 774bd481603SBarry Smith Notes: 775bd481603SBarry Smith See the chapter in the users manual on performance for more details 776bd481603SBarry Smith 7771d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 778bd481603SBarry Smith 779bd481603SBarry Smith Concepts: preallocation^Matrix 780bd481603SBarry Smith 781bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(), 782bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal() 783bd481603SBarry Smith M*/ 784c4f061fbSSatish Balay #define MatPreallocateSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\ 785c4f061fbSSatish Balay {\ 786c1ac3661SBarry Smith PetscInt __l;\ 787ef66eb69SBarry Smith _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\ 788ef66eb69SBarry Smith _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\ 789c4f061fbSSatish Balay for (__l=0;__l<nrows;__l++) {\ 790ef66eb69SBarry Smith _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 791c4f061fbSSatish Balay }\ 792c4f061fbSSatish Balay } 793c4f061fbSSatish Balay 794bd481603SBarry Smith /*MC 795bd481603SBarry Smith MatPreallocateSymmetricSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be 796bd481603SBarry Smith inserted using a local number of the rows and columns 797bd481603SBarry Smith 798bd481603SBarry Smith Synopsis: 799c1ac3661SBarry Smith PetscErrorCode MatPreallocateSymmetricSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 800bd481603SBarry Smith 801bd481603SBarry Smith Not Collective 802bd481603SBarry Smith 803bd481603SBarry Smith Input Parameters: 804bd481603SBarry Smith + map - the mapping between local numbering and global numbering 805bd481603SBarry Smith . nrows - the number of rows indicated 8061d73ed98SBarry Smith . rows - the indices of the rows 807bd481603SBarry Smith . ncols - the number of columns in the matrix 808bd481603SBarry Smith . cols - the columns indicated 809bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 810bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 811bd481603SBarry Smith 812bd481603SBarry Smith 813bd481603SBarry Smith Level: intermediate 814bd481603SBarry Smith 815bd481603SBarry Smith Notes: 816bd481603SBarry Smith See the chapter in the users manual on performance for more details 817bd481603SBarry Smith 818bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 819bd481603SBarry Smith 820bd481603SBarry Smith Concepts: preallocation^Matrix 821bd481603SBarry Smith 822bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(), 823bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal() 824bd481603SBarry Smith M*/ 825d3d32019SBarry Smith #define MatPreallocateSymmetricSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\ 826d3d32019SBarry Smith {\ 827c1ac3661SBarry Smith PetscInt __l;\ 828d3d32019SBarry Smith _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\ 829d3d32019SBarry Smith _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\ 830d3d32019SBarry Smith for (__l=0;__l<nrows;__l++) {\ 831d3d32019SBarry Smith _4_ierr = MatPreallocateSymmetricSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 832d3d32019SBarry Smith }\ 833d3d32019SBarry Smith } 834d3d32019SBarry Smith 835bd481603SBarry Smith /*MC 836bd481603SBarry Smith MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be 837bd481603SBarry Smith inserted using a local number of the rows and columns 838bd481603SBarry Smith 839bd481603SBarry Smith Synopsis: 840c1ac3661SBarry Smith PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 841bd481603SBarry Smith 842bd481603SBarry Smith Not Collective 843bd481603SBarry Smith 844bd481603SBarry Smith Input Parameters: 84564075487SBarry Smith + row - the row 846bd481603SBarry Smith . ncols - the number of columns in the matrix 847a50f8bf6SHong Zhang - cols - the columns indicated 848a50f8bf6SHong Zhang 849a50f8bf6SHong Zhang Output Parameters: 850a50f8bf6SHong Zhang + dnz - the array that will be passed to the matrix preallocation routines 851bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 852bd481603SBarry Smith 853bd481603SBarry Smith 854bd481603SBarry Smith Level: intermediate 855bd481603SBarry Smith 856bd481603SBarry Smith Notes: 857bd481603SBarry Smith See the chapter in the users manual on performance for more details 858bd481603SBarry Smith 859bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 860bd481603SBarry Smith 8611620fd73SBarry Smith This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize(). 8621620fd73SBarry Smith 863bd481603SBarry Smith Concepts: preallocation^Matrix 864bd481603SBarry Smith 865bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(), 866bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal() 867bd481603SBarry Smith M*/ 868c4f061fbSSatish Balay #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\ 869c1ac3661SBarry Smith { PetscInt __i; \ 870*e32f2f54SBarry Smith if (row < __rstart) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Trying to set preallocation for row %D less than first local row %D",row,__rstart);\ 871*e32f2f54SBarry Smith if (row >= __rstart+__nrows) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Trying to set preallocation for row %D greater than last local row %D",row,__rstart+__nrows-1);\ 8727c922b88SBarry Smith for (__i=0; __i<nc; __i++) {\ 87364075487SBarry Smith if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \ 8747cc688d7SBarry Smith else dnz[row - __rstart]++;\ 8757c922b88SBarry Smith }\ 8767c922b88SBarry Smith } 8777c922b88SBarry Smith 878bd481603SBarry Smith /*MC 879bd481603SBarry Smith MatPreallocateSymmetricSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be 880bd481603SBarry Smith inserted using a local number of the rows and columns 881bd481603SBarry Smith 882bd481603SBarry Smith Synopsis: 883c1ac3661SBarry Smith PetscErrorCode MatPreallocateSymmetricSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 884bd481603SBarry Smith 885bd481603SBarry Smith Not Collective 886bd481603SBarry Smith 887bd481603SBarry Smith Input Parameters: 888bd481603SBarry Smith + nrows - the number of rows indicated 8891d73ed98SBarry Smith . rows - the indices of the rows 890bd481603SBarry Smith . ncols - the number of columns in the matrix 891bd481603SBarry Smith . cols - the columns indicated 892bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 893bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 894bd481603SBarry Smith 895bd481603SBarry Smith 896bd481603SBarry Smith Level: intermediate 897bd481603SBarry Smith 898bd481603SBarry Smith Notes: 899bd481603SBarry Smith See the chapter in the users manual on performance for more details 900bd481603SBarry Smith 901bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 902bd481603SBarry Smith 9031620fd73SBarry Smith This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize(). 9041620fd73SBarry Smith 905bd481603SBarry Smith Concepts: preallocation^Matrix 906bd481603SBarry Smith 907bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(), 908bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal() 909bd481603SBarry Smith M*/ 910d3d32019SBarry Smith #define MatPreallocateSymmetricSet(row,nc,cols,dnz,onz) 0;\ 911c1ac3661SBarry Smith { PetscInt __i; \ 912d3d32019SBarry Smith for (__i=0; __i<nc; __i++) {\ 913d3d32019SBarry Smith if (cols[__i] >= __end) onz[row - __rstart]++; \ 914d3d32019SBarry Smith else if (cols[__i] >= row) dnz[row - __rstart]++;\ 915d3d32019SBarry Smith }\ 916d3d32019SBarry Smith } 917d3d32019SBarry Smith 918bd481603SBarry Smith /*MC 91916371a99SBarry Smith MatPreallocateLocation - An alternative to MatPreallocationSet() that puts the nonzero locations into the matrix if it exists 92016371a99SBarry Smith 92116371a99SBarry Smith Synopsis: 92216371a99SBarry Smith PetscErrorCode MatPreallocateLocations(Mat A,PetscInt row,PetscInt ncols,PetscInt *cols,PetscInt *dnz,PetscInt *onz) 92316371a99SBarry Smith 92416371a99SBarry Smith Not Collective 92516371a99SBarry Smith 92616371a99SBarry Smith Input Parameters: 92716371a99SBarry Smith . A - matrix 92816371a99SBarry Smith . row - row where values exist (must be local to this process) 92916371a99SBarry Smith . ncols - number of columns 93016371a99SBarry Smith . cols - columns with nonzeros 93116371a99SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 93216371a99SBarry Smith - ozn - the other array passed to the matrix preallocation routines 93316371a99SBarry Smith 93416371a99SBarry Smith 93516371a99SBarry Smith Level: intermediate 93616371a99SBarry Smith 93716371a99SBarry Smith Notes: 93816371a99SBarry Smith See the chapter in the users manual on performance for more details 93916371a99SBarry Smith 94016371a99SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 94116371a99SBarry Smith 94216371a99SBarry Smith This is a MACRO not a function because it uses a bunch of variables private to the MatPreallocation.... routines. 94316371a99SBarry Smith 94416371a99SBarry Smith Concepts: preallocation^Matrix 94516371a99SBarry Smith 94616371a99SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(), 94716371a99SBarry Smith MatPreallocateSymmetricInitialize(), MatPreallocateSymmetricSetLocal() 94816371a99SBarry Smith M*/ 94916371a99SBarry 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);} 95016371a99SBarry Smith 95116371a99SBarry Smith 95216371a99SBarry Smith /*MC 9530d2a7ff2SSatish Balay MatPreallocateFinalize - Ends the block of code that will count the number of nonzeros per 954bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 955bd481603SBarry Smith 956bd481603SBarry Smith Synopsis: 957c1ac3661SBarry Smith PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz) 958bd481603SBarry Smith 959bd481603SBarry Smith Collective on MPI_Comm 960bd481603SBarry Smith 961bd481603SBarry Smith Input Parameters: 96216371a99SBarry Smith + dnz - the array that was be passed to the matrix preallocation routines 963bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 964bd481603SBarry Smith 965bd481603SBarry Smith 966bd481603SBarry Smith Level: intermediate 967bd481603SBarry Smith 968bd481603SBarry Smith Notes: 969bd481603SBarry Smith See the chapter in the users manual on performance for more details 970bd481603SBarry Smith 971bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 972bd481603SBarry Smith 9731620fd73SBarry Smith This is a MACRO not a function because it closes the { started in MatPreallocateInitialize(). 9741620fd73SBarry Smith 975bd481603SBarry Smith Concepts: preallocation^Matrix 976bd481603SBarry Smith 977bd481603SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(), 978bd481603SBarry Smith MatPreallocateSymmetricInitialize(), MatPreallocateSymmetricSetLocal() 979bd481603SBarry Smith M*/ 980a89bc211SBarry Smith #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree2(dnz,onz);CHKERRQ(_4_ierr);} 9817c922b88SBarry Smith 982bd481603SBarry Smith 983bd481603SBarry Smith 9847b80b807SBarry Smith /* Routines unique to particular data structures */ 985be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellGetContext(Mat,void **); 986ba966639SSatish Balay PetscPolymorphicFunction(MatShellGetContext,(Mat A),(A,&t),void*,t) 987698d4c6aSKris Buschelman 988be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInodeAdjustForInodes(Mat,IS*,IS*); 989be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *); 990698d4c6aSKris Buschelman 991be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetColumnIndices(Mat,PetscInt[]); 992be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetColumnIndices(Mat,PetscInt[]); 993be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 994c75a6043SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 995c75a6043SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 9967b80b807SBarry Smith 997a96a251dSBarry Smith #define MAT_SKIP_ALLOCATION -4 998a96a251dSBarry Smith 999be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]); 1000ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz)) 1001be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]); 1002ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqSBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz)) 1003be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]); 1004ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqAIJSetPreallocation,(Mat A,const PetscInt nnz[]),(A,0,nnz)) 1005273d9f13SBarry Smith 1006be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 1007ba966639SSatish Balay PetscPolymorphicSubroutine(MatMPIBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[],const PetscInt onz[]),(A,bs,0,nnz,0,onz)) 1008be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 1009be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 101053803f26SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []); 1011725b52f3SLisandro Dalcin EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 1012be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]); 1013be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 1014be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]); 1015be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIDenseSetPreallocation(Mat,PetscScalar[]); 1016284134d9SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation(Mat,PetscScalar[]); 1017be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,PetscInt*[]); 1018be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,PetscInt*[]); 1019be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAdicSetLocalFunction(Mat,void (*)(void)); 1020273d9f13SBarry Smith 1021be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetLDA(Mat,PetscInt); 1022ab92ecdeSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDenseGetLocalMatrix(Mat,Mat*); 10231b807ce4Svictorle 1024be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStoreValues(Mat); 1025be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRetrieveValues(Mat); 10262e8a6d31SBarry Smith 1027be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDAADSetCtx(Mat,void*); 10283a7fca6bSBarry Smith 10297b80b807SBarry Smith /* 10307b80b807SBarry Smith These routines are not usually accessed directly, rather solving is 103194b7f48cSBarry Smith done through the KSP and PC interfaces. 10327b80b807SBarry Smith */ 10337b80b807SBarry Smith 1034d9274352SBarry Smith /*E 1035d9274352SBarry Smith MatOrderingType - String with the name of a PETSc matrix ordering or the creation function 1036d9274352SBarry Smith with an optional dynamic library name, for example 1037d9274352SBarry Smith http://www.mcs.anl.gov/petsc/lib.a:orderingcreate() 1038d9274352SBarry Smith 1039d9274352SBarry Smith Level: beginner 1040d9274352SBarry Smith 1041e5a9bf91SBarry Smith Cannot use const because the PC objects manipulate the string 1042e5a9bf91SBarry Smith 1043d9274352SBarry Smith .seealso: MatGetOrdering() 1044d9274352SBarry Smith E*/ 10453eaad2caSSatish Balay #define MatOrderingType char* 1046b12f92e5SBarry Smith #define MATORDERING_NATURAL "natural" 1047b12f92e5SBarry Smith #define MATORDERING_ND "nd" 1048b12f92e5SBarry Smith #define MATORDERING_1WD "1wd" 1049b12f92e5SBarry Smith #define MATORDERING_RCM "rcm" 1050b12f92e5SBarry Smith #define MATORDERING_QMD "qmd" 1051b12f92e5SBarry Smith #define MATORDERING_ROWLENGTH "rowlength" 1052e106eecfSBarry Smith #define MATORDERING_DSC_ND "dsc_nd" /* these three are only for DSCPACK, see its documentation for details */ 105362152c8bSBarry Smith #define MATORDERING_DSC_MMD "dsc_mmd" 105462152c8bSBarry Smith #define MATORDERING_DSC_MDF "dsc_mdf" 1055c06d978dSMatthew Knepley #define MATORDERING_CONSTRAINED "constrained" 1056c06d978dSMatthew Knepley #define MATORDERING_IDENTITY "identity" 1057c06d978dSMatthew Knepley #define MATORDERING_REVERSE "reverse" 10588fa24674SBarry Smith #define MATORDERING_FLOW "flow" 105999671722SJed Brown #define MATORDERING_AMD "amd" 1060b12f92e5SBarry Smith 1061f2416520SSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOrdering(Mat,const MatOrderingType,IS*,IS*); 1062be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOrderingList(PetscFList *list); 1063f2416520SSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,const MatOrderingType,IS*,IS*)); 106430de9b25SBarry Smith 106530de9b25SBarry Smith /*MC 10661890ba74SBarry Smith MatOrderingRegisterDynamic - Adds a new sparse matrix ordering to the matrix package. 106730de9b25SBarry Smith 106830de9b25SBarry Smith Synopsis: 10691890ba74SBarry Smith PetscErrorCode MatOrderingRegisterDynamic(const char *name_ordering,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatOrdering)) 107030de9b25SBarry Smith 107130de9b25SBarry Smith Not Collective 107230de9b25SBarry Smith 107330de9b25SBarry Smith Input Parameters: 107430de9b25SBarry Smith + sname - name of ordering (for example MATORDERING_ND) 107530de9b25SBarry Smith . path - location of library where creation routine is 107630de9b25SBarry Smith . name - name of function that creates the ordering type,a string 107730de9b25SBarry Smith - function - function pointer that creates the ordering 107830de9b25SBarry Smith 107930de9b25SBarry Smith Level: developer 108030de9b25SBarry Smith 108130de9b25SBarry Smith If dynamic libraries are used, then the fourth input argument (function) 108230de9b25SBarry Smith is ignored. 108330de9b25SBarry Smith 108430de9b25SBarry Smith Sample usage: 108530de9b25SBarry Smith .vb 108630de9b25SBarry Smith MatOrderingRegisterDynamic("my_order",/home/username/my_lib/lib/libO/solaris/mylib.a, 108730de9b25SBarry Smith "MyOrder",MyOrder); 108830de9b25SBarry Smith .ve 108930de9b25SBarry Smith 109030de9b25SBarry Smith Then, your partitioner can be chosen with the procedural interface via 109130de9b25SBarry Smith $ MatOrderingSetType(part,"my_order) 109230de9b25SBarry Smith or at runtime via the option 10932401956bSBarry Smith $ -pc_factor_mat_ordering_type my_order 109430de9b25SBarry Smith 1095ab901514SBarry Smith ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values. 109630de9b25SBarry Smith 109730de9b25SBarry Smith .keywords: matrix, ordering, register 109830de9b25SBarry Smith 109930de9b25SBarry Smith .seealso: MatOrderingRegisterDestroy(), MatOrderingRegisterAll() 110030de9b25SBarry Smith M*/ 1101aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 1102f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,0) 1103b12f92e5SBarry Smith #else 1104f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,d) 1105b12f92e5SBarry Smith #endif 110630de9b25SBarry Smith 1107be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegisterDestroy(void); 1108be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegisterAll(const char[]); 11092bad1931SBarry Smith extern PetscTruth MatOrderingRegisterAllCalled; 1110b0a32e0cSBarry Smith extern PetscFList MatOrderingList; 1111d4fbbf0eSBarry Smith 1112be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS); 1113a2ce50c7SBarry Smith 1114d91e6319SBarry Smith /*S 1115d90ac83dSHong Zhang MatFactorShiftType - Numeric Shift. 1116d90ac83dSHong Zhang 1117d90ac83dSHong Zhang Level: beginner 1118d90ac83dSHong Zhang 1119d90ac83dSHong Zhang S*/ 1120d90ac83dSHong Zhang typedef enum {MAT_SHIFT_NONE,MAT_SHIFT_NONZERO,MAT_SHIFT_POSITIVE_DEFINITE,MAT_SHIFT_INBLOCKS} MatFactorShiftType; 1121d90ac83dSHong Zhang extern const char *MatFactorShiftTypes[]; 1122d90ac83dSHong Zhang 1123d90ac83dSHong Zhang /*S 11242401956bSBarry Smith MatFactorInfo - Data passed into the matrix factorization routines 1125b00f7748SHong Zhang 112661cad860SBarry Smith In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE, that is use 112761cad860SBarry Smith $ MatFactorInfo info(MAT_FACTORINFO_SIZE) 1128b00f7748SHong Zhang 112915e8a5b3SHong Zhang Notes: These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC. 1130b00f7748SHong Zhang 113161cad860SBarry Smith You can use MatFactorInfoInitialize() to set default values. 113261cad860SBarry Smith 1133b00f7748SHong Zhang Level: developer 1134b00f7748SHong Zhang 1135d7d82daaSBarry Smith .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(), 1136d7d82daaSBarry Smith MatFactorInfoInitialize() 1137b00f7748SHong Zhang 1138b00f7748SHong Zhang S*/ 1139b00f7748SHong Zhang typedef struct { 114015e8a5b3SHong Zhang PetscReal diagonal_fill; /* force diagonal to fill in if initially not filled */ 114185317021SBarry Smith PetscReal usedt; 114215e8a5b3SHong Zhang PetscReal dt; /* drop tolerance */ 1143b00f7748SHong Zhang PetscReal dtcol; /* tolerance for pivoting */ 114415e8a5b3SHong Zhang PetscReal dtcount; /* maximum nonzeros to be allowed per row */ 114567e28bfeSBarry Smith PetscReal fill; /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */ 1146348344bbSBarry Smith PetscReal levels; /* ICC/ILU(levels) */ 1147bcd9e38bSBarry Smith PetscReal pivotinblocks; /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0 1148bcd9e38bSBarry Smith factorization may be faster if do not pivot */ 114915e8a5b3SHong Zhang PetscReal zeropivot; /* pivot is called zero if less than this */ 1150f4db908eSBarry Smith PetscReal shifttype; /* type of shift added to matrix factor to prevent zero pivots */ 1151f4db908eSBarry Smith PetscReal shiftamount; /* how large the shift is */ 115215e8a5b3SHong Zhang } MatFactorInfo; 1153ffa6d0a5SLois Curfman McInnes 1154be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFactorInfoInitialize(MatFactorInfo*); 11550481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactor(Mat,IS,const MatFactorInfo*); 11560481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*); 11570481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactorNumeric(Mat,Mat,const MatFactorInfo*); 11580481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactor(Mat,IS,IS,const MatFactorInfo*); 11590481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUFactor(Mat,IS,IS,const MatFactorInfo*); 11600481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*); 11610481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*); 11620481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatICCFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*); 11630481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatICCFactor(Mat,IS,const MatFactorInfo*); 11640481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactorNumeric(Mat,Mat,const MatFactorInfo*); 1165be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*); 1166be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolve(Mat,Vec,Vec); 1167be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatForwardSolve(Mat,Vec,Vec); 1168be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatBackwardSolve(Mat,Vec,Vec); 1169be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveAdd(Mat,Vec,Vec,Vec); 1170be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveTranspose(Mat,Vec,Vec); 1171be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveTransposeAdd(Mat,Vec,Vec,Vec); 1172be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolves(Mat,Vecs,Vecs); 11738ed539a5SBarry Smith 1174be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUnfactored(Mat); 1175bb5a7306SBarry Smith 1176d91e6319SBarry Smith /*E 1177d91e6319SBarry Smith MatSORType - What type of (S)SOR to perform 1178bb1eb677SSatish Balay 1179d91e6319SBarry Smith Level: beginner 1180d91e6319SBarry Smith 1181d9274352SBarry Smith May be bitwise ORd together 1182d9274352SBarry Smith 1183d91e6319SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 1184d91e6319SBarry Smith 11854e7234bfSBarry Smith MatSORType may be bitwise ORd together, so do not change the numbers 11864e7234bfSBarry Smith 118741f059aeSBarry Smith .seealso: MatSOR() 1188d91e6319SBarry Smith E*/ 1189ee50ffe9SBarry Smith typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3, 1190ee50ffe9SBarry Smith SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8, 1191ee50ffe9SBarry Smith SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16, 119284cb2905SBarry Smith SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType; 119341f059aeSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSOR(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec); 11948ed539a5SBarry Smith 1195d4fbbf0eSBarry Smith /* 1196639f9d9dSBarry Smith These routines are for efficiently computing Jacobians via finite differences. 1197639f9d9dSBarry Smith */ 1198b12f92e5SBarry Smith 1199d9274352SBarry Smith /*E 1200d9274352SBarry Smith MatColoringType - String with the name of a PETSc matrix coloring or the creation function 1201d9274352SBarry Smith with an optional dynamic library name, for example 1202d9274352SBarry Smith http://www.mcs.anl.gov/petsc/lib.a:coloringcreate() 1203d9274352SBarry Smith 1204d9274352SBarry Smith Level: beginner 1205d9274352SBarry Smith 1206d9274352SBarry Smith .seealso: MatGetColoring() 1207d9274352SBarry Smith E*/ 1208a313700dSBarry Smith #define MatColoringType char* 1209b12f92e5SBarry Smith #define MATCOLORING_NATURAL "natural" 1210b12f92e5SBarry Smith #define MATCOLORING_SL "sl" 1211b12f92e5SBarry Smith #define MATCOLORING_LF "lf" 1212b12f92e5SBarry Smith #define MATCOLORING_ID "id" 1213b12f92e5SBarry Smith 1214ab970d9bSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColoring(Mat,const MatColoringType,ISColoring*); 12152e90c967SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,MatColoringType,ISColoring *)); 121630de9b25SBarry Smith 121730de9b25SBarry Smith /*MC 121830de9b25SBarry Smith MatColoringRegisterDynamic - Adds a new sparse matrix coloring to the 121930de9b25SBarry Smith matrix package. 122030de9b25SBarry Smith 122130de9b25SBarry Smith Synopsis: 12221890ba74SBarry Smith PetscErrorCode MatColoringRegisterDynamic(const char *name_coloring,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatColoring)) 122330de9b25SBarry Smith 122430de9b25SBarry Smith Not Collective 122530de9b25SBarry Smith 122630de9b25SBarry Smith Input Parameters: 122730de9b25SBarry Smith + sname - name of Coloring (for example MATCOLORING_SL) 122830de9b25SBarry Smith . path - location of library where creation routine is 122930de9b25SBarry Smith . name - name of function that creates the Coloring type, a string 123030de9b25SBarry Smith - function - function pointer that creates the coloring 123130de9b25SBarry Smith 123230de9b25SBarry Smith Level: developer 123330de9b25SBarry Smith 123430de9b25SBarry Smith If dynamic libraries are used, then the fourth input argument (function) 123530de9b25SBarry Smith is ignored. 123630de9b25SBarry Smith 123730de9b25SBarry Smith Sample usage: 123830de9b25SBarry Smith .vb 123930de9b25SBarry Smith MatColoringRegisterDynamic("my_color",/home/username/my_lib/lib/libO/solaris/mylib.a, 124030de9b25SBarry Smith "MyColor",MyColor); 124130de9b25SBarry Smith .ve 124230de9b25SBarry Smith 124330de9b25SBarry Smith Then, your partitioner can be chosen with the procedural interface via 124430de9b25SBarry Smith $ MatColoringSetType(part,"my_color") 124530de9b25SBarry Smith or at runtime via the option 124630de9b25SBarry Smith $ -mat_coloring_type my_color 124730de9b25SBarry Smith 1248ab901514SBarry Smith $PETSC_ARCH occuring in pathname will be replaced with appropriate values. 124930de9b25SBarry Smith 125030de9b25SBarry Smith .keywords: matrix, Coloring, register 125130de9b25SBarry Smith 125230de9b25SBarry Smith .seealso: MatColoringRegisterDestroy(), MatColoringRegisterAll() 125330de9b25SBarry Smith M*/ 1254aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 1255f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,0) 1256b12f92e5SBarry Smith #else 1257f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,d) 1258b12f92e5SBarry Smith #endif 125930de9b25SBarry Smith 12602bad1931SBarry Smith extern PetscTruth MatColoringRegisterAllCalled; 1261f1144a30SSatish Balay 1262be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegisterAll(const char[]); 1263be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegisterDestroy(void); 1264be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*); 1265639f9d9dSBarry Smith 1266d9274352SBarry Smith /*S 1267d9274352SBarry Smith MatFDColoring - Object for computing a sparse Jacobian via finite differences 1268d9274352SBarry Smith and coloring 1269639f9d9dSBarry Smith 1270d9274352SBarry Smith Level: beginner 1271d9274352SBarry Smith 1272d9274352SBarry Smith Concepts: coloring, sparse Jacobian, finite differences 1273d9274352SBarry Smith 1274d9274352SBarry Smith .seealso: MatFDColoringCreate() 1275d9274352SBarry Smith S*/ 1276e2a1c21fSSatish Balay typedef struct _p_MatFDColoring* MatFDColoring; 1277639f9d9dSBarry Smith 1278be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringCreate(Mat,ISColoring,MatFDColoring *); 1279be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringDestroy(MatFDColoring); 1280be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringView(MatFDColoring,PetscViewer); 1281be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*); 12824a9d489dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**); 1283be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal); 1284be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFromOptions(MatFDColoring); 1285be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApply(Mat,MatFDColoring,Vec,MatStructure*,void *); 1286be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApplyTS(Mat,MatFDColoring,PetscReal,Vec,MatStructure*,void *); 1287be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetF(MatFDColoring,Vec); 1288be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,PetscInt*[]); 1289639f9d9dSBarry Smith /* 12900752156aSBarry Smith These routines are for partitioning matrices: currently used only 12913eda8832SBarry Smith for adjacency matrix, MatCreateMPIAdj(). 12920752156aSBarry Smith */ 1293ca161407SBarry Smith 1294d9274352SBarry Smith /*S 1295d9274352SBarry Smith MatPartitioning - Object for managing the partitioning of a matrix or graph 1296d9274352SBarry Smith 1297d9274352SBarry Smith Level: beginner 1298d9274352SBarry Smith 1299d9274352SBarry Smith Concepts: partitioning 1300d9274352SBarry Smith 1301743a1026Svictorle .seealso: MatPartitioningCreate(), MatPartitioningType 1302d9274352SBarry Smith S*/ 130391e9ee9fSBarry Smith typedef struct _p_MatPartitioning* MatPartitioning; 1304d9274352SBarry Smith 1305d9274352SBarry Smith /*E 13065bc6e7ccSvictorle MatPartitioningType - String with the name of a PETSc matrix partitioning or the creation function 1307d9274352SBarry Smith with an optional dynamic library name, for example 1308d9274352SBarry Smith http://www.mcs.anl.gov/petsc/lib.a:partitioningcreate() 1309d9274352SBarry Smith 1310d9274352SBarry Smith Level: beginner 1311d9274352SBarry Smith 1312b547a13bSvictorle .seealso: MatPartitioningCreate(), MatPartitioning 1313d9274352SBarry Smith E*/ 1314a313700dSBarry Smith #define MatPartitioningType char* 13158ba1e511SMatthew Knepley #define MAT_PARTITIONING_CURRENT "current" 1316c2ec2a73SHong Zhang #define MAT_PARTITIONING_SQUARE "square" 13178ba1e511SMatthew Knepley #define MAT_PARTITIONING_PARMETIS "parmetis" 131871306c60Sjroman #define MAT_PARTITIONING_CHACO "chaco" 131971306c60Sjroman #define MAT_PARTITIONING_JOSTLE "jostle" 132071306c60Sjroman #define MAT_PARTITIONING_PARTY "party" 132171306c60Sjroman #define MAT_PARTITIONING_SCOTCH "scotch" 132271306c60Sjroman 1323ca161407SBarry Smith 1324be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningCreate(MPI_Comm,MatPartitioning*); 1325a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetType(MatPartitioning,const MatPartitioningType); 1326be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetNParts(MatPartitioning,PetscInt); 1327be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetAdjacency(MatPartitioning,Mat); 1328be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]); 1329be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []); 1330be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningApply(MatPartitioning,IS*); 1331be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningDestroy(MatPartitioning); 13322aabb6bbSBarry Smith 1333be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatPartitioning)); 133430de9b25SBarry Smith 133530de9b25SBarry Smith /*MC 133630de9b25SBarry Smith MatPartitioningRegisterDynamic - Adds a new sparse matrix partitioning to the 133730de9b25SBarry Smith matrix package. 133830de9b25SBarry Smith 133930de9b25SBarry Smith Synopsis: 13401890ba74SBarry Smith PetscErrorCode MatPartitioningRegisterDynamic(const char *name_partitioning,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatPartitioning)) 134130de9b25SBarry Smith 134230de9b25SBarry Smith Not Collective 134330de9b25SBarry Smith 134430de9b25SBarry Smith Input Parameters: 134530de9b25SBarry Smith + sname - name of partitioning (for example MAT_PARTITIONING_CURRENT) or parmetis 134630de9b25SBarry Smith . path - location of library where creation routine is 134730de9b25SBarry Smith . name - name of function that creates the partitioning type, a string 134830de9b25SBarry Smith - function - function pointer that creates the partitioning type 134930de9b25SBarry Smith 135030de9b25SBarry Smith Level: developer 135130de9b25SBarry Smith 135230de9b25SBarry Smith If dynamic libraries are used, then the fourth input argument (function) 135330de9b25SBarry Smith is ignored. 135430de9b25SBarry Smith 135530de9b25SBarry Smith Sample usage: 135630de9b25SBarry Smith .vb 135730de9b25SBarry Smith MatPartitioningRegisterDynamic("my_part",/home/username/my_lib/lib/libO/solaris/mylib.a, 135830de9b25SBarry Smith "MyPartCreate",MyPartCreate); 135930de9b25SBarry Smith .ve 136030de9b25SBarry Smith 136130de9b25SBarry Smith Then, your partitioner can be chosen with the procedural interface via 136230de9b25SBarry Smith $ MatPartitioningSetType(part,"my_part") 136330de9b25SBarry Smith or at runtime via the option 136430de9b25SBarry Smith $ -mat_partitioning_type my_part 136530de9b25SBarry Smith 1366ab901514SBarry Smith $PETSC_ARCH occuring in pathname will be replaced with appropriate values. 136730de9b25SBarry Smith 136830de9b25SBarry Smith .keywords: matrix, partitioning, register 136930de9b25SBarry Smith 137030de9b25SBarry Smith .seealso: MatPartitioningRegisterDestroy(), MatPartitioningRegisterAll() 137130de9b25SBarry Smith M*/ 1372aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 1373f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,0) 13742aabb6bbSBarry Smith #else 1375f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,d) 13762aabb6bbSBarry Smith #endif 13772aabb6bbSBarry Smith 13782bad1931SBarry Smith extern PetscTruth MatPartitioningRegisterAllCalled; 1379f1144a30SSatish Balay 1380be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegisterAll(const char[]); 1381be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegisterDestroy(void); 13822bad1931SBarry Smith 1383be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningView(MatPartitioning,PetscViewer); 1384be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetFromOptions(MatPartitioning); 1385a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningGetType(MatPartitioning,const MatPartitioningType*); 1386ca161407SBarry Smith 1387be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningParmetisSetCoarseSequential(MatPartitioning); 13880e8af3caSMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningParmetisGetEdgeCut(MatPartitioning, PetscInt *); 13890752156aSBarry Smith 1390be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningJostleSetCoarseLevel(MatPartitioning,PetscReal); 1391be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningJostleSetCoarseSequential(MatPartitioning); 139271306c60Sjroman 1393954b3ec6SBarry Smith typedef enum {MP_CHACO_MULTILEVEL_KL,MP_CHACO_SPECTRAL,MP_CHACO_LINEAR,MP_CHACO_RANDOM, MP_CHACO_SCATTERED} MPChacoGlobalType; 1394be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetGlobal(MatPartitioning, MPChacoGlobalType); 139571306c60Sjroman typedef enum { MP_CHACO_KERNIGHAN_LIN, MP_CHACO_NONE } MPChacoLocalType; 1396be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetLocal(MatPartitioning, MPChacoLocalType); 1397be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal); 139871306c60Sjroman typedef enum { MP_CHACO_LANCZOS, MP_CHACO_RQI_SYMMLQ } MPChacoEigenType; 1399be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType); 1400be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenTol(MatPartitioning, PetscReal); 1401be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenNumber(MatPartitioning, PetscInt); 140271306c60Sjroman 140371306c60Sjroman #define MP_PARTY_OPT "opt" 140471306c60Sjroman #define MP_PARTY_LIN "lin" 140571306c60Sjroman #define MP_PARTY_SCA "sca" 140671306c60Sjroman #define MP_PARTY_RAN "ran" 140771306c60Sjroman #define MP_PARTY_GBF "gbf" 140871306c60Sjroman #define MP_PARTY_GCF "gcf" 140971306c60Sjroman #define MP_PARTY_BUB "bub" 141071306c60Sjroman #define MP_PARTY_DEF "def" 1411be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetGlobal(MatPartitioning, const char*); 141271306c60Sjroman #define MP_PARTY_HELPFUL_SETS "hs" 141371306c60Sjroman #define MP_PARTY_KERNIGHAN_LIN "kl" 141471306c60Sjroman #define MP_PARTY_NONE "no" 1415be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetLocal(MatPartitioning, const char*); 1416be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal); 1417be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetBipart(MatPartitioning,PetscTruth); 1418be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscTruth); 141971306c60Sjroman 142071306c60Sjroman typedef enum { MP_SCOTCH_GREEDY, MP_SCOTCH_GPS, MP_SCOTCH_GR_GPS } MPScotchGlobalType; 1421be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetArch(MatPartitioning,const char*); 1422be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetMultilevel(MatPartitioning); 1423be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetGlobal(MatPartitioning,MPScotchGlobalType); 1424be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetCoarseLevel(MatPartitioning,PetscReal); 1425be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetHostList(MatPartitioning,const char*); 142671306c60Sjroman typedef enum { MP_SCOTCH_KERNIGHAN_LIN, MP_SCOTCH_NONE } MPScotchLocalType; 1427be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetLocal(MatPartitioning,MPScotchLocalType); 1428be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetMapping(MatPartitioning); 1429be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetStrategy(MatPartitioning,char*); 143071306c60Sjroman 1431591294e4SBarry Smith EXTERN PetscErrorCode MatMeshToVertexGraph(Mat,PetscInt,Mat*); 143214196391SBarry Smith EXTERN PetscErrorCode MatMeshToCellGraph(Mat,PetscInt,Mat*); 1433591294e4SBarry Smith 14340752156aSBarry Smith /* 14350a835dfdSSatish Balay If you add entries here you must also add them to finclude/petscmat.h 1436d4fbbf0eSBarry Smith */ 14371c1c02c0SLois Curfman McInnes typedef enum { MATOP_SET_VALUES=0, 14381c1c02c0SLois Curfman McInnes MATOP_GET_ROW=1, 14391c1c02c0SLois Curfman McInnes MATOP_RESTORE_ROW=2, 14401c1c02c0SLois Curfman McInnes MATOP_MULT=3, 14411c1c02c0SLois Curfman McInnes MATOP_MULT_ADD=4, 14427c922b88SBarry Smith MATOP_MULT_TRANSPOSE=5, 14437c922b88SBarry Smith MATOP_MULT_TRANSPOSE_ADD=6, 14441c1c02c0SLois Curfman McInnes MATOP_SOLVE=7, 14451c1c02c0SLois Curfman McInnes MATOP_SOLVE_ADD=8, 14467c922b88SBarry Smith MATOP_SOLVE_TRANSPOSE=9, 14477c922b88SBarry Smith MATOP_SOLVE_TRANSPOSE_ADD=10, 14481c1c02c0SLois Curfman McInnes MATOP_LUFACTOR=11, 14491c1c02c0SLois Curfman McInnes MATOP_CHOLESKYFACTOR=12, 1450a714c06dSBarry Smith MATOP_SOR=13, 14511c1c02c0SLois Curfman McInnes MATOP_TRANSPOSE=14, 14521c1c02c0SLois Curfman McInnes MATOP_GETINFO=15, 14531c1c02c0SLois Curfman McInnes MATOP_EQUAL=16, 14541c1c02c0SLois Curfman McInnes MATOP_GET_DIAGONAL=17, 14551c1c02c0SLois Curfman McInnes MATOP_DIAGONAL_SCALE=18, 14561c1c02c0SLois Curfman McInnes MATOP_NORM=19, 14571c1c02c0SLois Curfman McInnes MATOP_ASSEMBLY_BEGIN=20, 14581c1c02c0SLois Curfman McInnes MATOP_ASSEMBLY_END=21, 1459d519adbfSMatthew Knepley MATOP_SET_OPTION=22, 1460d519adbfSMatthew Knepley MATOP_ZERO_ENTRIES=23, 1461d519adbfSMatthew Knepley MATOP_ZERO_ROWS=24, 1462d519adbfSMatthew Knepley MATOP_LUFACTOR_SYMBOLIC=25, 1463d519adbfSMatthew Knepley MATOP_LUFACTOR_NUMERIC=26, 1464d519adbfSMatthew Knepley MATOP_CHOLESKY_FACTOR_SYMBOLIC=27, 1465d519adbfSMatthew Knepley MATOP_CHOLESKY_FACTOR_NUMERIC=28, 1466d519adbfSMatthew Knepley MATOP_SETUP_PREALLOCATION=29, 1467d519adbfSMatthew Knepley MATOP_ILUFACTOR_SYMBOLIC=30, 1468d519adbfSMatthew Knepley MATOP_ICCFACTOR_SYMBOLIC=31, 1469d519adbfSMatthew Knepley MATOP_GET_ARRAY=32, 1470d519adbfSMatthew Knepley MATOP_RESTORE_ARRAY=33, 1471d519adbfSMatthew Knepley MATOP_DUPLICATE=34, 1472d519adbfSMatthew Knepley MATOP_FORWARD_SOLVE=35, 1473d519adbfSMatthew Knepley MATOP_BACKWARD_SOLVE=36, 1474d519adbfSMatthew Knepley MATOP_ILUFACTOR=37, 1475d519adbfSMatthew Knepley MATOP_ICCFACTOR=38, 1476d519adbfSMatthew Knepley MATOP_AXPY=39, 1477d519adbfSMatthew Knepley MATOP_GET_SUBMATRICES=40, 1478d519adbfSMatthew Knepley MATOP_INCREASE_OVERLAP=41, 1479d519adbfSMatthew Knepley MATOP_GET_VALUES=42, 1480d519adbfSMatthew Knepley MATOP_COPY=43, 1481d519adbfSMatthew Knepley MATOP_GET_ROW_MAX=44, 1482d519adbfSMatthew Knepley MATOP_SCALE=45, 1483d519adbfSMatthew Knepley MATOP_SHIFT=46, 148435153367SBarry Smith MATOP_DIAGONAL_SET=47, 1485d519adbfSMatthew Knepley MATOP_ILUDT_FACTOR=48, 1486d519adbfSMatthew Knepley MATOP_SET_BLOCK_SIZE=49, 1487d519adbfSMatthew Knepley MATOP_GET_ROW_IJ=50, 1488d519adbfSMatthew Knepley MATOP_RESTORE_ROW_IJ=51, 1489d519adbfSMatthew Knepley MATOP_GET_COLUMN_IJ=52, 1490d519adbfSMatthew Knepley MATOP_RESTORE_COLUMN_IJ=53, 1491d519adbfSMatthew Knepley MATOP_FDCOLORING_CREATE=54, 1492d519adbfSMatthew Knepley MATOP_COLORING_PATCH=55, 1493d519adbfSMatthew Knepley MATOP_SET_UNFACTORED=56, 1494d519adbfSMatthew Knepley MATOP_PERMUTE=57, 1495d519adbfSMatthew Knepley MATOP_SET_VALUES_BLOCKED=58, 1496d519adbfSMatthew Knepley MATOP_GET_SUBMATRIX=59, 1497d519adbfSMatthew Knepley MATOP_DESTROY=60, 1498d519adbfSMatthew Knepley MATOP_VIEW=61, 1499d519adbfSMatthew Knepley MATOP_CONVERT_FROM=62, 1500d519adbfSMatthew Knepley MATOP_USE_SCALED_FORM=63, 1501d519adbfSMatthew Knepley MATOP_SCALE_SYSTEM=64, 1502d519adbfSMatthew Knepley MATOP_UNSCALE_SYSTEM=65, 1503d519adbfSMatthew Knepley MATOP_SET_LOCAL_TO_GLOBAL_MAP=66, 1504d519adbfSMatthew Knepley MATOP_SET_VALUES_LOCAL=67, 1505d519adbfSMatthew Knepley MATOP_ZERO_ROWS_LOCAL=68, 1506d519adbfSMatthew Knepley MATOP_GET_ROW_MAX_ABS=69, 1507d519adbfSMatthew Knepley MATOP_GET_ROW_MIN_ABS=70, 1508d519adbfSMatthew Knepley MATOP_CONVERT=71, 1509d519adbfSMatthew Knepley MATOP_SET_COLORING=72, 1510d519adbfSMatthew Knepley MATOP_SET_VALUES_ADIC=73, 1511d519adbfSMatthew Knepley MATOP_SET_VALUES_ADIFOR=74, 1512d519adbfSMatthew Knepley MATOP_FD_COLORING_APPLY=75, 1513d519adbfSMatthew Knepley MATOP_SET_FROM_OPTIONS=76, 1514d519adbfSMatthew Knepley MATOP_MULT_CON=77, 1515d519adbfSMatthew Knepley MATOP_MULT_TRANSPOSE_CON=78, 1516d519adbfSMatthew Knepley MATOP_PERMUTE_SPARSIFY=79, 1517d519adbfSMatthew Knepley MATOP_MULT_MULTIPLE=80, 1518d519adbfSMatthew Knepley MATOP_SOLVE_MULTIPLE=81, 1519d519adbfSMatthew Knepley MATOP_GET_INERTIA=82, 1520d519adbfSMatthew Knepley MATOP_LOAD=83, 1521d519adbfSMatthew Knepley MATOP_IS_SYMMETRIC=84, 1522d519adbfSMatthew Knepley MATOP_IS_HERMITIAN=85, 1523d519adbfSMatthew Knepley MATOP_IS_STRUCTURALLY_SYMMETRIC=86, 152441f059aeSBarry Smith MATOP_DUMMY=87, 1525d519adbfSMatthew Knepley MATOP_GET_VECS=88, 1526d519adbfSMatthew Knepley MATOP_MAT_MULT=89, 1527d519adbfSMatthew Knepley MATOP_MAT_MULT_SYMBOLIC=90, 1528d519adbfSMatthew Knepley MATOP_MAT_MULT_NUMERIC=91, 1529d519adbfSMatthew Knepley MATOP_PTAP=92, 1530d519adbfSMatthew Knepley MATOP_PTAP_SYMBOLIC=93, 1531d519adbfSMatthew Knepley MATOP_PTAP_NUMERIC=94, 1532d519adbfSMatthew Knepley MATOP_MAT_MULTTRANSPOSE=95, 15331763ac4eSSatish Balay MATOP_MAT_MULTTRANSPOSE_SYM=96, 15341763ac4eSSatish Balay MATOP_MAT_MULTTRANSPOSE_NUM=97, 1535d519adbfSMatthew Knepley MATOP_PTAP_SYMBOLIC_SEQAIJ=98, 1536d519adbfSMatthew Knepley MATOP_PTAP_NUMERIC_SEQAIJ=99, 1537d519adbfSMatthew Knepley MATOP_PTAP_SYMBOLIC_MPIAIJ=100, 1538d519adbfSMatthew Knepley MATOP_PTAP_NUMERIC_MPIAIJ=101, 1539d519adbfSMatthew Knepley MATOP_CONJUGATE=102, 1540d519adbfSMatthew Knepley MATOP_SET_SIZES=103, 1541d519adbfSMatthew Knepley MATOP_SET_VALUES_ROW=104, 1542d519adbfSMatthew Knepley MATOP_REAL_PART=105, 1543d519adbfSMatthew Knepley MATOP_IMAG_PART=106, 1544d519adbfSMatthew Knepley MATOP_GET_ROW_UTRIANGULAR=107, 1545d519adbfSMatthew Knepley MATOP_RESTORE_ROW_UTRIANGULAR=108, 1546d519adbfSMatthew Knepley MATOP_MATSOLVE=109, 1547d519adbfSMatthew Knepley MATOP_GET_REDUNDANTMATRIX=110, 1548d519adbfSMatthew Knepley MATOP_GET_ROW_MIN=111, 1549d519adbfSMatthew Knepley MATOP_GET_COLUMN_VEC=112, 1550d519adbfSMatthew Knepley MATOP_MISSING_DIAGONAL=113, 1551d519adbfSMatthew Knepley MATOP_MATGETSEQNONZEROSTRUCTURE=114, 155289c6957cSBarry Smith MATOP_CREATE=115, 155389c6957cSBarry Smith MATOP_GET_GHOSTS=116, 1554eeffb40dSHong Zhang MATOP_MULT_DIAGONAL_BLOCK=119, 1555547795f9SHong Zhang MATOP_HERMITIANTRANSPOSE=120, 1556547795f9SHong Zhang MATOP_MULTHERMITIANTRANSPOSE=121, 1557547795f9SHong Zhang MATOP_MULTHERMITIANTRANSPOSEADD=122 1558fae171e0SBarry Smith } MatOperation; 1559be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatHasOperation(Mat,MatOperation,PetscTruth*); 1560be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellSetOperation(Mat,MatOperation,void(*)(void)); 1561be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellGetOperation(Mat,MatOperation,void(**)(void)); 1562be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellSetContext(Mat,void*); 1563112a2221SBarry Smith 156490ace30eSBarry Smith /* 156590ace30eSBarry Smith Codes for matrices stored on disk. By default they are 156690ace30eSBarry Smith stored in a universal format. By changing the format with 15677973ad04SBarry Smith PetscViewerSetFormat(viewer,PETSC_VIEWER_NATIVE); the matrices will 156890ace30eSBarry Smith be stored in a way natural for the matrix, for example dense matrices 156990ace30eSBarry Smith would be stored as dense. Matrices stored this way may only be 157090ace30eSBarry Smith read into matrices of the same time. 157190ace30eSBarry Smith */ 157290ace30eSBarry Smith #define MATRIX_BINARY_FORMAT_DENSE -1 157390ace30eSBarry Smith 1574be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetHashTableFactor(Mat,PetscReal); 1575be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatISGetLocalMat(Mat,Mat*); 15761f4e1ec7SBarry Smith 1577d9274352SBarry Smith /*S 1578d9274352SBarry Smith MatNullSpace - Object that removes a null space from a vector, i.e. 1579d9274352SBarry Smith orthogonalizes the vector to a subsapce 1580d9274352SBarry Smith 1581f7a9e4ceSBarry Smith Level: advanced 1582d9274352SBarry Smith 1583d9274352SBarry Smith Concepts: matrix; linear operator, null space 1584d9274352SBarry Smith 15856e1639daSBarry Smith Users manual sections: 15866e1639daSBarry Smith . sec_singular 15876e1639daSBarry Smith 1588d9274352SBarry Smith .seealso: MatNullSpaceCreate() 1589d9274352SBarry Smith S*/ 159074637425SBarry Smith typedef struct _p_MatNullSpace* MatNullSpace; 1591d9274352SBarry Smith 1592be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceCreate(MPI_Comm,PetscTruth,PetscInt,const Vec[],MatNullSpace*); 1593b22b330cSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(MatNullSpace,Vec,void*),void*); 1594be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceDestroy(MatNullSpace); 1595be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceRemove(MatNullSpace,Vec,Vec*); 1596be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceAttach(Mat,MatNullSpace); 159795902228SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceTest(MatNullSpace,Mat,PetscTruth *); 159874637425SBarry Smith 1599be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatReorderingSeqSBAIJ(Mat,IS); 1600be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPISBAIJSetHashTableFactor(Mat,PetscReal); 1601be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetColumnIndices(Mat,PetscInt *); 160246129b97SKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJInvertBlockDiagonal(Mat); 16033f1d51d7SBarry Smith 1604be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMAIJ(Mat,PetscInt,Mat*); 1605be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMAIJRedimension(Mat,PetscInt,Mat*); 1606be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMAIJGetAIJ(Mat,Mat*); 1607c4f061fbSSatish Balay 1608be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatComputeExplicitOperator(Mat,Mat*); 1609b0a32e0cSBarry Smith 1610be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalScaleLocal(Mat,Vec); 161104f1ad80SBarry Smith 1612fef1beadSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*); 16133ec795f1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetBase(Mat,Vec,Vec); 16143ec795f1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*); 1615e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*)); 1616e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec)); 1617e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDAddNullSpace(Mat,MatNullSpace); 1618e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt); 1619e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDResetHHistory(Mat); 1620e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctionError(Mat,PetscReal); 1621e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetPeriod(Mat,PetscInt); 1622e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDGetH(Mat,PetscScalar *); 16236aa9148fSLisandro Dalcin EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetOptionsPrefix(Mat,const char[]); 1624e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFromOptions(Mat); 1625e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*); 1626e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*); 1627e884886fSBarry Smith 16286370053bSBarry Smith /*S 16296370053bSBarry Smith MatMFFD - A data structured used to manage the computation of the h differencing parameter for matrix-free 16306370053bSBarry Smith Jacobian vector products 1631e884886fSBarry Smith 16326370053bSBarry Smith Notes: MATMFFD is a specific MatType which uses the MatMFFD data structure 16336370053bSBarry Smith 16346370053bSBarry Smith MatMFFD*() methods actually take the Mat as their first argument. Not a MatMFFD data structure 16356370053bSBarry Smith 16366370053bSBarry Smith Level: developer 16376370053bSBarry Smith 16386370053bSBarry Smith .seealso: MATMFFD, MatCreateMFFD(), MatMFFDSetFuction(), MatMFFDSetType(), MatMFFDRegister() 16396370053bSBarry Smith S*/ 1640e884886fSBarry Smith typedef struct _p_MatMFFD* MatMFFD; 1641e884886fSBarry Smith 1642e884886fSBarry Smith /*E 1643e884886fSBarry Smith MatMFFDType - algorithm used to compute the h used in computing matrix-vector products via differencing of the function 1644e884886fSBarry Smith 1645e884886fSBarry Smith Level: beginner 1646e884886fSBarry Smith 1647e884886fSBarry Smith .seealso: MatMFFDSetType(), MatMFFDRegister() 1648e884886fSBarry Smith E*/ 1649a313700dSBarry Smith #define MatMFFDType char* 1650e884886fSBarry Smith #define MATMFFD_DS "ds" 1651e884886fSBarry Smith #define MATMFFD_WP "wp" 1652e884886fSBarry Smith 1653a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetType(Mat,const MatMFFDType); 1654e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatMFFD)); 1655e884886fSBarry Smith 1656e884886fSBarry Smith /*MC 1657e884886fSBarry Smith MatMFFDRegisterDynamic - Adds a method to the MatMFFD registry. 1658e884886fSBarry Smith 1659e884886fSBarry Smith Synopsis: 16601890ba74SBarry Smith PetscErrorCode MatMFFDRegisterDynamic(const char *name_solver,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatMFFD)) 1661e884886fSBarry Smith 1662e884886fSBarry Smith Not Collective 1663e884886fSBarry Smith 1664e884886fSBarry Smith Input Parameters: 1665e884886fSBarry Smith + name_solver - name of a new user-defined compute-h module 1666e884886fSBarry Smith . path - path (either absolute or relative) the library containing this solver 1667e884886fSBarry Smith . name_create - name of routine to create method context 1668e884886fSBarry Smith - routine_create - routine to create method context 1669e884886fSBarry Smith 1670e884886fSBarry Smith Level: developer 1671e884886fSBarry Smith 1672e884886fSBarry Smith Notes: 1673e884886fSBarry Smith MatMFFDRegisterDynamic() may be called multiple times to add several user-defined solvers. 1674e884886fSBarry Smith 1675e884886fSBarry Smith If dynamic libraries are used, then the fourth input argument (routine_create) 1676e884886fSBarry Smith is ignored. 1677e884886fSBarry Smith 1678e884886fSBarry Smith Sample usage: 1679e884886fSBarry Smith .vb 1680e884886fSBarry Smith MatMFFDRegisterDynamic("my_h",/home/username/my_lib/lib/libO/solaris/mylib.a, 1681e884886fSBarry Smith "MyHCreate",MyHCreate); 1682e884886fSBarry Smith .ve 1683e884886fSBarry Smith 1684e884886fSBarry Smith Then, your solver can be chosen with the procedural interface via 1685e884886fSBarry Smith $ MatMFFDSetType(mfctx,"my_h") 1686e884886fSBarry Smith or at runtime via the option 1687e884886fSBarry Smith $ -snes_mf_type my_h 1688e884886fSBarry Smith 1689e884886fSBarry Smith .keywords: MatMFFD, register 1690e884886fSBarry Smith 1691e884886fSBarry Smith .seealso: MatMFFDRegisterAll(), MatMFFDRegisterDestroy() 1692e884886fSBarry Smith M*/ 1693e884886fSBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 1694e884886fSBarry Smith #define MatMFFDRegisterDynamic(a,b,c,d) MatMFFDRegister(a,b,c,0) 1695e884886fSBarry Smith #else 1696e884886fSBarry Smith #define MatMFFDRegisterDynamic(a,b,c,d) MatMFFDRegister(a,b,c,d) 1697e884886fSBarry Smith #endif 1698e884886fSBarry Smith 1699e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDRegisterAll(const char[]); 1700e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDRegisterDestroy(void); 17011d0fab5eSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDDSSetUmin(Mat,PetscReal); 1702e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDWPSetComputeNormU(Mat,PetscTruth); 1703e884886fSBarry Smith 1704e884886fSBarry Smith 1705be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *); 1706be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *); 17077dbadf16SMatthew Knepley 170897969023SHong Zhang /* 170997969023SHong Zhang PETSc interface to MUMPS 171097969023SHong Zhang */ 171197969023SHong Zhang #ifdef PETSC_HAVE_MUMPS 171297969023SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMumpsSetIcntl(Mat,PetscInt,PetscInt); 171397969023SHong Zhang #endif 171423a5497aSJed Brown 171523a5497aSJed Brown PETSC_EXTERN_CXX_END 171623a5497aSJed Brown #endif 1717