12eac72dbSBarry Smith /* 22eac72dbSBarry Smith Include file for the matrix component of PETSc 32eac72dbSBarry Smith */ 40a835dfdSSatish Balay #ifndef __PETSCMAT_H 50a835dfdSSatish Balay #define __PETSCMAT_H 60a835dfdSSatish Balay #include "petscvec.h" 7e9fa29b7SSatish Balay PETSC_EXTERN_CXX_BEGIN 82eac72dbSBarry Smith 9d9274352SBarry Smith /*S 10d9274352SBarry Smith Mat - Abstract PETSc matrix object 112eac72dbSBarry Smith 12d91e6319SBarry Smith Level: beginner 13d91e6319SBarry Smith 14d9274352SBarry Smith Concepts: matrix; linear operator 15d9274352SBarry Smith 16d9274352SBarry Smith .seealso: MatCreate(), MatType, MatSetType() 17d9274352SBarry Smith S*/ 18d9274352SBarry Smith typedef struct _p_Mat* Mat; 19d9274352SBarry Smith 20d9274352SBarry Smith /*E 21d9274352SBarry Smith MatType - String with the name of a PETSc matrix or the creation function 22d9274352SBarry Smith with an optional dynamic library name, for example 23d9274352SBarry Smith http://www.mcs.anl.gov/petsc/lib.a:mymatcreate() 24d9274352SBarry Smith 25d9274352SBarry Smith Level: beginner 26d9274352SBarry Smith 27c7393fdbSBarry Smith .seealso: MatSetType(), Mat, MatSolverPackage 28d91e6319SBarry Smith E*/ 29a313700dSBarry Smith #define MatType char* 30273d9f13SBarry Smith #define MATSAME "same" 31273d9f13SBarry Smith #define MATSEQMAIJ "seqmaij" 32273d9f13SBarry Smith #define MATMPIMAIJ "mpimaij" 33209238afSKris Buschelman #define MATMAIJ "maij" 34273d9f13SBarry Smith #define MATIS "is" 35273d9f13SBarry Smith #define MATSEQAIJ "seqaij" 36273d9f13SBarry Smith #define MATMPIAIJ "mpiaij" 37209238afSKris Buschelman #define MATAIJ "aij" 38273d9f13SBarry Smith #define MATSHELL "shell" 39209238afSKris Buschelman #define MATSEQDENSE "seqdense" 40273d9f13SBarry Smith #define MATMPIDENSE "mpidense" 41209238afSKris Buschelman #define MATDENSE "dense" 42273d9f13SBarry Smith #define MATSEQBAIJ "seqbaij" 43273d9f13SBarry Smith #define MATMPIBAIJ "mpibaij" 44209238afSKris Buschelman #define MATBAIJ "baij" 45273d9f13SBarry Smith #define MATMPIADJ "mpiadj" 46273d9f13SBarry Smith #define MATSEQSBAIJ "seqsbaij" 47273d9f13SBarry Smith #define MATMPISBAIJ "mpisbaij" 48209238afSKris Buschelman #define MATSBAIJ "sbaij" 49cebc7f6cSBarry Smith #define MATDAAD "daad" 50cebc7f6cSBarry Smith #define MATMFFD "mffd" 51c8a8475eSBarry Smith #define MATNORMAL "normal" 52ab92ecdeSBarry Smith #define MATLRC "lrc" 534b8d99adSRichard Tran Mills #define MATSEQCSRPERM "seqcsrperm" 5418def556SRichard Tran Mills #define MATMPICSRPERM "mpicsrperm" 55814655b8SBarry Smith #define MATCSRPERM "csrperm" 5681824310SBarry Smith #define MATSEQCRL "seqcrl" 57c02b24c5SBarry Smith #define MATMPICRL "mpicrl" 58c02b24c5SBarry Smith #define MATCRL "crl" 592a6744ebSBarry Smith #define MATSCATTER "scatter" 60421e10b8SBarry Smith #define MATBLOCKMAT "blockmat" 61793850ffSBarry Smith #define MATCOMPOSITE "composite" 625a7f1df3SHong Zhang #define MATSEQFFTW "seqfftw" 63557cca28SSatish Balay #define MATTRANSPOSEMAT "transpose" 6472ca8751SBarry Smith #define MATSCHURCOMPLEMENT "schurcomplement" 651d6018f0SLisandro Dalcin #define MATPYTHON "python" 66f91d8e95SBarry Smith #define MATHYPRESTRUCT "hyprestruct" 67a9e6138eSGlenn Hammond #define MATHYPRESSTRUCT "hypresstruct" 68ee1cef2cSJed Brown #define MATSUBMATRIX "submatrix" 69*c6570e9aSVictor Minden #define MATSEQAIJCUDA "seqaijcuda" 70d91e6319SBarry Smith 71164db6ccSDmitry Karpeev 72dbc6227fSDmitry Karpeev #define MATDD "matdd" 73164db6ccSDmitry Karpeev 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*); 117db4efbfdSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetFactorAvailable(Mat,const MatSolverPackage,MatFactorType,PetscTruth*); 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 214273d9f13SBarry Smith extern PetscTruth 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)) 281be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAdj(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscInt[],Mat*); 282be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 283ba966639SSatish 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) 284ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A) 285ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A) 286ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A) 287ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A)) 288ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A)) 289ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A)) 290be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPISBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 291ba966639SSatish 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) 292ba966639SSatish 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) 293ba966639SSatish 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) 294ba966639SSatish 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) 295ba966639SSatish 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)) 296ba966639SSatish 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)) 297ba966639SSatish 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)) 298ba966639SSatish 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) 299ba966639SSatish 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) 300ba966639SSatish 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) 301ba966639SSatish 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) 302ba966639SSatish 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)) 303ba966639SSatish 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)) 304ba966639SSatish 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)) 305be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateShell(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,void *,Mat*); 306ba966639SSatish 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) 307ba966639SSatish Balay PetscPolymorphicFunction(MatCreateShell,(PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx),(PETSC_COMM_WORLD,m,n,M,N,ctx,&A),Mat,A) 308be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateAdic(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,void (*)(void),Mat*); 309be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateNormal(Mat,Mat*); 310ba966639SSatish Balay PetscPolymorphicFunction(MatCreateNormal,(Mat mat),(mat,&A),Mat,A) 311ab92ecdeSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateLRC(Mat,Mat,Mat,Mat*); 312284134d9SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateIS(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,ISLocalToGlobalMapping,Mat*); 3131472f72bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 3141472f72bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPICRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 3152a6744ebSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateScatter(MPI_Comm,VecScatter,Mat*); 3162a6744ebSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScatterSetVecScatter(Mat,VecScatter); 3172a6744ebSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScatterGetVecScatter(Mat,VecScatter*); 3188cdf2d9bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateBlockMat(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,Mat*); 319793850ffSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompositeAddMat(Mat,Mat); 320b52f573bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompositeMerge(Mat); 321793850ffSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateComposite(MPI_Comm,PetscInt,const Mat*,Mat*); 3226d7c1e57SBarry Smith typedef enum {MAT_COMPOSITE_ADDITIVE,MAT_COMPOSITE_MULTIPLICATIVE} MatCompositeType; 3236d7c1e57SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompositeSetType(Mat,MatCompositeType); 3246d7c1e57SBarry Smith 325dbc6227fSDmitry Karpeev typedef enum {MATDD_BLOCK_COMM_DEFAULT = 0, MATDD_BLOCK_COMM_SELF = -1, MATDD_BLOCK_COMM_DETERMINE = -2} MatDDBlockCommType; 326dbc6227fSDmitry Karpeev EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDDAIJSetPreallocation(Mat A,PetscInt nz,PetscInt *nnz); 327dbc6227fSDmitry Karpeev EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDDGetDefaltBlockType(Mat A, const MatType *type); 328dbc6227fSDmitry Karpeev EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDDSetDefaltBlockType(Mat A, const MatType type); 329dbc6227fSDmitry Karpeev EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDDSetPreallocation_AIJ(Mat A,PetscInt nz,PetscInt *nnz); 330dbc6227fSDmitry Karpeev EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDDAddBlock(Mat A, PetscInt rowblock, PetscInt colblock, const MatType blockmattype, MatDDBlockCommType blockcommtype, Mat *block); 331dbc6227fSDmitry Karpeev EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDDSetBlock(Mat A, PetscInt rowblock, PetscInt colblock, Mat block); 332dbc6227fSDmitry Karpeev EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDDGetBlock(Mat A, PetscInt rowblock, PetscInt colblock, Mat *block); 333dbc6227fSDmitry Karpeev EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDDSetScatters(Mat A, PetscInt blockcount, Mat scatters[]); 334dbc6227fSDmitry Karpeev EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDDSetGathers(Mat A, PetscInt blockcount, Mat gathers[]); 33562ee7b53SDmitry Karpeev 336ba2f8784SDmitry Karpeev #if defined PETSC_HAVE_MATIM 337ba2f8784SDmitry Karpeev EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIMSetIS(Mat A, IS in, IS out); 338ba2f8784SDmitry Karpeev EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIMGetIS(Mat A, IS *in, IS *out); 339ba2f8784SDmitry Karpeev #endif 340ba2f8784SDmitry Karpeev 341ba2f8784SDmitry Karpeev 342a6a5cd3fSDmitry Karpeev 3435a7f1df3SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqFFTW(MPI_Comm,PetscInt,const PetscInt[],Mat*); 344f20085c4SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateTranspose(Mat,Mat*); 345ee1cef2cSJed Brown EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSubMatrix(Mat,IS,IS,Mat*); 346ee1cef2cSJed Brown EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSubMatrixUpdate(Mat,Mat,IS,IS); 3471472f72bSBarry Smith 3481d6018f0SLisandro Dalcin EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreatePython(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const char[],Mat*); 3491d6018f0SLisandro Dalcin EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPythonSetType(Mat,const char[]); 3501d6018f0SLisandro Dalcin 3511d6018f0SLisandro Dalcin 352f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUp(Mat); 353be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroy(Mat); 35421c89e3eSBarry Smith 3550c451bc4SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConjugate(Mat); 35699cafbc1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRealPart(Mat); 35799cafbc1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatImaginaryPart(Mat); 358bfc7d00aSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetDiagonalBlock(Mat,PetscTruth*,MatReuse,Mat*); 359fe97e370SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetTrace(Mat,PetscScalar*); 36099cafbc1SBarry Smith 3618ed539a5SBarry Smith /* ------------------------------------------------------------*/ 362be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 363be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 36487d4246cSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesRow(Mat,PetscInt,const PetscScalar[]); 365f259bd47SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesRowLocal(Mat,PetscInt,const PetscScalar[]); 36684cb2905SBarry Smith 3672ef4de8bSBarry Smith /*S 3682ef4de8bSBarry Smith MatStencil - Data structure (C struct) for storing information about a single row or 3692ef4de8bSBarry Smith column of a matrix as index on an associated grid. 3702ef4de8bSBarry Smith 3712ef4de8bSBarry Smith Level: beginner 3722ef4de8bSBarry Smith 3732ef4de8bSBarry Smith Concepts: matrix; linear operator 3742ef4de8bSBarry Smith 375d7bd93baSBarry Smith .seealso: MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockStencil() 3762ef4de8bSBarry Smith S*/ 377435da068SBarry Smith typedef struct { 378c1ac3661SBarry Smith PetscInt k,j,i,c; 379435da068SBarry Smith } MatStencil; 3802ef4de8bSBarry Smith 381be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode); 382be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode); 383be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt); 384435da068SBarry Smith 385be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetColoring(Mat,ISColoring); 386be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesAdic(Mat,void*); 387be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesAdifor(Mat,PetscInt,void*); 3883a7fca6bSBarry Smith 389d91e6319SBarry Smith /*E 390d91e6319SBarry Smith MatAssemblyType - Indicates if the matrix is now to be used, or if you plan 391d91e6319SBarry Smith to continue to add values to it 392d91e6319SBarry Smith 393d91e6319SBarry Smith Level: beginner 394d91e6319SBarry Smith 395d91e6319SBarry Smith .seealso: MatAssemblyBegin(), MatAssemblyEnd() 396d91e6319SBarry Smith E*/ 3976d4a8577SBarry Smith typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType; 398be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssemblyBegin(Mat,MatAssemblyType); 399be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssemblyEnd(Mat,MatAssemblyType); 400be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssembled(Mat,PetscTruth*); 4014f9c727eSBarry Smith 4021d73ed98SBarry Smith 40330de9b25SBarry Smith 404d91e6319SBarry Smith /*E 405d91e6319SBarry Smith MatOption - Options that may be set for a matrix and its behavior or storage 406d91e6319SBarry Smith 407d91e6319SBarry Smith Level: beginner 408d91e6319SBarry Smith 4090a835dfdSSatish Balay Any additions/changes here MUST also be made in include/finclude/petscmat.h 410d91e6319SBarry Smith 411d91e6319SBarry Smith .seealso: MatSetOption() 412d91e6319SBarry Smith E*/ 413512a5fc5SBarry Smith typedef enum {MAT_ROW_ORIENTED,MAT_NEW_NONZERO_LOCATIONS, 4144e0d8c25SBarry Smith MAT_SYMMETRIC, 4154e0d8c25SBarry Smith MAT_STRUCTURALLY_SYMMETRIC, 4168047e77bSBarry Smith MAT_NEW_DIAGONALS,MAT_IGNORE_OFF_PROC_ENTRIES, 4174e0d8c25SBarry Smith MAT_NEW_NONZERO_LOCATION_ERR, 4184e0d8c25SBarry Smith MAT_NEW_NONZERO_ALLOCATION_ERR,MAT_USE_HASH_TABLE, 419a9817697SBarry Smith MAT_KEEP_NONZERO_PATTERN,MAT_IGNORE_ZERO_ENTRIES, 4204e0d8c25SBarry Smith MAT_USE_INODES, 4214e0d8c25SBarry Smith MAT_HERMITIAN, 4224e0d8c25SBarry Smith MAT_SYMMETRY_ETERNAL, 4234e0d8c25SBarry Smith MAT_USE_COMPRESSEDROW, 4244e0d8c25SBarry Smith MAT_IGNORE_LOWER_TRIANGULAR,MAT_ERROR_LOWER_TRIANGULAR, 42528b2fa4aSMatthew Knepley MAT_GETROW_UPPERTRIANGULAR,MAT_UNUSED_NONZERO_LOCATION_ERR, 42628b2fa4aSMatthew Knepley NUM_MAT_OPTIONS} MatOption; 427290bbb0aSBarry Smith extern const char *MatOptions[]; 4284e0d8c25SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetOption(Mat,MatOption,PetscTruth); 429a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetType(Mat,const MatType*); 430a313700dSBarry Smith PetscPolymorphicFunction(MatGetType,(Mat mat),(mat,&t),const MatType,t) 43184cb2905SBarry Smith 432be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]); 433be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 434be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 435f5edf698SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowUpperTriangular(Mat); 436f5edf698SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRowUpperTriangular(Mat); 437be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 438be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 439be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumnVector(Mat,Vec,PetscInt); 440be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetArray(Mat,PetscScalar *[]); 441ba966639SSatish Balay PetscPolymorphicFunction(MatGetArray,(Mat mat),(mat,&a),PetscScalar*,a) 442be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreArray(Mat,PetscScalar *[]); 443be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBlockSize(Mat,PetscInt *); 444ba966639SSatish Balay PetscPolymorphicFunction(MatGetBlockSize,(Mat mat),(mat,&a),PetscInt,a) 445be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetBlockSize(Mat,PetscInt); 4467b80b807SBarry Smith 4471620fd73SBarry Smith 448be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMult(Mat,Vec,Vec); 44989c6957cSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultDiagonalBlock(Mat,Vec,Vec); 450be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultAdd(Mat,Vec,Vec,Vec); 451ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultAdd,(Mat A,Vec x,Vec y),(A,x,y,y)) 452be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTranspose(Mat,Vec,Vec); 453547795f9SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultHermitianTranspose(Mat,Vec,Vec); 454be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsTranspose(Mat,Mat,PetscReal,PetscTruth*); 455ba966639SSatish Balay PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B,PetscReal tol),(A,B,tol,&t),PetscTruth,t) 456e37f64cdSBarry Smith PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B),(A,B,0,&t),PetscTruth,t) 4571cbb95d3SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitianTranspose(Mat,Mat,PetscReal,PetscTruth*); 458be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeAdd(Mat,Vec,Vec,Vec); 4595cb05578SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultHermitianTransposeAdd(Mat,Vec,Vec,Vec); 460ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultTransposeAdd,(Mat A,Vec x,Vec y),(A,x,y,y)) 461be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultConstrained(Mat,Vec,Vec); 462be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeConstrained(Mat,Vec,Vec); 4632bebee5dSHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatSolve(Mat,Mat,Mat); 464f5edf698SHong Zhang 465d91e6319SBarry Smith /*E 466d91e6319SBarry Smith MatDuplicateOption - Indicates if a duplicated sparse matrix should have 467d91e6319SBarry Smith its numerical values copied over or just its nonzero structure. 468d91e6319SBarry Smith 469d91e6319SBarry Smith Level: beginner 470d91e6319SBarry Smith 471d91e6319SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 472d91e6319SBarry Smith 47370dcbbb9SBarry Smith $ MAT_SHARE_NONZERO_PATTERN - the i and j arrays in the new matrix will be shared with the original matrix 47470dcbbb9SBarry Smith $ this also triggers the MAT_DO_NOT_COPY_VALUES option. This is used when you 47570dcbbb9SBarry Smith $ have several matrices with the same nonzero pattern. 47670dcbbb9SBarry Smith 477d91e6319SBarry Smith .seealso: MatDuplicate() 478d91e6319SBarry Smith E*/ 47970dcbbb9SBarry Smith typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES,MAT_SHARE_NONZERO_PATTERN} MatDuplicateOption; 4802e8a6d31SBarry Smith 481a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert(Mat,const MatType,MatReuse,Mat*); 482a313700dSBarry Smith PetscPolymorphicFunction(MatConvert,(Mat A,const MatType t),(A,t,MAT_INITIAL_MATRIX,&a),Mat,a) 483be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDuplicate(Mat,MatDuplicateOption,Mat*); 484ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A,MatDuplicateOption o),(A,o,&a),Mat,a) 485ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A),(A,MAT_COPY_VALUES,&a),Mat,a) 48694a9d846SBarry Smith 487cb5b572fSBarry Smith 488be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCopy(Mat,Mat,MatStructure); 489be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatView(Mat,PetscViewer); 490be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsSymmetric(Mat,PetscReal,PetscTruth*); 491ba966639SSatish Balay PetscPolymorphicFunction(MatIsSymmetric,(Mat A,PetscReal tol),(A,tol,&t),PetscTruth,t) 492e37f64cdSBarry Smith PetscPolymorphicFunction(MatIsSymmetric,(Mat A),(A,0,&t),PetscTruth,t) 493be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsStructurallySymmetric(Mat,PetscTruth*); 494ba966639SSatish Balay PetscPolymorphicFunction(MatIsStructurallySymmetric,(Mat A),(A,&t),PetscTruth,t) 4951cbb95d3SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitian(Mat,PetscReal,PetscTruth*); 4961cbb95d3SBarry Smith PetscPolymorphicFunction(MatIsHermitian,(Mat A),(A,0,&t),PetscTruth,t) 497be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsSymmetricKnown(Mat,PetscTruth*,PetscTruth*); 498be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitianKnown(Mat,PetscTruth*,PetscTruth*); 49964c62002SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMissingDiagonal(Mat,PetscTruth *,PetscInt *); 500a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLoad(PetscViewer,const MatType,Mat*); 501fbdbba38SShri Abhyankar EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLoadnew(PetscViewer,Mat); 502a313700dSBarry Smith PetscPolymorphicFunction(MatLoad,(PetscViewer v,const MatType t),(v,t,&a),Mat,a) 5037b80b807SBarry Smith 5048f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *); 5058f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRowIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *); 5068f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumnIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *); 5078f7157efSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreColumnIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *); 508d4fbbf0eSBarry Smith 509d91e6319SBarry Smith /*S 510d91e6319SBarry Smith MatInfo - Context of matrix information, used with MatGetInfo() 511d91e6319SBarry Smith 512d91e6319SBarry Smith In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE 513d91e6319SBarry Smith 514d91e6319SBarry Smith Level: intermediate 515d91e6319SBarry Smith 516d91e6319SBarry Smith Concepts: matrix^nonzero information 517d91e6319SBarry Smith 518d9274352SBarry Smith .seealso: MatGetInfo(), MatInfoType 519d91e6319SBarry Smith S*/ 5204e220ebcSLois Curfman McInnes typedef struct { 521b0a32e0cSBarry Smith PetscLogDouble block_size; /* block size */ 522b0a32e0cSBarry Smith PetscLogDouble nz_allocated,nz_used,nz_unneeded; /* number of nonzeros */ 523b0a32e0cSBarry Smith PetscLogDouble memory; /* memory allocated */ 524b0a32e0cSBarry Smith PetscLogDouble assemblies; /* number of matrix assemblies called */ 525b0a32e0cSBarry Smith PetscLogDouble mallocs; /* number of mallocs during MatSetValues() */ 526b0a32e0cSBarry Smith PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */ 527b0a32e0cSBarry Smith PetscLogDouble factor_mallocs; /* number of mallocs during factorization */ 5284e220ebcSLois Curfman McInnes } MatInfo; 5294e220ebcSLois Curfman McInnes 530d9274352SBarry Smith /*E 531d9274352SBarry Smith MatInfoType - Indicates if you want information about the local part of the matrix, 532d9274352SBarry Smith the entire parallel matrix or the maximum over all the local parts. 533d9274352SBarry Smith 534d9274352SBarry Smith Level: beginner 535d9274352SBarry Smith 536d9274352SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 537d9274352SBarry Smith 538d9274352SBarry Smith .seealso: MatGetInfo(), MatInfo 539d9274352SBarry Smith E*/ 5407b80b807SBarry Smith typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType; 541be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetInfo(Mat,MatInfoType,MatInfo*); 542be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatValid(Mat,PetscTruth*); 543be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetDiagonal(Mat,Vec); 544985db425SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMax(Mat,Vec,PetscInt[]); 545985db425SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMin(Mat,Vec,PetscInt[]); 546985db425SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMaxAbs(Mat,Vec,PetscInt[]); 547c87e5d42SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMinAbs(Mat,Vec,PetscInt[]); 54886697f06SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowSum(Mat,Vec); 549fc4dec0aSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatTranspose(Mat,MatReuse,Mat*); 550fc4dec0aSBarry Smith PetscPolymorphicFunction(MatTranspose,(Mat A),(A,MAT_INITIAL_MATRIX,&t),Mat,t) 551eeffb40dSHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatHermitianTranspose(Mat,MatReuse,Mat*); 552be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPermute(Mat,IS,IS,Mat *); 553ba966639SSatish Balay PetscPolymorphicFunction(MatPermute,(Mat A,IS is1,IS is2),(A,is1,is2,&t),Mat,t) 554be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPermuteSparsify(Mat,PetscInt,PetscReal,PetscReal,IS,IS,Mat *); 555be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalScale(Mat,Vec,Vec); 556be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalSet(Mat,Vec,InsertMode); 557be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatEqual(Mat,Mat,PetscTruth*); 558ba966639SSatish Balay PetscPolymorphicFunction(MatEqual,(Mat A,Mat B),(A,B,&t),PetscTruth,t) 559be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultEqual(Mat,Mat,PetscInt,PetscTruth*); 560be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultAddEqual(Mat,Mat,PetscInt,PetscTruth*); 561be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeEqual(Mat,Mat,PetscInt,PetscTruth*); 562be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscTruth*); 5637b80b807SBarry Smith 564be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNorm(Mat,NormType,PetscReal *); 565ba966639SSatish Balay PetscPolymorphicFunction(MatNorm,(Mat A,NormType t),(A,t,&n),PetscReal,n) 566242f1d38SBarry Smith EXTERN PetscErrorCode MatGetColumnNorms(Mat,NormType,PetscReal *); 567be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroEntries(Mat); 568f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar); 569f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsIS(Mat,IS,PetscScalar); 570f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroColumns(Mat,PetscInt,const PetscInt [],const PetscScalar*); 571f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroColumnsIS(Mat,IS,const PetscScalar*); 5727b80b807SBarry Smith 573be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatUseScaledForm(Mat,PetscTruth); 574be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScaleSystem(Mat,Vec,Vec); 575be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatUnScaleSystem(Mat,Vec,Vec); 5765ef9f2a5SBarry Smith 577be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSize(Mat,PetscInt*,PetscInt*); 578be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalSize(Mat,PetscInt*,PetscInt*); 579be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRange(Mat,PetscInt*,PetscInt*); 5808ee2e534SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRanges(Mat,const PetscInt**); 581f4b8409dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRangeColumn(Mat,PetscInt*,PetscInt*); 58272ca8751SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRangesColumn(Mat,const PetscInt**); 5837b80b807SBarry Smith 584be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]); 585be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroyMatrices(PetscInt,Mat *[]); 5864aa3045dSJed Brown EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrix(Mat,IS,IS,MatReuse,Mat *); 587f6d58c54SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSeqNonzeroStructure(Mat,Mat*); 588f6d58c54SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroySeqNonzeroStructure(Mat*); 5895494a064SHong Zhang 590be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge(MPI_Comm,Mat,PetscInt,MatReuse,Mat*); 591be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPI(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*); 592be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPISymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*); 593be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPINumeric(Mat,Mat); 594be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroy_MPIAIJ_SeqsToMPI(Mat); 595be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalMat(Mat,MatReuse,Mat*); 596be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*); 597be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,PetscInt*,Mat*); 5981d79065fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBrowsOfAoCols(Mat,Mat,MatReuse,PetscInt**,PetscInt**,MatScalar**,Mat*); 59943eb5e2fSMatthew Knepley #if defined (PETSC_USE_CTABLE) 600cd116e26SSatish Balay #include "petscctable.h" 60143eb5e2fSMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetCommunicationStructs(Mat, Vec *, PetscTable *, VecScatter *); 60243eb5e2fSMatthew Knepley #else 60343eb5e2fSMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetCommunicationStructs(Mat, Vec *, PetscInt *[], VecScatter *); 60443eb5e2fSMatthew Knepley #endif 6058c7482ecSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetGhosts(Mat, PetscInt *,const PetscInt *[]); 6068efafbd8SBarry Smith 607be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt); 6087b80b807SBarry Smith 609be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*); 610be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultSymbolic(Mat,Mat,PetscReal,Mat*); 611be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultNumeric(Mat,Mat,Mat); 61222440eb1SKris Buschelman 613be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*); 614be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAPSymbolic(Mat,Mat,PetscReal,Mat*); 615be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAPNumeric(Mat,Mat,Mat); 61622440eb1SKris Buschelman 617be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTranspose(Mat,Mat,MatReuse,PetscReal,Mat*); 618be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTransposeSymbolic(Mat,Mat,PetscReal,Mat*); 619be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTransposeNumeric(Mat,Mat,Mat); 620bc011b1eSHong Zhang 621f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAXPY(Mat,PetscScalar,Mat,MatStructure); 62204aac2b0SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAYPX(Mat,PetscScalar,Mat,MatStructure); 6231c741599SHong Zhang 624f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScale(Mat,PetscScalar); 625f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShift(Mat,PetscScalar); 6267b80b807SBarry Smith 627be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping); 628be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetLocalToGlobalMappingBlock(Mat,ISLocalToGlobalMapping); 629f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar); 630f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsLocalIS(Mat,IS,PetscScalar); 631be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 632be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 633052efed2SBarry Smith 634be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStashSetInitialSize(Mat,PetscInt,PetscInt); 635be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*); 63690f02eecSBarry Smith 637be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInterpolate(Mat,Vec,Vec); 6380c451bc4SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInterpolateAdd(Mat,Vec,Vec,Vec); 639ba966639SSatish Balay PetscPolymorphicSubroutine(MatInterpolateAdd,(Mat A,Vec x,Vec y),(A,x,y,y)) 640be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestrict(Mat,Vec,Vec); 641be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetVecs(Mat,Vec*,Vec*); 64269db28dcSHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRedundantMatrix(Mat,PetscInt,MPI_Comm,PetscInt,MatReuse,Mat*); 643bd481603SBarry Smith 644bd481603SBarry Smith /*MC 6451620fd73SBarry Smith MatSetValue - Set a single entry into a matrix. 6461620fd73SBarry Smith 6471620fd73SBarry Smith Not collective 6481620fd73SBarry Smith 6491620fd73SBarry Smith Input Parameters: 6501620fd73SBarry Smith + m - the matrix 6511620fd73SBarry Smith . row - the row location of the entry 6521620fd73SBarry Smith . col - the column location of the entry 6531620fd73SBarry Smith . value - the value to insert 6541620fd73SBarry Smith - mode - either INSERT_VALUES or ADD_VALUES 6551620fd73SBarry Smith 6561620fd73SBarry Smith Notes: 6571620fd73SBarry Smith For efficiency one should use MatSetValues() and set several or many 6581620fd73SBarry Smith values simultaneously if possible. 6591620fd73SBarry Smith 6601620fd73SBarry Smith Level: beginner 6611620fd73SBarry Smith 6621620fd73SBarry Smith .seealso: MatSetValues(), MatSetValueLocal() 6631620fd73SBarry Smith M*/ 6641620fd73SBarry 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);} 6651620fd73SBarry Smith 6661620fd73SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatGetValue(Mat v,PetscInt i,PetscInt j,PetscScalar *va) {return MatGetValues(v,1,&i,1,&j,va);} 6671620fd73SBarry Smith 6681620fd73SBarry 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);} 6691620fd73SBarry Smith 6701620fd73SBarry Smith /*MC 6710d2a7ff2SSatish Balay MatPreallocateInitialize - Begins the block of code that will count the number of nonzeros per 672bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 673bd481603SBarry Smith 674bd481603SBarry Smith Synopsis: 675c1ac3661SBarry Smith PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz) 676bd481603SBarry Smith 677bd481603SBarry Smith Collective on MPI_Comm 678bd481603SBarry Smith 679bd481603SBarry Smith Input Parameters: 680bd481603SBarry Smith + comm - the communicator that will share the eventually allocated matrix 681859fcb39SBarry Smith . nrows - the number of LOCAL rows in the matrix 682859fcb39SBarry Smith - ncols - the number of LOCAL columns in the matrix 683bd481603SBarry Smith 684bd481603SBarry Smith Output Parameters: 685bd481603SBarry Smith + dnz - the array that will be passed to the matrix preallocation routines 686bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 687bd481603SBarry Smith 688bd481603SBarry Smith 689bd481603SBarry Smith Level: intermediate 690bd481603SBarry Smith 691bd481603SBarry Smith Notes: 692bd481603SBarry Smith See the chapter in the users manual on performance for more details 693bd481603SBarry Smith 6941d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 695bd481603SBarry Smith 696bd481603SBarry Smith Use MatPreallocateInitializeSymmetric() for symmetric matrices (MPISBAIJ matrices) 697bd481603SBarry Smith 6981620fd73SBarry Smith This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize(). 6991620fd73SBarry Smith 700bd481603SBarry Smith Concepts: preallocation^Matrix 701bd481603SBarry Smith 702bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(), 703bd481603SBarry Smith MatPreallocateInitializeSymmetric(), MatPreallocateSymmetricSetLocal() 704bd481603SBarry Smith M*/ 705c4f061fbSSatish Balay #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \ 7067c922b88SBarry Smith { \ 707a89bc211SBarry Smith PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__start,__end; \ 708a89bc211SBarry Smith _4_ierr = PetscMalloc2(__nrows,PetscInt,&dnz,__nrows,PetscInt,&onz);CHKERRQ(_4_ierr); \ 709a89bc211SBarry Smith _4_ierr = PetscMemzero(dnz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\ 710a89bc211SBarry Smith _4_ierr = PetscMemzero(onz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\ 711c1ac3661SBarry Smith _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ctmp;\ 712a89bc211SBarry Smith _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows; 7137c922b88SBarry Smith 714bd481603SBarry Smith /*MC 7150d2a7ff2SSatish Balay MatPreallocateSymmetricInitialize - Begins the block of code that will count the number of nonzeros per 716bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 717bd481603SBarry Smith 718bd481603SBarry Smith Synopsis: 719c1ac3661SBarry Smith PetscErrorCode MatPreallocateSymmetricInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz) 720bd481603SBarry Smith 721bd481603SBarry Smith Collective on MPI_Comm 722bd481603SBarry Smith 723bd481603SBarry Smith Input Parameters: 724bd481603SBarry Smith + comm - the communicator that will share the eventually allocated matrix 725859fcb39SBarry Smith . nrows - the number of LOCAL rows in the matrix 726859fcb39SBarry Smith - ncols - the number of LOCAL columns in the matrix 727bd481603SBarry Smith 728bd481603SBarry Smith Output Parameters: 729bd481603SBarry Smith + dnz - the array that will be passed to the matrix preallocation routines 730bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 731bd481603SBarry Smith 732bd481603SBarry Smith 733bd481603SBarry Smith Level: intermediate 734bd481603SBarry Smith 735bd481603SBarry Smith Notes: 736bd481603SBarry Smith See the chapter in the users manual on performance for more details 737bd481603SBarry Smith 7381d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 739bd481603SBarry Smith 7401620fd73SBarry Smith This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize(). 7411620fd73SBarry Smith 742bd481603SBarry Smith Concepts: preallocation^Matrix 743bd481603SBarry Smith 744bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(), 745bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal() 746bd481603SBarry Smith M*/ 747222b16d4SBarry Smith #define MatPreallocateSymmetricInitialize(comm,nrows,ncols,dnz,onz) 0; \ 748222b16d4SBarry Smith { \ 749a89bc211SBarry Smith PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__end; \ 750a89bc211SBarry Smith _4_ierr = PetscMalloc2(__nrows,PetscInt,&dnz,__nrows,PetscInt,&onz);CHKERRQ(_4_ierr); \ 751a89bc211SBarry Smith _4_ierr = PetscMemzero(dnz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\ 752a89bc211SBarry Smith _4_ierr = PetscMemzero(onz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\ 753c1ac3661SBarry Smith _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr);\ 754a89bc211SBarry Smith _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows; 755222b16d4SBarry Smith 756bd481603SBarry Smith /*MC 757bd481603SBarry Smith MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be 758bd481603SBarry Smith inserted using a local number of the rows and columns 759bd481603SBarry Smith 760bd481603SBarry Smith Synopsis: 761c1ac3661SBarry Smith PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 762bd481603SBarry Smith 763bd481603SBarry Smith Not Collective 764bd481603SBarry Smith 765bd481603SBarry Smith Input Parameters: 766bd481603SBarry Smith + map - the mapping between local numbering and global numbering 767bd481603SBarry Smith . nrows - the number of rows indicated 7681d73ed98SBarry Smith . rows - the indices of the rows 769bd481603SBarry Smith . ncols - the number of columns in the matrix 770bd481603SBarry Smith . cols - the columns indicated 771bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 772bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 773bd481603SBarry Smith 774bd481603SBarry Smith 775bd481603SBarry Smith Level: intermediate 776bd481603SBarry Smith 777bd481603SBarry Smith Notes: 778bd481603SBarry Smith See the chapter in the users manual on performance for more details 779bd481603SBarry Smith 7801d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 781bd481603SBarry Smith 782bd481603SBarry Smith Concepts: preallocation^Matrix 783bd481603SBarry Smith 784bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(), 785bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal() 786bd481603SBarry Smith M*/ 787c4f061fbSSatish Balay #define MatPreallocateSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\ 788c4f061fbSSatish Balay {\ 789c1ac3661SBarry Smith PetscInt __l;\ 790ef66eb69SBarry Smith _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\ 791ef66eb69SBarry Smith _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\ 792c4f061fbSSatish Balay for (__l=0;__l<nrows;__l++) {\ 793ef66eb69SBarry Smith _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 794c4f061fbSSatish Balay }\ 795c4f061fbSSatish Balay } 796c4f061fbSSatish Balay 797bd481603SBarry Smith /*MC 798bd481603SBarry Smith MatPreallocateSymmetricSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be 799bd481603SBarry Smith inserted using a local number of the rows and columns 800bd481603SBarry Smith 801bd481603SBarry Smith Synopsis: 802c1ac3661SBarry Smith PetscErrorCode MatPreallocateSymmetricSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 803bd481603SBarry Smith 804bd481603SBarry Smith Not Collective 805bd481603SBarry Smith 806bd481603SBarry Smith Input Parameters: 807bd481603SBarry Smith + map - the mapping between local numbering and global numbering 808bd481603SBarry Smith . nrows - the number of rows indicated 8091d73ed98SBarry Smith . rows - the indices of the rows 810bd481603SBarry Smith . ncols - the number of columns in the matrix 811bd481603SBarry Smith . cols - the columns indicated 812bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 813bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 814bd481603SBarry Smith 815bd481603SBarry Smith 816bd481603SBarry Smith Level: intermediate 817bd481603SBarry Smith 818bd481603SBarry Smith Notes: 819bd481603SBarry Smith See the chapter in the users manual on performance for more details 820bd481603SBarry Smith 821bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 822bd481603SBarry Smith 823bd481603SBarry Smith Concepts: preallocation^Matrix 824bd481603SBarry Smith 825bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(), 826bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal() 827bd481603SBarry Smith M*/ 828d3d32019SBarry Smith #define MatPreallocateSymmetricSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\ 829d3d32019SBarry Smith {\ 830c1ac3661SBarry Smith PetscInt __l;\ 831d3d32019SBarry Smith _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\ 832d3d32019SBarry Smith _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\ 833d3d32019SBarry Smith for (__l=0;__l<nrows;__l++) {\ 834d3d32019SBarry Smith _4_ierr = MatPreallocateSymmetricSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 835d3d32019SBarry Smith }\ 836d3d32019SBarry Smith } 837d3d32019SBarry Smith 838bd481603SBarry Smith /*MC 839bd481603SBarry Smith MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be 840bd481603SBarry Smith inserted using a local number of the rows and columns 841bd481603SBarry Smith 842bd481603SBarry Smith Synopsis: 843c1ac3661SBarry Smith PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 844bd481603SBarry Smith 845bd481603SBarry Smith Not Collective 846bd481603SBarry Smith 847bd481603SBarry Smith Input Parameters: 84864075487SBarry Smith + row - the row 849bd481603SBarry Smith . ncols - the number of columns in the matrix 850a50f8bf6SHong Zhang - cols - the columns indicated 851a50f8bf6SHong Zhang 852a50f8bf6SHong Zhang Output Parameters: 853a50f8bf6SHong Zhang + dnz - the array that will be passed to the matrix preallocation routines 854bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 855bd481603SBarry Smith 856bd481603SBarry Smith 857bd481603SBarry Smith Level: intermediate 858bd481603SBarry Smith 859bd481603SBarry Smith Notes: 860bd481603SBarry Smith See the chapter in the users manual on performance for more details 861bd481603SBarry Smith 862bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 863bd481603SBarry Smith 8641620fd73SBarry Smith This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize(). 8651620fd73SBarry Smith 866bd481603SBarry Smith Concepts: preallocation^Matrix 867bd481603SBarry Smith 868bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(), 869bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal() 870bd481603SBarry Smith M*/ 871c4f061fbSSatish Balay #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\ 872c1ac3661SBarry Smith { PetscInt __i; \ 873e32f2f54SBarry 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);\ 874e32f2f54SBarry 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);\ 8757c922b88SBarry Smith for (__i=0; __i<nc; __i++) {\ 87664075487SBarry Smith if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \ 8777cc688d7SBarry Smith else dnz[row - __rstart]++;\ 8787c922b88SBarry Smith }\ 8797c922b88SBarry Smith } 8807c922b88SBarry Smith 881bd481603SBarry Smith /*MC 882bd481603SBarry Smith MatPreallocateSymmetricSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be 883bd481603SBarry Smith inserted using a local number of the rows and columns 884bd481603SBarry Smith 885bd481603SBarry Smith Synopsis: 886c1ac3661SBarry Smith PetscErrorCode MatPreallocateSymmetricSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 887bd481603SBarry Smith 888bd481603SBarry Smith Not Collective 889bd481603SBarry Smith 890bd481603SBarry Smith Input Parameters: 891bd481603SBarry Smith + nrows - the number of rows indicated 8921d73ed98SBarry Smith . rows - the indices of the rows 893bd481603SBarry Smith . ncols - the number of columns in the matrix 894bd481603SBarry Smith . cols - the columns indicated 895bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 896bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 897bd481603SBarry Smith 898bd481603SBarry Smith 899bd481603SBarry Smith Level: intermediate 900bd481603SBarry Smith 901bd481603SBarry Smith Notes: 902bd481603SBarry Smith See the chapter in the users manual on performance for more details 903bd481603SBarry Smith 904bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 905bd481603SBarry Smith 9061620fd73SBarry Smith This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize(). 9071620fd73SBarry Smith 908bd481603SBarry Smith Concepts: preallocation^Matrix 909bd481603SBarry Smith 910bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(), 911bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal() 912bd481603SBarry Smith M*/ 913d3d32019SBarry Smith #define MatPreallocateSymmetricSet(row,nc,cols,dnz,onz) 0;\ 914c1ac3661SBarry Smith { PetscInt __i; \ 915d3d32019SBarry Smith for (__i=0; __i<nc; __i++) {\ 916d3d32019SBarry Smith if (cols[__i] >= __end) onz[row - __rstart]++; \ 917d3d32019SBarry Smith else if (cols[__i] >= row) dnz[row - __rstart]++;\ 918d3d32019SBarry Smith }\ 919d3d32019SBarry Smith } 920d3d32019SBarry Smith 921bd481603SBarry Smith /*MC 92216371a99SBarry Smith MatPreallocateLocation - An alternative to MatPreallocationSet() that puts the nonzero locations into the matrix if it exists 92316371a99SBarry Smith 92416371a99SBarry Smith Synopsis: 92516371a99SBarry Smith PetscErrorCode MatPreallocateLocations(Mat A,PetscInt row,PetscInt ncols,PetscInt *cols,PetscInt *dnz,PetscInt *onz) 92616371a99SBarry Smith 92716371a99SBarry Smith Not Collective 92816371a99SBarry Smith 92916371a99SBarry Smith Input Parameters: 93016371a99SBarry Smith . A - matrix 93116371a99SBarry Smith . row - row where values exist (must be local to this process) 93216371a99SBarry Smith . ncols - number of columns 93316371a99SBarry Smith . cols - columns with nonzeros 93416371a99SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 93516371a99SBarry Smith - ozn - the other array passed to the matrix preallocation routines 93616371a99SBarry Smith 93716371a99SBarry Smith 93816371a99SBarry Smith Level: intermediate 93916371a99SBarry Smith 94016371a99SBarry Smith Notes: 94116371a99SBarry Smith See the chapter in the users manual on performance for more details 94216371a99SBarry Smith 94316371a99SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 94416371a99SBarry Smith 94516371a99SBarry Smith This is a MACRO not a function because it uses a bunch of variables private to the MatPreallocation.... routines. 94616371a99SBarry Smith 94716371a99SBarry Smith Concepts: preallocation^Matrix 94816371a99SBarry Smith 94916371a99SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(), 95016371a99SBarry Smith MatPreallocateSymmetricInitialize(), MatPreallocateSymmetricSetLocal() 95116371a99SBarry Smith M*/ 95216371a99SBarry 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);} 95316371a99SBarry Smith 95416371a99SBarry Smith 95516371a99SBarry Smith /*MC 9560d2a7ff2SSatish Balay MatPreallocateFinalize - Ends the block of code that will count the number of nonzeros per 957bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 958bd481603SBarry Smith 959bd481603SBarry Smith Synopsis: 960c1ac3661SBarry Smith PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz) 961bd481603SBarry Smith 962bd481603SBarry Smith Collective on MPI_Comm 963bd481603SBarry Smith 964bd481603SBarry Smith Input Parameters: 96516371a99SBarry Smith + dnz - the array that was be passed to the matrix preallocation routines 966bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 967bd481603SBarry Smith 968bd481603SBarry Smith 969bd481603SBarry Smith Level: intermediate 970bd481603SBarry Smith 971bd481603SBarry Smith Notes: 972bd481603SBarry Smith See the chapter in the users manual on performance for more details 973bd481603SBarry Smith 974bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 975bd481603SBarry Smith 9761620fd73SBarry Smith This is a MACRO not a function because it closes the { started in MatPreallocateInitialize(). 9771620fd73SBarry Smith 978bd481603SBarry Smith Concepts: preallocation^Matrix 979bd481603SBarry Smith 980bd481603SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(), 981bd481603SBarry Smith MatPreallocateSymmetricInitialize(), MatPreallocateSymmetricSetLocal() 982bd481603SBarry Smith M*/ 983a89bc211SBarry Smith #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree2(dnz,onz);CHKERRQ(_4_ierr);} 9847c922b88SBarry Smith 985bd481603SBarry Smith 986bd481603SBarry Smith 9877b80b807SBarry Smith /* Routines unique to particular data structures */ 988be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellGetContext(Mat,void **); 989ba966639SSatish Balay PetscPolymorphicFunction(MatShellGetContext,(Mat A),(A,&t),void*,t) 990698d4c6aSKris Buschelman 991be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInodeAdjustForInodes(Mat,IS*,IS*); 992be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *); 993698d4c6aSKris Buschelman 994be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetColumnIndices(Mat,PetscInt[]); 995be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetColumnIndices(Mat,PetscInt[]); 996be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 997c75a6043SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 998c75a6043SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 9997b80b807SBarry Smith 1000a96a251dSBarry Smith #define MAT_SKIP_ALLOCATION -4 1001a96a251dSBarry Smith 1002be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]); 1003ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz)) 1004be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]); 1005ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqSBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz)) 1006be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]); 1007ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqAIJSetPreallocation,(Mat A,const PetscInt nnz[]),(A,0,nnz)) 1008273d9f13SBarry Smith 1009be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 1010ba966639SSatish Balay PetscPolymorphicSubroutine(MatMPIBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[],const PetscInt onz[]),(A,bs,0,nnz,0,onz)) 1011be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 1012be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 101353803f26SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []); 1014725b52f3SLisandro Dalcin EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 1015be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]); 1016be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 1017be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]); 1018be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIDenseSetPreallocation(Mat,PetscScalar[]); 1019284134d9SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation(Mat,PetscScalar[]); 1020be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,PetscInt*[]); 1021be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,PetscInt*[]); 1022be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAdicSetLocalFunction(Mat,void (*)(void)); 1023273d9f13SBarry Smith 1024be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetLDA(Mat,PetscInt); 1025ab92ecdeSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDenseGetLocalMatrix(Mat,Mat*); 10261b807ce4Svictorle 1027be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStoreValues(Mat); 1028be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRetrieveValues(Mat); 10292e8a6d31SBarry Smith 1030be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDAADSetCtx(Mat,void*); 10313a7fca6bSBarry Smith 10327b80b807SBarry Smith /* 10337b80b807SBarry Smith These routines are not usually accessed directly, rather solving is 103494b7f48cSBarry Smith done through the KSP and PC interfaces. 10357b80b807SBarry Smith */ 10367b80b807SBarry Smith 1037d9274352SBarry Smith /*E 1038d9274352SBarry Smith MatOrderingType - String with the name of a PETSc matrix ordering or the creation function 1039d9274352SBarry Smith with an optional dynamic library name, for example 1040d9274352SBarry Smith http://www.mcs.anl.gov/petsc/lib.a:orderingcreate() 1041d9274352SBarry Smith 1042d9274352SBarry Smith Level: beginner 1043d9274352SBarry Smith 1044e5a9bf91SBarry Smith Cannot use const because the PC objects manipulate the string 1045e5a9bf91SBarry Smith 1046d9274352SBarry Smith .seealso: MatGetOrdering() 1047d9274352SBarry Smith E*/ 10483eaad2caSSatish Balay #define MatOrderingType char* 10492692d6eeSBarry Smith #define MATORDERINGNATURAL "natural" 10502692d6eeSBarry Smith #define MATORDERINGND "nd" 10512692d6eeSBarry Smith #define MATORDERING1WD "1wd" 10522692d6eeSBarry Smith #define MATORDERINGRCM "rcm" 10532692d6eeSBarry Smith #define MATORDERINGQMD "qmd" 10542692d6eeSBarry Smith #define MATORDERINGROWLENGTH "rowlength" 10552692d6eeSBarry Smith #define MATORDERINGDSC_ND "dsc_nd" /* these three are only for DSCPACK, see its documentation for details */ 10562692d6eeSBarry Smith #define MATORDERINGDSC_MMD "dsc_mmd" 10572692d6eeSBarry Smith #define MATORDERINGDSC_MDF "dsc_mdf" 10582692d6eeSBarry Smith #define MATORDERINGCONSTRAINED "constrained" 10592692d6eeSBarry Smith #define MATORDERINGIDENTITY "identity" 10602692d6eeSBarry Smith #define MATORDERINGREVERSE "reverse" 10612692d6eeSBarry Smith #define MATORDERINGFLOW "flow" 10622692d6eeSBarry Smith #define MATORDERINGAMD "amd" 1063b12f92e5SBarry Smith 1064f2416520SSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOrdering(Mat,const MatOrderingType,IS*,IS*); 1065be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOrderingList(PetscFList *list); 1066f2416520SSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,const MatOrderingType,IS*,IS*)); 106730de9b25SBarry Smith 106830de9b25SBarry Smith /*MC 10691890ba74SBarry Smith MatOrderingRegisterDynamic - Adds a new sparse matrix ordering to the matrix package. 107030de9b25SBarry Smith 107130de9b25SBarry Smith Synopsis: 10721890ba74SBarry Smith PetscErrorCode MatOrderingRegisterDynamic(const char *name_ordering,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatOrdering)) 107330de9b25SBarry Smith 107430de9b25SBarry Smith Not Collective 107530de9b25SBarry Smith 107630de9b25SBarry Smith Input Parameters: 10772692d6eeSBarry Smith + sname - name of ordering (for example MATORDERINGND) 107830de9b25SBarry Smith . path - location of library where creation routine is 107930de9b25SBarry Smith . name - name of function that creates the ordering type,a string 108030de9b25SBarry Smith - function - function pointer that creates the ordering 108130de9b25SBarry Smith 108230de9b25SBarry Smith Level: developer 108330de9b25SBarry Smith 108430de9b25SBarry Smith If dynamic libraries are used, then the fourth input argument (function) 108530de9b25SBarry Smith is ignored. 108630de9b25SBarry Smith 108730de9b25SBarry Smith Sample usage: 108830de9b25SBarry Smith .vb 108930de9b25SBarry Smith MatOrderingRegisterDynamic("my_order",/home/username/my_lib/lib/libO/solaris/mylib.a, 109030de9b25SBarry Smith "MyOrder",MyOrder); 109130de9b25SBarry Smith .ve 109230de9b25SBarry Smith 109330de9b25SBarry Smith Then, your partitioner can be chosen with the procedural interface via 109430de9b25SBarry Smith $ MatOrderingSetType(part,"my_order) 109530de9b25SBarry Smith or at runtime via the option 10962401956bSBarry Smith $ -pc_factor_mat_ordering_type my_order 109730de9b25SBarry Smith 1098ab901514SBarry Smith ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values. 109930de9b25SBarry Smith 110030de9b25SBarry Smith .keywords: matrix, ordering, register 110130de9b25SBarry Smith 110230de9b25SBarry Smith .seealso: MatOrderingRegisterDestroy(), MatOrderingRegisterAll() 110330de9b25SBarry Smith M*/ 1104aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 1105f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,0) 1106b12f92e5SBarry Smith #else 1107f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,d) 1108b12f92e5SBarry Smith #endif 110930de9b25SBarry Smith 1110be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegisterDestroy(void); 1111be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegisterAll(const char[]); 11122bad1931SBarry Smith extern PetscTruth MatOrderingRegisterAllCalled; 1113b0a32e0cSBarry Smith extern PetscFList MatOrderingList; 1114d4fbbf0eSBarry Smith 1115be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS); 1116a2ce50c7SBarry Smith 1117d91e6319SBarry Smith /*S 1118d90ac83dSHong Zhang MatFactorShiftType - Numeric Shift. 1119d90ac83dSHong Zhang 1120d90ac83dSHong Zhang Level: beginner 1121d90ac83dSHong Zhang 1122d90ac83dSHong Zhang S*/ 1123d90ac83dSHong Zhang typedef enum {MAT_SHIFT_NONE,MAT_SHIFT_NONZERO,MAT_SHIFT_POSITIVE_DEFINITE,MAT_SHIFT_INBLOCKS} MatFactorShiftType; 1124d90ac83dSHong Zhang extern const char *MatFactorShiftTypes[]; 1125d90ac83dSHong Zhang 1126d90ac83dSHong Zhang /*S 11272401956bSBarry Smith MatFactorInfo - Data passed into the matrix factorization routines 1128b00f7748SHong Zhang 112961cad860SBarry Smith In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE, that is use 113061cad860SBarry Smith $ MatFactorInfo info(MAT_FACTORINFO_SIZE) 1131b00f7748SHong Zhang 113215e8a5b3SHong Zhang Notes: These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC. 1133b00f7748SHong Zhang 113461cad860SBarry Smith You can use MatFactorInfoInitialize() to set default values. 113561cad860SBarry Smith 1136b00f7748SHong Zhang Level: developer 1137b00f7748SHong Zhang 1138d7d82daaSBarry Smith .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(), 1139d7d82daaSBarry Smith MatFactorInfoInitialize() 1140b00f7748SHong Zhang 1141b00f7748SHong Zhang S*/ 1142b00f7748SHong Zhang typedef struct { 114315e8a5b3SHong Zhang PetscReal diagonal_fill; /* force diagonal to fill in if initially not filled */ 114485317021SBarry Smith PetscReal usedt; 114515e8a5b3SHong Zhang PetscReal dt; /* drop tolerance */ 1146b00f7748SHong Zhang PetscReal dtcol; /* tolerance for pivoting */ 114715e8a5b3SHong Zhang PetscReal dtcount; /* maximum nonzeros to be allowed per row */ 114867e28bfeSBarry Smith PetscReal fill; /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */ 1149348344bbSBarry Smith PetscReal levels; /* ICC/ILU(levels) */ 1150bcd9e38bSBarry Smith PetscReal pivotinblocks; /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0 1151bcd9e38bSBarry Smith factorization may be faster if do not pivot */ 115215e8a5b3SHong Zhang PetscReal zeropivot; /* pivot is called zero if less than this */ 1153f4db908eSBarry Smith PetscReal shifttype; /* type of shift added to matrix factor to prevent zero pivots */ 1154f4db908eSBarry Smith PetscReal shiftamount; /* how large the shift is */ 115515e8a5b3SHong Zhang } MatFactorInfo; 1156ffa6d0a5SLois Curfman McInnes 1157be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFactorInfoInitialize(MatFactorInfo*); 11580481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactor(Mat,IS,const MatFactorInfo*); 11590481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*); 11600481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactorNumeric(Mat,Mat,const MatFactorInfo*); 11610481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactor(Mat,IS,IS,const MatFactorInfo*); 11620481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUFactor(Mat,IS,IS,const MatFactorInfo*); 11630481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*); 11640481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*); 11650481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatICCFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*); 11660481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatICCFactor(Mat,IS,const MatFactorInfo*); 11670481f469SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactorNumeric(Mat,Mat,const MatFactorInfo*); 1168be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*); 1169be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolve(Mat,Vec,Vec); 1170be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatForwardSolve(Mat,Vec,Vec); 1171be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatBackwardSolve(Mat,Vec,Vec); 1172be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveAdd(Mat,Vec,Vec,Vec); 1173be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveTranspose(Mat,Vec,Vec); 1174be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveTransposeAdd(Mat,Vec,Vec,Vec); 1175be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolves(Mat,Vecs,Vecs); 11768ed539a5SBarry Smith 1177be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUnfactored(Mat); 1178bb5a7306SBarry Smith 1179d91e6319SBarry Smith /*E 1180d91e6319SBarry Smith MatSORType - What type of (S)SOR to perform 1181bb1eb677SSatish Balay 1182d91e6319SBarry Smith Level: beginner 1183d91e6319SBarry Smith 1184d9274352SBarry Smith May be bitwise ORd together 1185d9274352SBarry Smith 1186d91e6319SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 1187d91e6319SBarry Smith 11884e7234bfSBarry Smith MatSORType may be bitwise ORd together, so do not change the numbers 11894e7234bfSBarry Smith 119041f059aeSBarry Smith .seealso: MatSOR() 1191d91e6319SBarry Smith E*/ 1192ee50ffe9SBarry Smith typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3, 1193ee50ffe9SBarry Smith SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8, 1194ee50ffe9SBarry Smith SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16, 119584cb2905SBarry Smith SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType; 119641f059aeSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSOR(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec); 11978ed539a5SBarry Smith 1198d4fbbf0eSBarry Smith /* 1199639f9d9dSBarry Smith These routines are for efficiently computing Jacobians via finite differences. 1200639f9d9dSBarry Smith */ 1201b12f92e5SBarry Smith 1202d9274352SBarry Smith /*E 1203d9274352SBarry Smith MatColoringType - String with the name of a PETSc matrix coloring or the creation function 1204d9274352SBarry Smith with an optional dynamic library name, for example 1205d9274352SBarry Smith http://www.mcs.anl.gov/petsc/lib.a:coloringcreate() 1206d9274352SBarry Smith 1207d9274352SBarry Smith Level: beginner 1208d9274352SBarry Smith 1209d9274352SBarry Smith .seealso: MatGetColoring() 1210d9274352SBarry Smith E*/ 1211a313700dSBarry Smith #define MatColoringType char* 12122692d6eeSBarry Smith #define MATCOLORINGNATURAL "natural" 12132692d6eeSBarry Smith #define MATCOLORINGSL "sl" 12142692d6eeSBarry Smith #define MATCOLORINGLF "lf" 12152692d6eeSBarry Smith #define MATCOLORINGID "id" 1216b12f92e5SBarry Smith 1217ab970d9bSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColoring(Mat,const MatColoringType,ISColoring*); 12182e90c967SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,MatColoringType,ISColoring *)); 121930de9b25SBarry Smith 122030de9b25SBarry Smith /*MC 122130de9b25SBarry Smith MatColoringRegisterDynamic - Adds a new sparse matrix coloring to the 122230de9b25SBarry Smith matrix package. 122330de9b25SBarry Smith 122430de9b25SBarry Smith Synopsis: 12251890ba74SBarry Smith PetscErrorCode MatColoringRegisterDynamic(const char *name_coloring,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatColoring)) 122630de9b25SBarry Smith 122730de9b25SBarry Smith Not Collective 122830de9b25SBarry Smith 122930de9b25SBarry Smith Input Parameters: 12302692d6eeSBarry Smith + sname - name of Coloring (for example MATCOLORINGSL) 123130de9b25SBarry Smith . path - location of library where creation routine is 123230de9b25SBarry Smith . name - name of function that creates the Coloring type, a string 123330de9b25SBarry Smith - function - function pointer that creates the coloring 123430de9b25SBarry Smith 123530de9b25SBarry Smith Level: developer 123630de9b25SBarry Smith 123730de9b25SBarry Smith If dynamic libraries are used, then the fourth input argument (function) 123830de9b25SBarry Smith is ignored. 123930de9b25SBarry Smith 124030de9b25SBarry Smith Sample usage: 124130de9b25SBarry Smith .vb 124230de9b25SBarry Smith MatColoringRegisterDynamic("my_color",/home/username/my_lib/lib/libO/solaris/mylib.a, 124330de9b25SBarry Smith "MyColor",MyColor); 124430de9b25SBarry Smith .ve 124530de9b25SBarry Smith 124630de9b25SBarry Smith Then, your partitioner can be chosen with the procedural interface via 124730de9b25SBarry Smith $ MatColoringSetType(part,"my_color") 124830de9b25SBarry Smith or at runtime via the option 124930de9b25SBarry Smith $ -mat_coloring_type my_color 125030de9b25SBarry Smith 1251ab901514SBarry Smith $PETSC_ARCH occuring in pathname will be replaced with appropriate values. 125230de9b25SBarry Smith 125330de9b25SBarry Smith .keywords: matrix, Coloring, register 125430de9b25SBarry Smith 125530de9b25SBarry Smith .seealso: MatColoringRegisterDestroy(), MatColoringRegisterAll() 125630de9b25SBarry Smith M*/ 1257aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 1258f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,0) 1259b12f92e5SBarry Smith #else 1260f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,d) 1261b12f92e5SBarry Smith #endif 126230de9b25SBarry Smith 12632bad1931SBarry Smith extern PetscTruth MatColoringRegisterAllCalled; 1264f1144a30SSatish Balay 1265be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegisterAll(const char[]); 1266be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegisterDestroy(void); 1267be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*); 1268639f9d9dSBarry Smith 1269d9274352SBarry Smith /*S 1270d9274352SBarry Smith MatFDColoring - Object for computing a sparse Jacobian via finite differences 1271d9274352SBarry Smith and coloring 1272639f9d9dSBarry Smith 1273d9274352SBarry Smith Level: beginner 1274d9274352SBarry Smith 1275d9274352SBarry Smith Concepts: coloring, sparse Jacobian, finite differences 1276d9274352SBarry Smith 1277d9274352SBarry Smith .seealso: MatFDColoringCreate() 1278d9274352SBarry Smith S*/ 1279e2a1c21fSSatish Balay typedef struct _p_MatFDColoring* MatFDColoring; 1280639f9d9dSBarry Smith 1281be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringCreate(Mat,ISColoring,MatFDColoring *); 1282be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringDestroy(MatFDColoring); 1283be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringView(MatFDColoring,PetscViewer); 1284be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*); 12854a9d489dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**); 1286be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal); 1287be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFromOptions(MatFDColoring); 1288be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApply(Mat,MatFDColoring,Vec,MatStructure*,void *); 1289be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApplyTS(Mat,MatFDColoring,PetscReal,Vec,MatStructure*,void *); 1290be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetF(MatFDColoring,Vec); 1291be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,PetscInt*[]); 1292639f9d9dSBarry Smith /* 12930752156aSBarry Smith These routines are for partitioning matrices: currently used only 12943eda8832SBarry Smith for adjacency matrix, MatCreateMPIAdj(). 12950752156aSBarry Smith */ 1296ca161407SBarry Smith 1297d9274352SBarry Smith /*S 1298d9274352SBarry Smith MatPartitioning - Object for managing the partitioning of a matrix or graph 1299d9274352SBarry Smith 1300d9274352SBarry Smith Level: beginner 1301d9274352SBarry Smith 1302d9274352SBarry Smith Concepts: partitioning 1303d9274352SBarry Smith 1304743a1026Svictorle .seealso: MatPartitioningCreate(), MatPartitioningType 1305d9274352SBarry Smith S*/ 130691e9ee9fSBarry Smith typedef struct _p_MatPartitioning* MatPartitioning; 1307d9274352SBarry Smith 1308d9274352SBarry Smith /*E 13095bc6e7ccSvictorle MatPartitioningType - String with the name of a PETSc matrix partitioning or the creation function 1310d9274352SBarry Smith with an optional dynamic library name, for example 1311d9274352SBarry Smith http://www.mcs.anl.gov/petsc/lib.a:partitioningcreate() 1312d9274352SBarry Smith 1313d9274352SBarry Smith Level: beginner 1314d9274352SBarry Smith 1315b547a13bSvictorle .seealso: MatPartitioningCreate(), MatPartitioning 1316d9274352SBarry Smith E*/ 1317a313700dSBarry Smith #define MatPartitioningType char* 13182692d6eeSBarry Smith #define MATPARTITIONINGCURRENT "current" 13192692d6eeSBarry Smith #define MATPARTITIONINGSQUARE "square" 13202692d6eeSBarry Smith #define MATPARTITIONINGPARMETIS "parmetis" 13212692d6eeSBarry Smith #define MATPARTITIONINGCHACO "chaco" 13222692d6eeSBarry Smith #define MATPARTITIONINGJOSTLE "jostle" 13232692d6eeSBarry Smith #define MATPARTITIONINGPARTY "party" 13242692d6eeSBarry Smith #define MATPARTITIONINGSCOTCH "scotch" 132571306c60Sjroman 1326ca161407SBarry Smith 1327be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningCreate(MPI_Comm,MatPartitioning*); 1328a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetType(MatPartitioning,const MatPartitioningType); 1329be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetNParts(MatPartitioning,PetscInt); 1330be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetAdjacency(MatPartitioning,Mat); 1331be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]); 1332be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []); 1333be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningApply(MatPartitioning,IS*); 1334be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningDestroy(MatPartitioning); 13352aabb6bbSBarry Smith 1336be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatPartitioning)); 133730de9b25SBarry Smith 133830de9b25SBarry Smith /*MC 133930de9b25SBarry Smith MatPartitioningRegisterDynamic - Adds a new sparse matrix partitioning to the 134030de9b25SBarry Smith matrix package. 134130de9b25SBarry Smith 134230de9b25SBarry Smith Synopsis: 13431890ba74SBarry Smith PetscErrorCode MatPartitioningRegisterDynamic(const char *name_partitioning,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatPartitioning)) 134430de9b25SBarry Smith 134530de9b25SBarry Smith Not Collective 134630de9b25SBarry Smith 134730de9b25SBarry Smith Input Parameters: 13482692d6eeSBarry Smith + sname - name of partitioning (for example MATPARTITIONINGCURRENT) or parmetis 134930de9b25SBarry Smith . path - location of library where creation routine is 135030de9b25SBarry Smith . name - name of function that creates the partitioning type, a string 135130de9b25SBarry Smith - function - function pointer that creates the partitioning type 135230de9b25SBarry Smith 135330de9b25SBarry Smith Level: developer 135430de9b25SBarry Smith 135530de9b25SBarry Smith If dynamic libraries are used, then the fourth input argument (function) 135630de9b25SBarry Smith is ignored. 135730de9b25SBarry Smith 135830de9b25SBarry Smith Sample usage: 135930de9b25SBarry Smith .vb 136030de9b25SBarry Smith MatPartitioningRegisterDynamic("my_part",/home/username/my_lib/lib/libO/solaris/mylib.a, 136130de9b25SBarry Smith "MyPartCreate",MyPartCreate); 136230de9b25SBarry Smith .ve 136330de9b25SBarry Smith 136430de9b25SBarry Smith Then, your partitioner can be chosen with the procedural interface via 136530de9b25SBarry Smith $ MatPartitioningSetType(part,"my_part") 136630de9b25SBarry Smith or at runtime via the option 136730de9b25SBarry Smith $ -mat_partitioning_type my_part 136830de9b25SBarry Smith 1369ab901514SBarry Smith $PETSC_ARCH occuring in pathname will be replaced with appropriate values. 137030de9b25SBarry Smith 137130de9b25SBarry Smith .keywords: matrix, partitioning, register 137230de9b25SBarry Smith 137330de9b25SBarry Smith .seealso: MatPartitioningRegisterDestroy(), MatPartitioningRegisterAll() 137430de9b25SBarry Smith M*/ 1375aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 1376f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,0) 13772aabb6bbSBarry Smith #else 1378f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,d) 13792aabb6bbSBarry Smith #endif 13802aabb6bbSBarry Smith 13812bad1931SBarry Smith extern PetscTruth MatPartitioningRegisterAllCalled; 1382f1144a30SSatish Balay 1383be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegisterAll(const char[]); 1384be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegisterDestroy(void); 13852bad1931SBarry Smith 1386be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningView(MatPartitioning,PetscViewer); 1387be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetFromOptions(MatPartitioning); 1388a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningGetType(MatPartitioning,const MatPartitioningType*); 1389ca161407SBarry Smith 1390be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningParmetisSetCoarseSequential(MatPartitioning); 13910e8af3caSMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningParmetisGetEdgeCut(MatPartitioning, PetscInt *); 13920752156aSBarry Smith 1393be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningJostleSetCoarseLevel(MatPartitioning,PetscReal); 1394be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningJostleSetCoarseSequential(MatPartitioning); 139571306c60Sjroman 1396954b3ec6SBarry Smith typedef enum {MP_CHACO_MULTILEVEL_KL,MP_CHACO_SPECTRAL,MP_CHACO_LINEAR,MP_CHACO_RANDOM, MP_CHACO_SCATTERED} MPChacoGlobalType; 1397be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetGlobal(MatPartitioning, MPChacoGlobalType); 139871306c60Sjroman typedef enum { MP_CHACO_KERNIGHAN_LIN, MP_CHACO_NONE } MPChacoLocalType; 1399be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetLocal(MatPartitioning, MPChacoLocalType); 1400be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal); 140171306c60Sjroman typedef enum { MP_CHACO_LANCZOS, MP_CHACO_RQI_SYMMLQ } MPChacoEigenType; 1402be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType); 1403be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenTol(MatPartitioning, PetscReal); 1404be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenNumber(MatPartitioning, PetscInt); 140571306c60Sjroman 140671306c60Sjroman #define MP_PARTY_OPT "opt" 140771306c60Sjroman #define MP_PARTY_LIN "lin" 140871306c60Sjroman #define MP_PARTY_SCA "sca" 140971306c60Sjroman #define MP_PARTY_RAN "ran" 141071306c60Sjroman #define MP_PARTY_GBF "gbf" 141171306c60Sjroman #define MP_PARTY_GCF "gcf" 141271306c60Sjroman #define MP_PARTY_BUB "bub" 141371306c60Sjroman #define MP_PARTY_DEF "def" 1414be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetGlobal(MatPartitioning, const char*); 141571306c60Sjroman #define MP_PARTY_HELPFUL_SETS "hs" 141671306c60Sjroman #define MP_PARTY_KERNIGHAN_LIN "kl" 141771306c60Sjroman #define MP_PARTY_NONE "no" 1418be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetLocal(MatPartitioning, const char*); 1419be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal); 1420be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetBipart(MatPartitioning,PetscTruth); 1421be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscTruth); 142271306c60Sjroman 142371306c60Sjroman typedef enum { MP_SCOTCH_GREEDY, MP_SCOTCH_GPS, MP_SCOTCH_GR_GPS } MPScotchGlobalType; 1424be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetArch(MatPartitioning,const char*); 1425be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetMultilevel(MatPartitioning); 1426be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetGlobal(MatPartitioning,MPScotchGlobalType); 1427be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetCoarseLevel(MatPartitioning,PetscReal); 1428be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetHostList(MatPartitioning,const char*); 142971306c60Sjroman typedef enum { MP_SCOTCH_KERNIGHAN_LIN, MP_SCOTCH_NONE } MPScotchLocalType; 1430be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetLocal(MatPartitioning,MPScotchLocalType); 1431be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetMapping(MatPartitioning); 1432be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetStrategy(MatPartitioning,char*); 143371306c60Sjroman 1434591294e4SBarry Smith EXTERN PetscErrorCode MatMeshToVertexGraph(Mat,PetscInt,Mat*); 143514196391SBarry Smith EXTERN PetscErrorCode MatMeshToCellGraph(Mat,PetscInt,Mat*); 1436591294e4SBarry Smith 14370752156aSBarry Smith /* 14380a835dfdSSatish Balay If you add entries here you must also add them to finclude/petscmat.h 1439d4fbbf0eSBarry Smith */ 14401c1c02c0SLois Curfman McInnes typedef enum { MATOP_SET_VALUES=0, 14411c1c02c0SLois Curfman McInnes MATOP_GET_ROW=1, 14421c1c02c0SLois Curfman McInnes MATOP_RESTORE_ROW=2, 14431c1c02c0SLois Curfman McInnes MATOP_MULT=3, 14441c1c02c0SLois Curfman McInnes MATOP_MULT_ADD=4, 14457c922b88SBarry Smith MATOP_MULT_TRANSPOSE=5, 14467c922b88SBarry Smith MATOP_MULT_TRANSPOSE_ADD=6, 14471c1c02c0SLois Curfman McInnes MATOP_SOLVE=7, 14481c1c02c0SLois Curfman McInnes MATOP_SOLVE_ADD=8, 14497c922b88SBarry Smith MATOP_SOLVE_TRANSPOSE=9, 14507c922b88SBarry Smith MATOP_SOLVE_TRANSPOSE_ADD=10, 14511c1c02c0SLois Curfman McInnes MATOP_LUFACTOR=11, 14521c1c02c0SLois Curfman McInnes MATOP_CHOLESKYFACTOR=12, 1453a714c06dSBarry Smith MATOP_SOR=13, 14541c1c02c0SLois Curfman McInnes MATOP_TRANSPOSE=14, 14551c1c02c0SLois Curfman McInnes MATOP_GETINFO=15, 14561c1c02c0SLois Curfman McInnes MATOP_EQUAL=16, 14571c1c02c0SLois Curfman McInnes MATOP_GET_DIAGONAL=17, 14581c1c02c0SLois Curfman McInnes MATOP_DIAGONAL_SCALE=18, 14591c1c02c0SLois Curfman McInnes MATOP_NORM=19, 14601c1c02c0SLois Curfman McInnes MATOP_ASSEMBLY_BEGIN=20, 14611c1c02c0SLois Curfman McInnes MATOP_ASSEMBLY_END=21, 1462d519adbfSMatthew Knepley MATOP_SET_OPTION=22, 1463d519adbfSMatthew Knepley MATOP_ZERO_ENTRIES=23, 1464d519adbfSMatthew Knepley MATOP_ZERO_ROWS=24, 1465d519adbfSMatthew Knepley MATOP_LUFACTOR_SYMBOLIC=25, 1466d519adbfSMatthew Knepley MATOP_LUFACTOR_NUMERIC=26, 1467d519adbfSMatthew Knepley MATOP_CHOLESKY_FACTOR_SYMBOLIC=27, 1468d519adbfSMatthew Knepley MATOP_CHOLESKY_FACTOR_NUMERIC=28, 1469d519adbfSMatthew Knepley MATOP_SETUP_PREALLOCATION=29, 1470d519adbfSMatthew Knepley MATOP_ILUFACTOR_SYMBOLIC=30, 1471d519adbfSMatthew Knepley MATOP_ICCFACTOR_SYMBOLIC=31, 1472d519adbfSMatthew Knepley MATOP_GET_ARRAY=32, 1473d519adbfSMatthew Knepley MATOP_RESTORE_ARRAY=33, 1474d519adbfSMatthew Knepley MATOP_DUPLICATE=34, 1475d519adbfSMatthew Knepley MATOP_FORWARD_SOLVE=35, 1476d519adbfSMatthew Knepley MATOP_BACKWARD_SOLVE=36, 1477d519adbfSMatthew Knepley MATOP_ILUFACTOR=37, 1478d519adbfSMatthew Knepley MATOP_ICCFACTOR=38, 1479d519adbfSMatthew Knepley MATOP_AXPY=39, 1480d519adbfSMatthew Knepley MATOP_GET_SUBMATRICES=40, 1481d519adbfSMatthew Knepley MATOP_INCREASE_OVERLAP=41, 1482d519adbfSMatthew Knepley MATOP_GET_VALUES=42, 1483d519adbfSMatthew Knepley MATOP_COPY=43, 1484d519adbfSMatthew Knepley MATOP_GET_ROW_MAX=44, 1485d519adbfSMatthew Knepley MATOP_SCALE=45, 1486d519adbfSMatthew Knepley MATOP_SHIFT=46, 148735153367SBarry Smith MATOP_DIAGONAL_SET=47, 1488d519adbfSMatthew Knepley MATOP_ILUDT_FACTOR=48, 1489d519adbfSMatthew Knepley MATOP_SET_BLOCK_SIZE=49, 1490d519adbfSMatthew Knepley MATOP_GET_ROW_IJ=50, 1491d519adbfSMatthew Knepley MATOP_RESTORE_ROW_IJ=51, 1492d519adbfSMatthew Knepley MATOP_GET_COLUMN_IJ=52, 1493d519adbfSMatthew Knepley MATOP_RESTORE_COLUMN_IJ=53, 1494d519adbfSMatthew Knepley MATOP_FDCOLORING_CREATE=54, 1495d519adbfSMatthew Knepley MATOP_COLORING_PATCH=55, 1496d519adbfSMatthew Knepley MATOP_SET_UNFACTORED=56, 1497d519adbfSMatthew Knepley MATOP_PERMUTE=57, 1498d519adbfSMatthew Knepley MATOP_SET_VALUES_BLOCKED=58, 1499d519adbfSMatthew Knepley MATOP_GET_SUBMATRIX=59, 1500d519adbfSMatthew Knepley MATOP_DESTROY=60, 1501d519adbfSMatthew Knepley MATOP_VIEW=61, 1502d519adbfSMatthew Knepley MATOP_CONVERT_FROM=62, 1503d519adbfSMatthew Knepley MATOP_USE_SCALED_FORM=63, 1504d519adbfSMatthew Knepley MATOP_SCALE_SYSTEM=64, 1505d519adbfSMatthew Knepley MATOP_UNSCALE_SYSTEM=65, 1506d519adbfSMatthew Knepley MATOP_SET_LOCAL_TO_GLOBAL_MAP=66, 1507d519adbfSMatthew Knepley MATOP_SET_VALUES_LOCAL=67, 1508d519adbfSMatthew Knepley MATOP_ZERO_ROWS_LOCAL=68, 1509d519adbfSMatthew Knepley MATOP_GET_ROW_MAX_ABS=69, 1510d519adbfSMatthew Knepley MATOP_GET_ROW_MIN_ABS=70, 1511d519adbfSMatthew Knepley MATOP_CONVERT=71, 1512d519adbfSMatthew Knepley MATOP_SET_COLORING=72, 1513d519adbfSMatthew Knepley MATOP_SET_VALUES_ADIC=73, 1514d519adbfSMatthew Knepley MATOP_SET_VALUES_ADIFOR=74, 1515d519adbfSMatthew Knepley MATOP_FD_COLORING_APPLY=75, 1516d519adbfSMatthew Knepley MATOP_SET_FROM_OPTIONS=76, 1517d519adbfSMatthew Knepley MATOP_MULT_CON=77, 1518d519adbfSMatthew Knepley MATOP_MULT_TRANSPOSE_CON=78, 1519d519adbfSMatthew Knepley MATOP_PERMUTE_SPARSIFY=79, 1520d519adbfSMatthew Knepley MATOP_MULT_MULTIPLE=80, 1521d519adbfSMatthew Knepley MATOP_SOLVE_MULTIPLE=81, 1522d519adbfSMatthew Knepley MATOP_GET_INERTIA=82, 1523d519adbfSMatthew Knepley MATOP_LOAD=83, 1524d519adbfSMatthew Knepley MATOP_IS_SYMMETRIC=84, 1525d519adbfSMatthew Knepley MATOP_IS_HERMITIAN=85, 1526d519adbfSMatthew Knepley MATOP_IS_STRUCTURALLY_SYMMETRIC=86, 152741f059aeSBarry Smith MATOP_DUMMY=87, 1528d519adbfSMatthew Knepley MATOP_GET_VECS=88, 1529d519adbfSMatthew Knepley MATOP_MAT_MULT=89, 1530d519adbfSMatthew Knepley MATOP_MAT_MULT_SYMBOLIC=90, 1531d519adbfSMatthew Knepley MATOP_MAT_MULT_NUMERIC=91, 1532d519adbfSMatthew Knepley MATOP_PTAP=92, 1533d519adbfSMatthew Knepley MATOP_PTAP_SYMBOLIC=93, 1534d519adbfSMatthew Knepley MATOP_PTAP_NUMERIC=94, 1535d519adbfSMatthew Knepley MATOP_MAT_MULTTRANSPOSE=95, 15361763ac4eSSatish Balay MATOP_MAT_MULTTRANSPOSE_SYM=96, 15371763ac4eSSatish Balay MATOP_MAT_MULTTRANSPOSE_NUM=97, 1538d519adbfSMatthew Knepley MATOP_PTAP_SYMBOLIC_SEQAIJ=98, 1539d519adbfSMatthew Knepley MATOP_PTAP_NUMERIC_SEQAIJ=99, 1540d519adbfSMatthew Knepley MATOP_PTAP_SYMBOLIC_MPIAIJ=100, 1541d519adbfSMatthew Knepley MATOP_PTAP_NUMERIC_MPIAIJ=101, 1542d519adbfSMatthew Knepley MATOP_CONJUGATE=102, 1543d519adbfSMatthew Knepley MATOP_SET_SIZES=103, 1544d519adbfSMatthew Knepley MATOP_SET_VALUES_ROW=104, 1545d519adbfSMatthew Knepley MATOP_REAL_PART=105, 1546d519adbfSMatthew Knepley MATOP_IMAG_PART=106, 1547d519adbfSMatthew Knepley MATOP_GET_ROW_UTRIANGULAR=107, 1548d519adbfSMatthew Knepley MATOP_RESTORE_ROW_UTRIANGULAR=108, 1549d519adbfSMatthew Knepley MATOP_MATSOLVE=109, 1550d519adbfSMatthew Knepley MATOP_GET_REDUNDANTMATRIX=110, 1551d519adbfSMatthew Knepley MATOP_GET_ROW_MIN=111, 1552d519adbfSMatthew Knepley MATOP_GET_COLUMN_VEC=112, 1553d519adbfSMatthew Knepley MATOP_MISSING_DIAGONAL=113, 1554d519adbfSMatthew Knepley MATOP_MATGETSEQNONZEROSTRUCTURE=114, 155589c6957cSBarry Smith MATOP_CREATE=115, 155689c6957cSBarry Smith MATOP_GET_GHOSTS=116, 1557eeffb40dSHong Zhang MATOP_MULT_DIAGONAL_BLOCK=119, 1558547795f9SHong Zhang MATOP_HERMITIANTRANSPOSE=120, 1559547795f9SHong Zhang MATOP_MULTHERMITIANTRANSPOSE=121, 1560fbdbba38SShri Abhyankar MATOP_MULTHERMITIANTRANSPOSEADD=122, 1561fbdbba38SShri Abhyankar MATOP_LOADNEW=123, 1562fae171e0SBarry Smith } MatOperation; 1563be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatHasOperation(Mat,MatOperation,PetscTruth*); 1564be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellSetOperation(Mat,MatOperation,void(*)(void)); 1565be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellGetOperation(Mat,MatOperation,void(**)(void)); 1566be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellSetContext(Mat,void*); 1567112a2221SBarry Smith 156890ace30eSBarry Smith /* 156990ace30eSBarry Smith Codes for matrices stored on disk. By default they are 157090ace30eSBarry Smith stored in a universal format. By changing the format with 15717973ad04SBarry Smith PetscViewerSetFormat(viewer,PETSC_VIEWER_NATIVE); the matrices will 157290ace30eSBarry Smith be stored in a way natural for the matrix, for example dense matrices 157390ace30eSBarry Smith would be stored as dense. Matrices stored this way may only be 157490ace30eSBarry Smith read into matrices of the same time. 157590ace30eSBarry Smith */ 157690ace30eSBarry Smith #define MATRIX_BINARY_FORMAT_DENSE -1 157790ace30eSBarry Smith 1578be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetHashTableFactor(Mat,PetscReal); 1579be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatISGetLocalMat(Mat,Mat*); 15801f4e1ec7SBarry Smith 1581d9274352SBarry Smith /*S 1582d9274352SBarry Smith MatNullSpace - Object that removes a null space from a vector, i.e. 1583d9274352SBarry Smith orthogonalizes the vector to a subsapce 1584d9274352SBarry Smith 1585f7a9e4ceSBarry Smith Level: advanced 1586d9274352SBarry Smith 1587d9274352SBarry Smith Concepts: matrix; linear operator, null space 1588d9274352SBarry Smith 15896e1639daSBarry Smith Users manual sections: 15906e1639daSBarry Smith . sec_singular 15916e1639daSBarry Smith 1592d9274352SBarry Smith .seealso: MatNullSpaceCreate() 1593d9274352SBarry Smith S*/ 159474637425SBarry Smith typedef struct _p_MatNullSpace* MatNullSpace; 1595d9274352SBarry Smith 1596be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceCreate(MPI_Comm,PetscTruth,PetscInt,const Vec[],MatNullSpace*); 1597b22b330cSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(MatNullSpace,Vec,void*),void*); 1598be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceDestroy(MatNullSpace); 1599be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceRemove(MatNullSpace,Vec,Vec*); 1600be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceAttach(Mat,MatNullSpace); 160195902228SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceTest(MatNullSpace,Mat,PetscTruth *); 160274637425SBarry Smith 1603be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatReorderingSeqSBAIJ(Mat,IS); 1604be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPISBAIJSetHashTableFactor(Mat,PetscReal); 1605be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetColumnIndices(Mat,PetscInt *); 160646129b97SKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJInvertBlockDiagonal(Mat); 16073f1d51d7SBarry Smith 1608be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMAIJ(Mat,PetscInt,Mat*); 1609be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMAIJRedimension(Mat,PetscInt,Mat*); 1610be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMAIJGetAIJ(Mat,Mat*); 1611c4f061fbSSatish Balay 1612be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatComputeExplicitOperator(Mat,Mat*); 1613b0a32e0cSBarry Smith 1614be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalScaleLocal(Mat,Vec); 161504f1ad80SBarry Smith 1616fef1beadSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*); 16173ec795f1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetBase(Mat,Vec,Vec); 16183ec795f1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*); 1619e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*)); 1620e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec)); 1621e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDAddNullSpace(Mat,MatNullSpace); 1622e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt); 1623e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDResetHHistory(Mat); 1624e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctionError(Mat,PetscReal); 1625e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetPeriod(Mat,PetscInt); 1626e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDGetH(Mat,PetscScalar *); 16276aa9148fSLisandro Dalcin EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetOptionsPrefix(Mat,const char[]); 1628e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFromOptions(Mat); 1629e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*); 1630e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*); 1631e884886fSBarry Smith 16326370053bSBarry Smith /*S 16336370053bSBarry Smith MatMFFD - A data structured used to manage the computation of the h differencing parameter for matrix-free 16346370053bSBarry Smith Jacobian vector products 1635e884886fSBarry Smith 16366370053bSBarry Smith Notes: MATMFFD is a specific MatType which uses the MatMFFD data structure 16376370053bSBarry Smith 16386370053bSBarry Smith MatMFFD*() methods actually take the Mat as their first argument. Not a MatMFFD data structure 16396370053bSBarry Smith 16406370053bSBarry Smith Level: developer 16416370053bSBarry Smith 16426370053bSBarry Smith .seealso: MATMFFD, MatCreateMFFD(), MatMFFDSetFuction(), MatMFFDSetType(), MatMFFDRegister() 16436370053bSBarry Smith S*/ 1644e884886fSBarry Smith typedef struct _p_MatMFFD* MatMFFD; 1645e884886fSBarry Smith 1646e884886fSBarry Smith /*E 1647e884886fSBarry Smith MatMFFDType - algorithm used to compute the h used in computing matrix-vector products via differencing of the function 1648e884886fSBarry Smith 1649e884886fSBarry Smith Level: beginner 1650e884886fSBarry Smith 1651e884886fSBarry Smith .seealso: MatMFFDSetType(), MatMFFDRegister() 1652e884886fSBarry Smith E*/ 1653a313700dSBarry Smith #define MatMFFDType char* 1654e884886fSBarry Smith #define MATMFFD_DS "ds" 1655e884886fSBarry Smith #define MATMFFD_WP "wp" 1656e884886fSBarry Smith 1657a313700dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetType(Mat,const MatMFFDType); 1658e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatMFFD)); 1659e884886fSBarry Smith 1660e884886fSBarry Smith /*MC 1661e884886fSBarry Smith MatMFFDRegisterDynamic - Adds a method to the MatMFFD registry. 1662e884886fSBarry Smith 1663e884886fSBarry Smith Synopsis: 16641890ba74SBarry Smith PetscErrorCode MatMFFDRegisterDynamic(const char *name_solver,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatMFFD)) 1665e884886fSBarry Smith 1666e884886fSBarry Smith Not Collective 1667e884886fSBarry Smith 1668e884886fSBarry Smith Input Parameters: 1669e884886fSBarry Smith + name_solver - name of a new user-defined compute-h module 1670e884886fSBarry Smith . path - path (either absolute or relative) the library containing this solver 1671e884886fSBarry Smith . name_create - name of routine to create method context 1672e884886fSBarry Smith - routine_create - routine to create method context 1673e884886fSBarry Smith 1674e884886fSBarry Smith Level: developer 1675e884886fSBarry Smith 1676e884886fSBarry Smith Notes: 1677e884886fSBarry Smith MatMFFDRegisterDynamic() may be called multiple times to add several user-defined solvers. 1678e884886fSBarry Smith 1679e884886fSBarry Smith If dynamic libraries are used, then the fourth input argument (routine_create) 1680e884886fSBarry Smith is ignored. 1681e884886fSBarry Smith 1682e884886fSBarry Smith Sample usage: 1683e884886fSBarry Smith .vb 1684e884886fSBarry Smith MatMFFDRegisterDynamic("my_h",/home/username/my_lib/lib/libO/solaris/mylib.a, 1685e884886fSBarry Smith "MyHCreate",MyHCreate); 1686e884886fSBarry Smith .ve 1687e884886fSBarry Smith 1688e884886fSBarry Smith Then, your solver can be chosen with the procedural interface via 1689e884886fSBarry Smith $ MatMFFDSetType(mfctx,"my_h") 1690e884886fSBarry Smith or at runtime via the option 1691e884886fSBarry Smith $ -snes_mf_type my_h 1692e884886fSBarry Smith 1693e884886fSBarry Smith .keywords: MatMFFD, register 1694e884886fSBarry Smith 1695e884886fSBarry Smith .seealso: MatMFFDRegisterAll(), MatMFFDRegisterDestroy() 1696e884886fSBarry Smith M*/ 1697e884886fSBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 1698e884886fSBarry Smith #define MatMFFDRegisterDynamic(a,b,c,d) MatMFFDRegister(a,b,c,0) 1699e884886fSBarry Smith #else 1700e884886fSBarry Smith #define MatMFFDRegisterDynamic(a,b,c,d) MatMFFDRegister(a,b,c,d) 1701e884886fSBarry Smith #endif 1702e884886fSBarry Smith 1703e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDRegisterAll(const char[]); 1704e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDRegisterDestroy(void); 17051d0fab5eSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDDSSetUmin(Mat,PetscReal); 1706e884886fSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDWPSetComputeNormU(Mat,PetscTruth); 1707e884886fSBarry Smith 1708e884886fSBarry Smith 1709be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *); 1710be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *); 17117dbadf16SMatthew Knepley 171297969023SHong Zhang /* 171397969023SHong Zhang PETSc interface to MUMPS 171497969023SHong Zhang */ 171597969023SHong Zhang #ifdef PETSC_HAVE_MUMPS 1716a1f19f5aSHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetMumpsIcntl(Mat,PetscInt,PetscInt); 171797969023SHong Zhang #endif 171823a5497aSJed Brown 171923a5497aSJed Brown PETSC_EXTERN_CXX_END 172023a5497aSJed Brown #endif 1721