12eac72dbSBarry Smith /* 22eac72dbSBarry Smith Include file for the matrix component of PETSc 32eac72dbSBarry Smith */ 40a835dfdSSatish Balay #ifndef __PETSCMAT_H 50a835dfdSSatish Balay #define __PETSCMAT_H 60a835dfdSSatish Balay #include "petscvec.h" 7e9fa29b7SSatish Balay PETSC_EXTERN_CXX_BEGIN 82eac72dbSBarry Smith 9d9274352SBarry Smith /*S 10d9274352SBarry Smith Mat - Abstract PETSc matrix object 112eac72dbSBarry Smith 12d91e6319SBarry Smith Level: beginner 13d91e6319SBarry Smith 14d9274352SBarry Smith Concepts: matrix; linear operator 15d9274352SBarry Smith 16d9274352SBarry Smith .seealso: MatCreate(), MatType, MatSetType() 17d9274352SBarry Smith S*/ 18d9274352SBarry Smith typedef struct _p_Mat* Mat; 19d9274352SBarry Smith 20d9274352SBarry Smith /*E 21d9274352SBarry Smith MatType - String with the name of a PETSc matrix or the creation function 22d9274352SBarry Smith with an optional dynamic library name, for example 23d9274352SBarry Smith http://www.mcs.anl.gov/petsc/lib.a:mymatcreate() 24d9274352SBarry Smith 25d9274352SBarry Smith Level: beginner 26d9274352SBarry Smith 27c7393fdbSBarry Smith .seealso: MatSetType(), Mat, MatSolverPackage 28d91e6319SBarry Smith E*/ 29a313700dSBarry Smith #define MatType char* 30273d9f13SBarry Smith #define MATSAME "same" 315a11e1b2SBarry Smith #define MATMAIJ "maij" 32273d9f13SBarry Smith #define MATSEQMAIJ "seqmaij" 33273d9f13SBarry Smith #define MATMPIMAIJ "mpimaij" 34273d9f13SBarry Smith #define MATIS "is" 355a11e1b2SBarry Smith #define MATAIJ "aij" 36273d9f13SBarry Smith #define MATSEQAIJ "seqaij" 377d6a0e61SBarry Smith #define MATSEQAIJPTHREAD "seqaijpthread" 38bf2c1783SBarry Smith #define MATAIJPTHREAD "aijpthread" 39273d9f13SBarry Smith #define MATMPIAIJ "mpiaij" 405a11e1b2SBarry Smith #define MATAIJCRL "aijcrl" 415a11e1b2SBarry Smith #define MATSEQAIJCRL "seqaijcrl" 425a11e1b2SBarry Smith #define MATMPIAIJCRL "mpiaijcrl" 438154be41SBarry Smith #define MATAIJCUSP "aijcusp" 448154be41SBarry Smith #define MATSEQAIJCUSP "seqaijcusp" 458154be41SBarry Smith #define MATMPIAIJCUSP "mpiaijcusp" 465a11e1b2SBarry Smith #define MATAIJPERM "aijperm" 475a11e1b2SBarry Smith #define MATSEQAIJPERM "seqaijperm" 485a11e1b2SBarry Smith #define MATMPIAIJPERM "mpiaijperm" 49273d9f13SBarry Smith #define MATSHELL "shell" 505a11e1b2SBarry Smith #define MATDENSE "dense" 51209238afSKris Buschelman #define MATSEQDENSE "seqdense" 52273d9f13SBarry Smith #define MATMPIDENSE "mpidense" 535a11e1b2SBarry Smith #define MATBAIJ "baij" 54273d9f13SBarry Smith #define MATSEQBAIJ "seqbaij" 55273d9f13SBarry Smith #define MATMPIBAIJ "mpibaij" 56273d9f13SBarry Smith #define MATMPIADJ "mpiadj" 575a11e1b2SBarry Smith #define MATSBAIJ "sbaij" 58273d9f13SBarry Smith #define MATSEQSBAIJ "seqsbaij" 59273d9f13SBarry Smith #define MATMPISBAIJ "mpisbaij" 60c0cdd4a1SDahai Guo 61c0cdd4a1SDahai Guo #define MATSEQBSTRM "seqbstrm" 62c0cdd4a1SDahai Guo #define MATMPIBSTRM "mpibstrm" 63c0cdd4a1SDahai Guo #define MATBSTRM "bstrm" 64aa5a9175SDahai Guo #define MATSEQSBSTRM "seqsbstrm" 65aa5a9175SDahai Guo #define MATMPISBSTRM "mpisbstrm" 66aa5a9175SDahai Guo #define MATSBSTRM "sbstrm" 67c0cdd4a1SDahai Guo 68cebc7f6cSBarry Smith #define MATDAAD "daad" 69cebc7f6cSBarry Smith #define MATMFFD "mffd" 70c8a8475eSBarry Smith #define MATNORMAL "normal" 71ab92ecdeSBarry Smith #define MATLRC "lrc" 722a6744ebSBarry Smith #define MATSCATTER "scatter" 73421e10b8SBarry Smith #define MATBLOCKMAT "blockmat" 74793850ffSBarry Smith #define MATCOMPOSITE "composite" 751f8c7532SHong Zhang #define MATFFT "fft" 761f8c7532SHong Zhang #define MATFFTW "fftw" 77e133240eSMatthew G Knepley #define MATSEQCUFFT "seqcufft" 78557cca28SSatish Balay #define MATTRANSPOSEMAT "transpose" 7972ca8751SBarry Smith #define MATSCHURCOMPLEMENT "schurcomplement" 801d6018f0SLisandro Dalcin #define MATPYTHON "python" 81f91d8e95SBarry Smith #define MATHYPRESTRUCT "hyprestruct" 82a9e6138eSGlenn Hammond #define MATHYPRESSTRUCT "hypresstruct" 83ee1cef2cSJed Brown #define MATSUBMATRIX "submatrix" 842c0dbf93SJed Brown #define MATLOCALREF "localref" 85d8588912SDave May #define MATNEST "nest" 86773941d6SBarry Smith 879989ab13SBarry Smith /*E 88c7393fdbSBarry Smith MatSolverPackage - String with the name of a PETSc matrix solver type. 899989ab13SBarry Smith 909989ab13SBarry Smith For example: "petsc" indicates what PETSc provides, "superlu" indicates either 919989ab13SBarry Smith SuperLU or SuperLU_Dist etc. 929989ab13SBarry Smith 939989ab13SBarry Smith 949989ab13SBarry Smith Level: beginner 959989ab13SBarry Smith 965c9eb25fSBarry Smith .seealso: MatGetFactor(), Mat, MatSetType(), MatType 979989ab13SBarry Smith E*/ 98c7393fdbSBarry Smith #define MatSolverPackage char* 992692d6eeSBarry Smith #define MATSOLVERSPOOLES "spooles" 1002692d6eeSBarry Smith #define MATSOLVERSUPERLU "superlu" 1012692d6eeSBarry Smith #define MATSOLVERSUPERLU_DIST "superlu_dist" 1022692d6eeSBarry Smith #define MATSOLVERUMFPACK "umfpack" 10320db9a53SJed Brown #define MATSOLVERCHOLMOD "cholmod" 1042692d6eeSBarry Smith #define MATSOLVERESSL "essl" 1052692d6eeSBarry Smith #define MATSOLVERLUSOL "lusol" 1062692d6eeSBarry Smith #define MATSOLVERMUMPS "mumps" 1072692d6eeSBarry Smith #define MATSOLVERPASTIX "pastix" 1082692d6eeSBarry Smith #define MATSOLVERDSCPACK "dscpack" 1092692d6eeSBarry Smith #define MATSOLVERMATLAB "matlab" 1102692d6eeSBarry Smith #define MATSOLVERPETSC "petsc" 1112692d6eeSBarry Smith #define MATSOLVERPLAPACK "plapack" 1122692d6eeSBarry Smith #define MATSOLVERBAS "bas" 113773941d6SBarry Smith 114aa5a9175SDahai Guo #define MATSOLVERBSTRM "bstrm" 115aa5a9175SDahai Guo #define MATSOLVERSBSTRM "sbstrm" 116c0cdd4a1SDahai Guo 117b24902e0SBarry Smith /*E 118b24902e0SBarry Smith MatFactorType - indicates what type of factorization is requested 119b24902e0SBarry Smith 120b24902e0SBarry Smith Level: beginner 121b24902e0SBarry Smith 122b24902e0SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 123b24902e0SBarry Smith 124c7393fdbSBarry Smith .seealso: MatSolverPackage, MatGetFactor() 125b24902e0SBarry Smith E*/ 126599ef60dSHong Zhang typedef enum {MAT_FACTOR_NONE, MAT_FACTOR_LU, MAT_FACTOR_CHOLESKY, MAT_FACTOR_ILU, MAT_FACTOR_ICC,MAT_FACTOR_ILUDT} MatFactorType; 12734918c53SJed Brown extern const char *const MatFactorTypes[]; 128e92e720dSBarry Smith 1297087cfbeSBarry Smith extern PetscErrorCode MatGetFactor(Mat,const MatSolverPackage,MatFactorType,Mat*); 1307087cfbeSBarry Smith extern PetscErrorCode MatGetFactorAvailable(Mat,const MatSolverPackage,MatFactorType,PetscBool *); 1317087cfbeSBarry Smith extern PetscErrorCode MatFactorGetSolverPackage(Mat,const MatSolverPackage*); 1327087cfbeSBarry Smith extern PetscErrorCode MatGetFactorType(Mat,MatFactorType*); 1339989ab13SBarry Smith 134c06d978dSMatthew Knepley /* Logging support */ 1350700a824SBarry Smith #define MAT_FILE_CLASSID 1211216 /* used to indicate matrices in binary files */ 1367087cfbeSBarry Smith extern PetscClassId MAT_CLASSID; 1377087cfbeSBarry Smith extern PetscClassId MAT_FDCOLORING_CLASSID; 1387087cfbeSBarry Smith extern PetscClassId MAT_PARTITIONING_CLASSID; 1397087cfbeSBarry Smith extern PetscClassId MAT_NULLSPACE_CLASSID; 1407087cfbeSBarry Smith extern PetscClassId MATMFFD_CLASSID; 141c06d978dSMatthew Knepley 142ceb03754SKris Buschelman /*E 143ceb03754SKris Buschelman MatReuse - Indicates if matrices obtained from a previous call to MatGetSubMatrices() 144d6eff37eSBarry Smith or MatGetSubMatrix() are to be reused to store the new matrix values. For MatConvert() is used to indicate 145d6eff37eSBarry Smith that the input matrix is to be replaced with the converted matrix. 146ceb03754SKris Buschelman 147ceb03754SKris Buschelman Level: beginner 148ceb03754SKris Buschelman 149ceb03754SKris Buschelman Any additions/changes here MUST also be made in include/finclude/petscmat.h 150ceb03754SKris Buschelman 1510c451bc4SBarry Smith .seealso: MatGetSubMatrices(), MatGetSubMatrix(), MatDestroyMatrices(), MatConvert() 152ceb03754SKris Buschelman E*/ 153dfe085dbSJed Brown typedef enum {MAT_INITIAL_MATRIX,MAT_REUSE_MATRIX,MAT_IGNORE_MATRIX} MatReuse; 154ceb03754SKris Buschelman 1555494a064SHong Zhang /*E 1565494a064SHong Zhang MatGetSubMatrixOption - Indicates if matrices obtained from a call to MatGetSubMatrices() 157829201f2SHong Zhang include the matrix values. Currently it is only used by MatGetSeqNonzerostructure(). 1585494a064SHong Zhang 1595494a064SHong Zhang Level: beginner 1605494a064SHong Zhang 161829201f2SHong Zhang .seealso: MatGetSeqNonzerostructure() 1625494a064SHong Zhang E*/ 1635494a064SHong Zhang typedef enum {MAT_DO_NOT_GET_VALUES,MAT_GET_VALUES} MatGetSubMatrixOption; 1645494a064SHong Zhang 1657087cfbeSBarry Smith extern PetscErrorCode MatInitializePackage(const char[]); 166c06d978dSMatthew Knepley 1677087cfbeSBarry Smith extern PetscErrorCode MatCreate(MPI_Comm,Mat*); 168a69a960dSMatthew Knepley PetscPolymorphicFunction(MatCreate,(MPI_Comm comm),(comm,&A),Mat,A) 169a69a960dSMatthew Knepley PetscPolymorphicFunction(MatCreate,(),(PETSC_COMM_WORLD,&A),Mat,A) 1707087cfbeSBarry Smith extern PetscErrorCode MatSetSizes(Mat,PetscInt,PetscInt,PetscInt,PetscInt); 1717087cfbeSBarry Smith extern PetscErrorCode MatSetType(Mat,const MatType); 1727087cfbeSBarry Smith extern PetscErrorCode MatSetFromOptions(Mat); 1737087cfbeSBarry Smith extern PetscErrorCode MatSetUpPreallocation(Mat); 1747087cfbeSBarry Smith extern PetscErrorCode MatRegisterAll(const char[]); 1757087cfbeSBarry Smith extern PetscErrorCode MatRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat)); 17601bebe75SBarry Smith extern PetscErrorCode MatRegisterBaseName(const char[],const char[],const char[]); 1777087cfbeSBarry Smith extern PetscErrorCode MatSetOptionsPrefix(Mat,const char[]); 1787087cfbeSBarry Smith extern PetscErrorCode MatAppendOptionsPrefix(Mat,const char[]); 1797087cfbeSBarry Smith extern PetscErrorCode MatGetOptionsPrefix(Mat,const char*[]); 180f69a0ea3SMatthew Knepley 18130de9b25SBarry Smith /*MC 18230de9b25SBarry Smith MatRegisterDynamic - Adds a new matrix type 18330de9b25SBarry Smith 18430de9b25SBarry Smith Synopsis: 1851890ba74SBarry Smith PetscErrorCode MatRegisterDynamic(const char *name,const char *path,const char *name_create,PetscErrorCode (*routine_create)(Mat)) 18630de9b25SBarry Smith 18730de9b25SBarry Smith Not Collective 18830de9b25SBarry Smith 18930de9b25SBarry Smith Input Parameters: 19030de9b25SBarry Smith + name - name of a new user-defined matrix type 19130de9b25SBarry Smith . path - path (either absolute or relative) the library containing this solver 19230de9b25SBarry Smith . name_create - name of routine to create method context 19330de9b25SBarry Smith - routine_create - routine to create method context 19430de9b25SBarry Smith 19530de9b25SBarry Smith Notes: 19630de9b25SBarry Smith MatRegisterDynamic() may be called multiple times to add several user-defined solvers. 19730de9b25SBarry Smith 19830de9b25SBarry Smith If dynamic libraries are used, then the fourth input argument (routine_create) 19930de9b25SBarry Smith is ignored. 20030de9b25SBarry Smith 20130de9b25SBarry Smith Sample usage: 20230de9b25SBarry Smith .vb 20330de9b25SBarry Smith MatRegisterDynamic("my_mat",/home/username/my_lib/lib/libO/solaris/mylib.a, 20430de9b25SBarry Smith "MyMatCreate",MyMatCreate); 20530de9b25SBarry Smith .ve 20630de9b25SBarry Smith 20730de9b25SBarry Smith Then, your solver can be chosen with the procedural interface via 20830de9b25SBarry Smith $ MatSetType(Mat,"my_mat") 20930de9b25SBarry Smith or at runtime via the option 21030de9b25SBarry Smith $ -mat_type my_mat 21130de9b25SBarry Smith 21230de9b25SBarry Smith Level: advanced 21330de9b25SBarry Smith 214ab901514SBarry Smith Notes: ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values. 21530de9b25SBarry Smith If your function is not being put into a shared library then use VecRegister() instead 21630de9b25SBarry Smith 21730de9b25SBarry Smith .keywords: Mat, register 21830de9b25SBarry Smith 21930de9b25SBarry Smith .seealso: MatRegisterAll(), MatRegisterDestroy() 22030de9b25SBarry Smith 22130de9b25SBarry Smith M*/ 222273d9f13SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 223273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,0) 224273d9f13SBarry Smith #else 225273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,d) 22630de9b25SBarry Smith #endif 22730de9b25SBarry Smith 228ace3abfcSBarry Smith extern PetscBool MatRegisterAllCalled; 229b0a32e0cSBarry Smith extern PetscFList MatList; 230b022a5c1SBarry Smith extern PetscFList MatColoringList; 231b022a5c1SBarry Smith extern PetscFList MatPartitioningList; 23228988994SBarry Smith 2333b224e63SBarry Smith /*E 2343b224e63SBarry Smith MatStructure - Indicates if the matrix has the same nonzero structure 2353b224e63SBarry Smith 2363b224e63SBarry Smith Level: beginner 2373b224e63SBarry Smith 2383b224e63SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 2393b224e63SBarry Smith 2403b224e63SBarry Smith .seealso: MatCopy(), KSPSetOperators(), PCSetOperators() 2413b224e63SBarry Smith E*/ 242214bc6a2SJed Brown typedef enum {DIFFERENT_NONZERO_PATTERN,SUBSET_NONZERO_PATTERN,SAME_NONZERO_PATTERN,SAME_PRECONDITIONER} MatStructure; 2433b224e63SBarry Smith 2447087cfbeSBarry Smith extern PetscErrorCode MatCreateSeqDense(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*); 2457087cfbeSBarry Smith extern PetscErrorCode MatCreateMPIDense(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*); 2467087cfbeSBarry Smith extern PetscErrorCode MatCreateSeqAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 247ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,nz,nnz,&A),Mat,A) 248ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,&A),Mat,A) 249ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,0,nnz,&A),Mat,A) 250ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,&A),Mat,A) 251ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,A)) 252ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,m,n,0,nnz,A)) 253ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,A)) 2547087cfbeSBarry Smith extern PetscErrorCode MatCreateMPIAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 255ba966639SSatish 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) 256ba966639SSatish 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) 257ba966639SSatish 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) 258ba966639SSatish 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) 259ba966639SSatish 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)) 260ba966639SSatish 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)) 261ba966639SSatish 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)) 262ba966639SSatish 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) 263ba966639SSatish 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) 264ba966639SSatish 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) 265ba966639SSatish 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) 266ba966639SSatish 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)) 267ba966639SSatish 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)) 268ba966639SSatish 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)) 2697087cfbeSBarry Smith extern PetscErrorCode MatCreateMPIAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *); 2707087cfbeSBarry Smith extern PetscErrorCode MatCreateMPIAIJWithSplitArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],PetscInt[],PetscInt[],PetscScalar[],Mat*); 2718d7a6e47SBarry Smith 2727087cfbeSBarry Smith extern PetscErrorCode MatCreateSeqBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 273ba966639SSatish 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) 274ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A) 275ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A) 276ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A) 277ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A)) 278ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A)) 279ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A)) 2807087cfbeSBarry Smith extern PetscErrorCode MatCreateMPIBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 281ba966639SSatish 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) 282ba966639SSatish 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) 283ba966639SSatish 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) 284ba966639SSatish 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) 285ba966639SSatish 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)) 286ba966639SSatish 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)) 287ba966639SSatish 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)) 288ba966639SSatish 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) 289ba966639SSatish 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) 290ba966639SSatish 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) 291ba966639SSatish 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) 292ba966639SSatish 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)) 293ba966639SSatish 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)) 294ba966639SSatish 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)) 2957087cfbeSBarry Smith extern PetscErrorCode MatCreateMPIBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat*); 296d21a29f3SJed Brown 2977087cfbeSBarry Smith extern PetscErrorCode MatCreateMPIAdj(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscInt[],Mat*); 2987087cfbeSBarry Smith extern PetscErrorCode MatCreateSeqSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 299ba966639SSatish 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) 300ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A) 301ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A) 302ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A) 303ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A)) 304ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A)) 305ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A)) 306d21a29f3SJed Brown 3077087cfbeSBarry Smith extern PetscErrorCode MatCreateMPISBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 308ba966639SSatish 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) 309ba966639SSatish 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) 310ba966639SSatish 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) 311ba966639SSatish 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) 312ba966639SSatish 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)) 313ba966639SSatish 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)) 314ba966639SSatish 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)) 315ba966639SSatish 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) 316ba966639SSatish 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) 317ba966639SSatish 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) 318ba966639SSatish 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) 319ba966639SSatish 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)) 320ba966639SSatish 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)) 321ba966639SSatish 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)) 3227087cfbeSBarry Smith extern PetscErrorCode MatCreateMPISBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *); 3237087cfbeSBarry Smith extern PetscErrorCode MatMPISBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 324dfb205c3SBarry Smith 3257087cfbeSBarry Smith extern PetscErrorCode MatCreateShell(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,void *,Mat*); 326ba966639SSatish 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) 327ba966639SSatish Balay PetscPolymorphicFunction(MatCreateShell,(PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx),(PETSC_COMM_WORLD,m,n,M,N,ctx,&A),Mat,A) 3287087cfbeSBarry Smith extern PetscErrorCode MatCreateAdic(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,void (*)(void),Mat*); 3297087cfbeSBarry Smith extern PetscErrorCode MatCreateNormal(Mat,Mat*); 330ba966639SSatish Balay PetscPolymorphicFunction(MatCreateNormal,(Mat mat),(mat,&A),Mat,A) 3317087cfbeSBarry Smith extern PetscErrorCode MatCreateLRC(Mat,Mat,Mat,Mat*); 3327087cfbeSBarry Smith extern PetscErrorCode MatCreateIS(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,ISLocalToGlobalMapping,Mat*); 3337087cfbeSBarry Smith extern PetscErrorCode MatCreateSeqAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 3347087cfbeSBarry Smith extern PetscErrorCode MatCreateMPIAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 335c0cdd4a1SDahai Guo 336c0cdd4a1SDahai Guo extern PetscErrorCode MatCreateSeqBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 337c0cdd4a1SDahai Guo extern PetscErrorCode MatCreateMPIBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 338aa5a9175SDahai Guo extern PetscErrorCode MatCreateSeqSBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 339aa5a9175SDahai Guo extern PetscErrorCode MatCreateMPISBSTRM(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 340c0cdd4a1SDahai Guo 3417087cfbeSBarry Smith extern PetscErrorCode MatCreateScatter(MPI_Comm,VecScatter,Mat*); 3427087cfbeSBarry Smith extern PetscErrorCode MatScatterSetVecScatter(Mat,VecScatter); 3437087cfbeSBarry Smith extern PetscErrorCode MatScatterGetVecScatter(Mat,VecScatter*); 3447087cfbeSBarry Smith extern PetscErrorCode MatCreateBlockMat(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,Mat*); 3457087cfbeSBarry Smith extern PetscErrorCode MatCompositeAddMat(Mat,Mat); 3467087cfbeSBarry Smith extern PetscErrorCode MatCompositeMerge(Mat); 3477087cfbeSBarry Smith extern PetscErrorCode MatCreateComposite(MPI_Comm,PetscInt,const Mat*,Mat*); 3486d7c1e57SBarry Smith typedef enum {MAT_COMPOSITE_ADDITIVE,MAT_COMPOSITE_MULTIPLICATIVE} MatCompositeType; 3497087cfbeSBarry Smith extern PetscErrorCode MatCompositeSetType(Mat,MatCompositeType); 3506d7c1e57SBarry Smith 351dedccee8SHong Zhang extern PetscErrorCode MatCreateFFT(MPI_Comm,PetscInt,const PetscInt[],const MatType,Mat*); 3527087cfbeSBarry Smith extern PetscErrorCode MatCreateSeqCUFFT(MPI_Comm,PetscInt,const PetscInt[],Mat*); 353dedccee8SHong Zhang 3547087cfbeSBarry Smith extern PetscErrorCode MatCreateTranspose(Mat,Mat*); 3557087cfbeSBarry Smith extern PetscErrorCode MatCreateSubMatrix(Mat,IS,IS,Mat*); 3567087cfbeSBarry Smith extern PetscErrorCode MatSubMatrixUpdate(Mat,Mat,IS,IS); 3577087cfbeSBarry Smith extern PetscErrorCode MatCreateLocalRef(Mat,IS,IS,Mat*); 3581472f72bSBarry Smith 3597087cfbeSBarry Smith extern PetscErrorCode MatPythonSetType(Mat,const char[]); 3601d6018f0SLisandro Dalcin 3617087cfbeSBarry Smith extern PetscErrorCode MatSetUp(Mat); 362fcfd50ebSBarry Smith extern PetscErrorCode MatDestroy(Mat*); 36321c89e3eSBarry Smith 3647087cfbeSBarry Smith extern PetscErrorCode MatConjugate(Mat); 3657087cfbeSBarry Smith extern PetscErrorCode MatRealPart(Mat); 3667087cfbeSBarry Smith extern PetscErrorCode MatImaginaryPart(Mat); 36711bd1e4dSLisandro Dalcin extern PetscErrorCode MatGetDiagonalBlock(Mat,Mat*); 3687087cfbeSBarry Smith extern PetscErrorCode MatGetTrace(Mat,PetscScalar*); 369bbead8a2SBarry Smith extern PetscErrorCode MatInvertBlockDiagonal(Mat,PetscScalar **); 37099cafbc1SBarry Smith 3718ed539a5SBarry Smith /* ------------------------------------------------------------*/ 3727087cfbeSBarry Smith extern PetscErrorCode MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 3737087cfbeSBarry Smith extern PetscErrorCode MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 3747087cfbeSBarry Smith extern PetscErrorCode MatSetValuesRow(Mat,PetscInt,const PetscScalar[]); 3757087cfbeSBarry Smith extern PetscErrorCode MatSetValuesRowLocal(Mat,PetscInt,const PetscScalar[]); 376*2f35a324SMatthew G Knepley extern PetscErrorCode MatSetValuesBatch(Mat,PetscInt,PetscInt,PetscInt[],const PetscScalar[]); 37784cb2905SBarry Smith 3782ef4de8bSBarry Smith /*S 3792ef4de8bSBarry Smith MatStencil - Data structure (C struct) for storing information about a single row or 3802ef4de8bSBarry Smith column of a matrix as index on an associated grid. 3812ef4de8bSBarry Smith 3822ef4de8bSBarry Smith Level: beginner 3832ef4de8bSBarry Smith 3842ef4de8bSBarry Smith Concepts: matrix; linear operator 3852ef4de8bSBarry Smith 386d7bd93baSBarry Smith .seealso: MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockStencil() 3872ef4de8bSBarry Smith S*/ 388435da068SBarry Smith typedef struct { 389c1ac3661SBarry Smith PetscInt k,j,i,c; 390435da068SBarry Smith } MatStencil; 3912ef4de8bSBarry Smith 3927087cfbeSBarry Smith extern PetscErrorCode MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode); 3937087cfbeSBarry Smith extern PetscErrorCode MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode); 3947087cfbeSBarry Smith extern PetscErrorCode MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt); 395435da068SBarry Smith 3967087cfbeSBarry Smith extern PetscErrorCode MatSetColoring(Mat,ISColoring); 3977087cfbeSBarry Smith extern PetscErrorCode MatSetValuesAdic(Mat,void*); 3987087cfbeSBarry Smith extern PetscErrorCode MatSetValuesAdifor(Mat,PetscInt,void*); 3993a7fca6bSBarry Smith 400d91e6319SBarry Smith /*E 401d91e6319SBarry Smith MatAssemblyType - Indicates if the matrix is now to be used, or if you plan 402d91e6319SBarry Smith to continue to add values to it 403d91e6319SBarry Smith 404d91e6319SBarry Smith Level: beginner 405d91e6319SBarry Smith 406d91e6319SBarry Smith .seealso: MatAssemblyBegin(), MatAssemblyEnd() 407d91e6319SBarry Smith E*/ 4086d4a8577SBarry Smith typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType; 4097087cfbeSBarry Smith extern PetscErrorCode MatAssemblyBegin(Mat,MatAssemblyType); 4107087cfbeSBarry Smith extern PetscErrorCode MatAssemblyEnd(Mat,MatAssemblyType); 4117087cfbeSBarry Smith extern PetscErrorCode MatAssembled(Mat,PetscBool *); 4124f9c727eSBarry Smith 4131d73ed98SBarry Smith 41430de9b25SBarry Smith 415d91e6319SBarry Smith /*E 416d91e6319SBarry Smith MatOption - Options that may be set for a matrix and its behavior or storage 417d91e6319SBarry Smith 418d91e6319SBarry Smith Level: beginner 419d91e6319SBarry Smith 4200a835dfdSSatish Balay Any additions/changes here MUST also be made in include/finclude/petscmat.h 421d91e6319SBarry Smith 422d91e6319SBarry Smith .seealso: MatSetOption() 423d91e6319SBarry Smith E*/ 424512a5fc5SBarry Smith typedef enum {MAT_ROW_ORIENTED,MAT_NEW_NONZERO_LOCATIONS, 4254e0d8c25SBarry Smith MAT_SYMMETRIC, 4264e0d8c25SBarry Smith MAT_STRUCTURALLY_SYMMETRIC, 4278047e77bSBarry Smith MAT_NEW_DIAGONALS,MAT_IGNORE_OFF_PROC_ENTRIES, 4284e0d8c25SBarry Smith MAT_NEW_NONZERO_LOCATION_ERR, 4294e0d8c25SBarry Smith MAT_NEW_NONZERO_ALLOCATION_ERR,MAT_USE_HASH_TABLE, 430a9817697SBarry Smith MAT_KEEP_NONZERO_PATTERN,MAT_IGNORE_ZERO_ENTRIES, 4314e0d8c25SBarry Smith MAT_USE_INODES, 4324e0d8c25SBarry Smith MAT_HERMITIAN, 4334e0d8c25SBarry Smith MAT_SYMMETRY_ETERNAL, 434cd6b891eSBarry Smith MAT_CHECK_COMPRESSED_ROW, 4354e0d8c25SBarry Smith MAT_IGNORE_LOWER_TRIANGULAR,MAT_ERROR_LOWER_TRIANGULAR, 43628b2fa4aSMatthew Knepley MAT_GETROW_UPPERTRIANGULAR,MAT_UNUSED_NONZERO_LOCATION_ERR, 4374cb17eb5SBarry Smith MAT_SPD,MAT_NO_OFF_PROC_ENTRIES,MAT_NO_OFF_PROC_ZERO_ROWS, 43828b2fa4aSMatthew Knepley NUM_MAT_OPTIONS} MatOption; 439290bbb0aSBarry Smith extern const char *MatOptions[]; 4407087cfbeSBarry Smith extern PetscErrorCode MatSetOption(Mat,MatOption,PetscBool ); 4417087cfbeSBarry Smith extern PetscErrorCode MatGetType(Mat,const MatType*); 442a313700dSBarry Smith PetscPolymorphicFunction(MatGetType,(Mat mat),(mat,&t),const MatType,t) 44384cb2905SBarry Smith 4447087cfbeSBarry Smith extern PetscErrorCode MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]); 4457087cfbeSBarry Smith extern PetscErrorCode MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 4467087cfbeSBarry Smith extern PetscErrorCode MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 4477087cfbeSBarry Smith extern PetscErrorCode MatGetRowUpperTriangular(Mat); 4487087cfbeSBarry Smith extern PetscErrorCode MatRestoreRowUpperTriangular(Mat); 4497087cfbeSBarry Smith extern PetscErrorCode MatGetColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 4507087cfbeSBarry Smith extern PetscErrorCode MatRestoreColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 4517087cfbeSBarry Smith extern PetscErrorCode MatGetColumnVector(Mat,Vec,PetscInt); 4527087cfbeSBarry Smith extern PetscErrorCode MatGetArray(Mat,PetscScalar *[]); 453ba966639SSatish Balay PetscPolymorphicFunction(MatGetArray,(Mat mat),(mat,&a),PetscScalar*,a) 4547087cfbeSBarry Smith extern PetscErrorCode MatRestoreArray(Mat,PetscScalar *[]); 4557087cfbeSBarry Smith extern PetscErrorCode MatGetBlockSize(Mat,PetscInt *); 456ba966639SSatish Balay PetscPolymorphicFunction(MatGetBlockSize,(Mat mat),(mat,&a),PetscInt,a) 4577087cfbeSBarry Smith extern PetscErrorCode MatSetBlockSize(Mat,PetscInt); 4587b80b807SBarry Smith 4591620fd73SBarry Smith 4607087cfbeSBarry Smith extern PetscErrorCode MatMult(Mat,Vec,Vec); 4617087cfbeSBarry Smith extern PetscErrorCode MatMultDiagonalBlock(Mat,Vec,Vec); 4627087cfbeSBarry Smith extern PetscErrorCode MatMultAdd(Mat,Vec,Vec,Vec); 463ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultAdd,(Mat A,Vec x,Vec y),(A,x,y,y)) 4647087cfbeSBarry Smith extern PetscErrorCode MatMultTranspose(Mat,Vec,Vec); 4657087cfbeSBarry Smith extern PetscErrorCode MatMultHermitianTranspose(Mat,Vec,Vec); 4667087cfbeSBarry Smith extern PetscErrorCode MatIsTranspose(Mat,Mat,PetscReal,PetscBool *); 467ace3abfcSBarry Smith PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B,PetscReal tol),(A,B,tol,&t),PetscBool ,t) 468ace3abfcSBarry Smith PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B),(A,B,0,&t),PetscBool ,t) 4697087cfbeSBarry Smith extern PetscErrorCode MatIsHermitianTranspose(Mat,Mat,PetscReal,PetscBool *); 4707087cfbeSBarry Smith extern PetscErrorCode MatMultTransposeAdd(Mat,Vec,Vec,Vec); 4717087cfbeSBarry Smith extern PetscErrorCode MatMultHermitianTransposeAdd(Mat,Vec,Vec,Vec); 472ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultTransposeAdd,(Mat A,Vec x,Vec y),(A,x,y,y)) 4737087cfbeSBarry Smith extern PetscErrorCode MatMultConstrained(Mat,Vec,Vec); 4747087cfbeSBarry Smith extern PetscErrorCode MatMultTransposeConstrained(Mat,Vec,Vec); 4757087cfbeSBarry Smith extern PetscErrorCode MatMatSolve(Mat,Mat,Mat); 476f5edf698SHong Zhang 477d91e6319SBarry Smith /*E 478d91e6319SBarry Smith MatDuplicateOption - Indicates if a duplicated sparse matrix should have 479d91e6319SBarry Smith its numerical values copied over or just its nonzero structure. 480d91e6319SBarry Smith 481d91e6319SBarry Smith Level: beginner 482d91e6319SBarry Smith 483d91e6319SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 484d91e6319SBarry Smith 48570dcbbb9SBarry Smith $ MAT_SHARE_NONZERO_PATTERN - the i and j arrays in the new matrix will be shared with the original matrix 48670dcbbb9SBarry Smith $ this also triggers the MAT_DO_NOT_COPY_VALUES option. This is used when you 48770dcbbb9SBarry Smith $ have several matrices with the same nonzero pattern. 48870dcbbb9SBarry Smith 489d91e6319SBarry Smith .seealso: MatDuplicate() 490d91e6319SBarry Smith E*/ 49170dcbbb9SBarry Smith typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES,MAT_SHARE_NONZERO_PATTERN} MatDuplicateOption; 4922e8a6d31SBarry Smith 4937087cfbeSBarry Smith extern PetscErrorCode MatConvert(Mat,const MatType,MatReuse,Mat*); 494a313700dSBarry Smith PetscPolymorphicFunction(MatConvert,(Mat A,const MatType t),(A,t,MAT_INITIAL_MATRIX,&a),Mat,a) 4957087cfbeSBarry Smith extern PetscErrorCode MatDuplicate(Mat,MatDuplicateOption,Mat*); 496ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A,MatDuplicateOption o),(A,o,&a),Mat,a) 497ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A),(A,MAT_COPY_VALUES,&a),Mat,a) 49894a9d846SBarry Smith 499cb5b572fSBarry Smith 5007087cfbeSBarry Smith extern PetscErrorCode MatCopy(Mat,Mat,MatStructure); 5017087cfbeSBarry Smith extern PetscErrorCode MatView(Mat,PetscViewer); 5027087cfbeSBarry Smith extern PetscErrorCode MatIsSymmetric(Mat,PetscReal,PetscBool *); 503ace3abfcSBarry Smith PetscPolymorphicFunction(MatIsSymmetric,(Mat A,PetscReal tol),(A,tol,&t),PetscBool ,t) 504ace3abfcSBarry Smith PetscPolymorphicFunction(MatIsSymmetric,(Mat A),(A,0,&t),PetscBool ,t) 5057087cfbeSBarry Smith extern PetscErrorCode MatIsStructurallySymmetric(Mat,PetscBool *); 506ace3abfcSBarry Smith PetscPolymorphicFunction(MatIsStructurallySymmetric,(Mat A),(A,&t),PetscBool ,t) 5077087cfbeSBarry Smith extern PetscErrorCode MatIsHermitian(Mat,PetscReal,PetscBool *); 508ace3abfcSBarry Smith PetscPolymorphicFunction(MatIsHermitian,(Mat A),(A,0,&t),PetscBool ,t) 5097087cfbeSBarry Smith extern PetscErrorCode MatIsSymmetricKnown(Mat,PetscBool *,PetscBool *); 5107087cfbeSBarry Smith extern PetscErrorCode MatIsHermitianKnown(Mat,PetscBool *,PetscBool *); 5117087cfbeSBarry Smith extern PetscErrorCode MatMissingDiagonal(Mat,PetscBool *,PetscInt *); 5127087cfbeSBarry Smith extern PetscErrorCode MatLoad(Mat, PetscViewer); 5137b80b807SBarry Smith 5147087cfbeSBarry Smith extern PetscErrorCode MatGetRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,PetscInt *[],PetscInt *[],PetscBool *); 5157087cfbeSBarry Smith extern PetscErrorCode MatRestoreRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,PetscInt *[],PetscInt *[],PetscBool *); 5167087cfbeSBarry Smith extern PetscErrorCode MatGetColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,PetscInt *[],PetscInt *[],PetscBool *); 5177087cfbeSBarry Smith extern PetscErrorCode MatRestoreColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,PetscInt *[],PetscInt *[],PetscBool *); 518d4fbbf0eSBarry Smith 519d91e6319SBarry Smith /*S 520d91e6319SBarry Smith MatInfo - Context of matrix information, used with MatGetInfo() 521d91e6319SBarry Smith 522d91e6319SBarry Smith In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE 523d91e6319SBarry Smith 524d91e6319SBarry Smith Level: intermediate 525d91e6319SBarry Smith 526d91e6319SBarry Smith Concepts: matrix^nonzero information 527d91e6319SBarry Smith 528d9274352SBarry Smith .seealso: MatGetInfo(), MatInfoType 529d91e6319SBarry Smith S*/ 5304e220ebcSLois Curfman McInnes typedef struct { 531b0a32e0cSBarry Smith PetscLogDouble block_size; /* block size */ 532b0a32e0cSBarry Smith PetscLogDouble nz_allocated,nz_used,nz_unneeded; /* number of nonzeros */ 533b0a32e0cSBarry Smith PetscLogDouble memory; /* memory allocated */ 534b0a32e0cSBarry Smith PetscLogDouble assemblies; /* number of matrix assemblies called */ 535b0a32e0cSBarry Smith PetscLogDouble mallocs; /* number of mallocs during MatSetValues() */ 536b0a32e0cSBarry Smith PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */ 537b0a32e0cSBarry Smith PetscLogDouble factor_mallocs; /* number of mallocs during factorization */ 5384e220ebcSLois Curfman McInnes } MatInfo; 5394e220ebcSLois Curfman McInnes 540d9274352SBarry Smith /*E 541d9274352SBarry Smith MatInfoType - Indicates if you want information about the local part of the matrix, 542d9274352SBarry Smith the entire parallel matrix or the maximum over all the local parts. 543d9274352SBarry Smith 544d9274352SBarry Smith Level: beginner 545d9274352SBarry Smith 546d9274352SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 547d9274352SBarry Smith 548d9274352SBarry Smith .seealso: MatGetInfo(), MatInfo 549d9274352SBarry Smith E*/ 5507b80b807SBarry Smith typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType; 5517087cfbeSBarry Smith extern PetscErrorCode MatGetInfo(Mat,MatInfoType,MatInfo*); 5527087cfbeSBarry Smith extern PetscErrorCode MatGetDiagonal(Mat,Vec); 5537087cfbeSBarry Smith extern PetscErrorCode MatGetRowMax(Mat,Vec,PetscInt[]); 5547087cfbeSBarry Smith extern PetscErrorCode MatGetRowMin(Mat,Vec,PetscInt[]); 5557087cfbeSBarry Smith extern PetscErrorCode MatGetRowMaxAbs(Mat,Vec,PetscInt[]); 5567087cfbeSBarry Smith extern PetscErrorCode MatGetRowMinAbs(Mat,Vec,PetscInt[]); 5577087cfbeSBarry Smith extern PetscErrorCode MatGetRowSum(Mat,Vec); 5587087cfbeSBarry Smith extern PetscErrorCode MatTranspose(Mat,MatReuse,Mat*); 559fc4dec0aSBarry Smith PetscPolymorphicFunction(MatTranspose,(Mat A),(A,MAT_INITIAL_MATRIX,&t),Mat,t) 5607087cfbeSBarry Smith extern PetscErrorCode MatHermitianTranspose(Mat,MatReuse,Mat*); 5617087cfbeSBarry Smith extern PetscErrorCode MatPermute(Mat,IS,IS,Mat *); 562ba966639SSatish Balay PetscPolymorphicFunction(MatPermute,(Mat A,IS is1,IS is2),(A,is1,is2,&t),Mat,t) 5637087cfbeSBarry Smith extern PetscErrorCode MatDiagonalScale(Mat,Vec,Vec); 5647087cfbeSBarry Smith extern PetscErrorCode MatDiagonalSet(Mat,Vec,InsertMode); 5657087cfbeSBarry Smith extern PetscErrorCode MatEqual(Mat,Mat,PetscBool *); 566ace3abfcSBarry Smith PetscPolymorphicFunction(MatEqual,(Mat A,Mat B),(A,B,&t),PetscBool ,t) 5677087cfbeSBarry Smith extern PetscErrorCode MatMultEqual(Mat,Mat,PetscInt,PetscBool *); 5687087cfbeSBarry Smith extern PetscErrorCode MatMultAddEqual(Mat,Mat,PetscInt,PetscBool *); 5697087cfbeSBarry Smith extern PetscErrorCode MatMultTransposeEqual(Mat,Mat,PetscInt,PetscBool *); 5707087cfbeSBarry Smith extern PetscErrorCode MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscBool *); 5717b80b807SBarry Smith 5727087cfbeSBarry Smith extern PetscErrorCode MatNorm(Mat,NormType,PetscReal *); 573ba966639SSatish Balay PetscPolymorphicFunction(MatNorm,(Mat A,NormType t),(A,t,&n),PetscReal,n) 57409573ac7SBarry Smith extern PetscErrorCode MatGetColumnNorms(Mat,NormType,PetscReal *); 5757087cfbeSBarry Smith extern PetscErrorCode MatZeroEntries(Mat); 5767087cfbeSBarry Smith extern PetscErrorCode MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 5777087cfbeSBarry Smith extern PetscErrorCode MatZeroRowsIS(Mat,IS,PetscScalar,Vec,Vec); 5787087cfbeSBarry Smith extern PetscErrorCode MatZeroRowsStencil(Mat,PetscInt,const MatStencil [],PetscScalar,Vec,Vec); 5797087cfbeSBarry Smith extern PetscErrorCode MatZeroRowsColumns(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 5807087cfbeSBarry Smith extern PetscErrorCode MatZeroRowsColumnsIS(Mat,IS,PetscScalar,Vec,Vec); 5817b80b807SBarry Smith 5827087cfbeSBarry Smith extern PetscErrorCode MatUseScaledForm(Mat,PetscBool ); 5837087cfbeSBarry Smith extern PetscErrorCode MatScaleSystem(Mat,Vec,Vec); 5847087cfbeSBarry Smith extern PetscErrorCode MatUnScaleSystem(Mat,Vec,Vec); 5855ef9f2a5SBarry Smith 5867087cfbeSBarry Smith extern PetscErrorCode MatGetSize(Mat,PetscInt*,PetscInt*); 5877087cfbeSBarry Smith extern PetscErrorCode MatGetLocalSize(Mat,PetscInt*,PetscInt*); 5887087cfbeSBarry Smith extern PetscErrorCode MatGetOwnershipRange(Mat,PetscInt*,PetscInt*); 5897087cfbeSBarry Smith extern PetscErrorCode MatGetOwnershipRanges(Mat,const PetscInt**); 5907087cfbeSBarry Smith extern PetscErrorCode MatGetOwnershipRangeColumn(Mat,PetscInt*,PetscInt*); 5917087cfbeSBarry Smith extern PetscErrorCode MatGetOwnershipRangesColumn(Mat,const PetscInt**); 5927b80b807SBarry Smith 5937087cfbeSBarry Smith extern PetscErrorCode MatGetSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]); 5947087cfbeSBarry Smith extern PetscErrorCode MatGetSubMatricesParallel(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]); 5957087cfbeSBarry Smith extern PetscErrorCode MatDestroyMatrices(PetscInt,Mat *[]); 5967087cfbeSBarry Smith extern PetscErrorCode MatGetSubMatrix(Mat,IS,IS,MatReuse,Mat *); 5977087cfbeSBarry Smith extern PetscErrorCode MatGetLocalSubMatrix(Mat,IS,IS,Mat*); 5987087cfbeSBarry Smith extern PetscErrorCode MatRestoreLocalSubMatrix(Mat,IS,IS,Mat*); 5997087cfbeSBarry Smith extern PetscErrorCode MatGetSeqNonzeroStructure(Mat,Mat*); 6007087cfbeSBarry Smith extern PetscErrorCode MatDestroySeqNonzeroStructure(Mat*); 6015494a064SHong Zhang 6027087cfbeSBarry Smith extern PetscErrorCode MatMerge(MPI_Comm,Mat,PetscInt,MatReuse,Mat*); 6037087cfbeSBarry Smith extern PetscErrorCode MatMerge_SeqsToMPI(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*); 6047087cfbeSBarry Smith extern PetscErrorCode MatMerge_SeqsToMPISymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*); 6057087cfbeSBarry Smith extern PetscErrorCode MatMerge_SeqsToMPINumeric(Mat,Mat); 6064a2b5492SBarry Smith extern PetscErrorCode MatMPIAIJGetLocalMat(Mat,MatReuse,Mat*); 6074a2b5492SBarry Smith extern PetscErrorCode MatMPIAIJGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*); 6087087cfbeSBarry Smith extern PetscErrorCode MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,PetscInt*,Mat*); 6097087cfbeSBarry Smith extern PetscErrorCode MatGetBrowsOfAoCols(Mat,Mat,MatReuse,PetscInt**,PetscInt**,MatScalar**,Mat*); 61043eb5e2fSMatthew Knepley #if defined (PETSC_USE_CTABLE) 611cd116e26SSatish Balay #include "petscctable.h" 6127087cfbeSBarry Smith extern PetscErrorCode MatGetCommunicationStructs(Mat, Vec *, PetscTable *, VecScatter *); 61343eb5e2fSMatthew Knepley #else 6147087cfbeSBarry Smith extern PetscErrorCode MatGetCommunicationStructs(Mat, Vec *, PetscInt *[], VecScatter *); 61543eb5e2fSMatthew Knepley #endif 6167087cfbeSBarry Smith extern PetscErrorCode MatGetGhosts(Mat, PetscInt *,const PetscInt *[]); 6178efafbd8SBarry Smith 6187087cfbeSBarry Smith extern PetscErrorCode MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt); 6197b80b807SBarry Smith 6207087cfbeSBarry Smith extern PetscErrorCode MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*); 6217087cfbeSBarry Smith extern PetscErrorCode MatMatMultSymbolic(Mat,Mat,PetscReal,Mat*); 6227087cfbeSBarry Smith extern PetscErrorCode MatMatMultNumeric(Mat,Mat,Mat); 62322440eb1SKris Buschelman 6247087cfbeSBarry Smith extern PetscErrorCode MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*); 6257087cfbeSBarry Smith extern PetscErrorCode MatPtAPSymbolic(Mat,Mat,PetscReal,Mat*); 6267087cfbeSBarry Smith extern PetscErrorCode MatPtAPNumeric(Mat,Mat,Mat); 62722440eb1SKris Buschelman 6287087cfbeSBarry Smith extern PetscErrorCode MatMatMultTranspose(Mat,Mat,MatReuse,PetscReal,Mat*); 6297087cfbeSBarry Smith extern PetscErrorCode MatMatMultTransposeSymbolic(Mat,Mat,PetscReal,Mat*); 6307087cfbeSBarry Smith extern PetscErrorCode MatMatMultTransposeNumeric(Mat,Mat,Mat); 631bc011b1eSHong Zhang 6327087cfbeSBarry Smith extern PetscErrorCode MatAXPY(Mat,PetscScalar,Mat,MatStructure); 6337087cfbeSBarry Smith extern PetscErrorCode MatAYPX(Mat,PetscScalar,Mat,MatStructure); 6341c741599SHong Zhang 6357087cfbeSBarry Smith extern PetscErrorCode MatScale(Mat,PetscScalar); 6367087cfbeSBarry Smith extern PetscErrorCode MatShift(Mat,PetscScalar); 6377b80b807SBarry Smith 6387087cfbeSBarry Smith extern PetscErrorCode MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping); 6397087cfbeSBarry Smith extern PetscErrorCode MatSetLocalToGlobalMappingBlock(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping); 6407087cfbeSBarry Smith extern PetscErrorCode MatGetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*); 6417087cfbeSBarry Smith extern PetscErrorCode MatGetLocalToGlobalMappingBlock(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*); 6427087cfbeSBarry Smith extern PetscErrorCode MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 6437087cfbeSBarry Smith extern PetscErrorCode MatZeroRowsLocalIS(Mat,IS,PetscScalar,Vec,Vec); 6447087cfbeSBarry Smith extern PetscErrorCode MatZeroRowsColumnsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 6457087cfbeSBarry Smith extern PetscErrorCode MatZeroRowsColumnsLocalIS(Mat,IS,PetscScalar,Vec,Vec); 6467087cfbeSBarry Smith extern PetscErrorCode MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 6477087cfbeSBarry Smith extern PetscErrorCode MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 648052efed2SBarry Smith 6497087cfbeSBarry Smith extern PetscErrorCode MatStashSetInitialSize(Mat,PetscInt,PetscInt); 6507087cfbeSBarry Smith extern PetscErrorCode MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*); 65190f02eecSBarry Smith 6527087cfbeSBarry Smith extern PetscErrorCode MatInterpolate(Mat,Vec,Vec); 6537087cfbeSBarry Smith extern PetscErrorCode MatInterpolateAdd(Mat,Vec,Vec,Vec); 654ba966639SSatish Balay PetscPolymorphicSubroutine(MatInterpolateAdd,(Mat A,Vec x,Vec y),(A,x,y,y)) 6557087cfbeSBarry Smith extern PetscErrorCode MatRestrict(Mat,Vec,Vec); 6567087cfbeSBarry Smith extern PetscErrorCode MatGetVecs(Mat,Vec*,Vec*); 6577087cfbeSBarry Smith extern PetscErrorCode MatGetRedundantMatrix(Mat,PetscInt,MPI_Comm,PetscInt,MatReuse,Mat*); 6587087cfbeSBarry Smith extern PetscErrorCode MatGetMultiProcBlock(Mat,MPI_Comm,Mat*); 659fc9bc008SSatish Balay extern PetscErrorCode MatFindZeroDiagonals(Mat,IS*); 660bd481603SBarry Smith 661bd481603SBarry Smith /*MC 6621620fd73SBarry Smith MatSetValue - Set a single entry into a matrix. 6631620fd73SBarry Smith 6641620fd73SBarry Smith Not collective 6651620fd73SBarry Smith 6661620fd73SBarry Smith Input Parameters: 6671620fd73SBarry Smith + m - the matrix 6681620fd73SBarry Smith . row - the row location of the entry 6691620fd73SBarry Smith . col - the column location of the entry 6701620fd73SBarry Smith . value - the value to insert 6711620fd73SBarry Smith - mode - either INSERT_VALUES or ADD_VALUES 6721620fd73SBarry Smith 6731620fd73SBarry Smith Notes: 6741620fd73SBarry Smith For efficiency one should use MatSetValues() and set several or many 6751620fd73SBarry Smith values simultaneously if possible. 6761620fd73SBarry Smith 6771620fd73SBarry Smith Level: beginner 6781620fd73SBarry Smith 6791620fd73SBarry Smith .seealso: MatSetValues(), MatSetValueLocal() 6801620fd73SBarry Smith M*/ 6811620fd73SBarry 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);} 6821620fd73SBarry Smith 6831620fd73SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatGetValue(Mat v,PetscInt i,PetscInt j,PetscScalar *va) {return MatGetValues(v,1,&i,1,&j,va);} 6841620fd73SBarry Smith 6851620fd73SBarry 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);} 6861620fd73SBarry Smith 6871620fd73SBarry Smith /*MC 6880d2a7ff2SSatish Balay MatPreallocateInitialize - Begins the block of code that will count the number of nonzeros per 689bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 690bd481603SBarry Smith 691bd481603SBarry Smith Synopsis: 692c1ac3661SBarry Smith PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz) 693bd481603SBarry Smith 694bd481603SBarry Smith Collective on MPI_Comm 695bd481603SBarry Smith 696bd481603SBarry Smith Input Parameters: 697bd481603SBarry Smith + comm - the communicator that will share the eventually allocated matrix 698859fcb39SBarry Smith . nrows - the number of LOCAL rows in the matrix 699859fcb39SBarry Smith - ncols - the number of LOCAL columns in the matrix 700bd481603SBarry Smith 701bd481603SBarry Smith Output Parameters: 702bd481603SBarry Smith + dnz - the array that will be passed to the matrix preallocation routines 703bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 704bd481603SBarry Smith 705bd481603SBarry Smith 706bd481603SBarry Smith Level: intermediate 707bd481603SBarry Smith 708bd481603SBarry Smith Notes: 7090598bfebSBarry Smith See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details. 710bd481603SBarry Smith 7111d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 712bd481603SBarry Smith 713bd481603SBarry Smith Use MatPreallocateInitializeSymmetric() for symmetric matrices (MPISBAIJ matrices) 714bd481603SBarry Smith 7151620fd73SBarry Smith This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize(). 7161620fd73SBarry Smith 717bd481603SBarry Smith Concepts: preallocation^Matrix 718bd481603SBarry Smith 719bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(), 720bd481603SBarry Smith MatPreallocateInitializeSymmetric(), MatPreallocateSymmetricSetLocal() 721bd481603SBarry Smith M*/ 722c4f061fbSSatish Balay #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \ 7237c922b88SBarry Smith { \ 724a89bc211SBarry Smith PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__start,__end; \ 725a89bc211SBarry Smith _4_ierr = PetscMalloc2(__nrows,PetscInt,&dnz,__nrows,PetscInt,&onz);CHKERRQ(_4_ierr); \ 726a89bc211SBarry Smith _4_ierr = PetscMemzero(dnz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\ 727a89bc211SBarry Smith _4_ierr = PetscMemzero(onz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\ 728c1ac3661SBarry Smith _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ctmp;\ 729a89bc211SBarry Smith _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows; 7307c922b88SBarry Smith 731bd481603SBarry Smith /*MC 7320d2a7ff2SSatish Balay MatPreallocateSymmetricInitialize - Begins the block of code that will count the number of nonzeros per 733bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 734bd481603SBarry Smith 735bd481603SBarry Smith Synopsis: 736c1ac3661SBarry Smith PetscErrorCode MatPreallocateSymmetricInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz) 737bd481603SBarry Smith 738bd481603SBarry Smith Collective on MPI_Comm 739bd481603SBarry Smith 740bd481603SBarry Smith Input Parameters: 741bd481603SBarry Smith + comm - the communicator that will share the eventually allocated matrix 742859fcb39SBarry Smith . nrows - the number of LOCAL rows in the matrix 743859fcb39SBarry Smith - ncols - the number of LOCAL columns in the matrix 744bd481603SBarry Smith 745bd481603SBarry Smith Output Parameters: 746bd481603SBarry Smith + dnz - the array that will be passed to the matrix preallocation routines 747bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 748bd481603SBarry Smith 749bd481603SBarry Smith 750bd481603SBarry Smith Level: intermediate 751bd481603SBarry Smith 752bd481603SBarry Smith Notes: 7530598bfebSBarry Smith See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details. 754bd481603SBarry Smith 7551d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 756bd481603SBarry Smith 7571620fd73SBarry Smith This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize(). 7581620fd73SBarry Smith 759bd481603SBarry Smith Concepts: preallocation^Matrix 760bd481603SBarry Smith 761bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(), 762bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal() 763bd481603SBarry Smith M*/ 764222b16d4SBarry Smith #define MatPreallocateSymmetricInitialize(comm,nrows,ncols,dnz,onz) 0; \ 765222b16d4SBarry Smith { \ 766a89bc211SBarry Smith PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ctmp = (ncols),__rstart,__end; \ 767a89bc211SBarry Smith _4_ierr = PetscMalloc2(__nrows,PetscInt,&dnz,__nrows,PetscInt,&onz);CHKERRQ(_4_ierr); \ 768a89bc211SBarry Smith _4_ierr = PetscMemzero(dnz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\ 769a89bc211SBarry Smith _4_ierr = PetscMemzero(onz,__nrows*sizeof(PetscInt));CHKERRQ(_4_ierr);\ 770c1ac3661SBarry Smith _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr);\ 771a89bc211SBarry Smith _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows; 772222b16d4SBarry Smith 773bd481603SBarry Smith /*MC 774bd481603SBarry Smith MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be 775bd481603SBarry Smith inserted using a local number of the rows and columns 776bd481603SBarry Smith 777bd481603SBarry Smith Synopsis: 778c1ac3661SBarry Smith PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 779bd481603SBarry Smith 780bd481603SBarry Smith Not Collective 781bd481603SBarry Smith 782bd481603SBarry Smith Input Parameters: 783784ac674SJed Brown + map - the row mapping from local numbering to global numbering 784bd481603SBarry Smith . nrows - the number of rows indicated 7851d73ed98SBarry Smith . rows - the indices of the rows 786784ac674SJed Brown . cmap - the column mapping from local to global numbering 787bd481603SBarry Smith . ncols - the number of columns in the matrix 788bd481603SBarry Smith . cols - the columns indicated 789bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 790bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 791bd481603SBarry Smith 792bd481603SBarry Smith 793bd481603SBarry Smith Level: intermediate 794bd481603SBarry Smith 795bd481603SBarry Smith Notes: 7960598bfebSBarry Smith See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details. 797bd481603SBarry Smith 7981d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 799bd481603SBarry Smith 800bd481603SBarry Smith Concepts: preallocation^Matrix 801bd481603SBarry Smith 802bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(), 803bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal() 804bd481603SBarry Smith M*/ 805784ac674SJed Brown #define MatPreallocateSetLocal(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \ 806c4f061fbSSatish Balay {\ 807c1ac3661SBarry Smith PetscInt __l;\ 808784ac674SJed Brown _4_ierr = ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\ 809784ac674SJed Brown _4_ierr = ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\ 810c4f061fbSSatish Balay for (__l=0;__l<nrows;__l++) {\ 811ef66eb69SBarry Smith _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 812c4f061fbSSatish Balay }\ 813c4f061fbSSatish Balay } 814c4f061fbSSatish Balay 815bd481603SBarry Smith /*MC 816bd481603SBarry Smith MatPreallocateSymmetricSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be 817bd481603SBarry Smith inserted using a local number of the rows and columns 818bd481603SBarry Smith 819bd481603SBarry Smith Synopsis: 820c1ac3661SBarry Smith PetscErrorCode MatPreallocateSymmetricSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 821bd481603SBarry Smith 822bd481603SBarry Smith Not Collective 823bd481603SBarry Smith 824bd481603SBarry Smith Input Parameters: 825bd481603SBarry Smith + map - the mapping between local numbering and global numbering 826bd481603SBarry Smith . nrows - the number of rows indicated 8271d73ed98SBarry Smith . rows - the indices of the rows 828bd481603SBarry Smith . ncols - the number of columns in the matrix 829bd481603SBarry Smith . cols - the columns indicated 830bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 831bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 832bd481603SBarry Smith 833bd481603SBarry Smith 834bd481603SBarry Smith Level: intermediate 835bd481603SBarry Smith 836bd481603SBarry Smith Notes: 8370598bfebSBarry Smith See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details. 838bd481603SBarry Smith 839bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 840bd481603SBarry Smith 841bd481603SBarry Smith Concepts: preallocation^Matrix 842bd481603SBarry Smith 843bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(), 844bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal() 845bd481603SBarry Smith M*/ 846d3d32019SBarry Smith #define MatPreallocateSymmetricSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\ 847d3d32019SBarry Smith {\ 848c1ac3661SBarry Smith PetscInt __l;\ 849d3d32019SBarry Smith _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\ 850d3d32019SBarry Smith _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\ 851d3d32019SBarry Smith for (__l=0;__l<nrows;__l++) {\ 852d3d32019SBarry Smith _4_ierr = MatPreallocateSymmetricSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 853d3d32019SBarry Smith }\ 854d3d32019SBarry Smith } 855d3d32019SBarry Smith 856bd481603SBarry Smith /*MC 857bd481603SBarry Smith MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be 858bd481603SBarry Smith inserted using a local number of the rows and columns 859bd481603SBarry Smith 860bd481603SBarry Smith Synopsis: 861c1ac3661SBarry Smith PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 862bd481603SBarry Smith 863bd481603SBarry Smith Not Collective 864bd481603SBarry Smith 865bd481603SBarry Smith Input Parameters: 86664075487SBarry Smith + row - the row 867bd481603SBarry Smith . ncols - the number of columns in the matrix 868a50f8bf6SHong Zhang - cols - the columns indicated 869a50f8bf6SHong Zhang 870a50f8bf6SHong Zhang Output Parameters: 871a50f8bf6SHong Zhang + dnz - the array that will be passed to the matrix preallocation routines 872bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 873bd481603SBarry Smith 874bd481603SBarry Smith 875bd481603SBarry Smith Level: intermediate 876bd481603SBarry Smith 877bd481603SBarry Smith Notes: 8780598bfebSBarry Smith See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details. 879bd481603SBarry Smith 880bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 881bd481603SBarry Smith 8821620fd73SBarry Smith This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize(). 8831620fd73SBarry Smith 884bd481603SBarry Smith Concepts: preallocation^Matrix 885bd481603SBarry Smith 886bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(), 887bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal() 888bd481603SBarry Smith M*/ 889c4f061fbSSatish Balay #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\ 890c1ac3661SBarry Smith { PetscInt __i; \ 891e32f2f54SBarry 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);\ 892e32f2f54SBarry 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);\ 8937c922b88SBarry Smith for (__i=0; __i<nc; __i++) {\ 89464075487SBarry Smith if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \ 8957cc688d7SBarry Smith else dnz[row - __rstart]++;\ 8967c922b88SBarry Smith }\ 8977c922b88SBarry Smith } 8987c922b88SBarry Smith 899bd481603SBarry Smith /*MC 900bd481603SBarry Smith MatPreallocateSymmetricSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be 901bd481603SBarry Smith inserted using a local number of the rows and columns 902bd481603SBarry Smith 903bd481603SBarry Smith Synopsis: 904c1ac3661SBarry Smith PetscErrorCode MatPreallocateSymmetricSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 905bd481603SBarry Smith 906bd481603SBarry Smith Not Collective 907bd481603SBarry Smith 908bd481603SBarry Smith Input Parameters: 909bd481603SBarry Smith + nrows - the number of rows indicated 9101d73ed98SBarry Smith . rows - the indices of the rows 911bd481603SBarry Smith . ncols - the number of columns in the matrix 912bd481603SBarry Smith . cols - the columns indicated 913bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 914bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 915bd481603SBarry Smith 916bd481603SBarry Smith 917bd481603SBarry Smith Level: intermediate 918bd481603SBarry Smith 919bd481603SBarry Smith Notes: 9200598bfebSBarry Smith See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details. 921bd481603SBarry Smith 922bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 923bd481603SBarry Smith 9241620fd73SBarry Smith This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize(). 9251620fd73SBarry Smith 926bd481603SBarry Smith Concepts: preallocation^Matrix 927bd481603SBarry Smith 928bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(), 929bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal() 930bd481603SBarry Smith M*/ 931d3d32019SBarry Smith #define MatPreallocateSymmetricSet(row,nc,cols,dnz,onz) 0;\ 932c1ac3661SBarry Smith { PetscInt __i; \ 933d3d32019SBarry Smith for (__i=0; __i<nc; __i++) {\ 934d3d32019SBarry Smith if (cols[__i] >= __end) onz[row - __rstart]++; \ 935d3d32019SBarry Smith else if (cols[__i] >= row) dnz[row - __rstart]++;\ 936d3d32019SBarry Smith }\ 937d3d32019SBarry Smith } 938d3d32019SBarry Smith 939bd481603SBarry Smith /*MC 94016371a99SBarry Smith MatPreallocateLocation - An alternative to MatPreallocationSet() that puts the nonzero locations into the matrix if it exists 94116371a99SBarry Smith 94216371a99SBarry Smith Synopsis: 94316371a99SBarry Smith PetscErrorCode MatPreallocateLocations(Mat A,PetscInt row,PetscInt ncols,PetscInt *cols,PetscInt *dnz,PetscInt *onz) 94416371a99SBarry Smith 94516371a99SBarry Smith Not Collective 94616371a99SBarry Smith 94716371a99SBarry Smith Input Parameters: 94816371a99SBarry Smith . A - matrix 94916371a99SBarry Smith . row - row where values exist (must be local to this process) 95016371a99SBarry Smith . ncols - number of columns 95116371a99SBarry Smith . cols - columns with nonzeros 95216371a99SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 95316371a99SBarry Smith - ozn - the other array passed to the matrix preallocation routines 95416371a99SBarry Smith 95516371a99SBarry Smith 95616371a99SBarry Smith Level: intermediate 95716371a99SBarry Smith 95816371a99SBarry Smith Notes: 9590598bfebSBarry Smith See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details. 96016371a99SBarry Smith 96116371a99SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 96216371a99SBarry Smith 96316371a99SBarry Smith This is a MACRO not a function because it uses a bunch of variables private to the MatPreallocation.... routines. 96416371a99SBarry Smith 96516371a99SBarry Smith Concepts: preallocation^Matrix 96616371a99SBarry Smith 96716371a99SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(), 96816371a99SBarry Smith MatPreallocateSymmetricInitialize(), MatPreallocateSymmetricSetLocal() 96916371a99SBarry Smith M*/ 97016371a99SBarry 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);} 97116371a99SBarry Smith 97216371a99SBarry Smith 97316371a99SBarry Smith /*MC 9740d2a7ff2SSatish Balay MatPreallocateFinalize - Ends the block of code that will count the number of nonzeros per 975bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 976bd481603SBarry Smith 977bd481603SBarry Smith Synopsis: 978c1ac3661SBarry Smith PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz) 979bd481603SBarry Smith 980bd481603SBarry Smith Collective on MPI_Comm 981bd481603SBarry Smith 982bd481603SBarry Smith Input Parameters: 98316371a99SBarry Smith + dnz - the array that was be passed to the matrix preallocation routines 984bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 985bd481603SBarry Smith 986bd481603SBarry Smith 987bd481603SBarry Smith Level: intermediate 988bd481603SBarry Smith 989bd481603SBarry Smith Notes: 9900598bfebSBarry Smith See the <A href="../../docs/manual.pdf#nameddest=ch_performance">Hints for Performance Improvment</A> chapter in the users manual for more details. 991bd481603SBarry Smith 992bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 993bd481603SBarry Smith 9941620fd73SBarry Smith This is a MACRO not a function because it closes the { started in MatPreallocateInitialize(). 9951620fd73SBarry Smith 996bd481603SBarry Smith Concepts: preallocation^Matrix 997bd481603SBarry Smith 998bd481603SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(), 999bd481603SBarry Smith MatPreallocateSymmetricInitialize(), MatPreallocateSymmetricSetLocal() 1000bd481603SBarry Smith M*/ 1001a89bc211SBarry Smith #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree2(dnz,onz);CHKERRQ(_4_ierr);} 10027c922b88SBarry Smith 1003bd481603SBarry Smith 1004bd481603SBarry Smith 10057b80b807SBarry Smith /* Routines unique to particular data structures */ 10068fe8eb27SJed Brown extern PetscErrorCode MatShellGetContext(Mat,void *); 1007ba966639SSatish Balay PetscPolymorphicFunction(MatShellGetContext,(Mat A),(A,&t),void*,t) 1008698d4c6aSKris Buschelman 10097087cfbeSBarry Smith extern PetscErrorCode MatInodeAdjustForInodes(Mat,IS*,IS*); 10107087cfbeSBarry Smith extern PetscErrorCode MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *); 1011698d4c6aSKris Buschelman 10127087cfbeSBarry Smith extern PetscErrorCode MatSeqAIJSetColumnIndices(Mat,PetscInt[]); 10137087cfbeSBarry Smith extern PetscErrorCode MatSeqBAIJSetColumnIndices(Mat,PetscInt[]); 10147087cfbeSBarry Smith extern PetscErrorCode MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 10157087cfbeSBarry Smith extern PetscErrorCode MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 10167087cfbeSBarry Smith extern PetscErrorCode MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 10177b80b807SBarry Smith 1018a96a251dSBarry Smith #define MAT_SKIP_ALLOCATION -4 1019a96a251dSBarry Smith 10207087cfbeSBarry Smith extern PetscErrorCode MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]); 1021ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz)) 10227087cfbeSBarry Smith extern PetscErrorCode MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]); 1023ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqSBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz)) 10247087cfbeSBarry Smith extern PetscErrorCode MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]); 1025ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqAIJSetPreallocation,(Mat A,const PetscInt nnz[]),(A,0,nnz)) 1026273d9f13SBarry Smith 10277087cfbeSBarry Smith extern PetscErrorCode MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 1028ba966639SSatish Balay PetscPolymorphicSubroutine(MatMPIBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[],const PetscInt onz[]),(A,bs,0,nnz,0,onz)) 10297087cfbeSBarry Smith extern PetscErrorCode MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 10307087cfbeSBarry Smith extern PetscErrorCode MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 10317087cfbeSBarry Smith extern PetscErrorCode MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []); 10327087cfbeSBarry Smith extern PetscErrorCode MatSeqBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 10337087cfbeSBarry Smith extern PetscErrorCode MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]); 10347087cfbeSBarry Smith extern PetscErrorCode MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 10357087cfbeSBarry Smith extern PetscErrorCode MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]); 10367087cfbeSBarry Smith extern PetscErrorCode MatMPIDenseSetPreallocation(Mat,PetscScalar[]); 10377087cfbeSBarry Smith extern PetscErrorCode MatSeqDenseSetPreallocation(Mat,PetscScalar[]); 10387087cfbeSBarry Smith extern PetscErrorCode MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,PetscInt*[]); 10397087cfbeSBarry Smith extern PetscErrorCode MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,PetscInt*[]); 10407087cfbeSBarry Smith extern PetscErrorCode MatAdicSetLocalFunction(Mat,void (*)(void)); 1041273d9f13SBarry Smith 10427087cfbeSBarry Smith extern PetscErrorCode MatSeqDenseSetLDA(Mat,PetscInt); 10437087cfbeSBarry Smith extern PetscErrorCode MatDenseGetLocalMatrix(Mat,Mat*); 10441b807ce4Svictorle 10457087cfbeSBarry Smith extern PetscErrorCode MatStoreValues(Mat); 10467087cfbeSBarry Smith extern PetscErrorCode MatRetrieveValues(Mat); 10472e8a6d31SBarry Smith 10487087cfbeSBarry Smith extern PetscErrorCode MatDAADSetCtx(Mat,void*); 10493a7fca6bSBarry Smith 1050b3a44c85SBarry Smith extern PetscErrorCode MatFindNonzeroRows(Mat,IS*); 10517b80b807SBarry Smith /* 10527b80b807SBarry Smith These routines are not usually accessed directly, rather solving is 105394b7f48cSBarry Smith done through the KSP and PC interfaces. 10547b80b807SBarry Smith */ 10557b80b807SBarry Smith 1056d9274352SBarry Smith /*E 1057d9274352SBarry Smith MatOrderingType - String with the name of a PETSc matrix ordering or the creation function 1058d9274352SBarry Smith with an optional dynamic library name, for example 1059d9274352SBarry Smith http://www.mcs.anl.gov/petsc/lib.a:orderingcreate() 1060d9274352SBarry Smith 1061d9274352SBarry Smith Level: beginner 1062d9274352SBarry Smith 1063e5a9bf91SBarry Smith Cannot use const because the PC objects manipulate the string 1064e5a9bf91SBarry Smith 1065d9274352SBarry Smith .seealso: MatGetOrdering() 1066d9274352SBarry Smith E*/ 10673eaad2caSSatish Balay #define MatOrderingType char* 10682692d6eeSBarry Smith #define MATORDERINGNATURAL "natural" 10692692d6eeSBarry Smith #define MATORDERINGND "nd" 10702692d6eeSBarry Smith #define MATORDERING1WD "1wd" 10712692d6eeSBarry Smith #define MATORDERINGRCM "rcm" 10722692d6eeSBarry Smith #define MATORDERINGQMD "qmd" 10732692d6eeSBarry Smith #define MATORDERINGROWLENGTH "rowlength" 10742692d6eeSBarry Smith #define MATORDERINGDSC_ND "dsc_nd" /* these three are only for DSCPACK, see its documentation for details */ 10752692d6eeSBarry Smith #define MATORDERINGDSC_MMD "dsc_mmd" 10762692d6eeSBarry Smith #define MATORDERINGDSC_MDF "dsc_mdf" 1077312e372aSBarry Smith #define MATORDERINGAMD "amd" /* only works if UMFPACK is installed with PETSc */ 1078b12f92e5SBarry Smith 10797087cfbeSBarry Smith extern PetscErrorCode MatGetOrdering(Mat,const MatOrderingType,IS*,IS*); 10807087cfbeSBarry Smith extern PetscErrorCode MatGetOrderingList(PetscFList *list); 10817087cfbeSBarry Smith extern PetscErrorCode MatOrderingRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,const MatOrderingType,IS*,IS*)); 108230de9b25SBarry Smith 108330de9b25SBarry Smith /*MC 10841890ba74SBarry Smith MatOrderingRegisterDynamic - Adds a new sparse matrix ordering to the matrix package. 108530de9b25SBarry Smith 108630de9b25SBarry Smith Synopsis: 10871890ba74SBarry Smith PetscErrorCode MatOrderingRegisterDynamic(const char *name_ordering,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatOrdering)) 108830de9b25SBarry Smith 108930de9b25SBarry Smith Not Collective 109030de9b25SBarry Smith 109130de9b25SBarry Smith Input Parameters: 10922692d6eeSBarry Smith + sname - name of ordering (for example MATORDERINGND) 109330de9b25SBarry Smith . path - location of library where creation routine is 109430de9b25SBarry Smith . name - name of function that creates the ordering type,a string 109530de9b25SBarry Smith - function - function pointer that creates the ordering 109630de9b25SBarry Smith 109730de9b25SBarry Smith Level: developer 109830de9b25SBarry Smith 109930de9b25SBarry Smith If dynamic libraries are used, then the fourth input argument (function) 110030de9b25SBarry Smith is ignored. 110130de9b25SBarry Smith 110230de9b25SBarry Smith Sample usage: 110330de9b25SBarry Smith .vb 110430de9b25SBarry Smith MatOrderingRegisterDynamic("my_order",/home/username/my_lib/lib/libO/solaris/mylib.a, 110530de9b25SBarry Smith "MyOrder",MyOrder); 110630de9b25SBarry Smith .ve 110730de9b25SBarry Smith 110830de9b25SBarry Smith Then, your partitioner can be chosen with the procedural interface via 110930de9b25SBarry Smith $ MatOrderingSetType(part,"my_order) 111030de9b25SBarry Smith or at runtime via the option 11112401956bSBarry Smith $ -pc_factor_mat_ordering_type my_order 111230de9b25SBarry Smith 1113ab901514SBarry Smith ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values. 111430de9b25SBarry Smith 111530de9b25SBarry Smith .keywords: matrix, ordering, register 111630de9b25SBarry Smith 111730de9b25SBarry Smith .seealso: MatOrderingRegisterDestroy(), MatOrderingRegisterAll() 111830de9b25SBarry Smith M*/ 1119aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 1120f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,0) 1121b12f92e5SBarry Smith #else 1122f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,d) 1123b12f92e5SBarry Smith #endif 112430de9b25SBarry Smith 11257087cfbeSBarry Smith extern PetscErrorCode MatOrderingRegisterDestroy(void); 11267087cfbeSBarry Smith extern PetscErrorCode MatOrderingRegisterAll(const char[]); 1127ace3abfcSBarry Smith extern PetscBool MatOrderingRegisterAllCalled; 1128b0a32e0cSBarry Smith extern PetscFList MatOrderingList; 1129d4fbbf0eSBarry Smith 11307087cfbeSBarry Smith extern PetscErrorCode MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS); 1131a2ce50c7SBarry Smith 1132d91e6319SBarry Smith /*S 1133d90ac83dSHong Zhang MatFactorShiftType - Numeric Shift. 1134d90ac83dSHong Zhang 1135d90ac83dSHong Zhang Level: beginner 1136d90ac83dSHong Zhang 1137d90ac83dSHong Zhang S*/ 1138d90ac83dSHong Zhang typedef enum {MAT_SHIFT_NONE,MAT_SHIFT_NONZERO,MAT_SHIFT_POSITIVE_DEFINITE,MAT_SHIFT_INBLOCKS} MatFactorShiftType; 1139d90ac83dSHong Zhang extern const char *MatFactorShiftTypes[]; 1140d90ac83dSHong Zhang 1141d90ac83dSHong Zhang /*S 11422401956bSBarry Smith MatFactorInfo - Data passed into the matrix factorization routines 1143b00f7748SHong Zhang 114461cad860SBarry Smith In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE, that is use 114561cad860SBarry Smith $ MatFactorInfo info(MAT_FACTORINFO_SIZE) 1146b00f7748SHong Zhang 114715e8a5b3SHong Zhang Notes: These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC. 1148b00f7748SHong Zhang 114961cad860SBarry Smith You can use MatFactorInfoInitialize() to set default values. 115061cad860SBarry Smith 1151b00f7748SHong Zhang Level: developer 1152b00f7748SHong Zhang 1153d7d82daaSBarry Smith .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(), 1154d7d82daaSBarry Smith MatFactorInfoInitialize() 1155b00f7748SHong Zhang 1156b00f7748SHong Zhang S*/ 1157b00f7748SHong Zhang typedef struct { 115815e8a5b3SHong Zhang PetscReal diagonal_fill; /* force diagonal to fill in if initially not filled */ 115985317021SBarry Smith PetscReal usedt; 116015e8a5b3SHong Zhang PetscReal dt; /* drop tolerance */ 1161b00f7748SHong Zhang PetscReal dtcol; /* tolerance for pivoting */ 116215e8a5b3SHong Zhang PetscReal dtcount; /* maximum nonzeros to be allowed per row */ 116367e28bfeSBarry Smith PetscReal fill; /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */ 1164348344bbSBarry Smith PetscReal levels; /* ICC/ILU(levels) */ 1165bcd9e38bSBarry Smith PetscReal pivotinblocks; /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0 1166bcd9e38bSBarry Smith factorization may be faster if do not pivot */ 116715e8a5b3SHong Zhang PetscReal zeropivot; /* pivot is called zero if less than this */ 1168f4db908eSBarry Smith PetscReal shifttype; /* type of shift added to matrix factor to prevent zero pivots */ 1169f4db908eSBarry Smith PetscReal shiftamount; /* how large the shift is */ 117015e8a5b3SHong Zhang } MatFactorInfo; 1171ffa6d0a5SLois Curfman McInnes 11727087cfbeSBarry Smith extern PetscErrorCode MatFactorInfoInitialize(MatFactorInfo*); 11737087cfbeSBarry Smith extern PetscErrorCode MatCholeskyFactor(Mat,IS,const MatFactorInfo*); 11747087cfbeSBarry Smith extern PetscErrorCode MatCholeskyFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*); 11757087cfbeSBarry Smith extern PetscErrorCode MatCholeskyFactorNumeric(Mat,Mat,const MatFactorInfo*); 11767087cfbeSBarry Smith extern PetscErrorCode MatLUFactor(Mat,IS,IS,const MatFactorInfo*); 11777087cfbeSBarry Smith extern PetscErrorCode MatILUFactor(Mat,IS,IS,const MatFactorInfo*); 11787087cfbeSBarry Smith extern PetscErrorCode MatLUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*); 11797087cfbeSBarry Smith extern PetscErrorCode MatILUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*); 11807087cfbeSBarry Smith extern PetscErrorCode MatICCFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*); 11817087cfbeSBarry Smith extern PetscErrorCode MatICCFactor(Mat,IS,const MatFactorInfo*); 11827087cfbeSBarry Smith extern PetscErrorCode MatLUFactorNumeric(Mat,Mat,const MatFactorInfo*); 11837087cfbeSBarry Smith extern PetscErrorCode MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*); 11847087cfbeSBarry Smith extern PetscErrorCode MatSolve(Mat,Vec,Vec); 11857087cfbeSBarry Smith extern PetscErrorCode MatForwardSolve(Mat,Vec,Vec); 11867087cfbeSBarry Smith extern PetscErrorCode MatBackwardSolve(Mat,Vec,Vec); 11877087cfbeSBarry Smith extern PetscErrorCode MatSolveAdd(Mat,Vec,Vec,Vec); 11887087cfbeSBarry Smith extern PetscErrorCode MatSolveTranspose(Mat,Vec,Vec); 11897087cfbeSBarry Smith extern PetscErrorCode MatSolveTransposeAdd(Mat,Vec,Vec,Vec); 11907087cfbeSBarry Smith extern PetscErrorCode MatSolves(Mat,Vecs,Vecs); 11918ed539a5SBarry Smith 11927087cfbeSBarry Smith extern PetscErrorCode MatSetUnfactored(Mat); 1193bb5a7306SBarry Smith 1194d91e6319SBarry Smith /*E 1195d91e6319SBarry Smith MatSORType - What type of (S)SOR to perform 1196bb1eb677SSatish Balay 1197d91e6319SBarry Smith Level: beginner 1198d91e6319SBarry Smith 1199d9274352SBarry Smith May be bitwise ORd together 1200d9274352SBarry Smith 1201d91e6319SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 1202d91e6319SBarry Smith 12034e7234bfSBarry Smith MatSORType may be bitwise ORd together, so do not change the numbers 12044e7234bfSBarry Smith 120541f059aeSBarry Smith .seealso: MatSOR() 1206d91e6319SBarry Smith E*/ 1207ee50ffe9SBarry Smith typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3, 1208ee50ffe9SBarry Smith SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8, 1209ee50ffe9SBarry Smith SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16, 121084cb2905SBarry Smith SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType; 12117087cfbeSBarry Smith extern PetscErrorCode MatSOR(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec); 12128ed539a5SBarry Smith 1213d4fbbf0eSBarry Smith /* 1214639f9d9dSBarry Smith These routines are for efficiently computing Jacobians via finite differences. 1215639f9d9dSBarry Smith */ 1216b12f92e5SBarry Smith 1217d9274352SBarry Smith /*E 1218d9274352SBarry Smith MatColoringType - String with the name of a PETSc matrix coloring or the creation function 1219d9274352SBarry Smith with an optional dynamic library name, for example 1220d9274352SBarry Smith http://www.mcs.anl.gov/petsc/lib.a:coloringcreate() 1221d9274352SBarry Smith 1222d9274352SBarry Smith Level: beginner 1223d9274352SBarry Smith 1224d9274352SBarry Smith .seealso: MatGetColoring() 1225d9274352SBarry Smith E*/ 1226a313700dSBarry Smith #define MatColoringType char* 12272692d6eeSBarry Smith #define MATCOLORINGNATURAL "natural" 12282692d6eeSBarry Smith #define MATCOLORINGSL "sl" 12292692d6eeSBarry Smith #define MATCOLORINGLF "lf" 12302692d6eeSBarry Smith #define MATCOLORINGID "id" 1231b12f92e5SBarry Smith 12327087cfbeSBarry Smith extern PetscErrorCode MatGetColoring(Mat,const MatColoringType,ISColoring*); 12337087cfbeSBarry Smith extern PetscErrorCode MatColoringRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,MatColoringType,ISColoring *)); 123430de9b25SBarry Smith 123530de9b25SBarry Smith /*MC 123630de9b25SBarry Smith MatColoringRegisterDynamic - Adds a new sparse matrix coloring to the 123730de9b25SBarry Smith matrix package. 123830de9b25SBarry Smith 123930de9b25SBarry Smith Synopsis: 12401890ba74SBarry Smith PetscErrorCode MatColoringRegisterDynamic(const char *name_coloring,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatColoring)) 124130de9b25SBarry Smith 124230de9b25SBarry Smith Not Collective 124330de9b25SBarry Smith 124430de9b25SBarry Smith Input Parameters: 12452692d6eeSBarry Smith + sname - name of Coloring (for example MATCOLORINGSL) 124630de9b25SBarry Smith . path - location of library where creation routine is 124730de9b25SBarry Smith . name - name of function that creates the Coloring type, a string 124830de9b25SBarry Smith - function - function pointer that creates the coloring 124930de9b25SBarry Smith 125030de9b25SBarry Smith Level: developer 125130de9b25SBarry Smith 125230de9b25SBarry Smith If dynamic libraries are used, then the fourth input argument (function) 125330de9b25SBarry Smith is ignored. 125430de9b25SBarry Smith 125530de9b25SBarry Smith Sample usage: 125630de9b25SBarry Smith .vb 125730de9b25SBarry Smith MatColoringRegisterDynamic("my_color",/home/username/my_lib/lib/libO/solaris/mylib.a, 125830de9b25SBarry Smith "MyColor",MyColor); 125930de9b25SBarry Smith .ve 126030de9b25SBarry Smith 126130de9b25SBarry Smith Then, your partitioner can be chosen with the procedural interface via 126230de9b25SBarry Smith $ MatColoringSetType(part,"my_color") 126330de9b25SBarry Smith or at runtime via the option 126430de9b25SBarry Smith $ -mat_coloring_type my_color 126530de9b25SBarry Smith 1266ab901514SBarry Smith $PETSC_ARCH occuring in pathname will be replaced with appropriate values. 126730de9b25SBarry Smith 126830de9b25SBarry Smith .keywords: matrix, Coloring, register 126930de9b25SBarry Smith 127030de9b25SBarry Smith .seealso: MatColoringRegisterDestroy(), MatColoringRegisterAll() 127130de9b25SBarry Smith M*/ 1272aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 1273f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,0) 1274b12f92e5SBarry Smith #else 1275f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,d) 1276b12f92e5SBarry Smith #endif 127730de9b25SBarry Smith 1278ace3abfcSBarry Smith extern PetscBool MatColoringRegisterAllCalled; 1279f1144a30SSatish Balay 12807087cfbeSBarry Smith extern PetscErrorCode MatColoringRegisterAll(const char[]); 12817087cfbeSBarry Smith extern PetscErrorCode MatColoringRegisterDestroy(void); 12827087cfbeSBarry Smith extern PetscErrorCode MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*); 1283639f9d9dSBarry Smith 1284d9274352SBarry Smith /*S 1285d9274352SBarry Smith MatFDColoring - Object for computing a sparse Jacobian via finite differences 1286d9274352SBarry Smith and coloring 1287639f9d9dSBarry Smith 1288d9274352SBarry Smith Level: beginner 1289d9274352SBarry Smith 1290d9274352SBarry Smith Concepts: coloring, sparse Jacobian, finite differences 1291d9274352SBarry Smith 1292d9274352SBarry Smith .seealso: MatFDColoringCreate() 1293d9274352SBarry Smith S*/ 1294e2a1c21fSSatish Balay typedef struct _p_MatFDColoring* MatFDColoring; 1295639f9d9dSBarry Smith 12967087cfbeSBarry Smith extern PetscErrorCode MatFDColoringCreate(Mat,ISColoring,MatFDColoring *); 1297fcfd50ebSBarry Smith extern PetscErrorCode MatFDColoringDestroy(MatFDColoring*); 12987087cfbeSBarry Smith extern PetscErrorCode MatFDColoringView(MatFDColoring,PetscViewer); 12997087cfbeSBarry Smith extern PetscErrorCode MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*); 13007087cfbeSBarry Smith extern PetscErrorCode MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**); 13017087cfbeSBarry Smith extern PetscErrorCode MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal); 13027087cfbeSBarry Smith extern PetscErrorCode MatFDColoringSetFromOptions(MatFDColoring); 13037087cfbeSBarry Smith extern PetscErrorCode MatFDColoringApply(Mat,MatFDColoring,Vec,MatStructure*,void *); 13047087cfbeSBarry Smith extern PetscErrorCode MatFDColoringApplyTS(Mat,MatFDColoring,PetscReal,Vec,MatStructure*,void *); 13057087cfbeSBarry Smith extern PetscErrorCode MatFDColoringSetF(MatFDColoring,Vec); 13067087cfbeSBarry Smith extern PetscErrorCode MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,PetscInt*[]); 1307639f9d9dSBarry Smith /* 13080752156aSBarry Smith These routines are for partitioning matrices: currently used only 13093eda8832SBarry Smith for adjacency matrix, MatCreateMPIAdj(). 13100752156aSBarry Smith */ 1311ca161407SBarry Smith 1312d9274352SBarry Smith /*S 1313d9274352SBarry Smith MatPartitioning - Object for managing the partitioning of a matrix or graph 1314d9274352SBarry Smith 1315d9274352SBarry Smith Level: beginner 1316d9274352SBarry Smith 1317d9274352SBarry Smith Concepts: partitioning 1318d9274352SBarry Smith 1319743a1026Svictorle .seealso: MatPartitioningCreate(), MatPartitioningType 1320d9274352SBarry Smith S*/ 132191e9ee9fSBarry Smith typedef struct _p_MatPartitioning* MatPartitioning; 1322d9274352SBarry Smith 1323d9274352SBarry Smith /*E 13245bc6e7ccSvictorle MatPartitioningType - String with the name of a PETSc matrix partitioning or the creation function 1325d9274352SBarry Smith with an optional dynamic library name, for example 1326d9274352SBarry Smith http://www.mcs.anl.gov/petsc/lib.a:partitioningcreate() 1327d9274352SBarry Smith 1328d9274352SBarry Smith Level: beginner 1329d9274352SBarry Smith 1330b547a13bSvictorle .seealso: MatPartitioningCreate(), MatPartitioning 1331d9274352SBarry Smith E*/ 1332a313700dSBarry Smith #define MatPartitioningType char* 13332692d6eeSBarry Smith #define MATPARTITIONINGCURRENT "current" 13342692d6eeSBarry Smith #define MATPARTITIONINGSQUARE "square" 13352692d6eeSBarry Smith #define MATPARTITIONINGPARMETIS "parmetis" 13362692d6eeSBarry Smith #define MATPARTITIONINGCHACO "chaco" 13372692d6eeSBarry Smith #define MATPARTITIONINGJOSTLE "jostle" 13382692d6eeSBarry Smith #define MATPARTITIONINGPARTY "party" 13392692d6eeSBarry Smith #define MATPARTITIONINGSCOTCH "scotch" 134071306c60Sjroman 1341ca161407SBarry Smith 13427087cfbeSBarry Smith extern PetscErrorCode MatPartitioningCreate(MPI_Comm,MatPartitioning*); 13437087cfbeSBarry Smith extern PetscErrorCode MatPartitioningSetType(MatPartitioning,const MatPartitioningType); 13447087cfbeSBarry Smith extern PetscErrorCode MatPartitioningSetNParts(MatPartitioning,PetscInt); 13457087cfbeSBarry Smith extern PetscErrorCode MatPartitioningSetAdjacency(MatPartitioning,Mat); 13467087cfbeSBarry Smith extern PetscErrorCode MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]); 13477087cfbeSBarry Smith extern PetscErrorCode MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []); 13487087cfbeSBarry Smith extern PetscErrorCode MatPartitioningApply(MatPartitioning,IS*); 1349fcfd50ebSBarry Smith extern PetscErrorCode MatPartitioningDestroy(MatPartitioning*); 13502aabb6bbSBarry Smith 13517087cfbeSBarry Smith extern PetscErrorCode MatPartitioningRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatPartitioning)); 135230de9b25SBarry Smith 135330de9b25SBarry Smith /*MC 135430de9b25SBarry Smith MatPartitioningRegisterDynamic - Adds a new sparse matrix partitioning to the 135530de9b25SBarry Smith matrix package. 135630de9b25SBarry Smith 135730de9b25SBarry Smith Synopsis: 13581890ba74SBarry Smith PetscErrorCode MatPartitioningRegisterDynamic(const char *name_partitioning,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatPartitioning)) 135930de9b25SBarry Smith 136030de9b25SBarry Smith Not Collective 136130de9b25SBarry Smith 136230de9b25SBarry Smith Input Parameters: 13632692d6eeSBarry Smith + sname - name of partitioning (for example MATPARTITIONINGCURRENT) or parmetis 136430de9b25SBarry Smith . path - location of library where creation routine is 136530de9b25SBarry Smith . name - name of function that creates the partitioning type, a string 136630de9b25SBarry Smith - function - function pointer that creates the partitioning type 136730de9b25SBarry Smith 136830de9b25SBarry Smith Level: developer 136930de9b25SBarry Smith 137030de9b25SBarry Smith If dynamic libraries are used, then the fourth input argument (function) 137130de9b25SBarry Smith is ignored. 137230de9b25SBarry Smith 137330de9b25SBarry Smith Sample usage: 137430de9b25SBarry Smith .vb 137530de9b25SBarry Smith MatPartitioningRegisterDynamic("my_part",/home/username/my_lib/lib/libO/solaris/mylib.a, 137630de9b25SBarry Smith "MyPartCreate",MyPartCreate); 137730de9b25SBarry Smith .ve 137830de9b25SBarry Smith 137930de9b25SBarry Smith Then, your partitioner can be chosen with the procedural interface via 138030de9b25SBarry Smith $ MatPartitioningSetType(part,"my_part") 138130de9b25SBarry Smith or at runtime via the option 138230de9b25SBarry Smith $ -mat_partitioning_type my_part 138330de9b25SBarry Smith 1384ab901514SBarry Smith $PETSC_ARCH occuring in pathname will be replaced with appropriate values. 138530de9b25SBarry Smith 138630de9b25SBarry Smith .keywords: matrix, partitioning, register 138730de9b25SBarry Smith 138830de9b25SBarry Smith .seealso: MatPartitioningRegisterDestroy(), MatPartitioningRegisterAll() 138930de9b25SBarry Smith M*/ 1390aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 1391f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,0) 13922aabb6bbSBarry Smith #else 1393f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,d) 13942aabb6bbSBarry Smith #endif 13952aabb6bbSBarry Smith 1396ace3abfcSBarry Smith extern PetscBool MatPartitioningRegisterAllCalled; 1397f1144a30SSatish Balay 13987087cfbeSBarry Smith extern PetscErrorCode MatPartitioningRegisterAll(const char[]); 13997087cfbeSBarry Smith extern PetscErrorCode MatPartitioningRegisterDestroy(void); 14002bad1931SBarry Smith 14017087cfbeSBarry Smith extern PetscErrorCode MatPartitioningView(MatPartitioning,PetscViewer); 14027087cfbeSBarry Smith extern PetscErrorCode MatPartitioningSetFromOptions(MatPartitioning); 14037087cfbeSBarry Smith extern PetscErrorCode MatPartitioningGetType(MatPartitioning,const MatPartitioningType*); 1404ca161407SBarry Smith 14057087cfbeSBarry Smith extern PetscErrorCode MatPartitioningParmetisSetCoarseSequential(MatPartitioning); 14067087cfbeSBarry Smith extern PetscErrorCode MatPartitioningParmetisGetEdgeCut(MatPartitioning, PetscInt *); 14070752156aSBarry Smith 14087087cfbeSBarry Smith extern PetscErrorCode MatPartitioningJostleSetCoarseLevel(MatPartitioning,PetscReal); 14097087cfbeSBarry Smith extern PetscErrorCode MatPartitioningJostleSetCoarseSequential(MatPartitioning); 141071306c60Sjroman 1411954b3ec6SBarry Smith typedef enum {MP_CHACO_MULTILEVEL_KL,MP_CHACO_SPECTRAL,MP_CHACO_LINEAR,MP_CHACO_RANDOM, MP_CHACO_SCATTERED} MPChacoGlobalType; 14127087cfbeSBarry Smith extern PetscErrorCode MatPartitioningChacoSetGlobal(MatPartitioning, MPChacoGlobalType); 141371306c60Sjroman typedef enum { MP_CHACO_KERNIGHAN_LIN, MP_CHACO_NONE } MPChacoLocalType; 14147087cfbeSBarry Smith extern PetscErrorCode MatPartitioningChacoSetLocal(MatPartitioning, MPChacoLocalType); 14157087cfbeSBarry Smith extern PetscErrorCode MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal); 141671306c60Sjroman typedef enum { MP_CHACO_LANCZOS, MP_CHACO_RQI_SYMMLQ } MPChacoEigenType; 14177087cfbeSBarry Smith extern PetscErrorCode MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType); 14187087cfbeSBarry Smith extern PetscErrorCode MatPartitioningChacoSetEigenTol(MatPartitioning, PetscReal); 14197087cfbeSBarry Smith extern PetscErrorCode MatPartitioningChacoSetEigenNumber(MatPartitioning, PetscInt); 142071306c60Sjroman 142171306c60Sjroman #define MP_PARTY_OPT "opt" 142271306c60Sjroman #define MP_PARTY_LIN "lin" 142371306c60Sjroman #define MP_PARTY_SCA "sca" 142471306c60Sjroman #define MP_PARTY_RAN "ran" 142571306c60Sjroman #define MP_PARTY_GBF "gbf" 142671306c60Sjroman #define MP_PARTY_GCF "gcf" 142771306c60Sjroman #define MP_PARTY_BUB "bub" 142871306c60Sjroman #define MP_PARTY_DEF "def" 14297087cfbeSBarry Smith extern PetscErrorCode MatPartitioningPartySetGlobal(MatPartitioning, const char*); 143071306c60Sjroman #define MP_PARTY_HELPFUL_SETS "hs" 143171306c60Sjroman #define MP_PARTY_KERNIGHAN_LIN "kl" 143271306c60Sjroman #define MP_PARTY_NONE "no" 14337087cfbeSBarry Smith extern PetscErrorCode MatPartitioningPartySetLocal(MatPartitioning, const char*); 14347087cfbeSBarry Smith extern PetscErrorCode MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal); 14357087cfbeSBarry Smith extern PetscErrorCode MatPartitioningPartySetBipart(MatPartitioning,PetscBool ); 14367087cfbeSBarry Smith extern PetscErrorCode MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscBool ); 143771306c60Sjroman 143871306c60Sjroman typedef enum { MP_SCOTCH_GREEDY, MP_SCOTCH_GPS, MP_SCOTCH_GR_GPS } MPScotchGlobalType; 14397087cfbeSBarry Smith extern PetscErrorCode MatPartitioningScotchSetArch(MatPartitioning,const char*); 14407087cfbeSBarry Smith extern PetscErrorCode MatPartitioningScotchSetMultilevel(MatPartitioning); 14417087cfbeSBarry Smith extern PetscErrorCode MatPartitioningScotchSetGlobal(MatPartitioning,MPScotchGlobalType); 14427087cfbeSBarry Smith extern PetscErrorCode MatPartitioningScotchSetCoarseLevel(MatPartitioning,PetscReal); 14437087cfbeSBarry Smith extern PetscErrorCode MatPartitioningScotchSetHostList(MatPartitioning,const char*); 144471306c60Sjroman typedef enum { MP_SCOTCH_KERNIGHAN_LIN, MP_SCOTCH_NONE } MPScotchLocalType; 14457087cfbeSBarry Smith extern PetscErrorCode MatPartitioningScotchSetLocal(MatPartitioning,MPScotchLocalType); 14467087cfbeSBarry Smith extern PetscErrorCode MatPartitioningScotchSetMapping(MatPartitioning); 14477087cfbeSBarry Smith extern PetscErrorCode MatPartitioningScotchSetStrategy(MatPartitioning,char*); 144871306c60Sjroman 144909573ac7SBarry Smith extern PetscErrorCode MatMeshToVertexGraph(Mat,PetscInt,Mat*); 145009573ac7SBarry Smith extern PetscErrorCode MatMeshToCellGraph(Mat,PetscInt,Mat*); 1451591294e4SBarry Smith 14520752156aSBarry Smith /* 14530a835dfdSSatish Balay If you add entries here you must also add them to finclude/petscmat.h 1454d4fbbf0eSBarry Smith */ 14551c1c02c0SLois Curfman McInnes typedef enum { MATOP_SET_VALUES=0, 14561c1c02c0SLois Curfman McInnes MATOP_GET_ROW=1, 14571c1c02c0SLois Curfman McInnes MATOP_RESTORE_ROW=2, 14581c1c02c0SLois Curfman McInnes MATOP_MULT=3, 14591c1c02c0SLois Curfman McInnes MATOP_MULT_ADD=4, 14607c922b88SBarry Smith MATOP_MULT_TRANSPOSE=5, 14617c922b88SBarry Smith MATOP_MULT_TRANSPOSE_ADD=6, 14621c1c02c0SLois Curfman McInnes MATOP_SOLVE=7, 14631c1c02c0SLois Curfman McInnes MATOP_SOLVE_ADD=8, 14647c922b88SBarry Smith MATOP_SOLVE_TRANSPOSE=9, 14657c922b88SBarry Smith MATOP_SOLVE_TRANSPOSE_ADD=10, 14661c1c02c0SLois Curfman McInnes MATOP_LUFACTOR=11, 14671c1c02c0SLois Curfman McInnes MATOP_CHOLESKYFACTOR=12, 1468a714c06dSBarry Smith MATOP_SOR=13, 14691c1c02c0SLois Curfman McInnes MATOP_TRANSPOSE=14, 14701c1c02c0SLois Curfman McInnes MATOP_GETINFO=15, 14711c1c02c0SLois Curfman McInnes MATOP_EQUAL=16, 14721c1c02c0SLois Curfman McInnes MATOP_GET_DIAGONAL=17, 14731c1c02c0SLois Curfman McInnes MATOP_DIAGONAL_SCALE=18, 14741c1c02c0SLois Curfman McInnes MATOP_NORM=19, 14751c1c02c0SLois Curfman McInnes MATOP_ASSEMBLY_BEGIN=20, 14761c1c02c0SLois Curfman McInnes MATOP_ASSEMBLY_END=21, 1477d519adbfSMatthew Knepley MATOP_SET_OPTION=22, 1478d519adbfSMatthew Knepley MATOP_ZERO_ENTRIES=23, 1479d519adbfSMatthew Knepley MATOP_ZERO_ROWS=24, 1480d519adbfSMatthew Knepley MATOP_LUFACTOR_SYMBOLIC=25, 1481d519adbfSMatthew Knepley MATOP_LUFACTOR_NUMERIC=26, 1482d519adbfSMatthew Knepley MATOP_CHOLESKY_FACTOR_SYMBOLIC=27, 1483d519adbfSMatthew Knepley MATOP_CHOLESKY_FACTOR_NUMERIC=28, 1484d519adbfSMatthew Knepley MATOP_SETUP_PREALLOCATION=29, 1485d519adbfSMatthew Knepley MATOP_ILUFACTOR_SYMBOLIC=30, 1486d519adbfSMatthew Knepley MATOP_ICCFACTOR_SYMBOLIC=31, 1487d519adbfSMatthew Knepley MATOP_GET_ARRAY=32, 1488d519adbfSMatthew Knepley MATOP_RESTORE_ARRAY=33, 1489d519adbfSMatthew Knepley MATOP_DUPLICATE=34, 1490d519adbfSMatthew Knepley MATOP_FORWARD_SOLVE=35, 1491d519adbfSMatthew Knepley MATOP_BACKWARD_SOLVE=36, 1492d519adbfSMatthew Knepley MATOP_ILUFACTOR=37, 1493d519adbfSMatthew Knepley MATOP_ICCFACTOR=38, 1494d519adbfSMatthew Knepley MATOP_AXPY=39, 1495d519adbfSMatthew Knepley MATOP_GET_SUBMATRICES=40, 1496d519adbfSMatthew Knepley MATOP_INCREASE_OVERLAP=41, 1497d519adbfSMatthew Knepley MATOP_GET_VALUES=42, 1498d519adbfSMatthew Knepley MATOP_COPY=43, 1499d519adbfSMatthew Knepley MATOP_GET_ROW_MAX=44, 1500d519adbfSMatthew Knepley MATOP_SCALE=45, 1501d519adbfSMatthew Knepley MATOP_SHIFT=46, 150235153367SBarry Smith MATOP_DIAGONAL_SET=47, 1503d519adbfSMatthew Knepley MATOP_ILUDT_FACTOR=48, 1504d519adbfSMatthew Knepley MATOP_SET_BLOCK_SIZE=49, 1505d519adbfSMatthew Knepley MATOP_GET_ROW_IJ=50, 1506d519adbfSMatthew Knepley MATOP_RESTORE_ROW_IJ=51, 1507d519adbfSMatthew Knepley MATOP_GET_COLUMN_IJ=52, 1508d519adbfSMatthew Knepley MATOP_RESTORE_COLUMN_IJ=53, 1509d519adbfSMatthew Knepley MATOP_FDCOLORING_CREATE=54, 1510d519adbfSMatthew Knepley MATOP_COLORING_PATCH=55, 1511d519adbfSMatthew Knepley MATOP_SET_UNFACTORED=56, 1512d519adbfSMatthew Knepley MATOP_PERMUTE=57, 1513d519adbfSMatthew Knepley MATOP_SET_VALUES_BLOCKED=58, 1514d519adbfSMatthew Knepley MATOP_GET_SUBMATRIX=59, 1515d519adbfSMatthew Knepley MATOP_DESTROY=60, 1516d519adbfSMatthew Knepley MATOP_VIEW=61, 1517d519adbfSMatthew Knepley MATOP_CONVERT_FROM=62, 1518d519adbfSMatthew Knepley MATOP_USE_SCALED_FORM=63, 1519d519adbfSMatthew Knepley MATOP_SCALE_SYSTEM=64, 1520d519adbfSMatthew Knepley MATOP_UNSCALE_SYSTEM=65, 1521d519adbfSMatthew Knepley MATOP_SET_LOCAL_TO_GLOBAL_MAP=66, 1522d519adbfSMatthew Knepley MATOP_SET_VALUES_LOCAL=67, 1523d519adbfSMatthew Knepley MATOP_ZERO_ROWS_LOCAL=68, 1524d519adbfSMatthew Knepley MATOP_GET_ROW_MAX_ABS=69, 1525d519adbfSMatthew Knepley MATOP_GET_ROW_MIN_ABS=70, 1526d519adbfSMatthew Knepley MATOP_CONVERT=71, 1527d519adbfSMatthew Knepley MATOP_SET_COLORING=72, 1528d519adbfSMatthew Knepley MATOP_SET_VALUES_ADIC=73, 1529d519adbfSMatthew Knepley MATOP_SET_VALUES_ADIFOR=74, 1530d519adbfSMatthew Knepley MATOP_FD_COLORING_APPLY=75, 1531d519adbfSMatthew Knepley MATOP_SET_FROM_OPTIONS=76, 1532d519adbfSMatthew Knepley MATOP_MULT_CON=77, 1533d519adbfSMatthew Knepley MATOP_MULT_TRANSPOSE_CON=78, 1534d519adbfSMatthew Knepley MATOP_PERMUTE_SPARSIFY=79, 1535d519adbfSMatthew Knepley MATOP_MULT_MULTIPLE=80, 1536d519adbfSMatthew Knepley MATOP_SOLVE_MULTIPLE=81, 1537d519adbfSMatthew Knepley MATOP_GET_INERTIA=82, 1538d519adbfSMatthew Knepley MATOP_LOAD=83, 1539d519adbfSMatthew Knepley MATOP_IS_SYMMETRIC=84, 1540d519adbfSMatthew Knepley MATOP_IS_HERMITIAN=85, 1541d519adbfSMatthew Knepley MATOP_IS_STRUCTURALLY_SYMMETRIC=86, 154241f059aeSBarry Smith MATOP_DUMMY=87, 1543d519adbfSMatthew Knepley MATOP_GET_VECS=88, 1544d519adbfSMatthew Knepley MATOP_MAT_MULT=89, 1545d519adbfSMatthew Knepley MATOP_MAT_MULT_SYMBOLIC=90, 1546d519adbfSMatthew Knepley MATOP_MAT_MULT_NUMERIC=91, 1547d519adbfSMatthew Knepley MATOP_PTAP=92, 1548d519adbfSMatthew Knepley MATOP_PTAP_SYMBOLIC=93, 1549d519adbfSMatthew Knepley MATOP_PTAP_NUMERIC=94, 1550d519adbfSMatthew Knepley MATOP_MAT_MULTTRANSPOSE=95, 15511763ac4eSSatish Balay MATOP_MAT_MULTTRANSPOSE_SYM=96, 15521763ac4eSSatish Balay MATOP_MAT_MULTTRANSPOSE_NUM=97, 1553d519adbfSMatthew Knepley MATOP_PTAP_SYMBOLIC_SEQAIJ=98, 1554d519adbfSMatthew Knepley MATOP_PTAP_NUMERIC_SEQAIJ=99, 1555d519adbfSMatthew Knepley MATOP_PTAP_SYMBOLIC_MPIAIJ=100, 1556d519adbfSMatthew Knepley MATOP_PTAP_NUMERIC_MPIAIJ=101, 1557d519adbfSMatthew Knepley MATOP_CONJUGATE=102, 1558d519adbfSMatthew Knepley MATOP_SET_SIZES=103, 1559d519adbfSMatthew Knepley MATOP_SET_VALUES_ROW=104, 1560d519adbfSMatthew Knepley MATOP_REAL_PART=105, 1561d519adbfSMatthew Knepley MATOP_IMAG_PART=106, 1562d519adbfSMatthew Knepley MATOP_GET_ROW_UTRIANGULAR=107, 1563d519adbfSMatthew Knepley MATOP_RESTORE_ROW_UTRIANGULAR=108, 1564d519adbfSMatthew Knepley MATOP_MATSOLVE=109, 1565d519adbfSMatthew Knepley MATOP_GET_REDUNDANTMATRIX=110, 1566d519adbfSMatthew Knepley MATOP_GET_ROW_MIN=111, 1567d519adbfSMatthew Knepley MATOP_GET_COLUMN_VEC=112, 1568d519adbfSMatthew Knepley MATOP_MISSING_DIAGONAL=113, 1569d519adbfSMatthew Knepley MATOP_MATGETSEQNONZEROSTRUCTURE=114, 157089c6957cSBarry Smith MATOP_CREATE=115, 157189c6957cSBarry Smith MATOP_GET_GHOSTS=116, 1572ea38565cSJed Brown MATOP_GET_LOCALSUBMATRIX=117, 1573ea38565cSJed Brown MATOP_RESTORE_LOCALSUBMATRIX=118, 1574eeffb40dSHong Zhang MATOP_MULT_DIAGONAL_BLOCK=119, 1575547795f9SHong Zhang MATOP_HERMITIANTRANSPOSE=120, 1576547795f9SHong Zhang MATOP_MULTHERMITIANTRANSPOSE=121, 1577fbdbba38SShri Abhyankar MATOP_MULTHERMITIANTRANSPOSEADD=122, 1578b9614d88SDmitry Karpeev MATOP_GETMULTIPROCBLOCK=123, 15790716a85fSBarry Smith MATOP_GETCOLUMNNORMS=125, 158037868618SMatthew G Knepley MATOP_GET_SUBMATRICES_PARALLEL=128, 158137868618SMatthew G Knepley MATOP_SET_VALUES_BATCH=129 1582fae171e0SBarry Smith } MatOperation; 15837087cfbeSBarry Smith extern PetscErrorCode MatHasOperation(Mat,MatOperation,PetscBool *); 15847087cfbeSBarry Smith extern PetscErrorCode MatShellSetOperation(Mat,MatOperation,void(*)(void)); 15857087cfbeSBarry Smith extern PetscErrorCode MatShellGetOperation(Mat,MatOperation,void(**)(void)); 15867087cfbeSBarry Smith extern PetscErrorCode MatShellSetContext(Mat,void*); 1587112a2221SBarry Smith 158890ace30eSBarry Smith /* 158990ace30eSBarry Smith Codes for matrices stored on disk. By default they are 159090ace30eSBarry Smith stored in a universal format. By changing the format with 15917973ad04SBarry Smith PetscViewerSetFormat(viewer,PETSC_VIEWER_NATIVE); the matrices will 159290ace30eSBarry Smith be stored in a way natural for the matrix, for example dense matrices 159390ace30eSBarry Smith would be stored as dense. Matrices stored this way may only be 1594f4403165SShri Abhyankar read into matrices of the same type. 159590ace30eSBarry Smith */ 159690ace30eSBarry Smith #define MATRIX_BINARY_FORMAT_DENSE -1 159790ace30eSBarry Smith 15987087cfbeSBarry Smith extern PetscErrorCode MatMPIBAIJSetHashTableFactor(Mat,PetscReal); 15997087cfbeSBarry Smith extern PetscErrorCode MatISGetLocalMat(Mat,Mat*); 16001f4e1ec7SBarry Smith 1601d9274352SBarry Smith /*S 1602d9274352SBarry Smith MatNullSpace - Object that removes a null space from a vector, i.e. 1603d9274352SBarry Smith orthogonalizes the vector to a subsapce 1604d9274352SBarry Smith 1605f7a9e4ceSBarry Smith Level: advanced 1606d9274352SBarry Smith 1607d9274352SBarry Smith Concepts: matrix; linear operator, null space 1608d9274352SBarry Smith 16096e1639daSBarry Smith Users manual sections: 16106e1639daSBarry Smith . sec_singular 16116e1639daSBarry Smith 1612d9274352SBarry Smith .seealso: MatNullSpaceCreate() 1613d9274352SBarry Smith S*/ 161474637425SBarry Smith typedef struct _p_MatNullSpace* MatNullSpace; 1615d9274352SBarry Smith 16167087cfbeSBarry Smith extern PetscErrorCode MatNullSpaceCreate(MPI_Comm,PetscBool ,PetscInt,const Vec[],MatNullSpace*); 16177087cfbeSBarry Smith extern PetscErrorCode MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(MatNullSpace,Vec,void*),void*); 1618d34fcf5fSBarry Smith extern PetscErrorCode MatNullSpaceDestroy(MatNullSpace*); 16197087cfbeSBarry Smith extern PetscErrorCode MatNullSpaceRemove(MatNullSpace,Vec,Vec*); 16207087cfbeSBarry Smith extern PetscErrorCode MatNullSpaceAttach(Mat,MatNullSpace); 16217087cfbeSBarry Smith extern PetscErrorCode MatNullSpaceTest(MatNullSpace,Mat,PetscBool *); 1622b717e993SJed Brown extern PetscErrorCode MatNullSpaceView(MatNullSpace,PetscViewer); 162374637425SBarry Smith 16247087cfbeSBarry Smith extern PetscErrorCode MatReorderingSeqSBAIJ(Mat,IS); 16257087cfbeSBarry Smith extern PetscErrorCode MatMPISBAIJSetHashTableFactor(Mat,PetscReal); 16267087cfbeSBarry Smith extern PetscErrorCode MatSeqSBAIJSetColumnIndices(Mat,PetscInt *); 16277087cfbeSBarry Smith extern PetscErrorCode MatSeqBAIJInvertBlockDiagonal(Mat); 16283f1d51d7SBarry Smith 16297087cfbeSBarry Smith extern PetscErrorCode MatCreateMAIJ(Mat,PetscInt,Mat*); 16307087cfbeSBarry Smith extern PetscErrorCode MatMAIJRedimension(Mat,PetscInt,Mat*); 16317087cfbeSBarry Smith extern PetscErrorCode MatMAIJGetAIJ(Mat,Mat*); 1632c4f061fbSSatish Balay 16337087cfbeSBarry Smith extern PetscErrorCode MatComputeExplicitOperator(Mat,Mat*); 1634b0a32e0cSBarry Smith 16357087cfbeSBarry Smith extern PetscErrorCode MatDiagonalScaleLocal(Mat,Vec); 163604f1ad80SBarry Smith 16377087cfbeSBarry Smith extern PetscErrorCode MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*); 16387087cfbeSBarry Smith extern PetscErrorCode MatMFFDSetBase(Mat,Vec,Vec); 16397087cfbeSBarry Smith extern PetscErrorCode MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*); 16407087cfbeSBarry Smith extern PetscErrorCode MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*)); 16417087cfbeSBarry Smith extern PetscErrorCode MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec)); 16427087cfbeSBarry Smith extern PetscErrorCode MatMFFDAddNullSpace(Mat,MatNullSpace); 16437087cfbeSBarry Smith extern PetscErrorCode MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt); 16447087cfbeSBarry Smith extern PetscErrorCode MatMFFDResetHHistory(Mat); 16457087cfbeSBarry Smith extern PetscErrorCode MatMFFDSetFunctionError(Mat,PetscReal); 16467087cfbeSBarry Smith extern PetscErrorCode MatMFFDSetPeriod(Mat,PetscInt); 16477087cfbeSBarry Smith extern PetscErrorCode MatMFFDGetH(Mat,PetscScalar *); 16487087cfbeSBarry Smith extern PetscErrorCode MatMFFDSetOptionsPrefix(Mat,const char[]); 16497087cfbeSBarry Smith extern PetscErrorCode MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*); 16507087cfbeSBarry Smith extern PetscErrorCode MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*); 1651e884886fSBarry Smith 16526370053bSBarry Smith /*S 16536370053bSBarry Smith MatMFFD - A data structured used to manage the computation of the h differencing parameter for matrix-free 16546370053bSBarry Smith Jacobian vector products 1655e884886fSBarry Smith 16566370053bSBarry Smith Notes: MATMFFD is a specific MatType which uses the MatMFFD data structure 16576370053bSBarry Smith 16586370053bSBarry Smith MatMFFD*() methods actually take the Mat as their first argument. Not a MatMFFD data structure 16596370053bSBarry Smith 16606370053bSBarry Smith Level: developer 16616370053bSBarry Smith 16626370053bSBarry Smith .seealso: MATMFFD, MatCreateMFFD(), MatMFFDSetFuction(), MatMFFDSetType(), MatMFFDRegister() 16636370053bSBarry Smith S*/ 1664e884886fSBarry Smith typedef struct _p_MatMFFD* MatMFFD; 1665e884886fSBarry Smith 1666e884886fSBarry Smith /*E 1667e884886fSBarry Smith MatMFFDType - algorithm used to compute the h used in computing matrix-vector products via differencing of the function 1668e884886fSBarry Smith 1669e884886fSBarry Smith Level: beginner 1670e884886fSBarry Smith 1671e884886fSBarry Smith .seealso: MatMFFDSetType(), MatMFFDRegister() 1672e884886fSBarry Smith E*/ 1673a313700dSBarry Smith #define MatMFFDType char* 1674e884886fSBarry Smith #define MATMFFD_DS "ds" 1675e884886fSBarry Smith #define MATMFFD_WP "wp" 1676e884886fSBarry Smith 16777087cfbeSBarry Smith extern PetscErrorCode MatMFFDSetType(Mat,const MatMFFDType); 16787087cfbeSBarry Smith extern PetscErrorCode MatMFFDRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatMFFD)); 1679e884886fSBarry Smith 1680e884886fSBarry Smith /*MC 1681e884886fSBarry Smith MatMFFDRegisterDynamic - Adds a method to the MatMFFD registry. 1682e884886fSBarry Smith 1683e884886fSBarry Smith Synopsis: 16841890ba74SBarry Smith PetscErrorCode MatMFFDRegisterDynamic(const char *name_solver,const char *path,const char *name_create,PetscErrorCode (*routine_create)(MatMFFD)) 1685e884886fSBarry Smith 1686e884886fSBarry Smith Not Collective 1687e884886fSBarry Smith 1688e884886fSBarry Smith Input Parameters: 1689e884886fSBarry Smith + name_solver - name of a new user-defined compute-h module 1690e884886fSBarry Smith . path - path (either absolute or relative) the library containing this solver 1691e884886fSBarry Smith . name_create - name of routine to create method context 1692e884886fSBarry Smith - routine_create - routine to create method context 1693e884886fSBarry Smith 1694e884886fSBarry Smith Level: developer 1695e884886fSBarry Smith 1696e884886fSBarry Smith Notes: 1697e884886fSBarry Smith MatMFFDRegisterDynamic() may be called multiple times to add several user-defined solvers. 1698e884886fSBarry Smith 1699e884886fSBarry Smith If dynamic libraries are used, then the fourth input argument (routine_create) 1700e884886fSBarry Smith is ignored. 1701e884886fSBarry Smith 1702e884886fSBarry Smith Sample usage: 1703e884886fSBarry Smith .vb 1704e884886fSBarry Smith MatMFFDRegisterDynamic("my_h",/home/username/my_lib/lib/libO/solaris/mylib.a, 1705e884886fSBarry Smith "MyHCreate",MyHCreate); 1706e884886fSBarry Smith .ve 1707e884886fSBarry Smith 1708e884886fSBarry Smith Then, your solver can be chosen with the procedural interface via 1709e884886fSBarry Smith $ MatMFFDSetType(mfctx,"my_h") 1710e884886fSBarry Smith or at runtime via the option 1711e884886fSBarry Smith $ -snes_mf_type my_h 1712e884886fSBarry Smith 1713e884886fSBarry Smith .keywords: MatMFFD, register 1714e884886fSBarry Smith 1715e884886fSBarry Smith .seealso: MatMFFDRegisterAll(), MatMFFDRegisterDestroy() 1716e884886fSBarry Smith M*/ 1717e884886fSBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 1718e884886fSBarry Smith #define MatMFFDRegisterDynamic(a,b,c,d) MatMFFDRegister(a,b,c,0) 1719e884886fSBarry Smith #else 1720e884886fSBarry Smith #define MatMFFDRegisterDynamic(a,b,c,d) MatMFFDRegister(a,b,c,d) 1721e884886fSBarry Smith #endif 1722e884886fSBarry Smith 17237087cfbeSBarry Smith extern PetscErrorCode MatMFFDRegisterAll(const char[]); 17247087cfbeSBarry Smith extern PetscErrorCode MatMFFDRegisterDestroy(void); 17257087cfbeSBarry Smith extern PetscErrorCode MatMFFDDSSetUmin(Mat,PetscReal); 17267087cfbeSBarry Smith extern PetscErrorCode MatMFFDWPSetComputeNormU(Mat,PetscBool ); 1727e884886fSBarry Smith 1728e884886fSBarry Smith 17297087cfbeSBarry Smith extern PetscErrorCode PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *); 17307087cfbeSBarry Smith extern PetscErrorCode PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *); 17317dbadf16SMatthew Knepley 173297969023SHong Zhang /* 173397969023SHong Zhang PETSc interface to MUMPS 173497969023SHong Zhang */ 173597969023SHong Zhang #ifdef PETSC_HAVE_MUMPS 1736f250808bSBarry Smith extern PetscErrorCode MatMumpsSetIcntl(Mat,PetscInt,PetscInt); 173797969023SHong Zhang #endif 173823a5497aSJed Brown 1739d954db57SHong Zhang /* 1740d954db57SHong Zhang PETSc interface to SUPERLU 1741d954db57SHong Zhang */ 1742d954db57SHong Zhang #ifdef PETSC_HAVE_SUPERLU 1743f250808bSBarry Smith extern PetscErrorCode MatSuperluSetILUDropTol(Mat,PetscReal); 1744d954db57SHong Zhang #endif 1745d954db57SHong Zhang 174690c53902SBarry Smith #if defined(PETSC_HAVE_CUSP) 174790c53902SBarry Smith extern PetscErrorCode MatCreateSeqAIJCUSP(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 174890c53902SBarry Smith extern PetscErrorCode MatCreateMPIAIJCUSP(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 174990c53902SBarry Smith #endif 175090c53902SBarry Smith 175154efbe56SHong Zhang /* 175254efbe56SHong Zhang PETSc interface to FFTW 175354efbe56SHong Zhang */ 175454efbe56SHong Zhang #if defined(PETSC_HAVE_FFTW) 175574a26884SAmlan Barua extern PetscErrorCode VecScatterPetscToFFTW(Mat,Vec,Vec); 175674a26884SAmlan Barua extern PetscErrorCode VecScatterFFTWToPetsc(Mat,Vec,Vec); 17574be45526SHong Zhang extern PetscErrorCode MatGetVecsFFTW(Mat,Vec*,Vec*,Vec*); 175854efbe56SHong Zhang #endif 175954efbe56SHong Zhang 1760c8883902SJed Brown extern PetscErrorCode MatCreateNest(MPI_Comm,PetscInt,const IS[],PetscInt,const IS[],const Mat[],Mat*); 1761c8883902SJed Brown extern PetscErrorCode MatNestGetSize(Mat,PetscInt*,PetscInt*); 1762c8883902SJed Brown extern PetscErrorCode MatNestGetSubMats(Mat,PetscInt*,PetscInt*,Mat***); 1763c8883902SJed Brown extern PetscErrorCode MatNestGetSubMat(Mat,PetscInt,PetscInt,Mat*); 17647087cfbeSBarry Smith extern PetscErrorCode MatNestSetVecType(Mat,const VecType); 1765c8883902SJed Brown extern PetscErrorCode MatNestSetSubMats(Mat,PetscInt,const IS[],PetscInt,const IS[],const Mat[]); 17660782ca92SJed Brown extern PetscErrorCode MatNestSetSubMat(Mat,PetscInt,PetscInt,Mat); 1767d8588912SDave May 176823a5497aSJed Brown PETSC_EXTERN_CXX_END 176923a5497aSJed Brown #endif 1770