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" 315a11e1b2SBarry Smith #define MATMAIJ "maij" 32273d9f13SBarry Smith #define MATSEQMAIJ "seqmaij" 33273d9f13SBarry Smith #define MATMPIMAIJ "mpimaij" 34273d9f13SBarry Smith #define MATIS "is" 355a11e1b2SBarry Smith #define MATAIJ "aij" 36273d9f13SBarry Smith #define MATSEQAIJ "seqaij" 37273d9f13SBarry Smith #define MATMPIAIJ "mpiaij" 385a11e1b2SBarry Smith #define MATAIJCRL "aijcrl" 395a11e1b2SBarry Smith #define MATSEQAIJCRL "seqaijcrl" 405a11e1b2SBarry Smith #define MATMPIAIJCRL "mpiaijcrl" 415a11e1b2SBarry Smith #define MATAIJCUDA "aijcuda" 425a11e1b2SBarry Smith #define MATSEQAIJCUDA "seqaijcuda" 435a11e1b2SBarry Smith #define MATMPIAIJCUDA "mpiaijcuda" 445a11e1b2SBarry Smith #define MATAIJPERM "aijperm" 455a11e1b2SBarry Smith #define MATSEQAIJPERM "seqaijperm" 465a11e1b2SBarry Smith #define MATMPIAIJPERM "mpiaijperm" 47273d9f13SBarry Smith #define MATSHELL "shell" 485a11e1b2SBarry Smith #define MATDENSE "dense" 49209238afSKris Buschelman #define MATSEQDENSE "seqdense" 50273d9f13SBarry Smith #define MATMPIDENSE "mpidense" 515a11e1b2SBarry Smith #define MATBAIJ "baij" 52273d9f13SBarry Smith #define MATSEQBAIJ "seqbaij" 53273d9f13SBarry Smith #define MATMPIBAIJ "mpibaij" 54273d9f13SBarry Smith #define MATMPIADJ "mpiadj" 555a11e1b2SBarry Smith #define MATSBAIJ "sbaij" 56273d9f13SBarry Smith #define MATSEQSBAIJ "seqsbaij" 57273d9f13SBarry Smith #define MATMPISBAIJ "mpisbaij" 58cebc7f6cSBarry Smith #define MATDAAD "daad" 59cebc7f6cSBarry Smith #define MATMFFD "mffd" 60c8a8475eSBarry Smith #define MATNORMAL "normal" 61ab92ecdeSBarry Smith #define MATLRC "lrc" 622a6744ebSBarry Smith #define MATSCATTER "scatter" 63421e10b8SBarry Smith #define MATBLOCKMAT "blockmat" 64793850ffSBarry Smith #define MATCOMPOSITE "composite" 655a7f1df3SHong Zhang #define MATSEQFFTW "seqfftw" 66e133240eSMatthew G Knepley #define MATSEQCUFFT "seqcufft" 67557cca28SSatish Balay #define MATTRANSPOSEMAT "transpose" 6872ca8751SBarry Smith #define MATSCHURCOMPLEMENT "schurcomplement" 691d6018f0SLisandro Dalcin #define MATPYTHON "python" 70f91d8e95SBarry Smith #define MATHYPRESTRUCT "hyprestruct" 71a9e6138eSGlenn Hammond #define MATHYPRESSTRUCT "hypresstruct" 72ee1cef2cSJed Brown #define MATSUBMATRIX "submatrix" 73dbc6227fSDmitry Karpeev #define MATDD "matdd" 74ba2f8784SDmitry Karpeev #define MATIM "matim" 75ba2f8784SDmitry Karpeev 76773941d6SBarry Smith 779989ab13SBarry Smith /*E 78c7393fdbSBarry Smith MatSolverPackage - String with the name of a PETSc matrix solver type. 799989ab13SBarry Smith 809989ab13SBarry Smith For example: "petsc" indicates what PETSc provides, "superlu" indicates either 819989ab13SBarry Smith SuperLU or SuperLU_Dist etc. 829989ab13SBarry Smith 839989ab13SBarry Smith 849989ab13SBarry Smith Level: beginner 859989ab13SBarry Smith 865c9eb25fSBarry Smith .seealso: MatGetFactor(), Mat, MatSetType(), MatType 879989ab13SBarry Smith E*/ 88c7393fdbSBarry Smith #define MatSolverPackage char* 892692d6eeSBarry Smith #define MATSOLVERSPOOLES "spooles" 902692d6eeSBarry Smith #define MATSOLVERSUPERLU "superlu" 912692d6eeSBarry Smith #define MATSOLVERSUPERLU_DIST "superlu_dist" 922692d6eeSBarry Smith #define MATSOLVERUMFPACK "umfpack" 9320db9a53SJed Brown #define MATSOLVERCHOLMOD "cholmod" 942692d6eeSBarry Smith #define MATSOLVERESSL "essl" 952692d6eeSBarry Smith #define MATSOLVERLUSOL "lusol" 962692d6eeSBarry Smith #define MATSOLVERMUMPS "mumps" 972692d6eeSBarry Smith #define MATSOLVERPASTIX "pastix" 982692d6eeSBarry Smith #define MATSOLVERDSCPACK "dscpack" 992692d6eeSBarry Smith #define MATSOLVERMATLAB "matlab" 1002692d6eeSBarry Smith #define MATSOLVERPETSC "petsc" 1012692d6eeSBarry Smith #define MATSOLVERPLAPACK "plapack" 1022692d6eeSBarry Smith #define MATSOLVERBAS "bas" 103773941d6SBarry Smith 104b24902e0SBarry Smith /*E 105b24902e0SBarry Smith MatFactorType - indicates what type of factorization is requested 106b24902e0SBarry Smith 107b24902e0SBarry Smith Level: beginner 108b24902e0SBarry Smith 109b24902e0SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 110b24902e0SBarry Smith 111c7393fdbSBarry Smith .seealso: MatSolverPackage, MatGetFactor() 112b24902e0SBarry Smith E*/ 113599ef60dSHong Zhang typedef enum {MAT_FACTOR_NONE, MAT_FACTOR_LU, MAT_FACTOR_CHOLESKY, MAT_FACTOR_ILU, MAT_FACTOR_ICC,MAT_FACTOR_ILUDT} MatFactorType; 11434918c53SJed Brown extern const char *const MatFactorTypes[]; 115e92e720dSBarry Smith 116c7393fdbSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetFactor(Mat,const MatSolverPackage,MatFactorType,Mat*); 117ace3abfcSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetFactorAvailable(Mat,const MatSolverPackage,MatFactorType,PetscBool *); 11835bd34faSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFactorGetSolverPackage(Mat,const MatSolverPackage*); 1196335c336SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetFactorType(Mat,MatFactorType*); 1209989ab13SBarry Smith 121c06d978dSMatthew Knepley /* Logging support */ 1220700a824SBarry Smith #define MAT_FILE_CLASSID 1211216 /* used to indicate matrices in binary files */ 1230700a824SBarry Smith extern PetscClassId PETSCMAT_DLLEXPORT MAT_CLASSID; 1240700a824SBarry Smith extern PetscClassId PETSCMAT_DLLEXPORT MAT_FDCOLORING_CLASSID; 1250700a824SBarry Smith extern PetscClassId PETSCMAT_DLLEXPORT MAT_PARTITIONING_CLASSID; 1260700a824SBarry Smith extern PetscClassId PETSCMAT_DLLEXPORT MAT_NULLSPACE_CLASSID; 1270700a824SBarry Smith extern PetscClassId PETSCMAT_DLLEXPORT MATMFFD_CLASSID; 128c06d978dSMatthew Knepley 129ceb03754SKris Buschelman /*E 130ceb03754SKris Buschelman MatReuse - Indicates if matrices obtained from a previous call to MatGetSubMatrices() 131ceb03754SKris Buschelman or MatGetSubMatrix() are to be reused to store the new matrix values. 132ceb03754SKris Buschelman 133ceb03754SKris Buschelman Level: beginner 134ceb03754SKris Buschelman 135ceb03754SKris Buschelman Any additions/changes here MUST also be made in include/finclude/petscmat.h 136ceb03754SKris Buschelman 1370c451bc4SBarry Smith .seealso: MatGetSubMatrices(), MatGetSubMatrix(), MatDestroyMatrices(), MatConvert() 138ceb03754SKris Buschelman E*/ 139dfe085dbSJed Brown typedef enum {MAT_INITIAL_MATRIX,MAT_REUSE_MATRIX,MAT_IGNORE_MATRIX} MatReuse; 140ceb03754SKris Buschelman 1415494a064SHong Zhang /*E 1425494a064SHong Zhang MatGetSubMatrixOption - Indicates if matrices obtained from a call to MatGetSubMatrices() 143829201f2SHong Zhang include the matrix values. Currently it is only used by MatGetSeqNonzerostructure(). 1445494a064SHong Zhang 1455494a064SHong Zhang Level: beginner 1465494a064SHong Zhang 147829201f2SHong Zhang .seealso: MatGetSeqNonzerostructure() 1485494a064SHong Zhang E*/ 1495494a064SHong Zhang typedef enum {MAT_DO_NOT_GET_VALUES,MAT_GET_VALUES} MatGetSubMatrixOption; 1505494a064SHong Zhang 151e5bd5246SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInitializePackage(const char[]); 152c06d978dSMatthew Knepley 153f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreate(MPI_Comm,Mat*); 154a69a960dSMatthew Knepley PetscPolymorphicFunction(MatCreate,(MPI_Comm comm),(comm,&A),Mat,A) 155a69a960dSMatthew Knepley PetscPolymorphicFunction(MatCreate,(),(PETSC_COMM_WORLD,&A),Mat,A) 156f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetSizes(Mat,PetscInt,PetscInt,PetscInt,PetscInt); 157a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetType(Mat,const MatType); 158be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetFromOptions(Mat); 159be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUpPreallocation(Mat); 160be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRegisterAll(const char[]); 161be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat)); 16230de9b25SBarry Smith 163f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetOptionsPrefix(Mat,const char[]); 164f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAppendOptionsPrefix(Mat,const char[]); 165f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOptionsPrefix(Mat,const char*[]); 166f69a0ea3SMatthew Knepley 16730de9b25SBarry Smith /*MC 16830de9b25SBarry Smith MatRegisterDynamic - Adds a new matrix type 16930de9b25SBarry Smith 17030de9b25SBarry Smith Synopsis: 1711890ba74SBarry Smith PetscErrorCode MatRegisterDynamic(const char *name,const char *path,const char *name_create,PetscErrorCode (*routine_create)(Mat)) 17230de9b25SBarry Smith 17330de9b25SBarry Smith Not Collective 17430de9b25SBarry Smith 17530de9b25SBarry Smith Input Parameters: 17630de9b25SBarry Smith + name - name of a new user-defined matrix type 17730de9b25SBarry Smith . path - path (either absolute or relative) the library containing this solver 17830de9b25SBarry Smith . name_create - name of routine to create method context 17930de9b25SBarry Smith - routine_create - routine to create method context 18030de9b25SBarry Smith 18130de9b25SBarry Smith Notes: 18230de9b25SBarry Smith MatRegisterDynamic() may be called multiple times to add several user-defined solvers. 18330de9b25SBarry Smith 18430de9b25SBarry Smith If dynamic libraries are used, then the fourth input argument (routine_create) 18530de9b25SBarry Smith is ignored. 18630de9b25SBarry Smith 18730de9b25SBarry Smith Sample usage: 18830de9b25SBarry Smith .vb 18930de9b25SBarry Smith MatRegisterDynamic("my_mat",/home/username/my_lib/lib/libO/solaris/mylib.a, 19030de9b25SBarry Smith "MyMatCreate",MyMatCreate); 19130de9b25SBarry Smith .ve 19230de9b25SBarry Smith 19330de9b25SBarry Smith Then, your solver can be chosen with the procedural interface via 19430de9b25SBarry Smith $ MatSetType(Mat,"my_mat") 19530de9b25SBarry Smith or at runtime via the option 19630de9b25SBarry Smith $ -mat_type my_mat 19730de9b25SBarry Smith 19830de9b25SBarry Smith Level: advanced 19930de9b25SBarry Smith 200ab901514SBarry Smith Notes: ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values. 20130de9b25SBarry Smith If your function is not being put into a shared library then use VecRegister() instead 20230de9b25SBarry Smith 20330de9b25SBarry Smith .keywords: Mat, register 20430de9b25SBarry Smith 20530de9b25SBarry Smith .seealso: MatRegisterAll(), MatRegisterDestroy() 20630de9b25SBarry Smith 20730de9b25SBarry Smith M*/ 208273d9f13SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 209273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,0) 210273d9f13SBarry Smith #else 211273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,d) 21230de9b25SBarry Smith #endif 21330de9b25SBarry Smith 214ace3abfcSBarry Smith extern PetscBool MatRegisterAllCalled; 215b0a32e0cSBarry Smith extern PetscFList MatList; 216b022a5c1SBarry Smith extern PetscFList MatColoringList; 217b022a5c1SBarry Smith extern PetscFList MatPartitioningList; 21828988994SBarry Smith 2193b224e63SBarry Smith /*E 2203b224e63SBarry Smith MatStructure - Indicates if the matrix has the same nonzero structure 2213b224e63SBarry Smith 2223b224e63SBarry Smith Level: beginner 2233b224e63SBarry Smith 2243b224e63SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 2253b224e63SBarry Smith 2263b224e63SBarry Smith .seealso: MatCopy(), KSPSetOperators(), PCSetOperators() 2273b224e63SBarry Smith E*/ 2283b224e63SBarry Smith typedef enum {SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER,SUBSET_NONZERO_PATTERN} MatStructure; 2293b224e63SBarry Smith 230be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqDense(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*); 231be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIDense(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*); 232be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 233ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,nz,nnz,&A),Mat,A) 234ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,&A),Mat,A) 235ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,0,nnz,&A),Mat,A) 236ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,&A),Mat,A) 237ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,A)) 238ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,m,n,0,nnz,A)) 239ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,A)) 240be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 241ba966639SSatish 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) 242ba966639SSatish 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) 243ba966639SSatish 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) 244ba966639SSatish 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) 245ba966639SSatish 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)) 246ba966639SSatish 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)) 247ba966639SSatish 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)) 248ba966639SSatish 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) 249ba966639SSatish 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) 250ba966639SSatish 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) 251ba966639SSatish 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) 252ba966639SSatish 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)) 253ba966639SSatish 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)) 254ba966639SSatish 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)) 25563ab93bfSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *); 2568d7a6e47SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJWithSplitArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],PetscInt[],PetscInt[],PetscScalar[],Mat*); 2578d7a6e47SBarry Smith 258be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 259ba966639SSatish 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) 260ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A) 261ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A) 262ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A) 263ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A)) 264ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A)) 265ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A)) 266be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 267ba966639SSatish 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) 268ba966639SSatish 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) 269ba966639SSatish 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) 270ba966639SSatish 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) 271ba966639SSatish 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)) 272ba966639SSatish 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)) 273ba966639SSatish 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)) 274ba966639SSatish 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) 275ba966639SSatish 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) 276ba966639SSatish 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) 277ba966639SSatish 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) 278ba966639SSatish 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)) 279ba966639SSatish 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)) 280ba966639SSatish 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)) 281d21a29f3SJed Brown EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat*); 282d21a29f3SJed Brown 283be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAdj(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscInt[],Mat*); 284be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 285ba966639SSatish 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) 286ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A) 287ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A) 288ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A) 289ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A)) 290ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A)) 291ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A)) 292d21a29f3SJed Brown 293be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPISBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 294ba966639SSatish 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) 295ba966639SSatish 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) 296ba966639SSatish 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) 297ba966639SSatish 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) 298ba966639SSatish 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)) 299ba966639SSatish 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)) 300ba966639SSatish 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)) 301ba966639SSatish 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) 302ba966639SSatish 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) 303ba966639SSatish 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) 304ba966639SSatish 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) 305ba966639SSatish 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)) 306ba966639SSatish 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)) 307ba966639SSatish 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)) 308dfb205c3SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPISBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *); 309dfb205c3SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPISBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 310dfb205c3SBarry Smith 311be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateShell(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,void *,Mat*); 312ba966639SSatish 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) 313ba966639SSatish Balay PetscPolymorphicFunction(MatCreateShell,(PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx),(PETSC_COMM_WORLD,m,n,M,N,ctx,&A),Mat,A) 314be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateAdic(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,void (*)(void),Mat*); 315be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateNormal(Mat,Mat*); 316ba966639SSatish Balay PetscPolymorphicFunction(MatCreateNormal,(Mat mat),(mat,&A),Mat,A) 317ab92ecdeSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateLRC(Mat,Mat,Mat,Mat*); 318284134d9SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateIS(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,ISLocalToGlobalMapping,Mat*); 3195a11e1b2SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 3205a11e1b2SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 3212a6744ebSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateScatter(MPI_Comm,VecScatter,Mat*); 3222a6744ebSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScatterSetVecScatter(Mat,VecScatter); 3232a6744ebSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScatterGetVecScatter(Mat,VecScatter*); 3248cdf2d9bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateBlockMat(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,Mat*); 325793850ffSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompositeAddMat(Mat,Mat); 326b52f573bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompositeMerge(Mat); 327793850ffSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateComposite(MPI_Comm,PetscInt,const Mat*,Mat*); 3286d7c1e57SBarry Smith typedef enum {MAT_COMPOSITE_ADDITIVE,MAT_COMPOSITE_MULTIPLICATIVE} MatCompositeType; 3296d7c1e57SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompositeSetType(Mat,MatCompositeType); 3306d7c1e57SBarry Smith 33145b63f25SDmitry Karpeev #if defined PETSC_HAVE_MATDD 332d214e108SDmitry Karpeev EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDDCreate(Mat A); 333d214e108SDmitry Karpeev typedef enum {MATDD_DOMAINS_COLUMN, MATDD_DOMAINS_ROW} MatDDDomainType; 334d214e108SDmitry Karpeev EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDDSetDomainsLocal(Mat A, MatDDDomainType type, PetscInt domain_count, const PetscInt *supported_domains, const PetscInt *domain_limits, PetscBool covering); 335d214e108SDmitry Karpeev EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDDSetDomainsLocalIS(Mat A, MatDDDomainType type, IS supported_domains, IS domain_limits, PetscBool covering); 336d214e108SDmitry Karpeev EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDDSetScatter(Mat A, Mat S); 337d214e108SDmitry Karpeev EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDDSetGather(Mat A, Mat G); 338d214e108SDmitry Karpeev /**/ 339dbc6227fSDmitry Karpeev typedef enum {MATDD_BLOCK_COMM_DEFAULT = 0, MATDD_BLOCK_COMM_SELF = -1, MATDD_BLOCK_COMM_DETERMINE = -2} MatDDBlockCommType; 340dbc6227fSDmitry Karpeev EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDDGetDefaltBlockType(Mat A, const MatType *type); 341dbc6227fSDmitry Karpeev EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDDSetDefaltBlockType(Mat A, const MatType type); 342d214e108SDmitry Karpeev EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDDAddBlockLocal(Mat A, PetscInt rowblock, PetscInt colblock, const MatType blockmattype, MatDDBlockCommType blockcommtype, Mat *block); 343d214e108SDmitry Karpeev EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDDSetBlockLocal(Mat A, PetscInt rowblock, PetscInt colblock, Mat block); 344d214e108SDmitry Karpeev EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDDGetBlockLocal(Mat A, PetscInt rowblock, PetscInt colblock, Mat *block); 345d214e108SDmitry Karpeev /**/ 346d214e108SDmitry Karpeev EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDDAIJSetPreallocation(Mat A,PetscInt nz,PetscInt *nnz); 34745b63f25SDmitry Karpeev #endif 34862ee7b53SDmitry Karpeev 349ba2f8784SDmitry Karpeev #if defined PETSC_HAVE_MATIM 350ba2f8784SDmitry Karpeev EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIMSetIS(Mat A, IS in, IS out); 351ba2f8784SDmitry Karpeev EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIMGetIS(Mat A, IS *in, IS *out); 352ba2f8784SDmitry Karpeev #endif 353ba2f8784SDmitry Karpeev 354ba2f8784SDmitry Karpeev 355a6a5cd3fSDmitry Karpeev 3565a7f1df3SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqFFTW(MPI_Comm,PetscInt,const PetscInt[],Mat*); 357e133240eSMatthew G Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqCUFFT(MPI_Comm,PetscInt,const PetscInt[],Mat*); 358f20085c4SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateTranspose(Mat,Mat*); 359ee1cef2cSJed Brown EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSubMatrix(Mat,IS,IS,Mat*); 360ee1cef2cSJed Brown EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSubMatrixUpdate(Mat,Mat,IS,IS); 3611472f72bSBarry Smith 3621d6018f0SLisandro Dalcin EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreatePython(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const char[],Mat*); 3631d6018f0SLisandro Dalcin EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPythonSetType(Mat,const char[]); 3641d6018f0SLisandro Dalcin 3651d6018f0SLisandro Dalcin 366f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUp(Mat); 367be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroy(Mat); 36821c89e3eSBarry Smith 3690c451bc4SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConjugate(Mat); 37099cafbc1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRealPart(Mat); 37199cafbc1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatImaginaryPart(Mat); 372ace3abfcSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetDiagonalBlock(Mat,PetscBool *,MatReuse,Mat*); 373fe97e370SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetTrace(Mat,PetscScalar*); 37499cafbc1SBarry Smith 3758ed539a5SBarry Smith /* ------------------------------------------------------------*/ 376be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 377be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 37887d4246cSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesRow(Mat,PetscInt,const PetscScalar[]); 379f259bd47SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesRowLocal(Mat,PetscInt,const PetscScalar[]); 38084cb2905SBarry Smith 3812ef4de8bSBarry Smith /*S 3822ef4de8bSBarry Smith MatStencil - Data structure (C struct) for storing information about a single row or 3832ef4de8bSBarry Smith column of a matrix as index on an associated grid. 3842ef4de8bSBarry Smith 3852ef4de8bSBarry Smith Level: beginner 3862ef4de8bSBarry Smith 3872ef4de8bSBarry Smith Concepts: matrix; linear operator 3882ef4de8bSBarry Smith 389d7bd93baSBarry Smith .seealso: MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockStencil() 3902ef4de8bSBarry Smith S*/ 391435da068SBarry Smith typedef struct { 392c1ac3661SBarry Smith PetscInt k,j,i,c; 393435da068SBarry Smith } MatStencil; 3942ef4de8bSBarry Smith 395be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode); 396be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode); 397be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt); 398435da068SBarry Smith 399be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetColoring(Mat,ISColoring); 400be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesAdic(Mat,void*); 401be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesAdifor(Mat,PetscInt,void*); 4023a7fca6bSBarry Smith 403d91e6319SBarry Smith /*E 404d91e6319SBarry Smith MatAssemblyType - Indicates if the matrix is now to be used, or if you plan 405d91e6319SBarry Smith to continue to add values to it 406d91e6319SBarry Smith 407d91e6319SBarry Smith Level: beginner 408d91e6319SBarry Smith 409d91e6319SBarry Smith .seealso: MatAssemblyBegin(), MatAssemblyEnd() 410d91e6319SBarry Smith E*/ 4116d4a8577SBarry Smith typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType; 412be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssemblyBegin(Mat,MatAssemblyType); 413be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssemblyEnd(Mat,MatAssemblyType); 414ace3abfcSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssembled(Mat,PetscBool *); 4154f9c727eSBarry Smith 4161d73ed98SBarry Smith 41730de9b25SBarry Smith 418d91e6319SBarry Smith /*E 419d91e6319SBarry Smith MatOption - Options that may be set for a matrix and its behavior or storage 420d91e6319SBarry Smith 421d91e6319SBarry Smith Level: beginner 422d91e6319SBarry Smith 4230a835dfdSSatish Balay Any additions/changes here MUST also be made in include/finclude/petscmat.h 424d91e6319SBarry Smith 425d91e6319SBarry Smith .seealso: MatSetOption() 426d91e6319SBarry Smith E*/ 427512a5fc5SBarry Smith typedef enum {MAT_ROW_ORIENTED,MAT_NEW_NONZERO_LOCATIONS, 4284e0d8c25SBarry Smith MAT_SYMMETRIC, 4294e0d8c25SBarry Smith MAT_STRUCTURALLY_SYMMETRIC, 4308047e77bSBarry Smith MAT_NEW_DIAGONALS,MAT_IGNORE_OFF_PROC_ENTRIES, 4314e0d8c25SBarry Smith MAT_NEW_NONZERO_LOCATION_ERR, 4324e0d8c25SBarry Smith MAT_NEW_NONZERO_ALLOCATION_ERR,MAT_USE_HASH_TABLE, 433a9817697SBarry Smith MAT_KEEP_NONZERO_PATTERN,MAT_IGNORE_ZERO_ENTRIES, 4344e0d8c25SBarry Smith MAT_USE_INODES, 4354e0d8c25SBarry Smith MAT_HERMITIAN, 4364e0d8c25SBarry Smith MAT_SYMMETRY_ETERNAL, 4374e0d8c25SBarry Smith MAT_USE_COMPRESSEDROW, 4384e0d8c25SBarry Smith MAT_IGNORE_LOWER_TRIANGULAR,MAT_ERROR_LOWER_TRIANGULAR, 43928b2fa4aSMatthew Knepley MAT_GETROW_UPPERTRIANGULAR,MAT_UNUSED_NONZERO_LOCATION_ERR, 4404cb17eb5SBarry Smith MAT_SPD,MAT_NO_OFF_PROC_ENTRIES,MAT_NO_OFF_PROC_ZERO_ROWS, 44128b2fa4aSMatthew Knepley NUM_MAT_OPTIONS} MatOption; 442290bbb0aSBarry Smith extern const char *MatOptions[]; 443ace3abfcSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetOption(Mat,MatOption,PetscBool ); 444a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetType(Mat,const MatType*); 445a313700dSBarry Smith PetscPolymorphicFunction(MatGetType,(Mat mat),(mat,&t),const MatType,t) 44684cb2905SBarry Smith 447be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]); 448be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 449be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 450f5edf698SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowUpperTriangular(Mat); 451f5edf698SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRowUpperTriangular(Mat); 452be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 453be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 454be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumnVector(Mat,Vec,PetscInt); 455be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetArray(Mat,PetscScalar *[]); 456ba966639SSatish Balay PetscPolymorphicFunction(MatGetArray,(Mat mat),(mat,&a),PetscScalar*,a) 457be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreArray(Mat,PetscScalar *[]); 458be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBlockSize(Mat,PetscInt *); 459ba966639SSatish Balay PetscPolymorphicFunction(MatGetBlockSize,(Mat mat),(mat,&a),PetscInt,a) 460be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetBlockSize(Mat,PetscInt); 4617b80b807SBarry Smith 4621620fd73SBarry Smith 463be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMult(Mat,Vec,Vec); 46489c6957cSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultDiagonalBlock(Mat,Vec,Vec); 465be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultAdd(Mat,Vec,Vec,Vec); 466ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultAdd,(Mat A,Vec x,Vec y),(A,x,y,y)) 467be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTranspose(Mat,Vec,Vec); 468547795f9SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultHermitianTranspose(Mat,Vec,Vec); 469ace3abfcSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsTranspose(Mat,Mat,PetscReal,PetscBool *); 470ace3abfcSBarry Smith PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B,PetscReal tol),(A,B,tol,&t),PetscBool ,t) 471ace3abfcSBarry Smith PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B),(A,B,0,&t),PetscBool ,t) 472ace3abfcSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitianTranspose(Mat,Mat,PetscReal,PetscBool *); 473be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeAdd(Mat,Vec,Vec,Vec); 4745cb05578SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultHermitianTransposeAdd(Mat,Vec,Vec,Vec); 475ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultTransposeAdd,(Mat A,Vec x,Vec y),(A,x,y,y)) 476be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultConstrained(Mat,Vec,Vec); 477be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeConstrained(Mat,Vec,Vec); 4782bebee5dSHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatSolve(Mat,Mat,Mat); 479f5edf698SHong Zhang 480d91e6319SBarry Smith /*E 481d91e6319SBarry Smith MatDuplicateOption - Indicates if a duplicated sparse matrix should have 482d91e6319SBarry Smith its numerical values copied over or just its nonzero structure. 483d91e6319SBarry Smith 484d91e6319SBarry Smith Level: beginner 485d91e6319SBarry Smith 486d91e6319SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 487d91e6319SBarry Smith 48870dcbbb9SBarry Smith $ MAT_SHARE_NONZERO_PATTERN - the i and j arrays in the new matrix will be shared with the original matrix 48970dcbbb9SBarry Smith $ this also triggers the MAT_DO_NOT_COPY_VALUES option. This is used when you 49070dcbbb9SBarry Smith $ have several matrices with the same nonzero pattern. 49170dcbbb9SBarry Smith 492d91e6319SBarry Smith .seealso: MatDuplicate() 493d91e6319SBarry Smith E*/ 49470dcbbb9SBarry Smith typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES,MAT_SHARE_NONZERO_PATTERN} MatDuplicateOption; 4952e8a6d31SBarry Smith 496a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert(Mat,const MatType,MatReuse,Mat*); 497a313700dSBarry Smith PetscPolymorphicFunction(MatConvert,(Mat A,const MatType t),(A,t,MAT_INITIAL_MATRIX,&a),Mat,a) 498be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDuplicate(Mat,MatDuplicateOption,Mat*); 499ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A,MatDuplicateOption o),(A,o,&a),Mat,a) 500ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A),(A,MAT_COPY_VALUES,&a),Mat,a) 50194a9d846SBarry Smith 502cb5b572fSBarry Smith 503be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCopy(Mat,Mat,MatStructure); 504be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatView(Mat,PetscViewer); 505ace3abfcSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsSymmetric(Mat,PetscReal,PetscBool *); 506ace3abfcSBarry Smith PetscPolymorphicFunction(MatIsSymmetric,(Mat A,PetscReal tol),(A,tol,&t),PetscBool ,t) 507ace3abfcSBarry Smith PetscPolymorphicFunction(MatIsSymmetric,(Mat A),(A,0,&t),PetscBool ,t) 508ace3abfcSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsStructurallySymmetric(Mat,PetscBool *); 509ace3abfcSBarry Smith PetscPolymorphicFunction(MatIsStructurallySymmetric,(Mat A),(A,&t),PetscBool ,t) 510ace3abfcSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitian(Mat,PetscReal,PetscBool *); 511ace3abfcSBarry Smith PetscPolymorphicFunction(MatIsHermitian,(Mat A),(A,0,&t),PetscBool ,t) 512ace3abfcSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsSymmetricKnown(Mat,PetscBool *,PetscBool *); 513ace3abfcSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitianKnown(Mat,PetscBool *,PetscBool *); 514ace3abfcSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMissingDiagonal(Mat,PetscBool *,PetscInt *); 515112444f4SShri Abhyankar EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLoad(Mat, PetscViewer); 5167b80b807SBarry Smith 517ace3abfcSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,PetscInt *[],PetscInt *[],PetscBool *); 518ace3abfcSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,PetscInt *[],PetscInt *[],PetscBool *); 519ace3abfcSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,PetscInt *[],PetscInt *[],PetscBool *); 520ace3abfcSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,PetscInt *[],PetscInt *[],PetscBool *); 521d4fbbf0eSBarry Smith 522d91e6319SBarry Smith /*S 523d91e6319SBarry Smith MatInfo - Context of matrix information, used with MatGetInfo() 524d91e6319SBarry Smith 525d91e6319SBarry Smith In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE 526d91e6319SBarry Smith 527d91e6319SBarry Smith Level: intermediate 528d91e6319SBarry Smith 529d91e6319SBarry Smith Concepts: matrix^nonzero information 530d91e6319SBarry Smith 531d9274352SBarry Smith .seealso: MatGetInfo(), MatInfoType 532d91e6319SBarry Smith S*/ 5334e220ebcSLois Curfman McInnes typedef struct { 534b0a32e0cSBarry Smith PetscLogDouble block_size; /* block size */ 535b0a32e0cSBarry Smith PetscLogDouble nz_allocated,nz_used,nz_unneeded; /* number of nonzeros */ 536b0a32e0cSBarry Smith PetscLogDouble memory; /* memory allocated */ 537b0a32e0cSBarry Smith PetscLogDouble assemblies; /* number of matrix assemblies called */ 538b0a32e0cSBarry Smith PetscLogDouble mallocs; /* number of mallocs during MatSetValues() */ 539b0a32e0cSBarry Smith PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */ 540b0a32e0cSBarry Smith PetscLogDouble factor_mallocs; /* number of mallocs during factorization */ 5414e220ebcSLois Curfman McInnes } MatInfo; 5424e220ebcSLois Curfman McInnes 543d9274352SBarry Smith /*E 544d9274352SBarry Smith MatInfoType - Indicates if you want information about the local part of the matrix, 545d9274352SBarry Smith the entire parallel matrix or the maximum over all the local parts. 546d9274352SBarry Smith 547d9274352SBarry Smith Level: beginner 548d9274352SBarry Smith 549d9274352SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 550d9274352SBarry Smith 551d9274352SBarry Smith .seealso: MatGetInfo(), MatInfo 552d9274352SBarry Smith E*/ 5537b80b807SBarry Smith typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType; 554be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetInfo(Mat,MatInfoType,MatInfo*); 555be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetDiagonal(Mat,Vec); 556985db425SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMax(Mat,Vec,PetscInt[]); 557985db425SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMin(Mat,Vec,PetscInt[]); 558985db425SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMaxAbs(Mat,Vec,PetscInt[]); 559c87e5d42SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMinAbs(Mat,Vec,PetscInt[]); 56086697f06SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowSum(Mat,Vec); 561fc4dec0aSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatTranspose(Mat,MatReuse,Mat*); 562fc4dec0aSBarry Smith PetscPolymorphicFunction(MatTranspose,(Mat A),(A,MAT_INITIAL_MATRIX,&t),Mat,t) 563eeffb40dSHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatHermitianTranspose(Mat,MatReuse,Mat*); 564be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPermute(Mat,IS,IS,Mat *); 565ba966639SSatish Balay PetscPolymorphicFunction(MatPermute,(Mat A,IS is1,IS is2),(A,is1,is2,&t),Mat,t) 566be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPermuteSparsify(Mat,PetscInt,PetscReal,PetscReal,IS,IS,Mat *); 567be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalScale(Mat,Vec,Vec); 568be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalSet(Mat,Vec,InsertMode); 569ace3abfcSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatEqual(Mat,Mat,PetscBool *); 570ace3abfcSBarry Smith PetscPolymorphicFunction(MatEqual,(Mat A,Mat B),(A,B,&t),PetscBool ,t) 571ace3abfcSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultEqual(Mat,Mat,PetscInt,PetscBool *); 572ace3abfcSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultAddEqual(Mat,Mat,PetscInt,PetscBool *); 573ace3abfcSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeEqual(Mat,Mat,PetscInt,PetscBool *); 574ace3abfcSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscBool *); 5757b80b807SBarry Smith 576be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNorm(Mat,NormType,PetscReal *); 577ba966639SSatish Balay PetscPolymorphicFunction(MatNorm,(Mat A,NormType t),(A,t,&n),PetscReal,n) 578242f1d38SBarry Smith EXTERN PetscErrorCode MatGetColumnNorms(Mat,NormType,PetscReal *); 579be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroEntries(Mat); 5802b40b63fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 5812b40b63fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsIS(Mat,IS,PetscScalar,Vec,Vec); 5822b40b63fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsStencil(Mat,PetscInt,const MatStencil [],PetscScalar,Vec,Vec); 5836e169961SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsColumns(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 5846e169961SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsColumnsIS(Mat,IS,PetscScalar,Vec,Vec); 5857b80b807SBarry Smith 586ace3abfcSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatUseScaledForm(Mat,PetscBool ); 587be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScaleSystem(Mat,Vec,Vec); 588be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatUnScaleSystem(Mat,Vec,Vec); 5895ef9f2a5SBarry Smith 590be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSize(Mat,PetscInt*,PetscInt*); 591be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalSize(Mat,PetscInt*,PetscInt*); 592be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRange(Mat,PetscInt*,PetscInt*); 5938ee2e534SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRanges(Mat,const PetscInt**); 594f4b8409dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRangeColumn(Mat,PetscInt*,PetscInt*); 59572ca8751SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRangesColumn(Mat,const PetscInt**); 5967b80b807SBarry Smith 597be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]); 598be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroyMatrices(PetscInt,Mat *[]); 5994aa3045dSJed Brown EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrix(Mat,IS,IS,MatReuse,Mat *); 600f6d58c54SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSeqNonzeroStructure(Mat,Mat*); 601f6d58c54SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroySeqNonzeroStructure(Mat*); 6025494a064SHong Zhang 603be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge(MPI_Comm,Mat,PetscInt,MatReuse,Mat*); 604be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPI(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*); 605be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPISymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*); 606be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPINumeric(Mat,Mat); 607be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroy_MPIAIJ_SeqsToMPI(Mat); 608be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalMat(Mat,MatReuse,Mat*); 609be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*); 610be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,PetscInt*,Mat*); 6111d79065fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBrowsOfAoCols(Mat,Mat,MatReuse,PetscInt**,PetscInt**,MatScalar**,Mat*); 61243eb5e2fSMatthew Knepley #if defined (PETSC_USE_CTABLE) 613cd116e26SSatish Balay #include "petscctable.h" 61443eb5e2fSMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetCommunicationStructs(Mat, Vec *, PetscTable *, VecScatter *); 61543eb5e2fSMatthew Knepley #else 61643eb5e2fSMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetCommunicationStructs(Mat, Vec *, PetscInt *[], VecScatter *); 61743eb5e2fSMatthew Knepley #endif 6188c7482ecSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetGhosts(Mat, PetscInt *,const PetscInt *[]); 6198efafbd8SBarry Smith 620be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt); 6217b80b807SBarry Smith 622be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*); 623be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultSymbolic(Mat,Mat,PetscReal,Mat*); 624be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultNumeric(Mat,Mat,Mat); 62522440eb1SKris Buschelman 626be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*); 627be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAPSymbolic(Mat,Mat,PetscReal,Mat*); 628be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAPNumeric(Mat,Mat,Mat); 62922440eb1SKris Buschelman 630be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTranspose(Mat,Mat,MatReuse,PetscReal,Mat*); 631be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTransposeSymbolic(Mat,Mat,PetscReal,Mat*); 632be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTransposeNumeric(Mat,Mat,Mat); 633bc011b1eSHong Zhang 634f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAXPY(Mat,PetscScalar,Mat,MatStructure); 63504aac2b0SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAYPX(Mat,PetscScalar,Mat,MatStructure); 6361c741599SHong Zhang 637f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScale(Mat,PetscScalar); 638f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShift(Mat,PetscScalar); 6397b80b807SBarry Smith 640be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping); 641be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetLocalToGlobalMappingBlock(Mat,ISLocalToGlobalMapping); 6422b40b63fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 6432b40b63fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsLocalIS(Mat,IS,PetscScalar,Vec,Vec); 644*56ef4db7SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsColumnsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 645*56ef4db7SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsColumnsLocalIS(Mat,IS,PetscScalar,Vec,Vec); 646be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 647be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 648052efed2SBarry Smith 649be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStashSetInitialSize(Mat,PetscInt,PetscInt); 650be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*); 65190f02eecSBarry Smith 652be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInterpolate(Mat,Vec,Vec); 6530c451bc4SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInterpolateAdd(Mat,Vec,Vec,Vec); 654ba966639SSatish Balay PetscPolymorphicSubroutine(MatInterpolateAdd,(Mat A,Vec x,Vec y),(A,x,y,y)) 655be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestrict(Mat,Vec,Vec); 656be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetVecs(Mat,Vec*,Vec*); 65769db28dcSHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRedundantMatrix(Mat,PetscInt,MPI_Comm,PetscInt,MatReuse,Mat*); 658d6037b41SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetMultiProcBlock(Mat,MPI_Comm,Mat*); 659bd481603SBarry Smith 660bd481603SBarry Smith /*MC 6611620fd73SBarry Smith MatSetValue - Set a single entry into a matrix. 6621620fd73SBarry Smith 6631620fd73SBarry Smith Not collective 6641620fd73SBarry Smith 6651620fd73SBarry Smith Input Parameters: 6661620fd73SBarry Smith + m - the matrix 6671620fd73SBarry Smith . row - the row location of the entry 6681620fd73SBarry Smith . col - the column location of the entry 6691620fd73SBarry Smith . value - the value to insert 6701620fd73SBarry Smith - mode - either INSERT_VALUES or ADD_VALUES 6711620fd73SBarry Smith 6721620fd73SBarry Smith Notes: 6731620fd73SBarry Smith For efficiency one should use MatSetValues() and set several or many 6741620fd73SBarry Smith values simultaneously if possible. 6751620fd73SBarry Smith 6761620fd73SBarry Smith Level: beginner 6771620fd73SBarry Smith 6781620fd73SBarry Smith .seealso: MatSetValues(), MatSetValueLocal() 6791620fd73SBarry Smith M*/ 6801620fd73SBarry 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);} 6811620fd73SBarry Smith 6821620fd73SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatGetValue(Mat v,PetscInt i,PetscInt j,PetscScalar *va) {return MatGetValues(v,1,&i,1,&j,va);} 6831620fd73SBarry Smith 6841620fd73SBarry 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);} 6851620fd73SBarry Smith 6861620fd73SBarry Smith /*MC 6870d2a7ff2SSatish Balay MatPreallocateInitialize - Begins the block of code that will count the number of nonzeros per 688bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 689bd481603SBarry Smith 690bd481603SBarry Smith Synopsis: 691c1ac3661SBarry Smith PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz) 692bd481603SBarry Smith 693bd481603SBarry Smith Collective on MPI_Comm 694bd481603SBarry Smith 695bd481603SBarry Smith Input Parameters: 696bd481603SBarry Smith + comm - the communicator that will share the eventually allocated matrix 697859fcb39SBarry Smith . nrows - the number of LOCAL rows in the matrix 698859fcb39SBarry Smith - ncols - the number of LOCAL columns in the matrix 699bd481603SBarry Smith 700bd481603SBarry Smith Output Parameters: 701bd481603SBarry Smith + dnz - the array that will be passed to the matrix preallocation routines 702bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 703bd481603SBarry Smith 704bd481603SBarry Smith 705bd481603SBarry Smith Level: intermediate 706bd481603SBarry Smith 707bd481603SBarry Smith Notes: 7080598bfebSBarry Smith See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details. 709bd481603SBarry Smith 7101d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 711bd481603SBarry Smith 712bd481603SBarry Smith Use MatPreallocateInitializeSymmetric() for symmetric matrices (MPISBAIJ matrices) 713bd481603SBarry Smith 7141620fd73SBarry Smith This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize(). 7151620fd73SBarry Smith 716bd481603SBarry Smith Concepts: preallocation^Matrix 717bd481603SBarry Smith 718bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(), 719bd481603SBarry Smith MatPreallocateInitializeSymmetric(), MatPreallocateSymmetricSetLocal() 720bd481603SBarry Smith M*/ 721c4f061fbSSatish Balay #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \ 7227c922b88SBarry Smith { \ 723a89bc211SBarry Smith PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__start,__end; \ 724a89bc211SBarry Smith _4_ierr = PetscMalloc2(__nrows,PetscInt,&dnz,__nrows,PetscInt,&onz);CHKERRQ(_4_ierr); \ 725a89bc211SBarry Smith _4_ierr = PetscMemzero(dnz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\ 726a89bc211SBarry Smith _4_ierr = PetscMemzero(onz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\ 727c1ac3661SBarry Smith _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ctmp;\ 728a89bc211SBarry Smith _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows; 7297c922b88SBarry Smith 730bd481603SBarry Smith /*MC 7310d2a7ff2SSatish Balay MatPreallocateSymmetricInitialize - Begins the block of code that will count the number of nonzeros per 732bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 733bd481603SBarry Smith 734bd481603SBarry Smith Synopsis: 735c1ac3661SBarry Smith PetscErrorCode MatPreallocateSymmetricInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz) 736bd481603SBarry Smith 737bd481603SBarry Smith Collective on MPI_Comm 738bd481603SBarry Smith 739bd481603SBarry Smith Input Parameters: 740bd481603SBarry Smith + comm - the communicator that will share the eventually allocated matrix 741859fcb39SBarry Smith . nrows - the number of LOCAL rows in the matrix 742859fcb39SBarry Smith - ncols - the number of LOCAL columns in the matrix 743bd481603SBarry Smith 744bd481603SBarry Smith Output Parameters: 745bd481603SBarry Smith + dnz - the array that will be passed to the matrix preallocation routines 746bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 747bd481603SBarry Smith 748bd481603SBarry Smith 749bd481603SBarry Smith Level: intermediate 750bd481603SBarry Smith 751bd481603SBarry Smith Notes: 7520598bfebSBarry Smith See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details. 753bd481603SBarry Smith 7541d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 755bd481603SBarry Smith 7561620fd73SBarry Smith This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize(). 7571620fd73SBarry Smith 758bd481603SBarry Smith Concepts: preallocation^Matrix 759bd481603SBarry Smith 760bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(), 761bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal() 762bd481603SBarry Smith M*/ 763222b16d4SBarry Smith #define MatPreallocateSymmetricInitialize(comm,nrows,ncols,dnz,onz) 0; \ 764222b16d4SBarry Smith { \ 765a89bc211SBarry Smith PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__end; \ 766a89bc211SBarry Smith _4_ierr = PetscMalloc2(__nrows,PetscInt,&dnz,__nrows,PetscInt,&onz);CHKERRQ(_4_ierr); \ 767a89bc211SBarry Smith _4_ierr = PetscMemzero(dnz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\ 768a89bc211SBarry Smith _4_ierr = PetscMemzero(onz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\ 769c1ac3661SBarry Smith _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr);\ 770a89bc211SBarry Smith _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows; 771222b16d4SBarry Smith 772bd481603SBarry Smith /*MC 773bd481603SBarry Smith MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be 774bd481603SBarry Smith inserted using a local number of the rows and columns 775bd481603SBarry Smith 776bd481603SBarry Smith Synopsis: 777c1ac3661SBarry Smith PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 778bd481603SBarry Smith 779bd481603SBarry Smith Not Collective 780bd481603SBarry Smith 781bd481603SBarry Smith Input Parameters: 782bd481603SBarry Smith + map - the mapping between local numbering and global numbering 783bd481603SBarry Smith . nrows - the number of rows indicated 7841d73ed98SBarry Smith . rows - the indices of the rows 785bd481603SBarry Smith . ncols - the number of columns in the matrix 786bd481603SBarry Smith . cols - the columns indicated 787bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 788bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 789bd481603SBarry Smith 790bd481603SBarry Smith 791bd481603SBarry Smith Level: intermediate 792bd481603SBarry Smith 793bd481603SBarry Smith Notes: 7940598bfebSBarry Smith See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details. 795bd481603SBarry Smith 7961d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 797bd481603SBarry Smith 798bd481603SBarry Smith Concepts: preallocation^Matrix 799bd481603SBarry Smith 800bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(), 801bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal() 802bd481603SBarry Smith M*/ 803c4f061fbSSatish Balay #define MatPreallocateSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\ 804c4f061fbSSatish Balay {\ 805c1ac3661SBarry Smith PetscInt __l;\ 806ef66eb69SBarry Smith _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\ 807ef66eb69SBarry Smith _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\ 808c4f061fbSSatish Balay for (__l=0;__l<nrows;__l++) {\ 809ef66eb69SBarry Smith _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 810c4f061fbSSatish Balay }\ 811c4f061fbSSatish Balay } 812c4f061fbSSatish Balay 813bd481603SBarry Smith /*MC 814bd481603SBarry Smith MatPreallocateSymmetricSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be 815bd481603SBarry Smith inserted using a local number of the rows and columns 816bd481603SBarry Smith 817bd481603SBarry Smith Synopsis: 818c1ac3661SBarry Smith PetscErrorCode MatPreallocateSymmetricSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 819bd481603SBarry Smith 820bd481603SBarry Smith Not Collective 821bd481603SBarry Smith 822bd481603SBarry Smith Input Parameters: 823bd481603SBarry Smith + map - the mapping between local numbering and global numbering 824bd481603SBarry Smith . nrows - the number of rows indicated 8251d73ed98SBarry Smith . rows - the indices of the rows 826bd481603SBarry Smith . ncols - the number of columns in the matrix 827bd481603SBarry Smith . cols - the columns indicated 828bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 829bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 830bd481603SBarry Smith 831bd481603SBarry Smith 832bd481603SBarry Smith Level: intermediate 833bd481603SBarry Smith 834bd481603SBarry Smith Notes: 8350598bfebSBarry Smith See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details. 836bd481603SBarry Smith 837bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 838bd481603SBarry Smith 839bd481603SBarry Smith Concepts: preallocation^Matrix 840bd481603SBarry Smith 841bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(), 842bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal() 843bd481603SBarry Smith M*/ 844d3d32019SBarry Smith #define MatPreallocateSymmetricSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\ 845d3d32019SBarry Smith {\ 846c1ac3661SBarry Smith PetscInt __l;\ 847d3d32019SBarry Smith _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\ 848d3d32019SBarry Smith _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\ 849d3d32019SBarry Smith for (__l=0;__l<nrows;__l++) {\ 850d3d32019SBarry Smith _4_ierr = MatPreallocateSymmetricSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 851d3d32019SBarry Smith }\ 852d3d32019SBarry Smith } 853d3d32019SBarry Smith 854bd481603SBarry Smith /*MC 855bd481603SBarry Smith MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be 856bd481603SBarry Smith inserted using a local number of the rows and columns 857bd481603SBarry Smith 858bd481603SBarry Smith Synopsis: 859c1ac3661SBarry Smith PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 860bd481603SBarry Smith 861bd481603SBarry Smith Not Collective 862bd481603SBarry Smith 863bd481603SBarry Smith Input Parameters: 86464075487SBarry Smith + row - the row 865bd481603SBarry Smith . ncols - the number of columns in the matrix 866a50f8bf6SHong Zhang - cols - the columns indicated 867a50f8bf6SHong Zhang 868a50f8bf6SHong Zhang Output Parameters: 869a50f8bf6SHong Zhang + dnz - the array that will be passed to the matrix preallocation routines 870bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 871bd481603SBarry Smith 872bd481603SBarry Smith 873bd481603SBarry Smith Level: intermediate 874bd481603SBarry Smith 875bd481603SBarry Smith Notes: 8760598bfebSBarry Smith See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details. 877bd481603SBarry Smith 878bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 879bd481603SBarry Smith 8801620fd73SBarry Smith This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize(). 8811620fd73SBarry Smith 882bd481603SBarry Smith Concepts: preallocation^Matrix 883bd481603SBarry Smith 884bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(), 885bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal() 886bd481603SBarry Smith M*/ 887c4f061fbSSatish Balay #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\ 888c1ac3661SBarry Smith { PetscInt __i; \ 889e32f2f54SBarry 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);\ 890e32f2f54SBarry 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);\ 8917c922b88SBarry Smith for (__i=0; __i<nc; __i++) {\ 89264075487SBarry Smith if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \ 8937cc688d7SBarry Smith else dnz[row - __rstart]++;\ 8947c922b88SBarry Smith }\ 8957c922b88SBarry Smith } 8967c922b88SBarry Smith 897bd481603SBarry Smith /*MC 898bd481603SBarry Smith MatPreallocateSymmetricSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be 899bd481603SBarry Smith inserted using a local number of the rows and columns 900bd481603SBarry Smith 901bd481603SBarry Smith Synopsis: 902c1ac3661SBarry Smith PetscErrorCode MatPreallocateSymmetricSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 903bd481603SBarry Smith 904bd481603SBarry Smith Not Collective 905bd481603SBarry Smith 906bd481603SBarry Smith Input Parameters: 907bd481603SBarry Smith + nrows - the number of rows indicated 9081d73ed98SBarry Smith . rows - the indices of the rows 909bd481603SBarry Smith . ncols - the number of columns in the matrix 910bd481603SBarry Smith . cols - the columns indicated 911bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 912bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 913bd481603SBarry Smith 914bd481603SBarry Smith 915bd481603SBarry Smith Level: intermediate 916bd481603SBarry Smith 917bd481603SBarry Smith Notes: 9180598bfebSBarry Smith See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details. 919bd481603SBarry Smith 920bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 921bd481603SBarry Smith 9221620fd73SBarry Smith This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize(). 9231620fd73SBarry Smith 924bd481603SBarry Smith Concepts: preallocation^Matrix 925bd481603SBarry Smith 926bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(), 927bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal() 928bd481603SBarry Smith M*/ 929d3d32019SBarry Smith #define MatPreallocateSymmetricSet(row,nc,cols,dnz,onz) 0;\ 930c1ac3661SBarry Smith { PetscInt __i; \ 931d3d32019SBarry Smith for (__i=0; __i<nc; __i++) {\ 932d3d32019SBarry Smith if (cols[__i] >= __end) onz[row - __rstart]++; \ 933d3d32019SBarry Smith else if (cols[__i] >= row) dnz[row - __rstart]++;\ 934d3d32019SBarry Smith }\ 935d3d32019SBarry Smith } 936d3d32019SBarry Smith 937bd481603SBarry Smith /*MC 93816371a99SBarry Smith MatPreallocateLocation - An alternative to MatPreallocationSet() that puts the nonzero locations into the matrix if it exists 93916371a99SBarry Smith 94016371a99SBarry Smith Synopsis: 94116371a99SBarry Smith PetscErrorCode MatPreallocateLocations(Mat A,PetscInt row,PetscInt ncols,PetscInt *cols,PetscInt *dnz,PetscInt *onz) 94216371a99SBarry Smith 94316371a99SBarry Smith Not Collective 94416371a99SBarry Smith 94516371a99SBarry Smith Input Parameters: 94616371a99SBarry Smith . A - matrix 94716371a99SBarry Smith . row - row where values exist (must be local to this process) 94816371a99SBarry Smith . ncols - number of columns 94916371a99SBarry Smith . cols - columns with nonzeros 95016371a99SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 95116371a99SBarry Smith - ozn - the other array passed to the matrix preallocation routines 95216371a99SBarry Smith 95316371a99SBarry Smith 95416371a99SBarry Smith Level: intermediate 95516371a99SBarry Smith 95616371a99SBarry Smith Notes: 9570598bfebSBarry Smith See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details. 95816371a99SBarry Smith 95916371a99SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 96016371a99SBarry Smith 96116371a99SBarry Smith This is a MACRO not a function because it uses a bunch of variables private to the MatPreallocation.... routines. 96216371a99SBarry Smith 96316371a99SBarry Smith Concepts: preallocation^Matrix 96416371a99SBarry Smith 96516371a99SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(), 96616371a99SBarry Smith MatPreallocateSymmetricInitialize(), MatPreallocateSymmetricSetLocal() 96716371a99SBarry Smith M*/ 96816371a99SBarry 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);} 96916371a99SBarry Smith 97016371a99SBarry Smith 97116371a99SBarry Smith /*MC 9720d2a7ff2SSatish Balay MatPreallocateFinalize - Ends the block of code that will count the number of nonzeros per 973bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 974bd481603SBarry Smith 975bd481603SBarry Smith Synopsis: 976c1ac3661SBarry Smith PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz) 977bd481603SBarry Smith 978bd481603SBarry Smith Collective on MPI_Comm 979bd481603SBarry Smith 980bd481603SBarry Smith Input Parameters: 98116371a99SBarry Smith + dnz - the array that was be passed to the matrix preallocation routines 982bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 983bd481603SBarry Smith 984bd481603SBarry Smith 985bd481603SBarry Smith Level: intermediate 986bd481603SBarry Smith 987bd481603SBarry Smith Notes: 9880598bfebSBarry Smith See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details. 989bd481603SBarry Smith 990bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 991bd481603SBarry Smith 9921620fd73SBarry Smith This is a MACRO not a function because it closes the { started in MatPreallocateInitialize(). 9931620fd73SBarry Smith 994bd481603SBarry Smith Concepts: preallocation^Matrix 995bd481603SBarry Smith 996bd481603SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(), 997bd481603SBarry Smith MatPreallocateSymmetricInitialize(), MatPreallocateSymmetricSetLocal() 998bd481603SBarry Smith M*/ 999a89bc211SBarry Smith #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree2(dnz,onz);CHKERRQ(_4_ierr);} 10007c922b88SBarry Smith 1001bd481603SBarry Smith 1002bd481603SBarry Smith 10037b80b807SBarry Smith /* Routines unique to particular data structures */ 1004be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellGetContext(Mat,void **); 1005ba966639SSatish Balay PetscPolymorphicFunction(MatShellGetContext,(Mat A),(A,&t),void*,t) 1006698d4c6aSKris Buschelman 1007be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInodeAdjustForInodes(Mat,IS*,IS*); 1008be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *); 1009698d4c6aSKris Buschelman 1010be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetColumnIndices(Mat,PetscInt[]); 1011be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetColumnIndices(Mat,PetscInt[]); 1012be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 1013c75a6043SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 1014c75a6043SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 10157b80b807SBarry Smith 1016a96a251dSBarry Smith #define MAT_SKIP_ALLOCATION -4 1017a96a251dSBarry Smith 1018be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]); 1019ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz)) 1020be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]); 1021ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqSBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz)) 1022be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]); 1023ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqAIJSetPreallocation,(Mat A,const PetscInt nnz[]),(A,0,nnz)) 1024273d9f13SBarry Smith 1025be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 1026ba966639SSatish Balay PetscPolymorphicSubroutine(MatMPIBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[],const PetscInt onz[]),(A,bs,0,nnz,0,onz)) 1027be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 1028be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 102953803f26SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []); 1030725b52f3SLisandro Dalcin EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 1031be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]); 1032be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 1033be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]); 1034be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIDenseSetPreallocation(Mat,PetscScalar[]); 1035284134d9SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation(Mat,PetscScalar[]); 1036be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,PetscInt*[]); 1037be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,PetscInt*[]); 1038be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAdicSetLocalFunction(Mat,void (*)(void)); 1039273d9f13SBarry Smith 1040be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetLDA(Mat,PetscInt); 1041ab92ecdeSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDenseGetLocalMatrix(Mat,Mat*); 10421b807ce4Svictorle 1043be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStoreValues(Mat); 1044be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRetrieveValues(Mat); 10452e8a6d31SBarry Smith 1046be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDAADSetCtx(Mat,void*); 10473a7fca6bSBarry Smith 10487b80b807SBarry Smith /* 10497b80b807SBarry Smith These routines are not usually accessed directly, rather solving is 105094b7f48cSBarry Smith done through the KSP and PC interfaces. 10517b80b807SBarry Smith */ 10527b80b807SBarry Smith 1053d9274352SBarry Smith /*E 1054d9274352SBarry Smith MatOrderingType - String with the name of a PETSc matrix ordering or the creation function 1055d9274352SBarry Smith with an optional dynamic library name, for example 1056d9274352SBarry Smith http://www.mcs.anl.gov/petsc/lib.a:orderingcreate() 1057d9274352SBarry Smith 1058d9274352SBarry Smith Level: beginner 1059d9274352SBarry Smith 1060e5a9bf91SBarry Smith Cannot use const because the PC objects manipulate the string 1061e5a9bf91SBarry Smith 1062d9274352SBarry Smith .seealso: MatGetOrdering() 1063d9274352SBarry Smith E*/ 10643eaad2caSSatish Balay #define MatOrderingType char* 10652692d6eeSBarry Smith #define MATORDERINGNATURAL "natural" 10662692d6eeSBarry Smith #define MATORDERINGND "nd" 10672692d6eeSBarry Smith #define MATORDERING1WD "1wd" 10682692d6eeSBarry Smith #define MATORDERINGRCM "rcm" 10692692d6eeSBarry Smith #define MATORDERINGQMD "qmd" 10702692d6eeSBarry Smith #define MATORDERINGROWLENGTH "rowlength" 10712692d6eeSBarry Smith #define MATORDERINGDSC_ND "dsc_nd" /* these three are only for DSCPACK, see its documentation for details */ 10722692d6eeSBarry Smith #define MATORDERINGDSC_MMD "dsc_mmd" 10732692d6eeSBarry Smith #define MATORDERINGDSC_MDF "dsc_mdf" 1074312e372aSBarry Smith #define MATORDERINGAMD "amd" /* only works if UMFPACK is installed with PETSc */ 1075b12f92e5SBarry Smith 1076f2416520SSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOrdering(Mat,const MatOrderingType,IS*,IS*); 1077be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOrderingList(PetscFList *list); 1078f2416520SSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,const MatOrderingType,IS*,IS*)); 107930de9b25SBarry Smith 108030de9b25SBarry Smith /*MC 10811890ba74SBarry Smith MatOrderingRegisterDynamic - Adds a new sparse matrix ordering to the matrix package. 108230de9b25SBarry Smith 108330de9b25SBarry Smith Synopsis: 10841890ba74SBarry Smith PetscErrorCode MatOrderingRegisterDynamic(const char *name_ordering,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatOrdering)) 108530de9b25SBarry Smith 108630de9b25SBarry Smith Not Collective 108730de9b25SBarry Smith 108830de9b25SBarry Smith Input Parameters: 10892692d6eeSBarry Smith + sname - name of ordering (for example MATORDERINGND) 109030de9b25SBarry Smith . path - location of library where creation routine is 109130de9b25SBarry Smith . name - name of function that creates the ordering type,a string 109230de9b25SBarry Smith - function - function pointer that creates the ordering 109330de9b25SBarry Smith 109430de9b25SBarry Smith Level: developer 109530de9b25SBarry Smith 109630de9b25SBarry Smith If dynamic libraries are used, then the fourth input argument (function) 109730de9b25SBarry Smith is ignored. 109830de9b25SBarry Smith 109930de9b25SBarry Smith Sample usage: 110030de9b25SBarry Smith .vb 110130de9b25SBarry Smith MatOrderingRegisterDynamic("my_order",/home/username/my_lib/lib/libO/solaris/mylib.a, 110230de9b25SBarry Smith "MyOrder",MyOrder); 110330de9b25SBarry Smith .ve 110430de9b25SBarry Smith 110530de9b25SBarry Smith Then, your partitioner can be chosen with the procedural interface via 110630de9b25SBarry Smith $ MatOrderingSetType(part,"my_order) 110730de9b25SBarry Smith or at runtime via the option 11082401956bSBarry Smith $ -pc_factor_mat_ordering_type my_order 110930de9b25SBarry Smith 1110ab901514SBarry Smith ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values. 111130de9b25SBarry Smith 111230de9b25SBarry Smith .keywords: matrix, ordering, register 111330de9b25SBarry Smith 111430de9b25SBarry Smith .seealso: MatOrderingRegisterDestroy(), MatOrderingRegisterAll() 111530de9b25SBarry Smith M*/ 1116aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 1117f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,0) 1118b12f92e5SBarry Smith #else 1119f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,d) 1120b12f92e5SBarry Smith #endif 112130de9b25SBarry Smith 1122be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegisterDestroy(void); 1123be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegisterAll(const char[]); 1124ace3abfcSBarry Smith extern PetscBool MatOrderingRegisterAllCalled; 1125b0a32e0cSBarry Smith extern PetscFList MatOrderingList; 1126d4fbbf0eSBarry Smith 1127be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS); 1128a2ce50c7SBarry Smith 1129d91e6319SBarry Smith /*S 1130d90ac83dSHong Zhang MatFactorShiftType - Numeric Shift. 1131d90ac83dSHong Zhang 1132d90ac83dSHong Zhang Level: beginner 1133d90ac83dSHong Zhang 1134d90ac83dSHong Zhang S*/ 1135d90ac83dSHong Zhang typedef enum {MAT_SHIFT_NONE,MAT_SHIFT_NONZERO,MAT_SHIFT_POSITIVE_DEFINITE,MAT_SHIFT_INBLOCKS} MatFactorShiftType; 1136d90ac83dSHong Zhang extern const char *MatFactorShiftTypes[]; 1137d90ac83dSHong Zhang 1138d90ac83dSHong Zhang /*S 11392401956bSBarry Smith MatFactorInfo - Data passed into the matrix factorization routines 1140b00f7748SHong Zhang 114161cad860SBarry Smith In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE, that is use 114261cad860SBarry Smith $ MatFactorInfo info(MAT_FACTORINFO_SIZE) 1143b00f7748SHong Zhang 114415e8a5b3SHong Zhang Notes: These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC. 1145b00f7748SHong Zhang 114661cad860SBarry Smith You can use MatFactorInfoInitialize() to set default values. 114761cad860SBarry Smith 1148b00f7748SHong Zhang Level: developer 1149b00f7748SHong Zhang 1150d7d82daaSBarry Smith .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(), 1151d7d82daaSBarry Smith MatFactorInfoInitialize() 1152b00f7748SHong Zhang 1153b00f7748SHong Zhang S*/ 1154b00f7748SHong Zhang typedef struct { 115515e8a5b3SHong Zhang PetscReal diagonal_fill; /* force diagonal to fill in if initially not filled */ 115685317021SBarry Smith PetscReal usedt; 115715e8a5b3SHong Zhang PetscReal dt; /* drop tolerance */ 1158b00f7748SHong Zhang PetscReal dtcol; /* tolerance for pivoting */ 115915e8a5b3SHong Zhang PetscReal dtcount; /* maximum nonzeros to be allowed per row */ 116067e28bfeSBarry Smith PetscReal fill; /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */ 1161348344bbSBarry Smith PetscReal levels; /* ICC/ILU(levels) */ 1162bcd9e38bSBarry Smith PetscReal pivotinblocks; /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0 1163bcd9e38bSBarry Smith factorization may be faster if do not pivot */ 116415e8a5b3SHong Zhang PetscReal zeropivot; /* pivot is called zero if less than this */ 1165f4db908eSBarry Smith PetscReal shifttype; /* type of shift added to matrix factor to prevent zero pivots */ 1166f4db908eSBarry Smith PetscReal shiftamount; /* how large the shift is */ 116715e8a5b3SHong Zhang } MatFactorInfo; 1168ffa6d0a5SLois Curfman McInnes 1169be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFactorInfoInitialize(MatFactorInfo*); 11700481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactor(Mat,IS,const MatFactorInfo*); 11710481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*); 11720481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactorNumeric(Mat,Mat,const MatFactorInfo*); 11730481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactor(Mat,IS,IS,const MatFactorInfo*); 11740481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUFactor(Mat,IS,IS,const MatFactorInfo*); 11750481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*); 11760481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*); 11770481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatICCFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*); 11780481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatICCFactor(Mat,IS,const MatFactorInfo*); 11790481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactorNumeric(Mat,Mat,const MatFactorInfo*); 1180be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*); 1181be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolve(Mat,Vec,Vec); 1182be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatForwardSolve(Mat,Vec,Vec); 1183be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatBackwardSolve(Mat,Vec,Vec); 1184be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveAdd(Mat,Vec,Vec,Vec); 1185be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveTranspose(Mat,Vec,Vec); 1186be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveTransposeAdd(Mat,Vec,Vec,Vec); 1187be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolves(Mat,Vecs,Vecs); 11888ed539a5SBarry Smith 1189be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUnfactored(Mat); 1190bb5a7306SBarry Smith 1191d91e6319SBarry Smith /*E 1192d91e6319SBarry Smith MatSORType - What type of (S)SOR to perform 1193bb1eb677SSatish Balay 1194d91e6319SBarry Smith Level: beginner 1195d91e6319SBarry Smith 1196d9274352SBarry Smith May be bitwise ORd together 1197d9274352SBarry Smith 1198d91e6319SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 1199d91e6319SBarry Smith 12004e7234bfSBarry Smith MatSORType may be bitwise ORd together, so do not change the numbers 12014e7234bfSBarry Smith 120241f059aeSBarry Smith .seealso: MatSOR() 1203d91e6319SBarry Smith E*/ 1204ee50ffe9SBarry Smith typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3, 1205ee50ffe9SBarry Smith SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8, 1206ee50ffe9SBarry Smith SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16, 120784cb2905SBarry Smith SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType; 120841f059aeSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSOR(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec); 12098ed539a5SBarry Smith 1210d4fbbf0eSBarry Smith /* 1211639f9d9dSBarry Smith These routines are for efficiently computing Jacobians via finite differences. 1212639f9d9dSBarry Smith */ 1213b12f92e5SBarry Smith 1214d9274352SBarry Smith /*E 1215d9274352SBarry Smith MatColoringType - String with the name of a PETSc matrix coloring or the creation function 1216d9274352SBarry Smith with an optional dynamic library name, for example 1217d9274352SBarry Smith http://www.mcs.anl.gov/petsc/lib.a:coloringcreate() 1218d9274352SBarry Smith 1219d9274352SBarry Smith Level: beginner 1220d9274352SBarry Smith 1221d9274352SBarry Smith .seealso: MatGetColoring() 1222d9274352SBarry Smith E*/ 1223a313700dSBarry Smith #define MatColoringType char* 12242692d6eeSBarry Smith #define MATCOLORINGNATURAL "natural" 12252692d6eeSBarry Smith #define MATCOLORINGSL "sl" 12262692d6eeSBarry Smith #define MATCOLORINGLF "lf" 12272692d6eeSBarry Smith #define MATCOLORINGID "id" 1228b12f92e5SBarry Smith 1229ab970d9bSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColoring(Mat,const MatColoringType,ISColoring*); 12302e90c967SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,MatColoringType,ISColoring *)); 123130de9b25SBarry Smith 123230de9b25SBarry Smith /*MC 123330de9b25SBarry Smith MatColoringRegisterDynamic - Adds a new sparse matrix coloring to the 123430de9b25SBarry Smith matrix package. 123530de9b25SBarry Smith 123630de9b25SBarry Smith Synopsis: 12371890ba74SBarry Smith PetscErrorCode MatColoringRegisterDynamic(const char *name_coloring,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatColoring)) 123830de9b25SBarry Smith 123930de9b25SBarry Smith Not Collective 124030de9b25SBarry Smith 124130de9b25SBarry Smith Input Parameters: 12422692d6eeSBarry Smith + sname - name of Coloring (for example MATCOLORINGSL) 124330de9b25SBarry Smith . path - location of library where creation routine is 124430de9b25SBarry Smith . name - name of function that creates the Coloring type, a string 124530de9b25SBarry Smith - function - function pointer that creates the coloring 124630de9b25SBarry Smith 124730de9b25SBarry Smith Level: developer 124830de9b25SBarry Smith 124930de9b25SBarry Smith If dynamic libraries are used, then the fourth input argument (function) 125030de9b25SBarry Smith is ignored. 125130de9b25SBarry Smith 125230de9b25SBarry Smith Sample usage: 125330de9b25SBarry Smith .vb 125430de9b25SBarry Smith MatColoringRegisterDynamic("my_color",/home/username/my_lib/lib/libO/solaris/mylib.a, 125530de9b25SBarry Smith "MyColor",MyColor); 125630de9b25SBarry Smith .ve 125730de9b25SBarry Smith 125830de9b25SBarry Smith Then, your partitioner can be chosen with the procedural interface via 125930de9b25SBarry Smith $ MatColoringSetType(part,"my_color") 126030de9b25SBarry Smith or at runtime via the option 126130de9b25SBarry Smith $ -mat_coloring_type my_color 126230de9b25SBarry Smith 1263ab901514SBarry Smith $PETSC_ARCH occuring in pathname will be replaced with appropriate values. 126430de9b25SBarry Smith 126530de9b25SBarry Smith .keywords: matrix, Coloring, register 126630de9b25SBarry Smith 126730de9b25SBarry Smith .seealso: MatColoringRegisterDestroy(), MatColoringRegisterAll() 126830de9b25SBarry Smith M*/ 1269aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 1270f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,0) 1271b12f92e5SBarry Smith #else 1272f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,d) 1273b12f92e5SBarry Smith #endif 127430de9b25SBarry Smith 1275ace3abfcSBarry Smith extern PetscBool MatColoringRegisterAllCalled; 1276f1144a30SSatish Balay 1277be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegisterAll(const char[]); 1278be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegisterDestroy(void); 1279be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*); 1280639f9d9dSBarry Smith 1281d9274352SBarry Smith /*S 1282d9274352SBarry Smith MatFDColoring - Object for computing a sparse Jacobian via finite differences 1283d9274352SBarry Smith and coloring 1284639f9d9dSBarry Smith 1285d9274352SBarry Smith Level: beginner 1286d9274352SBarry Smith 1287d9274352SBarry Smith Concepts: coloring, sparse Jacobian, finite differences 1288d9274352SBarry Smith 1289d9274352SBarry Smith .seealso: MatFDColoringCreate() 1290d9274352SBarry Smith S*/ 1291e2a1c21fSSatish Balay typedef struct _p_MatFDColoring* MatFDColoring; 1292639f9d9dSBarry Smith 1293be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringCreate(Mat,ISColoring,MatFDColoring *); 1294be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringDestroy(MatFDColoring); 1295be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringView(MatFDColoring,PetscViewer); 1296be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*); 12974a9d489dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**); 1298be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal); 1299be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFromOptions(MatFDColoring); 1300be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApply(Mat,MatFDColoring,Vec,MatStructure*,void *); 1301be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApplyTS(Mat,MatFDColoring,PetscReal,Vec,MatStructure*,void *); 1302be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetF(MatFDColoring,Vec); 1303be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,PetscInt*[]); 1304639f9d9dSBarry Smith /* 13050752156aSBarry Smith These routines are for partitioning matrices: currently used only 13063eda8832SBarry Smith for adjacency matrix, MatCreateMPIAdj(). 13070752156aSBarry Smith */ 1308ca161407SBarry Smith 1309d9274352SBarry Smith /*S 1310d9274352SBarry Smith MatPartitioning - Object for managing the partitioning of a matrix or graph 1311d9274352SBarry Smith 1312d9274352SBarry Smith Level: beginner 1313d9274352SBarry Smith 1314d9274352SBarry Smith Concepts: partitioning 1315d9274352SBarry Smith 1316743a1026Svictorle .seealso: MatPartitioningCreate(), MatPartitioningType 1317d9274352SBarry Smith S*/ 131891e9ee9fSBarry Smith typedef struct _p_MatPartitioning* MatPartitioning; 1319d9274352SBarry Smith 1320d9274352SBarry Smith /*E 13215bc6e7ccSvictorle MatPartitioningType - String with the name of a PETSc matrix partitioning or the creation function 1322d9274352SBarry Smith with an optional dynamic library name, for example 1323d9274352SBarry Smith http://www.mcs.anl.gov/petsc/lib.a:partitioningcreate() 1324d9274352SBarry Smith 1325d9274352SBarry Smith Level: beginner 1326d9274352SBarry Smith 1327b547a13bSvictorle .seealso: MatPartitioningCreate(), MatPartitioning 1328d9274352SBarry Smith E*/ 1329a313700dSBarry Smith #define MatPartitioningType char* 13302692d6eeSBarry Smith #define MATPARTITIONINGCURRENT "current" 13312692d6eeSBarry Smith #define MATPARTITIONINGSQUARE "square" 13322692d6eeSBarry Smith #define MATPARTITIONINGPARMETIS "parmetis" 13332692d6eeSBarry Smith #define MATPARTITIONINGCHACO "chaco" 13342692d6eeSBarry Smith #define MATPARTITIONINGJOSTLE "jostle" 13352692d6eeSBarry Smith #define MATPARTITIONINGPARTY "party" 13362692d6eeSBarry Smith #define MATPARTITIONINGSCOTCH "scotch" 133771306c60Sjroman 1338ca161407SBarry Smith 1339be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningCreate(MPI_Comm,MatPartitioning*); 1340a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetType(MatPartitioning,const MatPartitioningType); 1341be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetNParts(MatPartitioning,PetscInt); 1342be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetAdjacency(MatPartitioning,Mat); 1343be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]); 1344be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []); 1345be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningApply(MatPartitioning,IS*); 1346be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningDestroy(MatPartitioning); 13472aabb6bbSBarry Smith 1348be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatPartitioning)); 134930de9b25SBarry Smith 135030de9b25SBarry Smith /*MC 135130de9b25SBarry Smith MatPartitioningRegisterDynamic - Adds a new sparse matrix partitioning to the 135230de9b25SBarry Smith matrix package. 135330de9b25SBarry Smith 135430de9b25SBarry Smith Synopsis: 13551890ba74SBarry Smith PetscErrorCode MatPartitioningRegisterDynamic(const char *name_partitioning,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatPartitioning)) 135630de9b25SBarry Smith 135730de9b25SBarry Smith Not Collective 135830de9b25SBarry Smith 135930de9b25SBarry Smith Input Parameters: 13602692d6eeSBarry Smith + sname - name of partitioning (for example MATPARTITIONINGCURRENT) or parmetis 136130de9b25SBarry Smith . path - location of library where creation routine is 136230de9b25SBarry Smith . name - name of function that creates the partitioning type, a string 136330de9b25SBarry Smith - function - function pointer that creates the partitioning type 136430de9b25SBarry Smith 136530de9b25SBarry Smith Level: developer 136630de9b25SBarry Smith 136730de9b25SBarry Smith If dynamic libraries are used, then the fourth input argument (function) 136830de9b25SBarry Smith is ignored. 136930de9b25SBarry Smith 137030de9b25SBarry Smith Sample usage: 137130de9b25SBarry Smith .vb 137230de9b25SBarry Smith MatPartitioningRegisterDynamic("my_part",/home/username/my_lib/lib/libO/solaris/mylib.a, 137330de9b25SBarry Smith "MyPartCreate",MyPartCreate); 137430de9b25SBarry Smith .ve 137530de9b25SBarry Smith 137630de9b25SBarry Smith Then, your partitioner can be chosen with the procedural interface via 137730de9b25SBarry Smith $ MatPartitioningSetType(part,"my_part") 137830de9b25SBarry Smith or at runtime via the option 137930de9b25SBarry Smith $ -mat_partitioning_type my_part 138030de9b25SBarry Smith 1381ab901514SBarry Smith $PETSC_ARCH occuring in pathname will be replaced with appropriate values. 138230de9b25SBarry Smith 138330de9b25SBarry Smith .keywords: matrix, partitioning, register 138430de9b25SBarry Smith 138530de9b25SBarry Smith .seealso: MatPartitioningRegisterDestroy(), MatPartitioningRegisterAll() 138630de9b25SBarry Smith M*/ 1387aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 1388f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,0) 13892aabb6bbSBarry Smith #else 1390f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,d) 13912aabb6bbSBarry Smith #endif 13922aabb6bbSBarry Smith 1393ace3abfcSBarry Smith extern PetscBool MatPartitioningRegisterAllCalled; 1394f1144a30SSatish Balay 1395be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegisterAll(const char[]); 1396be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegisterDestroy(void); 13972bad1931SBarry Smith 1398be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningView(MatPartitioning,PetscViewer); 1399be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetFromOptions(MatPartitioning); 1400a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningGetType(MatPartitioning,const MatPartitioningType*); 1401ca161407SBarry Smith 1402be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningParmetisSetCoarseSequential(MatPartitioning); 14030e8af3caSMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningParmetisGetEdgeCut(MatPartitioning, PetscInt *); 14040752156aSBarry Smith 1405be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningJostleSetCoarseLevel(MatPartitioning,PetscReal); 1406be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningJostleSetCoarseSequential(MatPartitioning); 140771306c60Sjroman 1408954b3ec6SBarry Smith typedef enum {MP_CHACO_MULTILEVEL_KL,MP_CHACO_SPECTRAL,MP_CHACO_LINEAR,MP_CHACO_RANDOM, MP_CHACO_SCATTERED} MPChacoGlobalType; 1409be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetGlobal(MatPartitioning, MPChacoGlobalType); 141071306c60Sjroman typedef enum { MP_CHACO_KERNIGHAN_LIN, MP_CHACO_NONE } MPChacoLocalType; 1411be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetLocal(MatPartitioning, MPChacoLocalType); 1412be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal); 141371306c60Sjroman typedef enum { MP_CHACO_LANCZOS, MP_CHACO_RQI_SYMMLQ } MPChacoEigenType; 1414be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType); 1415be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenTol(MatPartitioning, PetscReal); 1416be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenNumber(MatPartitioning, PetscInt); 141771306c60Sjroman 141871306c60Sjroman #define MP_PARTY_OPT "opt" 141971306c60Sjroman #define MP_PARTY_LIN "lin" 142071306c60Sjroman #define MP_PARTY_SCA "sca" 142171306c60Sjroman #define MP_PARTY_RAN "ran" 142271306c60Sjroman #define MP_PARTY_GBF "gbf" 142371306c60Sjroman #define MP_PARTY_GCF "gcf" 142471306c60Sjroman #define MP_PARTY_BUB "bub" 142571306c60Sjroman #define MP_PARTY_DEF "def" 1426be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetGlobal(MatPartitioning, const char*); 142771306c60Sjroman #define MP_PARTY_HELPFUL_SETS "hs" 142871306c60Sjroman #define MP_PARTY_KERNIGHAN_LIN "kl" 142971306c60Sjroman #define MP_PARTY_NONE "no" 1430be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetLocal(MatPartitioning, const char*); 1431be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal); 1432ace3abfcSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetBipart(MatPartitioning,PetscBool ); 1433ace3abfcSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscBool ); 143471306c60Sjroman 143571306c60Sjroman typedef enum { MP_SCOTCH_GREEDY, MP_SCOTCH_GPS, MP_SCOTCH_GR_GPS } MPScotchGlobalType; 1436be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetArch(MatPartitioning,const char*); 1437be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetMultilevel(MatPartitioning); 1438be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetGlobal(MatPartitioning,MPScotchGlobalType); 1439be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetCoarseLevel(MatPartitioning,PetscReal); 1440be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetHostList(MatPartitioning,const char*); 144171306c60Sjroman typedef enum { MP_SCOTCH_KERNIGHAN_LIN, MP_SCOTCH_NONE } MPScotchLocalType; 1442be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetLocal(MatPartitioning,MPScotchLocalType); 1443be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetMapping(MatPartitioning); 1444be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetStrategy(MatPartitioning,char*); 144571306c60Sjroman 1446591294e4SBarry Smith EXTERN PetscErrorCode MatMeshToVertexGraph(Mat,PetscInt,Mat*); 144714196391SBarry Smith EXTERN PetscErrorCode MatMeshToCellGraph(Mat,PetscInt,Mat*); 1448591294e4SBarry Smith 14490752156aSBarry Smith /* 14500a835dfdSSatish Balay If you add entries here you must also add them to finclude/petscmat.h 1451d4fbbf0eSBarry Smith */ 14521c1c02c0SLois Curfman McInnes typedef enum { MATOP_SET_VALUES=0, 14531c1c02c0SLois Curfman McInnes MATOP_GET_ROW=1, 14541c1c02c0SLois Curfman McInnes MATOP_RESTORE_ROW=2, 14551c1c02c0SLois Curfman McInnes MATOP_MULT=3, 14561c1c02c0SLois Curfman McInnes MATOP_MULT_ADD=4, 14577c922b88SBarry Smith MATOP_MULT_TRANSPOSE=5, 14587c922b88SBarry Smith MATOP_MULT_TRANSPOSE_ADD=6, 14591c1c02c0SLois Curfman McInnes MATOP_SOLVE=7, 14601c1c02c0SLois Curfman McInnes MATOP_SOLVE_ADD=8, 14617c922b88SBarry Smith MATOP_SOLVE_TRANSPOSE=9, 14627c922b88SBarry Smith MATOP_SOLVE_TRANSPOSE_ADD=10, 14631c1c02c0SLois Curfman McInnes MATOP_LUFACTOR=11, 14641c1c02c0SLois Curfman McInnes MATOP_CHOLESKYFACTOR=12, 1465a714c06dSBarry Smith MATOP_SOR=13, 14661c1c02c0SLois Curfman McInnes MATOP_TRANSPOSE=14, 14671c1c02c0SLois Curfman McInnes MATOP_GETINFO=15, 14681c1c02c0SLois Curfman McInnes MATOP_EQUAL=16, 14691c1c02c0SLois Curfman McInnes MATOP_GET_DIAGONAL=17, 14701c1c02c0SLois Curfman McInnes MATOP_DIAGONAL_SCALE=18, 14711c1c02c0SLois Curfman McInnes MATOP_NORM=19, 14721c1c02c0SLois Curfman McInnes MATOP_ASSEMBLY_BEGIN=20, 14731c1c02c0SLois Curfman McInnes MATOP_ASSEMBLY_END=21, 1474d519adbfSMatthew Knepley MATOP_SET_OPTION=22, 1475d519adbfSMatthew Knepley MATOP_ZERO_ENTRIES=23, 1476d519adbfSMatthew Knepley MATOP_ZERO_ROWS=24, 1477d519adbfSMatthew Knepley MATOP_LUFACTOR_SYMBOLIC=25, 1478d519adbfSMatthew Knepley MATOP_LUFACTOR_NUMERIC=26, 1479d519adbfSMatthew Knepley MATOP_CHOLESKY_FACTOR_SYMBOLIC=27, 1480d519adbfSMatthew Knepley MATOP_CHOLESKY_FACTOR_NUMERIC=28, 1481d519adbfSMatthew Knepley MATOP_SETUP_PREALLOCATION=29, 1482d519adbfSMatthew Knepley MATOP_ILUFACTOR_SYMBOLIC=30, 1483d519adbfSMatthew Knepley MATOP_ICCFACTOR_SYMBOLIC=31, 1484d519adbfSMatthew Knepley MATOP_GET_ARRAY=32, 1485d519adbfSMatthew Knepley MATOP_RESTORE_ARRAY=33, 1486d519adbfSMatthew Knepley MATOP_DUPLICATE=34, 1487d519adbfSMatthew Knepley MATOP_FORWARD_SOLVE=35, 1488d519adbfSMatthew Knepley MATOP_BACKWARD_SOLVE=36, 1489d519adbfSMatthew Knepley MATOP_ILUFACTOR=37, 1490d519adbfSMatthew Knepley MATOP_ICCFACTOR=38, 1491d519adbfSMatthew Knepley MATOP_AXPY=39, 1492d519adbfSMatthew Knepley MATOP_GET_SUBMATRICES=40, 1493d519adbfSMatthew Knepley MATOP_INCREASE_OVERLAP=41, 1494d519adbfSMatthew Knepley MATOP_GET_VALUES=42, 1495d519adbfSMatthew Knepley MATOP_COPY=43, 1496d519adbfSMatthew Knepley MATOP_GET_ROW_MAX=44, 1497d519adbfSMatthew Knepley MATOP_SCALE=45, 1498d519adbfSMatthew Knepley MATOP_SHIFT=46, 149935153367SBarry Smith MATOP_DIAGONAL_SET=47, 1500d519adbfSMatthew Knepley MATOP_ILUDT_FACTOR=48, 1501d519adbfSMatthew Knepley MATOP_SET_BLOCK_SIZE=49, 1502d519adbfSMatthew Knepley MATOP_GET_ROW_IJ=50, 1503d519adbfSMatthew Knepley MATOP_RESTORE_ROW_IJ=51, 1504d519adbfSMatthew Knepley MATOP_GET_COLUMN_IJ=52, 1505d519adbfSMatthew Knepley MATOP_RESTORE_COLUMN_IJ=53, 1506d519adbfSMatthew Knepley MATOP_FDCOLORING_CREATE=54, 1507d519adbfSMatthew Knepley MATOP_COLORING_PATCH=55, 1508d519adbfSMatthew Knepley MATOP_SET_UNFACTORED=56, 1509d519adbfSMatthew Knepley MATOP_PERMUTE=57, 1510d519adbfSMatthew Knepley MATOP_SET_VALUES_BLOCKED=58, 1511d519adbfSMatthew Knepley MATOP_GET_SUBMATRIX=59, 1512d519adbfSMatthew Knepley MATOP_DESTROY=60, 1513d519adbfSMatthew Knepley MATOP_VIEW=61, 1514d519adbfSMatthew Knepley MATOP_CONVERT_FROM=62, 1515d519adbfSMatthew Knepley MATOP_USE_SCALED_FORM=63, 1516d519adbfSMatthew Knepley MATOP_SCALE_SYSTEM=64, 1517d519adbfSMatthew Knepley MATOP_UNSCALE_SYSTEM=65, 1518d519adbfSMatthew Knepley MATOP_SET_LOCAL_TO_GLOBAL_MAP=66, 1519d519adbfSMatthew Knepley MATOP_SET_VALUES_LOCAL=67, 1520d519adbfSMatthew Knepley MATOP_ZERO_ROWS_LOCAL=68, 1521d519adbfSMatthew Knepley MATOP_GET_ROW_MAX_ABS=69, 1522d519adbfSMatthew Knepley MATOP_GET_ROW_MIN_ABS=70, 1523d519adbfSMatthew Knepley MATOP_CONVERT=71, 1524d519adbfSMatthew Knepley MATOP_SET_COLORING=72, 1525d519adbfSMatthew Knepley MATOP_SET_VALUES_ADIC=73, 1526d519adbfSMatthew Knepley MATOP_SET_VALUES_ADIFOR=74, 1527d519adbfSMatthew Knepley MATOP_FD_COLORING_APPLY=75, 1528d519adbfSMatthew Knepley MATOP_SET_FROM_OPTIONS=76, 1529d519adbfSMatthew Knepley MATOP_MULT_CON=77, 1530d519adbfSMatthew Knepley MATOP_MULT_TRANSPOSE_CON=78, 1531d519adbfSMatthew Knepley MATOP_PERMUTE_SPARSIFY=79, 1532d519adbfSMatthew Knepley MATOP_MULT_MULTIPLE=80, 1533d519adbfSMatthew Knepley MATOP_SOLVE_MULTIPLE=81, 1534d519adbfSMatthew Knepley MATOP_GET_INERTIA=82, 1535d519adbfSMatthew Knepley MATOP_LOAD=83, 1536d519adbfSMatthew Knepley MATOP_IS_SYMMETRIC=84, 1537d519adbfSMatthew Knepley MATOP_IS_HERMITIAN=85, 1538d519adbfSMatthew Knepley MATOP_IS_STRUCTURALLY_SYMMETRIC=86, 153941f059aeSBarry Smith MATOP_DUMMY=87, 1540d519adbfSMatthew Knepley MATOP_GET_VECS=88, 1541d519adbfSMatthew Knepley MATOP_MAT_MULT=89, 1542d519adbfSMatthew Knepley MATOP_MAT_MULT_SYMBOLIC=90, 1543d519adbfSMatthew Knepley MATOP_MAT_MULT_NUMERIC=91, 1544d519adbfSMatthew Knepley MATOP_PTAP=92, 1545d519adbfSMatthew Knepley MATOP_PTAP_SYMBOLIC=93, 1546d519adbfSMatthew Knepley MATOP_PTAP_NUMERIC=94, 1547d519adbfSMatthew Knepley MATOP_MAT_MULTTRANSPOSE=95, 15481763ac4eSSatish Balay MATOP_MAT_MULTTRANSPOSE_SYM=96, 15491763ac4eSSatish Balay MATOP_MAT_MULTTRANSPOSE_NUM=97, 1550d519adbfSMatthew Knepley MATOP_PTAP_SYMBOLIC_SEQAIJ=98, 1551d519adbfSMatthew Knepley MATOP_PTAP_NUMERIC_SEQAIJ=99, 1552d519adbfSMatthew Knepley MATOP_PTAP_SYMBOLIC_MPIAIJ=100, 1553d519adbfSMatthew Knepley MATOP_PTAP_NUMERIC_MPIAIJ=101, 1554d519adbfSMatthew Knepley MATOP_CONJUGATE=102, 1555d519adbfSMatthew Knepley MATOP_SET_SIZES=103, 1556d519adbfSMatthew Knepley MATOP_SET_VALUES_ROW=104, 1557d519adbfSMatthew Knepley MATOP_REAL_PART=105, 1558d519adbfSMatthew Knepley MATOP_IMAG_PART=106, 1559d519adbfSMatthew Knepley MATOP_GET_ROW_UTRIANGULAR=107, 1560d519adbfSMatthew Knepley MATOP_RESTORE_ROW_UTRIANGULAR=108, 1561d519adbfSMatthew Knepley MATOP_MATSOLVE=109, 1562d519adbfSMatthew Knepley MATOP_GET_REDUNDANTMATRIX=110, 1563d519adbfSMatthew Knepley MATOP_GET_ROW_MIN=111, 1564d519adbfSMatthew Knepley MATOP_GET_COLUMN_VEC=112, 1565d519adbfSMatthew Knepley MATOP_MISSING_DIAGONAL=113, 1566d519adbfSMatthew Knepley MATOP_MATGETSEQNONZEROSTRUCTURE=114, 156789c6957cSBarry Smith MATOP_CREATE=115, 156889c6957cSBarry Smith MATOP_GET_GHOSTS=116, 1569eeffb40dSHong Zhang MATOP_MULT_DIAGONAL_BLOCK=119, 1570547795f9SHong Zhang MATOP_HERMITIANTRANSPOSE=120, 1571547795f9SHong Zhang MATOP_MULTHERMITIANTRANSPOSE=121, 1572fbdbba38SShri Abhyankar MATOP_MULTHERMITIANTRANSPOSEADD=122, 1573503cfe85SJed Brown MATOP_GETMULTIPROCBLOCK=123 1574fae171e0SBarry Smith } MatOperation; 1575ace3abfcSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatHasOperation(Mat,MatOperation,PetscBool *); 1576be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellSetOperation(Mat,MatOperation,void(*)(void)); 1577be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellGetOperation(Mat,MatOperation,void(**)(void)); 1578be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellSetContext(Mat,void*); 1579112a2221SBarry Smith 158090ace30eSBarry Smith /* 158190ace30eSBarry Smith Codes for matrices stored on disk. By default they are 158290ace30eSBarry Smith stored in a universal format. By changing the format with 15837973ad04SBarry Smith PetscViewerSetFormat(viewer,PETSC_VIEWER_NATIVE); the matrices will 158490ace30eSBarry Smith be stored in a way natural for the matrix, for example dense matrices 158590ace30eSBarry Smith would be stored as dense. Matrices stored this way may only be 1586f4403165SShri Abhyankar read into matrices of the same type. 158790ace30eSBarry Smith */ 158890ace30eSBarry Smith #define MATRIX_BINARY_FORMAT_DENSE -1 158990ace30eSBarry Smith 1590be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetHashTableFactor(Mat,PetscReal); 1591be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatISGetLocalMat(Mat,Mat*); 15921f4e1ec7SBarry Smith 1593d9274352SBarry Smith /*S 1594d9274352SBarry Smith MatNullSpace - Object that removes a null space from a vector, i.e. 1595d9274352SBarry Smith orthogonalizes the vector to a subsapce 1596d9274352SBarry Smith 1597f7a9e4ceSBarry Smith Level: advanced 1598d9274352SBarry Smith 1599d9274352SBarry Smith Concepts: matrix; linear operator, null space 1600d9274352SBarry Smith 16016e1639daSBarry Smith Users manual sections: 16026e1639daSBarry Smith . sec_singular 16036e1639daSBarry Smith 1604d9274352SBarry Smith .seealso: MatNullSpaceCreate() 1605d9274352SBarry Smith S*/ 160674637425SBarry Smith typedef struct _p_MatNullSpace* MatNullSpace; 1607d9274352SBarry Smith 1608ace3abfcSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceCreate(MPI_Comm,PetscBool ,PetscInt,const Vec[],MatNullSpace*); 1609b22b330cSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(MatNullSpace,Vec,void*),void*); 1610be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceDestroy(MatNullSpace); 1611be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceRemove(MatNullSpace,Vec,Vec*); 1612be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceAttach(Mat,MatNullSpace); 1613ace3abfcSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceTest(MatNullSpace,Mat,PetscBool *); 161474637425SBarry Smith 1615be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatReorderingSeqSBAIJ(Mat,IS); 1616be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPISBAIJSetHashTableFactor(Mat,PetscReal); 1617be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetColumnIndices(Mat,PetscInt *); 161846129b97SKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJInvertBlockDiagonal(Mat); 16193f1d51d7SBarry Smith 1620be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMAIJ(Mat,PetscInt,Mat*); 1621be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMAIJRedimension(Mat,PetscInt,Mat*); 1622be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMAIJGetAIJ(Mat,Mat*); 1623c4f061fbSSatish Balay 1624be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatComputeExplicitOperator(Mat,Mat*); 1625b0a32e0cSBarry Smith 1626be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalScaleLocal(Mat,Vec); 162704f1ad80SBarry Smith 1628fef1beadSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*); 16293ec795f1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetBase(Mat,Vec,Vec); 16303ec795f1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*); 1631e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*)); 1632e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec)); 1633e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDAddNullSpace(Mat,MatNullSpace); 1634e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt); 1635e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDResetHHistory(Mat); 1636e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctionError(Mat,PetscReal); 1637e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetPeriod(Mat,PetscInt); 1638e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDGetH(Mat,PetscScalar *); 16396aa9148fSLisandro Dalcin EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetOptionsPrefix(Mat,const char[]); 1640e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFromOptions(Mat); 1641e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*); 1642e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*); 1643e884886fSBarry Smith 16446370053bSBarry Smith /*S 16456370053bSBarry Smith MatMFFD - A data structured used to manage the computation of the h differencing parameter for matrix-free 16466370053bSBarry Smith Jacobian vector products 1647e884886fSBarry Smith 16486370053bSBarry Smith Notes: MATMFFD is a specific MatType which uses the MatMFFD data structure 16496370053bSBarry Smith 16506370053bSBarry Smith MatMFFD*() methods actually take the Mat as their first argument. Not a MatMFFD data structure 16516370053bSBarry Smith 16526370053bSBarry Smith Level: developer 16536370053bSBarry Smith 16546370053bSBarry Smith .seealso: MATMFFD, MatCreateMFFD(), MatMFFDSetFuction(), MatMFFDSetType(), MatMFFDRegister() 16556370053bSBarry Smith S*/ 1656e884886fSBarry Smith typedef struct _p_MatMFFD* MatMFFD; 1657e884886fSBarry Smith 1658e884886fSBarry Smith /*E 1659e884886fSBarry Smith MatMFFDType - algorithm used to compute the h used in computing matrix-vector products via differencing of the function 1660e884886fSBarry Smith 1661e884886fSBarry Smith Level: beginner 1662e884886fSBarry Smith 1663e884886fSBarry Smith .seealso: MatMFFDSetType(), MatMFFDRegister() 1664e884886fSBarry Smith E*/ 1665a313700dSBarry Smith #define MatMFFDType char* 1666e884886fSBarry Smith #define MATMFFD_DS "ds" 1667e884886fSBarry Smith #define MATMFFD_WP "wp" 1668e884886fSBarry Smith 1669a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetType(Mat,const MatMFFDType); 1670e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatMFFD)); 1671e884886fSBarry Smith 1672e884886fSBarry Smith /*MC 1673e884886fSBarry Smith MatMFFDRegisterDynamic - Adds a method to the MatMFFD registry. 1674e884886fSBarry Smith 1675e884886fSBarry Smith Synopsis: 16761890ba74SBarry Smith PetscErrorCode MatMFFDRegisterDynamic(const char *name_solver,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatMFFD)) 1677e884886fSBarry Smith 1678e884886fSBarry Smith Not Collective 1679e884886fSBarry Smith 1680e884886fSBarry Smith Input Parameters: 1681e884886fSBarry Smith + name_solver - name of a new user-defined compute-h module 1682e884886fSBarry Smith . path - path (either absolute or relative) the library containing this solver 1683e884886fSBarry Smith . name_create - name of routine to create method context 1684e884886fSBarry Smith - routine_create - routine to create method context 1685e884886fSBarry Smith 1686e884886fSBarry Smith Level: developer 1687e884886fSBarry Smith 1688e884886fSBarry Smith Notes: 1689e884886fSBarry Smith MatMFFDRegisterDynamic() may be called multiple times to add several user-defined solvers. 1690e884886fSBarry Smith 1691e884886fSBarry Smith If dynamic libraries are used, then the fourth input argument (routine_create) 1692e884886fSBarry Smith is ignored. 1693e884886fSBarry Smith 1694e884886fSBarry Smith Sample usage: 1695e884886fSBarry Smith .vb 1696e884886fSBarry Smith MatMFFDRegisterDynamic("my_h",/home/username/my_lib/lib/libO/solaris/mylib.a, 1697e884886fSBarry Smith "MyHCreate",MyHCreate); 1698e884886fSBarry Smith .ve 1699e884886fSBarry Smith 1700e884886fSBarry Smith Then, your solver can be chosen with the procedural interface via 1701e884886fSBarry Smith $ MatMFFDSetType(mfctx,"my_h") 1702e884886fSBarry Smith or at runtime via the option 1703e884886fSBarry Smith $ -snes_mf_type my_h 1704e884886fSBarry Smith 1705e884886fSBarry Smith .keywords: MatMFFD, register 1706e884886fSBarry Smith 1707e884886fSBarry Smith .seealso: MatMFFDRegisterAll(), MatMFFDRegisterDestroy() 1708e884886fSBarry Smith M*/ 1709e884886fSBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 1710e884886fSBarry Smith #define MatMFFDRegisterDynamic(a,b,c,d) MatMFFDRegister(a,b,c,0) 1711e884886fSBarry Smith #else 1712e884886fSBarry Smith #define MatMFFDRegisterDynamic(a,b,c,d) MatMFFDRegister(a,b,c,d) 1713e884886fSBarry Smith #endif 1714e884886fSBarry Smith 1715e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDRegisterAll(const char[]); 1716e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDRegisterDestroy(void); 17171d0fab5eSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDDSSetUmin(Mat,PetscReal); 1718ace3abfcSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDWPSetComputeNormU(Mat,PetscBool ); 1719e884886fSBarry Smith 1720e884886fSBarry Smith 1721be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *); 1722be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *); 17237dbadf16SMatthew Knepley 172497969023SHong Zhang /* 172597969023SHong Zhang PETSc interface to MUMPS 172697969023SHong Zhang */ 172797969023SHong Zhang #ifdef PETSC_HAVE_MUMPS 1728a1f19f5aSHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetMumpsIcntl(Mat,PetscInt,PetscInt); 172997969023SHong Zhang #endif 173023a5497aSJed Brown 173123a5497aSJed Brown PETSC_EXTERN_CXX_END 173223a5497aSJed Brown #endif 1733