1 /* 2 Include file for the matrix component of PETSc 3 */ 4 #ifndef PETSCMAT_H 5 #define PETSCMAT_H 6 #include <petscvec.h> 7 8 /*S 9 Mat - Abstract PETSc matrix object used to manage all linear operators in PETSc, even those without 10 an explicit sparse representation (such as matrix-free operators) 11 12 Level: beginner 13 14 .seealso: MatCreate(), MatType, MatSetType(), MatDestroy() 15 S*/ 16 typedef struct _p_Mat* Mat; 17 18 /*J 19 MatType - String with the name of a PETSc matrix type 20 21 Level: beginner 22 23 .seealso: MatSetType(), Mat, MatSolverType, MatRegister() 24 J*/ 25 typedef const char* MatType; 26 #define MATSAME "same" 27 #define MATMAIJ "maij" 28 #define MATSEQMAIJ "seqmaij" 29 #define MATMPIMAIJ "mpimaij" 30 #define MATKAIJ "kaij" 31 #define MATSEQKAIJ "seqkaij" 32 #define MATMPIKAIJ "mpikaij" 33 #define MATIS "is" 34 #define MATAIJ "aij" 35 #define MATSEQAIJ "seqaij" 36 #define MATMPIAIJ "mpiaij" 37 #define MATAIJCRL "aijcrl" 38 #define MATSEQAIJCRL "seqaijcrl" 39 #define MATMPIAIJCRL "mpiaijcrl" 40 #define MATAIJCUSPARSE "aijcusparse" 41 #define MATSEQAIJCUSPARSE "seqaijcusparse" 42 #define MATMPIAIJCUSPARSE "mpiaijcusparse" 43 #define MATAIJVIENNACL "aijviennacl" 44 #define MATSEQAIJVIENNACL "seqaijviennacl" 45 #define MATMPIAIJVIENNACL "mpiaijviennacl" 46 #define MATAIJPERM "aijperm" 47 #define MATSEQAIJPERM "seqaijperm" 48 #define MATMPIAIJPERM "mpiaijperm" 49 #define MATAIJSELL "aijsell" 50 #define MATSEQAIJSELL "seqaijsell" 51 #define MATMPIAIJSELL "mpiaijsell" 52 #define MATAIJMKL "aijmkl" 53 #define MATSEQAIJMKL "seqaijmkl" 54 #define MATMPIAIJMKL "mpiaijmkl" 55 #define MATBAIJMKL "baijmkl" 56 #define MATSEQBAIJMKL "seqbaijmkl" 57 #define MATMPIBAIJMKL "mpibaijmkl" 58 #define MATSHELL "shell" 59 #define MATDENSE "dense" 60 #define MATDENSECUDA "densecuda" 61 #define MATSEQDENSE "seqdense" 62 #define MATSEQDENSECUDA "seqdensecuda" 63 #define MATMPIDENSE "mpidense" 64 #define MATMPIDENSECUDA "mpidensecuda" 65 #define MATELEMENTAL "elemental" 66 #define MATSCALAPACK "scalapack" 67 #define MATBAIJ "baij" 68 #define MATSEQBAIJ "seqbaij" 69 #define MATMPIBAIJ "mpibaij" 70 #define MATMPIADJ "mpiadj" 71 #define MATSBAIJ "sbaij" 72 #define MATSEQSBAIJ "seqsbaij" 73 #define MATMPISBAIJ "mpisbaij" 74 #define MATDAAD "daad" 75 #define MATMFFD "mffd" 76 #define MATNORMAL "normal" 77 #define MATNORMALHERMITIAN "normalh" 78 #define MATLRC "lrc" 79 #define MATSCATTER "scatter" 80 #define MATBLOCKMAT "blockmat" 81 #define MATCOMPOSITE "composite" 82 #define MATFFT "fft" 83 #define MATFFTW "fftw" 84 #define MATSEQCUFFT "seqcufft" 85 #define MATTRANSPOSEMAT "transpose" 86 #define MATSCHURCOMPLEMENT "schurcomplement" 87 #define MATPYTHON "python" 88 #define MATHYPRE "hypre" 89 #define MATHYPRESTRUCT "hyprestruct" 90 #define MATHYPRESSTRUCT "hypresstruct" 91 #define MATSUBMATRIX "submatrix" 92 #define MATLOCALREF "localref" 93 #define MATNEST "nest" 94 #define MATPREALLOCATOR "preallocator" 95 #define MATSELL "sell" 96 #define MATSEQSELL "seqsell" 97 #define MATMPISELL "mpisell" 98 #define MATDUMMY "dummy" 99 #define MATLMVM "lmvm" 100 #define MATLMVMDFP "lmvmdfp" 101 #define MATLMVMBFGS "lmvmbfgs" 102 #define MATLMVMSR1 "lmvmsr1" 103 #define MATLMVMBROYDEN "lmvmbroyden" 104 #define MATLMVMBADBROYDEN "lmvmbadbroyden" 105 #define MATLMVMSYMBROYDEN "lmvmsymbroyden" 106 #define MATLMVMSYMBADBROYDEN "lmvmsymbadbroyden" 107 #define MATLMVMDIAGBROYDEN "lmvmdiagbroyden" 108 #define MATCONSTANTDIAGONAL "constantdiagonal" 109 #define MATHARA "hara" 110 111 /*J 112 MatSolverType - String with the name of a PETSc matrix solver type. 113 114 For example: "petsc" indicates what PETSc provides, "superlu_dist" the parallel SuperLU_DIST package etc 115 116 Level: beginner 117 118 Notes: MATSOLVERUMFPACK, MATSOLVERCHOLMOD, MATSOLVERKLU form the SuiteSparse package for which you can use --download-suitesparse 119 120 .seealso: MatGetFactor(), PCFactorSetMatSolverType(), PCFactorGetMatSolverType() 121 J*/ 122 typedef const char* MatSolverType; 123 #define MATSOLVERSUPERLU "superlu" 124 #define MATSOLVERSUPERLU_DIST "superlu_dist" 125 #define MATSOLVERSTRUMPACK "strumpack" 126 #define MATSOLVERUMFPACK "umfpack" 127 #define MATSOLVERCHOLMOD "cholmod" 128 #define MATSOLVERKLU "klu" 129 #define MATSOLVERSPARSEELEMENTAL "sparseelemental" 130 #define MATSOLVERELEMENTAL "elemental" 131 #define MATSOLVERSCALAPACK "scalapack" 132 #define MATSOLVERESSL "essl" 133 #define MATSOLVERLUSOL "lusol" 134 #define MATSOLVERMUMPS "mumps" 135 #define MATSOLVERMKL_PARDISO "mkl_pardiso" 136 #define MATSOLVERMKL_CPARDISO "mkl_cpardiso" 137 #define MATSOLVERPASTIX "pastix" 138 #define MATSOLVERMATLAB "matlab" 139 #define MATSOLVERPETSC "petsc" 140 #define MATSOLVERBAS "bas" 141 #define MATSOLVERCUSPARSE "cusparse" 142 #define MATSOLVERCUDA "cuda" 143 144 /*E 145 MatFactorType - indicates what type of factorization is requested 146 147 Level: beginner 148 149 Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 150 151 .seealso: MatSolverType, MatGetFactor() 152 E*/ 153 typedef enum {MAT_FACTOR_NONE, MAT_FACTOR_LU, MAT_FACTOR_CHOLESKY, MAT_FACTOR_ILU, MAT_FACTOR_ICC,MAT_FACTOR_ILUDT} MatFactorType; 154 PETSC_EXTERN const char *const MatFactorTypes[]; 155 156 PETSC_EXTERN PetscErrorCode MatGetFactor(Mat,MatSolverType,MatFactorType,Mat*); 157 PETSC_EXTERN PetscErrorCode MatGetFactorAvailable(Mat,MatSolverType,MatFactorType,PetscBool *); 158 PETSC_EXTERN PetscErrorCode MatFactorGetSolverType(Mat,MatSolverType*); 159 PETSC_EXTERN PetscErrorCode MatGetFactorType(Mat,MatFactorType*); 160 PETSC_EXTERN PetscErrorCode MatSetFactorType(Mat,MatFactorType); 161 PETSC_EXTERN PetscErrorCode MatSolverTypeRegister(MatSolverType,MatType,MatFactorType,PetscErrorCode(*)(Mat,MatFactorType,Mat*)); 162 PETSC_EXTERN PetscErrorCode MatSolverTypeGet(MatSolverType,MatType,MatFactorType,PetscBool*,PetscBool*,PetscErrorCode (**)(Mat,MatFactorType,Mat*)); 163 typedef MatSolverType MatSolverPackage PETSC_DEPRECATED_TYPEDEF("Use MatSolverType (since version 3.9)"); 164 PETSC_DEPRECATED_FUNCTION("Use MatSolverTypeRegister() (since version 3.9)") PETSC_STATIC_INLINE PetscErrorCode MatSolverPackageRegister(MatSolverType stype,MatType mtype,MatFactorType ftype,PetscErrorCode(*f)(Mat,MatFactorType,Mat*)) 165 { return MatSolverTypeRegister(stype,mtype,ftype,f); } 166 PETSC_DEPRECATED_FUNCTION("Use MatSolverTypeGet() (since version 3.9)") PETSC_STATIC_INLINE PetscErrorCode MatSolverPackageGet(MatSolverType stype,MatType mtype,MatFactorType ftype,PetscBool *foundmtype,PetscBool *foundstype,PetscErrorCode(**f)(Mat,MatFactorType,Mat*)) 167 { return MatSolverTypeGet(stype,mtype,ftype,foundmtype,foundstype,f); } 168 169 /*E 170 MatProductType - indicates what type of matrix product is requested 171 172 Level: beginner 173 174 .seealso: MatSolverType, MatProductSetType() 175 E*/ 176 typedef enum {MATPRODUCT_UNSPECIFIED=0,MATPRODUCT_AB,MATPRODUCT_AtB,MATPRODUCT_ABt,MATPRODUCT_PtAP,MATPRODUCT_RARt,MATPRODUCT_ABC} MatProductType; 177 PETSC_EXTERN const char *const MatProductTypes[]; 178 179 /*J 180 MatProductAlgorithm - String with the name of an algorithm for a PETSc matrix product implementation 181 182 Level: beginner 183 184 .seealso: MatSetType(), Mat, MatSolverType, MatRegister(), MatProductSetAlgorithm(), MatProductType 185 J*/ 186 typedef const char* MatProductAlgorithm; 187 #define MATPRODUCTALGORITHM_DEFAULT "default" 188 189 PETSC_EXTERN PetscErrorCode MatProductCreate(Mat,Mat,Mat,Mat*); 190 PETSC_EXTERN PetscErrorCode MatProductCreateWithMat(Mat,Mat,Mat,Mat); 191 PETSC_EXTERN PetscErrorCode MatProductSetType(Mat,MatProductType); 192 PETSC_EXTERN PetscErrorCode MatProductSetAlgorithm(Mat,MatProductAlgorithm); 193 PETSC_EXTERN PetscErrorCode MatProductSetFill(Mat,PetscReal); 194 PETSC_EXTERN PetscErrorCode MatProductSetFromOptions(Mat); 195 PETSC_EXTERN PetscErrorCode MatProductSymbolic(Mat); 196 PETSC_EXTERN PetscErrorCode MatProductNumeric(Mat); 197 PETSC_EXTERN PetscErrorCode MatProductReplaceMats(Mat,Mat,Mat,Mat); 198 PETSC_EXTERN PetscErrorCode MatProductClear(Mat); 199 PETSC_EXTERN PetscErrorCode MatProductView(Mat,PetscViewer); 200 201 /* Logging support */ 202 #define MAT_FILE_CLASSID 1211216 /* used to indicate matrices in binary files */ 203 PETSC_EXTERN PetscClassId MAT_CLASSID; 204 PETSC_EXTERN PetscClassId MAT_COLORING_CLASSID; 205 PETSC_EXTERN PetscClassId MAT_FDCOLORING_CLASSID; 206 PETSC_EXTERN PetscClassId MAT_TRANSPOSECOLORING_CLASSID; 207 PETSC_EXTERN PetscClassId MAT_PARTITIONING_CLASSID; 208 PETSC_EXTERN PetscClassId MAT_COARSEN_CLASSID; 209 PETSC_EXTERN PetscClassId MAT_NULLSPACE_CLASSID; 210 PETSC_EXTERN PetscClassId MATMFFD_CLASSID; 211 212 /*E 213 MatReuse - Indicates if matrices obtained from a previous call to MatCreateSubMatrices(), MatCreateSubMatrix(), MatConvert() or several other functions 214 are to be reused to store the new matrix values. 215 216 $ MAT_INITIAL_MATRIX - create a new matrix 217 $ MAT_REUSE_MATRIX - reuse the matrix created with a previous call that used MAT_INITIAL_MATRIX 218 $ MAT_INPLACE_MATRIX - replace the first input matrix with the new matrix (not applicable to all functions) 219 $ MAT_IGNORE_MATRIX - do not create a new matrix or reuse a give matrix, just ignore that matrix argument (not applicable to all functions) 220 221 Level: beginner 222 223 Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 224 225 .seealso: MatCreateSubMatrices(), MatCreateSubMatrix(), MatDestroyMatrices(), MatConvert() 226 E*/ 227 typedef enum {MAT_INITIAL_MATRIX,MAT_REUSE_MATRIX,MAT_IGNORE_MATRIX,MAT_INPLACE_MATRIX} MatReuse; 228 229 /*E 230 MatCreateSubMatrixOption - Indicates if matrices obtained from a call to MatCreateSubMatrices() 231 include the matrix values. Currently it is only used by MatGetSeqNonzerostructure(). 232 233 Level: beginner 234 235 .seealso: MatGetSeqNonzerostructure() 236 E*/ 237 typedef enum {MAT_DO_NOT_GET_VALUES,MAT_GET_VALUES} MatCreateSubMatrixOption; 238 239 PETSC_EXTERN PetscErrorCode MatInitializePackage(void); 240 241 PETSC_EXTERN PetscErrorCode MatCreate(MPI_Comm,Mat*); 242 PETSC_EXTERN PetscErrorCode MatSetSizes(Mat,PetscInt,PetscInt,PetscInt,PetscInt); 243 PETSC_EXTERN PetscErrorCode MatSetType(Mat,MatType); 244 PETSC_EXTERN PetscErrorCode MatGetVecType(Mat,VecType*); 245 PETSC_EXTERN PetscErrorCode MatSetVecType(Mat,VecType); 246 PETSC_EXTERN PetscErrorCode MatSetFromOptions(Mat); 247 PETSC_EXTERN PetscErrorCode MatViewFromOptions(Mat,PetscObject,const char[]); 248 PETSC_EXTERN PetscErrorCode MatRegister(const char[],PetscErrorCode(*)(Mat)); 249 PETSC_EXTERN PetscErrorCode MatRegisterRootName(const char[],const char[],const char[]); 250 PETSC_EXTERN PetscErrorCode MatSetOptionsPrefix(Mat,const char[]); 251 PETSC_EXTERN PetscErrorCode MatAppendOptionsPrefix(Mat,const char[]); 252 PETSC_EXTERN PetscErrorCode MatGetOptionsPrefix(Mat,const char*[]); 253 PETSC_EXTERN PetscErrorCode MatSetErrorIfFailure(Mat,PetscBool); 254 255 PETSC_EXTERN PetscFunctionList MatList; 256 PETSC_EXTERN PetscFunctionList MatColoringList; 257 PETSC_EXTERN PetscFunctionList MatPartitioningList; 258 259 /*E 260 MatStructure - Indicates if two matrices have the same nonzero structure 261 262 Level: beginner 263 264 Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 265 266 .seealso: MatCopy(), MatAXPY() 267 E*/ 268 typedef enum {DIFFERENT_NONZERO_PATTERN,SUBSET_NONZERO_PATTERN,SAME_NONZERO_PATTERN} MatStructure; 269 270 #if defined PETSC_HAVE_MKL_SPARSE 271 PETSC_EXTERN PetscErrorCode MatCreateSeqAIJMKL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 272 PETSC_EXTERN PetscErrorCode MatCreateMPIAIJMKL(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 273 PETSC_EXTERN PetscErrorCode MatCreateBAIJMKL(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 274 PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJMKL(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 275 #endif 276 277 PETSC_EXTERN PetscErrorCode MatCreateSeqSELL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 278 PETSC_EXTERN PetscErrorCode MatCreateSELL(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 279 PETSC_EXTERN PetscErrorCode MatSeqSELLSetPreallocation(Mat,PetscInt,const PetscInt[]); 280 PETSC_EXTERN PetscErrorCode MatMPISELLSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 281 282 PETSC_EXTERN PetscErrorCode MatCreateSeqDense(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*); 283 PETSC_EXTERN PetscErrorCode MatCreateDense(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*); 284 PETSC_EXTERN PetscErrorCode MatCreateSeqAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 285 PETSC_EXTERN PetscErrorCode MatCreateAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 286 PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *); 287 PETSC_EXTERN PetscErrorCode MatUpdateMPIAIJWithArrays(Mat,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 288 PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithSplitArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],PetscInt[],PetscInt[],PetscScalar[],Mat*); 289 PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithSeqAIJ(MPI_Comm,Mat,Mat,const PetscInt[],Mat*); 290 291 PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 292 PETSC_EXTERN PetscErrorCode MatCreateBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 293 PETSC_EXTERN PetscErrorCode MatCreateMPIBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat*); 294 295 PETSC_EXTERN PetscErrorCode MatCreateMPIAdj(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscInt[],Mat*); 296 PETSC_EXTERN PetscErrorCode MatCreateSeqSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 297 298 PETSC_EXTERN PetscErrorCode MatCreateSBAIJ(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 299 PETSC_EXTERN PetscErrorCode MatCreateMPISBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[],Mat *); 300 PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 301 PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 302 PETSC_EXTERN PetscErrorCode MatXAIJSetPreallocation(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscInt[],const PetscInt[]); 303 304 PETSC_EXTERN PetscErrorCode MatCreateShell(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,void *,Mat*); 305 PETSC_EXTERN PetscErrorCode MatCreateNormal(Mat,Mat*); 306 PETSC_EXTERN PetscErrorCode MatCreateNormalHermitian(Mat,Mat*); 307 PETSC_EXTERN PetscErrorCode MatCreateLRC(Mat,Mat,Vec,Mat,Mat*); 308 PETSC_EXTERN PetscErrorCode MatLRCGetMats(Mat,Mat*,Mat*,Vec*,Mat*); 309 PETSC_EXTERN PetscErrorCode MatCreateIS(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,ISLocalToGlobalMapping,ISLocalToGlobalMapping,Mat*); 310 PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 311 PETSC_EXTERN PetscErrorCode MatCreateMPIAIJCRL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 312 313 PETSC_EXTERN PetscErrorCode MatCreateScatter(MPI_Comm,VecScatter,Mat*); 314 PETSC_EXTERN PetscErrorCode MatScatterSetVecScatter(Mat,VecScatter); 315 PETSC_EXTERN PetscErrorCode MatScatterGetVecScatter(Mat,VecScatter*); 316 PETSC_EXTERN PetscErrorCode MatCreateBlockMat(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,Mat*); 317 PETSC_EXTERN PetscErrorCode MatCompositeAddMat(Mat,Mat); 318 PETSC_EXTERN PetscErrorCode MatCompositeMerge(Mat); 319 typedef enum {MAT_COMPOSITE_MERGE_RIGHT,MAT_COMPOSITE_MERGE_LEFT} MatCompositeMergeType; 320 PETSC_EXTERN PetscErrorCode MatCompositeSetMergeType(Mat,MatCompositeMergeType); 321 PETSC_EXTERN PetscErrorCode MatCreateComposite(MPI_Comm,PetscInt,const Mat*,Mat*); 322 typedef enum {MAT_COMPOSITE_ADDITIVE,MAT_COMPOSITE_MULTIPLICATIVE} MatCompositeType; 323 PETSC_EXTERN PetscErrorCode MatCompositeSetType(Mat,MatCompositeType); 324 PETSC_EXTERN PetscErrorCode MatCompositeGetType(Mat,MatCompositeType*); 325 PETSC_EXTERN PetscErrorCode MatCompositeSetMatStructure(Mat,MatStructure); 326 PETSC_EXTERN PetscErrorCode MatCompositeGetMatStructure(Mat,MatStructure*); 327 PETSC_EXTERN PetscErrorCode MatCompositeGetNumberMat(Mat,PetscInt*); 328 PETSC_EXTERN PetscErrorCode MatCompositeGetMat(Mat,PetscInt,Mat*); 329 PETSC_EXTERN PetscErrorCode MatCompositeSetScalings(Mat,const PetscScalar*); 330 331 PETSC_EXTERN PetscErrorCode MatCreateFFT(MPI_Comm,PetscInt,const PetscInt[],MatType,Mat*); 332 PETSC_EXTERN PetscErrorCode MatCreateSeqCUFFT(MPI_Comm,PetscInt,const PetscInt[],Mat*); 333 334 PETSC_EXTERN PetscErrorCode MatCreateTranspose(Mat,Mat*); 335 PETSC_EXTERN PetscErrorCode MatTransposeGetMat(Mat,Mat*); 336 PETSC_EXTERN PetscErrorCode MatCreateHermitianTranspose(Mat,Mat*); 337 PETSC_EXTERN PetscErrorCode MatHermitianTransposeGetMat(Mat,Mat*); 338 PETSC_EXTERN PetscErrorCode MatCreateSubMatrixVirtual(Mat,IS,IS,Mat*); 339 PETSC_EXTERN PetscErrorCode MatSubMatrixVirtualUpdate(Mat,Mat,IS,IS); 340 PETSC_EXTERN PetscErrorCode MatCreateLocalRef(Mat,IS,IS,Mat*); 341 PETSC_EXTERN PetscErrorCode MatCreateConstantDiagonal(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar,Mat*); 342 343 #if defined(PETSC_HAVE_HYPRE) 344 PETSC_EXTERN PetscErrorCode MatHYPRESetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 345 #endif 346 347 PETSC_EXTERN PetscErrorCode MatPythonSetType(Mat,const char[]); 348 349 PETSC_EXTERN PetscErrorCode MatResetPreallocation(Mat); 350 PETSC_EXTERN PetscErrorCode MatSetUp(Mat); 351 PETSC_EXTERN PetscErrorCode MatDestroy(Mat*); 352 PETSC_EXTERN PetscErrorCode MatGetNonzeroState(Mat,PetscObjectState*); 353 354 PETSC_EXTERN PetscErrorCode MatConjugate(Mat); 355 PETSC_EXTERN PetscErrorCode MatRealPart(Mat); 356 PETSC_EXTERN PetscErrorCode MatImaginaryPart(Mat); 357 PETSC_EXTERN PetscErrorCode MatGetDiagonalBlock(Mat,Mat*); 358 PETSC_EXTERN PetscErrorCode MatGetTrace(Mat,PetscScalar*); 359 PETSC_EXTERN PetscErrorCode MatInvertBlockDiagonal(Mat,const PetscScalar **); 360 PETSC_EXTERN PetscErrorCode MatInvertVariableBlockDiagonal(Mat,PetscInt,const PetscInt*,PetscScalar*); 361 PETSC_EXTERN PetscErrorCode MatInvertBlockDiagonalMat(Mat,Mat); 362 363 /* ------------------------------------------------------------*/ 364 PETSC_EXTERN PetscErrorCode MatSetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 365 PETSC_EXTERN PetscErrorCode MatSetValuesBlocked(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 366 PETSC_EXTERN PetscErrorCode MatSetValuesRow(Mat,PetscInt,const PetscScalar[]); 367 PETSC_EXTERN PetscErrorCode MatSetValuesRowLocal(Mat,PetscInt,const PetscScalar[]); 368 PETSC_EXTERN PetscErrorCode MatSetValuesBatch(Mat,PetscInt,PetscInt,PetscInt[],const PetscScalar[]); 369 PETSC_EXTERN PetscErrorCode MatSetRandom(Mat,PetscRandom); 370 371 /*S 372 MatStencil - Data structure (C struct) for storing information about a single row or 373 column of a matrix as indexed on an associated grid. These are arguments to MatSetStencil() and MatSetBlockStencil() 374 375 The i,j, and k represent the logical coordinates over the entire grid (for 2 and 1 dimensional problems the k and j entries are ignored). 376 The c represents the the degrees of freedom at each grid point (the dof argument to DMDASetDOF()). If dof is 1 then this entry is ignored. 377 378 For stencil access to vectors see DMDAVecGetArray(), DMDAVecGetArrayF90(). 379 380 Fortran usage is different, see MatSetValuesStencil() for details. 381 382 Level: beginner 383 384 .seealso: MatSetValuesStencil(), MatSetStencil(), MatSetValuesBlockedStencil(), DMDAVecGetArray(), DMDAVecGetArrayF90() 385 S*/ 386 typedef struct { 387 PetscInt k,j,i,c; 388 } MatStencil; 389 390 PETSC_EXTERN PetscErrorCode MatSetValuesStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode); 391 PETSC_EXTERN PetscErrorCode MatSetValuesBlockedStencil(Mat,PetscInt,const MatStencil[],PetscInt,const MatStencil[],const PetscScalar[],InsertMode); 392 PETSC_EXTERN PetscErrorCode MatSetStencil(Mat,PetscInt,const PetscInt[],const PetscInt[],PetscInt); 393 394 /*E 395 MatAssemblyType - Indicates if the matrix is now to be used, or if you plan 396 to continue to add or insert values to it 397 398 Level: beginner 399 400 .seealso: MatAssemblyBegin(), MatAssemblyEnd() 401 E*/ 402 typedef enum {MAT_FLUSH_ASSEMBLY=1,MAT_FINAL_ASSEMBLY=0} MatAssemblyType; 403 PETSC_EXTERN PetscErrorCode MatAssemblyBegin(Mat,MatAssemblyType); 404 PETSC_EXTERN PetscErrorCode MatAssemblyEnd(Mat,MatAssemblyType); 405 PETSC_EXTERN PetscErrorCode MatAssembled(Mat,PetscBool *); 406 407 408 409 /*E 410 MatOption - Options that may be set for a matrix and its behavior or storage 411 412 Level: beginner 413 414 Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 415 Any additions/changes here must also be made in src/mat/interface/dlregismat.c in MatOptions[] 416 417 Developer Notes: 418 Entries that are negative need not be called collectively by all processes. 419 420 .seealso: MatSetOption() 421 E*/ 422 typedef enum {MAT_OPTION_MIN = -3, 423 MAT_UNUSED_NONZERO_LOCATION_ERR = -2, 424 MAT_ROW_ORIENTED = -1, 425 MAT_SYMMETRIC = 1, 426 MAT_STRUCTURALLY_SYMMETRIC = 2, 427 MAT_NEW_DIAGONALS = 3, 428 MAT_IGNORE_OFF_PROC_ENTRIES = 4, 429 MAT_USE_HASH_TABLE = 5, 430 MAT_KEEP_NONZERO_PATTERN = 6, 431 MAT_IGNORE_ZERO_ENTRIES = 7, 432 MAT_USE_INODES = 8, 433 MAT_HERMITIAN = 9, 434 MAT_SYMMETRY_ETERNAL = 10, 435 MAT_NEW_NONZERO_LOCATION_ERR = 11, 436 MAT_IGNORE_LOWER_TRIANGULAR = 12, 437 MAT_ERROR_LOWER_TRIANGULAR = 13, 438 MAT_GETROW_UPPERTRIANGULAR = 14, 439 MAT_SPD = 15, 440 MAT_NO_OFF_PROC_ZERO_ROWS = 16, 441 MAT_NO_OFF_PROC_ENTRIES = 17, 442 MAT_NEW_NONZERO_LOCATIONS = 18, 443 MAT_NEW_NONZERO_ALLOCATION_ERR = 19, 444 MAT_SUBSET_OFF_PROC_ENTRIES = 20, 445 MAT_SUBMAT_SINGLEIS = 21, 446 MAT_STRUCTURE_ONLY = 22, 447 MAT_SORTED_FULL = 23, 448 MAT_OPTION_MAX = 24} MatOption; 449 450 PETSC_EXTERN const char *const *MatOptions; 451 PETSC_EXTERN PetscErrorCode MatSetOption(Mat,MatOption,PetscBool); 452 PETSC_EXTERN PetscErrorCode MatGetOption(Mat,MatOption,PetscBool*); 453 PETSC_EXTERN PetscErrorCode MatPropagateSymmetryOptions(Mat,Mat); 454 PETSC_EXTERN PetscErrorCode MatGetType(Mat,MatType*); 455 456 PETSC_EXTERN PetscErrorCode MatGetValues(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]); 457 PETSC_EXTERN PetscErrorCode MatGetRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 458 PETSC_EXTERN PetscErrorCode MatRestoreRow(Mat,PetscInt,PetscInt *,const PetscInt *[],const PetscScalar*[]); 459 PETSC_EXTERN PetscErrorCode MatGetRowUpperTriangular(Mat); 460 PETSC_EXTERN PetscErrorCode MatRestoreRowUpperTriangular(Mat); 461 PETSC_EXTERN PetscErrorCode MatGetColumnVector(Mat,Vec,PetscInt); 462 PETSC_EXTERN PetscErrorCode MatSeqAIJGetArray(Mat,PetscScalar *[]); 463 PETSC_EXTERN PetscErrorCode MatSeqAIJGetArrayRead(Mat,const PetscScalar *[]); 464 PETSC_EXTERN PetscErrorCode MatSeqAIJRestoreArray(Mat,PetscScalar *[]); 465 PETSC_EXTERN PetscErrorCode MatSeqAIJRestoreArrayRead(Mat,const PetscScalar *[]); 466 PETSC_EXTERN PetscErrorCode MatSeqAIJGetMaxRowNonzeros(Mat,PetscInt*); 467 PETSC_EXTERN PetscErrorCode MatSeqAIJSetValuesLocalFast(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 468 PETSC_EXTERN PetscErrorCode MatSeqAIJSetType(Mat,MatType); 469 PETSC_EXTERN PetscErrorCode MatSeqAIJRegister(const char[],PetscErrorCode (*)(Mat,MatType,MatReuse,Mat *)); 470 PETSC_EXTERN PetscFunctionList MatSeqAIJList; 471 PETSC_EXTERN PetscErrorCode MatSeqBAIJGetArray(Mat,PetscScalar *[]); 472 PETSC_EXTERN PetscErrorCode MatSeqBAIJRestoreArray(Mat,PetscScalar *[]); 473 PETSC_EXTERN PetscErrorCode MatSeqSBAIJGetArray(Mat,PetscScalar *[]); 474 PETSC_EXTERN PetscErrorCode MatSeqSBAIJRestoreArray(Mat,PetscScalar *[]); 475 PETSC_EXTERN PetscErrorCode MatDenseGetArray(Mat,PetscScalar *[]); 476 PETSC_EXTERN PetscErrorCode MatDenseRestoreArray(Mat,PetscScalar *[]); 477 PETSC_EXTERN PetscErrorCode MatDensePlaceArray(Mat,const PetscScalar[]); 478 PETSC_EXTERN PetscErrorCode MatDenseReplaceArray(Mat,const PetscScalar[]); 479 PETSC_EXTERN PetscErrorCode MatDenseResetArray(Mat); 480 PETSC_EXTERN PetscErrorCode MatDenseGetArrayRead(Mat,const PetscScalar *[]); 481 PETSC_EXTERN PetscErrorCode MatDenseRestoreArrayRead(Mat,const PetscScalar *[]); 482 PETSC_EXTERN PetscErrorCode MatDenseGetArrayWrite(Mat,PetscScalar *[]); 483 PETSC_EXTERN PetscErrorCode MatDenseRestoreArrayWrite(Mat,PetscScalar *[]); 484 PETSC_EXTERN PetscErrorCode MatGetBlockSize(Mat,PetscInt *); 485 PETSC_EXTERN PetscErrorCode MatSetBlockSize(Mat,PetscInt); 486 PETSC_EXTERN PetscErrorCode MatGetBlockSizes(Mat,PetscInt *,PetscInt *); 487 PETSC_EXTERN PetscErrorCode MatSetBlockSizes(Mat,PetscInt,PetscInt); 488 PETSC_EXTERN PetscErrorCode MatSetBlockSizesFromMats(Mat,Mat,Mat); 489 PETSC_EXTERN PetscErrorCode MatSetVariableBlockSizes(Mat,PetscInt,PetscInt*); 490 PETSC_EXTERN PetscErrorCode MatGetVariableBlockSizes(Mat,PetscInt*,const PetscInt**); 491 492 PETSC_EXTERN PetscErrorCode MatDenseGetColumn(Mat,PetscInt,PetscScalar*[]); 493 PETSC_EXTERN PetscErrorCode MatDenseRestoreColumn(Mat,PetscScalar*[]); 494 PETSC_EXTERN PetscErrorCode MatDenseGetColumnVec(Mat,PetscInt,Vec*); 495 PETSC_EXTERN PetscErrorCode MatDenseRestoreColumnVec(Mat,PetscInt,Vec*); 496 PETSC_EXTERN PetscErrorCode MatDenseGetColumnVecRead(Mat,PetscInt,Vec*); 497 PETSC_EXTERN PetscErrorCode MatDenseRestoreColumnVecRead(Mat,PetscInt,Vec*); 498 PETSC_EXTERN PetscErrorCode MatDenseGetColumnVecWrite(Mat,PetscInt,Vec*); 499 PETSC_EXTERN PetscErrorCode MatDenseRestoreColumnVecWrite(Mat,PetscInt,Vec*); 500 PETSC_EXTERN PetscErrorCode MatDenseGetSubMatrix(Mat,PetscInt,PetscInt,Mat*); 501 PETSC_EXTERN PetscErrorCode MatDenseRestoreSubMatrix(Mat,Mat*); 502 503 PETSC_EXTERN PetscErrorCode MatMult(Mat,Vec,Vec); 504 PETSC_EXTERN PetscErrorCode MatMultDiagonalBlock(Mat,Vec,Vec); 505 PETSC_EXTERN PetscErrorCode MatMultAdd(Mat,Vec,Vec,Vec); 506 PETSC_EXTERN PetscErrorCode MatMultTranspose(Mat,Vec,Vec); 507 PETSC_EXTERN PetscErrorCode MatMultHermitianTranspose(Mat,Vec,Vec); 508 PETSC_EXTERN PetscErrorCode MatIsTranspose(Mat,Mat,PetscReal,PetscBool *); 509 PETSC_EXTERN PetscErrorCode MatIsHermitianTranspose(Mat,Mat,PetscReal,PetscBool *); 510 PETSC_EXTERN PetscErrorCode MatMultTransposeAdd(Mat,Vec,Vec,Vec); 511 PETSC_EXTERN PetscErrorCode MatMultHermitianTransposeAdd(Mat,Vec,Vec,Vec); 512 PETSC_EXTERN PetscErrorCode MatMultConstrained(Mat,Vec,Vec); 513 PETSC_EXTERN PetscErrorCode MatMultTransposeConstrained(Mat,Vec,Vec); 514 PETSC_EXTERN PetscErrorCode MatMatSolve(Mat,Mat,Mat); 515 PETSC_EXTERN PetscErrorCode MatMatSolveTranspose(Mat,Mat,Mat); 516 PETSC_EXTERN PetscErrorCode MatMatTransposeSolve(Mat,Mat,Mat); 517 PETSC_EXTERN PetscErrorCode MatResidual(Mat,Vec,Vec,Vec); 518 519 /*E 520 MatDuplicateOption - Indicates if a duplicated sparse matrix should have 521 its numerical values copied over or just its nonzero structure. 522 523 Level: beginner 524 525 Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 526 527 $ MAT_DO_NOT_COPY_VALUES - Create a matrix using the same nonzero pattern as the original matrix, 528 $ with zeros for the numerical values. 529 $ MAT_COPY_VALUES - Create a matrix with the same nonzero pattern as the original matrix 530 $ and with the same numerical values. 531 $ MAT_SHARE_NONZERO_PATTERN - Create a matrix that shares the nonzero structure with the previous matrix 532 $ and does not copy it, using zeros for the numerical values. The parent and 533 $ child matrices will share their index (i and j) arrays, and you cannot 534 $ insert new nonzero entries into either matrix. 535 536 Notes: 537 Many matrix types (including SeqAIJ) do not support the MAT_SHARE_NONZERO_PATTERN optimization; in 538 this case the behavior is as if MAT_DO_NOT_COPY_VALUES has been specified. 539 540 .seealso: MatDuplicate() 541 E*/ 542 typedef enum {MAT_DO_NOT_COPY_VALUES,MAT_COPY_VALUES,MAT_SHARE_NONZERO_PATTERN} MatDuplicateOption; 543 544 PETSC_EXTERN PetscErrorCode MatConvert(Mat,MatType,MatReuse,Mat*); 545 PETSC_EXTERN PetscErrorCode MatDuplicate(Mat,MatDuplicateOption,Mat*); 546 547 548 PETSC_EXTERN PetscErrorCode MatCopy(Mat,Mat,MatStructure); 549 PETSC_EXTERN PetscErrorCode MatView(Mat,PetscViewer); 550 PETSC_EXTERN PetscErrorCode MatIsSymmetric(Mat,PetscReal,PetscBool *); 551 PETSC_EXTERN PetscErrorCode MatIsStructurallySymmetric(Mat,PetscBool *); 552 PETSC_EXTERN PetscErrorCode MatIsHermitian(Mat,PetscReal,PetscBool *); 553 PETSC_EXTERN PetscErrorCode MatIsSymmetricKnown(Mat,PetscBool *,PetscBool *); 554 PETSC_EXTERN PetscErrorCode MatIsHermitianKnown(Mat,PetscBool *,PetscBool *); 555 PETSC_EXTERN PetscErrorCode MatMissingDiagonal(Mat,PetscBool *,PetscInt *); 556 PETSC_EXTERN PetscErrorCode MatLoad(Mat, PetscViewer); 557 558 PETSC_EXTERN PetscErrorCode MatGetRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool *); 559 PETSC_EXTERN PetscErrorCode MatRestoreRowIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool *); 560 PETSC_EXTERN PetscErrorCode MatGetColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool *); 561 PETSC_EXTERN PetscErrorCode MatRestoreColumnIJ(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool *); 562 563 /*S 564 MatInfo - Context of matrix information, used with MatGetInfo() 565 566 In Fortran this is simply a double precision array of dimension MAT_INFO_SIZE 567 568 Level: intermediate 569 570 .seealso: MatGetInfo(), MatInfoType 571 S*/ 572 typedef struct { 573 PetscLogDouble block_size; /* block size */ 574 PetscLogDouble nz_allocated,nz_used,nz_unneeded; /* number of nonzeros */ 575 PetscLogDouble memory; /* memory allocated */ 576 PetscLogDouble assemblies; /* number of matrix assemblies called */ 577 PetscLogDouble mallocs; /* number of mallocs during MatSetValues() */ 578 PetscLogDouble fill_ratio_given,fill_ratio_needed; /* fill ratio for LU/ILU */ 579 PetscLogDouble factor_mallocs; /* number of mallocs during factorization */ 580 } MatInfo; 581 582 /*E 583 MatInfoType - Indicates if you want information about the local part of the matrix, 584 the entire parallel matrix or the maximum over all the local parts. 585 586 Level: beginner 587 588 Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 589 590 .seealso: MatGetInfo(), MatInfo 591 E*/ 592 typedef enum {MAT_LOCAL=1,MAT_GLOBAL_MAX=2,MAT_GLOBAL_SUM=3} MatInfoType; 593 PETSC_EXTERN PetscErrorCode MatGetInfo(Mat,MatInfoType,MatInfo*); 594 PETSC_EXTERN PetscErrorCode MatGetDiagonal(Mat,Vec); 595 PETSC_EXTERN PetscErrorCode MatGetRowMax(Mat,Vec,PetscInt[]); 596 PETSC_EXTERN PetscErrorCode MatGetRowMin(Mat,Vec,PetscInt[]); 597 PETSC_EXTERN PetscErrorCode MatGetRowMaxAbs(Mat,Vec,PetscInt[]); 598 PETSC_EXTERN PetscErrorCode MatGetRowMinAbs(Mat,Vec,PetscInt[]); 599 PETSC_EXTERN PetscErrorCode MatGetRowSum(Mat,Vec); 600 PETSC_EXTERN PetscErrorCode MatTranspose(Mat,MatReuse,Mat*); 601 PETSC_EXTERN PetscErrorCode MatHermitianTranspose(Mat,MatReuse,Mat*); 602 PETSC_EXTERN PetscErrorCode MatPermute(Mat,IS,IS,Mat*); 603 PETSC_EXTERN PetscErrorCode MatDiagonalScale(Mat,Vec,Vec); 604 PETSC_EXTERN PetscErrorCode MatDiagonalSet(Mat,Vec,InsertMode); 605 606 PETSC_EXTERN PetscErrorCode MatEqual(Mat,Mat,PetscBool*); 607 PETSC_EXTERN PetscErrorCode MatMultEqual(Mat,Mat,PetscInt,PetscBool*); 608 PETSC_EXTERN PetscErrorCode MatMultAddEqual(Mat,Mat,PetscInt,PetscBool*); 609 PETSC_EXTERN PetscErrorCode MatMultTransposeEqual(Mat,Mat,PetscInt,PetscBool*); 610 PETSC_EXTERN PetscErrorCode MatMultTransposeAddEqual(Mat,Mat,PetscInt,PetscBool*); 611 PETSC_EXTERN PetscErrorCode MatMatMultEqual(Mat,Mat,Mat,PetscInt,PetscBool*); 612 PETSC_EXTERN PetscErrorCode MatTransposeMatMultEqual(Mat,Mat,Mat,PetscInt,PetscBool*); 613 PETSC_EXTERN PetscErrorCode MatMatTransposeMultEqual(Mat,Mat,Mat,PetscInt,PetscBool*); 614 PETSC_EXTERN PetscErrorCode MatPtAPMultEqual(Mat,Mat,Mat,PetscInt,PetscBool*); 615 PETSC_EXTERN PetscErrorCode MatRARtMultEqual(Mat,Mat,Mat,PetscInt,PetscBool*); 616 PETSC_EXTERN PetscErrorCode MatIsLinear(Mat,PetscInt,PetscBool*); 617 618 PETSC_EXTERN PetscErrorCode MatNorm(Mat,NormType,PetscReal*); 619 PETSC_EXTERN PetscErrorCode MatGetColumnNorms(Mat,NormType,PetscReal*); 620 PETSC_EXTERN PetscErrorCode MatZeroEntries(Mat); 621 PETSC_EXTERN PetscErrorCode MatZeroRows(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 622 PETSC_EXTERN PetscErrorCode MatZeroRowsIS(Mat,IS,PetscScalar,Vec,Vec); 623 PETSC_EXTERN PetscErrorCode MatZeroRowsStencil(Mat,PetscInt,const MatStencil [],PetscScalar,Vec,Vec); 624 PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsStencil(Mat,PetscInt,const MatStencil[],PetscScalar,Vec,Vec); 625 PETSC_EXTERN PetscErrorCode MatZeroRowsColumns(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 626 PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsIS(Mat,IS,PetscScalar,Vec,Vec); 627 628 PETSC_EXTERN PetscErrorCode MatGetSize(Mat,PetscInt*,PetscInt*); 629 PETSC_EXTERN PetscErrorCode MatGetLocalSize(Mat,PetscInt*,PetscInt*); 630 PETSC_EXTERN PetscErrorCode MatGetOwnershipRange(Mat,PetscInt*,PetscInt*); 631 PETSC_EXTERN PetscErrorCode MatGetOwnershipRanges(Mat,const PetscInt**); 632 PETSC_EXTERN PetscErrorCode MatGetOwnershipRangeColumn(Mat,PetscInt*,PetscInt*); 633 PETSC_EXTERN PetscErrorCode MatGetOwnershipRangesColumn(Mat,const PetscInt**); 634 PETSC_EXTERN PetscErrorCode MatGetOwnershipIS(Mat,IS*,IS*); 635 636 PETSC_EXTERN PetscErrorCode MatCreateSubMatrices(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]); 637 PETSC_DEPRECATED_FUNCTION("Use MatCreateSubMatrices() (since version 3.8)") PETSC_STATIC_INLINE PetscErrorCode MatGetSubMatrices(Mat mat,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *submat[]) {return MatCreateSubMatrices(mat,n,irow,icol,scall,submat);} 638 PETSC_EXTERN PetscErrorCode MatCreateSubMatricesMPI(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]); 639 PETSC_DEPRECATED_FUNCTION("Use MatCreateSubMatricesMPI() (since version 3.8)") PETSC_STATIC_INLINE PetscErrorCode MatGetSubMatricesMPI(Mat mat,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *submat[]) {return MatCreateSubMatricesMPI(mat,n,irow,icol,scall,submat);} 640 PETSC_EXTERN PetscErrorCode MatDestroyMatrices(PetscInt,Mat *[]); 641 PETSC_EXTERN PetscErrorCode MatDestroySubMatrices(PetscInt,Mat *[]); 642 PETSC_EXTERN PetscErrorCode MatCreateSubMatrix(Mat,IS,IS,MatReuse,Mat *); 643 PETSC_DEPRECATED_FUNCTION("Use MatCreateSubMatrix() (since version 3.8)") PETSC_STATIC_INLINE PetscErrorCode MatGetSubMatrix(Mat mat,IS isrow,IS iscol,MatReuse cll,Mat *newmat) {return MatCreateSubMatrix(mat,isrow,iscol,cll,newmat);} 644 PETSC_EXTERN PetscErrorCode MatGetLocalSubMatrix(Mat,IS,IS,Mat*); 645 PETSC_EXTERN PetscErrorCode MatRestoreLocalSubMatrix(Mat,IS,IS,Mat*); 646 PETSC_EXTERN PetscErrorCode MatGetSeqNonzeroStructure(Mat,Mat*); 647 PETSC_EXTERN PetscErrorCode MatDestroySeqNonzeroStructure(Mat*); 648 649 PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJ(MPI_Comm,Mat,PetscInt,PetscInt,MatReuse,Mat*); 650 PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJSymbolic(MPI_Comm,Mat,PetscInt,PetscInt,Mat*); 651 PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJNumeric(Mat,Mat); 652 PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMat(Mat,MatReuse,Mat*); 653 PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMatCondensed(Mat,MatReuse,IS*,IS*,Mat*); 654 PETSC_EXTERN PetscErrorCode MatGetBrowsOfAcols(Mat,Mat,MatReuse,IS*,IS*,Mat*); 655 PETSC_EXTERN PetscErrorCode MatGetGhosts(Mat, PetscInt *,const PetscInt *[]); 656 657 PETSC_EXTERN PetscErrorCode MatIncreaseOverlap(Mat,PetscInt,IS[],PetscInt); 658 PETSC_EXTERN PetscErrorCode MatIncreaseOverlapSplit(Mat mat,PetscInt n,IS is[],PetscInt ov); 659 PETSC_EXTERN PetscErrorCode MatMPIAIJSetUseScalableIncreaseOverlap(Mat,PetscBool); 660 661 PETSC_EXTERN PetscErrorCode MatMatMult(Mat,Mat,MatReuse,PetscReal,Mat*); 662 663 PETSC_EXTERN PetscErrorCode MatMatMatMult(Mat,Mat,Mat,MatReuse,PetscReal,Mat*); 664 PETSC_EXTERN PetscErrorCode MatGalerkin(Mat,Mat,Mat,MatReuse,PetscReal,Mat*); 665 666 PETSC_EXTERN PetscErrorCode MatPtAP(Mat,Mat,MatReuse,PetscReal,Mat*); 667 PETSC_EXTERN PetscErrorCode MatRARt(Mat,Mat,MatReuse,PetscReal,Mat*); 668 669 PETSC_EXTERN PetscErrorCode MatTransposeMatMult(Mat,Mat,MatReuse,PetscReal,Mat*); 670 PETSC_EXTERN PetscErrorCode MatMatTransposeMult(Mat,Mat,MatReuse,PetscReal,Mat*); 671 672 PETSC_EXTERN PetscErrorCode MatAXPY(Mat,PetscScalar,Mat,MatStructure); 673 PETSC_EXTERN PetscErrorCode MatAYPX(Mat,PetscScalar,Mat,MatStructure); 674 675 PETSC_EXTERN PetscErrorCode MatScale(Mat,PetscScalar); 676 PETSC_EXTERN PetscErrorCode MatShift(Mat,PetscScalar); 677 678 PETSC_EXTERN PetscErrorCode MatSetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping); 679 PETSC_EXTERN PetscErrorCode MatGetLocalToGlobalMapping(Mat,ISLocalToGlobalMapping*,ISLocalToGlobalMapping*); 680 PETSC_EXTERN PetscErrorCode MatGetLayouts(Mat,PetscLayout*,PetscLayout*); 681 PETSC_EXTERN PetscErrorCode MatZeroRowsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 682 PETSC_EXTERN PetscErrorCode MatZeroRowsLocalIS(Mat,IS,PetscScalar,Vec,Vec); 683 PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocal(Mat,PetscInt,const PetscInt [],PetscScalar,Vec,Vec); 684 PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocalIS(Mat,IS,PetscScalar,Vec,Vec); 685 PETSC_EXTERN PetscErrorCode MatGetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]); 686 PETSC_EXTERN PetscErrorCode MatSetValuesLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 687 PETSC_EXTERN PetscErrorCode MatSetValuesBlockedLocal(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 688 689 PETSC_EXTERN PetscErrorCode MatStashSetInitialSize(Mat,PetscInt,PetscInt); 690 PETSC_EXTERN PetscErrorCode MatStashGetInfo(Mat,PetscInt*,PetscInt*,PetscInt*,PetscInt*); 691 692 PETSC_EXTERN PetscErrorCode MatInterpolate(Mat,Vec,Vec); 693 PETSC_EXTERN PetscErrorCode MatInterpolateAdd(Mat,Vec,Vec,Vec); 694 PETSC_EXTERN PetscErrorCode MatRestrict(Mat,Vec,Vec); 695 PETSC_EXTERN PetscErrorCode MatCreateVecs(Mat,Vec*,Vec*); 696 PETSC_DEPRECATED_FUNCTION("Use MatCreateVecs() (since version 3.6)") PETSC_STATIC_INLINE PetscErrorCode MatGetVecs(Mat mat,Vec *x,Vec *y) {return MatCreateVecs(mat,x,y);} 697 PETSC_EXTERN PetscErrorCode MatCreateRedundantMatrix(Mat,PetscInt,MPI_Comm,MatReuse,Mat*); 698 PETSC_EXTERN PetscErrorCode MatGetMultiProcBlock(Mat,MPI_Comm,MatReuse,Mat*); 699 PETSC_EXTERN PetscErrorCode MatFindZeroDiagonals(Mat,IS*); 700 PETSC_EXTERN PetscErrorCode MatFindOffBlockDiagonalEntries(Mat,IS*); 701 PETSC_EXTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat(MPI_Comm,Mat,PetscInt,MatReuse,Mat*); 702 703 /*MC 704 MatSetValue - Set a single entry into a matrix. 705 706 Not collective 707 708 Synopsis: 709 #include <petscmat.h> 710 PetscErrorCode MatSetValue(Mat m,PetscInt row,PetscInt col,PetscScalar value,InsertMode mode) 711 712 Input Parameters: 713 + m - the matrix 714 . row - the row location of the entry 715 . col - the column location of the entry 716 . value - the value to insert 717 - mode - either INSERT_VALUES or ADD_VALUES 718 719 Notes: 720 For efficiency one should use MatSetValues() and set several or many 721 values simultaneously if possible. 722 723 Level: beginner 724 725 .seealso: MatSetValues(), MatSetValueLocal() 726 M*/ 727 PETSC_STATIC_INLINE PetscErrorCode MatSetValue(Mat v,PetscInt i,PetscInt j,PetscScalar va,InsertMode mode) {return MatSetValues(v,1,&i,1,&j,&va,mode);} 728 729 PETSC_STATIC_INLINE PetscErrorCode MatGetValue(Mat v,PetscInt i,PetscInt j,PetscScalar *va) {return MatGetValues(v,1,&i,1,&j,va);} 730 731 PETSC_STATIC_INLINE PetscErrorCode MatSetValueLocal(Mat v,PetscInt i,PetscInt j,PetscScalar va,InsertMode mode) {return MatSetValuesLocal(v,1,&i,1,&j,&va,mode);} 732 733 /*MC 734 MatPreallocateInitialize - Begins the block of code that will count the number of nonzeros per 735 row in a matrix providing the data that one can use to correctly preallocate the matrix. 736 737 Synopsis: 738 #include <petscmat.h> 739 PetscErrorCode MatPreallocateInitialize(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz) 740 741 Collective 742 743 Input Parameters: 744 + comm - the communicator that will share the eventually allocated matrix 745 . nrows - the number of LOCAL rows in the matrix 746 - ncols - the number of LOCAL columns in the matrix 747 748 Output Parameters: 749 + dnz - the array that will be passed to the matrix preallocation routines 750 - onz - the other array passed to the matrix preallocation routines 751 752 Level: intermediate 753 754 Notes: 755 See Users-Manual: ch_performance for more details. 756 757 Do not malloc or free dnz and onz, that is handled internally by these routines 758 759 This is a MACRO not a function because it has a leading { that is closed by PetscPreallocateFinalize(). 760 761 .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(), 762 MatPreallocateSymmetricSetLocalBlock() 763 M*/ 764 #define MatPreallocateInitialize(comm,nrows,ncols,dnz,onz) 0; \ 765 do { \ 766 PetscErrorCode _4_ierr; PetscInt __nrows = (nrows),__ncols = (ncols),__rstart,__start,__end; \ 767 _4_ierr = PetscCalloc2((size_t)__nrows,&dnz,(size_t)__nrows,&onz);CHKERRQ(_4_ierr); \ 768 __start = 0; __end = __start; \ 769 _4_ierr = MPI_Scan(&__ncols,&__end,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __start = __end - __ncols;\ 770 _4_ierr = MPI_Scan(&__nrows,&__rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(_4_ierr); __rstart = __rstart - __nrows; \ 771 do { } while(0) 772 773 /*MC 774 MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be 775 inserted using a local number of the rows and columns 776 777 Synopsis: 778 #include <petscmat.h> 779 PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 780 781 Not Collective 782 783 Input Parameters: 784 + map - the row mapping from local numbering to global numbering 785 . nrows - the number of rows indicated 786 . rows - the indices of the rows 787 . cmap - the column mapping from local to global numbering 788 . ncols - the number of columns in the matrix 789 . cols - the columns indicated 790 . dnz - the array that will be passed to the matrix preallocation routines 791 - onz - the other array passed to the matrix preallocation routines 792 793 Level: intermediate 794 795 Notes: 796 See Users-Manual: ch_performance for more details. 797 798 Do not malloc or free dnz and onz, that is handled internally by these routines 799 800 .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock() 801 MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocalRemoveDups() 802 M*/ 803 #define MatPreallocateSetLocal(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \ 804 do {\ 805 PetscInt __l;\ 806 _4_ierr = ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\ 807 _4_ierr = ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\ 808 for (__l=0;__l<nrows;__l++) {\ 809 _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 810 }\ 811 } while (0) 812 813 /*MC 814 MatPreallocateSetLocalRemoveDups - Indicates the locations (rows and columns) in the matrix where nonzeros will be 815 inserted using a local number of the rows and columns. This version removes any duplicate columns in cols 816 817 Synopsis: 818 #include <petscmat.h> 819 PetscErrorCode MatPreallocateSetLocalRemoveDups(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 820 821 Not Collective 822 823 Input Parameters: 824 + map - the row mapping from local numbering to global numbering 825 . nrows - the number of rows indicated 826 . rows - the indices of the rows (these values are mapped to the global values) 827 . cmap - the column mapping from local to global numbering 828 . ncols - the number of columns in the matrix (this value will be changed if duplicate columns are found) 829 . cols - the columns indicated (these values are mapped to the global values, they are then sorted and duplicates removed) 830 . dnz - the array that will be passed to the matrix preallocation routines 831 - onz - the other array passed to the matrix preallocation routines 832 833 Level: intermediate 834 835 Notes: 836 See Users-Manual: ch_performance for more details. 837 838 Do not malloc or free dnz and onz, that is handled internally by these routines 839 840 .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock() 841 MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocal() 842 M*/ 843 #define MatPreallocateSetLocalRemoveDups(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \ 844 do {\ 845 PetscInt __l;\ 846 _4_ierr = ISLocalToGlobalMappingApply(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\ 847 _4_ierr = ISLocalToGlobalMappingApply(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\ 848 _4_ierr = PetscSortRemoveDupsInt(&ncols,cols);CHKERRQ(_4_ierr);\ 849 for (__l=0;__l<nrows;__l++) {\ 850 _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 851 }\ 852 } while (0) 853 854 /*MC 855 MatPreallocateSetLocalBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be 856 inserted using a local number of the rows and columns 857 858 Synopsis: 859 #include <petscmat.h> 860 PetscErrorCode MatPreallocateSetLocalBlock(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 861 862 Not Collective 863 864 Input Parameters: 865 + map - the row mapping from local numbering to global numbering 866 . nrows - the number of rows indicated 867 . rows - the indices of the rows 868 . cmap - the column mapping from local to global numbering 869 . ncols - the number of columns in the matrix 870 . cols - the columns indicated 871 . dnz - the array that will be passed to the matrix preallocation routines 872 - onz - the other array passed to the matrix preallocation routines 873 874 Level: intermediate 875 876 Notes: 877 See Users-Manual: ch_performance for more details. 878 879 Do not malloc or free dnz and onz, that is handled internally by these routines 880 881 .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock() 882 MatPreallocateInitialize(), MatPreallocateSymmetricSetLocalBlock() 883 M*/ 884 #define MatPreallocateSetLocalBlock(rmap,nrows,rows,cmap,ncols,cols,dnz,onz) 0; \ 885 do {\ 886 PetscInt __l;\ 887 _4_ierr = ISLocalToGlobalMappingApplyBlock(rmap,nrows,rows,rows);CHKERRQ(_4_ierr);\ 888 _4_ierr = ISLocalToGlobalMappingApplyBlock(cmap,ncols,cols,cols);CHKERRQ(_4_ierr);\ 889 for (__l=0;__l<nrows;__l++) {\ 890 _4_ierr = MatPreallocateSet((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 891 }\ 892 } while (0) 893 894 /*MC 895 MatPreallocateSymmetricSetLocalBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be 896 inserted using a local number of the rows and columns 897 898 Synopsis: 899 #include <petscmat.h> 900 PetscErrorCode MatPreallocateSymmetricSetLocalBlock(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 901 902 Not Collective 903 904 Input Parameters: 905 + map - the mapping between local numbering and global numbering 906 . nrows - the number of rows indicated 907 . rows - the indices of the rows 908 . ncols - the number of columns in the matrix 909 . cols - the columns indicated 910 . dnz - the array that will be passed to the matrix preallocation routines 911 - onz - the other array passed to the matrix preallocation routines 912 913 Level: intermediate 914 915 Notes: 916 See Users-Manual: ch_performance for more details. 917 918 Do not malloc or free dnz and onz that is handled internally by these routines 919 920 .seealso: MatPreallocateFinalize(), MatPreallocateSet() 921 MatPreallocateInitialize(), MatPreallocateSetLocal() 922 M*/ 923 #define MatPreallocateSymmetricSetLocalBlock(map,nrows,rows,ncols,cols,dnz,onz) 0;\ 924 do {\ 925 PetscInt __l;\ 926 _4_ierr = ISLocalToGlobalMappingApplyBlock(map,nrows,rows,rows);CHKERRQ(_4_ierr);\ 927 _4_ierr = ISLocalToGlobalMappingApplyBlock(map,ncols,cols,cols);CHKERRQ(_4_ierr);\ 928 for (__l=0;__l<nrows;__l++) {\ 929 _4_ierr = MatPreallocateSymmetricSetBlock((rows)[__l],ncols,cols,dnz,onz);CHKERRQ(_4_ierr);\ 930 }\ 931 } while (0) 932 933 /*MC 934 MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be 935 inserted using a local number of the rows and columns 936 937 Synopsis: 938 #include <petscmat.h> 939 PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 940 941 Not Collective 942 943 Input Parameters: 944 + row - the row 945 . ncols - the number of columns in the matrix 946 - cols - the columns indicated 947 948 Output Parameters: 949 + dnz - the array that will be passed to the matrix preallocation routines 950 - onz - the other array passed to the matrix preallocation routines 951 952 Level: intermediate 953 954 Notes: 955 See Users-Manual: ch_performance for more details. 956 957 Do not malloc or free dnz and onz that is handled internally by these routines 958 959 This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize(). 960 961 .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock() 962 MatPreallocateInitialize(), MatPreallocateSetLocal() 963 M*/ 964 #define MatPreallocateSet(row,nc,cols,dnz,onz) 0;\ 965 do { PetscInt __i; \ 966 if (row < __rstart) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Trying to set preallocation for row %D less than first local row %D",row,__rstart);\ 967 if (row >= __rstart+__nrows) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Trying to set preallocation for row %D greater than last local row %D",row,__rstart+__nrows-1);\ 968 for (__i=0; __i<nc; __i++) {\ 969 if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \ 970 else if (dnz[row - __rstart] < __ncols) dnz[row - __rstart]++;\ 971 }\ 972 } while (0) 973 974 /*MC 975 MatPreallocateSymmetricSetBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be 976 inserted using a local number of the rows and columns 977 978 Synopsis: 979 #include <petscmat.h> 980 PetscErrorCode MatPreallocateSymmetricSetBlock(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 981 982 Not Collective 983 984 Input Parameters: 985 + nrows - the number of rows indicated 986 . rows - the indices of the rows 987 . ncols - the number of columns in the matrix 988 . cols - the columns indicated 989 . dnz - the array that will be passed to the matrix preallocation routines 990 - onz - the other array passed to the matrix preallocation routines 991 992 Level: intermediate 993 994 Notes: 995 See Users-Manual: ch_performance for more details. 996 997 Do not malloc or free dnz and onz that is handled internally by these routines 998 999 This is a MACRO not a function because it uses variables declared in MatPreallocateInitialize(). 1000 1001 .seealso: MatPreallocateFinalize(), MatPreallocateSet(), MatPreallocateInitialize(), 1002 MatPreallocateSymmetricSetLocalBlock(), MatPreallocateSetLocal() 1003 M*/ 1004 #define MatPreallocateSymmetricSetBlock(row,nc,cols,dnz,onz) 0;\ 1005 do { PetscInt __i; \ 1006 for (__i=0; __i<nc; __i++) {\ 1007 if (cols[__i] >= __end) onz[row - __rstart]++; \ 1008 else if (cols[__i] >= row && dnz[row - __rstart] < __ncols) dnz[row - __rstart]++;\ 1009 }\ 1010 } while (0) 1011 1012 /*MC 1013 MatPreallocateLocation - An alternative to MatPreallocateSet() that puts the nonzero locations into the matrix if it exists 1014 1015 Synopsis: 1016 #include <petscmat.h> 1017 PetscErrorCode MatPreallocateLocations(Mat A,PetscInt row,PetscInt ncols,PetscInt *cols,PetscInt *dnz,PetscInt *onz) 1018 1019 Not Collective 1020 1021 Input Parameters: 1022 + A - matrix 1023 . row - row where values exist (must be local to this process) 1024 . ncols - number of columns 1025 . cols - columns with nonzeros 1026 . dnz - the array that will be passed to the matrix preallocation routines 1027 - onz - the other array passed to the matrix preallocation routines 1028 1029 Level: intermediate 1030 1031 Notes: 1032 See Users-Manual: ch_performance for more details. 1033 1034 Do not malloc or free dnz and onz that is handled internally by these routines 1035 1036 This is a MACRO not a function because it uses a bunch of variables private to the MatPreallocation.... routines. 1037 1038 .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(), 1039 MatPreallocateSymmetricSetLocalBlock() 1040 M*/ 1041 #define MatPreallocateLocation(A,row,ncols,cols,dnz,onz) 0; do {if (A) {ierr = MatSetValues(A,1,&row,ncols,cols,NULL,INSERT_VALUES);CHKERRQ(ierr);} else {ierr = MatPreallocateSet(row,ncols,cols,dnz,onz);CHKERRQ(ierr);}} while(0) 1042 1043 1044 /*MC 1045 MatPreallocateFinalize - Ends the block of code that will count the number of nonzeros per 1046 row in a matrix providing the data that one can use to correctly preallocate the matrix. 1047 1048 Synopsis: 1049 #include <petscmat.h> 1050 PetscErrorCode MatPreallocateFinalize(PetscInt *dnz, PetscInt *onz) 1051 1052 Collective 1053 1054 Input Parameters: 1055 + dnz - the array that was be passed to the matrix preallocation routines 1056 - onz - the other array passed to the matrix preallocation routines 1057 1058 Level: intermediate 1059 1060 Notes: 1061 See Users-Manual: ch_performance for more details. 1062 1063 Do not malloc or free dnz and onz that is handled internally by these routines 1064 1065 This is a MACRO not a function because it closes the { started in MatPreallocateInitialize(). 1066 1067 .seealso: MatPreallocateInitialize(), MatPreallocateSet(), MatPreallocateSymmetricSetBlock(), MatPreallocateSetLocal(), 1068 MatPreallocateSymmetricSetLocalBlock() 1069 M*/ 1070 #define MatPreallocateFinalize(dnz,onz) 0;_4_ierr = PetscFree2(dnz,onz);CHKERRQ(_4_ierr);} while(0) 1071 1072 /* Routines unique to particular data structures */ 1073 PETSC_EXTERN PetscErrorCode MatShellGetContext(Mat,void *); 1074 1075 PETSC_EXTERN PetscErrorCode MatInodeAdjustForInodes(Mat,IS*,IS*); 1076 PETSC_EXTERN PetscErrorCode MatInodeGetInodeSizes(Mat,PetscInt *,PetscInt *[],PetscInt *); 1077 1078 PETSC_EXTERN PetscErrorCode MatSeqAIJSetColumnIndices(Mat,PetscInt[]); 1079 PETSC_EXTERN PetscErrorCode MatSeqBAIJSetColumnIndices(Mat,PetscInt[]); 1080 PETSC_EXTERN PetscErrorCode MatCreateSeqAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 1081 PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 1082 PETSC_EXTERN PetscErrorCode MatCreateSeqSBAIJWithArrays(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*); 1083 PETSC_EXTERN PetscErrorCode MatCreateSeqAIJFromTriple(MPI_Comm,PetscInt,PetscInt,PetscInt[],PetscInt[],PetscScalar[],Mat*,PetscInt,PetscBool); 1084 1085 #define MAT_SKIP_ALLOCATION -4 1086 1087 PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]); 1088 PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]); 1089 PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocation(Mat,PetscInt,const PetscInt[]); 1090 1091 PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 1092 PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 1093 PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 1094 PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocationCSR(Mat,const PetscInt [],const PetscInt [],const PetscScalar []); 1095 PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 1096 PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocationCSR(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]); 1097 PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocationCSR(Mat,PetscInt,const PetscInt[],const PetscInt[],const PetscScalar[]); 1098 PETSC_EXTERN PetscErrorCode MatMPIAdjSetPreallocation(Mat,PetscInt[],PetscInt[],PetscInt[]); 1099 PETSC_EXTERN PetscErrorCode MatMPIAdjToSeq(Mat,Mat*); 1100 PETSC_EXTERN PetscErrorCode MatMPIDenseSetPreallocation(Mat,PetscScalar[]); 1101 PETSC_EXTERN PetscErrorCode MatSeqDenseSetPreallocation(Mat,PetscScalar[]); 1102 PETSC_EXTERN PetscErrorCode MatMPIAIJGetSeqAIJ(Mat,Mat*,Mat*,const PetscInt*[]); 1103 PETSC_EXTERN PetscErrorCode MatMPIBAIJGetSeqBAIJ(Mat,Mat*,Mat*,const PetscInt*[]); 1104 PETSC_EXTERN PetscErrorCode MatMPIAdjCreateNonemptySubcommMat(Mat,Mat*); 1105 1106 PETSC_EXTERN PetscErrorCode MatDenseGetLDA(Mat,PetscInt*); 1107 PETSC_EXTERN PetscErrorCode MatDenseSetLDA(Mat,PetscInt); 1108 PETSC_DEPRECATED_FUNCTION("Use MatDenseSetLDA() (since version 3.14)") PETSC_STATIC_INLINE PetscErrorCode MatSeqDenseSetLDA(Mat A,PetscInt lda) {return MatDenseSetLDA(A,lda);} 1109 PETSC_EXTERN PetscErrorCode MatDenseGetLocalMatrix(Mat,Mat*); 1110 1111 PETSC_EXTERN PetscErrorCode MatBlockMatSetPreallocation(Mat,PetscInt,PetscInt,const PetscInt[]); 1112 1113 PETSC_EXTERN PetscErrorCode MatStoreValues(Mat); 1114 PETSC_EXTERN PetscErrorCode MatRetrieveValues(Mat); 1115 1116 PETSC_EXTERN PetscErrorCode MatDAADSetCtx(Mat,void*); 1117 1118 PETSC_EXTERN PetscErrorCode MatFindNonzeroRows(Mat,IS*); 1119 PETSC_EXTERN PetscErrorCode MatFindZeroRows(Mat,IS*); 1120 /* 1121 These routines are not usually accessed directly, rather solving is 1122 done through the KSP and PC interfaces. 1123 */ 1124 1125 /*J 1126 MatOrderingType - String with the name of a PETSc matrix ordering 1127 1128 Level: beginner 1129 1130 .seealso: MatGetOrdering() 1131 J*/ 1132 typedef const char* MatOrderingType; 1133 #define MATORDERINGNATURAL "natural" 1134 #define MATORDERINGND "nd" 1135 #define MATORDERING1WD "1wd" 1136 #define MATORDERINGRCM "rcm" 1137 #define MATORDERINGQMD "qmd" 1138 #define MATORDERINGROWLENGTH "rowlength" 1139 #define MATORDERINGWBM "wbm" 1140 #define MATORDERINGSPECTRAL "spectral" 1141 #define MATORDERINGAMD "amd" /* only works if UMFPACK is installed with PETSc */ 1142 #define MATORDERINGNATURAL_OR_ND "natural_or_nd" /* special coase used for Cholesky and ICC, allows ND when AIJ matrix is used but Natural when SBAIJ is used */ 1143 1144 PETSC_EXTERN PetscErrorCode MatGetOrdering(Mat,MatOrderingType,IS*,IS*); 1145 PETSC_EXTERN PetscErrorCode MatGetOrderingList(PetscFunctionList*); 1146 PETSC_EXTERN PetscErrorCode MatOrderingRegister(const char[],PetscErrorCode(*)(Mat,MatOrderingType,IS*,IS*)); 1147 PETSC_EXTERN PetscFunctionList MatOrderingList; 1148 1149 PETSC_EXTERN PetscErrorCode MatReorderForNonzeroDiagonal(Mat,PetscReal,IS,IS); 1150 PETSC_EXTERN PetscErrorCode MatCreateLaplacian(Mat,PetscReal,PetscBool,Mat*); 1151 1152 /*S 1153 MatFactorShiftType - Numeric Shift for factorizations 1154 1155 Level: beginner 1156 1157 .seealso: MatGetFactor() 1158 S*/ 1159 typedef enum {MAT_SHIFT_NONE,MAT_SHIFT_NONZERO,MAT_SHIFT_POSITIVE_DEFINITE,MAT_SHIFT_INBLOCKS} MatFactorShiftType; 1160 PETSC_EXTERN const char *const MatFactorShiftTypes[]; 1161 PETSC_EXTERN const char *const MatFactorShiftTypesDetail[]; 1162 1163 /*S 1164 MatFactorError - indicates what type of error was generated in a matrix factorization 1165 1166 Level: beginner 1167 1168 Developer Notes: 1169 Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 1170 1171 .seealso: MatGetFactor() 1172 S*/ 1173 typedef enum {MAT_FACTOR_NOERROR,MAT_FACTOR_STRUCT_ZEROPIVOT,MAT_FACTOR_NUMERIC_ZEROPIVOT,MAT_FACTOR_OUTMEMORY,MAT_FACTOR_OTHER} MatFactorError; 1174 1175 PETSC_EXTERN PetscErrorCode MatFactorGetError(Mat,MatFactorError*); 1176 PETSC_EXTERN PetscErrorCode MatFactorClearError(Mat); 1177 PETSC_EXTERN PetscErrorCode MatFactorGetErrorZeroPivot(Mat,PetscReal*,PetscInt*); 1178 1179 /*S 1180 MatFactorInfo - Data passed into the matrix factorization routines, and information about the resulting factorization 1181 1182 In Fortran these are simply double precision arrays of size MAT_FACTORINFO_SIZE, that is use 1183 $ MatFactorInfo info(MAT_FACTORINFO_SIZE) 1184 1185 Notes: 1186 These are not usually directly used by users, instead use PC type of LU, ILU, CHOLESKY or ICC. 1187 1188 You can use MatFactorInfoInitialize() to set default values. 1189 1190 Level: developer 1191 1192 .seealso: MatLUFactorSymbolic(), MatILUFactorSymbolic(), MatCholeskyFactorSymbolic(), MatICCFactorSymbolic(), MatICCFactor(), 1193 MatFactorInfoInitialize() 1194 1195 S*/ 1196 typedef struct { 1197 PetscReal diagonal_fill; /* force diagonal to fill in if initially not filled */ 1198 PetscReal usedt; 1199 PetscReal dt; /* drop tolerance */ 1200 PetscReal dtcol; /* tolerance for pivoting */ 1201 PetscReal dtcount; /* maximum nonzeros to be allowed per row */ 1202 PetscReal fill; /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */ 1203 PetscReal levels; /* ICC/ILU(levels) */ 1204 PetscReal pivotinblocks; /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0 1205 factorization may be faster if do not pivot */ 1206 PetscReal zeropivot; /* pivot is called zero if less than this */ 1207 PetscReal shifttype; /* type of shift added to matrix factor to prevent zero pivots */ 1208 PetscReal shiftamount; /* how large the shift is */ 1209 } MatFactorInfo; 1210 1211 PETSC_EXTERN PetscErrorCode MatFactorInfoInitialize(MatFactorInfo*); 1212 PETSC_EXTERN PetscErrorCode MatCholeskyFactor(Mat,IS,const MatFactorInfo*); 1213 PETSC_EXTERN PetscErrorCode MatCholeskyFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*); 1214 PETSC_EXTERN PetscErrorCode MatCholeskyFactorNumeric(Mat,Mat,const MatFactorInfo*); 1215 PETSC_EXTERN PetscErrorCode MatLUFactor(Mat,IS,IS,const MatFactorInfo*); 1216 PETSC_EXTERN PetscErrorCode MatILUFactor(Mat,IS,IS,const MatFactorInfo*); 1217 PETSC_EXTERN PetscErrorCode MatLUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*); 1218 PETSC_EXTERN PetscErrorCode MatILUFactorSymbolic(Mat,Mat,IS,IS,const MatFactorInfo*); 1219 PETSC_EXTERN PetscErrorCode MatICCFactorSymbolic(Mat,Mat,IS,const MatFactorInfo*); 1220 PETSC_EXTERN PetscErrorCode MatICCFactor(Mat,IS,const MatFactorInfo*); 1221 PETSC_EXTERN PetscErrorCode MatLUFactorNumeric(Mat,Mat,const MatFactorInfo*); 1222 PETSC_EXTERN PetscErrorCode MatGetInertia(Mat,PetscInt*,PetscInt*,PetscInt*); 1223 PETSC_EXTERN PetscErrorCode MatSolve(Mat,Vec,Vec); 1224 PETSC_EXTERN PetscErrorCode MatForwardSolve(Mat,Vec,Vec); 1225 PETSC_EXTERN PetscErrorCode MatBackwardSolve(Mat,Vec,Vec); 1226 PETSC_EXTERN PetscErrorCode MatSolveAdd(Mat,Vec,Vec,Vec); 1227 PETSC_EXTERN PetscErrorCode MatSolveTranspose(Mat,Vec,Vec); 1228 PETSC_EXTERN PetscErrorCode MatSolveTransposeAdd(Mat,Vec,Vec,Vec); 1229 PETSC_EXTERN PetscErrorCode MatSolves(Mat,Vecs,Vecs); 1230 PETSC_EXTERN PetscErrorCode MatSetUnfactored(Mat); 1231 1232 typedef enum {MAT_FACTOR_SCHUR_UNFACTORED, MAT_FACTOR_SCHUR_FACTORED, MAT_FACTOR_SCHUR_INVERTED} MatFactorSchurStatus; 1233 PETSC_EXTERN PetscErrorCode MatFactorSetSchurIS(Mat,IS); 1234 PETSC_EXTERN PetscErrorCode MatFactorGetSchurComplement(Mat,Mat*,MatFactorSchurStatus*); 1235 PETSC_EXTERN PetscErrorCode MatFactorRestoreSchurComplement(Mat,Mat*,MatFactorSchurStatus); 1236 PETSC_EXTERN PetscErrorCode MatFactorInvertSchurComplement(Mat); 1237 PETSC_EXTERN PetscErrorCode MatFactorCreateSchurComplement(Mat,Mat*,MatFactorSchurStatus*); 1238 PETSC_EXTERN PetscErrorCode MatFactorSolveSchurComplement(Mat,Vec,Vec); 1239 PETSC_EXTERN PetscErrorCode MatFactorSolveSchurComplementTranspose(Mat,Vec,Vec); 1240 PETSC_EXTERN PetscErrorCode MatFactorFactorizeSchurComplement(Mat); 1241 1242 /*E 1243 MatSORType - What type of (S)SOR to perform 1244 1245 Level: beginner 1246 1247 May be bitwise ORd together 1248 1249 Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 1250 1251 MatSORType may be bitwise ORd together, so do not change the numbers 1252 1253 .seealso: MatSOR() 1254 E*/ 1255 typedef enum {SOR_FORWARD_SWEEP=1,SOR_BACKWARD_SWEEP=2,SOR_SYMMETRIC_SWEEP=3, 1256 SOR_LOCAL_FORWARD_SWEEP=4,SOR_LOCAL_BACKWARD_SWEEP=8, 1257 SOR_LOCAL_SYMMETRIC_SWEEP=12,SOR_ZERO_INITIAL_GUESS=16, 1258 SOR_EISENSTAT=32,SOR_APPLY_UPPER=64,SOR_APPLY_LOWER=128} MatSORType; 1259 PETSC_EXTERN PetscErrorCode MatSOR(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec); 1260 1261 /* 1262 These routines are for efficiently computing Jacobians via finite differences. 1263 */ 1264 1265 /*S 1266 MatColoring - Object for managing the coloring of matrices. 1267 1268 Level: beginner 1269 1270 Notes: 1271 Coloring of matrices can be computed directly from the sparse matrix nonzero structure via the MatColoring object or from the mesh from which the 1272 matrix comes from via DMCreateColoring(). In general using the mesh produces a more optimal coloring (fewer colors). 1273 1274 Once a coloring is available MatFDColoringCreate() creates an object that can be used to efficiently compute Jacobians using that coloring. This 1275 same object can also be used to efficiently convert data created by Automatic Differentation tools to PETSc sparse matrices. 1276 1277 .seealso: MatFDColoringCreate(), MatColoringWeightType, ISColoring, MatFDColoring, DMCreateColoring(), MatColoringCreate(), MatOrdering, MatPartitioning 1278 S*/ 1279 typedef struct _p_MatColoring* MatColoring; 1280 1281 /*J 1282 MatColoringType - String with the name of a PETSc matrix coloring 1283 1284 Level: beginner 1285 1286 .seealso: MatColoringSetType(), MatColoring 1287 J*/ 1288 typedef const char* MatColoringType; 1289 #define MATCOLORINGJP "jp" 1290 #define MATCOLORINGPOWER "power" 1291 #define MATCOLORINGNATURAL "natural" 1292 #define MATCOLORINGSL "sl" 1293 #define MATCOLORINGLF "lf" 1294 #define MATCOLORINGID "id" 1295 #define MATCOLORINGGREEDY "greedy" 1296 1297 /*E 1298 MatColoringWeightType - Type of weight scheme 1299 1300 Not Collective 1301 1302 + MAT_COLORING_RANDOM - Random weights 1303 . MAT_COLORING_LEXICAL - Lexical weighting based upon global numbering. 1304 - MAT_COLORING_LF - Last-first weighting. 1305 1306 Level: intermediate 1307 1308 Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 1309 1310 .seealso: MatColoring, MatColoringCreate() 1311 E*/ 1312 typedef enum {MAT_COLORING_WEIGHT_RANDOM,MAT_COLORING_WEIGHT_LEXICAL,MAT_COLORING_WEIGHT_LF,MAT_COLORING_WEIGHT_SL} MatColoringWeightType; 1313 1314 PETSC_EXTERN PetscErrorCode MatColoringCreate(Mat,MatColoring*); 1315 PETSC_EXTERN PetscErrorCode MatColoringGetDegrees(Mat,PetscInt,PetscInt*); 1316 PETSC_EXTERN PetscErrorCode MatColoringDestroy(MatColoring*); 1317 PETSC_EXTERN PetscErrorCode MatColoringView(MatColoring,PetscViewer); 1318 PETSC_EXTERN PetscErrorCode MatColoringSetType(MatColoring,MatColoringType); 1319 PETSC_EXTERN PetscErrorCode MatColoringSetFromOptions(MatColoring); 1320 PETSC_EXTERN PetscErrorCode MatColoringSetDistance(MatColoring,PetscInt); 1321 PETSC_EXTERN PetscErrorCode MatColoringGetDistance(MatColoring,PetscInt*); 1322 PETSC_EXTERN PetscErrorCode MatColoringSetMaxColors(MatColoring,PetscInt); 1323 PETSC_EXTERN PetscErrorCode MatColoringGetMaxColors(MatColoring,PetscInt*); 1324 PETSC_EXTERN PetscErrorCode MatColoringApply(MatColoring,ISColoring*); 1325 PETSC_EXTERN PetscErrorCode MatColoringRegister(const char[],PetscErrorCode(*)(MatColoring)); 1326 PETSC_EXTERN PetscErrorCode MatColoringPatch(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*); 1327 PETSC_EXTERN PetscErrorCode MatColoringSetWeightType(MatColoring,MatColoringWeightType); 1328 PETSC_EXTERN PetscErrorCode MatColoringSetWeights(MatColoring,PetscReal*,PetscInt*); 1329 PETSC_EXTERN PetscErrorCode MatColoringCreateWeights(MatColoring,PetscReal **,PetscInt **lperm); 1330 PETSC_EXTERN PetscErrorCode MatColoringTest(MatColoring,ISColoring); 1331 PETSC_DEPRECATED_FUNCTION("Use MatColoringTest() (since version 3.10)") PETSC_STATIC_INLINE PetscErrorCode MatColoringTestValid(MatColoring matcoloring,ISColoring iscoloring) {return MatColoringTest(matcoloring,iscoloring);} 1332 PETSC_EXTERN PetscErrorCode MatISColoringTest(Mat,ISColoring); 1333 1334 /*S 1335 MatFDColoring - Object for computing a sparse Jacobian via finite differences 1336 and coloring 1337 1338 Level: beginner 1339 1340 Notes: 1341 This object is creating utilizing a coloring provided by the MatColoring object or DMCreateColoring() 1342 1343 .seealso: MatFDColoringCreate(), MatColoring, DMCreateColoring() 1344 S*/ 1345 typedef struct _p_MatFDColoring* MatFDColoring; 1346 1347 PETSC_EXTERN PetscErrorCode MatFDColoringCreate(Mat,ISColoring,MatFDColoring *); 1348 PETSC_EXTERN PetscErrorCode MatFDColoringDestroy(MatFDColoring*); 1349 PETSC_EXTERN PetscErrorCode MatFDColoringView(MatFDColoring,PetscViewer); 1350 PETSC_EXTERN PetscErrorCode MatFDColoringSetFunction(MatFDColoring,PetscErrorCode (*)(void),void*); 1351 PETSC_EXTERN PetscErrorCode MatFDColoringGetFunction(MatFDColoring,PetscErrorCode (**)(void),void**); 1352 PETSC_EXTERN PetscErrorCode MatFDColoringSetParameters(MatFDColoring,PetscReal,PetscReal); 1353 PETSC_EXTERN PetscErrorCode MatFDColoringSetFromOptions(MatFDColoring); 1354 PETSC_EXTERN PetscErrorCode MatFDColoringApply(Mat,MatFDColoring,Vec,void *); 1355 PETSC_EXTERN PetscErrorCode MatFDColoringSetF(MatFDColoring,Vec); 1356 PETSC_EXTERN PetscErrorCode MatFDColoringGetPerturbedColumns(MatFDColoring,PetscInt*,const PetscInt*[]); 1357 PETSC_EXTERN PetscErrorCode MatFDColoringSetUp(Mat,ISColoring,MatFDColoring); 1358 PETSC_EXTERN PetscErrorCode MatFDColoringSetBlockSize(MatFDColoring,PetscInt,PetscInt); 1359 PETSC_EXTERN PetscErrorCode MatFDColoringSetValues(Mat,MatFDColoring,const PetscScalar*); 1360 1361 /*S 1362 MatTransposeColoring - Object for computing a sparse matrix product C=A*B^T via coloring 1363 1364 Level: beginner 1365 1366 .seealso: MatTransposeColoringCreate() 1367 S*/ 1368 typedef struct _p_MatTransposeColoring* MatTransposeColoring; 1369 1370 PETSC_EXTERN PetscErrorCode MatTransposeColoringCreate(Mat,ISColoring,MatTransposeColoring *); 1371 PETSC_EXTERN PetscErrorCode MatTransColoringApplySpToDen(MatTransposeColoring,Mat,Mat); 1372 PETSC_EXTERN PetscErrorCode MatTransColoringApplyDenToSp(MatTransposeColoring,Mat,Mat); 1373 PETSC_EXTERN PetscErrorCode MatTransposeColoringDestroy(MatTransposeColoring*); 1374 1375 /* 1376 These routines are for partitioning matrices: currently used only 1377 for adjacency matrix, MatCreateMPIAdj(). 1378 */ 1379 1380 /*S 1381 MatPartitioning - Object for managing the partitioning of a matrix or graph 1382 1383 Level: beginner 1384 1385 Notes: 1386 There is also a PetscPartitioner object that provides the same functionality. It can utilize the MatPartitioning operations 1387 via PetscPartitionerSetType(p,PETSCPARTITIONERMATPARTITIONING) 1388 1389 Developers Note: 1390 It is an extra maintainance and documentation cost to have two objects with the same functionality. 1391 1392 .seealso: MatPartitioningCreate(), MatPartitioningType, MatColoring, MatGetOrdering() 1393 S*/ 1394 typedef struct _p_MatPartitioning* MatPartitioning; 1395 1396 /*J 1397 MatPartitioningType - String with the name of a PETSc matrix partitioning 1398 1399 Level: beginner 1400 dm 1401 .seealso: MatPartitioningCreate(), MatPartitioning 1402 J*/ 1403 typedef const char* MatPartitioningType; 1404 #define MATPARTITIONINGCURRENT "current" 1405 #define MATPARTITIONINGAVERAGE "average" 1406 #define MATPARTITIONINGSQUARE "square" 1407 #define MATPARTITIONINGPARMETIS "parmetis" 1408 #define MATPARTITIONINGCHACO "chaco" 1409 #define MATPARTITIONINGPARTY "party" 1410 #define MATPARTITIONINGPTSCOTCH "ptscotch" 1411 #define MATPARTITIONINGHIERARCH "hierarch" 1412 1413 PETSC_EXTERN PetscErrorCode MatPartitioningCreate(MPI_Comm,MatPartitioning*); 1414 PETSC_EXTERN PetscErrorCode MatPartitioningSetType(MatPartitioning,MatPartitioningType); 1415 PETSC_EXTERN PetscErrorCode MatPartitioningSetNParts(MatPartitioning,PetscInt); 1416 PETSC_EXTERN PetscErrorCode MatPartitioningSetAdjacency(MatPartitioning,Mat); 1417 PETSC_EXTERN PetscErrorCode MatPartitioningSetVertexWeights(MatPartitioning,const PetscInt[]); 1418 PETSC_EXTERN PetscErrorCode MatPartitioningSetPartitionWeights(MatPartitioning,const PetscReal []); 1419 PETSC_EXTERN PetscErrorCode MatPartitioningSetUseEdgeWeights(MatPartitioning,PetscBool); 1420 PETSC_EXTERN PetscErrorCode MatPartitioningGetUseEdgeWeights(MatPartitioning,PetscBool*); 1421 PETSC_EXTERN PetscErrorCode MatPartitioningApply(MatPartitioning,IS*); 1422 PETSC_EXTERN PetscErrorCode MatPartitioningImprove(MatPartitioning,IS*); 1423 PETSC_EXTERN PetscErrorCode MatPartitioningViewImbalance(MatPartitioning,IS); 1424 PETSC_EXTERN PetscErrorCode MatPartitioningApplyND(MatPartitioning,IS*); 1425 PETSC_EXTERN PetscErrorCode MatPartitioningDestroy(MatPartitioning*); 1426 PETSC_EXTERN PetscErrorCode MatPartitioningRegister(const char[],PetscErrorCode (*)(MatPartitioning)); 1427 PETSC_EXTERN PetscErrorCode MatPartitioningView(MatPartitioning,PetscViewer); 1428 PETSC_EXTERN PetscErrorCode MatPartitioningViewFromOptions(MatPartitioning,PetscObject,const char[]); 1429 PETSC_EXTERN PetscErrorCode MatPartitioningSetFromOptions(MatPartitioning); 1430 PETSC_EXTERN PetscErrorCode MatPartitioningGetType(MatPartitioning,MatPartitioningType*); 1431 1432 PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetRepartition(MatPartitioning part); 1433 PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetCoarseSequential(MatPartitioning); 1434 PETSC_EXTERN PetscErrorCode MatPartitioningParmetisGetEdgeCut(MatPartitioning, PetscInt *); 1435 1436 typedef enum { MP_CHACO_MULTILEVEL=1,MP_CHACO_SPECTRAL=2,MP_CHACO_LINEAR=4,MP_CHACO_RANDOM=5,MP_CHACO_SCATTERED=6 } MPChacoGlobalType; 1437 PETSC_EXTERN const char *const MPChacoGlobalTypes[]; 1438 typedef enum { MP_CHACO_KERNIGHAN=1,MP_CHACO_NONE=2 } MPChacoLocalType; 1439 PETSC_EXTERN const char *const MPChacoLocalTypes[]; 1440 typedef enum { MP_CHACO_LANCZOS=0,MP_CHACO_RQI=1 } MPChacoEigenType; 1441 PETSC_EXTERN const char *const MPChacoEigenTypes[]; 1442 1443 PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetGlobal(MatPartitioning,MPChacoGlobalType); 1444 PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetGlobal(MatPartitioning,MPChacoGlobalType*); 1445 PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetLocal(MatPartitioning,MPChacoLocalType); 1446 PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetLocal(MatPartitioning,MPChacoLocalType*); 1447 PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetCoarseLevel(MatPartitioning,PetscReal); 1448 PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenSolver(MatPartitioning,MPChacoEigenType); 1449 PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenSolver(MatPartitioning,MPChacoEigenType*); 1450 PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenTol(MatPartitioning,PetscReal); 1451 PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenTol(MatPartitioning,PetscReal*); 1452 PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenNumber(MatPartitioning,PetscInt); 1453 PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenNumber(MatPartitioning,PetscInt*); 1454 1455 #define MP_PARTY_OPT "opt" 1456 #define MP_PARTY_LIN "lin" 1457 #define MP_PARTY_SCA "sca" 1458 #define MP_PARTY_RAN "ran" 1459 #define MP_PARTY_GBF "gbf" 1460 #define MP_PARTY_GCF "gcf" 1461 #define MP_PARTY_BUB "bub" 1462 #define MP_PARTY_DEF "def" 1463 PETSC_EXTERN PetscErrorCode MatPartitioningPartySetGlobal(MatPartitioning,const char*); 1464 #define MP_PARTY_HELPFUL_SETS "hs" 1465 #define MP_PARTY_KERNIGHAN_LIN "kl" 1466 #define MP_PARTY_NONE "no" 1467 PETSC_EXTERN PetscErrorCode MatPartitioningPartySetLocal(MatPartitioning,const char*); 1468 PETSC_EXTERN PetscErrorCode MatPartitioningPartySetCoarseLevel(MatPartitioning,PetscReal); 1469 PETSC_EXTERN PetscErrorCode MatPartitioningPartySetBipart(MatPartitioning,PetscBool); 1470 PETSC_EXTERN PetscErrorCode MatPartitioningPartySetMatchOptimization(MatPartitioning,PetscBool); 1471 1472 typedef enum { MP_PTSCOTCH_DEFAULT,MP_PTSCOTCH_QUALITY,MP_PTSCOTCH_SPEED,MP_PTSCOTCH_BALANCE,MP_PTSCOTCH_SAFETY,MP_PTSCOTCH_SCALABILITY } MPPTScotchStrategyType; 1473 PETSC_EXTERN const char *const MPPTScotchStrategyTypes[]; 1474 1475 PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetImbalance(MatPartitioning,PetscReal); 1476 PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetImbalance(MatPartitioning,PetscReal*); 1477 PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetStrategy(MatPartitioning,MPPTScotchStrategyType); 1478 PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetStrategy(MatPartitioning,MPPTScotchStrategyType*); 1479 1480 /* 1481 * hierarchical partitioning 1482 */ 1483 PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalGetFineparts(MatPartitioning,IS*); 1484 PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalGetCoarseparts(MatPartitioning,IS*); 1485 PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalSetNcoarseparts(MatPartitioning,PetscInt); 1486 PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalSetNfineparts(MatPartitioning, PetscInt); 1487 1488 PETSC_EXTERN PetscErrorCode MatMeshToVertexGraph(Mat,PetscInt,Mat*); 1489 PETSC_EXTERN PetscErrorCode MatMeshToCellGraph(Mat,PetscInt,Mat*); 1490 1491 /* 1492 If you add entries here you must also add them to include/petscmat.h 1493 and src/mat/f90-mod/petscmat.h 1494 */ 1495 typedef enum { MATOP_SET_VALUES=0, 1496 MATOP_GET_ROW=1, 1497 MATOP_RESTORE_ROW=2, 1498 MATOP_MULT=3, 1499 MATOP_MULT_ADD=4, 1500 MATOP_MULT_TRANSPOSE=5, 1501 MATOP_MULT_TRANSPOSE_ADD=6, 1502 MATOP_SOLVE=7, 1503 MATOP_SOLVE_ADD=8, 1504 MATOP_SOLVE_TRANSPOSE=9, 1505 MATOP_SOLVE_TRANSPOSE_ADD=10, 1506 MATOP_LUFACTOR=11, 1507 MATOP_CHOLESKYFACTOR=12, 1508 MATOP_SOR=13, 1509 MATOP_TRANSPOSE=14, 1510 MATOP_GETINFO=15, 1511 MATOP_EQUAL=16, 1512 MATOP_GET_DIAGONAL=17, 1513 MATOP_DIAGONAL_SCALE=18, 1514 MATOP_NORM=19, 1515 MATOP_ASSEMBLY_BEGIN=20, 1516 MATOP_ASSEMBLY_END=21, 1517 MATOP_SET_OPTION=22, 1518 MATOP_ZERO_ENTRIES=23, 1519 MATOP_ZERO_ROWS=24, 1520 MATOP_LUFACTOR_SYMBOLIC=25, 1521 MATOP_LUFACTOR_NUMERIC=26, 1522 MATOP_CHOLESKY_FACTOR_SYMBOLIC=27, 1523 MATOP_CHOLESKY_FACTOR_NUMERIC=28, 1524 MATOP_SETUP_PREALLOCATION=29, 1525 MATOP_ILUFACTOR_SYMBOLIC=30, 1526 MATOP_ICCFACTOR_SYMBOLIC=31, 1527 MATOP_GET_DIAGONAL_BLOCK=32, 1528 MATOP_FREE_INTER_STRUCT=33, 1529 MATOP_DUPLICATE=34, 1530 MATOP_FORWARD_SOLVE=35, 1531 MATOP_BACKWARD_SOLVE=36, 1532 MATOP_ILUFACTOR=37, 1533 MATOP_ICCFACTOR=38, 1534 MATOP_AXPY=39, 1535 MATOP_CREATE_SUBMATRICES=40, 1536 MATOP_INCREASE_OVERLAP=41, 1537 MATOP_GET_VALUES=42, 1538 MATOP_COPY=43, 1539 MATOP_GET_ROW_MAX=44, 1540 MATOP_SCALE=45, 1541 MATOP_SHIFT=46, 1542 MATOP_DIAGONAL_SET=47, 1543 MATOP_ZERO_ROWS_COLUMNS=48, 1544 MATOP_SET_RANDOM=49, 1545 MATOP_GET_ROW_IJ=50, 1546 MATOP_RESTORE_ROW_IJ=51, 1547 MATOP_GET_COLUMN_IJ=52, 1548 MATOP_RESTORE_COLUMN_IJ=53, 1549 MATOP_FDCOLORING_CREATE=54, 1550 MATOP_COLORING_PATCH=55, 1551 MATOP_SET_UNFACTORED=56, 1552 MATOP_PERMUTE=57, 1553 MATOP_SET_VALUES_BLOCKED=58, 1554 MATOP_CREATE_SUBMATRIX=59, 1555 MATOP_DESTROY=60, 1556 MATOP_VIEW=61, 1557 MATOP_CONVERT_FROM=62, 1558 MATOP_MATMAT_MULT=63, 1559 MATOP_MATMAT_MULT_SYMBOLIC=64, 1560 MATOP_MATMAT_MULT_NUMERIC=65, 1561 MATOP_SET_LOCAL_TO_GLOBAL_MAP=66, 1562 MATOP_SET_VALUES_LOCAL=67, 1563 MATOP_ZERO_ROWS_LOCAL=68, 1564 MATOP_GET_ROW_MAX_ABS=69, 1565 MATOP_GET_ROW_MIN_ABS=70, 1566 MATOP_CONVERT=71, 1567 MATOP_SET_COLORING=72, 1568 /* MATOP_PLACEHOLDER_73=73, */ 1569 MATOP_SET_VALUES_ADIFOR=74, 1570 MATOP_FD_COLORING_APPLY=75, 1571 MATOP_SET_FROM_OPTIONS=76, 1572 MATOP_MULT_CONSTRAINED=77, 1573 MATOP_MULT_TRANSPOSE_CONSTRAIN=78, 1574 MATOP_FIND_ZERO_DIAGONALS=79, 1575 MATOP_MULT_MULTIPLE=80, 1576 MATOP_SOLVE_MULTIPLE=81, 1577 MATOP_GET_INERTIA=82, 1578 MATOP_LOAD=83, 1579 MATOP_IS_SYMMETRIC=84, 1580 MATOP_IS_HERMITIAN=85, 1581 MATOP_IS_STRUCTURALLY_SYMMETRIC=86, 1582 MATOP_SET_VALUES_BLOCKEDLOCAL=87, 1583 MATOP_CREATE_VECS=88, 1584 MATOP_MAT_MULT=89, 1585 MATOP_MAT_MULT_SYMBOLIC=90, 1586 MATOP_MAT_MULT_NUMERIC=91, 1587 MATOP_PTAP=92, 1588 MATOP_PTAP_SYMBOLIC=93, 1589 MATOP_PTAP_NUMERIC=94, 1590 MATOP_MAT_TRANSPOSE_MULT=95, 1591 MATOP_MAT_TRANSPOSE_MULT_SYMBO=96, 1592 MATOP_MAT_TRANSPOSE_MULT_NUMER=97, 1593 /* MATOP_PLACEHOLDER_98=98, */ 1594 MATOP_PRODUCTSETFROMOPTIONS=99, 1595 MATOP_PRODUCTSYMBOLIC=100, 1596 MATOP_PRODUCTNUMERIC=101, 1597 MATOP_CONJUGATE=102, 1598 /* MATOP_PLACEHOLDER_103=103, */ 1599 MATOP_SET_VALUES_ROW=104, 1600 MATOP_REAL_PART=105, 1601 MATOP_IMAGINARY_PART=106, 1602 MATOP_GET_ROW_UPPER_TRIANGULAR=107, 1603 MATOP_RESTORE_ROW_UPPER_TRIANG=108, 1604 MATOP_MAT_SOLVE=109, 1605 MATOP_MAT_SOLVE_TRANSPOSE=110, 1606 MATOP_GET_ROW_MIN=111, 1607 MATOP_GET_COLUMN_VECTOR=112, 1608 MATOP_MISSING_DIAGONAL=113, 1609 MATOP_GET_SEQ_NONZERO_STRUCTUR=114, 1610 MATOP_CREATE=115, 1611 MATOP_GET_GHOSTS=116, 1612 MATOP_GET_LOCAL_SUB_MATRIX=117, 1613 MATOP_RESTORE_LOCALSUB_MATRIX=118, 1614 MATOP_MULT_DIAGONAL_BLOCK=119, 1615 MATOP_HERMITIAN_TRANSPOSE=120, 1616 MATOP_MULT_HERMITIAN_TRANSPOSE=121, 1617 MATOP_MULT_HERMITIAN_TRANS_ADD=122, 1618 MATOP_GET_MULTI_PROC_BLOCK=123, 1619 MATOP_FIND_NONZERO_ROWS=124, 1620 MATOP_GET_COLUMN_NORMS=125, 1621 MATOP_INVERT_BLOCK_DIAGONAL=126, 1622 /* MATOP_PLACEHOLDER_127=127, */ 1623 MATOP_CREATE_SUB_MATRICES_MPI=128, 1624 MATOP_SET_VALUES_BATCH=129, 1625 MATOP_TRANSPOSE_MAT_MULT=130, 1626 MATOP_TRANSPOSE_MAT_MULT_SYMBO=131, 1627 MATOP_TRANSPOSE_MAT_MULT_NUMER=132, 1628 MATOP_TRANSPOSE_COLORING_CREAT=133, 1629 MATOP_TRANS_COLORING_APPLY_SPT=134, 1630 MATOP_TRANS_COLORING_APPLY_DEN=135, 1631 MATOP_RART=136, 1632 MATOP_RART_SYMBOLIC=137, 1633 MATOP_RART_NUMERIC=138, 1634 MATOP_SET_BLOCK_SIZES=139, 1635 MATOP_AYPX=140, 1636 MATOP_RESIDUAL=141, 1637 MATOP_FDCOLORING_SETUP=142, 1638 MATOP_MPICONCATENATESEQ=144, 1639 MATOP_DESTROYSUBMATRICES=145, 1640 MATOP_TRANSPOSE_SOLVE=146, 1641 MATOP_GET_VALUES_LOCAL=147 1642 } MatOperation; 1643 PETSC_EXTERN PetscErrorCode MatSetOperation(Mat,MatOperation,void(*)(void)); 1644 PETSC_EXTERN PetscErrorCode MatGetOperation(Mat,MatOperation,void(**)(void)); 1645 PETSC_EXTERN PetscErrorCode MatHasOperation(Mat,MatOperation,PetscBool*); 1646 PETSC_EXTERN PetscErrorCode MatHasCongruentLayouts(Mat,PetscBool*); 1647 PETSC_DEPRECATED_FUNCTION("Use MatProductClear() (since version 3.14)") PETSC_STATIC_INLINE PetscErrorCode MatFreeIntermediateDataStructures(Mat A) {return MatProductClear(A);} 1648 PETSC_EXTERN PetscErrorCode MatShellSetOperation(Mat,MatOperation,void(*)(void)); 1649 PETSC_EXTERN PetscErrorCode MatShellGetOperation(Mat,MatOperation,void(**)(void)); 1650 PETSC_EXTERN PetscErrorCode MatShellSetContext(Mat,void*); 1651 PETSC_EXTERN PetscErrorCode MatShellSetVecType(Mat,VecType); 1652 PETSC_EXTERN PetscErrorCode MatShellTestMult(Mat,PetscErrorCode (*)(void*,Vec,Vec),Vec,void*,PetscBool*); 1653 PETSC_EXTERN PetscErrorCode MatShellTestMultTranspose(Mat,PetscErrorCode (*)(void*,Vec,Vec),Vec,void*,PetscBool*); 1654 PETSC_EXTERN PetscErrorCode MatShellSetManageScalingShifts(Mat); 1655 PETSC_EXTERN PetscErrorCode MatShellSetMatProductOperation(Mat,MatProductType,PetscErrorCode (*)(Mat,Mat,Mat,void**),PetscErrorCode (*)(Mat,Mat,Mat,void*),PetscErrorCode (*)(void*),MatType,MatType); 1656 PETSC_EXTERN PetscErrorCode MatIsShell(Mat,PetscBool*); 1657 1658 /* 1659 Codes for matrices stored on disk. By default they are 1660 stored in a universal format. By changing the format with 1661 PetscViewerPushFormat(viewer,PETSC_VIEWER_NATIVE); the matrices will 1662 be stored in a way natural for the matrix, for example dense matrices 1663 would be stored as dense. Matrices stored this way may only be 1664 read into matrices of the same type. 1665 */ 1666 #define MATRIX_BINARY_FORMAT_DENSE -1 1667 1668 PETSC_EXTERN PetscErrorCode MatMPIBAIJSetHashTableFactor(Mat,PetscReal); 1669 1670 PETSC_EXTERN PetscErrorCode MatISSetLocalMatType(Mat,MatType); 1671 PETSC_EXTERN PetscErrorCode MatISSetPreallocation(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 1672 PETSC_EXTERN PetscErrorCode MatISStoreL2L(Mat,PetscBool); 1673 PETSC_EXTERN PetscErrorCode MatISFixLocalEmpty(Mat,PetscBool); 1674 PETSC_EXTERN PetscErrorCode MatISGetLocalMat(Mat,Mat*); 1675 PETSC_EXTERN PetscErrorCode MatISRestoreLocalMat(Mat,Mat*); 1676 PETSC_EXTERN PetscErrorCode MatISSetLocalMat(Mat,Mat); 1677 PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use the MatConvert() interface (since version 3.10)") PetscErrorCode MatISGetMPIXAIJ(Mat,MatReuse,Mat*); 1678 1679 /*S 1680 MatNullSpace - Object that removes a null space from a vector, i.e. 1681 orthogonalizes the vector to a subsapce 1682 1683 Level: advanced 1684 1685 Users manual sections: 1686 . sec_singular 1687 1688 .seealso: MatNullSpaceCreate() 1689 S*/ 1690 typedef struct _p_MatNullSpace* MatNullSpace; 1691 1692 PETSC_EXTERN PetscErrorCode MatNullSpaceCreate(MPI_Comm,PetscBool ,PetscInt,const Vec[],MatNullSpace*); 1693 PETSC_EXTERN PetscErrorCode MatNullSpaceSetFunction(MatNullSpace,PetscErrorCode (*)(MatNullSpace,Vec,void*),void*); 1694 PETSC_EXTERN PetscErrorCode MatNullSpaceDestroy(MatNullSpace*); 1695 PETSC_EXTERN PetscErrorCode MatNullSpaceRemove(MatNullSpace,Vec); 1696 PETSC_EXTERN PetscErrorCode MatGetNullSpace(Mat, MatNullSpace *); 1697 PETSC_EXTERN PetscErrorCode MatGetTransposeNullSpace(Mat, MatNullSpace *); 1698 PETSC_EXTERN PetscErrorCode MatSetTransposeNullSpace(Mat,MatNullSpace); 1699 PETSC_EXTERN PetscErrorCode MatSetNullSpace(Mat,MatNullSpace); 1700 PETSC_EXTERN PetscErrorCode MatSetNearNullSpace(Mat,MatNullSpace); 1701 PETSC_EXTERN PetscErrorCode MatGetNearNullSpace(Mat,MatNullSpace*); 1702 PETSC_EXTERN PetscErrorCode MatNullSpaceTest(MatNullSpace,Mat,PetscBool *); 1703 PETSC_EXTERN PetscErrorCode MatNullSpaceView(MatNullSpace,PetscViewer); 1704 PETSC_EXTERN PetscErrorCode MatNullSpaceGetVecs(MatNullSpace,PetscBool*,PetscInt*,const Vec**); 1705 PETSC_EXTERN PetscErrorCode MatNullSpaceCreateRigidBody(Vec,MatNullSpace*); 1706 1707 PETSC_EXTERN PetscErrorCode MatReorderingSeqSBAIJ(Mat,IS); 1708 PETSC_EXTERN PetscErrorCode MatMPISBAIJSetHashTableFactor(Mat,PetscReal); 1709 PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetColumnIndices(Mat,PetscInt *); 1710 PETSC_EXTERN PetscErrorCode MatSeqBAIJInvertBlockDiagonal(Mat); 1711 1712 PETSC_EXTERN PetscErrorCode MatCreateMAIJ(Mat,PetscInt,Mat*); 1713 PETSC_EXTERN PetscErrorCode MatMAIJRedimension(Mat,PetscInt,Mat*); 1714 PETSC_EXTERN PetscErrorCode MatMAIJGetAIJ(Mat,Mat*); 1715 1716 PETSC_EXTERN PetscErrorCode MatComputeOperator(Mat,MatType,Mat*); 1717 PETSC_EXTERN PetscErrorCode MatComputeOperatorTranspose(Mat,MatType,Mat*); 1718 1719 PETSC_DEPRECATED_FUNCTION("Use MatComputeOperator() (since version 3.12)") PETSC_STATIC_INLINE PetscErrorCode MatComputeExplicitOperator(Mat A,Mat* B) { return MatComputeOperator(A,NULL,B); } 1720 PETSC_DEPRECATED_FUNCTION("Use MatComputeOperatorTranspose() (since version 3.12)") PETSC_STATIC_INLINE PetscErrorCode MatComputeExplicitOperatorTranspose(Mat A,Mat* B) { return MatComputeOperatorTranspose(A,NULL,B); } 1721 1722 PETSC_EXTERN PetscErrorCode MatCreateKAIJ(Mat,PetscInt,PetscInt,const PetscScalar[],const PetscScalar[],Mat*); 1723 PETSC_EXTERN PetscErrorCode MatKAIJGetAIJ(Mat,Mat*); 1724 PETSC_EXTERN PetscErrorCode MatKAIJGetS(Mat,PetscInt*,PetscInt*,PetscScalar**); 1725 PETSC_EXTERN PetscErrorCode MatKAIJGetSRead(Mat,PetscInt*,PetscInt*,const PetscScalar**); 1726 PETSC_EXTERN PetscErrorCode MatKAIJRestoreS(Mat,PetscScalar**); 1727 PETSC_EXTERN PetscErrorCode MatKAIJRestoreSRead(Mat,const PetscScalar**); 1728 PETSC_EXTERN PetscErrorCode MatKAIJGetT(Mat,PetscInt*,PetscInt*,PetscScalar**); 1729 PETSC_EXTERN PetscErrorCode MatKAIJGetTRead(Mat,PetscInt*,PetscInt*,const PetscScalar**); 1730 PETSC_EXTERN PetscErrorCode MatKAIJRestoreT(Mat,PetscScalar**); 1731 PETSC_EXTERN PetscErrorCode MatKAIJRestoreTRead(Mat,const PetscScalar**); 1732 PETSC_EXTERN PetscErrorCode MatKAIJSetAIJ(Mat,Mat); 1733 PETSC_EXTERN PetscErrorCode MatKAIJSetS(Mat,PetscInt,PetscInt,const PetscScalar[]); 1734 PETSC_EXTERN PetscErrorCode MatKAIJSetT(Mat,PetscInt,PetscInt,const PetscScalar[]); 1735 PETSC_EXTERN PetscErrorCode MatKAIJGetScaledIdentity(Mat,PetscBool*); 1736 1737 PETSC_EXTERN PetscErrorCode MatDiagonalScaleLocal(Mat,Vec); 1738 1739 PETSC_EXTERN PetscErrorCode MatCreateMFFD(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,Mat*); 1740 PETSC_EXTERN PetscErrorCode MatMFFDSetBase(Mat,Vec,Vec); 1741 PETSC_EXTERN PetscErrorCode MatMFFDSetFunction(Mat,PetscErrorCode(*)(void*,Vec,Vec),void*); 1742 PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioni(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*)); 1743 PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioniBase(Mat,PetscErrorCode (*)(void*,Vec)); 1744 PETSC_EXTERN PetscErrorCode MatMFFDSetHHistory(Mat,PetscScalar[],PetscInt); 1745 PETSC_EXTERN PetscErrorCode MatMFFDResetHHistory(Mat); 1746 PETSC_EXTERN PetscErrorCode MatMFFDSetFunctionError(Mat,PetscReal); 1747 PETSC_EXTERN PetscErrorCode MatMFFDSetPeriod(Mat,PetscInt); 1748 PETSC_EXTERN PetscErrorCode MatMFFDGetH(Mat,PetscScalar *); 1749 PETSC_EXTERN PetscErrorCode MatMFFDSetOptionsPrefix(Mat,const char[]); 1750 PETSC_EXTERN PetscErrorCode MatMFFDCheckPositivity(void*,Vec,Vec,PetscScalar*); 1751 PETSC_EXTERN PetscErrorCode MatMFFDSetCheckh(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*); 1752 1753 /*S 1754 MatMFFD - A data structured used to manage the computation of the h differencing parameter for matrix-free 1755 Jacobian vector products 1756 1757 Notes: 1758 MATMFFD is a specific MatType which uses the MatMFFD data structure 1759 1760 MatMFFD*() methods actually take the Mat as their first argument. Not a MatMFFD data structure 1761 1762 Level: developer 1763 1764 .seealso: MATMFFD, MatCreateMFFD(), MatMFFDSetFuction(), MatMFFDSetType(), MatMFFDRegister() 1765 S*/ 1766 typedef struct _p_MatMFFD* MatMFFD; 1767 1768 /*J 1769 MatMFFDType - algorithm used to compute the h used in computing matrix-vector products via differencing of the function 1770 1771 Level: beginner 1772 1773 .seealso: MatMFFDSetType(), MatMFFDRegister() 1774 J*/ 1775 typedef const char* MatMFFDType; 1776 #define MATMFFD_DS "ds" 1777 #define MATMFFD_WP "wp" 1778 1779 PETSC_EXTERN PetscErrorCode MatMFFDSetType(Mat,MatMFFDType); 1780 PETSC_EXTERN PetscErrorCode MatMFFDRegister(const char[],PetscErrorCode (*)(MatMFFD)); 1781 1782 PETSC_EXTERN PetscErrorCode MatMFFDDSSetUmin(Mat,PetscReal); 1783 PETSC_EXTERN PetscErrorCode MatMFFDWPSetComputeNormU(Mat,PetscBool ); 1784 1785 PETSC_EXTERN PetscErrorCode MatFDColoringSetType(MatFDColoring,MatMFFDType); 1786 1787 PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *); 1788 PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *); 1789 1790 #ifdef PETSC_HAVE_HARA 1791 PETSC_EXTERN_TYPEDEF typedef PetscScalar (*MatHaraKernel)(PetscInt,PetscReal[],PetscReal[],void*); 1792 PETSC_EXTERN PetscErrorCode MatCreateHaraFromKernel(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscReal[],MatHaraKernel,void*,PetscReal,PetscInt,PetscInt,Mat*); 1793 PETSC_EXTERN PetscErrorCode MatCreateHaraFromMat(Mat,PetscInt,const PetscReal[],PetscReal,PetscInt,PetscInt,PetscInt,PetscReal,Mat*); 1794 PETSC_EXTERN PetscErrorCode MatHaraSetSamplingMat(Mat,Mat,PetscInt,PetscReal); 1795 #endif 1796 1797 /* 1798 PETSc interface to MUMPS 1799 */ 1800 #ifdef PETSC_HAVE_MUMPS 1801 PETSC_EXTERN PetscErrorCode MatMumpsSetIcntl(Mat,PetscInt,PetscInt); 1802 PETSC_EXTERN PetscErrorCode MatMumpsGetIcntl(Mat,PetscInt,PetscInt*); 1803 PETSC_EXTERN PetscErrorCode MatMumpsSetCntl(Mat,PetscInt,PetscReal); 1804 PETSC_EXTERN PetscErrorCode MatMumpsGetCntl(Mat,PetscInt,PetscReal*); 1805 1806 PETSC_EXTERN PetscErrorCode MatMumpsGetInfo(Mat,PetscInt,PetscInt*); 1807 PETSC_EXTERN PetscErrorCode MatMumpsGetInfog(Mat,PetscInt,PetscInt*); 1808 PETSC_EXTERN PetscErrorCode MatMumpsGetRinfo(Mat,PetscInt,PetscReal*); 1809 PETSC_EXTERN PetscErrorCode MatMumpsGetRinfog(Mat,PetscInt,PetscReal*); 1810 PETSC_EXTERN PetscErrorCode MatMumpsGetInverse(Mat,Mat); 1811 PETSC_EXTERN PetscErrorCode MatMumpsGetInverseTranspose(Mat,Mat); 1812 #endif 1813 1814 /* 1815 PETSc interface to Mkl_Pardiso 1816 */ 1817 #ifdef PETSC_HAVE_MKL_PARDISO 1818 PETSC_EXTERN PetscErrorCode MatMkl_PardisoSetCntl(Mat,PetscInt,PetscInt); 1819 #endif 1820 1821 /* 1822 PETSc interface to Mkl_CPardiso 1823 */ 1824 #ifdef PETSC_HAVE_MKL_CPARDISO 1825 PETSC_EXTERN PetscErrorCode MatMkl_CPardisoSetCntl(Mat,PetscInt,PetscInt); 1826 #endif 1827 1828 /* 1829 PETSc interface to SUPERLU 1830 */ 1831 #ifdef PETSC_HAVE_SUPERLU 1832 PETSC_EXTERN PetscErrorCode MatSuperluSetILUDropTol(Mat,PetscReal); 1833 #endif 1834 1835 /* 1836 PETSc interface to SUPERLU_DIST 1837 */ 1838 #ifdef PETSC_HAVE_SUPERLU_DIST 1839 PETSC_EXTERN PetscErrorCode MatSuperluDistGetDiagU(Mat,PetscScalar*); 1840 #endif 1841 1842 /* 1843 PETSc interface to STRUMPACK 1844 */ 1845 #ifdef PETSC_HAVE_STRUMPACK 1846 /*E 1847 MatSTRUMPACKReordering - sparsity reducing ordering to be used in STRUMPACK 1848 1849 Level: intermediate 1850 E*/ 1851 typedef enum {MAT_STRUMPACK_NATURAL=0, 1852 MAT_STRUMPACK_METIS=1, 1853 MAT_STRUMPACK_PARMETIS=2, 1854 MAT_STRUMPACK_SCOTCH=3, 1855 MAT_STRUMPACK_PTSCOTCH=4, 1856 MAT_STRUMPACK_RCM=5} MatSTRUMPACKReordering; 1857 PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetReordering(Mat,MatSTRUMPACKReordering); 1858 PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetColPerm(Mat,PetscBool); 1859 PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSRelTol(Mat,PetscReal); 1860 PETSC_DEPRECATED_FUNCTION("Use MatSTRUMPACKSetHSSRelTol() (since version 3.9)") PETSC_STATIC_INLINE PetscErrorCode MatSTRUMPACKSetHSSRelCompTol(Mat mat,PetscReal rtol) {return MatSTRUMPACKSetHSSRelTol(mat,rtol);} 1861 PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSAbsTol(Mat,PetscReal); 1862 PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSMinSepSize(Mat,PetscInt); 1863 PETSC_DEPRECATED_FUNCTION("Use MatSTRUMPACKSetHSSMinSepSize() (since version 3.9)") PETSC_STATIC_INLINE PetscErrorCode MatSTRUMPACKSetHSSMinSize(Mat mat,PetscInt hssminsize) {return MatSTRUMPACKSetHSSMinSepSize(mat,hssminsize);} 1864 PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSMaxRank(Mat,PetscInt); 1865 PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSLeafSize(Mat,PetscInt); 1866 #endif 1867 1868 1869 PETSC_EXTERN PetscErrorCode MatBindToCPU(Mat,PetscBool); 1870 PETSC_DEPRECATED_FUNCTION("Use MatBindToCPU (since v3.13)") PETSC_STATIC_INLINE PetscErrorCode MatPinToCPU(Mat A,PetscBool flg) {return MatBindToCPU(A,flg);} 1871 1872 #ifdef PETSC_HAVE_CUDA 1873 /*E 1874 MatCUSPARSEStorageFormat - indicates the storage format for CUSPARSE (GPU) 1875 matrices. 1876 1877 Not Collective 1878 1879 + MAT_CUSPARSE_CSR - Compressed Sparse Row 1880 . MAT_CUSPARSE_ELL - Ellpack (requires CUDA 4.2 or later). 1881 - MAT_CUSPARSE_HYB - Hybrid, a combination of Ellpack and Coordinate format (requires CUDA 4.2 or later). 1882 1883 Level: intermediate 1884 1885 Any additions/changes here MUST also be made in include/petsc/finclude/petscmat.h 1886 1887 .seealso: MatCUSPARSESetFormat(), MatCUSPARSEFormatOperation 1888 E*/ 1889 1890 typedef enum {MAT_CUSPARSE_CSR, MAT_CUSPARSE_ELL, MAT_CUSPARSE_HYB} MatCUSPARSEStorageFormat; 1891 1892 /* these will be strings associated with enumerated type defined above */ 1893 PETSC_EXTERN const char *const MatCUSPARSEStorageFormats[]; 1894 1895 /*E 1896 MatCUSPARSEFormatOperation - indicates the operation of CUSPARSE (GPU) 1897 matrices whose operation should use a particular storage format. 1898 1899 Not Collective 1900 1901 + MAT_CUSPARSE_MULT_DIAG - sets the storage format for the diagonal matrix in the parallel MatMult 1902 . MAT_CUSPARSE_MULT_OFFDIAG - sets the storage format for the offdiagonal matrix in the parallel MatMult 1903 . MAT_CUSPARSE_MULT - sets the storage format for the entire matrix in the serial (single GPU) MatMult 1904 - MAT_CUSPARSE_ALL - sets the storage format for all CUSPARSE (GPU) matrices 1905 1906 Level: intermediate 1907 1908 .seealso: MatCUSPARSESetFormat(), MatCUSPARSEStorageFormat 1909 E*/ 1910 typedef enum {MAT_CUSPARSE_MULT_DIAG, MAT_CUSPARSE_MULT_OFFDIAG, MAT_CUSPARSE_MULT, MAT_CUSPARSE_ALL} MatCUSPARSEFormatOperation; 1911 1912 PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 1913 PETSC_EXTERN PetscErrorCode MatCreateAIJCUSPARSE(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 1914 PETSC_EXTERN PetscErrorCode MatCUSPARSESetFormat(Mat,MatCUSPARSEFormatOperation,MatCUSPARSEStorageFormat); 1915 1916 PETSC_EXTERN PetscErrorCode MatCreateDenseCUDA(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscScalar[],Mat*); 1917 PETSC_EXTERN PetscErrorCode MatCreateSeqDenseCUDA(MPI_Comm,PetscInt,PetscInt,PetscScalar[],Mat*); 1918 PETSC_EXTERN PetscErrorCode MatMPIDenseCUDASetPreallocation(Mat,PetscScalar[]); 1919 PETSC_EXTERN PetscErrorCode MatSeqDenseCUDASetPreallocation(Mat,PetscScalar[]); 1920 PETSC_EXTERN PetscErrorCode MatDenseCUDAGetArrayWrite(Mat,PetscScalar**); 1921 PETSC_EXTERN PetscErrorCode MatDenseCUDAGetArrayRead(Mat,const PetscScalar**); 1922 PETSC_EXTERN PetscErrorCode MatDenseCUDAGetArray(Mat,PetscScalar**); 1923 PETSC_EXTERN PetscErrorCode MatDenseCUDARestoreArrayWrite(Mat,PetscScalar**); 1924 PETSC_EXTERN PetscErrorCode MatDenseCUDARestoreArrayRead(Mat,const PetscScalar**); 1925 PETSC_EXTERN PetscErrorCode MatDenseCUDARestoreArray(Mat,PetscScalar**); 1926 PETSC_EXTERN PetscErrorCode MatDenseCUDAPlaceArray(Mat,const PetscScalar*); 1927 PETSC_EXTERN PetscErrorCode MatDenseCUDAReplaceArray(Mat,const PetscScalar*); 1928 PETSC_EXTERN PetscErrorCode MatDenseCUDAResetArray(Mat); 1929 #endif 1930 1931 #if defined(PETSC_HAVE_VIENNACL) 1932 PETSC_EXTERN PetscErrorCode MatCreateSeqAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,const PetscInt[],Mat*); 1933 PETSC_EXTERN PetscErrorCode MatCreateAIJViennaCL(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,const PetscInt[],PetscInt,const PetscInt[],Mat*); 1934 #endif 1935 1936 /* 1937 PETSc interface to FFTW 1938 */ 1939 #if defined(PETSC_HAVE_FFTW) 1940 PETSC_EXTERN PetscErrorCode VecScatterPetscToFFTW(Mat,Vec,Vec); 1941 PETSC_EXTERN PetscErrorCode VecScatterFFTWToPetsc(Mat,Vec,Vec); 1942 PETSC_EXTERN PetscErrorCode MatCreateVecsFFTW(Mat,Vec*,Vec*,Vec*); 1943 #endif 1944 1945 #if defined(PETSC_HAVE_SCALAPACK) 1946 PETSC_EXTERN PetscErrorCode MatCreateScaLAPACK(MPI_Comm,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,Mat*); 1947 PETSC_EXTERN PetscErrorCode MatScaLAPACKSetBlockSizes(Mat,PetscInt,PetscInt); 1948 PETSC_EXTERN PetscErrorCode MatScaLAPACKGetBlockSizes(Mat,PetscInt*,PetscInt*); 1949 #endif 1950 1951 PETSC_EXTERN PetscErrorCode MatCreateNest(MPI_Comm,PetscInt,const IS[],PetscInt,const IS[],const Mat[],Mat*); 1952 PETSC_EXTERN PetscErrorCode MatNestGetSize(Mat,PetscInt*,PetscInt*); 1953 PETSC_EXTERN PetscErrorCode MatNestGetISs(Mat,IS[],IS[]); 1954 PETSC_EXTERN PetscErrorCode MatNestGetLocalISs(Mat,IS[],IS[]); 1955 PETSC_EXTERN PetscErrorCode MatNestGetSubMats(Mat,PetscInt*,PetscInt*,Mat***); 1956 PETSC_EXTERN PetscErrorCode MatNestGetSubMat(Mat,PetscInt,PetscInt,Mat*); 1957 PETSC_EXTERN PetscErrorCode MatNestSetVecType(Mat,VecType); 1958 PETSC_EXTERN PetscErrorCode MatNestSetSubMats(Mat,PetscInt,const IS[],PetscInt,const IS[],const Mat[]); 1959 PETSC_EXTERN PetscErrorCode MatNestSetSubMat(Mat,PetscInt,PetscInt,Mat); 1960 1961 PETSC_EXTERN PetscErrorCode MatChop(Mat,PetscReal); 1962 PETSC_EXTERN PetscErrorCode MatComputeBandwidth(Mat,PetscReal,PetscInt*); 1963 1964 PETSC_EXTERN PetscErrorCode MatSubdomainsCreateCoalesce(Mat,PetscInt,PetscInt*,IS**); 1965 1966 PETSC_EXTERN PetscErrorCode MatPreallocatorPreallocate(Mat,PetscBool,Mat); 1967 1968 PETSC_INTERN PetscErrorCode MatHeaderMerge(Mat,Mat*); 1969 PETSC_EXTERN PetscErrorCode MatHeaderReplace(Mat,Mat*); 1970 1971 #endif 1972