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 27d9274352SBarry Smith .seealso: MatSetType(), Mat 28d91e6319SBarry Smith E*/ 29273d9f13SBarry Smith #define MATSAME "same" 30273d9f13SBarry Smith #define MATSEQMAIJ "seqmaij" 31273d9f13SBarry Smith #define MATMPIMAIJ "mpimaij" 32209238afSKris Buschelman #define MATMAIJ "maij" 33273d9f13SBarry Smith #define MATIS "is" 34273d9f13SBarry Smith #define MATMPIROWBS "mpirowbs" 35273d9f13SBarry Smith #define MATSEQAIJ "seqaij" 36273d9f13SBarry Smith #define MATMPIAIJ "mpiaij" 37209238afSKris Buschelman #define MATAIJ "aij" 38273d9f13SBarry Smith #define MATSHELL "shell" 39273d9f13SBarry Smith #define MATSEQBDIAG "seqbdiag" 40273d9f13SBarry Smith #define MATMPIBDIAG "mpibdiag" 41209238afSKris Buschelman #define MATBDIAG "bdiag" 42209238afSKris Buschelman #define MATSEQDENSE "seqdense" 43273d9f13SBarry Smith #define MATMPIDENSE "mpidense" 44209238afSKris Buschelman #define MATDENSE "dense" 45273d9f13SBarry Smith #define MATSEQBAIJ "seqbaij" 46273d9f13SBarry Smith #define MATMPIBAIJ "mpibaij" 47209238afSKris Buschelman #define MATBAIJ "baij" 48273d9f13SBarry Smith #define MATMPIADJ "mpiadj" 49273d9f13SBarry Smith #define MATSEQSBAIJ "seqsbaij" 50273d9f13SBarry Smith #define MATMPISBAIJ "mpisbaij" 51209238afSKris Buschelman #define MATSBAIJ "sbaij" 52cebc7f6cSBarry Smith #define MATDAAD "daad" 53cebc7f6cSBarry Smith #define MATMFFD "mffd" 54c8a8475eSBarry Smith #define MATNORMAL "normal" 55ab92ecdeSBarry Smith #define MATLRC "lrc" 56b3a1e11cSKris Buschelman #define MATSEQAIJSPOOLES "seqaijspooles" 57d10c748bSKris Buschelman #define MATMPIAIJSPOOLES "mpiaijspooles" 589abb65ffSKris Buschelman #define MATSEQSBAIJSPOOLES "seqsbaijspooles" 5922191285SKris Buschelman #define MATMPISBAIJSPOOLES "mpisbaijspooles" 60bb4d4166SHong Zhang #define MATAIJSPOOLES "aijspooles" 61bb4d4166SHong Zhang #define MATSBAIJSPOOLES "sbaijspooles" 6214b396bbSKris Buschelman #define MATSUPERLU "superlu" 631d96aa28SKris Buschelman #define MATSUPERLU_DIST "superlu_dist" 641677d0b8SKris Buschelman #define MATUMFPACK "umfpack" 65e8aa55a4SKris Buschelman #define MATESSL "essl" 664eb8e494SKris Buschelman #define MATLUSOL "lusol" 67397b6df1SKris Buschelman #define MATAIJMUMPS "aijmumps" 68397b6df1SKris Buschelman #define MATSBAIJMUMPS "sbaijmumps" 698da957c5SKris Buschelman #define MATDSCPACK "dscpack" 707065e2aaSKris Buschelman #define MATMATLAB "matlab" 714b8d99adSRichard Tran Mills #define MATSEQCSRPERM "seqcsrperm" 7218def556SRichard Tran Mills #define MATMPICSRPERM "mpicsrperm" 73814655b8SBarry Smith #define MATCSRPERM "csrperm" 7481824310SBarry Smith #define MATSEQCRL "seqcrl" 75c02b24c5SBarry Smith #define MATMPICRL "mpicrl" 76c02b24c5SBarry Smith #define MATCRL "crl" 77c0aa2d19SHong Zhang #define MATPLAPACK "plapack" 78f69a0ea3SMatthew Knepley #define MatType const char* 79d91e6319SBarry Smith 80c06d978dSMatthew Knepley /* Logging support */ 81552e946dSBarry Smith #define MAT_FILE_COOKIE 1211216 /* used to indicate matrices in binary files */ 82be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_COOKIE; 83be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_FDCOLORING_COOKIE; 84be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_PARTITIONING_COOKIE; 85be1d678aSKris Buschelman extern PetscCookie PETSCMAT_DLLEXPORT MAT_NULLSPACE_COOKIE; 866849ba73SBarry Smith extern PetscEvent MAT_Mult, MAT_MultMatrixFree, MAT_Mults, MAT_MultConstrained, MAT_MultAdd, MAT_MultTranspose; 876849ba73SBarry Smith extern PetscEvent MAT_MultTransposeConstrained, MAT_MultTransposeAdd, MAT_Solve, MAT_Solves, MAT_SolveAdd, MAT_SolveTranspose; 886849ba73SBarry Smith extern PetscEvent MAT_SolveTransposeAdd, MAT_Relax, MAT_ForwardSolve, MAT_BackwardSolve, MAT_LUFactor, MAT_LUFactorSymbolic; 896849ba73SBarry Smith extern PetscEvent MAT_LUFactorNumeric, MAT_CholeskyFactor, MAT_CholeskyFactorSymbolic, MAT_CholeskyFactorNumeric, MAT_ILUFactor; 906849ba73SBarry Smith extern PetscEvent MAT_ILUFactorSymbolic, MAT_ICCFactorSymbolic, MAT_Copy, MAT_Convert, MAT_Scale, MAT_AssemblyBegin; 916849ba73SBarry Smith extern PetscEvent MAT_AssemblyEnd, MAT_SetValues, MAT_GetValues, MAT_GetRow, MAT_GetSubMatrices, MAT_GetColoring, MAT_GetOrdering; 926849ba73SBarry Smith extern PetscEvent MAT_IncreaseOverlap, MAT_Partitioning, MAT_ZeroEntries, MAT_Load, MAT_View, MAT_AXPY, MAT_FDColoringCreate; 936849ba73SBarry Smith extern PetscEvent MAT_FDColoringApply, MAT_Transpose, MAT_FDColoringFunction; 94cb00f407SKris Buschelman extern PetscEvent MAT_MatMult, MAT_MatMultSymbolic, MAT_MatMultNumeric; 95534c1384SKris Buschelman extern PetscEvent MAT_PtAP, MAT_PtAPSymbolic, MAT_PtAPNumeric; 96bc011b1eSHong Zhang extern PetscEvent MAT_MatMultTranspose, MAT_MatMultTransposeSymbolic, MAT_MatMultTransposeNumeric; 97c06d978dSMatthew Knepley 98ceb03754SKris Buschelman /*E 99ceb03754SKris Buschelman MatReuse - Indicates if matrices obtained from a previous call to MatGetSubMatrices() 100ceb03754SKris Buschelman or MatGetSubMatrix() are to be reused to store the new matrix values. 101ceb03754SKris Buschelman 102ceb03754SKris Buschelman Level: beginner 103ceb03754SKris Buschelman 104ceb03754SKris Buschelman Any additions/changes here MUST also be made in include/finclude/petscmat.h 105ceb03754SKris Buschelman 1060c451bc4SBarry Smith .seealso: MatGetSubMatrices(), MatGetSubMatrix(), MatDestroyMatrices(), MatConvert() 107ceb03754SKris Buschelman E*/ 108ceb03754SKris Buschelman typedef enum {MAT_INITIAL_MATRIX,MAT_REUSE_MATRIX} MatReuse; 109ceb03754SKris Buschelman 110be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInitializePackage(char *); 111c06d978dSMatthew Knepley 112f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreate(MPI_Comm,Mat*); 113a69a960dSMatthew Knepley PetscPolymorphicFunction(MatCreate,(MPI_Comm comm),(comm,&A),Mat,A) 114a69a960dSMatthew Knepley PetscPolymorphicFunction(MatCreate,(),(PETSC_COMM_WORLD,&A),Mat,A) 115f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetSizes(Mat,PetscInt,PetscInt,PetscInt,PetscInt); 116f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetType(Mat,MatType); 117be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetFromOptions(Mat); 118be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUpPreallocation(Mat); 119be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRegisterAll(const char[]); 120be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat)); 12130de9b25SBarry Smith 122f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetOptionsPrefix(Mat,const char[]); 123f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAppendOptionsPrefix(Mat,const char[]); 124f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOptionsPrefix(Mat,const char*[]); 125f69a0ea3SMatthew Knepley 12630de9b25SBarry Smith /*MC 12730de9b25SBarry Smith MatRegisterDynamic - Adds a new matrix type 12830de9b25SBarry Smith 12930de9b25SBarry Smith Synopsis: 130c1ac3661SBarry Smith PetscErrorCode MatRegisterDynamic(char *name,char *path,char *name_create,PetscErrorCode (*routine_create)(Mat)) 13130de9b25SBarry Smith 13230de9b25SBarry Smith Not Collective 13330de9b25SBarry Smith 13430de9b25SBarry Smith Input Parameters: 13530de9b25SBarry Smith + name - name of a new user-defined matrix type 13630de9b25SBarry Smith . path - path (either absolute or relative) the library containing this solver 13730de9b25SBarry Smith . name_create - name of routine to create method context 13830de9b25SBarry Smith - routine_create - routine to create method context 13930de9b25SBarry Smith 14030de9b25SBarry Smith Notes: 14130de9b25SBarry Smith MatRegisterDynamic() may be called multiple times to add several user-defined solvers. 14230de9b25SBarry Smith 14330de9b25SBarry Smith If dynamic libraries are used, then the fourth input argument (routine_create) 14430de9b25SBarry Smith is ignored. 14530de9b25SBarry Smith 14630de9b25SBarry Smith Sample usage: 14730de9b25SBarry Smith .vb 14830de9b25SBarry Smith MatRegisterDynamic("my_mat",/home/username/my_lib/lib/libO/solaris/mylib.a, 14930de9b25SBarry Smith "MyMatCreate",MyMatCreate); 15030de9b25SBarry Smith .ve 15130de9b25SBarry Smith 15230de9b25SBarry Smith Then, your solver can be chosen with the procedural interface via 15330de9b25SBarry Smith $ MatSetType(Mat,"my_mat") 15430de9b25SBarry Smith or at runtime via the option 15530de9b25SBarry Smith $ -mat_type my_mat 15630de9b25SBarry Smith 15730de9b25SBarry Smith Level: advanced 15830de9b25SBarry Smith 159ab901514SBarry Smith Notes: ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values. 16030de9b25SBarry Smith If your function is not being put into a shared library then use VecRegister() instead 16130de9b25SBarry Smith 16230de9b25SBarry Smith .keywords: Mat, register 16330de9b25SBarry Smith 16430de9b25SBarry Smith .seealso: MatRegisterAll(), MatRegisterDestroy() 16530de9b25SBarry Smith 16630de9b25SBarry Smith M*/ 167273d9f13SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 168273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,0) 169273d9f13SBarry Smith #else 170273d9f13SBarry Smith #define MatRegisterDynamic(a,b,c,d) MatRegister(a,b,c,d) 17130de9b25SBarry Smith #endif 17230de9b25SBarry Smith 173273d9f13SBarry Smith extern PetscTruth MatRegisterAllCalled; 174b0a32e0cSBarry Smith extern PetscFList MatList; 17528988994SBarry Smith 176be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqDense(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*); 177be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIDense(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*); 178be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 179ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,nz,nnz,&A),Mat,A) 180ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,&A),Mat,A) 181ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,m,n,0,nnz,&A),Mat,A) 182ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqAIJ,(PetscInt m,PetscInt n),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,&A),Mat,A) 183ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,m,n,nz,PETSC_NULL,A)) 184ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,m,n,0,nnz,A)) 185ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqAIJ,(PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,m,n,0,PETSC_NULL,A)) 186be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 187ba966639SSatish 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) 188ba966639SSatish 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) 189ba966639SSatish 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) 190ba966639SSatish 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) 191ba966639SSatish 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)) 192ba966639SSatish 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)) 193ba966639SSatish 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)) 194ba966639SSatish 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) 195ba966639SSatish 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) 196ba966639SSatish 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) 197ba966639SSatish 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) 198ba966639SSatish 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)) 199ba966639SSatish 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)) 200ba966639SSatish 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)) 201be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIRowbs(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 202be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqBDiag(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscScalar*[],Mat*); 203be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIBDiag(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscScalar*[],Mat*); 204be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 205ba966639SSatish 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) 206ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A) 207ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A) 208ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A) 209ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A)) 210ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A)) 211ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A)) 212be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 213ba966639SSatish 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) 214ba966639SSatish 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) 215ba966639SSatish 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) 216ba966639SSatish 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) 217ba966639SSatish 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)) 218ba966639SSatish 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)) 219ba966639SSatish 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)) 220ba966639SSatish 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) 221ba966639SSatish 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) 222ba966639SSatish 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) 223ba966639SSatish 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) 224ba966639SSatish 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)) 225ba966639SSatish 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)) 226ba966639SSatish 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)) 227be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPIAdj(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscInt[],Mat*); 228be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 229ba966639SSatish 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) 230ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,&A),Mat,A) 231ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[]),(PETSC_COMM_SELF,bs,m,n,0,nnz,&A),Mat,A) 232ba966639SSatish Balay PetscPolymorphicFunction(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,&A),Mat,A) 233ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,PetscInt nz,Mat *A),(PETSC_COMM_SELF,bs,m,n,nz,PETSC_NULL,A)) 234ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,const PetscInt nnz[],Mat *A),(PETSC_COMM_SELF,bs,m,n,0,nnz,A)) 235ba966639SSatish Balay PetscPolymorphicSubroutine(MatCreateSeqSBAIJ,(PetscInt bs,PetscInt m,PetscInt n,Mat *A),(PETSC_COMM_SELF,bs,m,n,0,PETSC_NULL,A)) 236be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPISBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 237ba966639SSatish 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) 238ba966639SSatish 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) 239ba966639SSatish 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) 240ba966639SSatish 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) 241ba966639SSatish 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)) 242ba966639SSatish 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)) 243ba966639SSatish 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)) 244ba966639SSatish 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) 245ba966639SSatish 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) 246ba966639SSatish 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) 247ba966639SSatish 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) 248ba966639SSatish 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)) 249ba966639SSatish 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)) 250ba966639SSatish 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)) 251be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateShell(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,void *,Mat*); 252ba966639SSatish 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) 253ba966639SSatish Balay PetscPolymorphicFunction(MatCreateShell,(PetscInt m,PetscInt n,PetscInt M,PetscInt N,void *ctx),(PETSC_COMM_WORLD,m,n,M,N,ctx,&A),Mat,A) 254be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateAdic(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,void (*)(void),Mat*); 255be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateNormal(Mat,Mat*); 256ba966639SSatish Balay PetscPolymorphicFunction(MatCreateNormal,(Mat mat),(mat,&A),Mat,A) 257ab92ecdeSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateLRC(Mat,Mat,Mat,Mat*); 258284134d9SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateIS(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,ISLocalToGlobalMapping,Mat*); 2591472f72bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 2601472f72bSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMPICRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 2611472f72bSBarry Smith 262f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUp(Mat); 263be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroy(Mat); 26421c89e3eSBarry Smith 265be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPrintHelp(Mat); 266be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetPetscMaps(Mat,PetscMap*,PetscMap*); 267ec0117caSBarry Smith 2680c451bc4SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConjugate(Mat); 26999cafbc1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRealPart(Mat); 27099cafbc1SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatImaginaryPart(Mat); 27199cafbc1SBarry Smith 2728ed539a5SBarry Smith /* ------------------------------------------------------------*/ 273be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 274be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 27587d4246cSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesRow(Mat,PetscInt,const PetscScalar[]); 276f259bd47SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesRowLocal(Mat,PetscInt,const PetscScalar[]); 27784cb2905SBarry Smith 2782ef4de8bSBarry Smith /*S 2792ef4de8bSBarry Smith MatStencil - Data structure (C struct) for storing information about a single row or 2802ef4de8bSBarry Smith column of a matrix as index on an associated grid. 2812ef4de8bSBarry Smith 2822ef4de8bSBarry Smith Level: beginner 2832ef4de8bSBarry Smith 2842ef4de8bSBarry Smith Concepts: matrix; linear operator 2852ef4de8bSBarry Smith 286d7bd93baSBarry Smith .seealso: MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockStencil() 2872ef4de8bSBarry Smith S*/ 288435da068SBarry Smith typedef struct { 289c1ac3661SBarry Smith PetscInt k,j,i,c; 290435da068SBarry Smith } MatStencil; 2912ef4de8bSBarry Smith 292be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode); 293be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode); 294be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt); 295435da068SBarry Smith 296be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetColoring(Mat,ISColoring); 297be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesAdic(Mat,void*); 298be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesAdifor(Mat,PetscInt,void*); 2993a7fca6bSBarry Smith 300d91e6319SBarry Smith /*E 301d91e6319SBarry Smith MatAssemblyType - Indicates if the matrix is now to be used, or if you plan 302d91e6319SBarry Smith to continue to add values to it 303d91e6319SBarry Smith 304d91e6319SBarry Smith Level: beginner 305d91e6319SBarry Smith 306d91e6319SBarry Smith .seealso: MatAssemblyBegin(), MatAssemblyEnd() 307d91e6319SBarry Smith E*/ 3086d4a8577SBarry Smith typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType; 309be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssemblyBegin(Mat,MatAssemblyType); 310be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssemblyEnd(Mat,MatAssemblyType); 311be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAssembled(Mat,PetscTruth*); 3124f9c727eSBarry Smith 313be1d678aSKris Buschelman extern PETSCMAT_DLLEXPORT PetscInt MatSetValue_Row; 314be1d678aSKris Buschelman extern PETSCMAT_DLLEXPORT PetscInt MatSetValue_Column; 315be1d678aSKris Buschelman extern PETSCMAT_DLLEXPORT PetscScalar MatSetValue_Value; 3161d73ed98SBarry Smith 31730de9b25SBarry Smith /*MC 31830de9b25SBarry Smith MatSetValue - Set a single entry into a matrix. 31930de9b25SBarry Smith 32030de9b25SBarry Smith Synopsis: 321c1ac3661SBarry Smith PetscErrorCode MatSetValue(Mat m,PetscInt row,PetscInt col,PetscScalar value,InsertMode mode); 32230de9b25SBarry Smith 32330de9b25SBarry Smith Not collective 32430de9b25SBarry Smith 32530de9b25SBarry Smith Input Parameters: 32630de9b25SBarry Smith + m - the matrix 32730de9b25SBarry Smith . row - the row location of the entry 32830de9b25SBarry Smith . col - the column location of the entry 32930de9b25SBarry Smith . value - the value to insert 33030de9b25SBarry Smith - mode - either INSERT_VALUES or ADD_VALUES 33130de9b25SBarry Smith 33230de9b25SBarry Smith Notes: 33330de9b25SBarry Smith For efficiency one should use MatSetValues() and set several or many 33430de9b25SBarry Smith values simultaneously if possible. 33530de9b25SBarry Smith 33630de9b25SBarry Smith Level: beginner 33730de9b25SBarry Smith 33830de9b25SBarry Smith .seealso: MatSetValues(), MatSetValueLocal() 33930de9b25SBarry Smith M*/ 340b951964fSBarry Smith #define MatSetValue(v,i,j,va,mode) \ 341f1144a30SSatish Balay ((MatSetValue_Row = i,MatSetValue_Column = j,MatSetValue_Value = va,0) || \ 3421d73ed98SBarry Smith MatSetValues(v,1,&MatSetValue_Row,1,&MatSetValue_Column,&MatSetValue_Value,mode)) 34330de9b25SBarry Smith 344ea06a074SBarry Smith #define MatGetValue(v,i,j,va) \ 345f1144a30SSatish Balay ((MatSetValue_Row = i,MatSetValue_Column = j,0) || \ 3461d73ed98SBarry Smith MatGetValues(v,1,&MatSetValue_Row,1,&MatSetValue_Column,&va)) 34730de9b25SBarry Smith 348d91e6319SBarry Smith #define MatSetValueLocal(v,i,j,va,mode) \ 349f1144a30SSatish Balay ((MatSetValue_Row = i,MatSetValue_Column = j,MatSetValue_Value = va,0) || \ 3501d73ed98SBarry Smith MatSetValuesLocal(v,1,&MatSetValue_Row,1,&MatSetValue_Column,&MatSetValue_Value,mode)) 35130de9b25SBarry Smith 352d91e6319SBarry Smith /*E 353d91e6319SBarry Smith MatOption - Options that may be set for a matrix and its behavior or storage 354d91e6319SBarry Smith 355d91e6319SBarry Smith Level: beginner 356d91e6319SBarry Smith 3570a835dfdSSatish Balay Any additions/changes here MUST also be made in include/finclude/petscmat.h 358d91e6319SBarry Smith 359d91e6319SBarry Smith .seealso: MatSetOption() 360d91e6319SBarry Smith E*/ 3616d4a8577SBarry Smith typedef enum {MAT_ROW_ORIENTED=1,MAT_COLUMN_ORIENTED=2,MAT_ROWS_SORTED=4, 3626d4a8577SBarry Smith MAT_COLUMNS_SORTED=8,MAT_NO_NEW_NONZERO_LOCATIONS=16, 3636d4a8577SBarry Smith MAT_YES_NEW_NONZERO_LOCATIONS=32,MAT_SYMMETRIC=64, 3646ca9ecd3SBarry Smith MAT_STRUCTURALLY_SYMMETRIC=65,MAT_NO_NEW_DIAGONALS=66, 3656ca9ecd3SBarry Smith MAT_YES_NEW_DIAGONALS=67,MAT_INODE_LIMIT_1=68,MAT_INODE_LIMIT_2=69, 3666ca9ecd3SBarry Smith MAT_INODE_LIMIT_3=70,MAT_INODE_LIMIT_4=71,MAT_INODE_LIMIT_5=72, 3676ca9ecd3SBarry Smith MAT_IGNORE_OFF_PROC_ENTRIES=73,MAT_ROWS_UNSORTED=74, 3684787f768SSatish Balay MAT_COLUMNS_UNSORTED=75,MAT_NEW_NONZERO_LOCATION_ERR=76, 3697c922b88SBarry Smith MAT_NEW_NONZERO_ALLOCATION_ERR=77,MAT_USE_HASH_TABLE=78, 3702bad1931SBarry Smith MAT_KEEP_ZEROED_ROWS=79,MAT_IGNORE_ZERO_ENTRIES=80,MAT_USE_INODES=81, 371efcf0fc3SBarry Smith MAT_DO_NOT_USE_INODES=82,MAT_NOT_SYMMETRIC=83,MAT_HERMITIAN=84, 3729a4540c5SBarry Smith MAT_NOT_STRUCTURALLY_SYMMETRIC=85,MAT_NOT_HERMITIAN=86, 373d487561eSHong Zhang MAT_SYMMETRY_ETERNAL=87,MAT_NOT_SYMMETRY_ETERNAL=88, 374941593c8SHong Zhang MAT_USE_COMPRESSEDROW=89,MAT_DO_NOT_USE_COMPRESSEDROW=90, 375f5edf698SHong Zhang MAT_IGNORE_LOWER_TRIANGULAR=91,MAT_ERROR_LOWER_TRIANGULAR=92,MAT_GETROW_UPPERTRIANGULAR=93} MatOption; 376be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetOption(Mat,MatOption); 377be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetType(Mat,MatType*); 378ba966639SSatish Balay PetscPolymorphicFunction(MatGetType,(Mat mat),(mat,&t),MatType,t) 37984cb2905SBarry Smith 380be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]); 381be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 382be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 383f5edf698SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowUpperTriangular(Mat); 384f5edf698SHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRowUpperTriangular(Mat); 385be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 386be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreColumn(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 387be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumnVector(Mat,Vec,PetscInt); 388be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetArray(Mat,PetscScalar *[]); 389ba966639SSatish Balay PetscPolymorphicFunction(MatGetArray,(Mat mat),(mat,&a),PetscScalar*,a) 390be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreArray(Mat,PetscScalar *[]); 391be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBlockSize(Mat,PetscInt *); 392ba966639SSatish Balay PetscPolymorphicFunction(MatGetBlockSize,(Mat mat),(mat,&a),PetscInt,a) 393be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetBlockSize(Mat,PetscInt); 3947b80b807SBarry Smith 395be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMult(Mat,Vec,Vec); 396be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultAdd(Mat,Vec,Vec,Vec); 397ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultAdd,(Mat A,Vec x,Vec y),(A,x,y,y)) 398be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTranspose(Mat,Vec,Vec); 399be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsTranspose(Mat,Mat,PetscReal,PetscTruth*); 400ba966639SSatish Balay PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B,PetscReal tol),(A,B,tol,&t),PetscTruth,t) 401ba966639SSatish Balay PetscPolymorphicFunction(MatIsTranspose,(Mat A,Mat B),(A,B,0.0,&t),PetscTruth,t) 402be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeAdd(Mat,Vec,Vec,Vec); 403ba966639SSatish Balay PetscPolymorphicSubroutine(MatMultTransposeAdd,(Mat A,Vec x,Vec y),(A,x,y,y)) 404be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultConstrained(Mat,Vec,Vec); 405be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeConstrained(Mat,Vec,Vec); 406*2bebee5dSHong Zhang EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatSolve(Mat,Mat,Mat); 407f5edf698SHong Zhang 408d91e6319SBarry Smith /*E 409d91e6319SBarry Smith MatDuplicateOption - Indicates if a duplicated sparse matrix should have 410d91e6319SBarry Smith its numerical values copied over or just its nonzero structure. 411d91e6319SBarry Smith 412d91e6319SBarry Smith Level: beginner 413d91e6319SBarry Smith 414d91e6319SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 415d91e6319SBarry Smith 416d91e6319SBarry Smith .seealso: MatDuplicate() 417d91e6319SBarry Smith E*/ 4182e8a6d31SBarry Smith typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES} MatDuplicateOption; 4192e8a6d31SBarry Smith 420be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvertRegister(const char[],const char[],const char[],PetscErrorCode (*)(Mat,MatType,MatReuse,Mat*)); 421273d9f13SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 422273d9f13SBarry Smith #define MatConvertRegisterDynamic(a,b,c,d) MatConvertRegister(a,b,c,0) 423273d9f13SBarry Smith #else 424273d9f13SBarry Smith #define MatConvertRegisterDynamic(a,b,c,d) MatConvertRegister(a,b,c,d) 425273d9f13SBarry Smith #endif 426be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvertRegisterAll(const char[]); 427be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvertRegisterDestroy(void); 428273d9f13SBarry Smith extern PetscTruth MatConvertRegisterAllCalled; 429b0a32e0cSBarry Smith extern PetscFList MatConvertList; 430f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert(Mat,MatType,MatReuse,Mat*); 431ba966639SSatish Balay PetscPolymorphicFunction(MatConvert,(Mat A,MatType t),(A,t,MAT_INITIAL_MATRIX,&a),Mat,a) 432be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDuplicate(Mat,MatDuplicateOption,Mat*); 433ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A,MatDuplicateOption o),(A,o,&a),Mat,a) 434ba966639SSatish Balay PetscPolymorphicFunction(MatDuplicate,(Mat A),(A,MAT_COPY_VALUES,&a),Mat,a) 43594a9d846SBarry Smith 436d91e6319SBarry Smith /*E 437d91e6319SBarry Smith MatStructure - Indicates if the matrix has the same nonzero structure 438d91e6319SBarry Smith 439d91e6319SBarry Smith Level: beginner 440d91e6319SBarry Smith 441d91e6319SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 442d91e6319SBarry Smith 44394b7f48cSBarry Smith .seealso: MatCopy(), KSPSetOperators(), PCSetOperators() 444d91e6319SBarry Smith E*/ 445c537a176SHong Zhang typedef enum {SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER,SUBSET_NONZERO_PATTERN} MatStructure; 446cb5b572fSBarry Smith 447be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCopy(Mat,Mat,MatStructure); 448be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatView(Mat,PetscViewer); 449be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsSymmetric(Mat,PetscReal,PetscTruth*); 450ba966639SSatish Balay PetscPolymorphicFunction(MatIsSymmetric,(Mat A,PetscReal tol),(A,tol,&t),PetscTruth,t) 451ba966639SSatish Balay PetscPolymorphicFunction(MatIsSymmetric,(Mat A),(A,0.0,&t),PetscTruth,t) 452be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsStructurallySymmetric(Mat,PetscTruth*); 453ba966639SSatish Balay PetscPolymorphicFunction(MatIsStructurallySymmetric,(Mat A),(A,&t),PetscTruth,t) 454be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitian(Mat,PetscTruth*); 455ba966639SSatish Balay PetscPolymorphicFunction(MatIsHermitian,(Mat A),(A,&t),PetscTruth,t) 456be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsSymmetricKnown(Mat,PetscTruth*,PetscTruth*); 457be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIsHermitianKnown(Mat,PetscTruth*,PetscTruth*); 458f69a0ea3SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLoad(PetscViewer,MatType,Mat*); 459ba966639SSatish Balay PetscPolymorphicFunction(MatLoad,(PetscViewer v,MatType t),(v,t,&a),Mat,a) 4607b80b807SBarry Smith 461be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowIJ(Mat,PetscInt,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *); 462be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreRowIJ(Mat,PetscInt,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *); 463be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColumnIJ(Mat,PetscInt,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *); 464be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestoreColumnIJ(Mat,PetscInt,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *); 465d4fbbf0eSBarry Smith 466d91e6319SBarry Smith /*S 467d91e6319SBarry Smith MatInfo - Context of matrix information, used with MatGetInfo() 468d91e6319SBarry Smith 469d91e6319SBarry Smith In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE 470d91e6319SBarry Smith 471d91e6319SBarry Smith Level: intermediate 472d91e6319SBarry Smith 473d91e6319SBarry Smith Concepts: matrix^nonzero information 474d91e6319SBarry Smith 475d9274352SBarry Smith .seealso: MatGetInfo(), MatInfoType 476d91e6319SBarry Smith S*/ 4774e220ebcSLois Curfman McInnes typedef struct { 478b0a32e0cSBarry Smith PetscLogDouble rows_global,columns_global; /* number of global rows and columns */ 479b0a32e0cSBarry Smith PetscLogDouble rows_local,columns_local; /* number of local rows and columns */ 480b0a32e0cSBarry Smith PetscLogDouble block_size; /* block size */ 481b0a32e0cSBarry Smith PetscLogDouble nz_allocated,nz_used,nz_unneeded; /* number of nonzeros */ 482b0a32e0cSBarry Smith PetscLogDouble memory; /* memory allocated */ 483b0a32e0cSBarry Smith PetscLogDouble assemblies; /* number of matrix assemblies called */ 484b0a32e0cSBarry Smith PetscLogDouble mallocs; /* number of mallocs during MatSetValues() */ 485b0a32e0cSBarry Smith PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */ 486b0a32e0cSBarry Smith PetscLogDouble factor_mallocs; /* number of mallocs during factorization */ 4874e220ebcSLois Curfman McInnes } MatInfo; 4884e220ebcSLois Curfman McInnes 489d9274352SBarry Smith /*E 490d9274352SBarry Smith MatInfoType - Indicates if you want information about the local part of the matrix, 491d9274352SBarry Smith the entire parallel matrix or the maximum over all the local parts. 492d9274352SBarry Smith 493d9274352SBarry Smith Level: beginner 494d9274352SBarry Smith 495d9274352SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 496d9274352SBarry Smith 497d9274352SBarry Smith .seealso: MatGetInfo(), MatInfo 498d9274352SBarry Smith E*/ 4997b80b807SBarry Smith typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType; 500be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetInfo(Mat,MatInfoType,MatInfo*); 501be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatValid(Mat,PetscTruth*); 502be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetDiagonal(Mat,Vec); 503be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetRowMax(Mat,Vec); 504be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatTranspose(Mat,Mat*); 505ba966639SSatish Balay PetscPolymorphicFunction(MatTranspose,(Mat A),(A,&t),Mat,t) 506be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPermute(Mat,IS,IS,Mat *); 507ba966639SSatish Balay PetscPolymorphicFunction(MatPermute,(Mat A,IS is1,IS is2),(A,is1,is2,&t),Mat,t) 508be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPermuteSparsify(Mat,PetscInt,PetscReal,PetscReal,IS,IS,Mat *); 509be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalScale(Mat,Vec,Vec); 510be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalSet(Mat,Vec,InsertMode); 511be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatEqual(Mat,Mat,PetscTruth*); 512ba966639SSatish Balay PetscPolymorphicFunction(MatEqual,(Mat A,Mat B),(A,B,&t),PetscTruth,t) 513be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultEqual(Mat,Mat,PetscInt,PetscTruth*); 514be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultAddEqual(Mat,Mat,PetscInt,PetscTruth*); 515be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeEqual(Mat,Mat,PetscInt,PetscTruth*); 516be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscTruth*); 5177b80b807SBarry Smith 518be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNorm(Mat,NormType,PetscReal *); 519ba966639SSatish Balay PetscPolymorphicFunction(MatNorm,(Mat A,NormType t),(A,t,&n),PetscReal,n) 520be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroEntries(Mat); 521f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar); 522f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsIS(Mat,IS,PetscScalar); 523f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroColumns(Mat,PetscInt,const PetscInt [],const PetscScalar*); 524f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroColumnsIS(Mat,IS,const PetscScalar*); 5257b80b807SBarry Smith 526be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatUseScaledForm(Mat,PetscTruth); 527be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScaleSystem(Mat,Vec,Vec); 528be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatUnScaleSystem(Mat,Vec,Vec); 5295ef9f2a5SBarry Smith 530be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSize(Mat,PetscInt*,PetscInt*); 531be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalSize(Mat,PetscInt*,PetscInt*); 532be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOwnershipRange(Mat,PetscInt*,PetscInt*); 5337b80b807SBarry Smith 534be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]); 535be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroyMatrices(PetscInt,Mat *[]); 536be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrix(Mat,IS,IS,PetscInt,MatReuse,Mat *); 537abf0e977SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetSubMatrixRaw(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscInt,MatReuse,Mat *); 538be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge(MPI_Comm,Mat,PetscInt,MatReuse,Mat*); 539be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPI(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*); 540be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPISymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*); 541be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMerge_SeqsToMPINumeric(Mat,Mat); 542be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDestroy_MPIAIJ_SeqsToMPI(Mat); 543be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalMat(Mat,MatReuse,Mat*); 544be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*); 545be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,PetscInt*,Mat*); 546be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetBrowsOfAoCols(Mat,Mat,MatReuse,PetscInt**,PetscScalar**,Mat*); 54743eb5e2fSMatthew Knepley #if defined (PETSC_USE_CTABLE) 548cd116e26SSatish Balay #include "petscctable.h" 54943eb5e2fSMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetCommunicationStructs(Mat, Vec *, PetscTable *, VecScatter *); 55043eb5e2fSMatthew Knepley #else 55143eb5e2fSMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetCommunicationStructs(Mat, Vec *, PetscInt *[], VecScatter *); 55243eb5e2fSMatthew Knepley #endif 5538efafbd8SBarry Smith 554be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt); 5557b80b807SBarry Smith 556be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*); 557be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultSymbolic(Mat,Mat,PetscReal,Mat*); 558be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultNumeric(Mat,Mat,Mat); 55922440eb1SKris Buschelman 560be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*); 561be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAPSymbolic(Mat,Mat,PetscReal,Mat*); 562be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPtAPNumeric(Mat,Mat,Mat); 56322440eb1SKris Buschelman 564be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTranspose(Mat,Mat,MatReuse,PetscReal,Mat*); 565be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTransposeSymbolic(Mat,Mat,PetscReal,Mat*); 566be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMatMultTransposeNumeric(Mat,Mat,Mat); 567bc011b1eSHong Zhang 568f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAXPY(Mat,PetscScalar,Mat,MatStructure); 569f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAYPX(Mat,PetscScalar,Mat); 570be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCompress(Mat); 5711c741599SHong Zhang 572f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatScale(Mat,PetscScalar); 573f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShift(Mat,PetscScalar); 5747b80b807SBarry Smith 575be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping); 576be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetLocalToGlobalMappingBlock(Mat,ISLocalToGlobalMapping); 577f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar); 578f4df32b1SMatthew Knepley EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatZeroRowsLocalIS(Mat,IS,PetscScalar); 579be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 580be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 581052efed2SBarry Smith 582be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStashSetInitialSize(Mat,PetscInt,PetscInt); 583be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*); 58490f02eecSBarry Smith 585be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInterpolate(Mat,Vec,Vec); 5860c451bc4SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInterpolateAdd(Mat,Vec,Vec,Vec); 587ba966639SSatish Balay PetscPolymorphicSubroutine(MatInterpolateAdd,(Mat A,Vec x,Vec y),(A,x,y,y)) 588be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRestrict(Mat,Vec,Vec); 589be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetVecs(Mat,Vec*,Vec*); 590bd481603SBarry Smith 591bd481603SBarry Smith /*MC 592bd481603SBarry Smith MatPreallocInitialize - Begins the block of code that will count the number of nonzeros per 593bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 594bd481603SBarry Smith 595bd481603SBarry Smith Synopsis: 596c1ac3661SBarry Smith PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz) 597bd481603SBarry Smith 598bd481603SBarry Smith Collective on MPI_Comm 599bd481603SBarry Smith 600bd481603SBarry Smith Input Parameters: 601bd481603SBarry Smith + comm - the communicator that will share the eventually allocated matrix 602859fcb39SBarry Smith . nrows - the number of LOCAL rows in the matrix 603859fcb39SBarry Smith - ncols - the number of LOCAL columns in the matrix 604bd481603SBarry Smith 605bd481603SBarry Smith Output Parameters: 606bd481603SBarry Smith + dnz - the array that will be passed to the matrix preallocation routines 607bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 608bd481603SBarry Smith 609bd481603SBarry Smith 610bd481603SBarry Smith Level: intermediate 611bd481603SBarry Smith 612bd481603SBarry Smith Notes: 613bd481603SBarry Smith See the chapter in the users manual on performance for more details 614bd481603SBarry Smith 6151d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 616bd481603SBarry Smith 617bd481603SBarry Smith Use MatPreallocateInitializeSymmetric() for symmetric matrices (MPISBAIJ matrices) 618bd481603SBarry Smith 619bd481603SBarry Smith Concepts: preallocation^Matrix 620bd481603SBarry Smith 621bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(), 622bd481603SBarry Smith MatPreallocateInitializeSymmetric(), MatPreallocateSymmetricSetLocal() 623bd481603SBarry Smith M*/ 624c4f061fbSSatish Balay #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \ 6257c922b88SBarry Smith { \ 626c1ac3661SBarry Smith PetscErrorCode _4_ierr; PetscInt __tmp = (nrows),__ctmp = (ncols),__rstart,__start,__end; \ 627c1ac3661SBarry Smith _4_ierr = PetscMalloc(2*__tmp*sizeof(PetscInt),&dnz);CHKERRQ(_4_ierr);onz = dnz + __tmp;\ 628c1ac3661SBarry Smith _4_ierr = PetscMemzero(dnz,2*__tmp*sizeof(PetscInt));CHKERRQ(_4_ierr);\ 629c1ac3661SBarry Smith _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ctmp;\ 630c1ac3661SBarry Smith _4_ierr = MPI_Scan(&__tmp,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __tmp; 6317c922b88SBarry Smith 632bd481603SBarry Smith /*MC 633bd481603SBarry Smith MatPreallocSymmetricInitialize - Begins the block of code that will count the number of nonzeros per 634bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 635bd481603SBarry Smith 636bd481603SBarry Smith Synopsis: 637c1ac3661SBarry Smith PetscErrorCode MatPreallocateSymmetricInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz) 638bd481603SBarry Smith 639bd481603SBarry Smith Collective on MPI_Comm 640bd481603SBarry Smith 641bd481603SBarry Smith Input Parameters: 642bd481603SBarry Smith + comm - the communicator that will share the eventually allocated matrix 643859fcb39SBarry Smith . nrows - the number of LOCAL rows in the matrix 644859fcb39SBarry Smith - ncols - the number of LOCAL columns in the matrix 645bd481603SBarry Smith 646bd481603SBarry Smith Output Parameters: 647bd481603SBarry Smith + dnz - the array that will be passed to the matrix preallocation routines 648bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 649bd481603SBarry Smith 650bd481603SBarry Smith 651bd481603SBarry Smith Level: intermediate 652bd481603SBarry Smith 653bd481603SBarry Smith Notes: 654bd481603SBarry Smith See the chapter in the users manual on performance for more details 655bd481603SBarry Smith 6561d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 657bd481603SBarry Smith 658bd481603SBarry Smith Concepts: preallocation^Matrix 659bd481603SBarry Smith 660bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(), 661bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal() 662bd481603SBarry Smith M*/ 663222b16d4SBarry Smith #define MatPreallocateSymmetricInitialize(comm,nrows,ncols,dnz,onz) 0; \ 664222b16d4SBarry Smith { \ 665c1ac3661SBarry Smith PetscErrorCode _4_ierr; PetscInt __tmp = (nrows),__ctmp = (ncols),__rstart,__end; \ 666c1ac3661SBarry Smith _4_ierr = PetscMalloc(2*__tmp*sizeof(PetscInt),&dnz);CHKERRQ(_4_ierr);onz = dnz + __tmp;\ 667c1ac3661SBarry Smith _4_ierr = PetscMemzero(dnz,2*__tmp*sizeof(PetscInt));CHKERRQ(_4_ierr);\ 668c1ac3661SBarry Smith _4_ierr = MPI_Scan(&__ctmp,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr);\ 669c1ac3661SBarry Smith _4_ierr = MPI_Scan(&__tmp,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __tmp; 670222b16d4SBarry Smith 671bd481603SBarry Smith /*MC 672bd481603SBarry Smith MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be 673bd481603SBarry Smith inserted using a local number of the rows and columns 674bd481603SBarry Smith 675bd481603SBarry Smith Synopsis: 676c1ac3661SBarry Smith PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 677bd481603SBarry Smith 678bd481603SBarry Smith Not Collective 679bd481603SBarry Smith 680bd481603SBarry Smith Input Parameters: 681bd481603SBarry Smith + map - the mapping between local numbering and global numbering 682bd481603SBarry Smith . nrows - the number of rows indicated 6831d73ed98SBarry Smith . rows - the indices of the rows 684bd481603SBarry Smith . ncols - the number of columns in the matrix 685bd481603SBarry Smith . cols - the columns indicated 686bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 687bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 688bd481603SBarry Smith 689bd481603SBarry Smith 690bd481603SBarry Smith Level: intermediate 691bd481603SBarry Smith 692bd481603SBarry Smith Notes: 693bd481603SBarry Smith See the chapter in the users manual on performance for more details 694bd481603SBarry Smith 6951d73ed98SBarry Smith Do not malloc or free dnz and onz, that is handled internally by these routines 696bd481603SBarry Smith 697bd481603SBarry Smith Concepts: preallocation^Matrix 698bd481603SBarry Smith 699bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(), 700bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal() 701bd481603SBarry Smith M*/ 702c4f061fbSSatish Balay #define MatPreallocateSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\ 703c4f061fbSSatish Balay {\ 704c1ac3661SBarry Smith PetscInt __l;\ 705ef66eb69SBarry Smith _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\ 706ef66eb69SBarry Smith _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\ 707c4f061fbSSatish Balay for (__l=0;__l<nrows;__l++) {\ 708ef66eb69SBarry Smith _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 709c4f061fbSSatish Balay }\ 710c4f061fbSSatish Balay } 711c4f061fbSSatish Balay 712bd481603SBarry Smith /*MC 713bd481603SBarry Smith MatPreallocateSymmetricSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be 714bd481603SBarry Smith inserted using a local number of the rows and columns 715bd481603SBarry Smith 716bd481603SBarry Smith Synopsis: 717c1ac3661SBarry Smith PetscErrorCode MatPreallocateSymmetricSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 718bd481603SBarry Smith 719bd481603SBarry Smith Not Collective 720bd481603SBarry Smith 721bd481603SBarry Smith Input Parameters: 722bd481603SBarry Smith + map - the mapping between local numbering and global numbering 723bd481603SBarry Smith . nrows - the number of rows indicated 7241d73ed98SBarry Smith . rows - the indices of the rows 725bd481603SBarry Smith . ncols - the number of columns in the matrix 726bd481603SBarry Smith . cols - the columns indicated 727bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 728bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 729bd481603SBarry Smith 730bd481603SBarry Smith 731bd481603SBarry Smith Level: intermediate 732bd481603SBarry Smith 733bd481603SBarry Smith Notes: 734bd481603SBarry Smith See the chapter in the users manual on performance for more details 735bd481603SBarry Smith 736bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 737bd481603SBarry Smith 738bd481603SBarry Smith Concepts: preallocation^Matrix 739bd481603SBarry Smith 740bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(), 741bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal() 742bd481603SBarry Smith M*/ 743d3d32019SBarry Smith #define MatPreallocateSymmetricSetLocal(map,nrows,rows,ncols,cols,dnz,onz) 0;\ 744d3d32019SBarry Smith {\ 745c1ac3661SBarry Smith PetscInt __l;\ 746d3d32019SBarry Smith _4_ierr = ISLocalToGlobalMappingApply(map,nrows,rows,rows);CHKERRQ(_4_ierr);\ 747d3d32019SBarry Smith _4_ierr = ISLocalToGlobalMappingApply(map,ncols,cols,cols);CHKERRQ(_4_ierr);\ 748d3d32019SBarry Smith for (__l=0;__l<nrows;__l++) {\ 749d3d32019SBarry Smith _4_ierr = MatPreallocateSymmetricSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 750d3d32019SBarry Smith }\ 751d3d32019SBarry Smith } 752d3d32019SBarry Smith 753bd481603SBarry Smith /*MC 754bd481603SBarry Smith MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be 755bd481603SBarry Smith inserted using a local number of the rows and columns 756bd481603SBarry Smith 757bd481603SBarry Smith Synopsis: 758c1ac3661SBarry Smith PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 759bd481603SBarry Smith 760bd481603SBarry Smith Not Collective 761bd481603SBarry Smith 762bd481603SBarry Smith Input Parameters: 763bd481603SBarry Smith + nrows - the number of rows indicated 7641d73ed98SBarry Smith . rows - the indices of the rows 765bd481603SBarry Smith . ncols - the number of columns in the matrix 766a50f8bf6SHong Zhang - cols - the columns indicated 767a50f8bf6SHong Zhang 768a50f8bf6SHong Zhang Output Parameters: 769a50f8bf6SHong Zhang + dnz - the array that will be passed to the matrix preallocation routines 770bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 771bd481603SBarry Smith 772bd481603SBarry Smith 773bd481603SBarry Smith Level: intermediate 774bd481603SBarry Smith 775bd481603SBarry Smith Notes: 776bd481603SBarry Smith See the chapter in the users manual on performance for more details 777bd481603SBarry Smith 778bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 779bd481603SBarry Smith 780bd481603SBarry Smith Concepts: preallocation^Matrix 781bd481603SBarry Smith 782bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(), 783bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal() 784bd481603SBarry Smith M*/ 785c4f061fbSSatish Balay #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\ 786c1ac3661SBarry Smith { PetscInt __i; \ 7877c922b88SBarry Smith for (__i=0; __i<nc; __i++) {\ 7887c922b88SBarry Smith if (cols[__i] < __start || cols[__i] >= __end) onz[row - __rstart]++; \ 7897c922b88SBarry Smith }\ 7907c922b88SBarry Smith dnz[row - __rstart] = nc - onz[row - __rstart];\ 7917c922b88SBarry Smith } 7927c922b88SBarry Smith 793bd481603SBarry Smith /*MC 794bd481603SBarry Smith MatPreallocateSymmetricSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be 795bd481603SBarry Smith inserted using a local number of the rows and columns 796bd481603SBarry Smith 797bd481603SBarry Smith Synopsis: 798c1ac3661SBarry Smith PetscErrorCode MatPreallocateSymmetricSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 799bd481603SBarry Smith 800bd481603SBarry Smith Not Collective 801bd481603SBarry Smith 802bd481603SBarry Smith Input Parameters: 803bd481603SBarry Smith + nrows - the number of rows indicated 8041d73ed98SBarry Smith . rows - the indices of the rows 805bd481603SBarry Smith . ncols - the number of columns in the matrix 806bd481603SBarry Smith . cols - the columns indicated 807bd481603SBarry Smith . dnz - the array that will be passed to the matrix preallocation routines 808bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 809bd481603SBarry Smith 810bd481603SBarry Smith 811bd481603SBarry Smith Level: intermediate 812bd481603SBarry Smith 813bd481603SBarry Smith Notes: 814bd481603SBarry Smith See the chapter in the users manual on performance for more details 815bd481603SBarry Smith 816bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 817bd481603SBarry Smith 818bd481603SBarry Smith Concepts: preallocation^Matrix 819bd481603SBarry Smith 820bd481603SBarry Smith .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateInitialize(), 821bd481603SBarry Smith MatPreallocateInitialize(), MatPreallocateSymmetricSetLocal(), MatPreallocateSetLocal() 822bd481603SBarry Smith M*/ 823d3d32019SBarry Smith #define MatPreallocateSymmetricSet(row,nc,cols,dnz,onz) 0;\ 824c1ac3661SBarry Smith { PetscInt __i; \ 825d3d32019SBarry Smith for (__i=0; __i<nc; __i++) {\ 826d3d32019SBarry Smith if (cols[__i] >= __end) onz[row - __rstart]++; \ 827d3d32019SBarry Smith else if (cols[__i] >= row) dnz[row - __rstart]++;\ 828d3d32019SBarry Smith }\ 829d3d32019SBarry Smith } 830d3d32019SBarry Smith 831bd481603SBarry Smith /*MC 832bd481603SBarry Smith MatPreallocFinalize - Ends the block of code that will count the number of nonzeros per 833bd481603SBarry Smith row in a matrix providing the data that one can use to correctly preallocate the matrix. 834bd481603SBarry Smith 835bd481603SBarry Smith Synopsis: 836c1ac3661SBarry Smith PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz) 837bd481603SBarry Smith 838bd481603SBarry Smith Collective on MPI_Comm 839bd481603SBarry Smith 840bd481603SBarry Smith Input Parameters: 841bd481603SBarry Smith + dnz - the array that will be passed to the matrix preallocation routines 842bd481603SBarry Smith - ozn - the other array passed to the matrix preallocation routines 843bd481603SBarry Smith 844bd481603SBarry Smith 845bd481603SBarry Smith Level: intermediate 846bd481603SBarry Smith 847bd481603SBarry Smith Notes: 848bd481603SBarry Smith See the chapter in the users manual on performance for more details 849bd481603SBarry Smith 850bd481603SBarry Smith Do not malloc or free dnz and onz that is handled internally by these routines 851bd481603SBarry Smith 852bd481603SBarry Smith Concepts: preallocation^Matrix 853bd481603SBarry Smith 854bd481603SBarry Smith .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSet(), MatPreallocateSetLocal(), 855bd481603SBarry Smith MatPreallocateSymmetricInitialize(), MatPreallocateSymmetricSetLocal() 856bd481603SBarry Smith M*/ 857ef66eb69SBarry Smith #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree(dnz);CHKERRQ(_4_ierr);} 8587c922b88SBarry Smith 859bd481603SBarry Smith 860bd481603SBarry Smith 8617b80b807SBarry Smith /* Routines unique to particular data structures */ 862be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellGetContext(Mat,void **); 863ba966639SSatish Balay PetscPolymorphicFunction(MatShellGetContext,(Mat A),(A,&t),void*,t) 864be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatBDiagGetData(Mat,PetscInt*,PetscInt*,PetscInt*[],PetscInt*[],PetscScalar***); 865698d4c6aSKris Buschelman 866be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInodeAdjustForInodes(Mat,IS*,IS*); 867be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *); 868698d4c6aSKris Buschelman 869be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetColumnIndices(Mat,PetscInt[]); 870be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetColumnIndices(Mat,PetscInt[]); 871be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 8727b80b807SBarry Smith 873a96a251dSBarry Smith #define MAT_SKIP_ALLOCATION -4 874a96a251dSBarry Smith 875be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]); 876ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz)) 877be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]); 878ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqSBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[]),(A,bs,0,nnz)) 879be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]); 880ba966639SSatish Balay PetscPolymorphicSubroutine(MatSeqAIJSetPreallocation,(Mat A,const PetscInt nnz[]),(A,0,nnz)) 881be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDensePreallocation(Mat,PetscScalar[]); 882be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBDiagSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscScalar*[]); 883273d9f13SBarry Smith 884be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 885ba966639SSatish Balay PetscPolymorphicSubroutine(MatMPIBAIJSetPreallocation,(Mat A,PetscInt bs,const PetscInt nnz[],const PetscInt onz[]),(A,bs,0,nnz,0,onz)) 886be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 887be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 888be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]); 889be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 890be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIDensePreallocation(Mat,PetscScalar[]); 891be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBDiagSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscScalar*[]); 892be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]); 893be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIDenseSetPreallocation(Mat,PetscScalar[]); 894284134d9SBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation(Mat,PetscScalar[]); 895be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIRowbsSetPreallocation(Mat,PetscInt,const PetscInt[]); 896be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,PetscInt*[]); 897be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,PetscInt*[]); 898be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatAdicSetLocalFunction(Mat,void (*)(void)); 899273d9f13SBarry Smith 900be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetLDA(Mat,PetscInt); 901ab92ecdeSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDenseGetLocalMatrix(Mat,Mat*); 9021b807ce4Svictorle 903be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatStoreValues(Mat); 904be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRetrieveValues(Mat); 9052e8a6d31SBarry Smith 906be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDAADSetCtx(Mat,void*); 9073a7fca6bSBarry Smith 9087b80b807SBarry Smith /* 9097b80b807SBarry Smith These routines are not usually accessed directly, rather solving is 91094b7f48cSBarry Smith done through the KSP and PC interfaces. 9117b80b807SBarry Smith */ 9127b80b807SBarry Smith 913d9274352SBarry Smith /*E 914d9274352SBarry Smith MatOrderingType - String with the name of a PETSc matrix ordering or the creation function 915d9274352SBarry Smith with an optional dynamic library name, for example 916d9274352SBarry Smith http://www.mcs.anl.gov/petsc/lib.a:orderingcreate() 917d9274352SBarry Smith 918d9274352SBarry Smith Level: beginner 919d9274352SBarry Smith 920d9274352SBarry Smith .seealso: MatGetOrdering() 921d9274352SBarry Smith E*/ 92249773a63SBarry Smith #define MatOrderingType char* 923b12f92e5SBarry Smith #define MATORDERING_NATURAL "natural" 924b12f92e5SBarry Smith #define MATORDERING_ND "nd" 925b12f92e5SBarry Smith #define MATORDERING_1WD "1wd" 926b12f92e5SBarry Smith #define MATORDERING_RCM "rcm" 927b12f92e5SBarry Smith #define MATORDERING_QMD "qmd" 928b12f92e5SBarry Smith #define MATORDERING_ROWLENGTH "rowlength" 92962152c8bSBarry Smith #define MATORDERING_DSC_ND "dsc_nd" 93062152c8bSBarry Smith #define MATORDERING_DSC_MMD "dsc_mmd" 93162152c8bSBarry Smith #define MATORDERING_DSC_MDF "dsc_mdf" 932c06d978dSMatthew Knepley #define MATORDERING_CONSTRAINED "constrained" 933c06d978dSMatthew Knepley #define MATORDERING_IDENTITY "identity" 934c06d978dSMatthew Knepley #define MATORDERING_REVERSE "reverse" 935b12f92e5SBarry Smith 936be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOrdering(Mat,const MatOrderingType,IS*,IS*); 937be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetOrderingList(PetscFList *list); 938be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,const MatOrderingType,IS*,IS*)); 93930de9b25SBarry Smith 94030de9b25SBarry Smith /*MC 94130de9b25SBarry Smith MatOrderingRegisterDynamic - Adds a new sparse matrix ordering to the 94230de9b25SBarry Smith matrix package. 94330de9b25SBarry Smith 94430de9b25SBarry Smith Synopsis: 945c1ac3661SBarry Smith PetscErrorCode MatOrderingRegisterDynamic(char *name_ordering,char *path,char *name_create,PetscErrorCode (*routine_create)(MatOrdering)) 94630de9b25SBarry Smith 94730de9b25SBarry Smith Not Collective 94830de9b25SBarry Smith 94930de9b25SBarry Smith Input Parameters: 95030de9b25SBarry Smith + sname - name of ordering (for example MATORDERING_ND) 95130de9b25SBarry Smith . path - location of library where creation routine is 95230de9b25SBarry Smith . name - name of function that creates the ordering type,a string 95330de9b25SBarry Smith - function - function pointer that creates the ordering 95430de9b25SBarry Smith 95530de9b25SBarry Smith Level: developer 95630de9b25SBarry Smith 95730de9b25SBarry Smith If dynamic libraries are used, then the fourth input argument (function) 95830de9b25SBarry Smith is ignored. 95930de9b25SBarry Smith 96030de9b25SBarry Smith Sample usage: 96130de9b25SBarry Smith .vb 96230de9b25SBarry Smith MatOrderingRegisterDynamic("my_order",/home/username/my_lib/lib/libO/solaris/mylib.a, 96330de9b25SBarry Smith "MyOrder",MyOrder); 96430de9b25SBarry Smith .ve 96530de9b25SBarry Smith 96630de9b25SBarry Smith Then, your partitioner can be chosen with the procedural interface via 96730de9b25SBarry Smith $ MatOrderingSetType(part,"my_order) 96830de9b25SBarry Smith or at runtime via the option 96930de9b25SBarry Smith $ -pc_ilu_mat_ordering_type my_order 97030de9b25SBarry Smith $ -pc_lu_mat_ordering_type my_order 97130de9b25SBarry Smith 972ab901514SBarry Smith ${PETSC_ARCH} occuring in pathname will be replaced with appropriate values. 97330de9b25SBarry Smith 97430de9b25SBarry Smith .keywords: matrix, ordering, register 97530de9b25SBarry Smith 97630de9b25SBarry Smith .seealso: MatOrderingRegisterDestroy(), MatOrderingRegisterAll() 97730de9b25SBarry Smith M*/ 978aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 979f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,0) 980b12f92e5SBarry Smith #else 981f1af5d2fSBarry Smith #define MatOrderingRegisterDynamic(a,b,c,d) MatOrderingRegister(a,b,c,d) 982b12f92e5SBarry Smith #endif 98330de9b25SBarry Smith 984be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegisterDestroy(void); 985be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatOrderingRegisterAll(const char[]); 9862bad1931SBarry Smith extern PetscTruth MatOrderingRegisterAllCalled; 987b0a32e0cSBarry Smith extern PetscFList MatOrderingList; 988d4fbbf0eSBarry Smith 989be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS); 990a2ce50c7SBarry Smith 991d91e6319SBarry Smith /*S 99215e8a5b3SHong Zhang MatFactorInfo - Data based into the matrix factorization routines 993b00f7748SHong Zhang 99415e8a5b3SHong Zhang In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE 995b00f7748SHong Zhang 99615e8a5b3SHong Zhang Notes: These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC. 997b00f7748SHong Zhang 998b00f7748SHong Zhang Level: developer 999b00f7748SHong Zhang 1000d7d82daaSBarry Smith .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(), 1001d7d82daaSBarry Smith MatFactorInfoInitialize() 1002b00f7748SHong Zhang 1003b00f7748SHong Zhang S*/ 1004b00f7748SHong Zhang typedef struct { 10050a29876aSHong Zhang PetscReal shiftnz; /* scaling of identity added to matrix to prevent zero pivots */ 1006fbf22428SSatish Balay PetscReal shiftpd; /* if true, shift until positive pivots */ 10072cea7109SBarry Smith PetscReal shift_fraction; /* record shift fraction taken */ 100815e8a5b3SHong Zhang PetscReal diagonal_fill; /* force diagonal to fill in if initially not filled */ 100915e8a5b3SHong Zhang PetscReal dt; /* drop tolerance */ 1010b00f7748SHong Zhang PetscReal dtcol; /* tolerance for pivoting */ 101115e8a5b3SHong Zhang PetscReal dtcount; /* maximum nonzeros to be allowed per row */ 1012f6275e2eSBarry Smith PetscReal fill; /* expected fill; nonzeros in factored matrix/nonzeros in original matrix*/ 1013348344bbSBarry Smith PetscReal levels; /* ICC/ILU(levels) */ 1014bcd9e38bSBarry Smith PetscReal pivotinblocks; /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0 1015bcd9e38bSBarry Smith factorization may be faster if do not pivot */ 101615e8a5b3SHong Zhang PetscReal zeropivot; /* pivot is called zero if less than this */ 101715e8a5b3SHong Zhang } MatFactorInfo; 1018ffa6d0a5SLois Curfman McInnes 1019be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFactorInfoInitialize(MatFactorInfo*); 1020be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactor(Mat,IS,MatFactorInfo*); 1021be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactorSymbolic(Mat,IS,MatFactorInfo*,Mat*); 1022be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCholeskyFactorNumeric(Mat,MatFactorInfo*,Mat*); 1023be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactor(Mat,IS,IS,MatFactorInfo*); 1024be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUFactor(Mat,IS,IS,MatFactorInfo*); 1025be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactorSymbolic(Mat,IS,IS,MatFactorInfo*,Mat*); 1026be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUFactorSymbolic(Mat,IS,IS,MatFactorInfo*,Mat*); 1027be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatICCFactorSymbolic(Mat,IS,MatFactorInfo*,Mat*); 1028be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatICCFactor(Mat,IS,MatFactorInfo*); 1029be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatLUFactorNumeric(Mat,MatFactorInfo*,Mat*); 1030be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatILUDTFactor(Mat,IS,IS,MatFactorInfo*,Mat *); 1031be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*); 1032be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolve(Mat,Vec,Vec); 1033be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatForwardSolve(Mat,Vec,Vec); 1034be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatBackwardSolve(Mat,Vec,Vec); 1035be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveAdd(Mat,Vec,Vec,Vec); 1036be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveTranspose(Mat,Vec,Vec); 1037be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolveTransposeAdd(Mat,Vec,Vec,Vec); 1038be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSolves(Mat,Vecs,Vecs); 10398ed539a5SBarry Smith 1040be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSetUnfactored(Mat); 1041bb5a7306SBarry Smith 1042d91e6319SBarry Smith /*E 1043d91e6319SBarry Smith MatSORType - What type of (S)SOR to perform 1044bb1eb677SSatish Balay 1045d91e6319SBarry Smith Level: beginner 1046d91e6319SBarry Smith 1047d9274352SBarry Smith May be bitwise ORd together 1048d9274352SBarry Smith 1049d91e6319SBarry Smith Any additions/changes here MUST also be made in include/finclude/petscmat.h 1050d91e6319SBarry Smith 10514e7234bfSBarry Smith MatSORType may be bitwise ORd together, so do not change the numbers 10524e7234bfSBarry Smith 1053a1d92eedSBarry Smith .seealso: MatRelax(), MatPBRelax() 1054d91e6319SBarry Smith E*/ 1055ee50ffe9SBarry Smith typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3, 1056ee50ffe9SBarry Smith SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8, 1057ee50ffe9SBarry Smith SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16, 105884cb2905SBarry Smith SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType; 1059be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatRelax(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec); 1060be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPBRelax(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec); 10618ed539a5SBarry Smith 1062d4fbbf0eSBarry Smith /* 1063639f9d9dSBarry Smith These routines are for efficiently computing Jacobians via finite differences. 1064639f9d9dSBarry Smith */ 1065b12f92e5SBarry Smith 1066d9274352SBarry Smith /*E 1067d9274352SBarry Smith MatColoringType - String with the name of a PETSc matrix coloring or the creation function 1068d9274352SBarry Smith with an optional dynamic library name, for example 1069d9274352SBarry Smith http://www.mcs.anl.gov/petsc/lib.a:coloringcreate() 1070d9274352SBarry Smith 1071d9274352SBarry Smith Level: beginner 1072d9274352SBarry Smith 1073d9274352SBarry Smith .seealso: MatGetColoring() 1074d9274352SBarry Smith E*/ 107549773a63SBarry Smith #define MatColoringType char* 1076b12f92e5SBarry Smith #define MATCOLORING_NATURAL "natural" 1077b12f92e5SBarry Smith #define MATCOLORING_SL "sl" 1078b12f92e5SBarry Smith #define MATCOLORING_LF "lf" 1079b12f92e5SBarry Smith #define MATCOLORING_ID "id" 1080b12f92e5SBarry Smith 1081be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetColoring(Mat,const MatColoringType,ISColoring*); 1082be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegister(const char[],const char[],const char[],PetscErrorCode(*)(Mat,const MatColoringType,ISColoring *)); 108330de9b25SBarry Smith 108430de9b25SBarry Smith /*MC 108530de9b25SBarry Smith MatColoringRegisterDynamic - Adds a new sparse matrix coloring to the 108630de9b25SBarry Smith matrix package. 108730de9b25SBarry Smith 108830de9b25SBarry Smith Synopsis: 1089c1ac3661SBarry Smith PetscErrorCode MatColoringRegisterDynamic(char *name_coloring,char *path,char *name_create,PetscErrorCode (*routine_create)(MatColoring)) 109030de9b25SBarry Smith 109130de9b25SBarry Smith Not Collective 109230de9b25SBarry Smith 109330de9b25SBarry Smith Input Parameters: 109430de9b25SBarry Smith + sname - name of Coloring (for example MATCOLORING_SL) 109530de9b25SBarry Smith . path - location of library where creation routine is 109630de9b25SBarry Smith . name - name of function that creates the Coloring type, a string 109730de9b25SBarry Smith - function - function pointer that creates the coloring 109830de9b25SBarry Smith 109930de9b25SBarry Smith Level: developer 110030de9b25SBarry Smith 110130de9b25SBarry Smith If dynamic libraries are used, then the fourth input argument (function) 110230de9b25SBarry Smith is ignored. 110330de9b25SBarry Smith 110430de9b25SBarry Smith Sample usage: 110530de9b25SBarry Smith .vb 110630de9b25SBarry Smith MatColoringRegisterDynamic("my_color",/home/username/my_lib/lib/libO/solaris/mylib.a, 110730de9b25SBarry Smith "MyColor",MyColor); 110830de9b25SBarry Smith .ve 110930de9b25SBarry Smith 111030de9b25SBarry Smith Then, your partitioner can be chosen with the procedural interface via 111130de9b25SBarry Smith $ MatColoringSetType(part,"my_color") 111230de9b25SBarry Smith or at runtime via the option 111330de9b25SBarry Smith $ -mat_coloring_type my_color 111430de9b25SBarry Smith 1115ab901514SBarry Smith $PETSC_ARCH occuring in pathname will be replaced with appropriate values. 111630de9b25SBarry Smith 111730de9b25SBarry Smith .keywords: matrix, Coloring, register 111830de9b25SBarry Smith 111930de9b25SBarry Smith .seealso: MatColoringRegisterDestroy(), MatColoringRegisterAll() 112030de9b25SBarry Smith M*/ 1121aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 1122f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,0) 1123b12f92e5SBarry Smith #else 1124f1af5d2fSBarry Smith #define MatColoringRegisterDynamic(a,b,c,d) MatColoringRegister(a,b,c,d) 1125b12f92e5SBarry Smith #endif 112630de9b25SBarry Smith 11272bad1931SBarry Smith extern PetscTruth MatColoringRegisterAllCalled; 1128f1144a30SSatish Balay 1129be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegisterAll(const char[]); 1130be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringRegisterDestroy(void); 1131be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*); 1132639f9d9dSBarry Smith 1133d9274352SBarry Smith /*S 1134d9274352SBarry Smith MatFDColoring - Object for computing a sparse Jacobian via finite differences 1135d9274352SBarry Smith and coloring 1136639f9d9dSBarry Smith 1137d9274352SBarry Smith Level: beginner 1138d9274352SBarry Smith 1139d9274352SBarry Smith Concepts: coloring, sparse Jacobian, finite differences 1140d9274352SBarry Smith 1141d9274352SBarry Smith .seealso: MatFDColoringCreate() 1142d9274352SBarry Smith S*/ 1143e2a1c21fSSatish Balay typedef struct _p_MatFDColoring* MatFDColoring; 1144639f9d9dSBarry Smith 1145be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringCreate(Mat,ISColoring,MatFDColoring *); 1146be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringDestroy(MatFDColoring); 1147be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringView(MatFDColoring,PetscViewer); 1148be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*); 11494a9d489dSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**); 1150be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal); 1151be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFrequency(MatFDColoring,PetscInt); 1152be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetFrequency(MatFDColoring,PetscInt*); 1153be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetFromOptions(MatFDColoring); 1154be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApply(Mat,MatFDColoring,Vec,MatStructure*,void *); 1155be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringApplyTS(Mat,MatFDColoring,PetscReal,Vec,MatStructure*,void *); 1156be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetRecompute(MatFDColoring); 1157be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringSetF(MatFDColoring,Vec); 1158be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,PetscInt*[]); 1159639f9d9dSBarry Smith /* 11600752156aSBarry Smith These routines are for partitioning matrices: currently used only 11613eda8832SBarry Smith for adjacency matrix, MatCreateMPIAdj(). 11620752156aSBarry Smith */ 1163ca161407SBarry Smith 1164d9274352SBarry Smith /*S 1165d9274352SBarry Smith MatPartitioning - Object for managing the partitioning of a matrix or graph 1166d9274352SBarry Smith 1167d9274352SBarry Smith Level: beginner 1168d9274352SBarry Smith 1169d9274352SBarry Smith Concepts: partitioning 1170d9274352SBarry Smith 1171743a1026Svictorle .seealso: MatPartitioningCreate(), MatPartitioningType 1172d9274352SBarry Smith S*/ 117391e9ee9fSBarry Smith typedef struct _p_MatPartitioning* MatPartitioning; 1174d9274352SBarry Smith 1175d9274352SBarry Smith /*E 11765bc6e7ccSvictorle MatPartitioningType - String with the name of a PETSc matrix partitioning or the creation function 1177d9274352SBarry Smith with an optional dynamic library name, for example 1178d9274352SBarry Smith http://www.mcs.anl.gov/petsc/lib.a:partitioningcreate() 1179d9274352SBarry Smith 1180d9274352SBarry Smith Level: beginner 1181d9274352SBarry Smith 1182b547a13bSvictorle .seealso: MatPartitioningCreate(), MatPartitioning 1183d9274352SBarry Smith E*/ 118449773a63SBarry Smith #define MatPartitioningType char* 11858ba1e511SMatthew Knepley #define MAT_PARTITIONING_CURRENT "current" 1186c2ec2a73SHong Zhang #define MAT_PARTITIONING_SQUARE "square" 11878ba1e511SMatthew Knepley #define MAT_PARTITIONING_PARMETIS "parmetis" 118871306c60Sjroman #define MAT_PARTITIONING_CHACO "chaco" 118971306c60Sjroman #define MAT_PARTITIONING_JOSTLE "jostle" 119071306c60Sjroman #define MAT_PARTITIONING_PARTY "party" 119171306c60Sjroman #define MAT_PARTITIONING_SCOTCH "scotch" 119271306c60Sjroman 1193ca161407SBarry Smith 1194be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningCreate(MPI_Comm,MatPartitioning*); 1195be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetType(MatPartitioning,const MatPartitioningType); 1196be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetNParts(MatPartitioning,PetscInt); 1197be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetAdjacency(MatPartitioning,Mat); 1198be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]); 1199be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []); 1200be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningApply(MatPartitioning,IS*); 1201be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningDestroy(MatPartitioning); 12022aabb6bbSBarry Smith 1203be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegister(const char[],const char[],const char[],PetscErrorCode (*)(MatPartitioning)); 120430de9b25SBarry Smith 120530de9b25SBarry Smith /*MC 120630de9b25SBarry Smith MatPartitioningRegisterDynamic - Adds a new sparse matrix partitioning to the 120730de9b25SBarry Smith matrix package. 120830de9b25SBarry Smith 120930de9b25SBarry Smith Synopsis: 1210c1ac3661SBarry Smith PetscErrorCode MatPartitioningRegisterDynamic(char *name_partitioning,char *path,char *name_create,PetscErrorCode (*routine_create)(MatPartitioning)) 121130de9b25SBarry Smith 121230de9b25SBarry Smith Not Collective 121330de9b25SBarry Smith 121430de9b25SBarry Smith Input Parameters: 121530de9b25SBarry Smith + sname - name of partitioning (for example MAT_PARTITIONING_CURRENT) or parmetis 121630de9b25SBarry Smith . path - location of library where creation routine is 121730de9b25SBarry Smith . name - name of function that creates the partitioning type, a string 121830de9b25SBarry Smith - function - function pointer that creates the partitioning type 121930de9b25SBarry Smith 122030de9b25SBarry Smith Level: developer 122130de9b25SBarry Smith 122230de9b25SBarry Smith If dynamic libraries are used, then the fourth input argument (function) 122330de9b25SBarry Smith is ignored. 122430de9b25SBarry Smith 122530de9b25SBarry Smith Sample usage: 122630de9b25SBarry Smith .vb 122730de9b25SBarry Smith MatPartitioningRegisterDynamic("my_part",/home/username/my_lib/lib/libO/solaris/mylib.a, 122830de9b25SBarry Smith "MyPartCreate",MyPartCreate); 122930de9b25SBarry Smith .ve 123030de9b25SBarry Smith 123130de9b25SBarry Smith Then, your partitioner can be chosen with the procedural interface via 123230de9b25SBarry Smith $ MatPartitioningSetType(part,"my_part") 123330de9b25SBarry Smith or at runtime via the option 123430de9b25SBarry Smith $ -mat_partitioning_type my_part 123530de9b25SBarry Smith 1236ab901514SBarry Smith $PETSC_ARCH occuring in pathname will be replaced with appropriate values. 123730de9b25SBarry Smith 123830de9b25SBarry Smith .keywords: matrix, partitioning, register 123930de9b25SBarry Smith 124030de9b25SBarry Smith .seealso: MatPartitioningRegisterDestroy(), MatPartitioningRegisterAll() 124130de9b25SBarry Smith M*/ 1242aa482453SBarry Smith #if defined(PETSC_USE_DYNAMIC_LIBRARIES) 1243f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,0) 12442aabb6bbSBarry Smith #else 1245f1af5d2fSBarry Smith #define MatPartitioningRegisterDynamic(a,b,c,d) MatPartitioningRegister(a,b,c,d) 12462aabb6bbSBarry Smith #endif 12472aabb6bbSBarry Smith 12482bad1931SBarry Smith extern PetscTruth MatPartitioningRegisterAllCalled; 1249f1144a30SSatish Balay 1250be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegisterAll(const char[]); 1251be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningRegisterDestroy(void); 12522bad1931SBarry Smith 1253be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningView(MatPartitioning,PetscViewer); 1254be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningSetFromOptions(MatPartitioning); 1255be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningGetType(MatPartitioning,MatPartitioningType*); 1256ca161407SBarry Smith 1257be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningParmetisSetCoarseSequential(MatPartitioning); 12580752156aSBarry Smith 1259be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningJostleSetCoarseLevel(MatPartitioning,PetscReal); 1260be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningJostleSetCoarseSequential(MatPartitioning); 126171306c60Sjroman 1262954b3ec6SBarry Smith typedef enum {MP_CHACO_MULTILEVEL_KL,MP_CHACO_SPECTRAL,MP_CHACO_LINEAR,MP_CHACO_RANDOM, MP_CHACO_SCATTERED} MPChacoGlobalType; 1263be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetGlobal(MatPartitioning, MPChacoGlobalType); 126471306c60Sjroman typedef enum { MP_CHACO_KERNIGHAN_LIN, MP_CHACO_NONE } MPChacoLocalType; 1265be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetLocal(MatPartitioning, MPChacoLocalType); 1266be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal); 126771306c60Sjroman typedef enum { MP_CHACO_LANCZOS, MP_CHACO_RQI_SYMMLQ } MPChacoEigenType; 1268be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType); 1269be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenTol(MatPartitioning, PetscReal); 1270be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningChacoSetEigenNumber(MatPartitioning, PetscInt); 127171306c60Sjroman 127271306c60Sjroman #define MP_PARTY_OPT "opt" 127371306c60Sjroman #define MP_PARTY_LIN "lin" 127471306c60Sjroman #define MP_PARTY_SCA "sca" 127571306c60Sjroman #define MP_PARTY_RAN "ran" 127671306c60Sjroman #define MP_PARTY_GBF "gbf" 127771306c60Sjroman #define MP_PARTY_GCF "gcf" 127871306c60Sjroman #define MP_PARTY_BUB "bub" 127971306c60Sjroman #define MP_PARTY_DEF "def" 1280be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetGlobal(MatPartitioning, const char*); 128171306c60Sjroman #define MP_PARTY_HELPFUL_SETS "hs" 128271306c60Sjroman #define MP_PARTY_KERNIGHAN_LIN "kl" 128371306c60Sjroman #define MP_PARTY_NONE "no" 1284be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetLocal(MatPartitioning, const char*); 1285be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal); 1286be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetBipart(MatPartitioning,PetscTruth); 1287be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscTruth); 128871306c60Sjroman 128971306c60Sjroman typedef enum { MP_SCOTCH_GREEDY, MP_SCOTCH_GPS, MP_SCOTCH_GR_GPS } MPScotchGlobalType; 1290be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetArch(MatPartitioning,const char*); 1291be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetMultilevel(MatPartitioning); 1292be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetGlobal(MatPartitioning,MPScotchGlobalType); 1293be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetCoarseLevel(MatPartitioning,PetscReal); 1294be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetHostList(MatPartitioning,const char*); 129571306c60Sjroman typedef enum { MP_SCOTCH_KERNIGHAN_LIN, MP_SCOTCH_NONE } MPScotchLocalType; 1296be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetLocal(MatPartitioning,MPScotchLocalType); 1297be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetMapping(MatPartitioning); 1298be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatPartitioningScotchSetStrategy(MatPartitioning,char*); 129971306c60Sjroman 13000752156aSBarry Smith /* 13010a835dfdSSatish Balay If you add entries here you must also add them to finclude/petscmat.h 1302d4fbbf0eSBarry Smith */ 13031c1c02c0SLois Curfman McInnes typedef enum { MATOP_SET_VALUES=0, 13041c1c02c0SLois Curfman McInnes MATOP_GET_ROW=1, 13051c1c02c0SLois Curfman McInnes MATOP_RESTORE_ROW=2, 13061c1c02c0SLois Curfman McInnes MATOP_MULT=3, 13071c1c02c0SLois Curfman McInnes MATOP_MULT_ADD=4, 13087c922b88SBarry Smith MATOP_MULT_TRANSPOSE=5, 13097c922b88SBarry Smith MATOP_MULT_TRANSPOSE_ADD=6, 13101c1c02c0SLois Curfman McInnes MATOP_SOLVE=7, 13111c1c02c0SLois Curfman McInnes MATOP_SOLVE_ADD=8, 13127c922b88SBarry Smith MATOP_SOLVE_TRANSPOSE=9, 13137c922b88SBarry Smith MATOP_SOLVE_TRANSPOSE_ADD=10, 13141c1c02c0SLois Curfman McInnes MATOP_LUFACTOR=11, 13151c1c02c0SLois Curfman McInnes MATOP_CHOLESKYFACTOR=12, 13161c1c02c0SLois Curfman McInnes MATOP_RELAX=13, 13171c1c02c0SLois Curfman McInnes MATOP_TRANSPOSE=14, 13181c1c02c0SLois Curfman McInnes MATOP_GETINFO=15, 13191c1c02c0SLois Curfman McInnes MATOP_EQUAL=16, 13201c1c02c0SLois Curfman McInnes MATOP_GET_DIAGONAL=17, 13211c1c02c0SLois Curfman McInnes MATOP_DIAGONAL_SCALE=18, 13221c1c02c0SLois Curfman McInnes MATOP_NORM=19, 13231c1c02c0SLois Curfman McInnes MATOP_ASSEMBLY_BEGIN=20, 13241c1c02c0SLois Curfman McInnes MATOP_ASSEMBLY_END=21, 13251c1c02c0SLois Curfman McInnes MATOP_COMPRESS=22, 13261c1c02c0SLois Curfman McInnes MATOP_SET_OPTION=23, 13271c1c02c0SLois Curfman McInnes MATOP_ZERO_ENTRIES=24, 13281c1c02c0SLois Curfman McInnes MATOP_ZERO_ROWS=25, 13291c1c02c0SLois Curfman McInnes MATOP_LUFACTOR_SYMBOLIC=26, 13301c1c02c0SLois Curfman McInnes MATOP_LUFACTOR_NUMERIC=27, 13311c1c02c0SLois Curfman McInnes MATOP_CHOLESKY_FACTOR_SYMBOLIC=28, 13321c1c02c0SLois Curfman McInnes MATOP_CHOLESKY_FACTOR_NUMERIC=29, 1333d643ce63SMatthew Knepley MATOP_SETUP_PREALLOCATION=30, 1334d643ce63SMatthew Knepley MATOP_ILUFACTOR_SYMBOLIC=31, 1335d643ce63SMatthew Knepley MATOP_ICCFACTOR_SYMBOLIC=32, 1336d643ce63SMatthew Knepley MATOP_GET_ARRAY=33, 1337d643ce63SMatthew Knepley MATOP_RESTORE_ARRAY=34, 1338e26f1c3cSBarry Smith MATOP_DUPLICATE=35, 1339d643ce63SMatthew Knepley MATOP_FORWARD_SOLVE=36, 1340d643ce63SMatthew Knepley MATOP_BACKWARD_SOLVE=37, 1341d643ce63SMatthew Knepley MATOP_ILUFACTOR=38, 1342d643ce63SMatthew Knepley MATOP_ICCFACTOR=39, 1343d643ce63SMatthew Knepley MATOP_AXPY=40, 1344d643ce63SMatthew Knepley MATOP_GET_SUBMATRICES=41, 1345d643ce63SMatthew Knepley MATOP_INCREASE_OVERLAP=42, 1346d643ce63SMatthew Knepley MATOP_GET_VALUES=43, 1347d643ce63SMatthew Knepley MATOP_COPY=44, 1348d643ce63SMatthew Knepley MATOP_PRINT_HELP=45, 1349d643ce63SMatthew Knepley MATOP_SCALE=46, 1350d643ce63SMatthew Knepley MATOP_SHIFT=47, 1351d643ce63SMatthew Knepley MATOP_DIAGONAL_SHIFT=48, 1352d643ce63SMatthew Knepley MATOP_ILUDT_FACTOR=49, 1353d643ce63SMatthew Knepley MATOP_GET_BLOCK_SIZE=50, 1354d643ce63SMatthew Knepley MATOP_GET_ROW_IJ=51, 1355d643ce63SMatthew Knepley MATOP_RESTORE_ROW_IJ=52, 1356d643ce63SMatthew Knepley MATOP_GET_COLUMN_IJ=53, 1357d643ce63SMatthew Knepley MATOP_RESTORE_COLUMN_IJ=54, 1358d643ce63SMatthew Knepley MATOP_FDCOLORING_CREATE=55, 1359d643ce63SMatthew Knepley MATOP_COLORING_PATCH=56, 1360d643ce63SMatthew Knepley MATOP_SET_UNFACTORED=57, 1361d643ce63SMatthew Knepley MATOP_PERMUTE=58, 1362d643ce63SMatthew Knepley MATOP_SET_VALUES_BLOCKED=59, 1363d643ce63SMatthew Knepley MATOP_GET_SUBMATRIX=60, 1364d643ce63SMatthew Knepley MATOP_DESTROY=61, 1365d643ce63SMatthew Knepley MATOP_VIEW=62, 1366d643ce63SMatthew Knepley MATOP_GET_MAPS=63, 1367d643ce63SMatthew Knepley MATOP_USE_SCALED_FORM=64, 1368d643ce63SMatthew Knepley MATOP_SCALE_SYSTEM=65, 1369d643ce63SMatthew Knepley MATOP_UNSCALE_SYSTEM=66, 1370ba16a8cbSKris Buschelman MATOP_SET_LOCAL_TO_GLOBAL_MAP=67, 1371d643ce63SMatthew Knepley MATOP_SET_VALUES_LOCAL=68, 1372d643ce63SMatthew Knepley MATOP_ZERO_ROWS_LOCAL=69, 1373d643ce63SMatthew Knepley MATOP_GET_ROW_MAX=70, 1374d643ce63SMatthew Knepley MATOP_CONVERT=71, 1375d643ce63SMatthew Knepley MATOP_SET_COLORING=72, 1376d643ce63SMatthew Knepley MATOP_SET_VALUES_ADIC=73, 1377d643ce63SMatthew Knepley MATOP_SET_VALUES_ADIFOR=74, 1378d643ce63SMatthew Knepley MATOP_FD_COLORING_APPLY=75, 1379d643ce63SMatthew Knepley MATOP_SET_FROM_OPTIONS=76, 1380ba16a8cbSKris Buschelman MATOP_MULT_CON=77, 1381ba16a8cbSKris Buschelman MATOP_MULT_TRANSPOSE_CON=78, 1382ba16a8cbSKris Buschelman MATOP_ILU_FACTOR_SYMBOLIC_CON=79, 1383d643ce63SMatthew Knepley MATOP_PERMUTE_SPARSIFY=80, 1384d643ce63SMatthew Knepley MATOP_MULT_MULTIPLE=81, 138541acf15aSKris Buschelman MATOP_SOLVE_MULTIPLE=82, 138641acf15aSKris Buschelman MATOP_GET_INERTIA=83, 138741acf15aSKris Buschelman MATOP_LOAD=84, 138841acf15aSKris Buschelman MATOP_IS_SYMMETRIC=85, 138941acf15aSKris Buschelman MATOP_IS_HERMITIAN=86, 139041acf15aSKris Buschelman MATOP_IS_STRUCTURALLY_SYMMETRIC=87, 139141acf15aSKris Buschelman MATOP_PB_RELAX=88, 139241acf15aSKris Buschelman MATOP_GET_VECS=89, 139341acf15aSKris Buschelman MATOP_MAT_MULT=90, 139441acf15aSKris Buschelman MATOP_MAT_MULT_SYMBOLIC=91, 139541acf15aSKris Buschelman MATOP_MAT_MULT_NUMERIC=92, 139641acf15aSKris Buschelman MATOP_PTAP=93, 139741acf15aSKris Buschelman MATOP_PTAP_SYMBOLIC=94, 1398bc011b1eSHong Zhang MATOP_PTAP_NUMERIC=95, 1399bc011b1eSHong Zhang MATOP_MAT_MULTTRANSPOSE=96, 1400bc011b1eSHong Zhang MATOP_MAT_MULTTRANSPOSE_SYMBOLIC=97, 14017ba1a0bfSKris Buschelman MATOP_MAT_MULTTRANSPOSE_NUMERIC=98, 14027ba1a0bfSKris Buschelman MATOP_PTAP_SYMBOLIC_SEQAIJ=99, 14037ba1a0bfSKris Buschelman MATOP_PTAP_NUMERIC_SEQAIJ=100, 14047ba1a0bfSKris Buschelman MATOP_PTAP_SYMBOLIC_MPIAIJ=101, 140587d4246cSBarry Smith MATOP_PTAP_NUMERIC_MPIAIJ=102, 1406f5edf698SHong Zhang MATOP_SET_VALUES_ROW = 105, 1407d5c63f83SSatish Balay MATOP_GET_ROW_UTRIANGULAR=108, 1408*2bebee5dSHong Zhang MATOP_RESTORE_ROW_UTRIANGULAR=109, 1409*2bebee5dSHong Zhang MATOP_MATSOLVE=110 1410fae171e0SBarry Smith } MatOperation; 1411be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatHasOperation(Mat,MatOperation,PetscTruth*); 1412be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellSetOperation(Mat,MatOperation,void(*)(void)); 1413be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellGetOperation(Mat,MatOperation,void(**)(void)); 1414be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatShellSetContext(Mat,void*); 1415112a2221SBarry Smith 141690ace30eSBarry Smith /* 141790ace30eSBarry Smith Codes for matrices stored on disk. By default they are 141890ace30eSBarry Smith stored in a universal format. By changing the format with 1419fb9695e5SSatish Balay PetscViewerSetFormat(viewer,PETSC_VIEWER_BINARY_NATIVE); the matrices will 142090ace30eSBarry Smith be stored in a way natural for the matrix, for example dense matrices 142190ace30eSBarry Smith would be stored as dense. Matrices stored this way may only be 142290ace30eSBarry Smith read into matrices of the same time. 142390ace30eSBarry Smith */ 142490ace30eSBarry Smith #define MATRIX_BINARY_FORMAT_DENSE -1 142590ace30eSBarry Smith 1426be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJSetHashTableFactor(Mat,PetscReal); 1427be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPIRowbsGetColor(Mat,ISColoring *); 1428860d1616SSatish Balay 1429be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatISGetLocalMat(Mat,Mat*); 14301f4e1ec7SBarry Smith 1431d9274352SBarry Smith /*S 1432d9274352SBarry Smith MatNullSpace - Object that removes a null space from a vector, i.e. 1433d9274352SBarry Smith orthogonalizes the vector to a subsapce 1434d9274352SBarry Smith 1435f7a9e4ceSBarry Smith Level: advanced 1436d9274352SBarry Smith 1437d9274352SBarry Smith Concepts: matrix; linear operator, null space 1438d9274352SBarry Smith 14396e1639daSBarry Smith Users manual sections: 14406e1639daSBarry Smith . sec_singular 14416e1639daSBarry Smith 1442d9274352SBarry Smith .seealso: MatNullSpaceCreate() 1443d9274352SBarry Smith S*/ 144474637425SBarry Smith typedef struct _p_MatNullSpace* MatNullSpace; 1445d9274352SBarry Smith 1446be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceCreate(MPI_Comm,PetscTruth,PetscInt,const Vec[],MatNullSpace*); 1447281d1b2eSBarry Smith EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(Vec,void*),void*); 1448be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceDestroy(MatNullSpace); 1449be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceRemove(MatNullSpace,Vec,Vec*); 1450be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceAttach(Mat,MatNullSpace); 1451be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatNullSpaceTest(MatNullSpace,Mat); 145274637425SBarry Smith 1453be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatReorderingSeqSBAIJ(Mat,IS); 1454be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMPISBAIJSetHashTableFactor(Mat,PetscReal); 1455be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqSBAIJSetColumnIndices(Mat,PetscInt *); 145646129b97SKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatSeqBAIJInvertBlockDiagonal(Mat); 14573f1d51d7SBarry Smith 1458be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMAIJ(Mat,PetscInt,Mat*); 1459be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMAIJRedimension(Mat,PetscInt,Mat*); 1460be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatMAIJGetAIJ(Mat,Mat*); 1461c4f061fbSSatish Balay 1462be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatComputeExplicitOperator(Mat,Mat*); 1463b0a32e0cSBarry Smith 1464be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalScaleLocal(Mat,Vec); 146504f1ad80SBarry Smith 1466be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *); 1467be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *); 14687dbadf16SMatthew Knepley 1469e9fa29b7SSatish Balay PETSC_EXTERN_CXX_END 14702eac72dbSBarry Smith #endif 1471