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