1 /* 2 Include file for the matrix component of PETSc 3 */ 4 #ifndef PETSCMAT_H 5 #define PETSCMAT_H 6 7 #include <petscvec.h> 8 9 /* SUBMANSEC = Mat */ 10 11 /*S 12 Mat - Abstract PETSc matrix object used to manage all linear operators in PETSc, even those without 13 an explicit sparse representation (such as matrix-free operators) 14 15 Level: beginner 16 17 Note: 18 See [](doc_matrix) and `MatType` for available matrix types 19 20 .seealso: [](chapter_matrices), [](doc_matrix), `MatCreate()`, `MatType`, `MatSetType()`, `MatDestroy()` 21 S*/ 22 typedef struct _p_Mat *Mat; 23 24 /*J 25 MatType - String with the name of a PETSc matrix type 26 27 Level: beginner 28 29 Note: 30 [](doc_matrix) for a table of available matrix types 31 32 .seealso: [](chapter_matrices), `MatSetType()`, `Mat`, `MatSolverType`, `MatRegister()` 33 J*/ 34 typedef const char *MatType; 35 #define MATSAME "same" 36 #define MATMAIJ "maij" 37 #define MATSEQMAIJ "seqmaij" 38 #define MATMPIMAIJ "mpimaij" 39 #define MATKAIJ "kaij" 40 #define MATSEQKAIJ "seqkaij" 41 #define MATMPIKAIJ "mpikaij" 42 #define MATIS "is" 43 #define MATAIJ "aij" 44 #define MATSEQAIJ "seqaij" 45 #define MATMPIAIJ "mpiaij" 46 #define MATAIJCRL "aijcrl" 47 #define MATSEQAIJCRL "seqaijcrl" 48 #define MATMPIAIJCRL "mpiaijcrl" 49 #define MATAIJCUSPARSE "aijcusparse" 50 #define MATSEQAIJCUSPARSE "seqaijcusparse" 51 #define MATMPIAIJCUSPARSE "mpiaijcusparse" 52 #define MATAIJHIPSPARSE "aijhipsparse" 53 #define MATSEQAIJHIPSPARSE "seqaijhipsparse" 54 #define MATMPIAIJHIPSPARSE "mpiaijhipsparse" 55 #define MATAIJKOKKOS "aijkokkos" 56 #define MATSEQAIJKOKKOS "seqaijkokkos" 57 #define MATMPIAIJKOKKOS "mpiaijkokkos" 58 #define MATAIJVIENNACL "aijviennacl" 59 #define MATSEQAIJVIENNACL "seqaijviennacl" 60 #define MATMPIAIJVIENNACL "mpiaijviennacl" 61 #define MATAIJPERM "aijperm" 62 #define MATSEQAIJPERM "seqaijperm" 63 #define MATMPIAIJPERM "mpiaijperm" 64 #define MATAIJSELL "aijsell" 65 #define MATSEQAIJSELL "seqaijsell" 66 #define MATMPIAIJSELL "mpiaijsell" 67 #define MATAIJMKL "aijmkl" 68 #define MATSEQAIJMKL "seqaijmkl" 69 #define MATMPIAIJMKL "mpiaijmkl" 70 #define MATBAIJMKL "baijmkl" 71 #define MATSEQBAIJMKL "seqbaijmkl" 72 #define MATMPIBAIJMKL "mpibaijmkl" 73 #define MATSHELL "shell" 74 #define MATCENTERING "centering" 75 #define MATDENSE "dense" 76 #define MATDENSECUDA "densecuda" 77 #define MATDENSEHIP "densehip" 78 #define MATSEQDENSE "seqdense" 79 #define MATSEQDENSECUDA "seqdensecuda" 80 #define MATSEQDENSEHIP "seqdensehip" 81 #define MATMPIDENSE "mpidense" 82 #define MATMPIDENSECUDA "mpidensecuda" 83 #define MATMPIDENSEHIP "mpidensehip" 84 #define MATELEMENTAL "elemental" 85 #define MATSCALAPACK "scalapack" 86 #define MATBAIJ "baij" 87 #define MATSEQBAIJ "seqbaij" 88 #define MATMPIBAIJ "mpibaij" 89 #define MATMPIADJ "mpiadj" 90 #define MATSBAIJ "sbaij" 91 #define MATSEQSBAIJ "seqsbaij" 92 #define MATMPISBAIJ "mpisbaij" 93 #define MATMFFD "mffd" 94 #define MATNORMAL "normal" 95 #define MATNORMALHERMITIAN "normalh" 96 #define MATLRC "lrc" 97 #define MATSCATTER "scatter" 98 #define MATBLOCKMAT "blockmat" 99 #define MATCOMPOSITE "composite" 100 #define MATFFT "fft" 101 #define MATFFTW "fftw" 102 #define MATSEQCUFFT "seqcufft" 103 #define MATSEQHIPFFT "seqhipfft" 104 #define MATTRANSPOSEMAT PETSC_DEPRECATED_MACRO("GCC warning \"MATTRANSPOSEMAT macro is deprecated use MATTRANSPOSEVIRTUAL (since version 3.18)\"") "transpose" 105 #define MATTRANSPOSEVIRTUAL "transpose" 106 #define MATHERMITIANTRANSPOSEVIRTUAL "hermitiantranspose" 107 #define MATSCHURCOMPLEMENT "schurcomplement" 108 #define MATPYTHON "python" 109 #define MATHYPRE "hypre" 110 #define MATHYPRESTRUCT "hyprestruct" 111 #define MATHYPRESSTRUCT "hypresstruct" 112 #define MATSUBMATRIX "submatrix" 113 #define MATLOCALREF "localref" 114 #define MATNEST "nest" 115 #define MATPREALLOCATOR "preallocator" 116 #define MATSELL "sell" 117 #define MATSEQSELL "seqsell" 118 #define MATMPISELL "mpisell" 119 #define MATDUMMY "dummy" 120 #define MATLMVM "lmvm" 121 #define MATLMVMDFP "lmvmdfp" 122 #define MATLMVMBFGS "lmvmbfgs" 123 #define MATLMVMSR1 "lmvmsr1" 124 #define MATLMVMBROYDEN "lmvmbroyden" 125 #define MATLMVMBADBROYDEN "lmvmbadbroyden" 126 #define MATLMVMSYMBROYDEN "lmvmsymbroyden" 127 #define MATLMVMSYMBADBROYDEN "lmvmsymbadbroyden" 128 #define MATLMVMDIAGBROYDEN "lmvmdiagbroyden" 129 #define MATCONSTANTDIAGONAL "constantdiagonal" 130 #define MATHTOOL "htool" 131 #define MATH2OPUS "h2opus" 132 133 /*J 134 MatSolverType - String with the name of a PETSc matrix solver type. 135 136 For example: "petsc" indicates what PETSc provides, "superlu_dist" the parallel SuperLU_DIST package etc 137 138 Level: beginner 139 140 Note: 141 `MATSOLVERUMFPACK`, `MATSOLVERCHOLMOD`, `MATSOLVERKLU`, `MATSOLVERSPQR` form the SuiteSparse package for which you can use --download-suitesparse 142 143 .seealso: [](sec_matfactor), [](chapter_matrices), `MatGetFactor()`, `PCFactorSetMatSolverType()`, `PCFactorGetMatSolverType()` 144 J*/ 145 typedef const char *MatSolverType; 146 #define MATSOLVERSUPERLU "superlu" 147 #define MATSOLVERSUPERLU_DIST "superlu_dist" 148 #define MATSOLVERSTRUMPACK "strumpack" 149 #define MATSOLVERUMFPACK "umfpack" 150 #define MATSOLVERCHOLMOD "cholmod" 151 #define MATSOLVERKLU "klu" 152 #define MATSOLVERELEMENTAL "elemental" 153 #define MATSOLVERSCALAPACK "scalapack" 154 #define MATSOLVERESSL "essl" 155 #define MATSOLVERLUSOL "lusol" 156 #define MATSOLVERMUMPS "mumps" 157 #define MATSOLVERMKL_PARDISO "mkl_pardiso" 158 #define MATSOLVERMKL_CPARDISO "mkl_cpardiso" 159 #define MATSOLVERPASTIX "pastix" 160 #define MATSOLVERMATLAB "matlab" 161 #define MATSOLVERPETSC "petsc" 162 #define MATSOLVERBAS "bas" 163 #define MATSOLVERCUSPARSE "cusparse" 164 #define MATSOLVERCUSPARSEBAND "cusparseband" 165 #define MATSOLVERCUDA "cuda" 166 #define MATSOLVERHIPSPARSE "hipsparse" 167 #define MATSOLVERHIPSPARSEBAND "hipsparseband" 168 #define MATSOLVERHIP "hip" 169 #define MATSOLVERKOKKOS "kokkos" 170 #define MATSOLVERKOKKOSDEVICE "kokkosdevice" 171 #define MATSOLVERSPQR "spqr" 172 173 /*E 174 MatFactorType - indicates what type of factorization is requested 175 176 Level: beginner 177 178 Values: 179 + `MAT_FACTOR_LU` - LU factorization 180 . `MAT_FACTOR_CHOLESKY` - Cholesky factorization 181 . `MAT_FACTOR_ILU` - ILU factorization 182 . `MAT_FACTOR_ICC` - incomplete Cholesky factorization 183 . `MAT_FACTOR_ILUDT` - ILU factorization with drop tolerance 184 - `MAT_FACTOR_QR` - QR factorization 185 186 .seealso: [](chapter_matrices), `MatSolverType`, `MatGetFactor()`, `MatGetFactorAvailable()`, `MatSolverTypeRegister()` 187 E*/ 188 typedef enum { 189 MAT_FACTOR_NONE, 190 MAT_FACTOR_LU, 191 MAT_FACTOR_CHOLESKY, 192 MAT_FACTOR_ILU, 193 MAT_FACTOR_ICC, 194 MAT_FACTOR_ILUDT, 195 MAT_FACTOR_QR, 196 MAT_FACTOR_NUM_TYPES 197 } MatFactorType; 198 PETSC_EXTERN const char *const MatFactorTypes[]; 199 200 PETSC_EXTERN PetscErrorCode MatGetFactor(Mat, MatSolverType, MatFactorType, Mat *); 201 PETSC_EXTERN PetscErrorCode MatGetFactorAvailable(Mat, MatSolverType, MatFactorType, PetscBool *); 202 PETSC_EXTERN PetscErrorCode MatFactorGetCanUseOrdering(Mat, PetscBool *); 203 PETSC_DEPRECATED_FUNCTION("Use MatFactorGetCanUseOrdering() (since version 3.15)") static inline PetscErrorCode MatFactorGetUseOrdering(Mat A, PetscBool *b) 204 { 205 return MatFactorGetCanUseOrdering(A, b); 206 } 207 PETSC_EXTERN PetscErrorCode MatFactorGetSolverType(Mat, MatSolverType *); 208 PETSC_EXTERN PetscErrorCode MatGetFactorType(Mat, MatFactorType *); 209 PETSC_EXTERN PetscErrorCode MatSetFactorType(Mat, MatFactorType); 210 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*MatSolverFunction)(Mat, MatFactorType, Mat *); 211 PETSC_EXTERN PetscErrorCode MatSolverTypeRegister(MatSolverType, MatType, MatFactorType, MatSolverFunction); 212 PETSC_EXTERN PetscErrorCode MatSolverTypeGet(MatSolverType, MatType, MatFactorType, PetscBool *, PetscBool *, MatSolverFunction *); 213 typedef MatSolverType MatSolverPackage PETSC_DEPRECATED_TYPEDEF("Use MatSolverType (since version 3.9)"); 214 PETSC_DEPRECATED_FUNCTION("Use MatSolverTypeRegister() (since version 3.9)") static inline PetscErrorCode MatSolverPackageRegister(MatSolverType stype, MatType mtype, MatFactorType ftype, MatSolverFunction f) 215 { 216 return MatSolverTypeRegister(stype, mtype, ftype, f); 217 } 218 PETSC_DEPRECATED_FUNCTION("Use MatSolverTypeGet() (since version 3.9)") static inline PetscErrorCode MatSolverPackageGet(MatSolverType stype, MatType mtype, MatFactorType ftype, PetscBool *foundmtype, PetscBool *foundstype, MatSolverFunction *f) 219 { 220 return MatSolverTypeGet(stype, mtype, ftype, foundmtype, foundstype, f); 221 } 222 223 /*E 224 MatProductType - indicates what type of matrix product is requested 225 226 Level: beginner 227 228 Values: 229 + `MATPRODUCT_AB` - product of two matrices 230 . `MATPRODUCT_AtB` - product of the transpose of a given matrix with a matrix 231 . `MATPRODUCT_ABt` - product of a matrix with the transpose of another given matrix 232 . `MATPRODUCT_PtAP` - the triple product of the transpose of a matrix with another matrix and itself 233 . `MATPRODUCT_RARt` - the triple product of a matrix, another matrix and the transpose of the first matrix 234 - `MATPRODUCT_ABC` - the product of three matrices 235 236 .seealso: [](sec_matmatproduct), [](chapter_matrices), `MatProductSetType()` 237 E*/ 238 typedef enum { 239 MATPRODUCT_UNSPECIFIED = 0, 240 MATPRODUCT_AB, 241 MATPRODUCT_AtB, 242 MATPRODUCT_ABt, 243 MATPRODUCT_PtAP, 244 MATPRODUCT_RARt, 245 MATPRODUCT_ABC 246 } MatProductType; 247 PETSC_EXTERN const char *const MatProductTypes[]; 248 249 /*J 250 MatProductAlgorithm - String with the name of an algorithm for a PETSc matrix product implementation 251 252 Level: beginner 253 254 .seealso: [](sec_matmatproduct), [](chapter_matrices), `MatSetType()`, `Mat`, `MatProductSetAlgorithm()`, `MatProductType` 255 J*/ 256 typedef const char *MatProductAlgorithm; 257 #define MATPRODUCTALGORITHMDEFAULT "default" 258 #define MATPRODUCTALGORITHMSORTED "sorted" 259 #define MATPRODUCTALGORITHMSCALABLE "scalable" 260 #define MATPRODUCTALGORITHMSCALABLEFAST "scalable_fast" 261 #define MATPRODUCTALGORITHMHEAP "heap" 262 #define MATPRODUCTALGORITHMBHEAP "btheap" 263 #define MATPRODUCTALGORITHMLLCONDENSED "llcondensed" 264 #define MATPRODUCTALGORITHMROWMERGE "rowmerge" 265 #define MATPRODUCTALGORITHMOUTERPRODUCT "outerproduct" 266 #define MATPRODUCTALGORITHMATB "at*b" 267 #define MATPRODUCTALGORITHMRAP "rap" 268 #define MATPRODUCTALGORITHMNONSCALABLE "nonscalable" 269 #define MATPRODUCTALGORITHMSEQMPI "seqmpi" 270 #define MATPRODUCTALGORITHMBACKEND "backend" 271 #define MATPRODUCTALGORITHMOVERLAPPING "overlapping" 272 #define MATPRODUCTALGORITHMMERGED "merged" 273 #define MATPRODUCTALGORITHMALLATONCE "allatonce" 274 #define MATPRODUCTALGORITHMALLATONCEMERGED "allatonce_merged" 275 #define MATPRODUCTALGORITHMALLGATHERV "allgatherv" 276 #define MATPRODUCTALGORITHMCYCLIC "cyclic" 277 #if defined(PETSC_HAVE_HYPRE) 278 #define MATPRODUCTALGORITHMHYPRE "hypre" 279 #endif 280 281 PETSC_EXTERN PetscErrorCode MatProductCreate(Mat, Mat, Mat, Mat *); 282 PETSC_EXTERN PetscErrorCode MatProductCreateWithMat(Mat, Mat, Mat, Mat); 283 PETSC_EXTERN PetscErrorCode MatProductSetType(Mat, MatProductType); 284 PETSC_EXTERN PetscErrorCode MatProductSetAlgorithm(Mat, MatProductAlgorithm); 285 PETSC_EXTERN PetscErrorCode MatProductSetFill(Mat, PetscReal); 286 PETSC_EXTERN PetscErrorCode MatProductSetFromOptions(Mat); 287 PETSC_EXTERN PetscErrorCode MatProductSymbolic(Mat); 288 PETSC_EXTERN PetscErrorCode MatProductNumeric(Mat); 289 PETSC_EXTERN PetscErrorCode MatProductReplaceMats(Mat, Mat, Mat, Mat); 290 PETSC_EXTERN PetscErrorCode MatProductClear(Mat); 291 PETSC_EXTERN PetscErrorCode MatProductView(Mat, PetscViewer); 292 PETSC_EXTERN PetscErrorCode MatProductGetType(Mat, MatProductType *); 293 PETSC_EXTERN PetscErrorCode MatProductGetMats(Mat, Mat *, Mat *, Mat *); 294 295 /* Logging support */ 296 #define MAT_FILE_CLASSID 1211216 /* used to indicate matrices in binary files */ 297 PETSC_EXTERN PetscClassId MAT_CLASSID; 298 PETSC_EXTERN PetscClassId MAT_COLORING_CLASSID; 299 PETSC_EXTERN PetscClassId MAT_FDCOLORING_CLASSID; 300 PETSC_EXTERN PetscClassId MAT_TRANSPOSECOLORING_CLASSID; 301 PETSC_EXTERN PetscClassId MAT_PARTITIONING_CLASSID; 302 PETSC_EXTERN PetscClassId MAT_COARSEN_CLASSID; 303 PETSC_EXTERN PetscClassId MAT_NULLSPACE_CLASSID; 304 PETSC_EXTERN PetscClassId MATMFFD_CLASSID; 305 306 /*E 307 MatReuse - Indicates if matrices obtained from a previous call to `MatCreateSubMatrices()`, `MatCreateSubMatrix()`, `MatConvert()` or several other functions 308 are to be reused to store the new matrix values. 309 310 Values: 311 + `MAT_INITIAL_MATRIX` - create a new matrix 312 . `MAT_REUSE_MATRIX` - reuse the matrix created with a previous call that used `MAT_INITIAL_MATRIX` 313 . `MAT_INPLACE_MATRIX` - replace the first input matrix with the new matrix (not applicable to all functions) 314 - `MAT_IGNORE_MATRIX` - do not create a new matrix or reuse a give matrix, just ignore that matrix argument (not applicable to all functions) 315 316 Level: beginner 317 318 .seealso: [](chapter_matrices), `Mat`, `MatCreateSubMatrices()`, `MatCreateSubMatrix()`, `MatDestroyMatrices()`, `MatConvert()` 319 E*/ 320 typedef enum { 321 MAT_INITIAL_MATRIX, 322 MAT_REUSE_MATRIX, 323 MAT_IGNORE_MATRIX, 324 MAT_INPLACE_MATRIX 325 } MatReuse; 326 327 /*E 328 MatCreateSubMatrixOption - Indicates if matrices obtained from a call to `MatCreateSubMatrices()` 329 include the matrix values. Currently it is only used by `MatGetSeqNonzeroStructure()`. 330 331 Level: developer 332 333 Values: 334 + `MAT_DO_NOT_GET_VALUES` - do not copy the matrix values 335 - `MAT_GET_VALUES` - copy the matrix values 336 337 Developer Note: 338 Why is not just a boolean used for this information? 339 340 .seealso: [](chapter_matrices), `Mat`, `MatDuplicateOption`, `PetscCopyMode`, `MatGetSeqNonzeroStructure()` 341 E*/ 342 typedef enum { 343 MAT_DO_NOT_GET_VALUES, 344 MAT_GET_VALUES 345 } MatCreateSubMatrixOption; 346 347 PETSC_EXTERN PetscErrorCode MatInitializePackage(void); 348 349 PETSC_EXTERN PetscErrorCode MatCreate(MPI_Comm, Mat *); 350 PETSC_EXTERN PetscErrorCode MatSetSizes(Mat, PetscInt, PetscInt, PetscInt, PetscInt); 351 PETSC_EXTERN PetscErrorCode MatSetType(Mat, MatType); 352 PETSC_EXTERN PetscErrorCode MatGetVecType(Mat, VecType *); 353 PETSC_EXTERN PetscErrorCode MatSetVecType(Mat, VecType); 354 PETSC_EXTERN PetscErrorCode MatSetFromOptions(Mat); 355 PETSC_EXTERN PetscErrorCode MatViewFromOptions(Mat, PetscObject, const char[]); 356 PETSC_EXTERN PetscErrorCode MatRegister(const char[], PetscErrorCode (*)(Mat)); 357 PETSC_EXTERN PetscErrorCode MatRegisterRootName(const char[], const char[], const char[]); 358 PETSC_EXTERN PetscErrorCode MatSetOptionsPrefix(Mat, const char[]); 359 PETSC_EXTERN PetscErrorCode MatSetOptionsPrefixFactor(Mat, const char[]); 360 PETSC_EXTERN PetscErrorCode MatAppendOptionsPrefixFactor(Mat, const char[]); 361 PETSC_EXTERN PetscErrorCode MatAppendOptionsPrefix(Mat, const char[]); 362 PETSC_EXTERN PetscErrorCode MatGetOptionsPrefix(Mat, const char *[]); 363 PETSC_EXTERN PetscErrorCode MatSetErrorIfFailure(Mat, PetscBool); 364 365 PETSC_EXTERN PetscFunctionList MatList; 366 PETSC_EXTERN PetscFunctionList MatColoringList; 367 PETSC_EXTERN PetscFunctionList MatPartitioningList; 368 369 /*E 370 MatStructure - Indicates if two matrices have the same nonzero structure 371 372 Level: beginner 373 374 Values: 375 + `SAME_NONZERO_PATTERN` - the two matrices have identical nonzero patterns 376 . `DIFFERENT_NONZERO_PATTERN` - the two matrices may have different nonzero patterns 377 . `SUBSET_NONZERO_PATTERN` - the nonzero pattern of the second matrix is a subset of the nonzero pattern of the first matrix 378 - `UNKNOWN_NONZERO_PATTERN` - there is no known relationship between the nonzero patterns. In this case the implementations 379 may try to detect a relationship to optimize the operation 380 381 .seealso: [](chapter_matrices), `Mat`, `MatCopy()`, `MatAXPY()`, `MatAYPX()` 382 E*/ 383 typedef enum { 384 DIFFERENT_NONZERO_PATTERN, 385 SUBSET_NONZERO_PATTERN, 386 SAME_NONZERO_PATTERN, 387 UNKNOWN_NONZERO_PATTERN 388 } MatStructure; 389 PETSC_EXTERN const char *const MatStructures[]; 390 391 #if defined PETSC_HAVE_MKL_SPARSE 392 PETSC_EXTERN PetscErrorCode MatCreateSeqAIJMKL(MPI_Comm, PetscInt, PetscInt, PetscInt, const PetscInt[], Mat *); 393 PETSC_EXTERN PetscErrorCode MatCreateMPIAIJMKL(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, const PetscInt[], PetscInt, const PetscInt[], Mat *); 394 PETSC_EXTERN PetscErrorCode MatCreateBAIJMKL(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, const PetscInt[], PetscInt, const PetscInt[], Mat *); 395 PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJMKL(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, const PetscInt[], Mat *); 396 #endif 397 398 PETSC_EXTERN PetscErrorCode MatCreateSeqSELL(MPI_Comm, PetscInt, PetscInt, PetscInt, const PetscInt[], Mat *); 399 PETSC_EXTERN PetscErrorCode MatCreateSELL(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, const PetscInt[], PetscInt, const PetscInt[], Mat *); 400 PETSC_EXTERN PetscErrorCode MatSeqSELLSetPreallocation(Mat, PetscInt, const PetscInt[]); 401 PETSC_EXTERN PetscErrorCode MatMPISELLSetPreallocation(Mat, PetscInt, const PetscInt[], PetscInt, const PetscInt[]); 402 403 PETSC_EXTERN PetscErrorCode MatCreateSeqDense(MPI_Comm, PetscInt, PetscInt, PetscScalar[], Mat *); 404 PETSC_EXTERN PetscErrorCode MatCreateDense(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar[], Mat *); 405 PETSC_EXTERN PetscErrorCode MatCreateSeqAIJ(MPI_Comm, PetscInt, PetscInt, PetscInt, const PetscInt[], Mat *); 406 PETSC_EXTERN PetscErrorCode MatCreateAIJ(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, const PetscInt[], PetscInt, const PetscInt[], Mat *); 407 PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithArrays(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], Mat *); 408 PETSC_EXTERN PetscErrorCode MatUpdateMPIAIJWithArrays(Mat, PetscInt, PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[]); 409 PETSC_EXTERN PetscErrorCode MatUpdateMPIAIJWithArray(Mat, const PetscScalar[]); 410 PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithSplitArrays(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt[], PetscInt[], PetscScalar[], PetscInt[], PetscInt[], PetscScalar[], Mat *); 411 PETSC_EXTERN PetscErrorCode MatCreateMPIAIJWithSeqAIJ(MPI_Comm, Mat, Mat, const PetscInt[], Mat *); 412 413 PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJ(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, const PetscInt[], Mat *); 414 PETSC_EXTERN PetscErrorCode MatCreateBAIJ(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, const PetscInt[], PetscInt, const PetscInt[], Mat *); 415 PETSC_EXTERN PetscErrorCode MatCreateMPIBAIJWithArrays(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], Mat *); 416 417 PETSC_EXTERN PetscErrorCode MatSetPreallocationCOO(Mat, PetscCount, PetscInt[], PetscInt[]); 418 PETSC_EXTERN PetscErrorCode MatSetPreallocationCOOLocal(Mat, PetscCount, PetscInt[], PetscInt[]); 419 PETSC_EXTERN PetscErrorCode MatSetValuesCOO(Mat, const PetscScalar[], InsertMode); 420 421 PETSC_EXTERN PetscErrorCode MatCreateMPIAdj(MPI_Comm, PetscInt, PetscInt, PetscInt[], PetscInt[], PetscInt[], Mat *); 422 PETSC_EXTERN PetscErrorCode MatCreateSeqSBAIJ(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, const PetscInt[], Mat *); 423 424 PETSC_EXTERN PetscErrorCode MatCreateSBAIJ(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, const PetscInt[], PetscInt, const PetscInt[], Mat *); 425 PETSC_EXTERN PetscErrorCode MatCreateMPISBAIJWithArrays(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], Mat *); 426 PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocationCSR(Mat, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[]); 427 PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocationCSR(Mat, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[]); 428 PETSC_EXTERN PetscErrorCode MatXAIJSetPreallocation(Mat, PetscInt, const PetscInt[], const PetscInt[], const PetscInt[], const PetscInt[]); 429 430 PETSC_EXTERN PetscErrorCode MatCreateShell(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, void *, Mat *); 431 PETSC_EXTERN PetscErrorCode MatCreateCentering(MPI_Comm, PetscInt, PetscInt, Mat *); 432 PETSC_EXTERN PetscErrorCode MatCreateNormal(Mat, Mat *); 433 PETSC_EXTERN PetscErrorCode MatCreateNormalHermitian(Mat, Mat *); 434 PETSC_EXTERN PetscErrorCode MatCreateLRC(Mat, Mat, Vec, Mat, Mat *); 435 PETSC_EXTERN PetscErrorCode MatLRCGetMats(Mat, Mat *, Mat *, Vec *, Mat *); 436 PETSC_EXTERN PetscErrorCode MatCreateIS(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, ISLocalToGlobalMapping, ISLocalToGlobalMapping, Mat *); 437 PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCRL(MPI_Comm, PetscInt, PetscInt, PetscInt, const PetscInt[], Mat *); 438 PETSC_EXTERN PetscErrorCode MatCreateMPIAIJCRL(MPI_Comm, PetscInt, PetscInt, PetscInt, const PetscInt[], PetscInt, const PetscInt[], Mat *); 439 440 PETSC_EXTERN PetscErrorCode MatCreateScatter(MPI_Comm, VecScatter, Mat *); 441 PETSC_EXTERN PetscErrorCode MatScatterSetVecScatter(Mat, VecScatter); 442 PETSC_EXTERN PetscErrorCode MatScatterGetVecScatter(Mat, VecScatter *); 443 PETSC_EXTERN PetscErrorCode MatCreateBlockMat(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt *, Mat *); 444 PETSC_EXTERN PetscErrorCode MatCompositeAddMat(Mat, Mat); 445 PETSC_EXTERN PetscErrorCode MatCompositeMerge(Mat); 446 typedef enum { 447 MAT_COMPOSITE_MERGE_RIGHT, 448 MAT_COMPOSITE_MERGE_LEFT 449 } MatCompositeMergeType; 450 PETSC_EXTERN PetscErrorCode MatCompositeSetMergeType(Mat, MatCompositeMergeType); 451 PETSC_EXTERN PetscErrorCode MatCreateComposite(MPI_Comm, PetscInt, const Mat *, Mat *); 452 typedef enum { 453 MAT_COMPOSITE_ADDITIVE, 454 MAT_COMPOSITE_MULTIPLICATIVE 455 } MatCompositeType; 456 PETSC_EXTERN PetscErrorCode MatCompositeSetType(Mat, MatCompositeType); 457 PETSC_EXTERN PetscErrorCode MatCompositeGetType(Mat, MatCompositeType *); 458 PETSC_EXTERN PetscErrorCode MatCompositeSetMatStructure(Mat, MatStructure); 459 PETSC_EXTERN PetscErrorCode MatCompositeGetMatStructure(Mat, MatStructure *); 460 PETSC_EXTERN PetscErrorCode MatCompositeGetNumberMat(Mat, PetscInt *); 461 PETSC_EXTERN PetscErrorCode MatCompositeGetMat(Mat, PetscInt, Mat *); 462 PETSC_EXTERN PetscErrorCode MatCompositeSetScalings(Mat, const PetscScalar *); 463 464 PETSC_EXTERN PetscErrorCode MatCreateFFT(MPI_Comm, PetscInt, const PetscInt[], MatType, Mat *); 465 PETSC_EXTERN PetscErrorCode MatCreateSeqCUFFT(MPI_Comm, PetscInt, const PetscInt[], Mat *); 466 467 PETSC_EXTERN PetscErrorCode MatCreateTranspose(Mat, Mat *); 468 PETSC_EXTERN PetscErrorCode MatTransposeGetMat(Mat, Mat *); 469 PETSC_EXTERN PetscErrorCode MatCreateHermitianTranspose(Mat, Mat *); 470 PETSC_EXTERN PetscErrorCode MatHermitianTransposeGetMat(Mat, Mat *); 471 PETSC_EXTERN PetscErrorCode MatNormalGetMat(Mat, Mat *); 472 PETSC_EXTERN PetscErrorCode MatNormalHermitianGetMat(Mat, Mat *); 473 PETSC_EXTERN PetscErrorCode MatCreateSubMatrixVirtual(Mat, IS, IS, Mat *); 474 PETSC_EXTERN PetscErrorCode MatSubMatrixVirtualUpdate(Mat, Mat, IS, IS); 475 PETSC_EXTERN PetscErrorCode MatCreateLocalRef(Mat, IS, IS, Mat *); 476 PETSC_EXTERN PetscErrorCode MatCreateConstantDiagonal(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar, Mat *); 477 478 #if defined(PETSC_HAVE_HYPRE) 479 PETSC_EXTERN PetscErrorCode MatHYPRESetPreallocation(Mat, PetscInt, const PetscInt[], PetscInt, const PetscInt[]); 480 #endif 481 482 PETSC_EXTERN PetscErrorCode MatPythonSetType(Mat, const char[]); 483 PETSC_EXTERN PetscErrorCode MatPythonGetType(Mat, const char *[]); 484 485 PETSC_EXTERN PetscErrorCode MatResetPreallocation(Mat); 486 PETSC_EXTERN PetscErrorCode MatSetUp(Mat); 487 PETSC_EXTERN PetscErrorCode MatDestroy(Mat *); 488 PETSC_EXTERN PetscErrorCode MatGetNonzeroState(Mat, PetscObjectState *); 489 490 PETSC_EXTERN PetscErrorCode MatConjugate(Mat); 491 PETSC_EXTERN PetscErrorCode MatRealPart(Mat); 492 PETSC_EXTERN PetscErrorCode MatImaginaryPart(Mat); 493 PETSC_EXTERN PetscErrorCode MatGetDiagonalBlock(Mat, Mat *); 494 PETSC_EXTERN PetscErrorCode MatGetTrace(Mat, PetscScalar *); 495 PETSC_EXTERN PetscErrorCode MatInvertBlockDiagonal(Mat, const PetscScalar **); 496 PETSC_EXTERN PetscErrorCode MatInvertVariableBlockDiagonal(Mat, PetscInt, const PetscInt *, PetscScalar *); 497 PETSC_EXTERN PetscErrorCode MatInvertBlockDiagonalMat(Mat, Mat); 498 PETSC_EXTERN PetscErrorCode MatInvertVariableBlockEnvelope(Mat, MatReuse, Mat *); 499 500 /* ------------------------------------------------------------*/ 501 PETSC_EXTERN PetscErrorCode MatSetValues(Mat, PetscInt, const PetscInt[], PetscInt, const PetscInt[], const PetscScalar[], InsertMode); 502 PETSC_EXTERN PetscErrorCode MatSetValuesIS(Mat, IS, IS, const PetscScalar[], InsertMode); 503 PETSC_EXTERN PetscErrorCode MatSetValuesBlocked(Mat, PetscInt, const PetscInt[], PetscInt, const PetscInt[], const PetscScalar[], InsertMode); 504 PETSC_EXTERN PetscErrorCode MatSetValuesRow(Mat, PetscInt, const PetscScalar[]); 505 PETSC_EXTERN PetscErrorCode MatSetValuesRowLocal(Mat, PetscInt, const PetscScalar[]); 506 PETSC_EXTERN PetscErrorCode MatSetValuesBatch(Mat, PetscInt, PetscInt, PetscInt[], const PetscScalar[]); 507 PETSC_EXTERN PetscErrorCode MatSetRandom(Mat, PetscRandom); 508 509 /*S 510 MatStencil - Data structure (C struct) for storing information about rows and 511 columns of a matrix as indexed on an associated grid. These are arguments to `MatSetStencil()` and `MatSetBlockStencil()` 512 513 Level: beginner 514 515 Notes: 516 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). 517 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. 518 519 For stencil access to vectors see `DMDAVecGetArray()`, `DMDAVecGetArrayF90()`. 520 521 For staggered grids, see `DMStagStencil` 522 523 Fortran Note: 524 See `MatSetValuesStencil()` for details. 525 526 .seealso: [](chapter_matrices), `Mat`, `MatSetValuesStencil()`, `MatSetStencil()`, `MatSetValuesBlockedStencil()`, `DMDAVecGetArray()`, `DMDAVecGetArrayF90()`, 527 `DMStagStencil` 528 S*/ 529 typedef struct { 530 PetscInt k, j, i, c; 531 } MatStencil; 532 533 PETSC_EXTERN PetscErrorCode MatSetValuesStencil(Mat, PetscInt, const MatStencil[], PetscInt, const MatStencil[], const PetscScalar[], InsertMode); 534 PETSC_EXTERN PetscErrorCode MatSetValuesBlockedStencil(Mat, PetscInt, const MatStencil[], PetscInt, const MatStencil[], const PetscScalar[], InsertMode); 535 PETSC_EXTERN PetscErrorCode MatSetStencil(Mat, PetscInt, const PetscInt[], const PetscInt[], PetscInt); 536 537 /*E 538 MatAssemblyType - Indicates if the process of setting values into the matrix is complete, and the matrix is ready for use 539 540 Level: beginner 541 542 Values: 543 + `MAT_FLUSH_ASSEMBLY` - you will continue to put values into the matrix 544 - `MAT_FINAL_ASSEMBLY` - you wish to use the matrix with the values currently inserted 545 546 .seealso: [](chapter_matrices), `Mat`, `MatSetValues`, `MatAssemblyBegin()`, `MatAssemblyEnd()` 547 E*/ 548 typedef enum { 549 MAT_FLUSH_ASSEMBLY = 1, 550 MAT_FINAL_ASSEMBLY = 0 551 } MatAssemblyType; 552 PETSC_EXTERN PetscErrorCode MatAssemblyBegin(Mat, MatAssemblyType); 553 PETSC_EXTERN PetscErrorCode MatAssemblyEnd(Mat, MatAssemblyType); 554 PETSC_EXTERN PetscErrorCode MatAssembled(Mat, PetscBool *); 555 556 /*E 557 MatOption - Options that may be set for a matrix that indicate properties of the matrix or affect its behavior or storage 558 559 Level: beginner 560 561 Note: 562 See `MatSetOption()` for the use of the options 563 564 Developer Note: 565 Entries that are negative need not be called collectively by all processes. 566 567 .seealso: [](chapter_matrices), `Mat`, `MatSetOption()` 568 E*/ 569 typedef enum { 570 MAT_OPTION_MIN = -3, 571 MAT_UNUSED_NONZERO_LOCATION_ERR = -2, 572 MAT_ROW_ORIENTED = -1, 573 MAT_SYMMETRIC = 1, 574 MAT_STRUCTURALLY_SYMMETRIC = 2, 575 MAT_FORCE_DIAGONAL_ENTRIES = 3, 576 MAT_IGNORE_OFF_PROC_ENTRIES = 4, 577 MAT_USE_HASH_TABLE = 5, 578 MAT_KEEP_NONZERO_PATTERN = 6, 579 MAT_IGNORE_ZERO_ENTRIES = 7, 580 MAT_USE_INODES = 8, 581 MAT_HERMITIAN = 9, 582 MAT_SYMMETRY_ETERNAL = 10, 583 MAT_NEW_NONZERO_LOCATION_ERR = 11, 584 MAT_IGNORE_LOWER_TRIANGULAR = 12, 585 MAT_ERROR_LOWER_TRIANGULAR = 13, 586 MAT_GETROW_UPPERTRIANGULAR = 14, 587 MAT_SPD = 15, 588 MAT_NO_OFF_PROC_ZERO_ROWS = 16, 589 MAT_NO_OFF_PROC_ENTRIES = 17, 590 MAT_NEW_NONZERO_LOCATIONS = 18, 591 MAT_NEW_NONZERO_ALLOCATION_ERR = 19, 592 MAT_SUBSET_OFF_PROC_ENTRIES = 20, 593 MAT_SUBMAT_SINGLEIS = 21, 594 MAT_STRUCTURE_ONLY = 22, 595 MAT_SORTED_FULL = 23, 596 MAT_FORM_EXPLICIT_TRANSPOSE = 24, 597 MAT_STRUCTURAL_SYMMETRY_ETERNAL = 25, 598 MAT_SPD_ETERNAL = 26, 599 MAT_OPTION_MAX = 27 600 } MatOption; 601 602 PETSC_EXTERN const char *const *MatOptions; 603 PETSC_EXTERN PetscErrorCode MatSetOption(Mat, MatOption, PetscBool); 604 PETSC_EXTERN PetscErrorCode MatGetOption(Mat, MatOption, PetscBool *); 605 PETSC_EXTERN PetscErrorCode MatPropagateSymmetryOptions(Mat, Mat); 606 PETSC_EXTERN PetscErrorCode MatGetType(Mat, MatType *); 607 608 PETSC_EXTERN PetscErrorCode MatGetValues(Mat, PetscInt, const PetscInt[], PetscInt, const PetscInt[], PetscScalar[]); 609 PETSC_EXTERN PetscErrorCode MatGetRow(Mat, PetscInt, PetscInt *, const PetscInt *[], const PetscScalar *[]); 610 PETSC_EXTERN PetscErrorCode MatRestoreRow(Mat, PetscInt, PetscInt *, const PetscInt *[], const PetscScalar *[]); 611 PETSC_EXTERN PetscErrorCode MatGetRowUpperTriangular(Mat); 612 PETSC_EXTERN PetscErrorCode MatRestoreRowUpperTriangular(Mat); 613 PETSC_EXTERN PetscErrorCode MatGetColumnVector(Mat, Vec, PetscInt); 614 PETSC_EXTERN PetscErrorCode MatSeqAIJGetArray(Mat, PetscScalar *[]); 615 PETSC_EXTERN PetscErrorCode MatSeqAIJGetArrayRead(Mat, const PetscScalar *[]); 616 PETSC_EXTERN PetscErrorCode MatSeqAIJGetArrayWrite(Mat, PetscScalar *[]); 617 PETSC_EXTERN PetscErrorCode MatSeqAIJRestoreArray(Mat, PetscScalar *[]); 618 PETSC_EXTERN PetscErrorCode MatSeqAIJRestoreArrayRead(Mat, const PetscScalar *[]); 619 PETSC_EXTERN PetscErrorCode MatSeqAIJRestoreArrayWrite(Mat, PetscScalar *[]); 620 PETSC_EXTERN PetscErrorCode MatSeqAIJGetMaxRowNonzeros(Mat, PetscInt *); 621 PETSC_EXTERN PetscErrorCode MatSeqAIJSetValuesLocalFast(Mat, PetscInt, const PetscInt[], PetscInt, const PetscInt[], const PetscScalar[], InsertMode); 622 PETSC_EXTERN PetscErrorCode MatSeqAIJSetType(Mat, MatType); 623 PETSC_EXTERN PetscErrorCode MatSeqAIJKron(Mat, Mat, MatReuse, Mat *); 624 PETSC_EXTERN PetscErrorCode MatSeqAIJRegister(const char[], PetscErrorCode (*)(Mat, MatType, MatReuse, Mat *)); 625 PETSC_EXTERN PetscFunctionList MatSeqAIJList; 626 PETSC_EXTERN PetscErrorCode MatSeqBAIJGetArray(Mat, PetscScalar *[]); 627 PETSC_EXTERN PetscErrorCode MatSeqBAIJRestoreArray(Mat, PetscScalar *[]); 628 PETSC_EXTERN PetscErrorCode MatSeqSBAIJGetArray(Mat, PetscScalar *[]); 629 PETSC_EXTERN PetscErrorCode MatSeqSBAIJRestoreArray(Mat, PetscScalar *[]); 630 PETSC_EXTERN PetscErrorCode MatDenseGetArray(Mat, PetscScalar *[]); 631 PETSC_EXTERN PetscErrorCode MatDenseRestoreArray(Mat, PetscScalar *[]); 632 PETSC_EXTERN PetscErrorCode MatDensePlaceArray(Mat, const PetscScalar[]); 633 PETSC_EXTERN PetscErrorCode MatDenseReplaceArray(Mat, const PetscScalar[]); 634 PETSC_EXTERN PetscErrorCode MatDenseResetArray(Mat); 635 PETSC_EXTERN PetscErrorCode MatDenseGetArrayRead(Mat, const PetscScalar *[]); 636 PETSC_EXTERN PetscErrorCode MatDenseRestoreArrayRead(Mat, const PetscScalar *[]); 637 PETSC_EXTERN PetscErrorCode MatDenseGetArrayWrite(Mat, PetscScalar *[]); 638 PETSC_EXTERN PetscErrorCode MatDenseRestoreArrayWrite(Mat, PetscScalar *[]); 639 PETSC_EXTERN PetscErrorCode MatDenseGetArrayAndMemType(Mat, PetscScalar *[], PetscMemType *); 640 PETSC_EXTERN PetscErrorCode MatDenseRestoreArrayAndMemType(Mat, PetscScalar *[]); 641 PETSC_EXTERN PetscErrorCode MatDenseGetArrayReadAndMemType(Mat, const PetscScalar *[], PetscMemType *); 642 PETSC_EXTERN PetscErrorCode MatDenseRestoreArrayReadAndMemType(Mat, const PetscScalar *[]); 643 PETSC_EXTERN PetscErrorCode MatDenseGetArrayWriteAndMemType(Mat, PetscScalar *[], PetscMemType *); 644 PETSC_EXTERN PetscErrorCode MatDenseRestoreArrayWriteAndMemType(Mat, PetscScalar *[]); 645 PETSC_EXTERN PetscErrorCode MatGetBlockSize(Mat, PetscInt *); 646 PETSC_EXTERN PetscErrorCode MatSetBlockSize(Mat, PetscInt); 647 PETSC_EXTERN PetscErrorCode MatGetBlockSizes(Mat, PetscInt *, PetscInt *); 648 PETSC_EXTERN PetscErrorCode MatSetBlockSizes(Mat, PetscInt, PetscInt); 649 PETSC_EXTERN PetscErrorCode MatSetBlockSizesFromMats(Mat, Mat, Mat); 650 PETSC_EXTERN PetscErrorCode MatSetVariableBlockSizes(Mat, PetscInt, PetscInt *); 651 PETSC_EXTERN PetscErrorCode MatGetVariableBlockSizes(Mat, PetscInt *, const PetscInt **); 652 653 PETSC_EXTERN PetscErrorCode MatDenseGetColumn(Mat, PetscInt, PetscScalar *[]); 654 PETSC_EXTERN PetscErrorCode MatDenseRestoreColumn(Mat, PetscScalar *[]); 655 PETSC_EXTERN PetscErrorCode MatDenseGetColumnVec(Mat, PetscInt, Vec *); 656 PETSC_EXTERN PetscErrorCode MatDenseRestoreColumnVec(Mat, PetscInt, Vec *); 657 PETSC_EXTERN PetscErrorCode MatDenseGetColumnVecRead(Mat, PetscInt, Vec *); 658 PETSC_EXTERN PetscErrorCode MatDenseRestoreColumnVecRead(Mat, PetscInt, Vec *); 659 PETSC_EXTERN PetscErrorCode MatDenseGetColumnVecWrite(Mat, PetscInt, Vec *); 660 PETSC_EXTERN PetscErrorCode MatDenseRestoreColumnVecWrite(Mat, PetscInt, Vec *); 661 PETSC_EXTERN PetscErrorCode MatDenseGetSubMatrix(Mat, PetscInt, PetscInt, PetscInt, PetscInt, Mat *); 662 PETSC_EXTERN PetscErrorCode MatDenseRestoreSubMatrix(Mat, Mat *); 663 664 PETSC_EXTERN PetscErrorCode MatMult(Mat, Vec, Vec); 665 PETSC_EXTERN PetscErrorCode MatMultDiagonalBlock(Mat, Vec, Vec); 666 PETSC_EXTERN PetscErrorCode MatMultAdd(Mat, Vec, Vec, Vec); 667 PETSC_EXTERN PetscErrorCode MatMultTranspose(Mat, Vec, Vec); 668 PETSC_EXTERN PetscErrorCode MatMultHermitianTranspose(Mat, Vec, Vec); 669 PETSC_EXTERN PetscErrorCode MatIsTranspose(Mat, Mat, PetscReal, PetscBool *); 670 PETSC_EXTERN PetscErrorCode MatIsHermitianTranspose(Mat, Mat, PetscReal, PetscBool *); 671 PETSC_EXTERN PetscErrorCode MatMultTransposeAdd(Mat, Vec, Vec, Vec); 672 PETSC_EXTERN PetscErrorCode MatMultHermitianTransposeAdd(Mat, Vec, Vec, Vec); 673 PETSC_EXTERN PetscErrorCode MatMatSolve(Mat, Mat, Mat); 674 PETSC_EXTERN PetscErrorCode MatMatSolveTranspose(Mat, Mat, Mat); 675 PETSC_EXTERN PetscErrorCode MatMatTransposeSolve(Mat, Mat, Mat); 676 PETSC_EXTERN PetscErrorCode MatResidual(Mat, Vec, Vec, Vec); 677 678 /*E 679 MatDuplicateOption - Indicates if a duplicated sparse matrix should have 680 its numerical values copied over or just its nonzero structure. 681 682 Level: beginner 683 684 Values: 685 + `MAT_DO_NOT_COPY_VALUES` - Create a matrix using the same nonzero pattern as the original matrix, 686 with zeros for the numerical values 687 . `MAT_COPY_VALUES` - Create a matrix with the same nonzero pattern as the original matrix 688 and with the same numerical values. 689 - `MAT_SHARE_NONZERO_PATTERN` - Create a matrix that shares the nonzero structure with the previous matrix 690 and does not copy it, using zeros for the numerical values. The parent and 691 child matrices will share their index (i and j) arrays, and you cannot 692 insert new nonzero entries into either matrix 693 694 Note: 695 Many matrix types (including `MATSEQAIJ`) do not support the `MAT_SHARE_NONZERO_PATTERN` optimization; in 696 this case the behavior is as if `MAT_DO_NOT_COPY_VALUES` has been specified. 697 698 .seealso: [](chapter_matrices), `Mat`, `MatDuplicate()` 699 E*/ 700 typedef enum { 701 MAT_DO_NOT_COPY_VALUES, 702 MAT_COPY_VALUES, 703 MAT_SHARE_NONZERO_PATTERN 704 } MatDuplicateOption; 705 706 PETSC_EXTERN PetscErrorCode MatConvert(Mat, MatType, MatReuse, Mat *); 707 PETSC_EXTERN PetscErrorCode MatDuplicate(Mat, MatDuplicateOption, Mat *); 708 709 PETSC_EXTERN PetscErrorCode MatCopy(Mat, Mat, MatStructure); 710 PETSC_EXTERN PetscErrorCode MatView(Mat, PetscViewer); 711 PETSC_EXTERN PetscErrorCode MatIsSymmetric(Mat, PetscReal, PetscBool *); 712 PETSC_EXTERN PetscErrorCode MatIsStructurallySymmetric(Mat, PetscBool *); 713 PETSC_EXTERN PetscErrorCode MatIsHermitian(Mat, PetscReal, PetscBool *); 714 PETSC_EXTERN PetscErrorCode MatIsSymmetricKnown(Mat, PetscBool *, PetscBool *); 715 PETSC_EXTERN PetscErrorCode MatIsHermitianKnown(Mat, PetscBool *, PetscBool *); 716 PETSC_EXTERN PetscErrorCode MatIsStructurallySymmetricKnown(Mat, PetscBool *, PetscBool *); 717 PETSC_EXTERN PetscErrorCode MatIsSPDKnown(Mat, PetscBool *, PetscBool *); 718 PETSC_EXTERN PetscErrorCode MatMissingDiagonal(Mat, PetscBool *, PetscInt *); 719 PETSC_EXTERN PetscErrorCode MatLoad(Mat, PetscViewer); 720 721 PETSC_EXTERN PetscErrorCode MatGetRowIJ(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscBool *); 722 PETSC_EXTERN PetscErrorCode MatRestoreRowIJ(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscBool *); 723 PETSC_EXTERN PetscErrorCode MatGetColumnIJ(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscBool *); 724 PETSC_EXTERN PetscErrorCode MatRestoreColumnIJ(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscBool *); 725 726 /*S 727 MatInfo - Context of matrix information, used with `MatGetInfo()` 728 729 Level: intermediate 730 731 Fortran Note: 732 Information is stored as a double-precision array of dimension `MAT_INFO_SIZE` 733 734 .seealso: [](chapter_matrices), `Mat`, `MatGetInfo()`, `MatInfoType` 735 S*/ 736 typedef struct { 737 PetscLogDouble block_size; /* block size */ 738 PetscLogDouble nz_allocated, nz_used, nz_unneeded; /* number of nonzeros */ 739 PetscLogDouble memory; /* memory allocated */ 740 PetscLogDouble assemblies; /* number of matrix assemblies called */ 741 PetscLogDouble mallocs; /* number of mallocs during MatSetValues() */ 742 PetscLogDouble fill_ratio_given, fill_ratio_needed; /* fill ratio for LU/ILU */ 743 PetscLogDouble factor_mallocs; /* number of mallocs during factorization */ 744 } MatInfo; 745 746 /*E 747 MatInfoType - Indicates if you want information about the local part of the matrix, 748 the entire parallel matrix or the maximum over all the local parts. 749 750 Level: beginner 751 752 Values: 753 + `MAT_LOCAL` - values for each MPI process part of the matrix 754 . `MAT_GLOBAL_MAX` - maximum of each value over all MPI processes 755 - `MAT_GLOBAL_SUM` - sum of each value over all MPI processes 756 757 .seealso: [](chapter_matrices), `MatGetInfo()`, `MatInfo` 758 E*/ 759 typedef enum { 760 MAT_LOCAL = 1, 761 MAT_GLOBAL_MAX = 2, 762 MAT_GLOBAL_SUM = 3 763 } MatInfoType; 764 PETSC_EXTERN PetscErrorCode MatGetInfo(Mat, MatInfoType, MatInfo *); 765 PETSC_EXTERN PetscErrorCode MatGetDiagonal(Mat, Vec); 766 PETSC_EXTERN PetscErrorCode MatGetRowMax(Mat, Vec, PetscInt[]); 767 PETSC_EXTERN PetscErrorCode MatGetRowMin(Mat, Vec, PetscInt[]); 768 PETSC_EXTERN PetscErrorCode MatGetRowMaxAbs(Mat, Vec, PetscInt[]); 769 PETSC_EXTERN PetscErrorCode MatGetRowMinAbs(Mat, Vec, PetscInt[]); 770 PETSC_EXTERN PetscErrorCode MatGetRowSum(Mat, Vec); 771 PETSC_EXTERN PetscErrorCode MatTranspose(Mat, MatReuse, Mat *); 772 PETSC_EXTERN PetscErrorCode MatTransposeSymbolic(Mat, Mat *); 773 PETSC_EXTERN PetscErrorCode MatTransposeSetPrecursor(Mat, Mat); 774 PETSC_EXTERN PetscErrorCode MatHermitianTranspose(Mat, MatReuse, Mat *); 775 PETSC_EXTERN PetscErrorCode MatPermute(Mat, IS, IS, Mat *); 776 PETSC_EXTERN PetscErrorCode MatDiagonalScale(Mat, Vec, Vec); 777 PETSC_EXTERN PetscErrorCode MatDiagonalSet(Mat, Vec, InsertMode); 778 779 PETSC_EXTERN PetscErrorCode MatEqual(Mat, Mat, PetscBool *); 780 PETSC_EXTERN PetscErrorCode MatMultEqual(Mat, Mat, PetscInt, PetscBool *); 781 PETSC_EXTERN PetscErrorCode MatMultAddEqual(Mat, Mat, PetscInt, PetscBool *); 782 PETSC_EXTERN PetscErrorCode MatMultTransposeEqual(Mat, Mat, PetscInt, PetscBool *); 783 PETSC_EXTERN PetscErrorCode MatMultTransposeAddEqual(Mat, Mat, PetscInt, PetscBool *); 784 PETSC_EXTERN PetscErrorCode MatMultHermitianTransposeEqual(Mat, Mat, PetscInt, PetscBool *); 785 PETSC_EXTERN PetscErrorCode MatMultHermitianTransposeAddEqual(Mat, Mat, PetscInt, PetscBool *); 786 PETSC_EXTERN PetscErrorCode MatMatMultEqual(Mat, Mat, Mat, PetscInt, PetscBool *); 787 PETSC_EXTERN PetscErrorCode MatTransposeMatMultEqual(Mat, Mat, Mat, PetscInt, PetscBool *); 788 PETSC_EXTERN PetscErrorCode MatMatTransposeMultEqual(Mat, Mat, Mat, PetscInt, PetscBool *); 789 PETSC_EXTERN PetscErrorCode MatPtAPMultEqual(Mat, Mat, Mat, PetscInt, PetscBool *); 790 PETSC_EXTERN PetscErrorCode MatRARtMultEqual(Mat, Mat, Mat, PetscInt, PetscBool *); 791 PETSC_EXTERN PetscErrorCode MatIsLinear(Mat, PetscInt, PetscBool *); 792 793 PETSC_EXTERN PetscErrorCode MatNorm(Mat, NormType, PetscReal *); 794 PETSC_EXTERN PetscErrorCode MatGetColumnNorms(Mat, NormType, PetscReal *); 795 PETSC_EXTERN PetscErrorCode MatGetColumnSums(Mat, PetscScalar *); 796 PETSC_EXTERN PetscErrorCode MatGetColumnSumsRealPart(Mat, PetscReal *); 797 PETSC_EXTERN PetscErrorCode MatGetColumnSumsImaginaryPart(Mat, PetscReal *); 798 PETSC_EXTERN PetscErrorCode MatGetColumnMeans(Mat, PetscScalar *); 799 PETSC_EXTERN PetscErrorCode MatGetColumnMeansRealPart(Mat, PetscReal *); 800 PETSC_EXTERN PetscErrorCode MatGetColumnMeansImaginaryPart(Mat, PetscReal *); 801 PETSC_EXTERN PetscErrorCode MatGetColumnReductions(Mat, PetscInt, PetscReal *); 802 PETSC_EXTERN PetscErrorCode MatZeroEntries(Mat); 803 PETSC_EXTERN PetscErrorCode MatSetInf(Mat); 804 PETSC_EXTERN PetscErrorCode MatZeroRows(Mat, PetscInt, const PetscInt[], PetscScalar, Vec, Vec); 805 PETSC_EXTERN PetscErrorCode MatZeroRowsIS(Mat, IS, PetscScalar, Vec, Vec); 806 PETSC_EXTERN PetscErrorCode MatZeroRowsStencil(Mat, PetscInt, const MatStencil[], PetscScalar, Vec, Vec); 807 PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsStencil(Mat, PetscInt, const MatStencil[], PetscScalar, Vec, Vec); 808 PETSC_EXTERN PetscErrorCode MatZeroRowsColumns(Mat, PetscInt, const PetscInt[], PetscScalar, Vec, Vec); 809 PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsIS(Mat, IS, PetscScalar, Vec, Vec); 810 811 PETSC_EXTERN PetscErrorCode MatGetSize(Mat, PetscInt *, PetscInt *); 812 PETSC_EXTERN PetscErrorCode MatGetLocalSize(Mat, PetscInt *, PetscInt *); 813 PETSC_EXTERN PetscErrorCode MatGetOwnershipRange(Mat, PetscInt *, PetscInt *); 814 PETSC_EXTERN PetscErrorCode MatGetOwnershipRanges(Mat, const PetscInt **); 815 PETSC_EXTERN PetscErrorCode MatGetOwnershipRangeColumn(Mat, PetscInt *, PetscInt *); 816 PETSC_EXTERN PetscErrorCode MatGetOwnershipRangesColumn(Mat, const PetscInt **); 817 PETSC_EXTERN PetscErrorCode MatGetOwnershipIS(Mat, IS *, IS *); 818 819 PETSC_EXTERN PetscErrorCode MatCreateSubMatrices(Mat, PetscInt, const IS[], const IS[], MatReuse, Mat *[]); 820 PETSC_DEPRECATED_FUNCTION("Use MatCreateSubMatrices() (since version 3.8)") static inline PetscErrorCode MatGetSubMatrices(Mat mat, PetscInt n, const IS irow[], const IS icol[], MatReuse scall, Mat *submat[]) 821 { 822 return MatCreateSubMatrices(mat, n, irow, icol, scall, submat); 823 } 824 PETSC_EXTERN PetscErrorCode MatCreateSubMatricesMPI(Mat, PetscInt, const IS[], const IS[], MatReuse, Mat *[]); 825 PETSC_DEPRECATED_FUNCTION("Use MatCreateSubMatricesMPI() (since version 3.8)") static inline PetscErrorCode MatGetSubMatricesMPI(Mat mat, PetscInt n, const IS irow[], const IS icol[], MatReuse scall, Mat *submat[]) 826 { 827 return MatCreateSubMatricesMPI(mat, n, irow, icol, scall, submat); 828 } 829 PETSC_EXTERN PetscErrorCode MatDestroyMatrices(PetscInt, Mat *[]); 830 PETSC_EXTERN PetscErrorCode MatDestroySubMatrices(PetscInt, Mat *[]); 831 PETSC_EXTERN PetscErrorCode MatCreateSubMatrix(Mat, IS, IS, MatReuse, Mat *); 832 PETSC_DEPRECATED_FUNCTION("Use MatCreateSubMatrix() (since version 3.8)") static inline PetscErrorCode MatGetSubMatrix(Mat mat, IS isrow, IS iscol, MatReuse cll, Mat *newmat) 833 { 834 return MatCreateSubMatrix(mat, isrow, iscol, cll, newmat); 835 } 836 PETSC_EXTERN PetscErrorCode MatGetLocalSubMatrix(Mat, IS, IS, Mat *); 837 PETSC_EXTERN PetscErrorCode MatRestoreLocalSubMatrix(Mat, IS, IS, Mat *); 838 PETSC_EXTERN PetscErrorCode MatGetSeqNonzeroStructure(Mat, Mat *); 839 PETSC_EXTERN PetscErrorCode MatDestroySeqNonzeroStructure(Mat *); 840 841 PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJ(MPI_Comm, Mat, PetscInt, PetscInt, MatReuse, Mat *); 842 PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJSymbolic(MPI_Comm, Mat, PetscInt, PetscInt, Mat *); 843 PETSC_EXTERN PetscErrorCode MatCreateMPIAIJSumSeqAIJNumeric(Mat, Mat); 844 PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMat(Mat, MatReuse, Mat *); 845 PETSC_EXTERN PetscErrorCode MatAIJGetLocalMat(Mat, Mat *); 846 PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMatCondensed(Mat, MatReuse, IS *, IS *, Mat *); 847 PETSC_EXTERN PetscErrorCode MatMPIAIJGetLocalMatMerge(Mat, MatReuse, IS *, Mat *); 848 PETSC_EXTERN PetscErrorCode MatMPIAIJGetNumberNonzeros(Mat, PetscCount *); 849 PETSC_EXTERN PetscErrorCode MatGetBrowsOfAcols(Mat, Mat, MatReuse, IS *, IS *, Mat *); 850 PETSC_EXTERN PetscErrorCode MatGetGhosts(Mat, PetscInt *, const PetscInt *[]); 851 852 PETSC_EXTERN PetscErrorCode MatIncreaseOverlap(Mat, PetscInt, IS[], PetscInt); 853 PETSC_EXTERN PetscErrorCode MatIncreaseOverlapSplit(Mat mat, PetscInt n, IS is[], PetscInt ov); 854 PETSC_EXTERN PetscErrorCode MatMPIAIJSetUseScalableIncreaseOverlap(Mat, PetscBool); 855 856 PETSC_EXTERN PetscErrorCode MatMatMult(Mat, Mat, MatReuse, PetscReal, Mat *); 857 858 PETSC_EXTERN PetscErrorCode MatMatMatMult(Mat, Mat, Mat, MatReuse, PetscReal, Mat *); 859 PETSC_EXTERN PetscErrorCode MatGalerkin(Mat, Mat, Mat, MatReuse, PetscReal, Mat *); 860 861 PETSC_EXTERN PetscErrorCode MatPtAP(Mat, Mat, MatReuse, PetscReal, Mat *); 862 PETSC_EXTERN PetscErrorCode MatRARt(Mat, Mat, MatReuse, PetscReal, Mat *); 863 864 PETSC_EXTERN PetscErrorCode MatTransposeMatMult(Mat, Mat, MatReuse, PetscReal, Mat *); 865 PETSC_EXTERN PetscErrorCode MatMatTransposeMult(Mat, Mat, MatReuse, PetscReal, Mat *); 866 867 PETSC_EXTERN PetscErrorCode MatAXPY(Mat, PetscScalar, Mat, MatStructure); 868 PETSC_EXTERN PetscErrorCode MatAYPX(Mat, PetscScalar, Mat, MatStructure); 869 870 PETSC_EXTERN PetscErrorCode MatScale(Mat, PetscScalar); 871 PETSC_EXTERN PetscErrorCode MatShift(Mat, PetscScalar); 872 873 PETSC_EXTERN PetscErrorCode MatSetLocalToGlobalMapping(Mat, ISLocalToGlobalMapping, ISLocalToGlobalMapping); 874 PETSC_EXTERN PetscErrorCode MatGetLocalToGlobalMapping(Mat, ISLocalToGlobalMapping *, ISLocalToGlobalMapping *); 875 PETSC_EXTERN PetscErrorCode MatGetLayouts(Mat, PetscLayout *, PetscLayout *); 876 PETSC_EXTERN PetscErrorCode MatSetLayouts(Mat, PetscLayout, PetscLayout); 877 PETSC_EXTERN PetscErrorCode MatZeroRowsLocal(Mat, PetscInt, const PetscInt[], PetscScalar, Vec, Vec); 878 PETSC_EXTERN PetscErrorCode MatZeroRowsLocalIS(Mat, IS, PetscScalar, Vec, Vec); 879 PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocal(Mat, PetscInt, const PetscInt[], PetscScalar, Vec, Vec); 880 PETSC_EXTERN PetscErrorCode MatZeroRowsColumnsLocalIS(Mat, IS, PetscScalar, Vec, Vec); 881 PETSC_EXTERN PetscErrorCode MatGetValuesLocal(Mat, PetscInt, const PetscInt[], PetscInt, const PetscInt[], PetscScalar[]); 882 PETSC_EXTERN PetscErrorCode MatSetValuesLocal(Mat, PetscInt, const PetscInt[], PetscInt, const PetscInt[], const PetscScalar[], InsertMode); 883 PETSC_EXTERN PetscErrorCode MatSetValuesBlockedLocal(Mat, PetscInt, const PetscInt[], PetscInt, const PetscInt[], const PetscScalar[], InsertMode); 884 885 PETSC_EXTERN PetscErrorCode MatStashSetInitialSize(Mat, PetscInt, PetscInt); 886 PETSC_EXTERN PetscErrorCode MatStashGetInfo(Mat, PetscInt *, PetscInt *, PetscInt *, PetscInt *); 887 888 PETSC_EXTERN PetscErrorCode MatInterpolate(Mat, Vec, Vec); 889 PETSC_EXTERN PetscErrorCode MatInterpolateAdd(Mat, Vec, Vec, Vec); 890 PETSC_EXTERN PetscErrorCode MatRestrict(Mat, Vec, Vec); 891 PETSC_EXTERN PetscErrorCode MatMatInterpolate(Mat, Mat, Mat *); 892 PETSC_EXTERN PetscErrorCode MatMatInterpolateAdd(Mat, Mat, Mat, Mat *); 893 PETSC_EXTERN PetscErrorCode MatMatRestrict(Mat, Mat, Mat *); 894 PETSC_EXTERN PetscErrorCode MatCreateVecs(Mat, Vec *, Vec *); 895 PETSC_DEPRECATED_FUNCTION("Use MatCreateVecs() (since version 3.6)") static inline PetscErrorCode MatGetVecs(Mat mat, Vec *x, Vec *y) 896 { 897 return MatCreateVecs(mat, x, y); 898 } 899 PETSC_EXTERN PetscErrorCode MatCreateRedundantMatrix(Mat, PetscInt, MPI_Comm, MatReuse, Mat *); 900 PETSC_EXTERN PetscErrorCode MatGetMultiProcBlock(Mat, MPI_Comm, MatReuse, Mat *); 901 PETSC_EXTERN PetscErrorCode MatFindZeroDiagonals(Mat, IS *); 902 PETSC_EXTERN PetscErrorCode MatFindOffBlockDiagonalEntries(Mat, IS *); 903 PETSC_EXTERN PetscErrorCode MatCreateMPIMatConcatenateSeqMat(MPI_Comm, Mat, PetscInt, MatReuse, Mat *); 904 905 /*MC 906 MatSetValue - Set a single entry into a matrix. 907 908 Not Collective 909 910 Synopsis: 911 #include <petscmat.h> 912 PetscErrorCode MatSetValue(Mat m,PetscInt row,PetscInt col,PetscScalar value,InsertMode mode) 913 914 Input Parameters: 915 + m - the matrix 916 . row - the row location of the entry 917 . col - the column location of the entry 918 . value - the value to insert 919 - mode - either `INSERT_VALUES` or `ADD_VALUES` 920 921 Level: beginner 922 923 Note: 924 For efficiency one should use `MatSetValues()` and set several values simultaneously. 925 926 .seealso: [](chapter_matrices), `Mat`, `MatAssemblyBegin()`, `MatAssemblyEnd`, `InsertMode`, `MatGetValue()`, `MatSetValues()`, 927 `MatSetValueLocal()`, `MatSetValuesLocal()` 928 M*/ 929 static inline PetscErrorCode MatSetValue(Mat v, PetscInt i, PetscInt j, PetscScalar va, InsertMode mode) 930 { 931 return MatSetValues(v, 1, &i, 1, &j, &va, mode); 932 } 933 934 /*@C 935 MatGetValue - Gets a single value from a matrix 936 937 Not Collective; can only return a value owned by the given process 938 939 Input Parameters: 940 + mat - the matrix 941 . row - the row location of the entry 942 - col - the column location of the entry 943 944 Output Parameter: 945 . va - the value 946 947 Level: advanced 948 949 Notes: 950 The matrix must have been assembled with `MatAssemblyBegin()` and `MatAssemblyEnd` before this call 951 952 For efficiency one should use `MatGetValues()` and get several values simultaneously. 953 954 See notes for `MatGetValues()`. 955 956 .seealso: [](chapter_matrices), `Mat`, `MatAssemblyBegin()`, `MatAssemblyEnd`, `MatSetValue()`, `MatGetValueLocal()`, `MatGetValues()` 957 @*/ 958 static inline PetscErrorCode MatGetValue(Mat mat, PetscInt row, PetscInt col, PetscScalar *va) 959 { 960 return MatGetValues(mat, 1, &row, 1, &col, va); 961 } 962 963 /*MC 964 MatSetValueLocal - Inserts or adds a single value into a matrix, using a local numbering of the nodes. 965 966 Not Collective 967 968 Input Parameters: 969 + m - the matrix 970 . row - the row location of the entry 971 . col - the column location of the entry 972 . value - the value to insert 973 - mode - either `INSERT_VALUES` or `ADD_VALUES` 974 975 Level: intermediate 976 977 Notes: 978 For efficiency one should use `MatSetValuesLocal()` and set several values simultaneously. 979 980 See notes for `MatSetValuesLocal()` for additional information on when and how this function can be used. 981 982 .seealso: [](chapter_matrices), `MatSetValue()`, `MatSetValuesLocal()` 983 M*/ 984 static inline PetscErrorCode MatSetValueLocal(Mat v, PetscInt i, PetscInt j, PetscScalar va, InsertMode mode) 985 { 986 return MatSetValuesLocal(v, 1, &i, 1, &j, &va, mode); 987 } 988 989 /*MC 990 MatPreallocateBegin - Begins the block of code that will count the number of nonzeros per 991 row in a matrix providing the data that one can use to correctly preallocate the matrix. 992 993 Synopsis: 994 #include <petscmat.h> 995 PetscErrorCode MatPreallocateBegin(MPI_Comm comm, PetscInt nrows, PetscInt ncols, PetscInt *dnz, PetscInt *onz) 996 997 Collective 998 999 Input Parameters: 1000 + comm - the communicator that will share the eventually allocated matrix 1001 . nrows - the number of LOCAL rows in the matrix 1002 - ncols - the number of LOCAL columns in the matrix 1003 1004 Output Parameters: 1005 + dnz - the array that will be passed to the matrix preallocation routines 1006 - onz - the other array passed to the matrix preallocation routines 1007 1008 Level: deprecated (since v3.19) 1009 1010 Notes: 1011 This routine is no longer needed since assembling matrices without explicit preallocation will not be slower than 1012 the use of this routine 1013 1014 This is a macro that handles its own error checking, it does not return an error code. 1015 1016 Do not malloc or free `dnz` and `onz`, that is handled internally by these routines 1017 1018 Developer Note: 1019 This is a MACRO not a function because it has a leading { that is closed by `PetscPreallocateFinalize()`. 1020 1021 .seealso: [](chapter_matrices), `MatPreallocateEnd()`, `MatPreallocateSet()`, `MatPreallocateSymmetricSetBlock()`, `MatPreallocateSetLocal()`, 1022 `MatPreallocateSymmetricSetLocalBlock()` 1023 M*/ 1024 #define MatPreallocateBegin(comm, nrows, ncols, dnz, onz) \ 1025 do { \ 1026 PetscInt __nrows = (nrows), __ncols = (ncols), __rstart, __end = 0; \ 1027 PetscInt PETSC_UNUSED __start; \ 1028 PetscCall(PetscCalloc2(__nrows, &(dnz), __nrows, &(onz))); \ 1029 PetscCallMPI(MPI_Scan(&__ncols, &__end, 1, MPIU_INT, MPI_SUM, comm)); \ 1030 __start = __end - __ncols; \ 1031 (void)__start; \ 1032 PetscCallMPI(MPI_Scan(&__nrows, &__rstart, 1, MPIU_INT, MPI_SUM, comm)); \ 1033 __rstart -= __nrows 1034 1035 #define MatPreallocateInitialize(...) PETSC_DEPRECATED_MACRO("GCC warning \"Use MatPreallocateBegin() (since version 3.18)\"") MatPreallocateBegin(__VA_ARGS__) 1036 1037 /*MC 1038 MatPreallocateSetLocal - Indicates the locations (rows and columns) in the matrix where nonzeros will be 1039 inserted using a local number of the rows and columns 1040 1041 Synopsis: 1042 #include <petscmat.h> 1043 PetscErrorCode MatPreallocateSetLocal(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 1044 1045 Not Collective 1046 1047 Input Parameters: 1048 + map - the row mapping from local numbering to global numbering 1049 . nrows - the number of rows indicated 1050 . rows - the indices of the rows 1051 . cmap - the column mapping from local to global numbering 1052 . ncols - the number of columns in the matrix 1053 . cols - the columns indicated 1054 . dnz - the array that will be passed to the matrix preallocation routines 1055 - onz - the other array passed to the matrix preallocation routines 1056 1057 Level: deprecated (since v3.19) 1058 1059 Notes: 1060 This routine is no longer needed since assembling matrices without explicit preallocation will not be slower than 1061 the use of this routine 1062 1063 Do not malloc or free `dnz` and `onz`, that is handled internally by these routines 1064 1065 .seealso: [](chapter_matrices), `MatPreallocateEnd()`, `MatPreallocateSet()`, `MatPreallocateSymmetricSetBlock()` 1066 `MatPreallocateBegin()`, `MatPreallocateSymmetricSetLocalBlock()`, `MatPreallocateSetLocalRemoveDups()` 1067 M*/ 1068 #define MatPreallocateSetLocal(rmap, nrows, rows, cmap, ncols, cols, dnz, onz) \ 1069 PetscMacroReturnStandard(PetscCall(ISLocalToGlobalMappingApply(rmap, nrows, rows, rows)); PetscCall(ISLocalToGlobalMappingApply(cmap, ncols, cols, cols)); for (PetscInt __l = 0; __l < nrows; __l++) PetscCall(MatPreallocateSet((rows)[__l], ncols, cols, dnz, onz));) 1070 1071 /*MC 1072 MatPreallocateSetLocalRemoveDups - Indicates the locations (rows and columns) in the matrix where nonzeros will be 1073 inserted using a local number of the rows and columns. This version removes any duplicate columns in cols 1074 1075 Synopsis: 1076 #include <petscmat.h> 1077 PetscErrorCode MatPreallocateSetLocalRemoveDups(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 1078 1079 Not Collective 1080 1081 Input Parameters: 1082 + map - the row mapping from local numbering to global numbering 1083 . nrows - the number of rows indicated 1084 . rows - the indices of the rows (these values are mapped to the global values) 1085 . cmap - the column mapping from local to global numbering 1086 . ncols - the number of columns in the matrix (this value will be changed if duplicate columns are found) 1087 . cols - the columns indicated (these values are mapped to the global values, they are then sorted and duplicates removed) 1088 . dnz - the array that will be passed to the matrix preallocation routines 1089 - onz - the other array passed to the matrix preallocation routines 1090 1091 Level: deprecated (since v3.19) 1092 1093 Notes: 1094 This routine is no longer needed since assembling matrices without explicit preallocation will not be slower than 1095 the use of this routine 1096 1097 Do not malloc or free `dnz` and `onz`, that is handled internally by these routines 1098 1099 .seealso: [](chapter_matrices), `MatPreallocateEnd()`, `MatPreallocateSet()`, `MatPreallocateSymmetricSetBlock()` 1100 `MatPreallocateBegin()`, `MatPreallocateSymmetricSetLocalBlock()`, `MatPreallocateSetLocal()` 1101 M*/ 1102 #define MatPreallocateSetLocalRemoveDups(rmap, nrows, rows, cmap, ncols, cols, dnz, onz) \ 1103 PetscMacroReturnStandard(PetscCall(ISLocalToGlobalMappingApply(rmap, nrows, rows, rows)); PetscCall(ISLocalToGlobalMappingApply(cmap, ncols, cols, cols)); PetscCall(PetscSortRemoveDupsInt(&ncols, cols)); for (PetscInt __l = 0; __l < nrows; __l++) PetscCall(MatPreallocateSet((rows)[__l], ncols, cols, dnz, onz));) 1104 1105 /*MC 1106 MatPreallocateSetLocalBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be 1107 inserted using a local number of the rows and columns 1108 1109 Synopsis: 1110 #include <petscmat.h> 1111 PetscErrorCode MatPreallocateSetLocalBlock(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 1112 1113 Not Collective 1114 1115 Input Parameters: 1116 + map - the row mapping from local numbering to global numbering 1117 . nrows - the number of rows indicated 1118 . rows - the indices of the rows 1119 . cmap - the column mapping from local to global numbering 1120 . ncols - the number of columns in the matrix 1121 . cols - the columns indicated 1122 . dnz - the array that will be passed to the matrix preallocation routines 1123 - onz - the other array passed to the matrix preallocation routines 1124 1125 Level: deprecated (since v3.19) 1126 1127 Notes: 1128 This routine is no longer needed since assembling matrices without explicit preallocation will not be slower than 1129 the use of this routine 1130 1131 Do not malloc or free `dnz` and `onz`, that is handled internally by these routines 1132 1133 .seealso: [](chapter_matrices), `MatPreallocateEnd()`, `MatPreallocateSet()`, `MatPreallocateSymmetricSetBlock()` 1134 `MatPreallocateBegin()`, `MatPreallocateSymmetricSetLocalBlock()` 1135 M*/ 1136 #define MatPreallocateSetLocalBlock(rmap, nrows, rows, cmap, ncols, cols, dnz, onz) \ 1137 PetscMacroReturnStandard(PetscCall(ISLocalToGlobalMappingApplyBlock(rmap, nrows, rows, rows)); PetscCall(ISLocalToGlobalMappingApplyBlock(cmap, ncols, cols, cols)); for (PetscInt __l = 0; __l < nrows; __l++) PetscCall(MatPreallocateSet((rows)[__l], ncols, cols, dnz, onz));) 1138 1139 /*MC 1140 MatPreallocateSymmetricSetLocalBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be 1141 inserted using a local number of the rows and columns 1142 1143 Synopsis: 1144 #include <petscmat.h> 1145 PetscErrorCode MatPreallocateSymmetricSetLocalBlock(ISLocalToGlobalMappping map,PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 1146 1147 Not Collective 1148 1149 Input Parameters: 1150 + map - the mapping between local numbering and global numbering 1151 . nrows - the number of rows indicated 1152 . rows - the indices of the rows 1153 . ncols - the number of columns in the matrix 1154 . cols - the columns indicated 1155 . dnz - the array that will be passed to the matrix preallocation routines 1156 - onz - the other array passed to the matrix preallocation routines 1157 1158 Level: deprecated (since v3.19) 1159 1160 Notes: 1161 This routine is no longer needed since assembling matrices without explicit preallocation will not be slower than 1162 the use of this routine 1163 1164 Do not malloc or free `dnz` and `onz` that is handled internally by these routines 1165 1166 .seealso: [](chapter_matrices), `MatPreallocateEnd()`, `MatPreallocateSet()` 1167 `MatPreallocateBegin()`, `MatPreallocateSetLocal()` 1168 M*/ 1169 #define MatPreallocateSymmetricSetLocalBlock(map, nrows, rows, ncols, cols, dnz, onz) \ 1170 PetscMacroReturnStandard(PetscCall(ISLocalToGlobalMappingApplyBlock(map, nrows, rows, rows)); PetscCall(ISLocalToGlobalMappingApplyBlock(map, ncols, cols, cols)); for (PetscInt __l = 0; __l < nrows; __l++) PetscCall(MatPreallocateSymmetricSetBlock((rows)[__l], ncols, cols, dnz, onz));) 1171 1172 /*MC 1173 MatPreallocateSet - Indicates the locations (rows and columns) in the matrix where nonzeros will be 1174 inserted using a local number of the rows and columns 1175 1176 Synopsis: 1177 #include <petscmat.h> 1178 PetscErrorCode MatPreallocateSet(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 1179 1180 Not Collective 1181 1182 Input Parameters: 1183 + row - the row 1184 . ncols - the number of columns in the matrix 1185 - cols - the columns indicated 1186 1187 Output Parameters: 1188 + dnz - the array that will be passed to the matrix preallocation routines 1189 - onz - the other array passed to the matrix preallocation routines 1190 1191 Level: deprecated (since v3.19) 1192 1193 Notes: 1194 This routine is no longer needed since assembling matrices without explicit preallocation will not be slower than 1195 the use of this routine 1196 1197 Do not malloc or free `dnz` and `onz` that is handled internally by these routines 1198 1199 This is a MACRO not a function because it uses variables declared in MatPreallocateBegin(). 1200 1201 .seealso: [](chapter_matrices), `MatPreallocateEnd()`, `MatPreallocateSet()`, `MatPreallocateSymmetricSetBlock()` 1202 `MatPreallocateBegin()`, `MatPreallocateSetLocal()` 1203 M*/ 1204 #define MatPreallocateSet(row, nc, cols, dnz, onz) \ 1205 PetscMacroReturnStandard(PetscCheck(row >= __rstart, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Trying to set preallocation for row %" PetscInt_FMT " less than first local row %" PetscInt_FMT, row, __rstart); PetscCheck(row < __rstart + __nrows, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Trying to set preallocation for row %" PetscInt_FMT " greater than last local row %" PetscInt_FMT, row, __rstart + __nrows - 1); for (PetscInt __i = 0; __i < nc; ++__i) { \ 1206 if ((cols)[__i] < __start || (cols)[__i] >= __end) onz[row - __rstart]++; \ 1207 else if (dnz[row - __rstart] < __ncols) dnz[row - __rstart]++; \ 1208 }) 1209 1210 /*MC 1211 MatPreallocateSymmetricSetBlock - Indicates the locations (rows and columns) in the matrix where nonzeros will be 1212 inserted using a local number of the rows and columns 1213 1214 Synopsis: 1215 #include <petscmat.h> 1216 PetscErrorCode MatPreallocateSymmetricSetBlock(PetscInt nrows, PetscInt *rows,PetscInt ncols, PetscInt *cols,PetscInt *dnz, PetscInt *onz) 1217 1218 Not Collective 1219 1220 Input Parameters: 1221 + nrows - the number of rows indicated 1222 . rows - the indices of the rows 1223 . ncols - the number of columns in the matrix 1224 . cols - the columns indicated 1225 . dnz - the array that will be passed to the matrix preallocation routines 1226 - onz - the other array passed to the matrix preallocation routines 1227 1228 Level: deprecated (since v3.19) 1229 1230 Notes: 1231 This routine is no longer needed since assembling matrices without explicit preallocation will not be slower than 1232 the use of this routine 1233 1234 Do not malloc or free `dnz` and `onz` that is handled internally by these routines 1235 1236 This is a MACRO not a function because it uses variables declared in MatPreallocateBegin(). 1237 1238 .seealso: [](chapter_matrices), `MatPreallocateEnd()`, `MatPreallocateSet()`, `MatPreallocateBegin()`, 1239 `MatPreallocateSymmetricSetLocalBlock()`, `MatPreallocateSetLocal()` 1240 M*/ 1241 #define MatPreallocateSymmetricSetBlock(row, nc, cols, dnz, onz) \ 1242 PetscMacroReturnStandard(for (PetscInt __i = 0; __i < nc; __i++) { \ 1243 if (cols[__i] >= __end) onz[row - __rstart]++; \ 1244 else if (cols[__i] >= row && dnz[row - __rstart] < __ncols) dnz[row - __rstart]++; \ 1245 }) 1246 1247 /*MC 1248 MatPreallocateLocation - An alternative to MatPreallocateSet() that puts the nonzero locations into the matrix if it exists 1249 1250 Synopsis: 1251 #include <petscmat.h> 1252 PetscErrorCode MatPreallocateLocations(Mat A,PetscInt row,PetscInt ncols,PetscInt *cols,PetscInt *dnz,PetscInt *onz) 1253 1254 Not Collective 1255 1256 Input Parameters: 1257 + A - matrix 1258 . row - row where values exist (must be local to this process) 1259 . ncols - number of columns 1260 . cols - columns with nonzeros 1261 . dnz - the array that will be passed to the matrix preallocation routines 1262 - onz - the other array passed to the matrix preallocation routines 1263 1264 Level: deprecated (since v3.19) 1265 1266 Notes: 1267 This routine is no longer needed since assembling matrices without explicit preallocation will not be slower than 1268 the use of this routine 1269 1270 Do not malloc or free `dnz` and `onz` that is handled internally by these routines 1271 1272 This is a MACRO not a function because it uses a bunch of variables private to the MatPreallocation.... routines. 1273 1274 .seealso: [](chapter_matrices), `MatPreallocateBegin()`, `MatPreallocateSet()`, `MatPreallocateSymmetricSetBlock()`, `MatPreallocateSetLocal()`, 1275 `MatPreallocateSymmetricSetLocalBlock()` 1276 M*/ 1277 #define MatPreallocateLocation(A, row, ncols, cols, dnz, onz) (A ? MatSetValues(A, 1, &row, ncols, cols, NULL, INSERT_VALUES) : MatPreallocateSet(row, ncols, cols, dnz, onz)) 1278 1279 /*MC 1280 MatPreallocateEnd - Ends the block of code that will count the number of nonzeros per 1281 row in a matrix providing the data that one can use to correctly preallocate the matrix. 1282 1283 Synopsis: 1284 #include <petscmat.h> 1285 PetscErrorCode MatPreallocateEnd(PetscInt *dnz, PetscInt *onz) 1286 1287 Collective 1288 1289 Input Parameters: 1290 + dnz - the array that was be passed to the matrix preallocation routines 1291 - onz - the other array passed to the matrix preallocation routines 1292 1293 Level: deprecated (since v3.19) 1294 1295 Notes: 1296 This routine is no longer needed since assembling matrices without explicit preallocation will not be slower than 1297 the use of this routine 1298 1299 Do not malloc or free `dnz` and `onz`, that is handled internally by these routines 1300 1301 Developer Note: 1302 This is a MACRO not a function because it closes the { started in MatPreallocateBegin(). 1303 1304 .seealso: [](chapter_matrices), `MatPreallocateBegin()`, `MatPreallocateSet()`, `MatPreallocateSymmetricSetBlock()`, `MatPreallocateSetLocal()`, 1305 `MatPreallocateSymmetricSetLocalBlock()` 1306 M*/ 1307 #define MatPreallocateEnd(dnz, onz) \ 1308 PetscCall(PetscFree2(dnz, onz)); \ 1309 } \ 1310 while (0) 1311 1312 #define MatPreallocateFinalize(...) PETSC_DEPRECATED_MACRO("GCC warning \"Use MatPreallocateEnd() (since version 3.18)\"") MatPreallocateEnd(__VA_ARGS__) 1313 1314 /* Routines unique to particular data structures */ 1315 PETSC_EXTERN PetscErrorCode MatShellGetContext(Mat, void *); 1316 1317 PETSC_EXTERN PetscErrorCode MatInodeAdjustForInodes(Mat, IS *, IS *); 1318 PETSC_EXTERN PetscErrorCode MatInodeGetInodeSizes(Mat, PetscInt *, PetscInt *[], PetscInt *); 1319 1320 PETSC_EXTERN PetscErrorCode MatSeqAIJSetColumnIndices(Mat, PetscInt[]); 1321 PETSC_EXTERN PetscErrorCode MatSeqBAIJSetColumnIndices(Mat, PetscInt[]); 1322 PETSC_EXTERN PetscErrorCode MatCreateSeqAIJWithArrays(MPI_Comm, PetscInt, PetscInt, PetscInt[], PetscInt[], PetscScalar[], Mat *); 1323 PETSC_EXTERN PetscErrorCode MatCreateSeqBAIJWithArrays(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt[], PetscInt[], PetscScalar[], Mat *); 1324 PETSC_EXTERN PetscErrorCode MatCreateSeqSBAIJWithArrays(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt[], PetscInt[], PetscScalar[], Mat *); 1325 PETSC_EXTERN PetscErrorCode MatCreateSeqAIJFromTriple(MPI_Comm, PetscInt, PetscInt, PetscInt[], PetscInt[], PetscScalar[], Mat *, PetscInt, PetscBool); 1326 1327 #define MAT_SKIP_ALLOCATION -4 1328 1329 PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocation(Mat, PetscInt, PetscInt, const PetscInt[]); 1330 PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetPreallocation(Mat, PetscInt, PetscInt, const PetscInt[]); 1331 PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocation(Mat, PetscInt, const PetscInt[]); 1332 PETSC_EXTERN PetscErrorCode MatSeqAIJSetTotalPreallocation(Mat, PetscInt); 1333 1334 PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocation(Mat, PetscInt, PetscInt, const PetscInt[], PetscInt, const PetscInt[]); 1335 PETSC_EXTERN PetscErrorCode MatMPISBAIJSetPreallocation(Mat, PetscInt, PetscInt, const PetscInt[], PetscInt, const PetscInt[]); 1336 PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocation(Mat, PetscInt, const PetscInt[], PetscInt, const PetscInt[]); 1337 PETSC_EXTERN PetscErrorCode MatSeqAIJSetPreallocationCSR(Mat, const PetscInt[], const PetscInt[], const PetscScalar[]); 1338 PETSC_EXTERN PetscErrorCode MatSeqBAIJSetPreallocationCSR(Mat, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[]); 1339 PETSC_EXTERN PetscErrorCode MatMPIAIJSetPreallocationCSR(Mat, const PetscInt[], const PetscInt[], const PetscScalar[]); 1340 PETSC_EXTERN PetscErrorCode MatMPIBAIJSetPreallocationCSR(Mat, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[]); 1341 PETSC_EXTERN PetscErrorCode MatMPIAdjSetPreallocation(Mat, PetscInt[], PetscInt[], PetscInt[]); 1342 PETSC_EXTERN PetscErrorCode MatMPIAdjToSeq(Mat, Mat *); 1343 PETSC_EXTERN PetscErrorCode MatMPIAdjToSeqRankZero(Mat, Mat *); 1344 PETSC_EXTERN PetscErrorCode MatMPIDenseSetPreallocation(Mat, PetscScalar[]); 1345 PETSC_EXTERN PetscErrorCode MatSeqDenseSetPreallocation(Mat, PetscScalar[]); 1346 PETSC_EXTERN PetscErrorCode MatMPIAIJGetSeqAIJ(Mat, Mat *, Mat *, const PetscInt *[]); 1347 PETSC_EXTERN PetscErrorCode MatMPIBAIJGetSeqBAIJ(Mat, Mat *, Mat *, const PetscInt *[]); 1348 PETSC_EXTERN PetscErrorCode MatMPIAdjCreateNonemptySubcommMat(Mat, Mat *); 1349 1350 PETSC_EXTERN PetscErrorCode MatDenseGetLDA(Mat, PetscInt *); 1351 PETSC_EXTERN PetscErrorCode MatDenseSetLDA(Mat, PetscInt); 1352 PETSC_DEPRECATED_FUNCTION("Use MatDenseSetLDA() (since version 3.14)") static inline PetscErrorCode MatSeqDenseSetLDA(Mat A, PetscInt lda) 1353 { 1354 return MatDenseSetLDA(A, lda); 1355 } 1356 PETSC_EXTERN PetscErrorCode MatDenseGetLocalMatrix(Mat, Mat *); 1357 1358 PETSC_EXTERN PetscErrorCode MatBlockMatSetPreallocation(Mat, PetscInt, PetscInt, const PetscInt[]); 1359 1360 PETSC_EXTERN PetscErrorCode MatStoreValues(Mat); 1361 PETSC_EXTERN PetscErrorCode MatRetrieveValues(Mat); 1362 1363 PETSC_EXTERN PetscErrorCode MatFindNonzeroRows(Mat, IS *); 1364 PETSC_EXTERN PetscErrorCode MatFindZeroRows(Mat, IS *); 1365 /* 1366 These routines are not usually accessed directly, rather solving is 1367 done through the KSP and PC interfaces. 1368 */ 1369 1370 /*J 1371 MatOrderingType - String with the name of a PETSc matrix ordering. These orderings are most commonly used 1372 to reduce fill in sparse factorizations. 1373 1374 Level: beginner 1375 1376 Notes: 1377 If `MATORDERINGEXTERNAL` is used then PETSc does not compute an ordering and instead the external package called utilizes one 1378 in the external the factorization package 1379 1380 There is no `MatOrdering` object, the ordering is obtained directly from the matrix with `MatGetOrdering()` 1381 1382 Developer Note: 1383 This API should be converted to an API similar to those for `MatColoring` and `MatPartitioning` 1384 1385 .seealso: [](chapter_matrices), [](sec_graph), `MatGetFactor()`, `MatGetOrdering()`, `MatColoringType`, `MatPartitioningType`, `MatCoarsenType`, `MatCoarsenType`, 1386 `PCFactorSetOrderingType()` 1387 J*/ 1388 typedef const char *MatOrderingType; 1389 #define MATORDERINGNATURAL "natural" 1390 #define MATORDERINGND "nd" 1391 #define MATORDERING1WD "1wd" 1392 #define MATORDERINGRCM "rcm" 1393 #define MATORDERINGQMD "qmd" 1394 #define MATORDERINGROWLENGTH "rowlength" 1395 #define MATORDERINGWBM "wbm" 1396 #define MATORDERINGSPECTRAL "spectral" 1397 #define MATORDERINGAMD "amd" /* only works if UMFPACK is installed with PETSc */ 1398 #define MATORDERINGMETISND "metisnd" /* only works if METIS is installed with PETSc */ 1399 #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 */ 1400 #define MATORDERINGEXTERNAL "external" /* uses an ordering type internal to the factorization package */ 1401 1402 PETSC_EXTERN PetscErrorCode MatGetOrdering(Mat, MatOrderingType, IS *, IS *); 1403 PETSC_EXTERN PetscErrorCode MatGetOrderingList(PetscFunctionList *); 1404 PETSC_EXTERN PetscErrorCode MatOrderingRegister(const char[], PetscErrorCode (*)(Mat, MatOrderingType, IS *, IS *)); 1405 PETSC_EXTERN PetscFunctionList MatOrderingList; 1406 1407 #include "petscmatcoarsen.h" 1408 1409 PETSC_EXTERN PetscErrorCode MatReorderForNonzeroDiagonal(Mat, PetscReal, IS, IS); 1410 PETSC_EXTERN PetscErrorCode MatCreateLaplacian(Mat, PetscReal, PetscBool, Mat *); 1411 1412 PETSC_EXTERN PetscErrorCode MatFactorGetPreferredOrdering(Mat, MatFactorType, MatOrderingType *); 1413 1414 /*S 1415 MatFactorShiftType - Type of numeric shift used for factorizations 1416 1417 Values: 1418 + `MAT_SHIFT_NONE` - do not shift the matrix diagonal entries 1419 . `MAT_SHIFT_NONZERO` - shift the entries to be non-zero 1420 . `MAT_SHIFT_POSITIVE_DEFINITE` - shift the entries to force the factorization to be positive definite 1421 - `MAT_SHIFT_INBLOCKS` - only shift the factors inside the small dense diagonal blocks of the matrix, for example with `MATBAIJ` 1422 1423 Level: intermediate 1424 1425 .seealso: [](chapter_matrices), `Mat`, `MatGetFactor()`, `PCFactorSetShiftType()` 1426 S*/ 1427 typedef enum { 1428 MAT_SHIFT_NONE, 1429 MAT_SHIFT_NONZERO, 1430 MAT_SHIFT_POSITIVE_DEFINITE, 1431 MAT_SHIFT_INBLOCKS 1432 } MatFactorShiftType; 1433 PETSC_EXTERN const char *const MatFactorShiftTypes[]; 1434 PETSC_EXTERN const char *const MatFactorShiftTypesDetail[]; 1435 1436 /*S 1437 MatFactorError - indicates what type of error was generated in a matrix factorization 1438 1439 Values: 1440 + `MAT_FACTOR_NOERROR` - there was no error during the factorization 1441 . `MAT_FACTOR_STRUCT_ZEROPIVOT` - there was a missing entry in a diagonal location of the matrix 1442 . `MAT_FACTOR_NUMERIC_ZEROPIVOT` - there was a (near) zero pivot during the factorization 1443 . `MAT_FACTOR_OUTMEMORY` - the factorization has run out of memory 1444 - `MAT_FACTOR_OTHER` - some other error has occurred. 1445 1446 Level: intermediate 1447 1448 Note: 1449 When a factorization is done in a preconditioner `PC` the error may be propagated up to a `PCFailedReason` or a `KSPConvergedReason` 1450 1451 .seealso: [](chapter_matrices), `Mat`, `MatGetFactor()`, `MatFactorGetError()`, `MatFactorGetErrorZeroPivot()`, `MatFactorClearError()`, 1452 `PCFailedReason`, `PCGetFailedReason()`, `KSPConvergedReason` 1453 S*/ 1454 typedef enum { 1455 MAT_FACTOR_NOERROR, 1456 MAT_FACTOR_STRUCT_ZEROPIVOT, 1457 MAT_FACTOR_NUMERIC_ZEROPIVOT, 1458 MAT_FACTOR_OUTMEMORY, 1459 MAT_FACTOR_OTHER 1460 } MatFactorError; 1461 1462 PETSC_EXTERN PetscErrorCode MatFactorGetError(Mat, MatFactorError *); 1463 PETSC_EXTERN PetscErrorCode MatFactorClearError(Mat); 1464 PETSC_EXTERN PetscErrorCode MatFactorGetErrorZeroPivot(Mat, PetscReal *, PetscInt *); 1465 1466 /*S 1467 MatFactorInfo - Data passed into the matrix factorization routines, and information about the resulting factorization 1468 1469 Level: developer 1470 1471 Note: 1472 You can use `MatFactorInfoInitialize()` to set default values. 1473 1474 Fortran Note: 1475 The data is available in a double precision arrays of size `MAT_FACTORINFO_SIZE` 1476 1477 .seealso: [](chapter_matrices), `Mat`, `MatInfo`, `MatGetFactor()`, `MatLUFactorSymbolic()`, `MatILUFactorSymbolic()`, `MatCholeskyFactorSymbolic()`, 1478 `MatICCFactorSymbolic()`, `MatICCFactor()`, `MatFactorInfoInitialize()` 1479 S*/ 1480 typedef struct { 1481 PetscReal diagonal_fill; /* force diagonal to fill in if initially not filled */ 1482 PetscReal usedt; 1483 PetscReal dt; /* drop tolerance */ 1484 PetscReal dtcol; /* tolerance for pivoting */ 1485 PetscReal dtcount; /* maximum nonzeros to be allowed per row */ 1486 PetscReal fill; /* expected fill, nonzeros in factored matrix/nonzeros in original matrix */ 1487 PetscReal levels; /* ICC/ILU(levels) */ 1488 PetscReal pivotinblocks; /* for BAIJ and SBAIJ matrices pivot in factorization on blocks, default 1.0 1489 factorization may be faster if do not pivot */ 1490 PetscReal zeropivot; /* pivot is called zero if less than this */ 1491 PetscReal shifttype; /* type of shift added to matrix factor to prevent zero pivots */ 1492 PetscReal shiftamount; /* how large the shift is */ 1493 } MatFactorInfo; 1494 1495 PETSC_EXTERN PetscErrorCode MatFactorInfoInitialize(MatFactorInfo *); 1496 PETSC_EXTERN PetscErrorCode MatCholeskyFactor(Mat, IS, const MatFactorInfo *); 1497 PETSC_EXTERN PetscErrorCode MatCholeskyFactorSymbolic(Mat, Mat, IS, const MatFactorInfo *); 1498 PETSC_EXTERN PetscErrorCode MatCholeskyFactorNumeric(Mat, Mat, const MatFactorInfo *); 1499 PETSC_EXTERN PetscErrorCode MatLUFactor(Mat, IS, IS, const MatFactorInfo *); 1500 PETSC_EXTERN PetscErrorCode MatILUFactor(Mat, IS, IS, const MatFactorInfo *); 1501 PETSC_EXTERN PetscErrorCode MatLUFactorSymbolic(Mat, Mat, IS, IS, const MatFactorInfo *); 1502 PETSC_EXTERN PetscErrorCode MatILUFactorSymbolic(Mat, Mat, IS, IS, const MatFactorInfo *); 1503 PETSC_EXTERN PetscErrorCode MatICCFactorSymbolic(Mat, Mat, IS, const MatFactorInfo *); 1504 PETSC_EXTERN PetscErrorCode MatICCFactor(Mat, IS, const MatFactorInfo *); 1505 PETSC_EXTERN PetscErrorCode MatLUFactorNumeric(Mat, Mat, const MatFactorInfo *); 1506 PETSC_EXTERN PetscErrorCode MatQRFactor(Mat, IS, const MatFactorInfo *); 1507 PETSC_EXTERN PetscErrorCode MatQRFactorSymbolic(Mat, Mat, IS, const MatFactorInfo *); 1508 PETSC_EXTERN PetscErrorCode MatQRFactorNumeric(Mat, Mat, const MatFactorInfo *); 1509 PETSC_EXTERN PetscErrorCode MatGetInertia(Mat, PetscInt *, PetscInt *, PetscInt *); 1510 PETSC_EXTERN PetscErrorCode MatSolve(Mat, Vec, Vec); 1511 PETSC_EXTERN PetscErrorCode MatForwardSolve(Mat, Vec, Vec); 1512 PETSC_EXTERN PetscErrorCode MatBackwardSolve(Mat, Vec, Vec); 1513 PETSC_EXTERN PetscErrorCode MatSolveAdd(Mat, Vec, Vec, Vec); 1514 PETSC_EXTERN PetscErrorCode MatSolveTranspose(Mat, Vec, Vec); 1515 PETSC_EXTERN PetscErrorCode MatSolveTransposeAdd(Mat, Vec, Vec, Vec); 1516 PETSC_EXTERN PetscErrorCode MatSolves(Mat, Vecs, Vecs); 1517 PETSC_EXTERN PetscErrorCode MatSetUnfactored(Mat); 1518 1519 typedef enum { 1520 MAT_FACTOR_SCHUR_UNFACTORED, 1521 MAT_FACTOR_SCHUR_FACTORED, 1522 MAT_FACTOR_SCHUR_INVERTED 1523 } MatFactorSchurStatus; 1524 PETSC_EXTERN PetscErrorCode MatFactorSetSchurIS(Mat, IS); 1525 PETSC_EXTERN PetscErrorCode MatFactorGetSchurComplement(Mat, Mat *, MatFactorSchurStatus *); 1526 PETSC_EXTERN PetscErrorCode MatFactorRestoreSchurComplement(Mat, Mat *, MatFactorSchurStatus); 1527 PETSC_EXTERN PetscErrorCode MatFactorInvertSchurComplement(Mat); 1528 PETSC_EXTERN PetscErrorCode MatFactorCreateSchurComplement(Mat, Mat *, MatFactorSchurStatus *); 1529 PETSC_EXTERN PetscErrorCode MatFactorSolveSchurComplement(Mat, Vec, Vec); 1530 PETSC_EXTERN PetscErrorCode MatFactorSolveSchurComplementTranspose(Mat, Vec, Vec); 1531 PETSC_EXTERN PetscErrorCode MatFactorFactorizeSchurComplement(Mat); 1532 1533 PETSC_EXTERN PetscErrorCode MatSeqDenseInvert(Mat); 1534 /*E 1535 MatSORType - What type of (S)SOR to perform 1536 1537 Values: 1538 + `SOR_FORWARD_SWEEP` - do a sweep from the first row of the matrix to the last 1539 . `SOR_BACKWARD_SWEEP` - do a sweep from the last row to the first 1540 . `SOR_SYMMETRIC_SWEEP` - do a sweep from the first row to the last and then back to the first 1541 . `SOR_LOCAL_FORWARD_SWEEP` - each MPI rank does its own forward sweep with no communication 1542 . `SOR_LOCAL_BACKWARD_SWEEP` - each MPI rank does its own backward sweep with no communication 1543 . `SOR_LOCAL_SYMMETRIC_SWEEP` - each MPI rank does its own symmetric sweep with no communication 1544 . `SOR_ZERO_INITIAL_GUESS` - indicates the initial solution is zero so the sweep can avoid unneeded computation 1545 . `SOR_EISENSTAT` - apply the Eisentat application of SOR, see `PCEISENSTAT` 1546 . `SOR_APPLY_UPPER` - multiply by the upper triangular portion of the matrix 1547 - `SOR_APPLY_LOWER` - multiply by the lower triangular portion of the matrix 1548 1549 Level: beginner 1550 1551 Note: 1552 These may be bitwise ORd together 1553 1554 Developer Note: 1555 Since `MatSORType` may be bitwise ORd together, so do not change the numerical values below 1556 1557 .seealso: [](chapter_matrices), `MatSOR()` 1558 E*/ 1559 typedef enum { 1560 SOR_FORWARD_SWEEP = 1, 1561 SOR_BACKWARD_SWEEP = 2, 1562 SOR_SYMMETRIC_SWEEP = 3, 1563 SOR_LOCAL_FORWARD_SWEEP = 4, 1564 SOR_LOCAL_BACKWARD_SWEEP = 8, 1565 SOR_LOCAL_SYMMETRIC_SWEEP = 12, 1566 SOR_ZERO_INITIAL_GUESS = 16, 1567 SOR_EISENSTAT = 32, 1568 SOR_APPLY_UPPER = 64, 1569 SOR_APPLY_LOWER = 128 1570 } MatSORType; 1571 PETSC_EXTERN PetscErrorCode MatSOR(Mat, Vec, PetscReal, MatSORType, PetscReal, PetscInt, PetscInt, Vec); 1572 1573 /*S 1574 MatColoring - Object for managing the coloring of matrices. 1575 1576 Level: beginner 1577 1578 Notes: 1579 Coloring of matrices can be computed directly from the sparse matrix nonzero structure via the `MatColoring` object or from the mesh from which the 1580 matrix comes from via `DMCreateColoring()`. In general using the mesh produces a more optimal coloring (fewer colors). 1581 1582 Once a coloring is available `MatFDColoringCreate()` creates an object that can be used to efficiently compute Jacobians using that coloring. This 1583 same object can also be used to efficiently convert data created by Automatic Differentiation tools to PETSc sparse matrices. 1584 1585 .seealso: [](chapter_matrices), [](sec_graph), `MatFDColoringCreate()`, `MatColoringWeightType`, `ISColoring`, `MatFDColoring`, `DMCreateColoring()`, `MatColoringCreate()`, 1586 `MatPartitioning`, `MatColoringType`, `MatPartitioningType`, `MatOrderingType`, `MatColoringSetWeightType()`, 1587 `MatColoringSetWeights()`, `MatCoarsenType`, `MatCoarsen` 1588 S*/ 1589 typedef struct _p_MatColoring *MatColoring; 1590 1591 /*J 1592 MatColoringType - String with the name of a PETSc matrix coloring 1593 1594 Level: beginner 1595 1596 .seealso: [](chapter_matrices), [](sec_graph), `Mat`, `MatFDColoringCreate()`, `MatColoringSetType()`, `MatColoring` 1597 J*/ 1598 typedef const char *MatColoringType; 1599 #define MATCOLORINGJP "jp" 1600 #define MATCOLORINGPOWER "power" 1601 #define MATCOLORINGNATURAL "natural" 1602 #define MATCOLORINGSL "sl" 1603 #define MATCOLORINGLF "lf" 1604 #define MATCOLORINGID "id" 1605 #define MATCOLORINGGREEDY "greedy" 1606 1607 /*E 1608 MatColoringWeightType - Type of weight scheme used for the coloring algorithm 1609 1610 Values: 1611 + `MAT_COLORING_RANDOM` - Random weights 1612 . `MAT_COLORING_LEXICAL` - Lexical weighting based upon global numbering. 1613 - `MAT_COLORING_LF` - Last-first weighting. 1614 1615 Level: intermediate 1616 1617 .seealso: [](chapter_matrices), `MatColoring`, `MatColoringCreate()`, `MatColoringSetWeightType()`, `MatColoringSetWeights()` 1618 E*/ 1619 typedef enum { 1620 MAT_COLORING_WEIGHT_RANDOM, 1621 MAT_COLORING_WEIGHT_LEXICAL, 1622 MAT_COLORING_WEIGHT_LF, 1623 MAT_COLORING_WEIGHT_SL 1624 } MatColoringWeightType; 1625 1626 PETSC_EXTERN PetscErrorCode MatColoringCreate(Mat, MatColoring *); 1627 PETSC_EXTERN PetscErrorCode MatColoringGetDegrees(Mat, PetscInt, PetscInt *); 1628 PETSC_EXTERN PetscErrorCode MatColoringDestroy(MatColoring *); 1629 PETSC_EXTERN PetscErrorCode MatColoringView(MatColoring, PetscViewer); 1630 PETSC_EXTERN PetscErrorCode MatColoringSetType(MatColoring, MatColoringType); 1631 PETSC_EXTERN PetscErrorCode MatColoringSetFromOptions(MatColoring); 1632 PETSC_EXTERN PetscErrorCode MatColoringSetDistance(MatColoring, PetscInt); 1633 PETSC_EXTERN PetscErrorCode MatColoringGetDistance(MatColoring, PetscInt *); 1634 PETSC_EXTERN PetscErrorCode MatColoringSetMaxColors(MatColoring, PetscInt); 1635 PETSC_EXTERN PetscErrorCode MatColoringGetMaxColors(MatColoring, PetscInt *); 1636 PETSC_EXTERN PetscErrorCode MatColoringApply(MatColoring, ISColoring *); 1637 PETSC_EXTERN PetscErrorCode MatColoringRegister(const char[], PetscErrorCode (*)(MatColoring)); 1638 PETSC_EXTERN PetscErrorCode MatColoringPatch(Mat, PetscInt, PetscInt, ISColoringValue[], ISColoring *); 1639 PETSC_EXTERN PetscErrorCode MatColoringSetWeightType(MatColoring, MatColoringWeightType); 1640 PETSC_EXTERN PetscErrorCode MatColoringSetWeights(MatColoring, PetscReal *, PetscInt *); 1641 PETSC_EXTERN PetscErrorCode MatColoringCreateWeights(MatColoring, PetscReal **, PetscInt **lperm); 1642 PETSC_EXTERN PetscErrorCode MatColoringTest(MatColoring, ISColoring); 1643 PETSC_DEPRECATED_FUNCTION("Use MatColoringTest() (since version 3.10)") static inline PetscErrorCode MatColoringTestValid(MatColoring matcoloring, ISColoring iscoloring) 1644 { 1645 return MatColoringTest(matcoloring, iscoloring); 1646 } 1647 PETSC_EXTERN PetscErrorCode MatISColoringTest(Mat, ISColoring); 1648 1649 /*S 1650 MatFDColoring - Object for computing a sparse Jacobian via finite differences with coloring 1651 1652 Level: beginner 1653 1654 Notes: 1655 This object is creating utilizing a coloring provided by the `MatColoring` object or `DMCreateColoring()` 1656 1657 The `SNES` option `-snes_fd_coloring` will cause the Jacobian needed by `SNES` to be computed via a use of this object 1658 1659 .seealso: [](chapter_matrices), `Mat`, `MatFDColoringCreate()`, `MatFDColoringSetFunction()`, `MatColoring`, `DMCreateColoring()` 1660 S*/ 1661 typedef struct _p_MatFDColoring *MatFDColoring; 1662 1663 PETSC_EXTERN PetscErrorCode MatFDColoringCreate(Mat, ISColoring, MatFDColoring *); 1664 PETSC_EXTERN PetscErrorCode MatFDColoringDestroy(MatFDColoring *); 1665 PETSC_EXTERN PetscErrorCode MatFDColoringView(MatFDColoring, PetscViewer); 1666 PETSC_EXTERN PetscErrorCode MatFDColoringSetFunction(MatFDColoring, PetscErrorCode (*)(void), void *); 1667 PETSC_EXTERN PetscErrorCode MatFDColoringGetFunction(MatFDColoring, PetscErrorCode (**)(void), void **); 1668 PETSC_EXTERN PetscErrorCode MatFDColoringSetParameters(MatFDColoring, PetscReal, PetscReal); 1669 PETSC_EXTERN PetscErrorCode MatFDColoringSetFromOptions(MatFDColoring); 1670 PETSC_EXTERN PetscErrorCode MatFDColoringApply(Mat, MatFDColoring, Vec, void *); 1671 PETSC_EXTERN PetscErrorCode MatFDColoringSetF(MatFDColoring, Vec); 1672 PETSC_EXTERN PetscErrorCode MatFDColoringGetPerturbedColumns(MatFDColoring, PetscInt *, const PetscInt *[]); 1673 PETSC_EXTERN PetscErrorCode MatFDColoringSetUp(Mat, ISColoring, MatFDColoring); 1674 PETSC_EXTERN PetscErrorCode MatFDColoringSetBlockSize(MatFDColoring, PetscInt, PetscInt); 1675 PETSC_EXTERN PetscErrorCode MatFDColoringSetValues(Mat, MatFDColoring, const PetscScalar *); 1676 1677 /*S 1678 MatTransposeColoring - Object for computing a sparse matrix product C=A*B^T via coloring 1679 1680 Level: developer 1681 1682 .seealso: [](chapter_matrices), `Mat`, `MatProductType`, `MatTransposeColoringCreate()` 1683 S*/ 1684 typedef struct _p_MatTransposeColoring *MatTransposeColoring; 1685 1686 PETSC_EXTERN PetscErrorCode MatTransposeColoringCreate(Mat, ISColoring, MatTransposeColoring *); 1687 PETSC_EXTERN PetscErrorCode MatTransColoringApplySpToDen(MatTransposeColoring, Mat, Mat); 1688 PETSC_EXTERN PetscErrorCode MatTransColoringApplyDenToSp(MatTransposeColoring, Mat, Mat); 1689 PETSC_EXTERN PetscErrorCode MatTransposeColoringDestroy(MatTransposeColoring *); 1690 1691 /*S 1692 MatPartitioning - Object for managing the partitioning of a matrix or graph 1693 1694 Level: beginner 1695 1696 Note: 1697 There is also a `PetscPartitioner` object that provides the same functionality. It can utilize the `MatPartitioning` operations 1698 via `PetscPartitionerSetType`(p,`PETSCPARTITIONERMATPARTITIONING`) 1699 1700 Developers Note: 1701 It is an extra maintenance and documentation cost to have two objects with the same functionality. `PetscPartitioner` should be removed 1702 1703 .seealso: [](chapter_matrices), [](sec_graph), `Mat`, `MatPartitioningCreate()`, `MatPartitioningType`, `MatColoring`, `MatGetOrdering()`, `MatOrderingType`, 1704 `MatCoarsenType`, `MatCoarsenType` 1705 S*/ 1706 typedef struct _p_MatPartitioning *MatPartitioning; 1707 1708 /*J 1709 MatPartitioningType - String with the name of a PETSc matrix partitioning 1710 1711 Level: beginner 1712 dm 1713 .seealso: [](chapter_matrices), [](sec_graph), `Mat`, `MatPartitioningCreate()`, `MatPartitioning`, `MatPartitioningSetType()`, `MatColoringType`, `MatOrderingType`, 1714 `MatCoarsenType`, `MatCoarsenType` 1715 J*/ 1716 typedef const char *MatPartitioningType; 1717 #define MATPARTITIONINGCURRENT "current" 1718 #define MATPARTITIONINGAVERAGE "average" 1719 #define MATPARTITIONINGSQUARE "square" 1720 #define MATPARTITIONINGPARMETIS "parmetis" 1721 #define MATPARTITIONINGCHACO "chaco" 1722 #define MATPARTITIONINGPARTY "party" 1723 #define MATPARTITIONINGPTSCOTCH "ptscotch" 1724 #define MATPARTITIONINGHIERARCH "hierarch" 1725 1726 PETSC_EXTERN PetscErrorCode MatPartitioningCreate(MPI_Comm, MatPartitioning *); 1727 PETSC_EXTERN PetscErrorCode MatPartitioningSetType(MatPartitioning, MatPartitioningType); 1728 PETSC_EXTERN PetscErrorCode MatPartitioningSetNParts(MatPartitioning, PetscInt); 1729 PETSC_EXTERN PetscErrorCode MatPartitioningSetAdjacency(MatPartitioning, Mat); 1730 PETSC_EXTERN PetscErrorCode MatPartitioningSetNumberVertexWeights(MatPartitioning, PetscInt); 1731 PETSC_EXTERN PetscErrorCode MatPartitioningSetVertexWeights(MatPartitioning, const PetscInt[]); 1732 PETSC_EXTERN PetscErrorCode MatPartitioningSetPartitionWeights(MatPartitioning, const PetscReal[]); 1733 PETSC_EXTERN PetscErrorCode MatPartitioningSetUseEdgeWeights(MatPartitioning, PetscBool); 1734 PETSC_EXTERN PetscErrorCode MatPartitioningGetUseEdgeWeights(MatPartitioning, PetscBool *); 1735 PETSC_EXTERN PetscErrorCode MatPartitioningApply(MatPartitioning, IS *); 1736 PETSC_EXTERN PetscErrorCode MatPartitioningImprove(MatPartitioning, IS *); 1737 PETSC_EXTERN PetscErrorCode MatPartitioningViewImbalance(MatPartitioning, IS); 1738 PETSC_EXTERN PetscErrorCode MatPartitioningApplyND(MatPartitioning, IS *); 1739 PETSC_EXTERN PetscErrorCode MatPartitioningDestroy(MatPartitioning *); 1740 PETSC_EXTERN PetscErrorCode MatPartitioningRegister(const char[], PetscErrorCode (*)(MatPartitioning)); 1741 PETSC_EXTERN PetscErrorCode MatPartitioningView(MatPartitioning, PetscViewer); 1742 PETSC_EXTERN PetscErrorCode MatPartitioningViewFromOptions(MatPartitioning, PetscObject, const char[]); 1743 PETSC_EXTERN PetscErrorCode MatPartitioningSetFromOptions(MatPartitioning); 1744 PETSC_EXTERN PetscErrorCode MatPartitioningGetType(MatPartitioning, MatPartitioningType *); 1745 1746 PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetRepartition(MatPartitioning part); 1747 PETSC_EXTERN PetscErrorCode MatPartitioningParmetisSetCoarseSequential(MatPartitioning); 1748 PETSC_EXTERN PetscErrorCode MatPartitioningParmetisGetEdgeCut(MatPartitioning, PetscInt *); 1749 1750 typedef enum { 1751 MP_CHACO_MULTILEVEL = 1, 1752 MP_CHACO_SPECTRAL = 2, 1753 MP_CHACO_LINEAR = 4, 1754 MP_CHACO_RANDOM = 5, 1755 MP_CHACO_SCATTERED = 6 1756 } MPChacoGlobalType; 1757 PETSC_EXTERN const char *const MPChacoGlobalTypes[]; 1758 typedef enum { 1759 MP_CHACO_KERNIGHAN = 1, 1760 MP_CHACO_NONE = 2 1761 } MPChacoLocalType; 1762 PETSC_EXTERN const char *const MPChacoLocalTypes[]; 1763 typedef enum { 1764 MP_CHACO_LANCZOS = 0, 1765 MP_CHACO_RQI = 1 1766 } MPChacoEigenType; 1767 PETSC_EXTERN const char *const MPChacoEigenTypes[]; 1768 1769 PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetGlobal(MatPartitioning, MPChacoGlobalType); 1770 PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetGlobal(MatPartitioning, MPChacoGlobalType *); 1771 PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetLocal(MatPartitioning, MPChacoLocalType); 1772 PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetLocal(MatPartitioning, MPChacoLocalType *); 1773 PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetCoarseLevel(MatPartitioning, PetscReal); 1774 PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenSolver(MatPartitioning, MPChacoEigenType); 1775 PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenSolver(MatPartitioning, MPChacoEigenType *); 1776 PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenTol(MatPartitioning, PetscReal); 1777 PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenTol(MatPartitioning, PetscReal *); 1778 PETSC_EXTERN PetscErrorCode MatPartitioningChacoSetEigenNumber(MatPartitioning, PetscInt); 1779 PETSC_EXTERN PetscErrorCode MatPartitioningChacoGetEigenNumber(MatPartitioning, PetscInt *); 1780 1781 #define MP_PARTY_OPT "opt" 1782 #define MP_PARTY_LIN "lin" 1783 #define MP_PARTY_SCA "sca" 1784 #define MP_PARTY_RAN "ran" 1785 #define MP_PARTY_GBF "gbf" 1786 #define MP_PARTY_GCF "gcf" 1787 #define MP_PARTY_BUB "bub" 1788 #define MP_PARTY_DEF "def" 1789 PETSC_EXTERN PetscErrorCode MatPartitioningPartySetGlobal(MatPartitioning, const char *); 1790 #define MP_PARTY_HELPFUL_SETS "hs" 1791 #define MP_PARTY_KERNIGHAN_LIN "kl" 1792 #define MP_PARTY_NONE "no" 1793 PETSC_EXTERN PetscErrorCode MatPartitioningPartySetLocal(MatPartitioning, const char *); 1794 PETSC_EXTERN PetscErrorCode MatPartitioningPartySetCoarseLevel(MatPartitioning, PetscReal); 1795 PETSC_EXTERN PetscErrorCode MatPartitioningPartySetBipart(MatPartitioning, PetscBool); 1796 PETSC_EXTERN PetscErrorCode MatPartitioningPartySetMatchOptimization(MatPartitioning, PetscBool); 1797 1798 typedef enum { 1799 MP_PTSCOTCH_DEFAULT, 1800 MP_PTSCOTCH_QUALITY, 1801 MP_PTSCOTCH_SPEED, 1802 MP_PTSCOTCH_BALANCE, 1803 MP_PTSCOTCH_SAFETY, 1804 MP_PTSCOTCH_SCALABILITY 1805 } MPPTScotchStrategyType; 1806 PETSC_EXTERN const char *const MPPTScotchStrategyTypes[]; 1807 1808 PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetImbalance(MatPartitioning, PetscReal); 1809 PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetImbalance(MatPartitioning, PetscReal *); 1810 PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchSetStrategy(MatPartitioning, MPPTScotchStrategyType); 1811 PETSC_EXTERN PetscErrorCode MatPartitioningPTScotchGetStrategy(MatPartitioning, MPPTScotchStrategyType *); 1812 1813 /* 1814 * hierarchical partitioning 1815 */ 1816 PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalGetFineparts(MatPartitioning, IS *); 1817 PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalGetCoarseparts(MatPartitioning, IS *); 1818 PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalSetNcoarseparts(MatPartitioning, PetscInt); 1819 PETSC_EXTERN PetscErrorCode MatPartitioningHierarchicalSetNfineparts(MatPartitioning, PetscInt); 1820 1821 PETSC_EXTERN PetscErrorCode MatMeshToCellGraph(Mat, PetscInt, Mat *); 1822 1823 /* 1824 If you add entries here you must also add them to src/mat/f90-mod/petscmat.h 1825 If any of the enum values are changed, also update dMatOps dict at src/binding/petsc4py/src/libpetsc4py/libpetsc4py.pyx 1826 */ 1827 typedef enum { 1828 MATOP_SET_VALUES = 0, 1829 MATOP_GET_ROW = 1, 1830 MATOP_RESTORE_ROW = 2, 1831 MATOP_MULT = 3, 1832 MATOP_MULT_ADD = 4, 1833 MATOP_MULT_TRANSPOSE = 5, 1834 MATOP_MULT_TRANSPOSE_ADD = 6, 1835 MATOP_SOLVE = 7, 1836 MATOP_SOLVE_ADD = 8, 1837 MATOP_SOLVE_TRANSPOSE = 9, 1838 MATOP_SOLVE_TRANSPOSE_ADD = 10, 1839 MATOP_LUFACTOR = 11, 1840 MATOP_CHOLESKYFACTOR = 12, 1841 MATOP_SOR = 13, 1842 MATOP_TRANSPOSE = 14, 1843 MATOP_GETINFO = 15, 1844 MATOP_EQUAL = 16, 1845 MATOP_GET_DIAGONAL = 17, 1846 MATOP_DIAGONAL_SCALE = 18, 1847 MATOP_NORM = 19, 1848 MATOP_ASSEMBLY_BEGIN = 20, 1849 MATOP_ASSEMBLY_END = 21, 1850 MATOP_SET_OPTION = 22, 1851 MATOP_ZERO_ENTRIES = 23, 1852 MATOP_ZERO_ROWS = 24, 1853 MATOP_LUFACTOR_SYMBOLIC = 25, 1854 MATOP_LUFACTOR_NUMERIC = 26, 1855 MATOP_CHOLESKY_FACTOR_SYMBOLIC = 27, 1856 MATOP_CHOLESKY_FACTOR_NUMERIC = 28, 1857 MATOP_SETUP = 29, 1858 MATOP_ILUFACTOR_SYMBOLIC = 30, 1859 MATOP_ICCFACTOR_SYMBOLIC = 31, 1860 MATOP_GET_DIAGONAL_BLOCK = 32, 1861 MATOP_SET_INF = 33, 1862 MATOP_DUPLICATE = 34, 1863 MATOP_FORWARD_SOLVE = 35, 1864 MATOP_BACKWARD_SOLVE = 36, 1865 MATOP_ILUFACTOR = 37, 1866 MATOP_ICCFACTOR = 38, 1867 MATOP_AXPY = 39, 1868 MATOP_CREATE_SUBMATRICES = 40, 1869 MATOP_INCREASE_OVERLAP = 41, 1870 MATOP_GET_VALUES = 42, 1871 MATOP_COPY = 43, 1872 MATOP_GET_ROW_MAX = 44, 1873 MATOP_SCALE = 45, 1874 MATOP_SHIFT = 46, 1875 MATOP_DIAGONAL_SET = 47, 1876 MATOP_ZERO_ROWS_COLUMNS = 48, 1877 MATOP_SET_RANDOM = 49, 1878 MATOP_GET_ROW_IJ = 50, 1879 MATOP_RESTORE_ROW_IJ = 51, 1880 MATOP_GET_COLUMN_IJ = 52, 1881 MATOP_RESTORE_COLUMN_IJ = 53, 1882 MATOP_FDCOLORING_CREATE = 54, 1883 MATOP_COLORING_PATCH = 55, 1884 MATOP_SET_UNFACTORED = 56, 1885 MATOP_PERMUTE = 57, 1886 MATOP_SET_VALUES_BLOCKED = 58, 1887 MATOP_CREATE_SUBMATRIX = 59, 1888 MATOP_DESTROY = 60, 1889 MATOP_VIEW = 61, 1890 MATOP_CONVERT_FROM = 62, 1891 /* MATOP_PLACEHOLDER_63=63 */ 1892 MATOP_MATMAT_MULT_SYMBOLIC = 64, 1893 MATOP_MATMAT_MULT_NUMERIC = 65, 1894 MATOP_SET_LOCAL_TO_GLOBAL_MAP = 66, 1895 MATOP_SET_VALUES_LOCAL = 67, 1896 MATOP_ZERO_ROWS_LOCAL = 68, 1897 MATOP_GET_ROW_MAX_ABS = 69, 1898 MATOP_GET_ROW_MIN_ABS = 70, 1899 MATOP_CONVERT = 71, 1900 MATOP_HAS_OPERATION = 72, 1901 /* MATOP_PLACEHOLDER_73=73, */ 1902 MATOP_SET_VALUES_ADIFOR = 74, 1903 MATOP_FD_COLORING_APPLY = 75, 1904 MATOP_SET_FROM_OPTIONS = 76, 1905 /* MATOP_PLACEHOLDER_77=77, */ 1906 /* MATOP_PLACEHOLDER_78=78, */ 1907 MATOP_FIND_ZERO_DIAGONALS = 79, 1908 MATOP_MULT_MULTIPLE = 80, 1909 MATOP_SOLVE_MULTIPLE = 81, 1910 MATOP_GET_INERTIA = 82, 1911 MATOP_LOAD = 83, 1912 MATOP_IS_SYMMETRIC = 84, 1913 MATOP_IS_HERMITIAN = 85, 1914 MATOP_IS_STRUCTURALLY_SYMMETRIC = 86, 1915 MATOP_SET_VALUES_BLOCKEDLOCAL = 87, 1916 MATOP_CREATE_VECS = 88, 1917 /* MATOP_PLACEHOLDER_89=89, */ 1918 MATOP_MAT_MULT_SYMBOLIC = 90, 1919 MATOP_MAT_MULT_NUMERIC = 91, 1920 /* MATOP_PLACEHOLDER_92=92, */ 1921 MATOP_PTAP_SYMBOLIC = 93, 1922 MATOP_PTAP_NUMERIC = 94, 1923 /* MATOP_PLACEHOLDER_95=95, */ 1924 MATOP_MAT_TRANSPOSE_MULT_SYMBO = 96, 1925 MATOP_MAT_TRANSPOSE_MULT_NUMER = 97, 1926 MATOP_BIND_TO_CPU = 98, 1927 MATOP_PRODUCTSETFROMOPTIONS = 99, 1928 MATOP_PRODUCTSYMBOLIC = 100, 1929 MATOP_PRODUCTNUMERIC = 101, 1930 MATOP_CONJUGATE = 102, 1931 MATOP_VIEW_NATIVE = 103, 1932 MATOP_SET_VALUES_ROW = 104, 1933 MATOP_REAL_PART = 105, 1934 MATOP_IMAGINARY_PART = 106, 1935 MATOP_GET_ROW_UPPER_TRIANGULAR = 107, 1936 MATOP_RESTORE_ROW_UPPER_TRIANG = 108, 1937 MATOP_MAT_SOLVE = 109, 1938 MATOP_MAT_SOLVE_TRANSPOSE = 110, 1939 MATOP_GET_ROW_MIN = 111, 1940 MATOP_GET_COLUMN_VECTOR = 112, 1941 MATOP_MISSING_DIAGONAL = 113, 1942 MATOP_GET_SEQ_NONZERO_STRUCTUR = 114, 1943 MATOP_CREATE = 115, 1944 MATOP_GET_GHOSTS = 116, 1945 MATOP_GET_LOCAL_SUB_MATRIX = 117, 1946 MATOP_RESTORE_LOCALSUB_MATRIX = 118, 1947 MATOP_MULT_DIAGONAL_BLOCK = 119, 1948 MATOP_HERMITIAN_TRANSPOSE = 120, 1949 MATOP_MULT_HERMITIAN_TRANSPOSE = 121, 1950 MATOP_MULT_HERMITIAN_TRANS_ADD = 122, 1951 MATOP_GET_MULTI_PROC_BLOCK = 123, 1952 MATOP_FIND_NONZERO_ROWS = 124, 1953 MATOP_GET_COLUMN_NORMS = 125, 1954 MATOP_INVERT_BLOCK_DIAGONAL = 126, 1955 MATOP_INVERT_VBLOCK_DIAGONAL = 127, 1956 MATOP_CREATE_SUB_MATRICES_MPI = 128, 1957 MATOP_SET_VALUES_BATCH = 129, 1958 /* MATOP_PLACEHOLDER_130=130, */ 1959 MATOP_TRANSPOSE_MAT_MULT_SYMBO = 131, 1960 MATOP_TRANSPOSE_MAT_MULT_NUMER = 132, 1961 MATOP_TRANSPOSE_COLORING_CREAT = 133, 1962 MATOP_TRANS_COLORING_APPLY_SPT = 134, 1963 MATOP_TRANS_COLORING_APPLY_DEN = 135, 1964 /* MATOP_PLACEHOLDER_136=136, */ 1965 MATOP_RART_SYMBOLIC = 137, 1966 MATOP_RART_NUMERIC = 138, 1967 MATOP_SET_BLOCK_SIZES = 139, 1968 MATOP_AYPX = 140, 1969 MATOP_RESIDUAL = 141, 1970 MATOP_FDCOLORING_SETUP = 142, 1971 MATOP_FIND_OFFBLOCK_ENTRIES = 143, 1972 MATOP_MPICONCATENATESEQ = 144, 1973 MATOP_DESTROYSUBMATRICES = 145, 1974 MATOP_TRANSPOSE_SOLVE = 146, 1975 MATOP_GET_VALUES_LOCAL = 147 1976 } MatOperation; 1977 PETSC_EXTERN PetscErrorCode MatSetOperation(Mat, MatOperation, void (*)(void)); 1978 PETSC_EXTERN PetscErrorCode MatGetOperation(Mat, MatOperation, void (**)(void)); 1979 PETSC_EXTERN PetscErrorCode MatHasOperation(Mat, MatOperation, PetscBool *); 1980 PETSC_EXTERN PetscErrorCode MatHasCongruentLayouts(Mat, PetscBool *); 1981 PETSC_DEPRECATED_FUNCTION("Use MatProductClear() (since version 3.14)") static inline PetscErrorCode MatFreeIntermediateDataStructures(Mat A) 1982 { 1983 return MatProductClear(A); 1984 } 1985 PETSC_EXTERN PetscErrorCode MatShellSetOperation(Mat, MatOperation, void (*)(void)); 1986 PETSC_EXTERN PetscErrorCode MatShellGetOperation(Mat, MatOperation, void (**)(void)); 1987 PETSC_EXTERN PetscErrorCode MatShellSetContext(Mat, void *); 1988 PETSC_EXTERN PetscErrorCode MatShellSetContextDestroy(Mat, PetscErrorCode (*)(void *)); 1989 PETSC_EXTERN PetscErrorCode MatShellSetVecType(Mat, VecType); 1990 PETSC_EXTERN PetscErrorCode MatShellTestMult(Mat, PetscErrorCode (*)(void *, Vec, Vec), Vec, void *, PetscBool *); 1991 PETSC_EXTERN PetscErrorCode MatShellTestMultTranspose(Mat, PetscErrorCode (*)(void *, Vec, Vec), Vec, void *, PetscBool *); 1992 PETSC_EXTERN PetscErrorCode MatShellSetManageScalingShifts(Mat); 1993 PETSC_EXTERN PetscErrorCode MatShellSetMatProductOperation(Mat, MatProductType, PetscErrorCode (*)(Mat, Mat, Mat, void **), PetscErrorCode (*)(Mat, Mat, Mat, void *), PetscErrorCode (*)(void *), MatType, MatType); 1994 PETSC_EXTERN PetscErrorCode MatIsShell(Mat, PetscBool *); 1995 1996 /* 1997 Codes for matrices stored on disk. By default they are 1998 stored in a universal format. By changing the format with 1999 PetscViewerPushFormat(viewer,PETSC_VIEWER_NATIVE); the matrices will 2000 be stored in a way natural for the matrix, for example dense matrices 2001 would be stored as dense. Matrices stored this way may only be 2002 read into matrices of the same type. 2003 */ 2004 #define MATRIX_BINARY_FORMAT_DENSE -1 2005 2006 PETSC_EXTERN PetscErrorCode MatMPIBAIJSetHashTableFactor(Mat, PetscReal); 2007 2008 PETSC_EXTERN PetscErrorCode MatISSetLocalMatType(Mat, MatType); 2009 PETSC_EXTERN PetscErrorCode MatISSetPreallocation(Mat, PetscInt, const PetscInt[], PetscInt, const PetscInt[]); 2010 PETSC_EXTERN PetscErrorCode MatISStoreL2L(Mat, PetscBool); 2011 PETSC_EXTERN PetscErrorCode MatISFixLocalEmpty(Mat, PetscBool); 2012 PETSC_EXTERN PetscErrorCode MatISGetLocalMat(Mat, Mat *); 2013 PETSC_EXTERN PetscErrorCode MatISRestoreLocalMat(Mat, Mat *); 2014 PETSC_EXTERN PetscErrorCode MatISSetLocalMat(Mat, Mat); 2015 PETSC_EXTERN PetscErrorCode MatISGetLocalToGlobalMapping(Mat, ISLocalToGlobalMapping *, ISLocalToGlobalMapping *); 2016 PETSC_EXTERN PETSC_DEPRECATED_FUNCTION("Use the MatConvert() interface (since version 3.10)") PetscErrorCode MatISGetMPIXAIJ(Mat, MatReuse, Mat *); 2017 2018 /*S 2019 MatNullSpace - Object that removes a null space from a vector, i.e. 2020 orthogonalizes the vector to a subspace 2021 2022 Level: advanced 2023 2024 .seealso: [](chapter_matrices), `Mat`, `MatNullSpaceCreate()`, `MatNullSpaceSetFunction()`, `MatGetNullSpace()`, `MatSetNullSpace()` 2025 S*/ 2026 typedef struct _p_MatNullSpace *MatNullSpace; 2027 2028 PETSC_EXTERN PetscErrorCode MatNullSpaceCreate(MPI_Comm, PetscBool, PetscInt, const Vec[], MatNullSpace *); 2029 PETSC_EXTERN PetscErrorCode MatNullSpaceSetFunction(MatNullSpace, PetscErrorCode (*)(MatNullSpace, Vec, void *), void *); 2030 PETSC_EXTERN PetscErrorCode MatNullSpaceDestroy(MatNullSpace *); 2031 PETSC_EXTERN PetscErrorCode MatNullSpaceRemove(MatNullSpace, Vec); 2032 PETSC_EXTERN PetscErrorCode MatGetNullSpace(Mat, MatNullSpace *); 2033 PETSC_EXTERN PetscErrorCode MatGetTransposeNullSpace(Mat, MatNullSpace *); 2034 PETSC_EXTERN PetscErrorCode MatSetTransposeNullSpace(Mat, MatNullSpace); 2035 PETSC_EXTERN PetscErrorCode MatSetNullSpace(Mat, MatNullSpace); 2036 PETSC_EXTERN PetscErrorCode MatSetNearNullSpace(Mat, MatNullSpace); 2037 PETSC_EXTERN PetscErrorCode MatGetNearNullSpace(Mat, MatNullSpace *); 2038 PETSC_EXTERN PetscErrorCode MatNullSpaceTest(MatNullSpace, Mat, PetscBool *); 2039 PETSC_EXTERN PetscErrorCode MatNullSpaceView(MatNullSpace, PetscViewer); 2040 PETSC_EXTERN PetscErrorCode MatNullSpaceGetVecs(MatNullSpace, PetscBool *, PetscInt *, const Vec **); 2041 PETSC_EXTERN PetscErrorCode MatNullSpaceCreateRigidBody(Vec, MatNullSpace *); 2042 2043 PETSC_EXTERN PetscErrorCode MatReorderingSeqSBAIJ(Mat, IS); 2044 PETSC_EXTERN PetscErrorCode MatMPISBAIJSetHashTableFactor(Mat, PetscReal); 2045 PETSC_EXTERN PetscErrorCode MatSeqSBAIJSetColumnIndices(Mat, PetscInt *); 2046 2047 PETSC_EXTERN PetscErrorCode MatCreateMAIJ(Mat, PetscInt, Mat *); 2048 PETSC_EXTERN PetscErrorCode MatMAIJRedimension(Mat, PetscInt, Mat *); 2049 PETSC_EXTERN PetscErrorCode MatMAIJGetAIJ(Mat, Mat *); 2050 2051 PETSC_EXTERN PetscErrorCode MatComputeOperator(Mat, MatType, Mat *); 2052 PETSC_EXTERN PetscErrorCode MatComputeOperatorTranspose(Mat, MatType, Mat *); 2053 2054 PETSC_DEPRECATED_FUNCTION("Use MatComputeOperator() (since version 3.12)") static inline PetscErrorCode MatComputeExplicitOperator(Mat A, Mat *B) 2055 { 2056 return MatComputeOperator(A, PETSC_NULLPTR, B); 2057 } 2058 PETSC_DEPRECATED_FUNCTION("Use MatComputeOperatorTranspose() (since version 3.12)") static inline PetscErrorCode MatComputeExplicitOperatorTranspose(Mat A, Mat *B) 2059 { 2060 return MatComputeOperatorTranspose(A, PETSC_NULLPTR, B); 2061 } 2062 2063 PETSC_EXTERN PetscErrorCode MatCreateKAIJ(Mat, PetscInt, PetscInt, const PetscScalar[], const PetscScalar[], Mat *); 2064 PETSC_EXTERN PetscErrorCode MatKAIJGetAIJ(Mat, Mat *); 2065 PETSC_EXTERN PetscErrorCode MatKAIJGetS(Mat, PetscInt *, PetscInt *, PetscScalar **); 2066 PETSC_EXTERN PetscErrorCode MatKAIJGetSRead(Mat, PetscInt *, PetscInt *, const PetscScalar **); 2067 PETSC_EXTERN PetscErrorCode MatKAIJRestoreS(Mat, PetscScalar **); 2068 PETSC_EXTERN PetscErrorCode MatKAIJRestoreSRead(Mat, const PetscScalar **); 2069 PETSC_EXTERN PetscErrorCode MatKAIJGetT(Mat, PetscInt *, PetscInt *, PetscScalar **); 2070 PETSC_EXTERN PetscErrorCode MatKAIJGetTRead(Mat, PetscInt *, PetscInt *, const PetscScalar **); 2071 PETSC_EXTERN PetscErrorCode MatKAIJRestoreT(Mat, PetscScalar **); 2072 PETSC_EXTERN PetscErrorCode MatKAIJRestoreTRead(Mat, const PetscScalar **); 2073 PETSC_EXTERN PetscErrorCode MatKAIJSetAIJ(Mat, Mat); 2074 PETSC_EXTERN PetscErrorCode MatKAIJSetS(Mat, PetscInt, PetscInt, const PetscScalar[]); 2075 PETSC_EXTERN PetscErrorCode MatKAIJSetT(Mat, PetscInt, PetscInt, const PetscScalar[]); 2076 PETSC_EXTERN PetscErrorCode MatKAIJGetScaledIdentity(Mat, PetscBool *); 2077 2078 PETSC_EXTERN PetscErrorCode MatDiagonalScaleLocal(Mat, Vec); 2079 2080 PETSC_EXTERN PetscErrorCode MatCreateMFFD(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, Mat *); 2081 PETSC_EXTERN PetscErrorCode MatMFFDSetBase(Mat, Vec, Vec); 2082 PETSC_EXTERN PetscErrorCode MatMFFDSetFunction(Mat, PetscErrorCode (*)(void *, Vec, Vec), void *); 2083 PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioni(Mat, PetscErrorCode (*)(void *, PetscInt, Vec, PetscScalar *)); 2084 PETSC_EXTERN PetscErrorCode MatMFFDSetFunctioniBase(Mat, PetscErrorCode (*)(void *, Vec)); 2085 PETSC_EXTERN PetscErrorCode MatMFFDSetHHistory(Mat, PetscScalar[], PetscInt); 2086 PETSC_EXTERN PetscErrorCode MatMFFDResetHHistory(Mat); 2087 PETSC_EXTERN PetscErrorCode MatMFFDSetFunctionError(Mat, PetscReal); 2088 PETSC_EXTERN PetscErrorCode MatMFFDSetPeriod(Mat, PetscInt); 2089 PETSC_EXTERN PetscErrorCode MatMFFDGetH(Mat, PetscScalar *); 2090 PETSC_EXTERN PetscErrorCode MatMFFDSetOptionsPrefix(Mat, const char[]); 2091 PETSC_EXTERN PetscErrorCode MatMFFDCheckPositivity(void *, Vec, Vec, PetscScalar *); 2092 PETSC_EXTERN PetscErrorCode MatMFFDSetCheckh(Mat, PetscErrorCode (*)(void *, Vec, Vec, PetscScalar *), void *); 2093 2094 /*S 2095 MatMFFD - A data structured used to manage the computation of the h differencing parameter for matrix-free 2096 Jacobian vector products 2097 2098 Level: developer 2099 2100 Notes: 2101 `MATMFFD` is a specific `MatType` which uses the `MatMFFD` data structure 2102 2103 MatMFFD*() methods actually take the `Mat` as their first argument. Not a `MatMFFD` data structure 2104 2105 This functionality is often obtained using `MatCreateSNESMF()` or with `SNES` solvers using `-snes_mf` or `-snes_mf_operator` 2106 2107 .seealso: [](chapter_matrices), `MatMFFDType`, `MATMFFD`, `MatCreateMFFD()`, `MatMFFDSetFuction()`, `MatMFFDSetType()`, `MatMFFDRegister()`, 2108 `MatCreateSNESMF()`, `SNES`, `-snes_mf`, `-snes_mf_operator` 2109 S*/ 2110 typedef struct _p_MatMFFD *MatMFFD; 2111 2112 /*J 2113 MatMFFDType - algorithm used to compute the `h` used in computing matrix-vector products via differencing of a function 2114 2115 Level: beginner 2116 2117 .seealso: [](chapter_matrices), `MatMFFDSetType()`, `MatMFFDRegister()`, `MatMFFDSetFunction()`, `MatCreateMFFD()` 2118 J*/ 2119 typedef const char *MatMFFDType; 2120 #define MATMFFD_DS "ds" 2121 #define MATMFFD_WP "wp" 2122 2123 PETSC_EXTERN PetscErrorCode MatMFFDSetType(Mat, MatMFFDType); 2124 PETSC_EXTERN PetscErrorCode MatMFFDRegister(const char[], PetscErrorCode (*)(MatMFFD)); 2125 2126 PETSC_EXTERN PetscErrorCode MatMFFDDSSetUmin(Mat, PetscReal); 2127 PETSC_EXTERN PetscErrorCode MatMFFDWPSetComputeNormU(Mat, PetscBool); 2128 2129 PETSC_EXTERN PetscErrorCode MatFDColoringSetType(MatFDColoring, MatMFFDType); 2130 2131 PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutMatrix(PetscViewer, PetscInt, PetscInt, PetscReal *); 2132 PETSC_EXTERN PetscErrorCode PetscViewerMathematicaPutCSRMatrix(PetscViewer, PetscInt, PetscInt, PetscInt *, PetscInt *, PetscReal *); 2133 2134 #ifdef PETSC_HAVE_H2OPUS 2135 PETSC_EXTERN_TYPEDEF typedef PetscScalar (*MatH2OpusKernel)(PetscInt, PetscReal[], PetscReal[], void *); 2136 PETSC_EXTERN PetscErrorCode MatCreateH2OpusFromKernel(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, const PetscReal[], PetscBool, MatH2OpusKernel, void *, PetscReal, PetscInt, PetscInt, Mat *); 2137 PETSC_EXTERN PetscErrorCode MatCreateH2OpusFromMat(Mat, PetscInt, const PetscReal[], PetscBool, PetscReal, PetscInt, PetscInt, PetscInt, PetscReal, Mat *); 2138 PETSC_EXTERN PetscErrorCode MatH2OpusSetSamplingMat(Mat, Mat, PetscInt, PetscReal); 2139 PETSC_EXTERN PetscErrorCode MatH2OpusOrthogonalize(Mat); 2140 PETSC_EXTERN PetscErrorCode MatH2OpusCompress(Mat, PetscReal); 2141 PETSC_EXTERN PetscErrorCode MatH2OpusSetNativeMult(Mat, PetscBool); 2142 PETSC_EXTERN PetscErrorCode MatH2OpusGetNativeMult(Mat, PetscBool *); 2143 PETSC_EXTERN PetscErrorCode MatH2OpusGetIndexMap(Mat, IS *); 2144 PETSC_EXTERN PetscErrorCode MatH2OpusMapVec(Mat, PetscBool, Vec, Vec *); 2145 PETSC_EXTERN PetscErrorCode MatH2OpusLowRankUpdate(Mat, Mat, Mat, PetscScalar); 2146 #endif 2147 2148 #ifdef PETSC_HAVE_HTOOL 2149 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode (*MatHtoolKernel)(PetscInt, PetscInt, PetscInt, const PetscInt *, const PetscInt *, PetscScalar *, void *); 2150 PETSC_EXTERN PetscErrorCode MatCreateHtoolFromKernel(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, const PetscReal[], const PetscReal[], MatHtoolKernel, void *, Mat *); 2151 PETSC_EXTERN PetscErrorCode MatHtoolSetKernel(Mat, MatHtoolKernel, void *); 2152 PETSC_EXTERN PetscErrorCode MatHtoolGetPermutationSource(Mat, IS *); 2153 PETSC_EXTERN PetscErrorCode MatHtoolGetPermutationTarget(Mat, IS *); 2154 PETSC_EXTERN PetscErrorCode MatHtoolUsePermutation(Mat, PetscBool); 2155 2156 /*E 2157 MatHtoolCompressorType - Indicates the type of compressor used by a `MATHTOOL` 2158 2159 Level: intermediate 2160 2161 Values: 2162 + `MAT_HTOOL_COMPRESSOR_SYMPARTIAL_ACA` (default) - symmetric partial adaptive cross approximation 2163 . `MAT_HTOOL_COMPRESSOR_FULL_ACA` - full adaptive cross approximation 2164 - `MAT_HTOOL_COMPRESSOR_SVD` - singular value decomposition 2165 2166 .seealso: [](chapter_matrices), `Mat`, `MatCreateHtoolFromKernel()`, `MATHTOOL`, `MatHtoolClusteringType` 2167 E*/ 2168 typedef enum { 2169 MAT_HTOOL_COMPRESSOR_SYMPARTIAL_ACA, 2170 MAT_HTOOL_COMPRESSOR_FULL_ACA, 2171 MAT_HTOOL_COMPRESSOR_SVD 2172 } MatHtoolCompressorType; 2173 2174 /*E 2175 MatHtoolClusteringType - Indicates the type of clustering used by a `MATHTOOL` 2176 2177 Level: intermediate 2178 2179 Values: 2180 + `MAT_HTOOL_CLUSTERING_PCA_REGULAR` (default) - axis computed via principle component analysis, split uniformly 2181 . `MAT_HTOOL_CLUSTERING_PCA_GEOMETRIC` - axis computed via principle component analysis, split barycentrically 2182 . `MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_REGULAR` - axis along the largest extent of the bounding box, split uniformly 2183 - `MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_GEOMETRIC` - axis along the largest extent of the bounding box, split barycentrically 2184 2185 Note: 2186 Higher-dimensional clustering is not yet supported in Htool, but once it is, one should add BOUNDING_BOX_{2,3} types 2187 2188 .seealso: [](chapter_matrices), `Mat`, `MatCreateHtoolFromKernel()`, `MATHTOOL`, `MatHtoolCompressorType` 2189 E*/ 2190 typedef enum { 2191 MAT_HTOOL_CLUSTERING_PCA_REGULAR, 2192 MAT_HTOOL_CLUSTERING_PCA_GEOMETRIC, 2193 MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_REGULAR, 2194 MAT_HTOOL_CLUSTERING_BOUNDING_BOX_1_GEOMETRIC 2195 } MatHtoolClusteringType; 2196 #endif 2197 2198 #ifdef PETSC_HAVE_MUMPS 2199 PETSC_EXTERN PetscErrorCode MatMumpsSetIcntl(Mat, PetscInt, PetscInt); 2200 PETSC_EXTERN PetscErrorCode MatMumpsGetIcntl(Mat, PetscInt, PetscInt *); 2201 PETSC_EXTERN PetscErrorCode MatMumpsSetCntl(Mat, PetscInt, PetscReal); 2202 PETSC_EXTERN PetscErrorCode MatMumpsGetCntl(Mat, PetscInt, PetscReal *); 2203 2204 PETSC_EXTERN PetscErrorCode MatMumpsGetInfo(Mat, PetscInt, PetscInt *); 2205 PETSC_EXTERN PetscErrorCode MatMumpsGetInfog(Mat, PetscInt, PetscInt *); 2206 PETSC_EXTERN PetscErrorCode MatMumpsGetRinfo(Mat, PetscInt, PetscReal *); 2207 PETSC_EXTERN PetscErrorCode MatMumpsGetRinfog(Mat, PetscInt, PetscReal *); 2208 PETSC_EXTERN PetscErrorCode MatMumpsGetNullPivots(Mat, PetscInt *, PetscInt **); 2209 PETSC_EXTERN PetscErrorCode MatMumpsGetInverse(Mat, Mat); 2210 PETSC_EXTERN PetscErrorCode MatMumpsGetInverseTranspose(Mat, Mat); 2211 #endif 2212 2213 #ifdef PETSC_HAVE_MKL_PARDISO 2214 PETSC_EXTERN PetscErrorCode MatMkl_PardisoSetCntl(Mat, PetscInt, PetscInt); 2215 #endif 2216 2217 #ifdef PETSC_HAVE_MKL_CPARDISO 2218 PETSC_EXTERN PetscErrorCode MatMkl_CPardisoSetCntl(Mat, PetscInt, PetscInt); 2219 #endif 2220 2221 #ifdef PETSC_HAVE_SUPERLU 2222 PETSC_EXTERN PetscErrorCode MatSuperluSetILUDropTol(Mat, PetscReal); 2223 #endif 2224 2225 #ifdef PETSC_HAVE_SUPERLU_DIST 2226 PETSC_EXTERN PetscErrorCode MatSuperluDistGetDiagU(Mat, PetscScalar *); 2227 #endif 2228 2229 #ifdef PETSC_HAVE_STRUMPACK 2230 /*E 2231 MatSTRUMPACKReordering - sparsity reducing ordering to be used in `MATSOLVERSTRUMPACK`` 2232 2233 Level: intermediate 2234 2235 Values: 2236 + `MAT_STRUMPACK_NATURAL` - use the current ordering 2237 . `MAT_STRUMPACK_METIS` - use MeTis to compute an ordering 2238 . `MAT_STRUMPACK_PARMETIS` - use ParMeTis to compute an ordering 2239 . `MAT_STRUMPACK_SCOTCH` - use Scotch to compute an ordering 2240 . `MAT_STRUMPACK_PTSCOTCH` - use parallel Scotch to compute an ordering 2241 - `MAT_STRUMPACK_RCM` - use an RCM ordering 2242 2243 Developer Note: 2244 Should be called `MatSTRUMPACKReorderingType` 2245 2246 .seealso: `Mat`, `MATSOLVERSTRUMPACK`, `MatGetFactor()`, `MatSTRUMPACKSetHSSRelTol()`, `MatSTRUMPACKSetReordering()` 2247 E*/ 2248 typedef enum { 2249 MAT_STRUMPACK_NATURAL = 0, 2250 MAT_STRUMPACK_METIS = 1, 2251 MAT_STRUMPACK_PARMETIS = 2, 2252 MAT_STRUMPACK_SCOTCH = 3, 2253 MAT_STRUMPACK_PTSCOTCH = 4, 2254 MAT_STRUMPACK_RCM = 5 2255 } MatSTRUMPACKReordering; 2256 2257 PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetReordering(Mat, MatSTRUMPACKReordering); 2258 PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetColPerm(Mat, PetscBool); 2259 PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSRelTol(Mat, PetscReal); 2260 PETSC_DEPRECATED_FUNCTION("Use MatSTRUMPACKSetHSSRelTol() (since version 3.9)") static inline PetscErrorCode MatSTRUMPACKSetHSSRelCompTol(Mat mat, PetscReal rtol) 2261 { 2262 return MatSTRUMPACKSetHSSRelTol(mat, rtol); 2263 } 2264 PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSAbsTol(Mat, PetscReal); 2265 PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSMinSepSize(Mat, PetscInt); 2266 PETSC_DEPRECATED_FUNCTION("Use MatSTRUMPACKSetHSSMinSepSize() (since version 3.9)") static inline PetscErrorCode MatSTRUMPACKSetHSSMinSize(Mat mat, PetscInt hssminsize) 2267 { 2268 return MatSTRUMPACKSetHSSMinSepSize(mat, hssminsize); 2269 } 2270 PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSMaxRank(Mat, PetscInt); 2271 PETSC_EXTERN PetscErrorCode MatSTRUMPACKSetHSSLeafSize(Mat, PetscInt); 2272 #endif 2273 2274 PETSC_EXTERN PetscErrorCode MatBindToCPU(Mat, PetscBool); 2275 PETSC_EXTERN PetscErrorCode MatBoundToCPU(Mat, PetscBool *); 2276 PETSC_DEPRECATED_FUNCTION("Use MatBindToCPU (since v3.13)") static inline PetscErrorCode MatPinToCPU(Mat A, PetscBool flg) 2277 { 2278 return MatBindToCPU(A, flg); 2279 } 2280 PETSC_EXTERN PetscErrorCode MatSetBindingPropagates(Mat, PetscBool); 2281 PETSC_EXTERN PetscErrorCode MatGetBindingPropagates(Mat, PetscBool *); 2282 2283 typedef struct _n_SplitCSRMat *PetscSplitCSRDataStructure; 2284 PETSC_EXTERN PetscErrorCode MatCUSPARSEGetDeviceMatWrite(Mat, PetscSplitCSRDataStructure *); 2285 2286 #ifdef PETSC_HAVE_KOKKOS_KERNELS 2287 PETSC_EXTERN PetscErrorCode MatKokkosGetDeviceMatWrite(Mat, PetscSplitCSRDataStructure *); 2288 PETSC_EXTERN PetscErrorCode MatSeqAIJKokkosSetDeviceMat(Mat, PetscSplitCSRDataStructure); 2289 PETSC_EXTERN PetscErrorCode MatSeqAIJKokkosGetDeviceMat(Mat, PetscSplitCSRDataStructure *); 2290 #endif 2291 2292 #ifdef PETSC_HAVE_CUDA 2293 /*E 2294 MatCUSPARSEStorageFormat - indicates the storage format for `MATAIJCUSPARSE` (GPU) 2295 matrices. 2296 2297 Values: 2298 + `MAT_CUSPARSE_CSR` - Compressed Sparse Row 2299 . `MAT_CUSPARSE_ELL` - Ellpack (requires CUDA 4.2 or later). 2300 - `MAT_CUSPARSE_HYB` - Hybrid, a combination of Ellpack and Coordinate format (requires CUDA 4.2 or later). 2301 2302 Level: intermediate 2303 2304 .seealso: [](chapter_matrices), `MATAIJCUSPARSE`, `MatCUSPARSESetFormat()`, `MatCUSPARSEFormatOperation` 2305 E*/ 2306 2307 typedef enum { 2308 MAT_CUSPARSE_CSR, 2309 MAT_CUSPARSE_ELL, 2310 MAT_CUSPARSE_HYB 2311 } MatCUSPARSEStorageFormat; 2312 2313 /* these will be strings associated with enumerated type defined above */ 2314 PETSC_EXTERN const char *const MatCUSPARSEStorageFormats[]; 2315 2316 /*E 2317 MatCUSPARSEFormatOperation - indicates the operation of `MATAIJCUSPARSE` (GPU) 2318 matrices whose operation should use a particular storage format. 2319 2320 Values: 2321 + `MAT_CUSPARSE_MULT_DIAG` - sets the storage format for the diagonal matrix in the parallel MatMult 2322 . `MAT_CUSPARSE_MULT_OFFDIAG` - sets the storage format for the offdiagonal matrix in the parallel MatMult 2323 . `MAT_CUSPARSE_MULT` - sets the storage format for the entire matrix in the serial (single GPU) MatMult 2324 - `MAT_CUSPARSE_ALL` - sets the storage format for all `MATAIJCUSPARSE` (GPU) matrices 2325 2326 Level: intermediate 2327 2328 .seealso: [](chapter_matrices), `MatCUSPARSESetFormat()`, `MatCUSPARSEStorageFormat` 2329 E*/ 2330 typedef enum { 2331 MAT_CUSPARSE_MULT_DIAG, 2332 MAT_CUSPARSE_MULT_OFFDIAG, 2333 MAT_CUSPARSE_MULT, 2334 MAT_CUSPARSE_ALL 2335 } MatCUSPARSEFormatOperation; 2336 2337 PETSC_EXTERN PetscErrorCode MatCreateSeqAIJCUSPARSE(MPI_Comm, PetscInt, PetscInt, PetscInt, const PetscInt[], Mat *); 2338 PETSC_EXTERN PetscErrorCode MatCreateAIJCUSPARSE(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, const PetscInt[], PetscInt, const PetscInt[], Mat *); 2339 PETSC_EXTERN PetscErrorCode MatCUSPARSESetFormat(Mat, MatCUSPARSEFormatOperation, MatCUSPARSEStorageFormat); 2340 PETSC_EXTERN PetscErrorCode MatCUSPARSESetUseCPUSolve(Mat, PetscBool); 2341 PETSC_EXTERN PetscErrorCode MatSeqAIJCUSPARSEGetIJ(Mat, PetscBool, const int **, const int **); 2342 PETSC_EXTERN PetscErrorCode MatSeqAIJCUSPARSERestoreIJ(Mat, PetscBool, const int **, const int **); 2343 PETSC_EXTERN PetscErrorCode MatSeqAIJCUSPARSEGetArrayRead(Mat, const PetscScalar **); 2344 PETSC_EXTERN PetscErrorCode MatSeqAIJCUSPARSERestoreArrayRead(Mat, const PetscScalar **); 2345 PETSC_EXTERN PetscErrorCode MatSeqAIJCUSPARSEGetArrayWrite(Mat, PetscScalar **); 2346 PETSC_EXTERN PetscErrorCode MatSeqAIJCUSPARSERestoreArrayWrite(Mat, PetscScalar **); 2347 PETSC_EXTERN PetscErrorCode MatSeqAIJCUSPARSEGetArray(Mat, PetscScalar **); 2348 PETSC_EXTERN PetscErrorCode MatSeqAIJCUSPARSERestoreArray(Mat, PetscScalar **); 2349 2350 PETSC_EXTERN PetscErrorCode MatCreateDenseCUDA(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar[], Mat *); 2351 PETSC_EXTERN PetscErrorCode MatCreateSeqDenseCUDA(MPI_Comm, PetscInt, PetscInt, PetscScalar[], Mat *); 2352 PETSC_EXTERN PetscErrorCode MatDenseCUDAGetArrayWrite(Mat, PetscScalar **); 2353 PETSC_EXTERN PetscErrorCode MatDenseCUDAGetArrayRead(Mat, const PetscScalar **); 2354 PETSC_EXTERN PetscErrorCode MatDenseCUDAGetArray(Mat, PetscScalar **); 2355 PETSC_EXTERN PetscErrorCode MatDenseCUDARestoreArrayWrite(Mat, PetscScalar **); 2356 PETSC_EXTERN PetscErrorCode MatDenseCUDARestoreArrayRead(Mat, const PetscScalar **); 2357 PETSC_EXTERN PetscErrorCode MatDenseCUDARestoreArray(Mat, PetscScalar **); 2358 PETSC_EXTERN PetscErrorCode MatDenseCUDAPlaceArray(Mat, const PetscScalar *); 2359 PETSC_EXTERN PetscErrorCode MatDenseCUDAReplaceArray(Mat, const PetscScalar *); 2360 PETSC_EXTERN PetscErrorCode MatDenseCUDAResetArray(Mat); 2361 #endif 2362 2363 #ifdef PETSC_HAVE_HIP 2364 /*E 2365 MatHIPSPARSEStorageFormat - indicates the storage format for `MATAIJHIPSPARSE` (GPU) 2366 matrices. 2367 2368 Level: intermediate 2369 2370 Values: 2371 + `MAT_HIPSPARSE_CSR` - Compressed Sparse Row 2372 . `MAT_HIPSPARSE_ELL` - Ellpack 2373 - `MAT_HIPSPARSE_HYB` - Hybrid, a combination of Ellpack and Coordinate format 2374 2375 .seealso: [](chapter_matrices), `MatHIPSPARSESetFormat()`, `MatHIPSPARSEFormatOperation` 2376 E*/ 2377 2378 typedef enum { 2379 MAT_HIPSPARSE_CSR, 2380 MAT_HIPSPARSE_ELL, 2381 MAT_HIPSPARSE_HYB 2382 } MatHIPSPARSEStorageFormat; 2383 2384 /* these will be strings associated with enumerated type defined above */ 2385 PETSC_EXTERN const char *const MatHIPSPARSEStorageFormats[]; 2386 2387 /*E 2388 MatHIPSPARSEFormatOperation - indicates the operation of `MATAIJHIPSPARSE` (GPU) 2389 matrices whose operation should use a particular storage format. 2390 2391 Level: intermediate 2392 2393 Values: 2394 + `MAT_HIPSPARSE_MULT_DIAG` - sets the storage format for the diagonal matrix in the parallel `MatMult()` 2395 . `MAT_HIPSPARSE_MULT_OFFDIAG` - sets the storage format for the offdiagonal matrix in the parallel `MatMul()t` 2396 . `MAT_HIPSPARSE_MULT` - sets the storage format for the entire matrix in the serial (single GPU) `MatMult()` 2397 - `MAT_HIPSPARSE_ALL` - sets the storage format for all HIPSPARSE (GPU) matrices 2398 2399 .seealso: [](chapter_matrices), `MatHIPSPARSESetFormat()`, `MatHIPSPARSEStorageFormat` 2400 E*/ 2401 typedef enum { 2402 MAT_HIPSPARSE_MULT_DIAG, 2403 MAT_HIPSPARSE_MULT_OFFDIAG, 2404 MAT_HIPSPARSE_MULT, 2405 MAT_HIPSPARSE_ALL 2406 } MatHIPSPARSEFormatOperation; 2407 2408 PETSC_EXTERN PetscErrorCode MatCreateSeqAIJHIPSPARSE(MPI_Comm, PetscInt, PetscInt, PetscInt, const PetscInt[], Mat *); 2409 PETSC_EXTERN PetscErrorCode MatCreateAIJHIPSPARSE(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, const PetscInt[], PetscInt, const PetscInt[], Mat *); 2410 PETSC_EXTERN PetscErrorCode MatHIPSPARSESetFormat(Mat, MatHIPSPARSEFormatOperation, MatHIPSPARSEStorageFormat); 2411 PETSC_EXTERN PetscErrorCode MatHIPSPARSESetUseCPUSolve(Mat, PetscBool); 2412 PETSC_EXTERN PetscErrorCode MatSeqAIJHIPSPARSEGetIJ(Mat, PetscBool, const int **, const int **); 2413 PETSC_EXTERN PetscErrorCode MatSeqAIJHIPSPARSERestoreIJ(Mat, PetscBool, const int **, const int **); 2414 PETSC_EXTERN PetscErrorCode MatSeqAIJHIPSPARSEGetArrayRead(Mat, const PetscScalar **); 2415 PETSC_EXTERN PetscErrorCode MatSeqAIJHIPSPARSERestoreArrayRead(Mat, const PetscScalar **); 2416 PETSC_EXTERN PetscErrorCode MatSeqAIJHIPSPARSEGetArrayWrite(Mat, PetscScalar **); 2417 PETSC_EXTERN PetscErrorCode MatSeqAIJHIPSPARSERestoreArrayWrite(Mat, PetscScalar **); 2418 PETSC_EXTERN PetscErrorCode MatSeqAIJHIPSPARSEGetArray(Mat, PetscScalar **); 2419 PETSC_EXTERN PetscErrorCode MatSeqAIJHIPSPARSERestoreArray(Mat, PetscScalar **); 2420 2421 PETSC_EXTERN PetscErrorCode MatCreateDenseHIP(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar[], Mat *); 2422 PETSC_EXTERN PetscErrorCode MatCreateSeqDenseHIP(MPI_Comm, PetscInt, PetscInt, PetscScalar[], Mat *); 2423 PETSC_EXTERN PetscErrorCode MatDenseHIPGetArrayWrite(Mat, PetscScalar **); 2424 PETSC_EXTERN PetscErrorCode MatDenseHIPGetArrayRead(Mat, const PetscScalar **); 2425 PETSC_EXTERN PetscErrorCode MatDenseHIPGetArray(Mat, PetscScalar **); 2426 PETSC_EXTERN PetscErrorCode MatDenseHIPRestoreArrayWrite(Mat, PetscScalar **); 2427 PETSC_EXTERN PetscErrorCode MatDenseHIPRestoreArrayRead(Mat, const PetscScalar **); 2428 PETSC_EXTERN PetscErrorCode MatDenseHIPRestoreArray(Mat, PetscScalar **); 2429 PETSC_EXTERN PetscErrorCode MatDenseHIPPlaceArray(Mat, const PetscScalar *); 2430 PETSC_EXTERN PetscErrorCode MatDenseHIPReplaceArray(Mat, const PetscScalar *); 2431 PETSC_EXTERN PetscErrorCode MatDenseHIPResetArray(Mat); 2432 #endif 2433 2434 #if defined(PETSC_HAVE_VIENNACL) 2435 PETSC_EXTERN PetscErrorCode MatCreateSeqAIJViennaCL(MPI_Comm, PetscInt, PetscInt, PetscInt, const PetscInt[], Mat *); 2436 PETSC_EXTERN PetscErrorCode MatCreateAIJViennaCL(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, const PetscInt[], PetscInt, const PetscInt[], Mat *); 2437 #endif 2438 2439 #if defined(PETSC_HAVE_FFTW) 2440 PETSC_EXTERN PetscErrorCode VecScatterPetscToFFTW(Mat, Vec, Vec); 2441 PETSC_EXTERN PetscErrorCode VecScatterFFTWToPetsc(Mat, Vec, Vec); 2442 PETSC_EXTERN PetscErrorCode MatCreateVecsFFTW(Mat, Vec *, Vec *, Vec *); 2443 #endif 2444 2445 #if defined(PETSC_HAVE_SCALAPACK) 2446 PETSC_EXTERN PetscErrorCode MatCreateScaLAPACK(MPI_Comm, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, PetscInt, Mat *); 2447 PETSC_EXTERN PetscErrorCode MatScaLAPACKSetBlockSizes(Mat, PetscInt, PetscInt); 2448 PETSC_EXTERN PetscErrorCode MatScaLAPACKGetBlockSizes(Mat, PetscInt *, PetscInt *); 2449 #endif 2450 2451 PETSC_EXTERN PetscErrorCode MatCreateNest(MPI_Comm, PetscInt, const IS[], PetscInt, const IS[], const Mat[], Mat *); 2452 PETSC_EXTERN PetscErrorCode MatNestGetSize(Mat, PetscInt *, PetscInt *); 2453 PETSC_EXTERN PetscErrorCode MatNestGetISs(Mat, IS[], IS[]); 2454 PETSC_EXTERN PetscErrorCode MatNestGetLocalISs(Mat, IS[], IS[]); 2455 PETSC_EXTERN PetscErrorCode MatNestGetSubMats(Mat, PetscInt *, PetscInt *, Mat ***); 2456 PETSC_EXTERN PetscErrorCode MatNestGetSubMat(Mat, PetscInt, PetscInt, Mat *); 2457 PETSC_EXTERN PetscErrorCode MatNestSetVecType(Mat, VecType); 2458 PETSC_EXTERN PetscErrorCode MatNestSetSubMats(Mat, PetscInt, const IS[], PetscInt, const IS[], const Mat[]); 2459 PETSC_EXTERN PetscErrorCode MatNestSetSubMat(Mat, PetscInt, PetscInt, Mat); 2460 2461 PETSC_EXTERN PetscErrorCode MatChop(Mat, PetscReal); 2462 PETSC_EXTERN PetscErrorCode MatComputeBandwidth(Mat, PetscReal, PetscInt *); 2463 2464 PETSC_EXTERN PetscErrorCode MatSubdomainsCreateCoalesce(Mat, PetscInt, PetscInt *, IS **); 2465 2466 PETSC_EXTERN PetscErrorCode MatPreallocatorPreallocate(Mat, PetscBool, Mat); 2467 2468 PETSC_INTERN PetscErrorCode MatHeaderMerge(Mat, Mat *); 2469 PETSC_EXTERN PetscErrorCode MatHeaderReplace(Mat, Mat *); 2470 2471 PETSC_EXTERN PetscErrorCode MatSeqAIJGetCSRAndMemType(Mat, const PetscInt **, const PetscInt **, PetscScalar **, PetscMemType *); 2472 2473 PETSC_EXTERN PetscErrorCode MatCreateGraph(Mat, PetscBool, PetscBool, PetscReal, Mat *); 2474 PETSC_EXTERN PetscErrorCode MatEliminateZeros(Mat); 2475 2476 PETSC_EXTERN PetscErrorCode VecCreateMatDense(Vec, PetscInt, PetscInt, PetscInt, PetscInt, PetscScalar *, Mat *); 2477 #endif 2478