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