1 /* 2 This is where the abstract matrix operations are defined 3 Portions of this code are under: 4 Copyright (c) 2022 Advanced Micro Devices, Inc. All rights reserved. 5 */ 6 7 #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/ 8 #include <petsc/private/isimpl.h> 9 #include <petsc/private/vecimpl.h> 10 11 /* Logging support */ 12 PetscClassId MAT_CLASSID; 13 PetscClassId MAT_COLORING_CLASSID; 14 PetscClassId MAT_FDCOLORING_CLASSID; 15 PetscClassId MAT_TRANSPOSECOLORING_CLASSID; 16 17 PetscLogEvent MAT_Mult, MAT_Mults, MAT_MultAdd, MAT_MultTranspose; 18 PetscLogEvent MAT_MultTransposeAdd, MAT_Solve, MAT_Solves, MAT_SolveAdd, MAT_SolveTranspose, MAT_MatSolve, MAT_MatTrSolve; 19 PetscLogEvent MAT_SolveTransposeAdd, MAT_SOR, MAT_ForwardSolve, MAT_BackwardSolve, MAT_LUFactor, MAT_LUFactorSymbolic; 20 PetscLogEvent MAT_LUFactorNumeric, MAT_CholeskyFactor, MAT_CholeskyFactorSymbolic, MAT_CholeskyFactorNumeric, MAT_ILUFactor; 21 PetscLogEvent MAT_ILUFactorSymbolic, MAT_ICCFactorSymbolic, MAT_Copy, MAT_Convert, MAT_Scale, MAT_AssemblyBegin; 22 PetscLogEvent MAT_QRFactorNumeric, MAT_QRFactorSymbolic, MAT_QRFactor; 23 PetscLogEvent MAT_AssemblyEnd, MAT_SetValues, MAT_GetValues, MAT_GetRow, MAT_GetRowIJ, MAT_CreateSubMats, MAT_GetOrdering, MAT_RedundantMat, MAT_GetSeqNonzeroStructure; 24 PetscLogEvent MAT_IncreaseOverlap, MAT_Partitioning, MAT_PartitioningND, MAT_Coarsen, MAT_ZeroEntries, MAT_Load, MAT_View, MAT_AXPY, MAT_FDColoringCreate; 25 PetscLogEvent MAT_FDColoringSetUp, MAT_FDColoringApply, MAT_Transpose, MAT_FDColoringFunction, MAT_CreateSubMat; 26 PetscLogEvent MAT_TransposeColoringCreate; 27 PetscLogEvent MAT_MatMult, MAT_MatMultSymbolic, MAT_MatMultNumeric; 28 PetscLogEvent MAT_PtAP, MAT_PtAPSymbolic, MAT_PtAPNumeric, MAT_RARt, MAT_RARtSymbolic, MAT_RARtNumeric; 29 PetscLogEvent MAT_MatTransposeMult, MAT_MatTransposeMultSymbolic, MAT_MatTransposeMultNumeric; 30 PetscLogEvent MAT_TransposeMatMult, MAT_TransposeMatMultSymbolic, MAT_TransposeMatMultNumeric; 31 PetscLogEvent MAT_MatMatMult, MAT_MatMatMultSymbolic, MAT_MatMatMultNumeric; 32 PetscLogEvent MAT_MultHermitianTranspose, MAT_MultHermitianTransposeAdd; 33 PetscLogEvent MAT_Getsymtranspose, MAT_Getsymtransreduced, MAT_GetBrowsOfAcols; 34 PetscLogEvent MAT_GetBrowsOfAocols, MAT_Getlocalmat, MAT_Getlocalmatcondensed, MAT_Seqstompi, MAT_Seqstompinum, MAT_Seqstompisym; 35 PetscLogEvent MAT_Applypapt, MAT_Applypapt_numeric, MAT_Applypapt_symbolic, MAT_GetSequentialNonzeroStructure; 36 PetscLogEvent MAT_GetMultiProcBlock; 37 PetscLogEvent MAT_CUSPARSECopyToGPU, MAT_CUSPARSECopyFromGPU, MAT_CUSPARSEGenerateTranspose, MAT_CUSPARSESolveAnalysis; 38 PetscLogEvent MAT_HIPSPARSECopyToGPU, MAT_HIPSPARSECopyFromGPU, MAT_HIPSPARSEGenerateTranspose, MAT_HIPSPARSESolveAnalysis; 39 PetscLogEvent MAT_PreallCOO, MAT_SetVCOO; 40 PetscLogEvent MAT_SetValuesBatch; 41 PetscLogEvent MAT_ViennaCLCopyToGPU; 42 PetscLogEvent MAT_DenseCopyToGPU, MAT_DenseCopyFromGPU; 43 PetscLogEvent MAT_Merge, MAT_Residual, MAT_SetRandom; 44 PetscLogEvent MAT_FactorFactS, MAT_FactorInvS; 45 PetscLogEvent MATCOLORING_Apply, MATCOLORING_Comm, MATCOLORING_Local, MATCOLORING_ISCreate, MATCOLORING_SetUp, MATCOLORING_Weights; 46 PetscLogEvent MAT_H2Opus_Build, MAT_H2Opus_Compress, MAT_H2Opus_Orthog, MAT_H2Opus_LR; 47 48 const char *const MatFactorTypes[] = {"NONE", "LU", "CHOLESKY", "ILU", "ICC", "ILUDT", "QR", "MatFactorType", "MAT_FACTOR_", NULL}; 49 50 /*@ 51 MatSetRandom - Sets all components of a matrix to random numbers. 52 53 Logically Collective 54 55 Input Parameters: 56 + x - the matrix 57 - rctx - the `PetscRandom` object, formed by `PetscRandomCreate()`, or `NULL` and 58 it will create one internally. 59 60 Output Parameter: 61 . x - the matrix 62 63 Example of Usage: 64 .vb 65 PetscRandomCreate(PETSC_COMM_WORLD,&rctx); 66 MatSetRandom(x,rctx); 67 PetscRandomDestroy(rctx); 68 .ve 69 70 Level: intermediate 71 72 Note: 73 For sparse matrices that have been preallocated but not been assembled it randomly selects appropriate locations, 74 for sparse matrices that already have locations it fills the locations with random numbers. It generates an error if used on sparse matrices that have 75 not been preallocated. 76 77 .seealso: `Mat`, `PetscRandom`, `PetscRandomCreate()`, `MatZeroEntries()`, `MatSetValues()`, `PetscRandomCreate()`, `PetscRandomDestroy()` 78 @*/ 79 PetscErrorCode MatSetRandom(Mat x, PetscRandom rctx) 80 { 81 PetscRandom randObj = NULL; 82 83 PetscFunctionBegin; 84 PetscValidHeaderSpecific(x, MAT_CLASSID, 1); 85 if (rctx) PetscValidHeaderSpecific(rctx, PETSC_RANDOM_CLASSID, 2); 86 PetscValidType(x, 1); 87 MatCheckPreallocated(x, 1); 88 89 if (!rctx) { 90 MPI_Comm comm; 91 PetscCall(PetscObjectGetComm((PetscObject)x, &comm)); 92 PetscCall(PetscRandomCreate(comm, &randObj)); 93 PetscCall(PetscRandomSetType(randObj, x->defaultrandtype)); 94 PetscCall(PetscRandomSetFromOptions(randObj)); 95 rctx = randObj; 96 } 97 PetscCall(PetscLogEventBegin(MAT_SetRandom, x, rctx, 0, 0)); 98 PetscUseTypeMethod(x, setrandom, rctx); 99 PetscCall(PetscLogEventEnd(MAT_SetRandom, x, rctx, 0, 0)); 100 101 PetscCall(MatAssemblyBegin(x, MAT_FINAL_ASSEMBLY)); 102 PetscCall(MatAssemblyEnd(x, MAT_FINAL_ASSEMBLY)); 103 PetscCall(PetscRandomDestroy(&randObj)); 104 PetscFunctionReturn(PETSC_SUCCESS); 105 } 106 107 /*@ 108 MatFactorGetErrorZeroPivot - returns the pivot value that was determined to be zero and the row it occurred in 109 110 Logically Collective 111 112 Input Parameter: 113 . mat - the factored matrix 114 115 Output Parameters: 116 + pivot - the pivot value computed 117 - row - the row that the zero pivot occurred. Note that this row must be interpreted carefully due to row reorderings and which processes 118 the share the matrix 119 120 Level: advanced 121 122 Notes: 123 This routine does not work for factorizations done with external packages. 124 125 This routine should only be called if `MatGetFactorError()` returns a value of `MAT_FACTOR_NUMERIC_ZEROPIVOT` 126 127 This can be called on non-factored matrices that come from, for example, matrices used in SOR. 128 129 .seealso: `Mat`, `MatZeroEntries()`, `MatFactor()`, `MatGetFactor()`, `MatLUFactorSymbolic()`, `MatCholeskyFactorSymbolic()`, `MatFactorClearError()`, `MatFactorGetErrorZeroPivot()`, 130 `MAT_FACTOR_NUMERIC_ZEROPIVOT` 131 @*/ 132 PetscErrorCode MatFactorGetErrorZeroPivot(Mat mat, PetscReal *pivot, PetscInt *row) 133 { 134 PetscFunctionBegin; 135 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 136 PetscValidRealPointer(pivot, 2); 137 PetscValidIntPointer(row, 3); 138 *pivot = mat->factorerror_zeropivot_value; 139 *row = mat->factorerror_zeropivot_row; 140 PetscFunctionReturn(PETSC_SUCCESS); 141 } 142 143 /*@ 144 MatFactorGetError - gets the error code from a factorization 145 146 Logically Collective 147 148 Input Parameters: 149 . mat - the factored matrix 150 151 Output Parameter: 152 . err - the error code 153 154 Level: advanced 155 156 Note: 157 This can also be called on non-factored matrices that come from, for example, matrices used in SOR. 158 159 .seealso: `Mat`, `MatZeroEntries()`, `MatFactor()`, `MatGetFactor()`, `MatLUFactorSymbolic()`, `MatCholeskyFactorSymbolic()`, `MatFactorClearError()`, `MatFactorGetErrorZeroPivot()`, 160 `MatFactorError` 161 @*/ 162 PetscErrorCode MatFactorGetError(Mat mat, MatFactorError *err) 163 { 164 PetscFunctionBegin; 165 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 166 PetscValidPointer(err, 2); 167 *err = mat->factorerrortype; 168 PetscFunctionReturn(PETSC_SUCCESS); 169 } 170 171 /*@ 172 MatFactorClearError - clears the error code in a factorization 173 174 Logically Collective 175 176 Input Parameter: 177 . mat - the factored matrix 178 179 Level: developer 180 181 Note: 182 This can also be called on non-factored matrices that come from, for example, matrices used in SOR. 183 184 .seealso: `Mat`, `MatZeroEntries()`, `MatFactor()`, `MatGetFactor()`, `MatLUFactorSymbolic()`, `MatCholeskyFactorSymbolic()`, `MatFactorGetError()`, `MatFactorGetErrorZeroPivot()`, 185 `MatGetErrorCode()`, `MatFactorError` 186 @*/ 187 PetscErrorCode MatFactorClearError(Mat mat) 188 { 189 PetscFunctionBegin; 190 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 191 mat->factorerrortype = MAT_FACTOR_NOERROR; 192 mat->factorerror_zeropivot_value = 0.0; 193 mat->factorerror_zeropivot_row = 0; 194 PetscFunctionReturn(PETSC_SUCCESS); 195 } 196 197 PETSC_INTERN PetscErrorCode MatFindNonzeroRowsOrCols_Basic(Mat mat, PetscBool cols, PetscReal tol, IS *nonzero) 198 { 199 Vec r, l; 200 const PetscScalar *al; 201 PetscInt i, nz, gnz, N, n; 202 203 PetscFunctionBegin; 204 PetscCall(MatCreateVecs(mat, &r, &l)); 205 if (!cols) { /* nonzero rows */ 206 PetscCall(MatGetSize(mat, &N, NULL)); 207 PetscCall(MatGetLocalSize(mat, &n, NULL)); 208 PetscCall(VecSet(l, 0.0)); 209 PetscCall(VecSetRandom(r, NULL)); 210 PetscCall(MatMult(mat, r, l)); 211 PetscCall(VecGetArrayRead(l, &al)); 212 } else { /* nonzero columns */ 213 PetscCall(MatGetSize(mat, NULL, &N)); 214 PetscCall(MatGetLocalSize(mat, NULL, &n)); 215 PetscCall(VecSet(r, 0.0)); 216 PetscCall(VecSetRandom(l, NULL)); 217 PetscCall(MatMultTranspose(mat, l, r)); 218 PetscCall(VecGetArrayRead(r, &al)); 219 } 220 if (tol <= 0.0) { 221 for (i = 0, nz = 0; i < n; i++) 222 if (al[i] != 0.0) nz++; 223 } else { 224 for (i = 0, nz = 0; i < n; i++) 225 if (PetscAbsScalar(al[i]) > tol) nz++; 226 } 227 PetscCall(MPIU_Allreduce(&nz, &gnz, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)mat))); 228 if (gnz != N) { 229 PetscInt *nzr; 230 PetscCall(PetscMalloc1(nz, &nzr)); 231 if (nz) { 232 if (tol < 0) { 233 for (i = 0, nz = 0; i < n; i++) 234 if (al[i] != 0.0) nzr[nz++] = i; 235 } else { 236 for (i = 0, nz = 0; i < n; i++) 237 if (PetscAbsScalar(al[i]) > tol) nzr[nz++] = i; 238 } 239 } 240 PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)mat), nz, nzr, PETSC_OWN_POINTER, nonzero)); 241 } else *nonzero = NULL; 242 if (!cols) { /* nonzero rows */ 243 PetscCall(VecRestoreArrayRead(l, &al)); 244 } else { 245 PetscCall(VecRestoreArrayRead(r, &al)); 246 } 247 PetscCall(VecDestroy(&l)); 248 PetscCall(VecDestroy(&r)); 249 PetscFunctionReturn(PETSC_SUCCESS); 250 } 251 252 /*@ 253 MatFindNonzeroRows - Locate all rows that are not completely zero in the matrix 254 255 Input Parameter: 256 . A - the matrix 257 258 Output Parameter: 259 . keptrows - the rows that are not completely zero 260 261 Level: intermediate 262 263 Note: 264 `keptrows` is set to `NULL` if all rows are nonzero. 265 266 .seealso: `Mat`, `MatFindZeroRows()` 267 @*/ 268 PetscErrorCode MatFindNonzeroRows(Mat mat, IS *keptrows) 269 { 270 PetscFunctionBegin; 271 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 272 PetscValidType(mat, 1); 273 PetscValidPointer(keptrows, 2); 274 PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 275 PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 276 if (mat->ops->findnonzerorows) PetscUseTypeMethod(mat, findnonzerorows, keptrows); 277 else PetscCall(MatFindNonzeroRowsOrCols_Basic(mat, PETSC_FALSE, 0.0, keptrows)); 278 PetscFunctionReturn(PETSC_SUCCESS); 279 } 280 281 /*@ 282 MatFindZeroRows - Locate all rows that are completely zero in the matrix 283 284 Input Parameter: 285 . A - the matrix 286 287 Output Parameter: 288 . zerorows - the rows that are completely zero 289 290 Level: intermediate 291 292 Note: 293 `zerorows` is set to `NULL` if no rows are zero. 294 295 .seealso: `Mat`, `MatFindNonzeroRows()` 296 @*/ 297 PetscErrorCode MatFindZeroRows(Mat mat, IS *zerorows) 298 { 299 IS keptrows; 300 PetscInt m, n; 301 302 PetscFunctionBegin; 303 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 304 PetscValidType(mat, 1); 305 PetscValidPointer(zerorows, 2); 306 PetscCall(MatFindNonzeroRows(mat, &keptrows)); 307 /* MatFindNonzeroRows sets keptrows to NULL if there are no zero rows. 308 In keeping with this convention, we set zerorows to NULL if there are no zero 309 rows. */ 310 if (keptrows == NULL) { 311 *zerorows = NULL; 312 } else { 313 PetscCall(MatGetOwnershipRange(mat, &m, &n)); 314 PetscCall(ISComplement(keptrows, m, n, zerorows)); 315 PetscCall(ISDestroy(&keptrows)); 316 } 317 PetscFunctionReturn(PETSC_SUCCESS); 318 } 319 320 /*@ 321 MatGetDiagonalBlock - Returns the part of the matrix associated with the on-process coupling 322 323 Not Collective 324 325 Input Parameters: 326 . A - the matrix 327 328 Output Parameters: 329 . a - the diagonal part (which is a SEQUENTIAL matrix) 330 331 Level: advanced 332 333 Notes: 334 See the manual page for `MatCreateAIJ()` for more information on the "diagonal part" of the matrix. 335 336 Use caution, as the reference count on the returned matrix is not incremented and it is used as part of `A`'s normal operation. 337 338 .seealso: `Mat`, `MatCreateAIJ()`, `MATAIJ`, `MATBAIJ`, `MATSBAIJ` 339 @*/ 340 PetscErrorCode MatGetDiagonalBlock(Mat A, Mat *a) 341 { 342 PetscFunctionBegin; 343 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 344 PetscValidType(A, 1); 345 PetscValidPointer(a, 2); 346 PetscCheck(!A->factortype, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 347 if (A->ops->getdiagonalblock) PetscUseTypeMethod(A, getdiagonalblock, a); 348 else { 349 PetscMPIInt size; 350 351 PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size)); 352 PetscCheck(size == 1, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Not for parallel matrix type %s", ((PetscObject)A)->type_name); 353 *a = A; 354 } 355 PetscFunctionReturn(PETSC_SUCCESS); 356 } 357 358 /*@ 359 MatGetTrace - Gets the trace of a matrix. The sum of the diagonal entries. 360 361 Collective 362 363 Input Parameters: 364 . mat - the matrix 365 366 Output Parameter: 367 . trace - the sum of the diagonal entries 368 369 Level: advanced 370 371 .seealso: `Mat` 372 @*/ 373 PetscErrorCode MatGetTrace(Mat mat, PetscScalar *trace) 374 { 375 Vec diag; 376 377 PetscFunctionBegin; 378 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 379 PetscValidScalarPointer(trace, 2); 380 PetscCall(MatCreateVecs(mat, &diag, NULL)); 381 PetscCall(MatGetDiagonal(mat, diag)); 382 PetscCall(VecSum(diag, trace)); 383 PetscCall(VecDestroy(&diag)); 384 PetscFunctionReturn(PETSC_SUCCESS); 385 } 386 387 /*@ 388 MatRealPart - Zeros out the imaginary part of the matrix 389 390 Logically Collective 391 392 Input Parameters: 393 . mat - the matrix 394 395 Level: advanced 396 397 .seealso: `Mat`, `MatImaginaryPart()` 398 @*/ 399 PetscErrorCode MatRealPart(Mat mat) 400 { 401 PetscFunctionBegin; 402 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 403 PetscValidType(mat, 1); 404 PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 405 PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 406 MatCheckPreallocated(mat, 1); 407 PetscUseTypeMethod(mat, realpart); 408 PetscFunctionReturn(PETSC_SUCCESS); 409 } 410 411 /*@C 412 MatGetGhosts - Get the global indices of all ghost nodes defined by the sparse matrix 413 414 Collective 415 416 Input Parameter: 417 . mat - the matrix 418 419 Output Parameters: 420 + nghosts - number of ghosts (note for `MATBAIJ` and `MATSBAIJ` matrices there is one ghost for each block) 421 - ghosts - the global indices of the ghost points 422 423 Level: advanced 424 425 Note: 426 `nghosts` and `ghosts` are suitable to pass into `VecCreateGhost()` 427 428 .seealso: `Mat`, `VecCreateGhost()` 429 @*/ 430 PetscErrorCode MatGetGhosts(Mat mat, PetscInt *nghosts, const PetscInt *ghosts[]) 431 { 432 PetscFunctionBegin; 433 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 434 PetscValidType(mat, 1); 435 PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 436 PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 437 if (mat->ops->getghosts) PetscUseTypeMethod(mat, getghosts, nghosts, ghosts); 438 else { 439 if (nghosts) *nghosts = 0; 440 if (ghosts) *ghosts = NULL; 441 } 442 PetscFunctionReturn(PETSC_SUCCESS); 443 } 444 445 /*@ 446 MatImaginaryPart - Moves the imaginary part of the matrix to the real part and zeros the imaginary part 447 448 Logically Collective 449 450 Input Parameters: 451 . mat - the matrix 452 453 Level: advanced 454 455 .seealso: `Mat`, `MatRealPart()` 456 @*/ 457 PetscErrorCode MatImaginaryPart(Mat mat) 458 { 459 PetscFunctionBegin; 460 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 461 PetscValidType(mat, 1); 462 PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 463 PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 464 MatCheckPreallocated(mat, 1); 465 PetscUseTypeMethod(mat, imaginarypart); 466 PetscFunctionReturn(PETSC_SUCCESS); 467 } 468 469 /*@ 470 MatMissingDiagonal - Determine if sparse matrix is missing a diagonal entry (or block entry for `MATBAIJ` and `MATSBAIJ` matrices) 471 472 Not Collective 473 474 Input Parameter: 475 . mat - the matrix 476 477 Output Parameters: 478 + missing - is any diagonal missing 479 - dd - first diagonal entry that is missing (optional) on this process 480 481 Level: advanced 482 483 .seealso: `Mat` 484 @*/ 485 PetscErrorCode MatMissingDiagonal(Mat mat, PetscBool *missing, PetscInt *dd) 486 { 487 PetscFunctionBegin; 488 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 489 PetscValidType(mat, 1); 490 PetscValidBoolPointer(missing, 2); 491 PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix %s", ((PetscObject)mat)->type_name); 492 PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 493 PetscUseTypeMethod(mat, missingdiagonal, missing, dd); 494 PetscFunctionReturn(PETSC_SUCCESS); 495 } 496 497 /*@C 498 MatGetRow - Gets a row of a matrix. You MUST call `MatRestoreRow()` 499 for each row that you get to ensure that your application does 500 not bleed memory. 501 502 Not Collective 503 504 Input Parameters: 505 + mat - the matrix 506 - row - the row to get 507 508 Output Parameters: 509 + ncols - if not `NULL`, the number of nonzeros in the row 510 . cols - if not `NULL`, the column numbers 511 - vals - if not `NULL`, the values 512 513 Level: advanced 514 515 Notes: 516 This routine is provided for people who need to have direct access 517 to the structure of a matrix. We hope that we provide enough 518 high-level matrix routines that few users will need it. 519 520 `MatGetRow()` always returns 0-based column indices, regardless of 521 whether the internal representation is 0-based (default) or 1-based. 522 523 For better efficiency, set cols and/or vals to `NULL` if you do 524 not wish to extract these quantities. 525 526 The user can only examine the values extracted with `MatGetRow()`; 527 the values cannot be altered. To change the matrix entries, one 528 must use `MatSetValues()`. 529 530 You can only have one call to `MatGetRow()` outstanding for a particular 531 matrix at a time, per processor. `MatGetRow()` can only obtain rows 532 associated with the given processor, it cannot get rows from the 533 other processors; for that we suggest using `MatCreateSubMatrices()`, then 534 MatGetRow() on the submatrix. The row index passed to `MatGetRow()` 535 is in the global number of rows. 536 537 Use `MatGetRowIJ()` and `MatRestoreRowIJ()` to access all the local indices of the sparse matrix. 538 539 Use `MatSeqAIJGetArray()` and similar functions to access the numerical values for certain matrix types directly. 540 541 Fortran Note: 542 The calling sequence is 543 .vb 544 MatGetRow(matrix,row,ncols,cols,values,ierr) 545 Mat matrix (input) 546 integer row (input) 547 integer ncols (output) 548 integer cols(maxcols) (output) 549 double precision (or double complex) values(maxcols) output 550 .ve 551 where maxcols >= maximum nonzeros in any row of the matrix. 552 553 Caution: 554 Do not try to change the contents of the output arrays (cols and vals). 555 In some cases, this may corrupt the matrix. 556 557 .seealso: `Mat`, `MatRestoreRow()`, `MatSetValues()`, `MatGetValues()`, `MatCreateSubMatrices()`, `MatGetDiagonal()`, `MatGetRowIJ()`, `MatRestoreRowIJ()` 558 @*/ 559 PetscErrorCode MatGetRow(Mat mat, PetscInt row, PetscInt *ncols, const PetscInt *cols[], const PetscScalar *vals[]) 560 { 561 PetscInt incols; 562 563 PetscFunctionBegin; 564 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 565 PetscValidType(mat, 1); 566 PetscCheck(mat->assembled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 567 PetscCheck(!mat->factortype, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 568 MatCheckPreallocated(mat, 1); 569 PetscCheck(row >= mat->rmap->rstart && row < mat->rmap->rend, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Only for local rows, %" PetscInt_FMT " not in [%" PetscInt_FMT ",%" PetscInt_FMT ")", row, mat->rmap->rstart, mat->rmap->rend); 570 PetscCall(PetscLogEventBegin(MAT_GetRow, mat, 0, 0, 0)); 571 PetscCall((*mat->ops->getrow)(mat, row, &incols, (PetscInt **)cols, (PetscScalar **)vals)); 572 if (ncols) *ncols = incols; 573 PetscCall(PetscLogEventEnd(MAT_GetRow, mat, 0, 0, 0)); 574 PetscFunctionReturn(PETSC_SUCCESS); 575 } 576 577 /*@ 578 MatConjugate - replaces the matrix values with their complex conjugates 579 580 Logically Collective 581 582 Input Parameters: 583 . mat - the matrix 584 585 Level: advanced 586 587 .seealso: `Mat`, `MatRealPart()`, `MatImaginaryPart()`, `VecConjugate()`, `MatTranspose()` 588 @*/ 589 PetscErrorCode MatConjugate(Mat mat) 590 { 591 PetscFunctionBegin; 592 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 593 PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 594 if (PetscDefined(USE_COMPLEX) && mat->hermitian != PETSC_BOOL3_TRUE) { 595 PetscUseTypeMethod(mat, conjugate); 596 PetscCall(PetscObjectStateIncrease((PetscObject)mat)); 597 } 598 PetscFunctionReturn(PETSC_SUCCESS); 599 } 600 601 /*@C 602 MatRestoreRow - Frees any temporary space allocated by `MatGetRow()`. 603 604 Not Collective 605 606 Input Parameters: 607 + mat - the matrix 608 . row - the row to get 609 . ncols, cols - the number of nonzeros and their columns 610 - vals - if nonzero the column values 611 612 Level: advanced 613 614 Notes: 615 This routine should be called after you have finished examining the entries. 616 617 This routine zeros out `ncols`, `cols`, and `vals`. This is to prevent accidental 618 us of the array after it has been restored. If you pass `NULL`, it will 619 not zero the pointers. Use of `cols` or `vals` after `MatRestoreRow()` is invalid. 620 621 Fortran Notes: 622 The calling sequence is 623 .vb 624 MatRestoreRow(matrix,row,ncols,cols,values,ierr) 625 Mat matrix (input) 626 integer row (input) 627 integer ncols (output) 628 integer cols(maxcols) (output) 629 double precision (or double complex) values(maxcols) output 630 .ve 631 Where maxcols >= maximum nonzeros in any row of the matrix. 632 633 In Fortran `MatRestoreRow()` MUST be called after `MatGetRow()` 634 before another call to `MatGetRow()` can be made. 635 636 .seealso: `Mat`, `MatGetRow()` 637 @*/ 638 PetscErrorCode MatRestoreRow(Mat mat, PetscInt row, PetscInt *ncols, const PetscInt *cols[], const PetscScalar *vals[]) 639 { 640 PetscFunctionBegin; 641 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 642 if (ncols) PetscValidIntPointer(ncols, 3); 643 PetscCheck(mat->assembled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 644 if (!mat->ops->restorerow) PetscFunctionReturn(PETSC_SUCCESS); 645 PetscCall((*mat->ops->restorerow)(mat, row, ncols, (PetscInt **)cols, (PetscScalar **)vals)); 646 if (ncols) *ncols = 0; 647 if (cols) *cols = NULL; 648 if (vals) *vals = NULL; 649 PetscFunctionReturn(PETSC_SUCCESS); 650 } 651 652 /*@ 653 MatGetRowUpperTriangular - Sets a flag to enable calls to `MatGetRow()` for matrix in `MATSBAIJ` format. 654 You should call `MatRestoreRowUpperTriangular()` after calling` MatGetRow()` and `MatRestoreRow()` to disable the flag. 655 656 Not Collective 657 658 Input Parameters: 659 . mat - the matrix 660 661 Level: advanced 662 663 Note: 664 The flag is to ensure that users are aware that `MatGetRow()` only provides the upper triangular part of the row for the matrices in `MATSBAIJ` format. 665 666 .seealso: `Mat`, `MATSBAIJ`, `MatRestoreRowUpperTriangular()` 667 @*/ 668 PetscErrorCode MatGetRowUpperTriangular(Mat mat) 669 { 670 PetscFunctionBegin; 671 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 672 PetscValidType(mat, 1); 673 PetscCheck(mat->assembled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 674 PetscCheck(!mat->factortype, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 675 MatCheckPreallocated(mat, 1); 676 if (!mat->ops->getrowuppertriangular) PetscFunctionReturn(PETSC_SUCCESS); 677 PetscUseTypeMethod(mat, getrowuppertriangular); 678 PetscFunctionReturn(PETSC_SUCCESS); 679 } 680 681 /*@ 682 MatRestoreRowUpperTriangular - Disable calls to `MatGetRow()` for matrix in `MATSBAIJ` format. 683 684 Not Collective 685 686 Input Parameters: 687 . mat - the matrix 688 689 Level: advanced 690 691 Note: 692 This routine should be called after you have finished calls to `MatGetRow()` and `MatRestoreRow()`. 693 694 .seealso: `Mat`, `MATSBAIJ`, `MatGetRowUpperTriangular()` 695 @*/ 696 PetscErrorCode MatRestoreRowUpperTriangular(Mat mat) 697 { 698 PetscFunctionBegin; 699 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 700 PetscValidType(mat, 1); 701 PetscCheck(mat->assembled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 702 PetscCheck(!mat->factortype, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 703 MatCheckPreallocated(mat, 1); 704 if (!mat->ops->restorerowuppertriangular) PetscFunctionReturn(PETSC_SUCCESS); 705 PetscUseTypeMethod(mat, restorerowuppertriangular); 706 PetscFunctionReturn(PETSC_SUCCESS); 707 } 708 709 /*@C 710 MatSetOptionsPrefix - Sets the prefix used for searching for all 711 `Mat` options in the database. 712 713 Logically Collective 714 715 Input Parameters: 716 + A - the matrix 717 - prefix - the prefix to prepend to all option names 718 719 Level: advanced 720 721 Notes: 722 A hyphen (-) must NOT be given at the beginning of the prefix name. 723 The first character of all runtime options is AUTOMATICALLY the hyphen. 724 725 This is NOT used for options for the factorization of the matrix. Normally the 726 prefix is automatically passed in from the PC calling the factorization. To set 727 it directly use `MatSetOptionsPrefixFactor()` 728 729 .seealso: `Mat`, `MatSetFromOptions()`, `MatSetOptionsPrefixFactor()` 730 @*/ 731 PetscErrorCode MatSetOptionsPrefix(Mat A, const char prefix[]) 732 { 733 PetscFunctionBegin; 734 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 735 PetscCall(PetscObjectSetOptionsPrefix((PetscObject)A, prefix)); 736 PetscFunctionReturn(PETSC_SUCCESS); 737 } 738 739 /*@C 740 MatSetOptionsPrefixFactor - Sets the prefix used for searching for all matrix factor options in the database for 741 for matrices created with `MatGetFactor()` 742 743 Logically Collective 744 745 Input Parameters: 746 + A - the matrix 747 - prefix - the prefix to prepend to all option names for the factored matrix 748 749 Level: developer 750 751 Notes: 752 A hyphen (-) must NOT be given at the beginning of the prefix name. 753 The first character of all runtime options is AUTOMATICALLY the hyphen. 754 755 Normally the prefix is automatically passed in from the `PC` calling the factorization. To set 756 it directly when not using `KSP`/`PC` use `MatSetOptionsPrefixFactor()` 757 758 .seealso: `Mat`, [Matrix Factorization](sec_matfactor), `MatGetFactor()`, `MatSetFromOptions()`, `MatSetOptionsPrefix()`, `MatAppendOptionsPrefixFactor()` 759 @*/ 760 PetscErrorCode MatSetOptionsPrefixFactor(Mat A, const char prefix[]) 761 { 762 PetscFunctionBegin; 763 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 764 if (prefix) { 765 PetscValidCharPointer(prefix, 2); 766 PetscCheck(prefix[0] != '-', PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Options prefix should not begin with a hyphen"); 767 if (prefix != A->factorprefix) { 768 PetscCall(PetscFree(A->factorprefix)); 769 PetscCall(PetscStrallocpy(prefix, &A->factorprefix)); 770 } 771 } else PetscCall(PetscFree(A->factorprefix)); 772 PetscFunctionReturn(PETSC_SUCCESS); 773 } 774 775 /*@C 776 MatAppendOptionsPrefixFactor - Appends to the prefix used for searching for all matrix factor options in the database for 777 for matrices created with `MatGetFactor()` 778 779 Logically Collective 780 781 Input Parameters: 782 + A - the matrix 783 - prefix - the prefix to prepend to all option names for the factored matrix 784 785 Level: developer 786 787 Notes: 788 A hyphen (-) must NOT be given at the beginning of the prefix name. 789 The first character of all runtime options is AUTOMATICALLY the hyphen. 790 791 Normally the prefix is automatically passed in from the `PC` calling the factorization. To set 792 it directly when not using `KSP`/`PC` use `MatAppendOptionsPrefixFactor()` 793 794 .seealso: `Mat`, [Matrix Factorization](sec_matfactor), `MatGetFactor()`, PetscOptionsCreate()`, `PetscOptionsDestroy()`, `PetscObjectSetOptionsPrefix()`, `PetscObjectPrependOptionsPrefix()`, 795 `PetscObjectGetOptionsPrefix()`, `TSAppendOptionsPrefix()`, `SNESAppendOptionsPrefix()`, `KSPAppendOptionsPrefix()`, `MatSetOptionsPrefixFactor()`, 796 `MatSetOptionsPrefix()` 797 @*/ 798 PetscErrorCode MatAppendOptionsPrefixFactor(Mat A, const char prefix[]) 799 { 800 char *buf = A->factorprefix; 801 size_t len1, len2; 802 803 PetscFunctionBegin; 804 PetscValidHeader(A, 1); 805 if (!prefix) PetscFunctionReturn(PETSC_SUCCESS); 806 if (!buf) { 807 PetscCall(MatSetOptionsPrefixFactor(A, prefix)); 808 PetscFunctionReturn(PETSC_SUCCESS); 809 } 810 PetscCheck(prefix[0] != '-', PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Options prefix should not begin with a hyphen"); 811 812 PetscCall(PetscStrlen(prefix, &len1)); 813 PetscCall(PetscStrlen(buf, &len2)); 814 PetscCall(PetscMalloc1(1 + len1 + len2, &A->factorprefix)); 815 PetscCall(PetscStrcpy(A->factorprefix, buf)); 816 PetscCall(PetscStrcat(A->factorprefix, prefix)); 817 PetscCall(PetscFree(buf)); 818 PetscFunctionReturn(PETSC_SUCCESS); 819 } 820 821 /*@C 822 MatAppendOptionsPrefix - Appends to the prefix used for searching for all 823 matrix options in the database. 824 825 Logically Collective 826 827 Input Parameters: 828 + A - the matrix 829 - prefix - the prefix to prepend to all option names 830 831 Level: advanced 832 833 Note: 834 A hyphen (-) must NOT be given at the beginning of the prefix name. 835 The first character of all runtime options is AUTOMATICALLY the hyphen. 836 837 .seealso: `Mat`, `MatGetOptionsPrefix()`, `MatAppendOptionsPrefixFactor()`, `MatSetOptionsPrefix()` 838 @*/ 839 PetscErrorCode MatAppendOptionsPrefix(Mat A, const char prefix[]) 840 { 841 PetscFunctionBegin; 842 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 843 PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)A, prefix)); 844 PetscFunctionReturn(PETSC_SUCCESS); 845 } 846 847 /*@C 848 MatGetOptionsPrefix - Gets the prefix used for searching for all 849 matrix options in the database. 850 851 Not Collective 852 853 Input Parameter: 854 . A - the matrix 855 856 Output Parameter: 857 . prefix - pointer to the prefix string used 858 859 Level: advanced 860 861 Fortran Note: 862 On the fortran side, the user should pass in a string `prefix` of 863 sufficient length to hold the prefix. 864 865 .seealso: `Mat`, `MatAppendOptionsPrefix()`, `MatSetOptionsPrefix()`, `MatAppendOptionsPrefixFactor()`, `MatSetOptionsPrefixFactor()` 866 @*/ 867 PetscErrorCode MatGetOptionsPrefix(Mat A, const char *prefix[]) 868 { 869 PetscFunctionBegin; 870 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 871 PetscValidPointer(prefix, 2); 872 PetscCall(PetscObjectGetOptionsPrefix((PetscObject)A, prefix)); 873 PetscFunctionReturn(PETSC_SUCCESS); 874 } 875 876 /*@ 877 MatResetPreallocation - Reset matrix to use the original nonzero pattern provided by users. 878 879 Collective 880 881 Input Parameters: 882 . A - the matrix 883 884 Level: beginner 885 886 Notes: 887 The allocated memory will be shrunk after calling `MatAssemblyBegin()` and `MatAssemblyEnd()` with `MAT_FINAL_ASSEMBLY`. 888 889 Users can reset the preallocation to access the original memory. 890 891 Currently only supported for `MATMPIAIJ` and `MATSEQAIJ` matrices. 892 893 .seealso: `Mat`, `MatSeqAIJSetPreallocation()`, `MatMPIAIJSetPreallocation()`, `MatXAIJSetPreallocation()` 894 @*/ 895 PetscErrorCode MatResetPreallocation(Mat A) 896 { 897 PetscFunctionBegin; 898 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 899 PetscValidType(A, 1); 900 PetscUseMethod(A, "MatResetPreallocation_C", (Mat), (A)); 901 PetscFunctionReturn(PETSC_SUCCESS); 902 } 903 904 /*@ 905 MatSetUp - Sets up the internal matrix data structures for later use. 906 907 Collective 908 909 Input Parameters: 910 . A - the matrix 911 912 Level: intermediate 913 914 Notes: 915 If the user has not set preallocation for this matrix then a default preallocation that is likely to be inefficient is used. 916 917 If a suitable preallocation routine is used, this function does not need to be called. 918 919 See the Performance chapter of the PETSc users manual for how to preallocate matrices 920 921 This routine is called internally by other matrix functions when needed so rarely needs to be called by users 922 923 .seealso: `Mat`, `MatMult()`, `MatCreate()`, `MatDestroy()` 924 @*/ 925 PetscErrorCode MatSetUp(Mat A) 926 { 927 PetscFunctionBegin; 928 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 929 if (!((PetscObject)A)->type_name) { 930 PetscMPIInt size; 931 932 PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size)); 933 PetscCall(MatSetType(A, size == 1 ? MATSEQAIJ : MATMPIAIJ)); 934 } 935 if (!A->preallocated && A->ops->setup) { 936 PetscCall(PetscInfo(A, "Warning not preallocating matrix storage\n")); 937 PetscUseTypeMethod(A, setup); 938 } 939 PetscCall(PetscLayoutSetUp(A->rmap)); 940 PetscCall(PetscLayoutSetUp(A->cmap)); 941 A->preallocated = PETSC_TRUE; 942 PetscFunctionReturn(PETSC_SUCCESS); 943 } 944 945 #if defined(PETSC_HAVE_SAWS) 946 #include <petscviewersaws.h> 947 #endif 948 949 /*@C 950 MatViewFromOptions - View properties of the matrix from the options database 951 952 Collective 953 954 Input Parameters: 955 + A - the matrix 956 . obj - optional additional object that provides the options prefix to use 957 - name - command line option 958 959 Options Database Key: 960 . -mat_view [viewertype]:... - the viewer and its options 961 962 Level: intermediate 963 964 Notes: 965 .vb 966 If no value is provided ascii:stdout is used 967 ascii[:[filename][:[format][:append]]] defaults to stdout - format can be one of ascii_info, ascii_info_detail, or ascii_matlab, 968 for example ascii::ascii_info prints just the information about the object not all details 969 unless :append is given filename opens in write mode, overwriting what was already there 970 binary[:[filename][:[format][:append]]] defaults to the file binaryoutput 971 draw[:drawtype[:filename]] for example, draw:tikz, draw:tikz:figure.tex or draw:x 972 socket[:port] defaults to the standard output port 973 saws[:communicatorname] publishes object to the Scientific Application Webserver (SAWs) 974 .ve 975 976 .seealso: `Mat`, `MatView()`, `PetscObjectViewFromOptions()`, `MatCreate()` 977 @*/ 978 PetscErrorCode MatViewFromOptions(Mat A, PetscObject obj, const char name[]) 979 { 980 PetscFunctionBegin; 981 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 982 PetscCall(PetscObjectViewFromOptions((PetscObject)A, obj, name)); 983 PetscFunctionReturn(PETSC_SUCCESS); 984 } 985 986 /*@C 987 MatView - display information about a matrix in a variety ways 988 989 Collective 990 991 Input Parameters: 992 + mat - the matrix 993 - viewer - visualization context 994 995 Notes: 996 The available visualization contexts include 997 + `PETSC_VIEWER_STDOUT_SELF` - for sequential matrices 998 . `PETSC_VIEWER_STDOUT_WORLD` - for parallel matrices created on `PETSC_COMM_WORLD` 999 . `PETSC_VIEWER_STDOUT_`(comm) - for matrices created on MPI communicator comm 1000 - `PETSC_VIEWER_DRAW_WORLD` - graphical display of nonzero structure 1001 1002 The user can open alternative visualization contexts with 1003 + `PetscViewerASCIIOpen()` - Outputs matrix to a specified file 1004 . `PetscViewerBinaryOpen()` - Outputs matrix in binary to a 1005 specified file; corresponding input uses MatLoad() 1006 . `PetscViewerDrawOpen()` - Outputs nonzero matrix structure to 1007 an X window display 1008 - `PetscViewerSocketOpen()` - Outputs matrix to Socket viewer. 1009 Currently only the sequential dense and AIJ 1010 matrix types support the Socket viewer. 1011 1012 The user can call `PetscViewerPushFormat()` to specify the output 1013 format of ASCII printed objects (when using `PETSC_VIEWER_STDOUT_SELF`, 1014 `PETSC_VIEWER_STDOUT_WORLD` and `PetscViewerASCIIOpen()`). Available formats include 1015 + `PETSC_VIEWER_DEFAULT` - default, prints matrix contents 1016 . `PETSC_VIEWER_ASCII_MATLAB` - prints matrix contents in Matlab format 1017 . `PETSC_VIEWER_ASCII_DENSE` - prints entire matrix including zeros 1018 . `PETSC_VIEWER_ASCII_COMMON` - prints matrix contents, using a sparse 1019 format common among all matrix types 1020 . `PETSC_VIEWER_ASCII_IMPL` - prints matrix contents, using an implementation-specific 1021 format (which is in many cases the same as the default) 1022 . `PETSC_VIEWER_ASCII_INFO` - prints basic information about the matrix 1023 size and structure (not the matrix entries) 1024 - `PETSC_VIEWER_ASCII_INFO_DETAIL` - prints more detailed information about 1025 the matrix structure 1026 1027 Options Database Keys: 1028 + -mat_view ::ascii_info - Prints info on matrix at conclusion of `MatAssemblyEnd()` 1029 . -mat_view ::ascii_info_detail - Prints more detailed info 1030 . -mat_view - Prints matrix in ASCII format 1031 . -mat_view ::ascii_matlab - Prints matrix in Matlab format 1032 . -mat_view draw - PetscDraws nonzero structure of matrix, using `MatView()` and `PetscDrawOpenX()`. 1033 . -display <name> - Sets display name (default is host) 1034 . -draw_pause <sec> - Sets number of seconds to pause after display 1035 . -mat_view socket - Sends matrix to socket, can be accessed from Matlab (see Users-Manual: ch_matlab for details) 1036 . -viewer_socket_machine <machine> - 1037 . -viewer_socket_port <port> - 1038 . -mat_view binary - save matrix to file in binary format 1039 - -viewer_binary_filename <name> - 1040 1041 Level: beginner 1042 1043 Notes: 1044 The ASCII viewers are only recommended for small matrices on at most a moderate number of processes, 1045 the program will seemingly hang and take hours for larger matrices, for larger matrices one should use the binary format. 1046 1047 In the debugger you can do "call MatView(mat,0)" to display the matrix. (The same holds for any PETSc object viewer). 1048 1049 See the manual page for `MatLoad()` for the exact format of the binary file when the binary 1050 viewer is used. 1051 1052 See share/petsc/matlab/PetscBinaryRead.m for a Matlab code that can read in the binary file when the binary 1053 viewer is used and lib/petsc/bin/PetscBinaryIO.py for loading them into Python. 1054 1055 One can use '-mat_view draw -draw_pause -1' to pause the graphical display of matrix nonzero structure, 1056 and then use the following mouse functions. 1057 .vb 1058 left mouse: zoom in 1059 middle mouse: zoom out 1060 right mouse: continue with the simulation 1061 .ve 1062 1063 .seealso: `Mat`, `PetscViewerPushFormat()`, `PetscViewerASCIIOpen()`, `PetscViewerDrawOpen()`, `PetscViewer`, 1064 `PetscViewerSocketOpen()`, `PetscViewerBinaryOpen()`, `MatLoad()`, `MatViewFromOptions()` 1065 @*/ 1066 PetscErrorCode MatView(Mat mat, PetscViewer viewer) 1067 { 1068 PetscInt rows, cols, rbs, cbs; 1069 PetscBool isascii, isstring, issaws; 1070 PetscViewerFormat format; 1071 PetscMPIInt size; 1072 1073 PetscFunctionBegin; 1074 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 1075 PetscValidType(mat, 1); 1076 if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)mat), &viewer)); 1077 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 1078 PetscCheckSameComm(mat, 1, viewer, 2); 1079 1080 PetscCall(PetscViewerGetFormat(viewer, &format)); 1081 PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)mat), &size)); 1082 if (size == 1 && format == PETSC_VIEWER_LOAD_BALANCE) PetscFunctionReturn(PETSC_SUCCESS); 1083 1084 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring)); 1085 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 1086 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSAWS, &issaws)); 1087 PetscCheck((isascii && (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL)) || !mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "No viewers for factored matrix except ASCII, info, or info_detail"); 1088 1089 PetscCall(PetscLogEventBegin(MAT_View, mat, viewer, 0, 0)); 1090 if (isascii) { 1091 if (!mat->preallocated) { 1092 PetscCall(PetscViewerASCIIPrintf(viewer, "Matrix has not been preallocated yet\n")); 1093 PetscFunctionReturn(PETSC_SUCCESS); 1094 } 1095 if (!mat->assembled) { 1096 PetscCall(PetscViewerASCIIPrintf(viewer, "Matrix has not been assembled yet\n")); 1097 PetscFunctionReturn(PETSC_SUCCESS); 1098 } 1099 PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)mat, viewer)); 1100 if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 1101 MatNullSpace nullsp, transnullsp; 1102 1103 PetscCall(PetscViewerASCIIPushTab(viewer)); 1104 PetscCall(MatGetSize(mat, &rows, &cols)); 1105 PetscCall(MatGetBlockSizes(mat, &rbs, &cbs)); 1106 if (rbs != 1 || cbs != 1) { 1107 if (rbs != cbs) PetscCall(PetscViewerASCIIPrintf(viewer, "rows=%" PetscInt_FMT ", cols=%" PetscInt_FMT ", rbs=%" PetscInt_FMT ", cbs=%" PetscInt_FMT "\n", rows, cols, rbs, cbs)); 1108 else PetscCall(PetscViewerASCIIPrintf(viewer, "rows=%" PetscInt_FMT ", cols=%" PetscInt_FMT ", bs=%" PetscInt_FMT "\n", rows, cols, rbs)); 1109 } else PetscCall(PetscViewerASCIIPrintf(viewer, "rows=%" PetscInt_FMT ", cols=%" PetscInt_FMT "\n", rows, cols)); 1110 if (mat->factortype) { 1111 MatSolverType solver; 1112 PetscCall(MatFactorGetSolverType(mat, &solver)); 1113 PetscCall(PetscViewerASCIIPrintf(viewer, "package used to perform factorization: %s\n", solver)); 1114 } 1115 if (mat->ops->getinfo) { 1116 MatInfo info; 1117 PetscCall(MatGetInfo(mat, MAT_GLOBAL_SUM, &info)); 1118 PetscCall(PetscViewerASCIIPrintf(viewer, "total: nonzeros=%.f, allocated nonzeros=%.f\n", info.nz_used, info.nz_allocated)); 1119 if (!mat->factortype) PetscCall(PetscViewerASCIIPrintf(viewer, "total number of mallocs used during MatSetValues calls=%" PetscInt_FMT "\n", (PetscInt)info.mallocs)); 1120 } 1121 PetscCall(MatGetNullSpace(mat, &nullsp)); 1122 PetscCall(MatGetTransposeNullSpace(mat, &transnullsp)); 1123 if (nullsp) PetscCall(PetscViewerASCIIPrintf(viewer, " has attached null space\n")); 1124 if (transnullsp && transnullsp != nullsp) PetscCall(PetscViewerASCIIPrintf(viewer, " has attached transposed null space\n")); 1125 PetscCall(MatGetNearNullSpace(mat, &nullsp)); 1126 if (nullsp) PetscCall(PetscViewerASCIIPrintf(viewer, " has attached near null space\n")); 1127 PetscCall(PetscViewerASCIIPushTab(viewer)); 1128 PetscCall(MatProductView(mat, viewer)); 1129 PetscCall(PetscViewerASCIIPopTab(viewer)); 1130 } 1131 } else if (issaws) { 1132 #if defined(PETSC_HAVE_SAWS) 1133 PetscMPIInt rank; 1134 1135 PetscCall(PetscObjectName((PetscObject)mat)); 1136 PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 1137 if (!((PetscObject)mat)->amsmem && rank == 0) PetscCall(PetscObjectViewSAWs((PetscObject)mat, viewer)); 1138 #endif 1139 } else if (isstring) { 1140 const char *type; 1141 PetscCall(MatGetType(mat, &type)); 1142 PetscCall(PetscViewerStringSPrintf(viewer, " MatType: %-7.7s", type)); 1143 PetscTryTypeMethod(mat, view, viewer); 1144 } 1145 if ((format == PETSC_VIEWER_NATIVE || format == PETSC_VIEWER_LOAD_BALANCE) && mat->ops->viewnative) { 1146 PetscCall(PetscViewerASCIIPushTab(viewer)); 1147 PetscUseTypeMethod(mat, viewnative, viewer); 1148 PetscCall(PetscViewerASCIIPopTab(viewer)); 1149 } else if (mat->ops->view) { 1150 PetscCall(PetscViewerASCIIPushTab(viewer)); 1151 PetscUseTypeMethod(mat, view, viewer); 1152 PetscCall(PetscViewerASCIIPopTab(viewer)); 1153 } 1154 if (isascii) { 1155 PetscCall(PetscViewerGetFormat(viewer, &format)); 1156 if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) PetscCall(PetscViewerASCIIPopTab(viewer)); 1157 } 1158 PetscCall(PetscLogEventEnd(MAT_View, mat, viewer, 0, 0)); 1159 PetscFunctionReturn(PETSC_SUCCESS); 1160 } 1161 1162 #if defined(PETSC_USE_DEBUG) 1163 #include <../src/sys/totalview/tv_data_display.h> 1164 PETSC_UNUSED static int TV_display_type(const struct _p_Mat *mat) 1165 { 1166 TV_add_row("Local rows", "int", &mat->rmap->n); 1167 TV_add_row("Local columns", "int", &mat->cmap->n); 1168 TV_add_row("Global rows", "int", &mat->rmap->N); 1169 TV_add_row("Global columns", "int", &mat->cmap->N); 1170 TV_add_row("Typename", TV_ascii_string_type, ((PetscObject)mat)->type_name); 1171 return TV_format_OK; 1172 } 1173 #endif 1174 1175 /*@C 1176 MatLoad - Loads a matrix that has been stored in binary/HDF5 format 1177 with `MatView()`. The matrix format is determined from the options database. 1178 Generates a parallel MPI matrix if the communicator has more than one 1179 processor. The default matrix type is `MATAIJ`. 1180 1181 Collective 1182 1183 Input Parameters: 1184 + mat - the newly loaded matrix, this needs to have been created with `MatCreate()` 1185 or some related function before a call to `MatLoad()` 1186 - viewer - binary/HDF5 file viewer 1187 1188 Options Database Keys: 1189 Used with block matrix formats (`MATSEQBAIJ`, ...) to specify 1190 block size 1191 . -matload_block_size <bs> - set block size 1192 1193 Level: beginner 1194 1195 Notes: 1196 If the `Mat` type has not yet been given then `MATAIJ` is used, call `MatSetFromOptions()` on the 1197 `Mat` before calling this routine if you wish to set it from the options database. 1198 1199 `MatLoad()` automatically loads into the options database any options 1200 given in the file filename.info where filename is the name of the file 1201 that was passed to the `PetscViewerBinaryOpen()`. The options in the info 1202 file will be ignored if you use the -viewer_binary_skip_info option. 1203 1204 If the type or size of mat is not set before a call to `MatLoad()`, PETSc 1205 sets the default matrix type AIJ and sets the local and global sizes. 1206 If type and/or size is already set, then the same are used. 1207 1208 In parallel, each processor can load a subset of rows (or the 1209 entire matrix). This routine is especially useful when a large 1210 matrix is stored on disk and only part of it is desired on each 1211 processor. For example, a parallel solver may access only some of 1212 the rows from each processor. The algorithm used here reads 1213 relatively small blocks of data rather than reading the entire 1214 matrix and then subsetting it. 1215 1216 Viewer's `PetscViewerType` must be either `PETSCVIEWERBINARY` or `PETSCVIEWERHDF5`. 1217 Such viewer can be created using `PetscViewerBinaryOpen()` or `PetscViewerHDF5Open()`, 1218 or the sequence like 1219 .vb 1220 `PetscViewer` v; 1221 `PetscViewerCreate`(`PETSC_COMM_WORLD`,&v); 1222 `PetscViewerSetType`(v,`PETSCVIEWERBINARY`); 1223 `PetscViewerSetFromOptions`(v); 1224 `PetscViewerFileSetMode`(v,`FILE_MODE_READ`); 1225 `PetscViewerFileSetName`(v,"datafile"); 1226 .ve 1227 The optional `PetscViewerSetFromOptions()` call allows overriding `PetscViewerSetType()` using the option 1228 $ -viewer_type {binary,hdf5} 1229 1230 See the example src/ksp/ksp/tutorials/ex27.c with the first approach, 1231 and src/mat/tutorials/ex10.c with the second approach. 1232 1233 Notes about the PETSc binary format: 1234 In case of `PETSCVIEWERBINARY`, a native PETSc binary format is used. Each of the blocks 1235 is read onto rank 0 and then shipped to its destination rank, one after another. 1236 Multiple objects, both matrices and vectors, can be stored within the same file. 1237 Their PetscObject name is ignored; they are loaded in the order of their storage. 1238 1239 Most users should not need to know the details of the binary storage 1240 format, since `MatLoad()` and `MatView()` completely hide these details. 1241 But for anyone who's interested, the standard binary matrix storage 1242 format is 1243 1244 $ PetscInt MAT_FILE_CLASSID 1245 $ PetscInt number of rows 1246 $ PetscInt number of columns 1247 $ PetscInt total number of nonzeros 1248 $ PetscInt *number nonzeros in each row 1249 $ PetscInt *column indices of all nonzeros (starting index is zero) 1250 $ PetscScalar *values of all nonzeros 1251 1252 PETSc automatically does the byte swapping for 1253 machines that store the bytes reversed, e.g. DEC alpha, freebsd, 1254 Linux, Microsoft Windows and the Intel Paragon; thus if you write your own binary 1255 read/write routines you have to swap the bytes; see `PetscBinaryRead()` 1256 and `PetscBinaryWrite()` to see how this may be done. 1257 1258 Notes about the HDF5 (MATLAB MAT-File Version 7.3) format: 1259 In case of `PETSCVIEWERHDF5`, a parallel HDF5 reader is used. 1260 Each processor's chunk is loaded independently by its owning rank. 1261 Multiple objects, both matrices and vectors, can be stored within the same file. 1262 They are looked up by their PetscObject name. 1263 1264 As the MATLAB MAT-File Version 7.3 format is also a HDF5 flavor, we decided to use 1265 by default the same structure and naming of the AIJ arrays and column count 1266 within the HDF5 file. This means that a MAT file saved with -v7.3 flag, e.g. 1267 $ save example.mat A b -v7.3 1268 can be directly read by this routine (see Reference 1 for details). 1269 Note that depending on your MATLAB version, this format might be a default, 1270 otherwise you can set it as default in Preferences. 1271 1272 Unless -nocompression flag is used to save the file in MATLAB, 1273 PETSc must be configured with ZLIB package. 1274 1275 See also examples src/mat/tutorials/ex10.c and src/ksp/ksp/tutorials/ex27.c 1276 1277 Current HDF5 (MAT-File) limitations: 1278 This reader currently supports only real `MATSEQAIJ`, `MATMPIAIJ`, `MATSEQDENSE` and `MATMPIDENSE` matrices. 1279 1280 Corresponding `MatView()` is not yet implemented. 1281 1282 The loaded matrix is actually a transpose of the original one in MATLAB, 1283 unless you push `PETSC_VIEWER_HDF5_MAT` format (see examples above). 1284 With this format, matrix is automatically transposed by PETSc, 1285 unless the matrix is marked as SPD or symmetric 1286 (see `MatSetOption()`, `MAT_SPD`, `MAT_SYMMETRIC`). 1287 1288 References: 1289 . * - MATLAB(R) Documentation, manual page of save(), https://www.mathworks.com/help/matlab/ref/save.html#btox10b-1-version 1290 1291 .seealso: `Mat`, `PetscViewerBinaryOpen()`, `PetscViewerSetType()`, `MatView()`, `VecLoad()` 1292 @*/ 1293 PetscErrorCode MatLoad(Mat mat, PetscViewer viewer) 1294 { 1295 PetscBool flg; 1296 1297 PetscFunctionBegin; 1298 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 1299 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 1300 1301 if (!((PetscObject)mat)->type_name) PetscCall(MatSetType(mat, MATAIJ)); 1302 1303 flg = PETSC_FALSE; 1304 PetscCall(PetscOptionsGetBool(((PetscObject)mat)->options, ((PetscObject)mat)->prefix, "-matload_symmetric", &flg, NULL)); 1305 if (flg) { 1306 PetscCall(MatSetOption(mat, MAT_SYMMETRIC, PETSC_TRUE)); 1307 PetscCall(MatSetOption(mat, MAT_SYMMETRY_ETERNAL, PETSC_TRUE)); 1308 } 1309 flg = PETSC_FALSE; 1310 PetscCall(PetscOptionsGetBool(((PetscObject)mat)->options, ((PetscObject)mat)->prefix, "-matload_spd", &flg, NULL)); 1311 if (flg) PetscCall(MatSetOption(mat, MAT_SPD, PETSC_TRUE)); 1312 1313 PetscCall(PetscLogEventBegin(MAT_Load, mat, viewer, 0, 0)); 1314 PetscUseTypeMethod(mat, load, viewer); 1315 PetscCall(PetscLogEventEnd(MAT_Load, mat, viewer, 0, 0)); 1316 PetscFunctionReturn(PETSC_SUCCESS); 1317 } 1318 1319 static PetscErrorCode MatDestroy_Redundant(Mat_Redundant **redundant) 1320 { 1321 Mat_Redundant *redund = *redundant; 1322 1323 PetscFunctionBegin; 1324 if (redund) { 1325 if (redund->matseq) { /* via MatCreateSubMatrices() */ 1326 PetscCall(ISDestroy(&redund->isrow)); 1327 PetscCall(ISDestroy(&redund->iscol)); 1328 PetscCall(MatDestroySubMatrices(1, &redund->matseq)); 1329 } else { 1330 PetscCall(PetscFree2(redund->send_rank, redund->recv_rank)); 1331 PetscCall(PetscFree(redund->sbuf_j)); 1332 PetscCall(PetscFree(redund->sbuf_a)); 1333 for (PetscInt i = 0; i < redund->nrecvs; i++) { 1334 PetscCall(PetscFree(redund->rbuf_j[i])); 1335 PetscCall(PetscFree(redund->rbuf_a[i])); 1336 } 1337 PetscCall(PetscFree4(redund->sbuf_nz, redund->rbuf_nz, redund->rbuf_j, redund->rbuf_a)); 1338 } 1339 1340 if (redund->subcomm) PetscCall(PetscCommDestroy(&redund->subcomm)); 1341 PetscCall(PetscFree(redund)); 1342 } 1343 PetscFunctionReturn(PETSC_SUCCESS); 1344 } 1345 1346 /*@C 1347 MatDestroy - Frees space taken by a matrix. 1348 1349 Collective 1350 1351 Input Parameter: 1352 . A - the matrix 1353 1354 Level: beginner 1355 1356 Developer Note: 1357 Some special arrays of matrices are not destroyed in this routine but instead by the routines called by 1358 `MatDestroySubMatrices()`. Thus one must be sure that any changes here must also be made in those routines. 1359 `MatHeaderMerge()` and `MatHeaderReplace()` also manipulate the data in the `Mat` object and likely need changes 1360 if changes are needed here. 1361 1362 .seealso: `Mat`, `MatCreate()` 1363 @*/ 1364 PetscErrorCode MatDestroy(Mat *A) 1365 { 1366 PetscFunctionBegin; 1367 if (!*A) PetscFunctionReturn(PETSC_SUCCESS); 1368 PetscValidHeaderSpecific(*A, MAT_CLASSID, 1); 1369 if (--((PetscObject)(*A))->refct > 0) { 1370 *A = NULL; 1371 PetscFunctionReturn(PETSC_SUCCESS); 1372 } 1373 1374 /* if memory was published with SAWs then destroy it */ 1375 PetscCall(PetscObjectSAWsViewOff((PetscObject)*A)); 1376 PetscTryTypeMethod((*A), destroy); 1377 1378 PetscCall(PetscFree((*A)->factorprefix)); 1379 PetscCall(PetscFree((*A)->defaultvectype)); 1380 PetscCall(PetscFree((*A)->defaultrandtype)); 1381 PetscCall(PetscFree((*A)->bsizes)); 1382 PetscCall(PetscFree((*A)->solvertype)); 1383 for (PetscInt i = 0; i < MAT_FACTOR_NUM_TYPES; i++) PetscCall(PetscFree((*A)->preferredordering[i])); 1384 if ((*A)->redundant && (*A)->redundant->matseq[0] == *A) (*A)->redundant->matseq[0] = NULL; 1385 PetscCall(MatDestroy_Redundant(&(*A)->redundant)); 1386 PetscCall(MatProductClear(*A)); 1387 PetscCall(MatNullSpaceDestroy(&(*A)->nullsp)); 1388 PetscCall(MatNullSpaceDestroy(&(*A)->transnullsp)); 1389 PetscCall(MatNullSpaceDestroy(&(*A)->nearnullsp)); 1390 PetscCall(MatDestroy(&(*A)->schur)); 1391 PetscCall(PetscLayoutDestroy(&(*A)->rmap)); 1392 PetscCall(PetscLayoutDestroy(&(*A)->cmap)); 1393 PetscCall(PetscHeaderDestroy(A)); 1394 PetscFunctionReturn(PETSC_SUCCESS); 1395 } 1396 1397 /*@C 1398 MatSetValues - Inserts or adds a block of values into a matrix. 1399 These values may be cached, so `MatAssemblyBegin()` and `MatAssemblyEnd()` 1400 MUST be called after all calls to `MatSetValues()` have been completed. 1401 1402 Not Collective 1403 1404 Input Parameters: 1405 + mat - the matrix 1406 . v - a logically two-dimensional array of values 1407 . m, idxm - the number of rows and their global indices 1408 . n, idxn - the number of columns and their global indices 1409 - addv - either `ADD_VALUES` to add values to any existing entries, or `INSERT_VALUES` to replace existing entries with new values 1410 1411 Level: beginner 1412 1413 Notes: 1414 If you create the matrix yourself (that is not with a call to `DMCreateMatrix()`) then you MUST call MatXXXXSetPreallocation() or 1415 `MatSetUp()` before using this routine 1416 1417 By default the values, v, are row-oriented. See `MatSetOption()` for other options. 1418 1419 Calls to `MatSetValues()` with the `INSERT_VALUES` and `ADD_VALUES` 1420 options cannot be mixed without intervening calls to the assembly 1421 routines. 1422 1423 `MatSetValues()` uses 0-based row and column numbers in Fortran 1424 as well as in C. 1425 1426 Negative indices may be passed in idxm and idxn, these rows and columns are 1427 simply ignored. This allows easily inserting element stiffness matrices 1428 with homogeneous Dirchlet boundary conditions that you don't want represented 1429 in the matrix. 1430 1431 Efficiency Alert: 1432 The routine `MatSetValuesBlocked()` may offer much better efficiency 1433 for users of block sparse formats (`MATSEQBAIJ` and `MATMPIBAIJ`). 1434 1435 Developer Note: 1436 This is labeled with C so does not automatically generate Fortran stubs and interfaces 1437 because it requires multiple Fortran interfaces depending on which arguments are scalar or arrays. 1438 1439 .seealso: `Mat`, `MatSetOption()`, `MatAssemblyBegin()`, `MatAssemblyEnd()`, `MatSetValuesBlocked()`, `MatSetValuesLocal()`, 1440 `InsertMode`, `INSERT_VALUES`, `ADD_VALUES` 1441 @*/ 1442 PetscErrorCode MatSetValues(Mat mat, PetscInt m, const PetscInt idxm[], PetscInt n, const PetscInt idxn[], const PetscScalar v[], InsertMode addv) 1443 { 1444 PetscFunctionBeginHot; 1445 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 1446 PetscValidType(mat, 1); 1447 if (!m || !n) PetscFunctionReturn(PETSC_SUCCESS); /* no values to insert */ 1448 PetscValidIntPointer(idxm, 3); 1449 PetscValidIntPointer(idxn, 5); 1450 MatCheckPreallocated(mat, 1); 1451 1452 if (mat->insertmode == NOT_SET_VALUES) mat->insertmode = addv; 1453 else PetscCheck(mat->insertmode == addv, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Cannot mix add values and insert values"); 1454 1455 if (PetscDefined(USE_DEBUG)) { 1456 PetscInt i, j; 1457 1458 PetscCheck(!mat->factortype, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 1459 for (i = 0; i < m; i++) { 1460 for (j = 0; j < n; j++) { 1461 if (mat->erroriffailure && PetscIsInfOrNanScalar(v[i * n + j])) 1462 #if defined(PETSC_USE_COMPLEX) 1463 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FP, "Inserting %g+i%g at matrix entry (%" PetscInt_FMT ",%" PetscInt_FMT ")", (double)PetscRealPart(v[i * n + j]), (double)PetscImaginaryPart(v[i * n + j]), idxm[i], idxn[j]); 1464 #else 1465 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FP, "Inserting %g at matrix entry (%" PetscInt_FMT ",%" PetscInt_FMT ")", (double)v[i * n + j], idxm[i], idxn[j]); 1466 #endif 1467 } 1468 } 1469 for (i = 0; i < m; i++) PetscCheck(idxm[i] < mat->rmap->N, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cannot insert in row %" PetscInt_FMT ", maximum is %" PetscInt_FMT, idxm[i], mat->rmap->N - 1); 1470 for (i = 0; i < n; i++) PetscCheck(idxn[i] < mat->cmap->N, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cannot insert in column %" PetscInt_FMT ", maximum is %" PetscInt_FMT, idxn[i], mat->cmap->N - 1); 1471 } 1472 1473 if (mat->assembled) { 1474 mat->was_assembled = PETSC_TRUE; 1475 mat->assembled = PETSC_FALSE; 1476 } 1477 PetscCall(PetscLogEventBegin(MAT_SetValues, mat, 0, 0, 0)); 1478 PetscUseTypeMethod(mat, setvalues, m, idxm, n, idxn, v, addv); 1479 PetscCall(PetscLogEventEnd(MAT_SetValues, mat, 0, 0, 0)); 1480 PetscFunctionReturn(PETSC_SUCCESS); 1481 } 1482 1483 /*@C 1484 MatSetValuesIS - Inserts or adds a block of values into a matrix using an `IS` to indicate the rows and columns 1485 These values may be cached, so `MatAssemblyBegin()` and `MatAssemblyEnd()` 1486 MUST be called after all calls to `MatSetValues()` have been completed. 1487 1488 Not Collective 1489 1490 Input Parameters: 1491 + mat - the matrix 1492 . v - a logically two-dimensional array of values 1493 . ism - the rows to provide 1494 . isn - the columns to provide 1495 - addv - either `ADD_VALUES` to add values to any existing entries, or `INSERT_VALUES` to replace existing entries with new values 1496 1497 Level: beginner 1498 1499 Notes: 1500 If you create the matrix yourself (that is not with a call to `DMCreateMatrix()`) then you MUST call MatXXXXSetPreallocation() or 1501 `MatSetUp()` before using this routine 1502 1503 By default the values, v, are row-oriented. See `MatSetOption()` for other options. 1504 1505 Calls to `MatSetValues()` with the `INSERT_VALUES` and `ADD_VALUES` 1506 options cannot be mixed without intervening calls to the assembly 1507 routines. 1508 1509 `MatSetValues()` uses 0-based row and column numbers in Fortran 1510 as well as in C. 1511 1512 Negative indices may be passed in ism and isn, these rows and columns are 1513 simply ignored. This allows easily inserting element stiffness matrices 1514 with homogeneous Dirchlet boundary conditions that you don't want represented 1515 in the matrix. 1516 1517 Efficiency Alert: 1518 The routine `MatSetValuesBlocked()` may offer much better efficiency 1519 for users of block sparse formats (`MATSEQBAIJ` and `MATMPIBAIJ`). 1520 1521 This is currently not optimized for any particular `ISType` 1522 1523 Developer Notes: 1524 This is labeled with C so does not automatically generate Fortran stubs and interfaces 1525 because it requires multiple Fortran interfaces depending on which arguments are scalar or arrays. 1526 1527 .seealso: `Mat`, `MatSetOption()`, `MatSetValues()`, `MatAssemblyBegin()`, `MatAssemblyEnd()`, `MatSetValuesBlocked()`, `MatSetValuesLocal()`, 1528 `InsertMode`, `INSERT_VALUES`, `ADD_VALUES`, `MatSetValues()` 1529 @*/ 1530 PetscErrorCode MatSetValuesIS(Mat mat, IS ism, IS isn, const PetscScalar v[], InsertMode addv) 1531 { 1532 PetscInt m, n; 1533 const PetscInt *rows, *cols; 1534 1535 PetscFunctionBeginHot; 1536 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 1537 PetscCall(ISGetIndices(ism, &rows)); 1538 PetscCall(ISGetIndices(isn, &cols)); 1539 PetscCall(ISGetLocalSize(ism, &m)); 1540 PetscCall(ISGetLocalSize(isn, &n)); 1541 PetscCall(MatSetValues(mat, m, rows, n, cols, v, addv)); 1542 PetscCall(ISRestoreIndices(ism, &rows)); 1543 PetscCall(ISRestoreIndices(isn, &cols)); 1544 PetscFunctionReturn(PETSC_SUCCESS); 1545 } 1546 1547 /*@ 1548 MatSetValuesRowLocal - Inserts a row (block row for `MATBAIJ` matrices) of nonzero 1549 values into a matrix 1550 1551 Not Collective 1552 1553 Input Parameters: 1554 + mat - the matrix 1555 . row - the (block) row to set 1556 - v - a logically two-dimensional array of values 1557 1558 Level: intermediate 1559 1560 Notes: 1561 The values, `v`, are column-oriented (for the block version) and sorted 1562 1563 All the nonzeros in the row must be provided 1564 1565 The matrix must have previously had its column indices set 1566 1567 The row must belong to this process 1568 1569 .seealso: `Mat`, `MatSetOption()`, `MatAssemblyBegin()`, `MatAssemblyEnd()`, `MatSetValuesBlocked()`, `MatSetValuesLocal()`, 1570 `InsertMode`, `INSERT_VALUES`, `ADD_VALUES`, `MatSetValues()`, `MatSetValuesRow()`, `MatSetLocalToGlobalMapping()` 1571 @*/ 1572 PetscErrorCode MatSetValuesRowLocal(Mat mat, PetscInt row, const PetscScalar v[]) 1573 { 1574 PetscInt globalrow; 1575 1576 PetscFunctionBegin; 1577 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 1578 PetscValidType(mat, 1); 1579 PetscValidScalarPointer(v, 3); 1580 PetscCall(ISLocalToGlobalMappingApply(mat->rmap->mapping, 1, &row, &globalrow)); 1581 PetscCall(MatSetValuesRow(mat, globalrow, v)); 1582 PetscFunctionReturn(PETSC_SUCCESS); 1583 } 1584 1585 /*@ 1586 MatSetValuesRow - Inserts a row (block row for `MATBAIJ` matrices) of nonzero 1587 values into a matrix 1588 1589 Not Collective 1590 1591 Input Parameters: 1592 + mat - the matrix 1593 . row - the (block) row to set 1594 - v - a logically two-dimensional (column major) array of values for block matrices with blocksize larger than one, otherwise a one dimensional array of values 1595 1596 Level: advanced 1597 1598 Notes: 1599 The values, `v`, are column-oriented for the block version. 1600 1601 All the nonzeros in the row must be provided 1602 1603 THE MATRIX MUST HAVE PREVIOUSLY HAD ITS COLUMN INDICES SET. IT IS RARE THAT THIS ROUTINE IS USED, usually `MatSetValues()` is used. 1604 1605 The row must belong to this process 1606 1607 .seealso: `Mat`, `MatSetValues()`, `MatSetOption()`, `MatAssemblyBegin()`, `MatAssemblyEnd()`, `MatSetValuesBlocked()`, `MatSetValuesLocal()`, 1608 `InsertMode`, `INSERT_VALUES`, `ADD_VALUES`, `MatSetValues()` 1609 @*/ 1610 PetscErrorCode MatSetValuesRow(Mat mat, PetscInt row, const PetscScalar v[]) 1611 { 1612 PetscFunctionBeginHot; 1613 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 1614 PetscValidType(mat, 1); 1615 MatCheckPreallocated(mat, 1); 1616 PetscValidScalarPointer(v, 3); 1617 PetscCheck(mat->insertmode != ADD_VALUES, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Cannot mix add and insert values"); 1618 PetscCheck(!mat->factortype, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 1619 mat->insertmode = INSERT_VALUES; 1620 1621 if (mat->assembled) { 1622 mat->was_assembled = PETSC_TRUE; 1623 mat->assembled = PETSC_FALSE; 1624 } 1625 PetscCall(PetscLogEventBegin(MAT_SetValues, mat, 0, 0, 0)); 1626 PetscUseTypeMethod(mat, setvaluesrow, row, v); 1627 PetscCall(PetscLogEventEnd(MAT_SetValues, mat, 0, 0, 0)); 1628 PetscFunctionReturn(PETSC_SUCCESS); 1629 } 1630 1631 /*@ 1632 MatSetValuesStencil - Inserts or adds a block of values into a matrix. 1633 Using structured grid indexing 1634 1635 Not Collective 1636 1637 Input Parameters: 1638 + mat - the matrix 1639 . m - number of rows being entered 1640 . idxm - grid coordinates (and component number when dof > 1) for matrix rows being entered 1641 . n - number of columns being entered 1642 . idxn - grid coordinates (and component number when dof > 1) for matrix columns being entered 1643 . v - a logically two-dimensional array of values 1644 - addv - either `ADD_VALUES` to add to existing entries at that location or `INSERT_VALUES` to replace existing entries with new values 1645 1646 Level: beginner 1647 1648 Notes: 1649 By default the values, `v`, are row-oriented. See `MatSetOption()` for other options. 1650 1651 Calls to `MatSetValuesStencil()` with the `INSERT_VALUES` and `ADD_VALUES` 1652 options cannot be mixed without intervening calls to the assembly 1653 routines. 1654 1655 The grid coordinates are across the entire grid, not just the local portion 1656 1657 `MatSetValuesStencil()` uses 0-based row and column numbers in Fortran 1658 as well as in C. 1659 1660 For setting/accessing vector values via array coordinates you can use the `DMDAVecGetArray()` routine 1661 1662 In order to use this routine you must either obtain the matrix with `DMCreateMatrix()` 1663 or call `MatSetLocalToGlobalMapping()` and `MatSetStencil()` first. 1664 1665 The columns and rows in the stencil passed in MUST be contained within the 1666 ghost region of the given process as set with DMDACreateXXX() or `MatSetStencil()`. For example, 1667 if you create a `DMDA` with an overlap of one grid level and on a particular process its first 1668 local nonghost x logical coordinate is 6 (so its first ghost x logical coordinate is 5) the 1669 first i index you can use in your column and row indices in `MatSetStencil()` is 5. 1670 1671 For periodic boundary conditions use negative indices for values to the left (below 0; that are to be 1672 obtained by wrapping values from right edge). For values to the right of the last entry using that index plus one 1673 etc to obtain values that obtained by wrapping the values from the left edge. This does not work for anything but the 1674 `DM_BOUNDARY_PERIODIC` boundary type. 1675 1676 For indices that don't mean anything for your case (like the k index when working in 2d) or the c index when you have 1677 a single value per point) you can skip filling those indices. 1678 1679 Inspired by the structured grid interface to the HYPRE package 1680 (https://computation.llnl.gov/projects/hypre-scalable-linear-solvers-multigrid-methods) 1681 1682 Efficiency Alert: 1683 The routine `MatSetValuesBlockedStencil()` may offer much better efficiency 1684 for users of block sparse formats (`MATSEQBAIJ` and `MATMPIBAIJ`). 1685 1686 Fortran Note: 1687 In Fortran idxm and idxn should be declared as 1688 $ MatStencil idxm(4,m),idxn(4,n) 1689 and the values inserted using 1690 .vb 1691 idxm(MatStencil_i,1) = i 1692 idxm(MatStencil_j,1) = j 1693 idxm(MatStencil_k,1) = k 1694 idxm(MatStencil_c,1) = c 1695 etc 1696 .ve 1697 1698 .seealso: `Mat`, `DMDA`, `MatSetOption()`, `MatAssemblyBegin()`, `MatAssemblyEnd()`, `MatSetValuesBlocked()`, `MatSetValuesLocal()` 1699 `MatSetValues()`, `MatSetValuesBlockedStencil()`, `MatSetStencil()`, `DMCreateMatrix()`, `DMDAVecGetArray()`, `MatStencil` 1700 @*/ 1701 PetscErrorCode MatSetValuesStencil(Mat mat, PetscInt m, const MatStencil idxm[], PetscInt n, const MatStencil idxn[], const PetscScalar v[], InsertMode addv) 1702 { 1703 PetscInt buf[8192], *bufm = NULL, *bufn = NULL, *jdxm, *jdxn; 1704 PetscInt j, i, dim = mat->stencil.dim, *dims = mat->stencil.dims + 1, tmp; 1705 PetscInt *starts = mat->stencil.starts, *dxm = (PetscInt *)idxm, *dxn = (PetscInt *)idxn, sdim = dim - (1 - (PetscInt)mat->stencil.noc); 1706 1707 PetscFunctionBegin; 1708 if (!m || !n) PetscFunctionReturn(PETSC_SUCCESS); /* no values to insert */ 1709 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 1710 PetscValidType(mat, 1); 1711 PetscValidPointer(idxm, 3); 1712 PetscValidPointer(idxn, 5); 1713 1714 if ((m + n) <= (PetscInt)(sizeof(buf) / sizeof(PetscInt))) { 1715 jdxm = buf; 1716 jdxn = buf + m; 1717 } else { 1718 PetscCall(PetscMalloc2(m, &bufm, n, &bufn)); 1719 jdxm = bufm; 1720 jdxn = bufn; 1721 } 1722 for (i = 0; i < m; i++) { 1723 for (j = 0; j < 3 - sdim; j++) dxm++; 1724 tmp = *dxm++ - starts[0]; 1725 for (j = 0; j < dim - 1; j++) { 1726 if ((*dxm++ - starts[j + 1]) < 0 || tmp < 0) tmp = -1; 1727 else tmp = tmp * dims[j] + *(dxm - 1) - starts[j + 1]; 1728 } 1729 if (mat->stencil.noc) dxm++; 1730 jdxm[i] = tmp; 1731 } 1732 for (i = 0; i < n; i++) { 1733 for (j = 0; j < 3 - sdim; j++) dxn++; 1734 tmp = *dxn++ - starts[0]; 1735 for (j = 0; j < dim - 1; j++) { 1736 if ((*dxn++ - starts[j + 1]) < 0 || tmp < 0) tmp = -1; 1737 else tmp = tmp * dims[j] + *(dxn - 1) - starts[j + 1]; 1738 } 1739 if (mat->stencil.noc) dxn++; 1740 jdxn[i] = tmp; 1741 } 1742 PetscCall(MatSetValuesLocal(mat, m, jdxm, n, jdxn, v, addv)); 1743 PetscCall(PetscFree2(bufm, bufn)); 1744 PetscFunctionReturn(PETSC_SUCCESS); 1745 } 1746 1747 /*@ 1748 MatSetValuesBlockedStencil - Inserts or adds a block of values into a matrix. 1749 Using structured grid indexing 1750 1751 Not Collective 1752 1753 Input Parameters: 1754 + mat - the matrix 1755 . m - number of rows being entered 1756 . idxm - grid coordinates for matrix rows being entered 1757 . n - number of columns being entered 1758 . idxn - grid coordinates for matrix columns being entered 1759 . v - a logically two-dimensional array of values 1760 - addv - either `ADD_VALUES` to add to existing entries or `INSERT_VALUES` to replace existing entries with new values 1761 1762 Level: beginner 1763 1764 Notes: 1765 By default the values, `v`, are row-oriented and unsorted. 1766 See `MatSetOption()` for other options. 1767 1768 Calls to `MatSetValuesBlockedStencil()` with the `INSERT_VALUES` and `ADD_VALUES` 1769 options cannot be mixed without intervening calls to the assembly 1770 routines. 1771 1772 The grid coordinates are across the entire grid, not just the local portion 1773 1774 `MatSetValuesBlockedStencil()` uses 0-based row and column numbers in Fortran 1775 as well as in C. 1776 1777 For setting/accessing vector values via array coordinates you can use the `DMDAVecGetArray()` routine 1778 1779 In order to use this routine you must either obtain the matrix with `DMCreateMatrix()` 1780 or call `MatSetBlockSize()`, `MatSetLocalToGlobalMapping()` and `MatSetStencil()` first. 1781 1782 The columns and rows in the stencil passed in MUST be contained within the 1783 ghost region of the given process as set with DMDACreateXXX() or `MatSetStencil()`. For example, 1784 if you create a `DMDA` with an overlap of one grid level and on a particular process its first 1785 local nonghost x logical coordinate is 6 (so its first ghost x logical coordinate is 5) the 1786 first i index you can use in your column and row indices in `MatSetStencil()` is 5. 1787 1788 Negative indices may be passed in idxm and idxn, these rows and columns are 1789 simply ignored. This allows easily inserting element stiffness matrices 1790 with homogeneous Dirchlet boundary conditions that you don't want represented 1791 in the matrix. 1792 1793 Inspired by the structured grid interface to the HYPRE package 1794 (https://computation.llnl.gov/projects/hypre-scalable-linear-solvers-multigrid-methods) 1795 1796 Fortran Note: 1797 In Fortran idxm and idxn should be declared as 1798 $ MatStencil idxm(4,m),idxn(4,n) 1799 and the values inserted using 1800 .vb 1801 idxm(MatStencil_i,1) = i 1802 idxm(MatStencil_j,1) = j 1803 idxm(MatStencil_k,1) = k 1804 etc 1805 .ve 1806 1807 .seealso: `Mat`, `DMDA`, `MatSetOption()`, `MatAssemblyBegin()`, `MatAssemblyEnd()`, `MatSetValuesBlocked()`, `MatSetValuesLocal()` 1808 `MatSetValues()`, `MatSetValuesStencil()`, `MatSetStencil()`, `DMCreateMatrix()`, `DMDAVecGetArray()`, `MatStencil`, 1809 `MatSetBlockSize()`, `MatSetLocalToGlobalMapping()` 1810 @*/ 1811 PetscErrorCode MatSetValuesBlockedStencil(Mat mat, PetscInt m, const MatStencil idxm[], PetscInt n, const MatStencil idxn[], const PetscScalar v[], InsertMode addv) 1812 { 1813 PetscInt buf[8192], *bufm = NULL, *bufn = NULL, *jdxm, *jdxn; 1814 PetscInt j, i, dim = mat->stencil.dim, *dims = mat->stencil.dims + 1, tmp; 1815 PetscInt *starts = mat->stencil.starts, *dxm = (PetscInt *)idxm, *dxn = (PetscInt *)idxn, sdim = dim - (1 - (PetscInt)mat->stencil.noc); 1816 1817 PetscFunctionBegin; 1818 if (!m || !n) PetscFunctionReturn(PETSC_SUCCESS); /* no values to insert */ 1819 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 1820 PetscValidType(mat, 1); 1821 PetscValidPointer(idxm, 3); 1822 PetscValidPointer(idxn, 5); 1823 PetscValidScalarPointer(v, 6); 1824 1825 if ((m + n) <= (PetscInt)(sizeof(buf) / sizeof(PetscInt))) { 1826 jdxm = buf; 1827 jdxn = buf + m; 1828 } else { 1829 PetscCall(PetscMalloc2(m, &bufm, n, &bufn)); 1830 jdxm = bufm; 1831 jdxn = bufn; 1832 } 1833 for (i = 0; i < m; i++) { 1834 for (j = 0; j < 3 - sdim; j++) dxm++; 1835 tmp = *dxm++ - starts[0]; 1836 for (j = 0; j < sdim - 1; j++) { 1837 if ((*dxm++ - starts[j + 1]) < 0 || tmp < 0) tmp = -1; 1838 else tmp = tmp * dims[j] + *(dxm - 1) - starts[j + 1]; 1839 } 1840 dxm++; 1841 jdxm[i] = tmp; 1842 } 1843 for (i = 0; i < n; i++) { 1844 for (j = 0; j < 3 - sdim; j++) dxn++; 1845 tmp = *dxn++ - starts[0]; 1846 for (j = 0; j < sdim - 1; j++) { 1847 if ((*dxn++ - starts[j + 1]) < 0 || tmp < 0) tmp = -1; 1848 else tmp = tmp * dims[j] + *(dxn - 1) - starts[j + 1]; 1849 } 1850 dxn++; 1851 jdxn[i] = tmp; 1852 } 1853 PetscCall(MatSetValuesBlockedLocal(mat, m, jdxm, n, jdxn, v, addv)); 1854 PetscCall(PetscFree2(bufm, bufn)); 1855 PetscFunctionReturn(PETSC_SUCCESS); 1856 } 1857 1858 /*@ 1859 MatSetStencil - Sets the grid information for setting values into a matrix via 1860 `MatSetValuesStencil()` 1861 1862 Not Collective 1863 1864 Input Parameters: 1865 + mat - the matrix 1866 . dim - dimension of the grid 1, 2, or 3 1867 . dims - number of grid points in x, y, and z direction, including ghost points on your processor 1868 . starts - starting point of ghost nodes on your processor in x, y, and z direction 1869 - dof - number of degrees of freedom per node 1870 1871 Level: beginner 1872 1873 Notes: 1874 Inspired by the structured grid interface to the HYPRE package 1875 (www.llnl.gov/CASC/hyper) 1876 1877 For matrices generated with `DMCreateMatrix()` this routine is automatically called and so not needed by the 1878 user. 1879 1880 .seealso: `Mat`, `MatStencil`, `MatSetOption()`, `MatAssemblyBegin()`, `MatAssemblyEnd()`, `MatSetValuesBlocked()`, `MatSetValuesLocal()` 1881 `MatSetValues()`, `MatSetValuesBlockedStencil()`, `MatSetValuesStencil()` 1882 @*/ 1883 PetscErrorCode MatSetStencil(Mat mat, PetscInt dim, const PetscInt dims[], const PetscInt starts[], PetscInt dof) 1884 { 1885 PetscFunctionBegin; 1886 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 1887 PetscValidIntPointer(dims, 3); 1888 PetscValidIntPointer(starts, 4); 1889 1890 mat->stencil.dim = dim + (dof > 1); 1891 for (PetscInt i = 0; i < dim; i++) { 1892 mat->stencil.dims[i] = dims[dim - i - 1]; /* copy the values in backwards */ 1893 mat->stencil.starts[i] = starts[dim - i - 1]; 1894 } 1895 mat->stencil.dims[dim] = dof; 1896 mat->stencil.starts[dim] = 0; 1897 mat->stencil.noc = (PetscBool)(dof == 1); 1898 PetscFunctionReturn(PETSC_SUCCESS); 1899 } 1900 1901 /*@C 1902 MatSetValuesBlocked - Inserts or adds a block of values into a matrix. 1903 1904 Not Collective 1905 1906 Input Parameters: 1907 + mat - the matrix 1908 . v - a logically two-dimensional array of values 1909 . m - the number of block rows 1910 . idxm - the global block indices 1911 . n - the number of block columns 1912 . idxn - the global block indices 1913 - addv - either `ADD_VALUES` to add values to any existing entries, or `INSERT_VALUES` replaces existing entries with new values 1914 1915 Level: intermediate 1916 1917 Notes: 1918 If you create the matrix yourself (that is not with a call to `DMCreateMatrix()`) then you MUST call 1919 MatXXXXSetPreallocation() or `MatSetUp()` before using this routine. 1920 1921 The `m` and `n` count the NUMBER of blocks in the row direction and column direction, 1922 NOT the total number of rows/columns; for example, if the block size is 2 and 1923 you are passing in values for rows 2,3,4,5 then m would be 2 (not 4). 1924 The values in idxm would be 1 2; that is the first index for each block divided by 1925 the block size. 1926 1927 Note that you must call `MatSetBlockSize()` when constructing this matrix (before 1928 preallocating it). 1929 1930 By default the values, `v`, are row-oriented, so the layout of 1931 `v` is the same as for `MatSetValues()`. See `MatSetOption()` for other options. 1932 1933 Calls to `MatSetValuesBlocked()` with the `INSERT_VALUES` and `ADD_VALUES` 1934 options cannot be mixed without intervening calls to the assembly 1935 routines. 1936 1937 `MatSetValuesBlocked()` uses 0-based row and column numbers in Fortran 1938 as well as in C. 1939 1940 Negative indices may be passed in `idxm` and `idxn`, these rows and columns are 1941 simply ignored. This allows easily inserting element stiffness matrices 1942 with homogeneous Dirchlet boundary conditions that you don't want represented 1943 in the matrix. 1944 1945 Each time an entry is set within a sparse matrix via `MatSetValues()`, 1946 internal searching must be done to determine where to place the 1947 data in the matrix storage space. By instead inserting blocks of 1948 entries via `MatSetValuesBlocked()`, the overhead of matrix assembly is 1949 reduced. 1950 1951 Example: 1952 .vb 1953 Suppose m=n=2 and block size(bs) = 2 The array is 1954 1955 1 2 | 3 4 1956 5 6 | 7 8 1957 - - - | - - - 1958 9 10 | 11 12 1959 13 14 | 15 16 1960 1961 v[] should be passed in like 1962 v[] = [1,2,3,4,5,6,7,8,9,10,11,12,13,14,15,16] 1963 1964 If you are not using row oriented storage of v (that is you called MatSetOption(mat,MAT_ROW_ORIENTED,PETSC_FALSE)) then 1965 v[] = [1,5,9,13,2,6,10,14,3,7,11,15,4,8,12,16] 1966 .ve 1967 1968 .seealso: `Mat`, `MatSetBlockSize()`, `MatSetOption()`, `MatAssemblyBegin()`, `MatAssemblyEnd()`, `MatSetValues()`, `MatSetValuesBlockedLocal()` 1969 @*/ 1970 PetscErrorCode MatSetValuesBlocked(Mat mat, PetscInt m, const PetscInt idxm[], PetscInt n, const PetscInt idxn[], const PetscScalar v[], InsertMode addv) 1971 { 1972 PetscFunctionBeginHot; 1973 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 1974 PetscValidType(mat, 1); 1975 if (!m || !n) PetscFunctionReturn(PETSC_SUCCESS); /* no values to insert */ 1976 PetscValidIntPointer(idxm, 3); 1977 PetscValidIntPointer(idxn, 5); 1978 MatCheckPreallocated(mat, 1); 1979 if (mat->insertmode == NOT_SET_VALUES) mat->insertmode = addv; 1980 else PetscCheck(mat->insertmode == addv, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Cannot mix add values and insert values"); 1981 if (PetscDefined(USE_DEBUG)) { 1982 PetscCheck(!mat->factortype, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 1983 PetscCheck(mat->ops->setvaluesblocked || mat->ops->setvalues, PETSC_COMM_SELF, PETSC_ERR_SUP, "Mat type %s", ((PetscObject)mat)->type_name); 1984 } 1985 if (PetscDefined(USE_DEBUG)) { 1986 PetscInt rbs, cbs, M, N, i; 1987 PetscCall(MatGetBlockSizes(mat, &rbs, &cbs)); 1988 PetscCall(MatGetSize(mat, &M, &N)); 1989 for (i = 0; i < m; i++) PetscCheck(idxm[i] * rbs < M, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Row block index %" PetscInt_FMT " (index %" PetscInt_FMT ") greater than row length %" PetscInt_FMT, i, idxm[i], M); 1990 for (i = 0; i < n; i++) PetscCheck(idxn[i] * cbs < N, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Column block index %" PetscInt_FMT " (index %" PetscInt_FMT ") great than column length %" PetscInt_FMT, i, idxn[i], N); 1991 } 1992 if (mat->assembled) { 1993 mat->was_assembled = PETSC_TRUE; 1994 mat->assembled = PETSC_FALSE; 1995 } 1996 PetscCall(PetscLogEventBegin(MAT_SetValues, mat, 0, 0, 0)); 1997 if (mat->ops->setvaluesblocked) { 1998 PetscUseTypeMethod(mat, setvaluesblocked, m, idxm, n, idxn, v, addv); 1999 } else { 2000 PetscInt buf[8192], *bufr = NULL, *bufc = NULL, *iidxm, *iidxn; 2001 PetscInt i, j, bs, cbs; 2002 2003 PetscCall(MatGetBlockSizes(mat, &bs, &cbs)); 2004 if (m * bs + n * cbs <= (PetscInt)(sizeof(buf) / sizeof(PetscInt))) { 2005 iidxm = buf; 2006 iidxn = buf + m * bs; 2007 } else { 2008 PetscCall(PetscMalloc2(m * bs, &bufr, n * cbs, &bufc)); 2009 iidxm = bufr; 2010 iidxn = bufc; 2011 } 2012 for (i = 0; i < m; i++) { 2013 for (j = 0; j < bs; j++) iidxm[i * bs + j] = bs * idxm[i] + j; 2014 } 2015 if (m != n || bs != cbs || idxm != idxn) { 2016 for (i = 0; i < n; i++) { 2017 for (j = 0; j < cbs; j++) iidxn[i * cbs + j] = cbs * idxn[i] + j; 2018 } 2019 } else iidxn = iidxm; 2020 PetscCall(MatSetValues(mat, m * bs, iidxm, n * cbs, iidxn, v, addv)); 2021 PetscCall(PetscFree2(bufr, bufc)); 2022 } 2023 PetscCall(PetscLogEventEnd(MAT_SetValues, mat, 0, 0, 0)); 2024 PetscFunctionReturn(PETSC_SUCCESS); 2025 } 2026 2027 /*@C 2028 MatGetValues - Gets a block of local values from a matrix. 2029 2030 Not Collective; can only return values that are owned by the give process 2031 2032 Input Parameters: 2033 + mat - the matrix 2034 . v - a logically two-dimensional array for storing the values 2035 . m - the number of rows 2036 . idxm - the global indices of the rows 2037 . n - the number of columns 2038 - idxn - the global indices of the columns 2039 2040 Level: advanced 2041 2042 Notes: 2043 The user must allocate space (m*n `PetscScalar`s) for the values, `v`. 2044 The values, `v`, are then returned in a row-oriented format, 2045 analogous to that used by default in `MatSetValues()`. 2046 2047 `MatGetValues()` uses 0-based row and column numbers in 2048 Fortran as well as in C. 2049 2050 `MatGetValues()` requires that the matrix has been assembled 2051 with `MatAssemblyBegin()`/`MatAssemblyEnd()`. Thus, calls to 2052 `MatSetValues()` and `MatGetValues()` CANNOT be made in succession 2053 without intermediate matrix assembly. 2054 2055 Negative row or column indices will be ignored and those locations in `v` will be 2056 left unchanged. 2057 2058 For the standard row-based matrix formats, `idxm` can only contain rows owned by the requesting MPI rank. 2059 That is, rows with global index greater than or equal to rstart and less than rend where rstart and rend are obtainable 2060 from `MatGetOwnershipRange`(mat,&rstart,&rend). 2061 2062 .seealso: `Mat`, `MatGetRow()`, `MatCreateSubMatrices()`, `MatSetValues()`, `MatGetOwnershipRange()`, `MatGetValuesLocal()`, `MatGetValue()` 2063 @*/ 2064 PetscErrorCode MatGetValues(Mat mat, PetscInt m, const PetscInt idxm[], PetscInt n, const PetscInt idxn[], PetscScalar v[]) 2065 { 2066 PetscFunctionBegin; 2067 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 2068 PetscValidType(mat, 1); 2069 if (!m || !n) PetscFunctionReturn(PETSC_SUCCESS); 2070 PetscValidIntPointer(idxm, 3); 2071 PetscValidIntPointer(idxn, 5); 2072 PetscValidScalarPointer(v, 6); 2073 PetscCheck(mat->assembled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 2074 PetscCheck(!mat->factortype, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 2075 MatCheckPreallocated(mat, 1); 2076 2077 PetscCall(PetscLogEventBegin(MAT_GetValues, mat, 0, 0, 0)); 2078 PetscUseTypeMethod(mat, getvalues, m, idxm, n, idxn, v); 2079 PetscCall(PetscLogEventEnd(MAT_GetValues, mat, 0, 0, 0)); 2080 PetscFunctionReturn(PETSC_SUCCESS); 2081 } 2082 2083 /*@C 2084 MatGetValuesLocal - retrieves values from certain locations in a matrix using the local numbering of the indices 2085 defined previously by `MatSetLocalToGlobalMapping()` 2086 2087 Not Collective 2088 2089 Input Parameters: 2090 + mat - the matrix 2091 . nrow - number of rows 2092 . irow - the row local indices 2093 . ncol - number of columns 2094 - icol - the column local indices 2095 2096 Output Parameter: 2097 . y - a logically two-dimensional array of values 2098 2099 Level: advanced 2100 2101 Notes: 2102 If you create the matrix yourself (that is not with a call to `DMCreateMatrix()`) then you MUST call `MatSetLocalToGlobalMapping()` before using this routine. 2103 2104 This routine can only return values that are owned by the requesting MPI rank. That is, for standard matrix formats, rows that, in the global numbering, 2105 are greater than or equal to rstart and less than rend where rstart and rend are obtainable from `MatGetOwnershipRange`(mat,&rstart,&rend). One can 2106 determine if the resulting global row associated with the local row r is owned by the requesting MPI rank by applying the `ISLocalToGlobalMapping` set 2107 with `MatSetLocalToGlobalMapping()`. 2108 2109 Developer Note: 2110 This is labelled with C so does not automatically generate Fortran stubs and interfaces 2111 because it requires multiple Fortran interfaces depending on which arguments are scalar or arrays. 2112 2113 .seealso: `Mat`, `MatAssemblyBegin()`, `MatAssemblyEnd()`, `MatSetValues()`, `MatSetLocalToGlobalMapping()`, 2114 `MatSetValuesLocal()`, `MatGetValues()` 2115 @*/ 2116 PetscErrorCode MatGetValuesLocal(Mat mat, PetscInt nrow, const PetscInt irow[], PetscInt ncol, const PetscInt icol[], PetscScalar y[]) 2117 { 2118 PetscFunctionBeginHot; 2119 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 2120 PetscValidType(mat, 1); 2121 MatCheckPreallocated(mat, 1); 2122 if (!nrow || !ncol) PetscFunctionReturn(PETSC_SUCCESS); /* no values to retrieve */ 2123 PetscValidIntPointer(irow, 3); 2124 PetscValidIntPointer(icol, 5); 2125 if (PetscDefined(USE_DEBUG)) { 2126 PetscCheck(!mat->factortype, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 2127 PetscCheck(mat->ops->getvalueslocal || mat->ops->getvalues, PETSC_COMM_SELF, PETSC_ERR_SUP, "Mat type %s", ((PetscObject)mat)->type_name); 2128 } 2129 PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 2130 PetscCall(PetscLogEventBegin(MAT_GetValues, mat, 0, 0, 0)); 2131 if (mat->ops->getvalueslocal) PetscUseTypeMethod(mat, getvalueslocal, nrow, irow, ncol, icol, y); 2132 else { 2133 PetscInt buf[8192], *bufr = NULL, *bufc = NULL, *irowm, *icolm; 2134 if ((nrow + ncol) <= (PetscInt)(sizeof(buf) / sizeof(PetscInt))) { 2135 irowm = buf; 2136 icolm = buf + nrow; 2137 } else { 2138 PetscCall(PetscMalloc2(nrow, &bufr, ncol, &bufc)); 2139 irowm = bufr; 2140 icolm = bufc; 2141 } 2142 PetscCheck(mat->rmap->mapping, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "MatGetValuesLocal() cannot proceed without local-to-global row mapping (See MatSetLocalToGlobalMapping())."); 2143 PetscCheck(mat->cmap->mapping, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "MatGetValuesLocal() cannot proceed without local-to-global column mapping (See MatSetLocalToGlobalMapping())."); 2144 PetscCall(ISLocalToGlobalMappingApply(mat->rmap->mapping, nrow, irow, irowm)); 2145 PetscCall(ISLocalToGlobalMappingApply(mat->cmap->mapping, ncol, icol, icolm)); 2146 PetscCall(MatGetValues(mat, nrow, irowm, ncol, icolm, y)); 2147 PetscCall(PetscFree2(bufr, bufc)); 2148 } 2149 PetscCall(PetscLogEventEnd(MAT_GetValues, mat, 0, 0, 0)); 2150 PetscFunctionReturn(PETSC_SUCCESS); 2151 } 2152 2153 /*@ 2154 MatSetValuesBatch - Adds (`ADD_VALUES`) many blocks of values into a matrix at once. The blocks must all be square and 2155 the same size. Currently, this can only be called once and creates the given matrix. 2156 2157 Not Collective 2158 2159 Input Parameters: 2160 + mat - the matrix 2161 . nb - the number of blocks 2162 . bs - the number of rows (and columns) in each block 2163 . rows - a concatenation of the rows for each block 2164 - v - a concatenation of logically two-dimensional arrays of values 2165 2166 Level: advanced 2167 2168 Note: 2169 `MatSetPreallocationCOO()` and `MatSetValuesCOO()` may be a better way to provide the values 2170 2171 In the future, we will extend this routine to handle rectangular blocks, and to allow multiple calls for a given matrix. 2172 2173 .seealso: `Mat`, `Mat`MatSetOption()`, `MatAssemblyBegin()`, `MatAssemblyEnd()`, `MatSetValuesBlocked()`, `MatSetValuesLocal()`, 2174 `InsertMode`, `INSERT_VALUES`, `ADD_VALUES`, `MatSetValues()`, `MatSetPreallocationCOO()`, `MatSetValuesCOO()` 2175 @*/ 2176 PetscErrorCode MatSetValuesBatch(Mat mat, PetscInt nb, PetscInt bs, PetscInt rows[], const PetscScalar v[]) 2177 { 2178 PetscFunctionBegin; 2179 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 2180 PetscValidType(mat, 1); 2181 PetscValidIntPointer(rows, 4); 2182 PetscValidScalarPointer(v, 5); 2183 PetscAssert(!mat->factortype, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 2184 2185 PetscCall(PetscLogEventBegin(MAT_SetValuesBatch, mat, 0, 0, 0)); 2186 if (mat->ops->setvaluesbatch) PetscUseTypeMethod(mat, setvaluesbatch, nb, bs, rows, v); 2187 else { 2188 for (PetscInt b = 0; b < nb; ++b) PetscCall(MatSetValues(mat, bs, &rows[b * bs], bs, &rows[b * bs], &v[b * bs * bs], ADD_VALUES)); 2189 } 2190 PetscCall(PetscLogEventEnd(MAT_SetValuesBatch, mat, 0, 0, 0)); 2191 PetscFunctionReturn(PETSC_SUCCESS); 2192 } 2193 2194 /*@ 2195 MatSetLocalToGlobalMapping - Sets a local-to-global numbering for use by 2196 the routine `MatSetValuesLocal()` to allow users to insert matrix entries 2197 using a local (per-processor) numbering. 2198 2199 Not Collective 2200 2201 Input Parameters: 2202 + x - the matrix 2203 . rmapping - row mapping created with `ISLocalToGlobalMappingCreate()` or `ISLocalToGlobalMappingCreateIS()` 2204 - cmapping - column mapping 2205 2206 Level: intermediate 2207 2208 Note: 2209 If the matrix is obtained with `DMCreateMatrix()` then this may already have been called on the matrix 2210 2211 .seealso: `Mat`, `DM`, `DMCreateMatrix()`, `MatGetLocalToGlobalMapping()`, `MatAssemblyBegin()`, `MatAssemblyEnd()`, `MatSetValues()`, `MatSetValuesLocal()`, `MatGetValuesLocal()` 2212 @*/ 2213 PetscErrorCode MatSetLocalToGlobalMapping(Mat x, ISLocalToGlobalMapping rmapping, ISLocalToGlobalMapping cmapping) 2214 { 2215 PetscFunctionBegin; 2216 PetscValidHeaderSpecific(x, MAT_CLASSID, 1); 2217 PetscValidType(x, 1); 2218 if (rmapping) PetscValidHeaderSpecific(rmapping, IS_LTOGM_CLASSID, 2); 2219 if (cmapping) PetscValidHeaderSpecific(cmapping, IS_LTOGM_CLASSID, 3); 2220 if (x->ops->setlocaltoglobalmapping) PetscUseTypeMethod(x, setlocaltoglobalmapping, rmapping, cmapping); 2221 else { 2222 PetscCall(PetscLayoutSetISLocalToGlobalMapping(x->rmap, rmapping)); 2223 PetscCall(PetscLayoutSetISLocalToGlobalMapping(x->cmap, cmapping)); 2224 } 2225 PetscFunctionReturn(PETSC_SUCCESS); 2226 } 2227 2228 /*@ 2229 MatGetLocalToGlobalMapping - Gets the local-to-global numbering set by `MatSetLocalToGlobalMapping()` 2230 2231 Not Collective 2232 2233 Input Parameter: 2234 . A - the matrix 2235 2236 Output Parameters: 2237 + rmapping - row mapping 2238 - cmapping - column mapping 2239 2240 Level: advanced 2241 2242 .seealso: `Mat`, `MatSetLocalToGlobalMapping()`, `MatSetValuesLocal()` 2243 @*/ 2244 PetscErrorCode MatGetLocalToGlobalMapping(Mat A, ISLocalToGlobalMapping *rmapping, ISLocalToGlobalMapping *cmapping) 2245 { 2246 PetscFunctionBegin; 2247 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 2248 PetscValidType(A, 1); 2249 if (rmapping) { 2250 PetscValidPointer(rmapping, 2); 2251 *rmapping = A->rmap->mapping; 2252 } 2253 if (cmapping) { 2254 PetscValidPointer(cmapping, 3); 2255 *cmapping = A->cmap->mapping; 2256 } 2257 PetscFunctionReturn(PETSC_SUCCESS); 2258 } 2259 2260 /*@ 2261 MatSetLayouts - Sets the `PetscLayout` objects for rows and columns of a matrix 2262 2263 Logically Collective 2264 2265 Input Parameters: 2266 + A - the matrix 2267 . rmap - row layout 2268 - cmap - column layout 2269 2270 Level: advanced 2271 2272 Note: 2273 The `PetscLayout` objects are usually created automatically for the matrix so this routine rarely needs to be called. 2274 2275 .seealso: `Mat`, `PetscLayout`, `MatCreateVecs()`, `MatGetLocalToGlobalMapping()`, `MatGetLayouts()` 2276 @*/ 2277 PetscErrorCode MatSetLayouts(Mat A, PetscLayout rmap, PetscLayout cmap) 2278 { 2279 PetscFunctionBegin; 2280 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 2281 PetscCall(PetscLayoutReference(rmap, &A->rmap)); 2282 PetscCall(PetscLayoutReference(cmap, &A->cmap)); 2283 PetscFunctionReturn(PETSC_SUCCESS); 2284 } 2285 2286 /*@ 2287 MatGetLayouts - Gets the `PetscLayout` objects for rows and columns 2288 2289 Not Collective 2290 2291 Input Parameter: 2292 . A - the matrix 2293 2294 Output Parameters: 2295 + rmap - row layout 2296 - cmap - column layout 2297 2298 Level: advanced 2299 2300 .seealso: `Mat`, [Matrix Layouts](sec_matlayout), `PetscLayout`, `MatCreateVecs()`, `MatGetLocalToGlobalMapping()`, `MatSetLayouts()` 2301 @*/ 2302 PetscErrorCode MatGetLayouts(Mat A, PetscLayout *rmap, PetscLayout *cmap) 2303 { 2304 PetscFunctionBegin; 2305 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 2306 PetscValidType(A, 1); 2307 if (rmap) { 2308 PetscValidPointer(rmap, 2); 2309 *rmap = A->rmap; 2310 } 2311 if (cmap) { 2312 PetscValidPointer(cmap, 3); 2313 *cmap = A->cmap; 2314 } 2315 PetscFunctionReturn(PETSC_SUCCESS); 2316 } 2317 2318 /*@C 2319 MatSetValuesLocal - Inserts or adds values into certain locations of a matrix, 2320 using a local numbering of the nodes. 2321 2322 Not Collective 2323 2324 Input Parameters: 2325 + mat - the matrix 2326 . nrow - number of rows 2327 . irow - the row local indices 2328 . ncol - number of columns 2329 . icol - the column local indices 2330 . y - a logically two-dimensional array of values 2331 - addv - either `INSERT_VALUES` to add values to any existing entries, or `INSERT_VALUES` to replace existing entries with new values 2332 2333 Level: intermediate 2334 2335 Notes: 2336 If you create the matrix yourself (that is not with a call to `DMCreateMatrix()`) then you MUST call MatXXXXSetPreallocation() or 2337 `MatSetUp()` before using this routine 2338 2339 If you create the matrix yourself (that is not with a call to `DMCreateMatrix()`) then you MUST call `MatSetLocalToGlobalMapping()` before using this routine 2340 2341 Calls to `MatSetValuesLocal()` with the `INSERT_VALUES` and `ADD_VALUES` 2342 options cannot be mixed without intervening calls to the assembly 2343 routines. 2344 2345 These values may be cached, so `MatAssemblyBegin()` and `MatAssemblyEnd()` 2346 MUST be called after all calls to `MatSetValuesLocal()` have been completed. 2347 2348 Developer Note: 2349 This is labeled with C so does not automatically generate Fortran stubs and interfaces 2350 because it requires multiple Fortran interfaces depending on which arguments are scalar or arrays. 2351 2352 .seealso: `Mat`, `MatAssemblyBegin()`, `MatAssemblyEnd()`, `MatSetValues()`, `MatSetLocalToGlobalMapping()`, 2353 `MatGetValuesLocal()` 2354 @*/ 2355 PetscErrorCode MatSetValuesLocal(Mat mat, PetscInt nrow, const PetscInt irow[], PetscInt ncol, const PetscInt icol[], const PetscScalar y[], InsertMode addv) 2356 { 2357 PetscFunctionBeginHot; 2358 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 2359 PetscValidType(mat, 1); 2360 MatCheckPreallocated(mat, 1); 2361 if (!nrow || !ncol) PetscFunctionReturn(PETSC_SUCCESS); /* no values to insert */ 2362 PetscValidIntPointer(irow, 3); 2363 PetscValidIntPointer(icol, 5); 2364 if (mat->insertmode == NOT_SET_VALUES) mat->insertmode = addv; 2365 else PetscCheck(mat->insertmode == addv, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Cannot mix add values and insert values"); 2366 if (PetscDefined(USE_DEBUG)) { 2367 PetscCheck(!mat->factortype, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 2368 PetscCheck(mat->ops->setvalueslocal || mat->ops->setvalues, PETSC_COMM_SELF, PETSC_ERR_SUP, "Mat type %s", ((PetscObject)mat)->type_name); 2369 } 2370 2371 if (mat->assembled) { 2372 mat->was_assembled = PETSC_TRUE; 2373 mat->assembled = PETSC_FALSE; 2374 } 2375 PetscCall(PetscLogEventBegin(MAT_SetValues, mat, 0, 0, 0)); 2376 if (mat->ops->setvalueslocal) PetscUseTypeMethod(mat, setvalueslocal, nrow, irow, ncol, icol, y, addv); 2377 else { 2378 PetscInt buf[8192], *bufr = NULL, *bufc = NULL; 2379 const PetscInt *irowm, *icolm; 2380 2381 if ((!mat->rmap->mapping && !mat->cmap->mapping) || (nrow + ncol) <= (PetscInt)(sizeof(buf) / sizeof(PetscInt))) { 2382 bufr = buf; 2383 bufc = buf + nrow; 2384 irowm = bufr; 2385 icolm = bufc; 2386 } else { 2387 PetscCall(PetscMalloc2(nrow, &bufr, ncol, &bufc)); 2388 irowm = bufr; 2389 icolm = bufc; 2390 } 2391 if (mat->rmap->mapping) PetscCall(ISLocalToGlobalMappingApply(mat->rmap->mapping, nrow, irow, bufr)); 2392 else irowm = irow; 2393 if (mat->cmap->mapping) { 2394 if (mat->cmap->mapping != mat->rmap->mapping || ncol != nrow || icol != irow) { 2395 PetscCall(ISLocalToGlobalMappingApply(mat->cmap->mapping, ncol, icol, bufc)); 2396 } else icolm = irowm; 2397 } else icolm = icol; 2398 PetscCall(MatSetValues(mat, nrow, irowm, ncol, icolm, y, addv)); 2399 if (bufr != buf) PetscCall(PetscFree2(bufr, bufc)); 2400 } 2401 PetscCall(PetscLogEventEnd(MAT_SetValues, mat, 0, 0, 0)); 2402 PetscFunctionReturn(PETSC_SUCCESS); 2403 } 2404 2405 /*@C 2406 MatSetValuesBlockedLocal - Inserts or adds values into certain locations of a matrix, 2407 using a local ordering of the nodes a block at a time. 2408 2409 Not Collective 2410 2411 Input Parameters: 2412 + x - the matrix 2413 . nrow - number of rows 2414 . irow - the row local indices 2415 . ncol - number of columns 2416 . icol - the column local indices 2417 . y - a logically two-dimensional array of values 2418 - addv - either `ADD_VALUES` to add values to any existing entries, or `INSERT_VALUES` to replace existing entries with new values 2419 2420 Level: intermediate 2421 2422 Notes: 2423 If you create the matrix yourself (that is not with a call to `DMCreateMatrix()`) then you MUST call MatXXXXSetPreallocation() or 2424 `MatSetUp()` before using this routine 2425 2426 If you create the matrix yourself (that is not with a call to `DMCreateMatrix()`) then you MUST call `MatSetBlockSize()` and `MatSetLocalToGlobalMapping()` 2427 before using this routineBefore calling `MatSetValuesLocal()`, the user must first set the 2428 2429 Calls to `MatSetValuesBlockedLocal()` with the `INSERT_VALUES` and `ADD_VALUES` 2430 options cannot be mixed without intervening calls to the assembly 2431 routines. 2432 2433 These values may be cached, so `MatAssemblyBegin()` and `MatAssemblyEnd()` 2434 MUST be called after all calls to `MatSetValuesBlockedLocal()` have been completed. 2435 2436 Developer Note: 2437 This is labeled with C so does not automatically generate Fortran stubs and interfaces 2438 because it requires multiple Fortran interfaces depending on which arguments are scalar or arrays. 2439 2440 .seealso: `Mat`, `MatSetBlockSize()`, `MatSetLocalToGlobalMapping()`, `MatAssemblyBegin()`, `MatAssemblyEnd()`, 2441 `MatSetValuesLocal()`, `MatSetValuesBlocked()` 2442 @*/ 2443 PetscErrorCode MatSetValuesBlockedLocal(Mat mat, PetscInt nrow, const PetscInt irow[], PetscInt ncol, const PetscInt icol[], const PetscScalar y[], InsertMode addv) 2444 { 2445 PetscFunctionBeginHot; 2446 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 2447 PetscValidType(mat, 1); 2448 MatCheckPreallocated(mat, 1); 2449 if (!nrow || !ncol) PetscFunctionReturn(PETSC_SUCCESS); /* no values to insert */ 2450 PetscValidIntPointer(irow, 3); 2451 PetscValidIntPointer(icol, 5); 2452 if (mat->insertmode == NOT_SET_VALUES) mat->insertmode = addv; 2453 else PetscCheck(mat->insertmode == addv, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Cannot mix add values and insert values"); 2454 if (PetscDefined(USE_DEBUG)) { 2455 PetscCheck(!mat->factortype, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 2456 PetscCheck(mat->ops->setvaluesblockedlocal || mat->ops->setvaluesblocked || mat->ops->setvalueslocal || mat->ops->setvalues, PETSC_COMM_SELF, PETSC_ERR_SUP, "Mat type %s", ((PetscObject)mat)->type_name); 2457 } 2458 2459 if (mat->assembled) { 2460 mat->was_assembled = PETSC_TRUE; 2461 mat->assembled = PETSC_FALSE; 2462 } 2463 if (PetscUnlikelyDebug(mat->rmap->mapping)) { /* Condition on the mapping existing, because MatSetValuesBlockedLocal_IS does not require it to be set. */ 2464 PetscInt irbs, rbs; 2465 PetscCall(MatGetBlockSizes(mat, &rbs, NULL)); 2466 PetscCall(ISLocalToGlobalMappingGetBlockSize(mat->rmap->mapping, &irbs)); 2467 PetscCheck(rbs == irbs, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Different row block sizes! mat %" PetscInt_FMT ", row l2g map %" PetscInt_FMT, rbs, irbs); 2468 } 2469 if (PetscUnlikelyDebug(mat->cmap->mapping)) { 2470 PetscInt icbs, cbs; 2471 PetscCall(MatGetBlockSizes(mat, NULL, &cbs)); 2472 PetscCall(ISLocalToGlobalMappingGetBlockSize(mat->cmap->mapping, &icbs)); 2473 PetscCheck(cbs == icbs, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Different col block sizes! mat %" PetscInt_FMT ", col l2g map %" PetscInt_FMT, cbs, icbs); 2474 } 2475 PetscCall(PetscLogEventBegin(MAT_SetValues, mat, 0, 0, 0)); 2476 if (mat->ops->setvaluesblockedlocal) PetscUseTypeMethod(mat, setvaluesblockedlocal, nrow, irow, ncol, icol, y, addv); 2477 else { 2478 PetscInt buf[8192], *bufr = NULL, *bufc = NULL; 2479 const PetscInt *irowm, *icolm; 2480 2481 if ((!mat->rmap->mapping && !mat->cmap->mapping) || (nrow + ncol) <= (PetscInt)(sizeof(buf) / sizeof(PetscInt))) { 2482 bufr = buf; 2483 bufc = buf + nrow; 2484 irowm = bufr; 2485 icolm = bufc; 2486 } else { 2487 PetscCall(PetscMalloc2(nrow, &bufr, ncol, &bufc)); 2488 irowm = bufr; 2489 icolm = bufc; 2490 } 2491 if (mat->rmap->mapping) PetscCall(ISLocalToGlobalMappingApplyBlock(mat->rmap->mapping, nrow, irow, bufr)); 2492 else irowm = irow; 2493 if (mat->cmap->mapping) { 2494 if (mat->cmap->mapping != mat->rmap->mapping || ncol != nrow || icol != irow) { 2495 PetscCall(ISLocalToGlobalMappingApplyBlock(mat->cmap->mapping, ncol, icol, bufc)); 2496 } else icolm = irowm; 2497 } else icolm = icol; 2498 PetscCall(MatSetValuesBlocked(mat, nrow, irowm, ncol, icolm, y, addv)); 2499 if (bufr != buf) PetscCall(PetscFree2(bufr, bufc)); 2500 } 2501 PetscCall(PetscLogEventEnd(MAT_SetValues, mat, 0, 0, 0)); 2502 PetscFunctionReturn(PETSC_SUCCESS); 2503 } 2504 2505 /*@ 2506 MatMultDiagonalBlock - Computes the matrix-vector product, y = Dx. Where D is defined by the inode or block structure of the diagonal 2507 2508 Collective 2509 2510 Input Parameters: 2511 + mat - the matrix 2512 - x - the vector to be multiplied 2513 2514 Output Parameters: 2515 . y - the result 2516 2517 Level: developer 2518 2519 Note: 2520 The vectors `x` and `y` cannot be the same. I.e., one cannot 2521 call `MatMultDiagonalBlock`(A,y,y). 2522 2523 .seealso: `Mat`, `MatMult()`, `MatMultTranspose()`, `MatMultAdd()`, `MatMultTransposeAdd()` 2524 @*/ 2525 PetscErrorCode MatMultDiagonalBlock(Mat mat, Vec x, Vec y) 2526 { 2527 PetscFunctionBegin; 2528 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 2529 PetscValidType(mat, 1); 2530 PetscValidHeaderSpecific(x, VEC_CLASSID, 2); 2531 PetscValidHeaderSpecific(y, VEC_CLASSID, 3); 2532 2533 PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 2534 PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 2535 PetscCheck(x != y, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "x and y must be different vectors"); 2536 MatCheckPreallocated(mat, 1); 2537 2538 PetscUseTypeMethod(mat, multdiagonalblock, x, y); 2539 PetscCall(PetscObjectStateIncrease((PetscObject)y)); 2540 PetscFunctionReturn(PETSC_SUCCESS); 2541 } 2542 2543 /*@ 2544 MatMult - Computes the matrix-vector product, y = Ax. 2545 2546 Neighbor-wise Collective 2547 2548 Input Parameters: 2549 + mat - the matrix 2550 - x - the vector to be multiplied 2551 2552 Output Parameters: 2553 . y - the result 2554 2555 Level: beginner 2556 2557 Note: 2558 The vectors `x` and `y` cannot be the same. I.e., one cannot 2559 call `MatMult`(A,y,y). 2560 2561 .seealso: `Mat`, `MatMultTranspose()`, `MatMultAdd()`, `MatMultTransposeAdd()` 2562 @*/ 2563 PetscErrorCode MatMult(Mat mat, Vec x, Vec y) 2564 { 2565 PetscFunctionBegin; 2566 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 2567 PetscValidType(mat, 1); 2568 PetscValidHeaderSpecific(x, VEC_CLASSID, 2); 2569 VecCheckAssembled(x); 2570 PetscValidHeaderSpecific(y, VEC_CLASSID, 3); 2571 PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 2572 PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 2573 PetscCheck(x != y, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "x and y must be different vectors"); 2574 PetscCheck(mat->cmap->N == x->map->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_SIZ, "Mat mat,Vec x: global dim %" PetscInt_FMT " %" PetscInt_FMT, mat->cmap->N, x->map->N); 2575 PetscCheck(mat->rmap->N == y->map->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_SIZ, "Mat mat,Vec y: global dim %" PetscInt_FMT " %" PetscInt_FMT, mat->rmap->N, y->map->N); 2576 PetscCheck(mat->cmap->n == x->map->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Mat mat,Vec x: local dim %" PetscInt_FMT " %" PetscInt_FMT, mat->cmap->n, x->map->n); 2577 PetscCheck(mat->rmap->n == y->map->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Mat mat,Vec y: local dim %" PetscInt_FMT " %" PetscInt_FMT, mat->rmap->n, y->map->n); 2578 PetscCall(VecSetErrorIfLocked(y, 3)); 2579 if (mat->erroriffailure) PetscCall(VecValidValues_Internal(x, 2, PETSC_TRUE)); 2580 MatCheckPreallocated(mat, 1); 2581 2582 PetscCall(VecLockReadPush(x)); 2583 PetscCall(PetscLogEventBegin(MAT_Mult, mat, x, y, 0)); 2584 PetscUseTypeMethod(mat, mult, x, y); 2585 PetscCall(PetscLogEventEnd(MAT_Mult, mat, x, y, 0)); 2586 if (mat->erroriffailure) PetscCall(VecValidValues_Internal(y, 3, PETSC_FALSE)); 2587 PetscCall(VecLockReadPop(x)); 2588 PetscFunctionReturn(PETSC_SUCCESS); 2589 } 2590 2591 /*@ 2592 MatMultTranspose - Computes matrix transpose times a vector y = A^T * x. 2593 2594 Neighbor-wise Collective 2595 2596 Input Parameters: 2597 + mat - the matrix 2598 - x - the vector to be multiplied 2599 2600 Output Parameters: 2601 . y - the result 2602 2603 Level: beginner 2604 2605 Notes: 2606 The vectors `x` and `y` cannot be the same. I.e., one cannot 2607 call `MatMultTranspose`(A,y,y). 2608 2609 For complex numbers this does NOT compute the Hermitian (complex conjugate) transpose multiple, 2610 use `MatMultHermitianTranspose()` 2611 2612 .seealso: `Mat`, `MatMult()`, `MatMultAdd()`, `MatMultTransposeAdd()`, `MatMultHermitianTranspose()`, `MatTranspose()` 2613 @*/ 2614 PetscErrorCode MatMultTranspose(Mat mat, Vec x, Vec y) 2615 { 2616 PetscErrorCode (*op)(Mat, Vec, Vec) = NULL; 2617 2618 PetscFunctionBegin; 2619 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 2620 PetscValidType(mat, 1); 2621 PetscValidHeaderSpecific(x, VEC_CLASSID, 2); 2622 VecCheckAssembled(x); 2623 PetscValidHeaderSpecific(y, VEC_CLASSID, 3); 2624 2625 PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 2626 PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 2627 PetscCheck(x != y, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "x and y must be different vectors"); 2628 PetscCheck(mat->cmap->N == y->map->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_SIZ, "Mat mat,Vec y: global dim %" PetscInt_FMT " %" PetscInt_FMT, mat->cmap->N, y->map->N); 2629 PetscCheck(mat->rmap->N == x->map->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_SIZ, "Mat mat,Vec x: global dim %" PetscInt_FMT " %" PetscInt_FMT, mat->rmap->N, x->map->N); 2630 PetscCheck(mat->cmap->n == y->map->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Mat mat,Vec y: local dim %" PetscInt_FMT " %" PetscInt_FMT, mat->cmap->n, y->map->n); 2631 PetscCheck(mat->rmap->n == x->map->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Mat mat,Vec x: local dim %" PetscInt_FMT " %" PetscInt_FMT, mat->rmap->n, x->map->n); 2632 if (mat->erroriffailure) PetscCall(VecValidValues_Internal(x, 2, PETSC_TRUE)); 2633 MatCheckPreallocated(mat, 1); 2634 2635 if (!mat->ops->multtranspose) { 2636 if (mat->symmetric == PETSC_BOOL3_TRUE && mat->ops->mult) op = mat->ops->mult; 2637 PetscCheck(op, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Matrix type %s does not have a multiply transpose defined or is symmetric and does not have a multiply defined", ((PetscObject)mat)->type_name); 2638 } else op = mat->ops->multtranspose; 2639 PetscCall(PetscLogEventBegin(MAT_MultTranspose, mat, x, y, 0)); 2640 PetscCall(VecLockReadPush(x)); 2641 PetscCall((*op)(mat, x, y)); 2642 PetscCall(VecLockReadPop(x)); 2643 PetscCall(PetscLogEventEnd(MAT_MultTranspose, mat, x, y, 0)); 2644 PetscCall(PetscObjectStateIncrease((PetscObject)y)); 2645 if (mat->erroriffailure) PetscCall(VecValidValues_Internal(y, 3, PETSC_FALSE)); 2646 PetscFunctionReturn(PETSC_SUCCESS); 2647 } 2648 2649 /*@ 2650 MatMultHermitianTranspose - Computes matrix Hermitian transpose times a vector. 2651 2652 Neighbor-wise Collective 2653 2654 Input Parameters: 2655 + mat - the matrix 2656 - x - the vector to be multilplied 2657 2658 Output Parameters: 2659 . y - the result 2660 2661 Level: beginner 2662 2663 Notes: 2664 The vectors `x` and `y` cannot be the same. I.e., one cannot 2665 call `MatMultHermitianTranspose`(A,y,y). 2666 2667 Also called the conjugate transpose, complex conjugate transpose, or adjoint. 2668 2669 For real numbers `MatMultTranspose()` and `MatMultHermitianTranspose()` are identical. 2670 2671 .seealso: `Mat`, `MatMult()`, `MatMultAdd()`, `MatMultHermitianTransposeAdd()`, `MatMultTranspose()` 2672 @*/ 2673 PetscErrorCode MatMultHermitianTranspose(Mat mat, Vec x, Vec y) 2674 { 2675 PetscFunctionBegin; 2676 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 2677 PetscValidType(mat, 1); 2678 PetscValidHeaderSpecific(x, VEC_CLASSID, 2); 2679 PetscValidHeaderSpecific(y, VEC_CLASSID, 3); 2680 2681 PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 2682 PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 2683 PetscCheck(x != y, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "x and y must be different vectors"); 2684 PetscCheck(mat->cmap->N == y->map->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_SIZ, "Mat mat,Vec y: global dim %" PetscInt_FMT " %" PetscInt_FMT, mat->cmap->N, y->map->N); 2685 PetscCheck(mat->rmap->N == x->map->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_SIZ, "Mat mat,Vec x: global dim %" PetscInt_FMT " %" PetscInt_FMT, mat->rmap->N, x->map->N); 2686 PetscCheck(mat->cmap->n == y->map->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Mat mat,Vec y: local dim %" PetscInt_FMT " %" PetscInt_FMT, mat->cmap->n, y->map->n); 2687 PetscCheck(mat->rmap->n == x->map->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Mat mat,Vec x: local dim %" PetscInt_FMT " %" PetscInt_FMT, mat->rmap->n, x->map->n); 2688 MatCheckPreallocated(mat, 1); 2689 2690 PetscCall(PetscLogEventBegin(MAT_MultHermitianTranspose, mat, x, y, 0)); 2691 #if defined(PETSC_USE_COMPLEX) 2692 if (mat->ops->multhermitiantranspose || (mat->hermitian == PETSC_BOOL3_TRUE && mat->ops->mult)) { 2693 PetscCall(VecLockReadPush(x)); 2694 if (mat->ops->multhermitiantranspose) PetscUseTypeMethod(mat, multhermitiantranspose, x, y); 2695 else PetscUseTypeMethod(mat, mult, x, y); 2696 PetscCall(VecLockReadPop(x)); 2697 } else { 2698 Vec w; 2699 PetscCall(VecDuplicate(x, &w)); 2700 PetscCall(VecCopy(x, w)); 2701 PetscCall(VecConjugate(w)); 2702 PetscCall(MatMultTranspose(mat, w, y)); 2703 PetscCall(VecDestroy(&w)); 2704 PetscCall(VecConjugate(y)); 2705 } 2706 PetscCall(PetscObjectStateIncrease((PetscObject)y)); 2707 #else 2708 PetscCall(MatMultTranspose(mat, x, y)); 2709 #endif 2710 PetscCall(PetscLogEventEnd(MAT_MultHermitianTranspose, mat, x, y, 0)); 2711 PetscFunctionReturn(PETSC_SUCCESS); 2712 } 2713 2714 /*@ 2715 MatMultAdd - Computes v3 = v2 + A * v1. 2716 2717 Neighbor-wise Collective 2718 2719 Input Parameters: 2720 + mat - the matrix 2721 - v1, v2 - the vectors 2722 2723 Output Parameters: 2724 . v3 - the result 2725 2726 Level: beginner 2727 2728 Note: 2729 The vectors `v1` and `v3` cannot be the same. I.e., one cannot 2730 call `MatMultAdd`(A,v1,v2,v1). 2731 2732 .seealso: `Mat`, `MatMultTranspose()`, `MatMult()`, `MatMultTransposeAdd()` 2733 @*/ 2734 PetscErrorCode MatMultAdd(Mat mat, Vec v1, Vec v2, Vec v3) 2735 { 2736 PetscFunctionBegin; 2737 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 2738 PetscValidType(mat, 1); 2739 PetscValidHeaderSpecific(v1, VEC_CLASSID, 2); 2740 PetscValidHeaderSpecific(v2, VEC_CLASSID, 3); 2741 PetscValidHeaderSpecific(v3, VEC_CLASSID, 4); 2742 2743 PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 2744 PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 2745 PetscCheck(mat->cmap->N == v1->map->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_SIZ, "Mat mat,Vec v1: global dim %" PetscInt_FMT " %" PetscInt_FMT, mat->cmap->N, v1->map->N); 2746 /* PetscCheck(mat->rmap->N == v2->map->N,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Mat mat,Vec v2: global dim %" PetscInt_FMT " %" PetscInt_FMT,mat->rmap->N,v2->map->N); 2747 PetscCheck(mat->rmap->N == v3->map->N,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Mat mat,Vec v3: global dim %" PetscInt_FMT " %" PetscInt_FMT,mat->rmap->N,v3->map->N); */ 2748 PetscCheck(mat->rmap->n == v3->map->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Mat mat,Vec v3: local dim %" PetscInt_FMT " %" PetscInt_FMT, mat->rmap->n, v3->map->n); 2749 PetscCheck(mat->rmap->n == v2->map->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Mat mat,Vec v2: local dim %" PetscInt_FMT " %" PetscInt_FMT, mat->rmap->n, v2->map->n); 2750 PetscCheck(v1 != v3, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_IDN, "v1 and v3 must be different vectors"); 2751 MatCheckPreallocated(mat, 1); 2752 2753 PetscCall(PetscLogEventBegin(MAT_MultAdd, mat, v1, v2, v3)); 2754 PetscCall(VecLockReadPush(v1)); 2755 PetscUseTypeMethod(mat, multadd, v1, v2, v3); 2756 PetscCall(VecLockReadPop(v1)); 2757 PetscCall(PetscLogEventEnd(MAT_MultAdd, mat, v1, v2, v3)); 2758 PetscCall(PetscObjectStateIncrease((PetscObject)v3)); 2759 PetscFunctionReturn(PETSC_SUCCESS); 2760 } 2761 2762 /*@ 2763 MatMultTransposeAdd - Computes v3 = v2 + A' * v1. 2764 2765 Neighbor-wise Collective 2766 2767 Input Parameters: 2768 + mat - the matrix 2769 - v1, v2 - the vectors 2770 2771 Output Parameters: 2772 . v3 - the result 2773 2774 Level: beginner 2775 2776 Note: 2777 The vectors `v1` and `v3` cannot be the same. I.e., one cannot 2778 call `MatMultTransposeAdd`(A,v1,v2,v1). 2779 2780 .seealso: `Mat`, `MatMultTranspose()`, `MatMultAdd()`, `MatMult()` 2781 @*/ 2782 PetscErrorCode MatMultTransposeAdd(Mat mat, Vec v1, Vec v2, Vec v3) 2783 { 2784 PetscErrorCode (*op)(Mat, Vec, Vec, Vec) = (!mat->ops->multtransposeadd && mat->symmetric) ? mat->ops->multadd : mat->ops->multtransposeadd; 2785 2786 PetscFunctionBegin; 2787 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 2788 PetscValidType(mat, 1); 2789 PetscValidHeaderSpecific(v1, VEC_CLASSID, 2); 2790 PetscValidHeaderSpecific(v2, VEC_CLASSID, 3); 2791 PetscValidHeaderSpecific(v3, VEC_CLASSID, 4); 2792 2793 PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 2794 PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 2795 PetscCheck(mat->rmap->N == v1->map->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_SIZ, "Mat mat,Vec v1: global dim %" PetscInt_FMT " %" PetscInt_FMT, mat->rmap->N, v1->map->N); 2796 PetscCheck(mat->cmap->N == v2->map->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_SIZ, "Mat mat,Vec v2: global dim %" PetscInt_FMT " %" PetscInt_FMT, mat->cmap->N, v2->map->N); 2797 PetscCheck(mat->cmap->N == v3->map->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_SIZ, "Mat mat,Vec v3: global dim %" PetscInt_FMT " %" PetscInt_FMT, mat->cmap->N, v3->map->N); 2798 PetscCheck(v1 != v3, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_IDN, "v1 and v3 must be different vectors"); 2799 PetscCheck(op, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Mat type %s", ((PetscObject)mat)->type_name); 2800 MatCheckPreallocated(mat, 1); 2801 2802 PetscCall(PetscLogEventBegin(MAT_MultTransposeAdd, mat, v1, v2, v3)); 2803 PetscCall(VecLockReadPush(v1)); 2804 PetscCall((*op)(mat, v1, v2, v3)); 2805 PetscCall(VecLockReadPop(v1)); 2806 PetscCall(PetscLogEventEnd(MAT_MultTransposeAdd, mat, v1, v2, v3)); 2807 PetscCall(PetscObjectStateIncrease((PetscObject)v3)); 2808 PetscFunctionReturn(PETSC_SUCCESS); 2809 } 2810 2811 /*@ 2812 MatMultHermitianTransposeAdd - Computes v3 = v2 + A^H * v1. 2813 2814 Neighbor-wise Collective 2815 2816 Input Parameters: 2817 + mat - the matrix 2818 - v1, v2 - the vectors 2819 2820 Output Parameters: 2821 . v3 - the result 2822 2823 Level: beginner 2824 2825 Note: 2826 The vectors `v1` and `v3` cannot be the same. I.e., one cannot 2827 call `MatMultHermitianTransposeAdd`(A,v1,v2,v1). 2828 2829 .seealso: `Mat`, `MatMultHermitianTranspose()`, `MatMultTranspose()`, `MatMultAdd()`, `MatMult()` 2830 @*/ 2831 PetscErrorCode MatMultHermitianTransposeAdd(Mat mat, Vec v1, Vec v2, Vec v3) 2832 { 2833 PetscFunctionBegin; 2834 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 2835 PetscValidType(mat, 1); 2836 PetscValidHeaderSpecific(v1, VEC_CLASSID, 2); 2837 PetscValidHeaderSpecific(v2, VEC_CLASSID, 3); 2838 PetscValidHeaderSpecific(v3, VEC_CLASSID, 4); 2839 2840 PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 2841 PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 2842 PetscCheck(v1 != v3, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_IDN, "v1 and v3 must be different vectors"); 2843 PetscCheck(mat->rmap->N == v1->map->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_SIZ, "Mat mat,Vec v1: global dim %" PetscInt_FMT " %" PetscInt_FMT, mat->rmap->N, v1->map->N); 2844 PetscCheck(mat->cmap->N == v2->map->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_SIZ, "Mat mat,Vec v2: global dim %" PetscInt_FMT " %" PetscInt_FMT, mat->cmap->N, v2->map->N); 2845 PetscCheck(mat->cmap->N == v3->map->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_SIZ, "Mat mat,Vec v3: global dim %" PetscInt_FMT " %" PetscInt_FMT, mat->cmap->N, v3->map->N); 2846 MatCheckPreallocated(mat, 1); 2847 2848 PetscCall(PetscLogEventBegin(MAT_MultHermitianTransposeAdd, mat, v1, v2, v3)); 2849 PetscCall(VecLockReadPush(v1)); 2850 if (mat->ops->multhermitiantransposeadd) PetscUseTypeMethod(mat, multhermitiantransposeadd, v1, v2, v3); 2851 else { 2852 Vec w, z; 2853 PetscCall(VecDuplicate(v1, &w)); 2854 PetscCall(VecCopy(v1, w)); 2855 PetscCall(VecConjugate(w)); 2856 PetscCall(VecDuplicate(v3, &z)); 2857 PetscCall(MatMultTranspose(mat, w, z)); 2858 PetscCall(VecDestroy(&w)); 2859 PetscCall(VecConjugate(z)); 2860 if (v2 != v3) { 2861 PetscCall(VecWAXPY(v3, 1.0, v2, z)); 2862 } else { 2863 PetscCall(VecAXPY(v3, 1.0, z)); 2864 } 2865 PetscCall(VecDestroy(&z)); 2866 } 2867 PetscCall(VecLockReadPop(v1)); 2868 PetscCall(PetscLogEventEnd(MAT_MultHermitianTransposeAdd, mat, v1, v2, v3)); 2869 PetscCall(PetscObjectStateIncrease((PetscObject)v3)); 2870 PetscFunctionReturn(PETSC_SUCCESS); 2871 } 2872 2873 /*@C 2874 MatGetFactorType - gets the type of factorization it is 2875 2876 Not Collective 2877 2878 Input Parameters: 2879 . mat - the matrix 2880 2881 Output Parameters: 2882 . t - the type, one of `MAT_FACTOR_NONE`, `MAT_FACTOR_LU`, `MAT_FACTOR_CHOLESKY`, `MAT_FACTOR_ILU`, `MAT_FACTOR_ICC,MAT_FACTOR_ILUDT`, `MAT_FACTOR_QR` 2883 2884 Level: intermediate 2885 2886 .seealso: `Mat`, [Matrix Factorization](sec_matfactor), `MatFactorType`, `MatGetFactor()`, `MatSetFactorType()`, `MAT_FACTOR_NONE`, `MAT_FACTOR_LU`, `MAT_FACTOR_CHOLESKY`, `MAT_FACTOR_ILU`, 2887 `MAT_FACTOR_ICC,MAT_FACTOR_ILUDT`, `MAT_FACTOR_QR` 2888 @*/ 2889 PetscErrorCode MatGetFactorType(Mat mat, MatFactorType *t) 2890 { 2891 PetscFunctionBegin; 2892 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 2893 PetscValidType(mat, 1); 2894 PetscValidPointer(t, 2); 2895 *t = mat->factortype; 2896 PetscFunctionReturn(PETSC_SUCCESS); 2897 } 2898 2899 /*@C 2900 MatSetFactorType - sets the type of factorization it is 2901 2902 Logically Collective 2903 2904 Input Parameters: 2905 + mat - the matrix 2906 - t - the type, one of `MAT_FACTOR_NONE`, `MAT_FACTOR_LU`, `MAT_FACTOR_CHOLESKY`, `MAT_FACTOR_ILU`, `MAT_FACTOR_ICC,MAT_FACTOR_ILUDT`, `MAT_FACTOR_QR` 2907 2908 Level: intermediate 2909 2910 .seealso: `Mat`, [Matrix Factorization](sec_matfactor), `MatFactorType`, `MatGetFactor()`, `MatGetFactorType()`, `MAT_FACTOR_NONE`, `MAT_FACTOR_LU`, `MAT_FACTOR_CHOLESKY`, `MAT_FACTOR_ILU`, 2911 `MAT_FACTOR_ICC,MAT_FACTOR_ILUDT`, `MAT_FACTOR_QR` 2912 @*/ 2913 PetscErrorCode MatSetFactorType(Mat mat, MatFactorType t) 2914 { 2915 PetscFunctionBegin; 2916 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 2917 PetscValidType(mat, 1); 2918 mat->factortype = t; 2919 PetscFunctionReturn(PETSC_SUCCESS); 2920 } 2921 2922 /*@C 2923 MatGetInfo - Returns information about matrix storage (number of 2924 nonzeros, memory, etc.). 2925 2926 Collective if `MAT_GLOBAL_MAX` or `MAT_GLOBAL_SUM` is used as the flag 2927 2928 Input Parameters: 2929 + mat - the matrix 2930 - flag - flag indicating the type of parameters to be returned (`MAT_LOCAL` - local matrix, `MAT_GLOBAL_MAX` - maximum over all processors, `MAT_GLOBAL_SUM` - sum over all processors) 2931 2932 Output Parameters: 2933 . info - matrix information context 2934 2935 Notes: 2936 The `MatInfo` context contains a variety of matrix data, including 2937 number of nonzeros allocated and used, number of mallocs during 2938 matrix assembly, etc. Additional information for factored matrices 2939 is provided (such as the fill ratio, number of mallocs during 2940 factorization, etc.). Much of this info is printed to `PETSC_STDOUT` 2941 when using the runtime options 2942 $ -info -mat_view ::ascii_info 2943 2944 Example for C/C++ Users: 2945 See the file ${PETSC_DIR}/include/petscmat.h for a complete list of 2946 data within the MatInfo context. For example, 2947 .vb 2948 MatInfo info; 2949 Mat A; 2950 double mal, nz_a, nz_u; 2951 2952 MatGetInfo(A,MAT_LOCAL,&info); 2953 mal = info.mallocs; 2954 nz_a = info.nz_allocated; 2955 .ve 2956 2957 Example for Fortran Users: 2958 Fortran users should declare info as a double precision 2959 array of dimension `MAT_INFO_SIZE`, and then extract the parameters 2960 of interest. See the file ${PETSC_DIR}/include/petsc/finclude/petscmat.h 2961 a complete list of parameter names. 2962 .vb 2963 double precision info(MAT_INFO_SIZE) 2964 double precision mal, nz_a 2965 Mat A 2966 integer ierr 2967 2968 call MatGetInfo(A,MAT_LOCAL,info,ierr) 2969 mal = info(MAT_INFO_MALLOCS) 2970 nz_a = info(MAT_INFO_NZ_ALLOCATED) 2971 .ve 2972 2973 Level: intermediate 2974 2975 Developer Note: 2976 Fortran interface is not autogenerated as the 2977 interface definition cannot be generated correctly [due to MatInfo] 2978 2979 .seealso: `Mat`, `MatInfo`, `MatStashGetInfo()` 2980 @*/ 2981 PetscErrorCode MatGetInfo(Mat mat, MatInfoType flag, MatInfo *info) 2982 { 2983 PetscFunctionBegin; 2984 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 2985 PetscValidType(mat, 1); 2986 PetscValidPointer(info, 3); 2987 MatCheckPreallocated(mat, 1); 2988 PetscUseTypeMethod(mat, getinfo, flag, info); 2989 PetscFunctionReturn(PETSC_SUCCESS); 2990 } 2991 2992 /* 2993 This is used by external packages where it is not easy to get the info from the actual 2994 matrix factorization. 2995 */ 2996 PetscErrorCode MatGetInfo_External(Mat A, MatInfoType flag, MatInfo *info) 2997 { 2998 PetscFunctionBegin; 2999 PetscCall(PetscMemzero(info, sizeof(MatInfo))); 3000 PetscFunctionReturn(PETSC_SUCCESS); 3001 } 3002 3003 /*@C 3004 MatLUFactor - Performs in-place LU factorization of matrix. 3005 3006 Collective 3007 3008 Input Parameters: 3009 + mat - the matrix 3010 . row - row permutation 3011 . col - column permutation 3012 - info - options for factorization, includes 3013 .vb 3014 fill - expected fill as ratio of original fill. 3015 dtcol - pivot tolerance (0 no pivot, 1 full column pivoting) 3016 Run with the option -info to determine an optimal value to use 3017 .ve 3018 Level: developer 3019 3020 Notes: 3021 Most users should employ the `KSP` interface for linear solvers 3022 instead of working directly with matrix algebra routines such as this. 3023 See, e.g., `KSPCreate()`. 3024 3025 This changes the state of the matrix to a factored matrix; it cannot be used 3026 for example with `MatSetValues()` unless one first calls `MatSetUnfactored()`. 3027 3028 This is really in-place only for dense matrices, the preferred approach is to use `MatGetFactor()`, `MatLUFactorSymbolic()`, and `MatLUFactorNumeric()` 3029 when not using `KSP`. 3030 3031 Developer Note: 3032 The Fortran interface is not autogenerated as the 3033 interface definition cannot be generated correctly [due to `MatFactorInfo`] 3034 3035 .seealso: [Matrix Factorization](sec_matfactor), `Mat`, `MatFactorType`, `MatLUFactorSymbolic()`, `MatLUFactorNumeric()`, `MatCholeskyFactor()`, 3036 `MatGetOrdering()`, `MatSetUnfactored()`, `MatFactorInfo`, `MatGetFactor()` 3037 @*/ 3038 PetscErrorCode MatLUFactor(Mat mat, IS row, IS col, const MatFactorInfo *info) 3039 { 3040 MatFactorInfo tinfo; 3041 3042 PetscFunctionBegin; 3043 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 3044 if (row) PetscValidHeaderSpecific(row, IS_CLASSID, 2); 3045 if (col) PetscValidHeaderSpecific(col, IS_CLASSID, 3); 3046 if (info) PetscValidPointer(info, 4); 3047 PetscValidType(mat, 1); 3048 PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 3049 PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 3050 MatCheckPreallocated(mat, 1); 3051 if (!info) { 3052 PetscCall(MatFactorInfoInitialize(&tinfo)); 3053 info = &tinfo; 3054 } 3055 3056 PetscCall(PetscLogEventBegin(MAT_LUFactor, mat, row, col, 0)); 3057 PetscUseTypeMethod(mat, lufactor, row, col, info); 3058 PetscCall(PetscLogEventEnd(MAT_LUFactor, mat, row, col, 0)); 3059 PetscCall(PetscObjectStateIncrease((PetscObject)mat)); 3060 PetscFunctionReturn(PETSC_SUCCESS); 3061 } 3062 3063 /*@C 3064 MatILUFactor - Performs in-place ILU factorization of matrix. 3065 3066 Collective 3067 3068 Input Parameters: 3069 + mat - the matrix 3070 . row - row permutation 3071 . col - column permutation 3072 - info - structure containing 3073 .vb 3074 levels - number of levels of fill. 3075 expected fill - as ratio of original fill. 3076 1 or 0 - indicating force fill on diagonal (improves robustness for matrices 3077 missing diagonal entries) 3078 .ve 3079 3080 Level: developer 3081 3082 Notes: 3083 Most users should employ the `KSP` interface for linear solvers 3084 instead of working directly with matrix algebra routines such as this. 3085 See, e.g., `KSPCreate()`. 3086 3087 Probably really in-place only when level of fill is zero, otherwise allocates 3088 new space to store factored matrix and deletes previous memory. The preferred approach is to use `MatGetFactor()`, `MatILUFactorSymbolic()`, and `MatILUFactorNumeric()` 3089 when not using `KSP`. 3090 3091 Developer Note: 3092 The Fortran interface is not autogenerated as the 3093 interface definition cannot be generated correctly [due to MatFactorInfo] 3094 3095 .seealso: `Mat`, [Matrix Factorization](sec_matfactor), `MatILUFactorSymbolic()`, `MatLUFactorNumeric()`, `MatCholeskyFactor()`, `MatFactorInfo` 3096 @*/ 3097 PetscErrorCode MatILUFactor(Mat mat, IS row, IS col, const MatFactorInfo *info) 3098 { 3099 PetscFunctionBegin; 3100 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 3101 if (row) PetscValidHeaderSpecific(row, IS_CLASSID, 2); 3102 if (col) PetscValidHeaderSpecific(col, IS_CLASSID, 3); 3103 PetscValidPointer(info, 4); 3104 PetscValidType(mat, 1); 3105 PetscCheck(mat->rmap->N == mat->cmap->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONG, "matrix must be square"); 3106 PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 3107 PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 3108 MatCheckPreallocated(mat, 1); 3109 3110 PetscCall(PetscLogEventBegin(MAT_ILUFactor, mat, row, col, 0)); 3111 PetscUseTypeMethod(mat, ilufactor, row, col, info); 3112 PetscCall(PetscLogEventEnd(MAT_ILUFactor, mat, row, col, 0)); 3113 PetscCall(PetscObjectStateIncrease((PetscObject)mat)); 3114 PetscFunctionReturn(PETSC_SUCCESS); 3115 } 3116 3117 /*@C 3118 MatLUFactorSymbolic - Performs symbolic LU factorization of matrix. 3119 Call this routine before calling `MatLUFactorNumeric()` and after `MatGetFactor()`. 3120 3121 Collective on fact 3122 3123 Input Parameters: 3124 + fact - the factor matrix obtained with `MatGetFactor()` 3125 . mat - the matrix 3126 . row, col - row and column permutations 3127 - info - options for factorization, includes 3128 .vb 3129 fill - expected fill as ratio of original fill. Run with the option -info to determine an optimal value to use 3130 dtcol - pivot tolerance (0 no pivot, 1 full column pivoting) 3131 .ve 3132 3133 Level: developer 3134 3135 Notes: 3136 See [Matrix Factorization](sec_matfactor) for additional information about factorizations 3137 3138 Most users should employ the simplified `KSP` interface for linear solvers 3139 instead of working directly with matrix algebra routines such as this. 3140 See, e.g., `KSPCreate()`. 3141 3142 Developer Note: 3143 The Fortran interface is not autogenerated as the 3144 interface definition cannot be generated correctly [due to `MatFactorInfo`] 3145 3146 .seealso: `Mat`, [Matrix Factorization](sec_matfactor), `MatGetFactor()`, `MatLUFactor()`, `MatLUFactorNumeric()`, `MatCholeskyFactor()`, `MatFactorInfo`, `MatFactorInfoInitialize()` 3147 @*/ 3148 PetscErrorCode MatLUFactorSymbolic(Mat fact, Mat mat, IS row, IS col, const MatFactorInfo *info) 3149 { 3150 MatFactorInfo tinfo; 3151 3152 PetscFunctionBegin; 3153 PetscValidHeaderSpecific(mat, MAT_CLASSID, 2); 3154 if (row) PetscValidHeaderSpecific(row, IS_CLASSID, 3); 3155 if (col) PetscValidHeaderSpecific(col, IS_CLASSID, 4); 3156 if (info) PetscValidPointer(info, 5); 3157 PetscValidType(mat, 2); 3158 PetscValidPointer(fact, 1); 3159 PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 3160 PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 3161 if (!(fact)->ops->lufactorsymbolic) { 3162 MatSolverType stype; 3163 PetscCall(MatFactorGetSolverType(fact, &stype)); 3164 SETERRQ(PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Matrix type %s symbolic LU using solver package %s", ((PetscObject)mat)->type_name, stype); 3165 } 3166 MatCheckPreallocated(mat, 2); 3167 if (!info) { 3168 PetscCall(MatFactorInfoInitialize(&tinfo)); 3169 info = &tinfo; 3170 } 3171 3172 if (!fact->trivialsymbolic) PetscCall(PetscLogEventBegin(MAT_LUFactorSymbolic, mat, row, col, 0)); 3173 PetscCall((fact->ops->lufactorsymbolic)(fact, mat, row, col, info)); 3174 if (!fact->trivialsymbolic) PetscCall(PetscLogEventEnd(MAT_LUFactorSymbolic, mat, row, col, 0)); 3175 PetscCall(PetscObjectStateIncrease((PetscObject)fact)); 3176 PetscFunctionReturn(PETSC_SUCCESS); 3177 } 3178 3179 /*@C 3180 MatLUFactorNumeric - Performs numeric LU factorization of a matrix. 3181 Call this routine after first calling `MatLUFactorSymbolic()` and `MatGetFactor()`. 3182 3183 Collective on fact 3184 3185 Input Parameters: 3186 + fact - the factor matrix obtained with `MatGetFactor()` 3187 . mat - the matrix 3188 - info - options for factorization 3189 3190 Level: developer 3191 3192 Notes: 3193 See `MatLUFactor()` for in-place factorization. See 3194 `MatCholeskyFactorNumeric()` for the symmetric, positive definite case. 3195 3196 Most users should employ the `KSP` interface for linear solvers 3197 instead of working directly with matrix algebra routines such as this. 3198 See, e.g., `KSPCreate()`. 3199 3200 Developer Note: 3201 The Fortran interface is not autogenerated as the 3202 interface definition cannot be generated correctly [due to `MatFactorInfo`] 3203 3204 .seealso: `Mat`, [Matrix Factorization](sec_matfactor), `MatGetFactor()`, `MatFactorInfo`, `MatLUFactorSymbolic()`, `MatLUFactor()`, `MatCholeskyFactor()` 3205 @*/ 3206 PetscErrorCode MatLUFactorNumeric(Mat fact, Mat mat, const MatFactorInfo *info) 3207 { 3208 MatFactorInfo tinfo; 3209 3210 PetscFunctionBegin; 3211 PetscValidHeaderSpecific(mat, MAT_CLASSID, 2); 3212 PetscValidType(mat, 2); 3213 PetscValidPointer(fact, 1); 3214 PetscValidHeaderSpecific(fact, MAT_CLASSID, 1); 3215 PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 3216 PetscCheck(mat->rmap->N == (fact)->rmap->N && mat->cmap->N == (fact)->cmap->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_SIZ, "Mat mat,Mat fact: global dimensions are different %" PetscInt_FMT " should = %" PetscInt_FMT " %" PetscInt_FMT " should = %" PetscInt_FMT, 3217 mat->rmap->N, (fact)->rmap->N, mat->cmap->N, (fact)->cmap->N); 3218 3219 PetscCheck((fact)->ops->lufactornumeric, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Mat type %s numeric LU", ((PetscObject)mat)->type_name); 3220 MatCheckPreallocated(mat, 2); 3221 if (!info) { 3222 PetscCall(MatFactorInfoInitialize(&tinfo)); 3223 info = &tinfo; 3224 } 3225 3226 if (!fact->trivialsymbolic) PetscCall(PetscLogEventBegin(MAT_LUFactorNumeric, mat, fact, 0, 0)); 3227 else PetscCall(PetscLogEventBegin(MAT_LUFactor, mat, fact, 0, 0)); 3228 PetscCall((fact->ops->lufactornumeric)(fact, mat, info)); 3229 if (!fact->trivialsymbolic) PetscCall(PetscLogEventEnd(MAT_LUFactorNumeric, mat, fact, 0, 0)); 3230 else PetscCall(PetscLogEventEnd(MAT_LUFactor, mat, fact, 0, 0)); 3231 PetscCall(MatViewFromOptions(fact, NULL, "-mat_factor_view")); 3232 PetscCall(PetscObjectStateIncrease((PetscObject)fact)); 3233 PetscFunctionReturn(PETSC_SUCCESS); 3234 } 3235 3236 /*@C 3237 MatCholeskyFactor - Performs in-place Cholesky factorization of a 3238 symmetric matrix. 3239 3240 Collective 3241 3242 Input Parameters: 3243 + mat - the matrix 3244 . perm - row and column permutations 3245 - f - expected fill as ratio of original fill 3246 3247 Level: developer 3248 3249 Notes: 3250 See `MatLUFactor()` for the nonsymmetric case. See also `MatGetFactor()`, 3251 `MatCholeskyFactorSymbolic()`, and `MatCholeskyFactorNumeric()`. 3252 3253 Most users should employ the `KSP` interface for linear solvers 3254 instead of working directly with matrix algebra routines such as this. 3255 See, e.g., `KSPCreate()`. 3256 3257 Developer Note: 3258 The Fortran interface is not autogenerated as the 3259 interface definition cannot be generated correctly [due to `MatFactorInfo`] 3260 3261 .seealso: `Mat`, [Matrix Factorization](sec_matfactor), `MatGetFactor()`, `MatFactorInfo`, `MatLUFactor()`, `MatCholeskyFactorSymbolic()`, `MatCholeskyFactorNumeric()` 3262 `MatGetOrdering()` 3263 @*/ 3264 PetscErrorCode MatCholeskyFactor(Mat mat, IS perm, const MatFactorInfo *info) 3265 { 3266 MatFactorInfo tinfo; 3267 3268 PetscFunctionBegin; 3269 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 3270 PetscValidType(mat, 1); 3271 if (perm) PetscValidHeaderSpecific(perm, IS_CLASSID, 2); 3272 if (info) PetscValidPointer(info, 3); 3273 PetscCheck(mat->rmap->N == mat->cmap->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONG, "Matrix must be square"); 3274 PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 3275 PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 3276 MatCheckPreallocated(mat, 1); 3277 if (!info) { 3278 PetscCall(MatFactorInfoInitialize(&tinfo)); 3279 info = &tinfo; 3280 } 3281 3282 PetscCall(PetscLogEventBegin(MAT_CholeskyFactor, mat, perm, 0, 0)); 3283 PetscUseTypeMethod(mat, choleskyfactor, perm, info); 3284 PetscCall(PetscLogEventEnd(MAT_CholeskyFactor, mat, perm, 0, 0)); 3285 PetscCall(PetscObjectStateIncrease((PetscObject)mat)); 3286 PetscFunctionReturn(PETSC_SUCCESS); 3287 } 3288 3289 /*@C 3290 MatCholeskyFactorSymbolic - Performs symbolic Cholesky factorization 3291 of a symmetric matrix. 3292 3293 Collective on fact 3294 3295 Input Parameters: 3296 + fact - the factor matrix obtained with `MatGetFactor()` 3297 . mat - the matrix 3298 . perm - row and column permutations 3299 - info - options for factorization, includes 3300 .vb 3301 fill - expected fill as ratio of original fill. 3302 dtcol - pivot tolerance (0 no pivot, 1 full column pivoting) 3303 Run with the option -info to determine an optimal value to use 3304 .ve 3305 3306 Level: developer 3307 3308 Notes: 3309 See `MatLUFactorSymbolic()` for the nonsymmetric case. See also 3310 `MatCholeskyFactor()` and `MatCholeskyFactorNumeric()`. 3311 3312 Most users should employ the `KSP` interface for linear solvers 3313 instead of working directly with matrix algebra routines such as this. 3314 See, e.g., `KSPCreate()`. 3315 3316 Developer Note: 3317 The Fortran interface is not autogenerated as the 3318 interface definition cannot be generated correctly [due to `MatFactorInfo`] 3319 3320 .seealso: `Mat`, [Matrix Factorization](sec_matfactor), `MatFactorInfo`, `MatGetFactor()`, `MatLUFactorSymbolic()`, `MatCholeskyFactor()`, `MatCholeskyFactorNumeric()` 3321 `MatGetOrdering()` 3322 @*/ 3323 PetscErrorCode MatCholeskyFactorSymbolic(Mat fact, Mat mat, IS perm, const MatFactorInfo *info) 3324 { 3325 MatFactorInfo tinfo; 3326 3327 PetscFunctionBegin; 3328 PetscValidHeaderSpecific(mat, MAT_CLASSID, 2); 3329 PetscValidType(mat, 2); 3330 if (perm) PetscValidHeaderSpecific(perm, IS_CLASSID, 3); 3331 if (info) PetscValidPointer(info, 4); 3332 PetscValidPointer(fact, 1); 3333 PetscCheck(mat->rmap->N == mat->cmap->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONG, "Matrix must be square"); 3334 PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 3335 PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 3336 if (!(fact)->ops->choleskyfactorsymbolic) { 3337 MatSolverType stype; 3338 PetscCall(MatFactorGetSolverType(fact, &stype)); 3339 SETERRQ(PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Mat type %s symbolic factor Cholesky using solver package %s", ((PetscObject)mat)->type_name, stype); 3340 } 3341 MatCheckPreallocated(mat, 2); 3342 if (!info) { 3343 PetscCall(MatFactorInfoInitialize(&tinfo)); 3344 info = &tinfo; 3345 } 3346 3347 if (!fact->trivialsymbolic) PetscCall(PetscLogEventBegin(MAT_CholeskyFactorSymbolic, mat, perm, 0, 0)); 3348 PetscCall((fact->ops->choleskyfactorsymbolic)(fact, mat, perm, info)); 3349 if (!fact->trivialsymbolic) PetscCall(PetscLogEventEnd(MAT_CholeskyFactorSymbolic, mat, perm, 0, 0)); 3350 PetscCall(PetscObjectStateIncrease((PetscObject)fact)); 3351 PetscFunctionReturn(PETSC_SUCCESS); 3352 } 3353 3354 /*@C 3355 MatCholeskyFactorNumeric - Performs numeric Cholesky factorization 3356 of a symmetric matrix. Call this routine after first calling `MatGetFactor()` and 3357 `MatCholeskyFactorSymbolic()`. 3358 3359 Collective on fact 3360 3361 Input Parameters: 3362 + fact - the factor matrix obtained with `MatGetFactor()` 3363 . mat - the initial matrix 3364 . info - options for factorization 3365 - fact - the symbolic factor of mat 3366 3367 Level: developer 3368 3369 Note: 3370 Most users should employ the `KSP` interface for linear solvers 3371 instead of working directly with matrix algebra routines such as this. 3372 See, e.g., `KSPCreate()`. 3373 3374 Developer Note: 3375 The Fortran interface is not autogenerated as the 3376 interface definition cannot be generated correctly [due to `MatFactorInfo`] 3377 3378 .seealso: `Mat`, [Matrix Factorization](sec_matfactor), `MatFactorInfo`, `MatGetFactor()`, `MatCholeskyFactorSymbolic()`, `MatCholeskyFactor()`, `MatLUFactorNumeric()` 3379 @*/ 3380 PetscErrorCode MatCholeskyFactorNumeric(Mat fact, Mat mat, const MatFactorInfo *info) 3381 { 3382 MatFactorInfo tinfo; 3383 3384 PetscFunctionBegin; 3385 PetscValidHeaderSpecific(mat, MAT_CLASSID, 2); 3386 PetscValidType(mat, 2); 3387 PetscValidPointer(fact, 1); 3388 PetscValidHeaderSpecific(fact, MAT_CLASSID, 1); 3389 PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 3390 PetscCheck((fact)->ops->choleskyfactornumeric, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Mat type %s numeric factor Cholesky", ((PetscObject)mat)->type_name); 3391 PetscCheck(mat->rmap->N == (fact)->rmap->N && mat->cmap->N == (fact)->cmap->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_SIZ, "Mat mat,Mat fact: global dim %" PetscInt_FMT " should = %" PetscInt_FMT " %" PetscInt_FMT " should = %" PetscInt_FMT, 3392 mat->rmap->N, (fact)->rmap->N, mat->cmap->N, (fact)->cmap->N); 3393 MatCheckPreallocated(mat, 2); 3394 if (!info) { 3395 PetscCall(MatFactorInfoInitialize(&tinfo)); 3396 info = &tinfo; 3397 } 3398 3399 if (!fact->trivialsymbolic) PetscCall(PetscLogEventBegin(MAT_CholeskyFactorNumeric, mat, fact, 0, 0)); 3400 else PetscCall(PetscLogEventBegin(MAT_CholeskyFactor, mat, fact, 0, 0)); 3401 PetscCall((fact->ops->choleskyfactornumeric)(fact, mat, info)); 3402 if (!fact->trivialsymbolic) PetscCall(PetscLogEventEnd(MAT_CholeskyFactorNumeric, mat, fact, 0, 0)); 3403 else PetscCall(PetscLogEventEnd(MAT_CholeskyFactor, mat, fact, 0, 0)); 3404 PetscCall(MatViewFromOptions(fact, NULL, "-mat_factor_view")); 3405 PetscCall(PetscObjectStateIncrease((PetscObject)fact)); 3406 PetscFunctionReturn(PETSC_SUCCESS); 3407 } 3408 3409 /*@ 3410 MatQRFactor - Performs in-place QR factorization of matrix. 3411 3412 Collective 3413 3414 Input Parameters: 3415 + mat - the matrix 3416 . col - column permutation 3417 - info - options for factorization, includes 3418 .vb 3419 fill - expected fill as ratio of original fill. 3420 dtcol - pivot tolerance (0 no pivot, 1 full column pivoting) 3421 Run with the option -info to determine an optimal value to use 3422 .ve 3423 3424 Level: developer 3425 3426 Notes: 3427 Most users should employ the `KSP` interface for linear solvers 3428 instead of working directly with matrix algebra routines such as this. 3429 See, e.g., `KSPCreate()`. 3430 3431 This changes the state of the matrix to a factored matrix; it cannot be used 3432 for example with `MatSetValues()` unless one first calls `MatSetUnfactored()`. 3433 3434 Developer Note: 3435 The Fortran interface is not autogenerated as the 3436 interface definition cannot be generated correctly [due to MatFactorInfo] 3437 3438 .seealso: `Mat`, [Matrix Factorization](sec_matfactor), `MatFactorInfo`, `MatGetFactor()`, `MatQRFactorSymbolic()`, `MatQRFactorNumeric()`, `MatLUFactor()`, 3439 `MatSetUnfactored()`, `MatFactorInfo`, `MatGetFactor()` 3440 @*/ 3441 PetscErrorCode MatQRFactor(Mat mat, IS col, const MatFactorInfo *info) 3442 { 3443 PetscFunctionBegin; 3444 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 3445 if (col) PetscValidHeaderSpecific(col, IS_CLASSID, 2); 3446 if (info) PetscValidPointer(info, 3); 3447 PetscValidType(mat, 1); 3448 PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 3449 PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 3450 MatCheckPreallocated(mat, 1); 3451 PetscCall(PetscLogEventBegin(MAT_QRFactor, mat, col, 0, 0)); 3452 PetscUseMethod(mat, "MatQRFactor_C", (Mat, IS, const MatFactorInfo *), (mat, col, info)); 3453 PetscCall(PetscLogEventEnd(MAT_QRFactor, mat, col, 0, 0)); 3454 PetscCall(PetscObjectStateIncrease((PetscObject)mat)); 3455 PetscFunctionReturn(PETSC_SUCCESS); 3456 } 3457 3458 /*@ 3459 MatQRFactorSymbolic - Performs symbolic QR factorization of matrix. 3460 Call this routine after `MatGetFactor()` but before calling `MatQRFactorNumeric()`. 3461 3462 Collective on fact 3463 3464 Input Parameters: 3465 + fact - the factor matrix obtained with `MatGetFactor()` 3466 . mat - the matrix 3467 . col - column permutation 3468 - info - options for factorization, includes 3469 .vb 3470 fill - expected fill as ratio of original fill. 3471 dtcol - pivot tolerance (0 no pivot, 1 full column pivoting) 3472 Run with the option -info to determine an optimal value to use 3473 .ve 3474 3475 Level: developer 3476 3477 Note: 3478 Most users should employ the `KSP` interface for linear solvers 3479 instead of working directly with matrix algebra routines such as this. 3480 See, e.g., `KSPCreate()`. 3481 3482 Developer Note: 3483 The Fortran interface is not autogenerated as the 3484 interface definition cannot be generated correctly [due to `MatFactorInfo`] 3485 3486 .seealso: `Mat`, [Matrix Factorization](sec_matfactor), `MatGetFactor()`, `MatFactorInfo`, `MatQRFactor()`, `MatQRFactorNumeric()`, `MatLUFactor()`, `MatFactorInfo`, `MatFactorInfoInitialize()` 3487 @*/ 3488 PetscErrorCode MatQRFactorSymbolic(Mat fact, Mat mat, IS col, const MatFactorInfo *info) 3489 { 3490 MatFactorInfo tinfo; 3491 3492 PetscFunctionBegin; 3493 PetscValidHeaderSpecific(mat, MAT_CLASSID, 2); 3494 if (col) PetscValidHeaderSpecific(col, IS_CLASSID, 3); 3495 if (info) PetscValidPointer(info, 4); 3496 PetscValidType(mat, 2); 3497 PetscValidPointer(fact, 1); 3498 PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 3499 PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 3500 MatCheckPreallocated(mat, 2); 3501 if (!info) { 3502 PetscCall(MatFactorInfoInitialize(&tinfo)); 3503 info = &tinfo; 3504 } 3505 3506 if (!fact->trivialsymbolic) PetscCall(PetscLogEventBegin(MAT_QRFactorSymbolic, fact, mat, col, 0)); 3507 PetscUseMethod(fact, "MatQRFactorSymbolic_C", (Mat, Mat, IS, const MatFactorInfo *), (fact, mat, col, info)); 3508 if (!fact->trivialsymbolic) PetscCall(PetscLogEventEnd(MAT_QRFactorSymbolic, fact, mat, col, 0)); 3509 PetscCall(PetscObjectStateIncrease((PetscObject)fact)); 3510 PetscFunctionReturn(PETSC_SUCCESS); 3511 } 3512 3513 /*@ 3514 MatQRFactorNumeric - Performs numeric QR factorization of a matrix. 3515 Call this routine after first calling `MatGetFactor()`, and `MatQRFactorSymbolic()`. 3516 3517 Collective on fact 3518 3519 Input Parameters: 3520 + fact - the factor matrix obtained with `MatGetFactor()` 3521 . mat - the matrix 3522 - info - options for factorization 3523 3524 Level: developer 3525 3526 Notes: 3527 See `MatQRFactor()` for in-place factorization. 3528 3529 Most users should employ the `KSP` interface for linear solvers 3530 instead of working directly with matrix algebra routines such as this. 3531 See, e.g., `KSPCreate()`. 3532 3533 Developer Note: 3534 The Fortran interface is not autogenerated as the 3535 interface definition cannot be generated correctly [due to `MatFactorInfo`] 3536 3537 .seealso: `Mat`, [Matrix Factorization](sec_matfactor), `MatFactorInfo`, `MatGetFactor()`, `MatQRFactor()`, `MatQRFactorSymbolic()`, `MatLUFactor()` 3538 @*/ 3539 PetscErrorCode MatQRFactorNumeric(Mat fact, Mat mat, const MatFactorInfo *info) 3540 { 3541 MatFactorInfo tinfo; 3542 3543 PetscFunctionBegin; 3544 PetscValidHeaderSpecific(mat, MAT_CLASSID, 2); 3545 PetscValidType(mat, 2); 3546 PetscValidPointer(fact, 1); 3547 PetscValidHeaderSpecific(fact, MAT_CLASSID, 1); 3548 PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 3549 PetscCheck(mat->rmap->N == fact->rmap->N && mat->cmap->N == fact->cmap->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_SIZ, "Mat mat,Mat fact: global dimensions are different %" PetscInt_FMT " should = %" PetscInt_FMT " %" PetscInt_FMT " should = %" PetscInt_FMT, 3550 mat->rmap->N, (fact)->rmap->N, mat->cmap->N, (fact)->cmap->N); 3551 3552 MatCheckPreallocated(mat, 2); 3553 if (!info) { 3554 PetscCall(MatFactorInfoInitialize(&tinfo)); 3555 info = &tinfo; 3556 } 3557 3558 if (!fact->trivialsymbolic) PetscCall(PetscLogEventBegin(MAT_QRFactorNumeric, mat, fact, 0, 0)); 3559 else PetscCall(PetscLogEventBegin(MAT_QRFactor, mat, fact, 0, 0)); 3560 PetscUseMethod(fact, "MatQRFactorNumeric_C", (Mat, Mat, const MatFactorInfo *), (fact, mat, info)); 3561 if (!fact->trivialsymbolic) PetscCall(PetscLogEventEnd(MAT_QRFactorNumeric, mat, fact, 0, 0)); 3562 else PetscCall(PetscLogEventEnd(MAT_QRFactor, mat, fact, 0, 0)); 3563 PetscCall(MatViewFromOptions(fact, NULL, "-mat_factor_view")); 3564 PetscCall(PetscObjectStateIncrease((PetscObject)fact)); 3565 PetscFunctionReturn(PETSC_SUCCESS); 3566 } 3567 3568 /*@ 3569 MatSolve - Solves A x = b, given a factored matrix. 3570 3571 Neighbor-wise Collective 3572 3573 Input Parameters: 3574 + mat - the factored matrix 3575 - b - the right-hand-side vector 3576 3577 Output Parameter: 3578 . x - the result vector 3579 3580 Level: developer 3581 3582 Notes: 3583 The vectors `b` and `x` cannot be the same. I.e., one cannot 3584 call `MatSolve`(A,x,x). 3585 3586 Most users should employ the `KSP` interface for linear solvers 3587 instead of working directly with matrix algebra routines such as this. 3588 See, e.g., `KSPCreate()`. 3589 3590 .seealso: `Mat`, [Matrix Factorization](sec_matfactor), `MatGetFactor()`, `MatLUFactor()`, `MatSolveAdd()`, `MatSolveTranspose()`, `MatSolveTransposeAdd()` 3591 @*/ 3592 PetscErrorCode MatSolve(Mat mat, Vec b, Vec x) 3593 { 3594 PetscFunctionBegin; 3595 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 3596 PetscValidType(mat, 1); 3597 PetscValidHeaderSpecific(b, VEC_CLASSID, 2); 3598 PetscValidHeaderSpecific(x, VEC_CLASSID, 3); 3599 PetscCheckSameComm(mat, 1, b, 2); 3600 PetscCheckSameComm(mat, 1, x, 3); 3601 PetscCheck(x != b, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_IDN, "x and b must be different vectors"); 3602 PetscCheck(mat->cmap->N == x->map->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_SIZ, "Mat mat,Vec x: global dim %" PetscInt_FMT " %" PetscInt_FMT, mat->cmap->N, x->map->N); 3603 PetscCheck(mat->rmap->N == b->map->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_SIZ, "Mat mat,Vec b: global dim %" PetscInt_FMT " %" PetscInt_FMT, mat->rmap->N, b->map->N); 3604 PetscCheck(mat->rmap->n == b->map->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Mat mat,Vec b: local dim %" PetscInt_FMT " %" PetscInt_FMT, mat->rmap->n, b->map->n); 3605 if (!mat->rmap->N && !mat->cmap->N) PetscFunctionReturn(PETSC_SUCCESS); 3606 MatCheckPreallocated(mat, 1); 3607 3608 PetscCall(PetscLogEventBegin(MAT_Solve, mat, b, x, 0)); 3609 if (mat->factorerrortype) { 3610 PetscCall(PetscInfo(mat, "MatFactorError %d\n", mat->factorerrortype)); 3611 PetscCall(VecSetInf(x)); 3612 } else PetscUseTypeMethod(mat, solve, b, x); 3613 PetscCall(PetscLogEventEnd(MAT_Solve, mat, b, x, 0)); 3614 PetscCall(PetscObjectStateIncrease((PetscObject)x)); 3615 PetscFunctionReturn(PETSC_SUCCESS); 3616 } 3617 3618 static PetscErrorCode MatMatSolve_Basic(Mat A, Mat B, Mat X, PetscBool trans) 3619 { 3620 Vec b, x; 3621 PetscInt N, i; 3622 PetscErrorCode (*f)(Mat, Vec, Vec); 3623 PetscBool Abound, Bneedconv = PETSC_FALSE, Xneedconv = PETSC_FALSE; 3624 3625 PetscFunctionBegin; 3626 if (A->factorerrortype) { 3627 PetscCall(PetscInfo(A, "MatFactorError %d\n", A->factorerrortype)); 3628 PetscCall(MatSetInf(X)); 3629 PetscFunctionReturn(PETSC_SUCCESS); 3630 } 3631 f = (!trans || (!A->ops->solvetranspose && A->symmetric)) ? A->ops->solve : A->ops->solvetranspose; 3632 PetscCheck(f, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Mat type %s", ((PetscObject)A)->type_name); 3633 PetscCall(MatBoundToCPU(A, &Abound)); 3634 if (!Abound) { 3635 PetscCall(PetscObjectTypeCompareAny((PetscObject)B, &Bneedconv, MATSEQDENSE, MATMPIDENSE, "")); 3636 PetscCall(PetscObjectTypeCompareAny((PetscObject)X, &Xneedconv, MATSEQDENSE, MATMPIDENSE, "")); 3637 } 3638 #if defined(PETSC_HAVE_CUDA) 3639 if (Bneedconv) PetscCall(MatConvert(B, MATDENSECUDA, MAT_INPLACE_MATRIX, &B)); 3640 if (Xneedconv) PetscCall(MatConvert(X, MATDENSECUDA, MAT_INPLACE_MATRIX, &X)); 3641 #elif (PETSC_HAVE_HIP) 3642 if (Bneedconv) PetscCall(MatConvert(B, MATDENSEHIP, MAT_INPLACE_MATRIX, &B)); 3643 if (Xneedconv) PetscCall(MatConvert(X, MATDENSEHIP, MAT_INPLACE_MATRIX, &X)); 3644 #endif 3645 PetscCall(MatGetSize(B, NULL, &N)); 3646 for (i = 0; i < N; i++) { 3647 PetscCall(MatDenseGetColumnVecRead(B, i, &b)); 3648 PetscCall(MatDenseGetColumnVecWrite(X, i, &x)); 3649 PetscCall((*f)(A, b, x)); 3650 PetscCall(MatDenseRestoreColumnVecWrite(X, i, &x)); 3651 PetscCall(MatDenseRestoreColumnVecRead(B, i, &b)); 3652 } 3653 if (Bneedconv) PetscCall(MatConvert(B, MATDENSE, MAT_INPLACE_MATRIX, &B)); 3654 if (Xneedconv) PetscCall(MatConvert(X, MATDENSE, MAT_INPLACE_MATRIX, &X)); 3655 PetscFunctionReturn(PETSC_SUCCESS); 3656 } 3657 3658 /*@ 3659 MatMatSolve - Solves A X = B, given a factored matrix. 3660 3661 Neighbor-wise Collective 3662 3663 Input Parameters: 3664 + A - the factored matrix 3665 - B - the right-hand-side matrix `MATDENSE` (or sparse `MATAIJ`-- when using MUMPS) 3666 3667 Output Parameter: 3668 . X - the result matrix (dense matrix) 3669 3670 Level: developer 3671 3672 Note: 3673 If `B` is a `MATDENSE` matrix then one can call `MatMatSolve`(A,B,B) except with `MATSOLVERMKL_CPARDISO`; 3674 otherwise, `B` and `X` cannot be the same. 3675 3676 .seealso: `Mat`, [Matrix Factorization](sec_matfactor), `MatGetFactor()`, `MatSolve()`, `MatMatSolveTranspose()`, `MatLUFactor()`, `MatCholeskyFactor()` 3677 @*/ 3678 PetscErrorCode MatMatSolve(Mat A, Mat B, Mat X) 3679 { 3680 PetscFunctionBegin; 3681 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 3682 PetscValidType(A, 1); 3683 PetscValidHeaderSpecific(B, MAT_CLASSID, 2); 3684 PetscValidHeaderSpecific(X, MAT_CLASSID, 3); 3685 PetscCheckSameComm(A, 1, B, 2); 3686 PetscCheckSameComm(A, 1, X, 3); 3687 PetscCheck(A->cmap->N == X->rmap->N, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_SIZ, "Mat A,Mat X: global dim %" PetscInt_FMT " %" PetscInt_FMT, A->cmap->N, X->rmap->N); 3688 PetscCheck(A->rmap->N == B->rmap->N, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_SIZ, "Mat A,Mat B: global dim %" PetscInt_FMT " %" PetscInt_FMT, A->rmap->N, B->rmap->N); 3689 PetscCheck(X->cmap->N == B->cmap->N, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_SIZ, "Solution matrix must have same number of columns as rhs matrix"); 3690 if (!A->rmap->N && !A->cmap->N) PetscFunctionReturn(PETSC_SUCCESS); 3691 PetscCheck(A->factortype, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "Unfactored matrix"); 3692 MatCheckPreallocated(A, 1); 3693 3694 PetscCall(PetscLogEventBegin(MAT_MatSolve, A, B, X, 0)); 3695 if (!A->ops->matsolve) { 3696 PetscCall(PetscInfo(A, "Mat type %s using basic MatMatSolve\n", ((PetscObject)A)->type_name)); 3697 PetscCall(MatMatSolve_Basic(A, B, X, PETSC_FALSE)); 3698 } else PetscUseTypeMethod(A, matsolve, B, X); 3699 PetscCall(PetscLogEventEnd(MAT_MatSolve, A, B, X, 0)); 3700 PetscCall(PetscObjectStateIncrease((PetscObject)X)); 3701 PetscFunctionReturn(PETSC_SUCCESS); 3702 } 3703 3704 /*@ 3705 MatMatSolveTranspose - Solves A^T X = B, given a factored matrix. 3706 3707 Neighbor-wise Collective 3708 3709 Input Parameters: 3710 + A - the factored matrix 3711 - B - the right-hand-side matrix (`MATDENSE` matrix) 3712 3713 Output Parameter: 3714 . X - the result matrix (dense matrix) 3715 3716 Level: developer 3717 3718 Note: 3719 The matrices `B` and `X` cannot be the same. I.e., one cannot 3720 call `MatMatSolveTranspose`(A,X,X). 3721 3722 .seealso: `Mat`, [Matrix Factorization](sec_matfactor), `MatGetFactor()`, `MatSolveTranspose()`, `MatMatSolve()`, `MatLUFactor()`, `MatCholeskyFactor()` 3723 @*/ 3724 PetscErrorCode MatMatSolveTranspose(Mat A, Mat B, Mat X) 3725 { 3726 PetscFunctionBegin; 3727 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 3728 PetscValidType(A, 1); 3729 PetscValidHeaderSpecific(B, MAT_CLASSID, 2); 3730 PetscValidHeaderSpecific(X, MAT_CLASSID, 3); 3731 PetscCheckSameComm(A, 1, B, 2); 3732 PetscCheckSameComm(A, 1, X, 3); 3733 PetscCheck(X != B, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_IDN, "X and B must be different matrices"); 3734 PetscCheck(A->cmap->N == X->rmap->N, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_SIZ, "Mat A,Mat X: global dim %" PetscInt_FMT " %" PetscInt_FMT, A->cmap->N, X->rmap->N); 3735 PetscCheck(A->rmap->N == B->rmap->N, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_SIZ, "Mat A,Mat B: global dim %" PetscInt_FMT " %" PetscInt_FMT, A->rmap->N, B->rmap->N); 3736 PetscCheck(A->rmap->n == B->rmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Mat A,Mat B: local dim %" PetscInt_FMT " %" PetscInt_FMT, A->rmap->n, B->rmap->n); 3737 PetscCheck(X->cmap->N >= B->cmap->N, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Solution matrix must have same number of columns as rhs matrix"); 3738 if (!A->rmap->N && !A->cmap->N) PetscFunctionReturn(PETSC_SUCCESS); 3739 PetscCheck(A->factortype, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "Unfactored matrix"); 3740 MatCheckPreallocated(A, 1); 3741 3742 PetscCall(PetscLogEventBegin(MAT_MatSolve, A, B, X, 0)); 3743 if (!A->ops->matsolvetranspose) { 3744 PetscCall(PetscInfo(A, "Mat type %s using basic MatMatSolveTranspose\n", ((PetscObject)A)->type_name)); 3745 PetscCall(MatMatSolve_Basic(A, B, X, PETSC_TRUE)); 3746 } else PetscUseTypeMethod(A, matsolvetranspose, B, X); 3747 PetscCall(PetscLogEventEnd(MAT_MatSolve, A, B, X, 0)); 3748 PetscCall(PetscObjectStateIncrease((PetscObject)X)); 3749 PetscFunctionReturn(PETSC_SUCCESS); 3750 } 3751 3752 /*@ 3753 MatMatTransposeSolve - Solves A X = B^T, given a factored matrix. 3754 3755 Neighbor-wise Collective 3756 3757 Input Parameters: 3758 + A - the factored matrix 3759 - Bt - the transpose of right-hand-side matrix as a `MATDENSE` 3760 3761 Output Parameter: 3762 . X - the result matrix (dense matrix) 3763 3764 Level: developer 3765 3766 Note: 3767 For MUMPS, it only supports centralized sparse compressed column format on the host processor for right hand side matrix. User must create B^T in sparse compressed row 3768 format on the host processor and call `MatMatTransposeSolve()` to implement MUMPS' `MatMatSolve()`. 3769 3770 .seealso: `Mat`, [Matrix Factorization](sec_matfactor), `MatMatSolve()`, `MatMatSolveTranspose()`, `MatLUFactor()`, `MatCholeskyFactor()` 3771 @*/ 3772 PetscErrorCode MatMatTransposeSolve(Mat A, Mat Bt, Mat X) 3773 { 3774 PetscFunctionBegin; 3775 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 3776 PetscValidType(A, 1); 3777 PetscValidHeaderSpecific(Bt, MAT_CLASSID, 2); 3778 PetscValidHeaderSpecific(X, MAT_CLASSID, 3); 3779 PetscCheckSameComm(A, 1, Bt, 2); 3780 PetscCheckSameComm(A, 1, X, 3); 3781 3782 PetscCheck(X != Bt, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_IDN, "X and B must be different matrices"); 3783 PetscCheck(A->cmap->N == X->rmap->N, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_SIZ, "Mat A,Mat X: global dim %" PetscInt_FMT " %" PetscInt_FMT, A->cmap->N, X->rmap->N); 3784 PetscCheck(A->rmap->N == Bt->cmap->N, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_SIZ, "Mat A,Mat Bt: global dim %" PetscInt_FMT " %" PetscInt_FMT, A->rmap->N, Bt->cmap->N); 3785 PetscCheck(X->cmap->N >= Bt->rmap->N, PetscObjectComm((PetscObject)X), PETSC_ERR_ARG_SIZ, "Solution matrix must have same number of columns as row number of the rhs matrix"); 3786 if (!A->rmap->N && !A->cmap->N) PetscFunctionReturn(PETSC_SUCCESS); 3787 PetscCheck(A->factortype, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "Unfactored matrix"); 3788 MatCheckPreallocated(A, 1); 3789 3790 PetscCall(PetscLogEventBegin(MAT_MatTrSolve, A, Bt, X, 0)); 3791 PetscUseTypeMethod(A, mattransposesolve, Bt, X); 3792 PetscCall(PetscLogEventEnd(MAT_MatTrSolve, A, Bt, X, 0)); 3793 PetscCall(PetscObjectStateIncrease((PetscObject)X)); 3794 PetscFunctionReturn(PETSC_SUCCESS); 3795 } 3796 3797 /*@ 3798 MatForwardSolve - Solves L x = b, given a factored matrix, A = LU, or 3799 U^T*D^(1/2) x = b, given a factored symmetric matrix, A = U^T*D*U, 3800 3801 Neighbor-wise Collective 3802 3803 Input Parameters: 3804 + mat - the factored matrix 3805 - b - the right-hand-side vector 3806 3807 Output Parameter: 3808 . x - the result vector 3809 3810 Level: developer 3811 3812 Notes: 3813 `MatSolve()` should be used for most applications, as it performs 3814 a forward solve followed by a backward solve. 3815 3816 The vectors `b` and `x` cannot be the same, i.e., one cannot 3817 call `MatForwardSolve`(A,x,x). 3818 3819 For matrix in `MATSEQBAIJ` format with block size larger than 1, 3820 the diagonal blocks are not implemented as D = D^(1/2) * D^(1/2) yet. 3821 `MatForwardSolve()` solves U^T*D y = b, and 3822 `MatBackwardSolve()` solves U x = y. 3823 Thus they do not provide a symmetric preconditioner. 3824 3825 .seealso: `Mat`, `MatBackwardSolve()`, `MatGetFactor()`, `MatSolve()`, `MatBackwardSolve()` 3826 @*/ 3827 PetscErrorCode MatForwardSolve(Mat mat, Vec b, Vec x) 3828 { 3829 PetscFunctionBegin; 3830 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 3831 PetscValidType(mat, 1); 3832 PetscValidHeaderSpecific(b, VEC_CLASSID, 2); 3833 PetscValidHeaderSpecific(x, VEC_CLASSID, 3); 3834 PetscCheckSameComm(mat, 1, b, 2); 3835 PetscCheckSameComm(mat, 1, x, 3); 3836 PetscCheck(x != b, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_IDN, "x and b must be different vectors"); 3837 PetscCheck(mat->cmap->N == x->map->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_SIZ, "Mat mat,Vec x: global dim %" PetscInt_FMT " %" PetscInt_FMT, mat->cmap->N, x->map->N); 3838 PetscCheck(mat->rmap->N == b->map->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_SIZ, "Mat mat,Vec b: global dim %" PetscInt_FMT " %" PetscInt_FMT, mat->rmap->N, b->map->N); 3839 PetscCheck(mat->rmap->n == b->map->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Mat mat,Vec b: local dim %" PetscInt_FMT " %" PetscInt_FMT, mat->rmap->n, b->map->n); 3840 if (!mat->rmap->N && !mat->cmap->N) PetscFunctionReturn(PETSC_SUCCESS); 3841 MatCheckPreallocated(mat, 1); 3842 3843 PetscCall(PetscLogEventBegin(MAT_ForwardSolve, mat, b, x, 0)); 3844 PetscUseTypeMethod(mat, forwardsolve, b, x); 3845 PetscCall(PetscLogEventEnd(MAT_ForwardSolve, mat, b, x, 0)); 3846 PetscCall(PetscObjectStateIncrease((PetscObject)x)); 3847 PetscFunctionReturn(PETSC_SUCCESS); 3848 } 3849 3850 /*@ 3851 MatBackwardSolve - Solves U x = b, given a factored matrix, A = LU. 3852 D^(1/2) U x = b, given a factored symmetric matrix, A = U^T*D*U, 3853 3854 Neighbor-wise Collective 3855 3856 Input Parameters: 3857 + mat - the factored matrix 3858 - b - the right-hand-side vector 3859 3860 Output Parameter: 3861 . x - the result vector 3862 3863 Level: developer 3864 3865 Notes: 3866 `MatSolve()` should be used for most applications, as it performs 3867 a forward solve followed by a backward solve. 3868 3869 The vectors `b` and `x` cannot be the same. I.e., one cannot 3870 call `MatBackwardSolve`(A,x,x). 3871 3872 For matrix in `MATSEQBAIJ` format with block size larger than 1, 3873 the diagonal blocks are not implemented as D = D^(1/2) * D^(1/2) yet. 3874 `MatForwardSolve()` solves U^T*D y = b, and 3875 `MatBackwardSolve()` solves U x = y. 3876 Thus they do not provide a symmetric preconditioner. 3877 3878 .seealso: `Mat`, `MatForwardSolve()`, `MatGetFactor()`, `MatSolve()`, `MatForwardSolve()` 3879 @*/ 3880 PetscErrorCode MatBackwardSolve(Mat mat, Vec b, Vec x) 3881 { 3882 PetscFunctionBegin; 3883 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 3884 PetscValidType(mat, 1); 3885 PetscValidHeaderSpecific(b, VEC_CLASSID, 2); 3886 PetscValidHeaderSpecific(x, VEC_CLASSID, 3); 3887 PetscCheckSameComm(mat, 1, b, 2); 3888 PetscCheckSameComm(mat, 1, x, 3); 3889 PetscCheck(x != b, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_IDN, "x and b must be different vectors"); 3890 PetscCheck(mat->cmap->N == x->map->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_SIZ, "Mat mat,Vec x: global dim %" PetscInt_FMT " %" PetscInt_FMT, mat->cmap->N, x->map->N); 3891 PetscCheck(mat->rmap->N == b->map->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_SIZ, "Mat mat,Vec b: global dim %" PetscInt_FMT " %" PetscInt_FMT, mat->rmap->N, b->map->N); 3892 PetscCheck(mat->rmap->n == b->map->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Mat mat,Vec b: local dim %" PetscInt_FMT " %" PetscInt_FMT, mat->rmap->n, b->map->n); 3893 if (!mat->rmap->N && !mat->cmap->N) PetscFunctionReturn(PETSC_SUCCESS); 3894 MatCheckPreallocated(mat, 1); 3895 3896 PetscCall(PetscLogEventBegin(MAT_BackwardSolve, mat, b, x, 0)); 3897 PetscUseTypeMethod(mat, backwardsolve, b, x); 3898 PetscCall(PetscLogEventEnd(MAT_BackwardSolve, mat, b, x, 0)); 3899 PetscCall(PetscObjectStateIncrease((PetscObject)x)); 3900 PetscFunctionReturn(PETSC_SUCCESS); 3901 } 3902 3903 /*@ 3904 MatSolveAdd - Computes x = y + inv(A)*b, given a factored matrix. 3905 3906 Neighbor-wise Collective 3907 3908 Input Parameters: 3909 + mat - the factored matrix 3910 . b - the right-hand-side vector 3911 - y - the vector to be added to 3912 3913 Output Parameter: 3914 . x - the result vector 3915 3916 Level: developer 3917 3918 Note: 3919 The vectors `b` and `x` cannot be the same. I.e., one cannot 3920 call `MatSolveAdd`(A,x,y,x). 3921 3922 .seealso: `Mat`, [Matrix Factorization](sec_matfactor), `MatSolve()`, `MatGetFactor()`, `MatSolveTranspose()`, `MatSolveTransposeAdd()` 3923 @*/ 3924 PetscErrorCode MatSolveAdd(Mat mat, Vec b, Vec y, Vec x) 3925 { 3926 PetscScalar one = 1.0; 3927 Vec tmp; 3928 3929 PetscFunctionBegin; 3930 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 3931 PetscValidType(mat, 1); 3932 PetscValidHeaderSpecific(y, VEC_CLASSID, 3); 3933 PetscValidHeaderSpecific(b, VEC_CLASSID, 2); 3934 PetscValidHeaderSpecific(x, VEC_CLASSID, 4); 3935 PetscCheckSameComm(mat, 1, b, 2); 3936 PetscCheckSameComm(mat, 1, y, 3); 3937 PetscCheckSameComm(mat, 1, x, 4); 3938 PetscCheck(x != b, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_IDN, "x and b must be different vectors"); 3939 PetscCheck(mat->cmap->N == x->map->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_SIZ, "Mat mat,Vec x: global dim %" PetscInt_FMT " %" PetscInt_FMT, mat->cmap->N, x->map->N); 3940 PetscCheck(mat->rmap->N == b->map->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_SIZ, "Mat mat,Vec b: global dim %" PetscInt_FMT " %" PetscInt_FMT, mat->rmap->N, b->map->N); 3941 PetscCheck(mat->rmap->N == y->map->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_SIZ, "Mat mat,Vec y: global dim %" PetscInt_FMT " %" PetscInt_FMT, mat->rmap->N, y->map->N); 3942 PetscCheck(mat->rmap->n == b->map->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Mat mat,Vec b: local dim %" PetscInt_FMT " %" PetscInt_FMT, mat->rmap->n, b->map->n); 3943 PetscCheck(x->map->n == y->map->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Vec x,Vec y: local dim %" PetscInt_FMT " %" PetscInt_FMT, x->map->n, y->map->n); 3944 if (!mat->rmap->N && !mat->cmap->N) PetscFunctionReturn(PETSC_SUCCESS); 3945 MatCheckPreallocated(mat, 1); 3946 3947 PetscCall(PetscLogEventBegin(MAT_SolveAdd, mat, b, x, y)); 3948 if (mat->factorerrortype) { 3949 PetscCall(PetscInfo(mat, "MatFactorError %d\n", mat->factorerrortype)); 3950 PetscCall(VecSetInf(x)); 3951 } else if (mat->ops->solveadd) { 3952 PetscUseTypeMethod(mat, solveadd, b, y, x); 3953 } else { 3954 /* do the solve then the add manually */ 3955 if (x != y) { 3956 PetscCall(MatSolve(mat, b, x)); 3957 PetscCall(VecAXPY(x, one, y)); 3958 } else { 3959 PetscCall(VecDuplicate(x, &tmp)); 3960 PetscCall(VecCopy(x, tmp)); 3961 PetscCall(MatSolve(mat, b, x)); 3962 PetscCall(VecAXPY(x, one, tmp)); 3963 PetscCall(VecDestroy(&tmp)); 3964 } 3965 } 3966 PetscCall(PetscLogEventEnd(MAT_SolveAdd, mat, b, x, y)); 3967 PetscCall(PetscObjectStateIncrease((PetscObject)x)); 3968 PetscFunctionReturn(PETSC_SUCCESS); 3969 } 3970 3971 /*@ 3972 MatSolveTranspose - Solves A' x = b, given a factored matrix. 3973 3974 Neighbor-wise Collective 3975 3976 Input Parameters: 3977 + mat - the factored matrix 3978 - b - the right-hand-side vector 3979 3980 Output Parameter: 3981 . x - the result vector 3982 3983 Level: developer 3984 3985 Notes: 3986 The vectors `b` and `x` cannot be the same. I.e., one cannot 3987 call `MatSolveTranspose`(A,x,x). 3988 3989 Most users should employ the `KSP` interface for linear solvers 3990 instead of working directly with matrix algebra routines such as this. 3991 See, e.g., `KSPCreate()`. 3992 3993 .seealso: `Mat`, `MatGetFactor()`, `KSP`, `MatSolve()`, `MatSolveAdd()`, `MatSolveTransposeAdd()` 3994 @*/ 3995 PetscErrorCode MatSolveTranspose(Mat mat, Vec b, Vec x) 3996 { 3997 PetscErrorCode (*f)(Mat, Vec, Vec) = (!mat->ops->solvetranspose && mat->symmetric) ? mat->ops->solve : mat->ops->solvetranspose; 3998 3999 PetscFunctionBegin; 4000 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 4001 PetscValidType(mat, 1); 4002 PetscValidHeaderSpecific(b, VEC_CLASSID, 2); 4003 PetscValidHeaderSpecific(x, VEC_CLASSID, 3); 4004 PetscCheckSameComm(mat, 1, b, 2); 4005 PetscCheckSameComm(mat, 1, x, 3); 4006 PetscCheck(x != b, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_IDN, "x and b must be different vectors"); 4007 PetscCheck(mat->rmap->N == x->map->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_SIZ, "Mat mat,Vec x: global dim %" PetscInt_FMT " %" PetscInt_FMT, mat->rmap->N, x->map->N); 4008 PetscCheck(mat->cmap->N == b->map->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_SIZ, "Mat mat,Vec b: global dim %" PetscInt_FMT " %" PetscInt_FMT, mat->cmap->N, b->map->N); 4009 if (!mat->rmap->N && !mat->cmap->N) PetscFunctionReturn(PETSC_SUCCESS); 4010 MatCheckPreallocated(mat, 1); 4011 PetscCall(PetscLogEventBegin(MAT_SolveTranspose, mat, b, x, 0)); 4012 if (mat->factorerrortype) { 4013 PetscCall(PetscInfo(mat, "MatFactorError %d\n", mat->factorerrortype)); 4014 PetscCall(VecSetInf(x)); 4015 } else { 4016 PetscCheck(f, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Matrix type %s", ((PetscObject)mat)->type_name); 4017 PetscCall((*f)(mat, b, x)); 4018 } 4019 PetscCall(PetscLogEventEnd(MAT_SolveTranspose, mat, b, x, 0)); 4020 PetscCall(PetscObjectStateIncrease((PetscObject)x)); 4021 PetscFunctionReturn(PETSC_SUCCESS); 4022 } 4023 4024 /*@ 4025 MatSolveTransposeAdd - Computes x = y + inv(Transpose(A)) b, given a 4026 factored matrix. 4027 4028 Neighbor-wise Collective 4029 4030 Input Parameters: 4031 + mat - the factored matrix 4032 . b - the right-hand-side vector 4033 - y - the vector to be added to 4034 4035 Output Parameter: 4036 . x - the result vector 4037 4038 Level: developer 4039 4040 Note: 4041 The vectors `b` and `x` cannot be the same. I.e., one cannot 4042 call `MatSolveTransposeAdd`(A,x,y,x). 4043 4044 .seealso: `Mat`, `MatGetFactor()`, `MatSolve()`, `MatSolveAdd()`, `MatSolveTranspose()` 4045 @*/ 4046 PetscErrorCode MatSolveTransposeAdd(Mat mat, Vec b, Vec y, Vec x) 4047 { 4048 PetscScalar one = 1.0; 4049 Vec tmp; 4050 PetscErrorCode (*f)(Mat, Vec, Vec, Vec) = (!mat->ops->solvetransposeadd && mat->symmetric) ? mat->ops->solveadd : mat->ops->solvetransposeadd; 4051 4052 PetscFunctionBegin; 4053 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 4054 PetscValidType(mat, 1); 4055 PetscValidHeaderSpecific(y, VEC_CLASSID, 3); 4056 PetscValidHeaderSpecific(b, VEC_CLASSID, 2); 4057 PetscValidHeaderSpecific(x, VEC_CLASSID, 4); 4058 PetscCheckSameComm(mat, 1, b, 2); 4059 PetscCheckSameComm(mat, 1, y, 3); 4060 PetscCheckSameComm(mat, 1, x, 4); 4061 PetscCheck(x != b, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_IDN, "x and b must be different vectors"); 4062 PetscCheck(mat->rmap->N == x->map->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_SIZ, "Mat mat,Vec x: global dim %" PetscInt_FMT " %" PetscInt_FMT, mat->rmap->N, x->map->N); 4063 PetscCheck(mat->cmap->N == b->map->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_SIZ, "Mat mat,Vec b: global dim %" PetscInt_FMT " %" PetscInt_FMT, mat->cmap->N, b->map->N); 4064 PetscCheck(mat->cmap->N == y->map->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_SIZ, "Mat mat,Vec y: global dim %" PetscInt_FMT " %" PetscInt_FMT, mat->cmap->N, y->map->N); 4065 PetscCheck(x->map->n == y->map->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Vec x,Vec y: local dim %" PetscInt_FMT " %" PetscInt_FMT, x->map->n, y->map->n); 4066 if (!mat->rmap->N && !mat->cmap->N) PetscFunctionReturn(PETSC_SUCCESS); 4067 MatCheckPreallocated(mat, 1); 4068 4069 PetscCall(PetscLogEventBegin(MAT_SolveTransposeAdd, mat, b, x, y)); 4070 if (mat->factorerrortype) { 4071 PetscCall(PetscInfo(mat, "MatFactorError %d\n", mat->factorerrortype)); 4072 PetscCall(VecSetInf(x)); 4073 } else if (f) { 4074 PetscCall((*f)(mat, b, y, x)); 4075 } else { 4076 /* do the solve then the add manually */ 4077 if (x != y) { 4078 PetscCall(MatSolveTranspose(mat, b, x)); 4079 PetscCall(VecAXPY(x, one, y)); 4080 } else { 4081 PetscCall(VecDuplicate(x, &tmp)); 4082 PetscCall(VecCopy(x, tmp)); 4083 PetscCall(MatSolveTranspose(mat, b, x)); 4084 PetscCall(VecAXPY(x, one, tmp)); 4085 PetscCall(VecDestroy(&tmp)); 4086 } 4087 } 4088 PetscCall(PetscLogEventEnd(MAT_SolveTransposeAdd, mat, b, x, y)); 4089 PetscCall(PetscObjectStateIncrease((PetscObject)x)); 4090 PetscFunctionReturn(PETSC_SUCCESS); 4091 } 4092 4093 /*@ 4094 MatSOR - Computes relaxation (SOR, Gauss-Seidel) sweeps. 4095 4096 Neighbor-wise Collective 4097 4098 Input Parameters: 4099 + mat - the matrix 4100 . b - the right hand side 4101 . omega - the relaxation factor 4102 . flag - flag indicating the type of SOR (see below) 4103 . shift - diagonal shift 4104 . its - the number of iterations 4105 - lits - the number of local iterations 4106 4107 Output Parameter: 4108 . x - the solution (can contain an initial guess, use option `SOR_ZERO_INITIAL_GUESS` to indicate no guess) 4109 4110 SOR Flags: 4111 + `SOR_FORWARD_SWEEP` - forward SOR 4112 . `SOR_BACKWARD_SWEEP` - backward SOR 4113 . `SOR_SYMMETRIC_SWEEP` - SSOR (symmetric SOR) 4114 . `SOR_LOCAL_FORWARD_SWEEP` - local forward SOR 4115 . `SOR_LOCAL_BACKWARD_SWEEP` - local forward SOR 4116 . `SOR_LOCAL_SYMMETRIC_SWEEP` - local SSOR 4117 . `SOR_EISENSTAT` - SOR with Eisenstat trick 4118 . `SOR_APPLY_UPPER`, `SOR_APPLY_LOWER` - applies 4119 upper/lower triangular part of matrix to 4120 vector (with omega) 4121 - `SOR_ZERO_INITIAL_GUESS` - zero initial guess 4122 4123 Level: developer 4124 4125 Notes: 4126 `SOR_LOCAL_FORWARD_SWEEP`, `SOR_LOCAL_BACKWARD_SWEEP`, and 4127 `SOR_LOCAL_SYMMETRIC_SWEEP` perform separate independent smoothings 4128 on each processor. 4129 4130 Application programmers will not generally use `MatSOR()` directly, 4131 but instead will employ the `KSP`/`PC` interface. 4132 4133 For `MATBAIJ`, `MATSBAIJ`, and `MATAIJ` matrices with Inodes this does a block SOR smoothing, otherwise it does a pointwise smoothing 4134 4135 Most users should employ the `KSP` interface for linear solvers 4136 instead of working directly with matrix algebra routines such as this. 4137 See, e.g., `KSPCreate()`. 4138 4139 Vectors `x` and `b` CANNOT be the same 4140 4141 Notes for Advanced Users: 4142 The flags are implemented as bitwise inclusive or operations. 4143 For example, use (`SOR_ZERO_INITIAL_GUESS` | `SOR_SYMMETRIC_SWEEP`) 4144 to specify a zero initial guess for SSOR. 4145 4146 Developer Note: 4147 We should add block SOR support for `MATAIJ` matrices with block size set to great than one and no inodes 4148 4149 .seealso: `Mat`, `MatMult()`, `KSP`, `PC`, `MatGetFactor()` 4150 @*/ 4151 PetscErrorCode MatSOR(Mat mat, Vec b, PetscReal omega, MatSORType flag, PetscReal shift, PetscInt its, PetscInt lits, Vec x) 4152 { 4153 PetscFunctionBegin; 4154 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 4155 PetscValidType(mat, 1); 4156 PetscValidHeaderSpecific(b, VEC_CLASSID, 2); 4157 PetscValidHeaderSpecific(x, VEC_CLASSID, 8); 4158 PetscCheckSameComm(mat, 1, b, 2); 4159 PetscCheckSameComm(mat, 1, x, 8); 4160 PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 4161 PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 4162 PetscCheck(mat->cmap->N == x->map->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_SIZ, "Mat mat,Vec x: global dim %" PetscInt_FMT " %" PetscInt_FMT, mat->cmap->N, x->map->N); 4163 PetscCheck(mat->rmap->N == b->map->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_SIZ, "Mat mat,Vec b: global dim %" PetscInt_FMT " %" PetscInt_FMT, mat->rmap->N, b->map->N); 4164 PetscCheck(mat->rmap->n == b->map->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Mat mat,Vec b: local dim %" PetscInt_FMT " %" PetscInt_FMT, mat->rmap->n, b->map->n); 4165 PetscCheck(its > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Relaxation requires global its %" PetscInt_FMT " positive", its); 4166 PetscCheck(lits > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Relaxation requires local its %" PetscInt_FMT " positive", lits); 4167 PetscCheck(b != x, PETSC_COMM_SELF, PETSC_ERR_ARG_IDN, "b and x vector cannot be the same"); 4168 4169 MatCheckPreallocated(mat, 1); 4170 PetscCall(PetscLogEventBegin(MAT_SOR, mat, b, x, 0)); 4171 PetscUseTypeMethod(mat, sor, b, omega, flag, shift, its, lits, x); 4172 PetscCall(PetscLogEventEnd(MAT_SOR, mat, b, x, 0)); 4173 PetscCall(PetscObjectStateIncrease((PetscObject)x)); 4174 PetscFunctionReturn(PETSC_SUCCESS); 4175 } 4176 4177 /* 4178 Default matrix copy routine. 4179 */ 4180 PetscErrorCode MatCopy_Basic(Mat A, Mat B, MatStructure str) 4181 { 4182 PetscInt i, rstart = 0, rend = 0, nz; 4183 const PetscInt *cwork; 4184 const PetscScalar *vwork; 4185 4186 PetscFunctionBegin; 4187 if (B->assembled) PetscCall(MatZeroEntries(B)); 4188 if (str == SAME_NONZERO_PATTERN) { 4189 PetscCall(MatGetOwnershipRange(A, &rstart, &rend)); 4190 for (i = rstart; i < rend; i++) { 4191 PetscCall(MatGetRow(A, i, &nz, &cwork, &vwork)); 4192 PetscCall(MatSetValues(B, 1, &i, nz, cwork, vwork, INSERT_VALUES)); 4193 PetscCall(MatRestoreRow(A, i, &nz, &cwork, &vwork)); 4194 } 4195 } else { 4196 PetscCall(MatAYPX(B, 0.0, A, str)); 4197 } 4198 PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 4199 PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 4200 PetscFunctionReturn(PETSC_SUCCESS); 4201 } 4202 4203 /*@ 4204 MatCopy - Copies a matrix to another matrix. 4205 4206 Collective 4207 4208 Input Parameters: 4209 + A - the matrix 4210 - str - `SAME_NONZERO_PATTERN` or `DIFFERENT_NONZERO_PATTERN` 4211 4212 Output Parameter: 4213 . B - where the copy is put 4214 4215 Level: intermediate 4216 4217 Notes: 4218 If you use `SAME_NONZERO_PATTERN `then the two matrices must have the same nonzero pattern or the routine will crash. 4219 4220 `MatCopy()` copies the matrix entries of a matrix to another existing 4221 matrix (after first zeroing the second matrix). A related routine is 4222 `MatConvert()`, which first creates a new matrix and then copies the data. 4223 4224 .seealso: `Mat`, `MatConvert()`, `MatDuplicate()` 4225 @*/ 4226 PetscErrorCode MatCopy(Mat A, Mat B, MatStructure str) 4227 { 4228 PetscInt i; 4229 4230 PetscFunctionBegin; 4231 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 4232 PetscValidHeaderSpecific(B, MAT_CLASSID, 2); 4233 PetscValidType(A, 1); 4234 PetscValidType(B, 2); 4235 PetscCheckSameComm(A, 1, B, 2); 4236 MatCheckPreallocated(B, 2); 4237 PetscCheck(A->assembled, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 4238 PetscCheck(!A->factortype, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 4239 PetscCheck(A->rmap->N == B->rmap->N && A->cmap->N == B->cmap->N, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_SIZ, "Mat A,Mat B: global dim (%" PetscInt_FMT ",%" PetscInt_FMT ") (%" PetscInt_FMT ",%" PetscInt_FMT ")", A->rmap->N, B->rmap->N, 4240 A->cmap->N, B->cmap->N); 4241 MatCheckPreallocated(A, 1); 4242 if (A == B) PetscFunctionReturn(PETSC_SUCCESS); 4243 4244 PetscCall(PetscLogEventBegin(MAT_Copy, A, B, 0, 0)); 4245 if (A->ops->copy) PetscUseTypeMethod(A, copy, B, str); 4246 else PetscCall(MatCopy_Basic(A, B, str)); 4247 4248 B->stencil.dim = A->stencil.dim; 4249 B->stencil.noc = A->stencil.noc; 4250 for (i = 0; i <= A->stencil.dim; i++) { 4251 B->stencil.dims[i] = A->stencil.dims[i]; 4252 B->stencil.starts[i] = A->stencil.starts[i]; 4253 } 4254 4255 PetscCall(PetscLogEventEnd(MAT_Copy, A, B, 0, 0)); 4256 PetscCall(PetscObjectStateIncrease((PetscObject)B)); 4257 PetscFunctionReturn(PETSC_SUCCESS); 4258 } 4259 4260 /*@C 4261 MatConvert - Converts a matrix to another matrix, either of the same 4262 or different type. 4263 4264 Collective 4265 4266 Input Parameters: 4267 + mat - the matrix 4268 . newtype - new matrix type. Use `MATSAME` to create a new matrix of the 4269 same type as the original matrix. 4270 - reuse - denotes if the destination matrix is to be created or reused. 4271 Use `MAT_INPLACE_MATRIX` for inplace conversion (that is when you want the input mat to be changed to contain the matrix in the new format), otherwise use 4272 `MAT_INITIAL_MATRIX` or `MAT_REUSE_MATRIX` (can only be used after the first call was made with `MAT_INITIAL_MATRIX`, causes the matrix space in M to be reused). 4273 4274 Output Parameter: 4275 . M - pointer to place new matrix 4276 4277 Level: intermediate 4278 4279 Notes: 4280 `MatConvert()` first creates a new matrix and then copies the data from 4281 the first matrix. A related routine is `MatCopy()`, which copies the matrix 4282 entries of one matrix to another already existing matrix context. 4283 4284 Cannot be used to convert a sequential matrix to parallel or parallel to sequential, 4285 the MPI communicator of the generated matrix is always the same as the communicator 4286 of the input matrix. 4287 4288 .seealso: `Mat`, `MatCopy()`, `MatDuplicate()`, `MAT_INITIAL_MATRIX`, `MAT_REUSE_MATRIX`, `MAT_INPLACE_MATRIX` 4289 @*/ 4290 PetscErrorCode MatConvert(Mat mat, MatType newtype, MatReuse reuse, Mat *M) 4291 { 4292 PetscBool sametype, issame, flg; 4293 PetscBool3 issymmetric, ishermitian; 4294 char convname[256], mtype[256]; 4295 Mat B; 4296 4297 PetscFunctionBegin; 4298 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 4299 PetscValidType(mat, 1); 4300 PetscValidPointer(M, 4); 4301 PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 4302 PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 4303 MatCheckPreallocated(mat, 1); 4304 4305 PetscCall(PetscOptionsGetString(((PetscObject)mat)->options, ((PetscObject)mat)->prefix, "-matconvert_type", mtype, sizeof(mtype), &flg)); 4306 if (flg) newtype = mtype; 4307 4308 PetscCall(PetscObjectTypeCompare((PetscObject)mat, newtype, &sametype)); 4309 PetscCall(PetscStrcmp(newtype, "same", &issame)); 4310 PetscCheck(!(reuse == MAT_INPLACE_MATRIX) || !(mat != *M), PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "MAT_INPLACE_MATRIX requires same input and output matrix"); 4311 PetscCheck(!(reuse == MAT_REUSE_MATRIX) || !(mat == *M), PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "MAT_REUSE_MATRIX means reuse matrix in final argument, perhaps you mean MAT_INPLACE_MATRIX"); 4312 4313 if ((reuse == MAT_INPLACE_MATRIX) && (issame || sametype)) { 4314 PetscCall(PetscInfo(mat, "Early return for inplace %s %d %d\n", ((PetscObject)mat)->type_name, sametype, issame)); 4315 PetscFunctionReturn(PETSC_SUCCESS); 4316 } 4317 4318 /* Cache Mat options because some converters use MatHeaderReplace */ 4319 issymmetric = mat->symmetric; 4320 ishermitian = mat->hermitian; 4321 4322 if ((sametype || issame) && (reuse == MAT_INITIAL_MATRIX) && mat->ops->duplicate) { 4323 PetscCall(PetscInfo(mat, "Calling duplicate for initial matrix %s %d %d\n", ((PetscObject)mat)->type_name, sametype, issame)); 4324 PetscUseTypeMethod(mat, duplicate, MAT_COPY_VALUES, M); 4325 } else { 4326 PetscErrorCode (*conv)(Mat, MatType, MatReuse, Mat *) = NULL; 4327 const char *prefix[3] = {"seq", "mpi", ""}; 4328 PetscInt i; 4329 /* 4330 Order of precedence: 4331 0) See if newtype is a superclass of the current matrix. 4332 1) See if a specialized converter is known to the current matrix. 4333 2) See if a specialized converter is known to the desired matrix class. 4334 3) See if a good general converter is registered for the desired class 4335 (as of 6/27/03 only MATMPIADJ falls into this category). 4336 4) See if a good general converter is known for the current matrix. 4337 5) Use a really basic converter. 4338 */ 4339 4340 /* 0) See if newtype is a superclass of the current matrix. 4341 i.e mat is mpiaij and newtype is aij */ 4342 for (i = 0; i < 2; i++) { 4343 PetscCall(PetscStrncpy(convname, prefix[i], sizeof(convname))); 4344 PetscCall(PetscStrlcat(convname, newtype, sizeof(convname))); 4345 PetscCall(PetscStrcmp(convname, ((PetscObject)mat)->type_name, &flg)); 4346 PetscCall(PetscInfo(mat, "Check superclass %s %s -> %d\n", convname, ((PetscObject)mat)->type_name, flg)); 4347 if (flg) { 4348 if (reuse == MAT_INPLACE_MATRIX) { 4349 PetscCall(PetscInfo(mat, "Early return\n")); 4350 PetscFunctionReturn(PETSC_SUCCESS); 4351 } else if (reuse == MAT_INITIAL_MATRIX && mat->ops->duplicate) { 4352 PetscCall(PetscInfo(mat, "Calling MatDuplicate\n")); 4353 PetscUseTypeMethod(mat, duplicate, MAT_COPY_VALUES, M); 4354 PetscFunctionReturn(PETSC_SUCCESS); 4355 } else if (reuse == MAT_REUSE_MATRIX && mat->ops->copy) { 4356 PetscCall(PetscInfo(mat, "Calling MatCopy\n")); 4357 PetscCall(MatCopy(mat, *M, SAME_NONZERO_PATTERN)); 4358 PetscFunctionReturn(PETSC_SUCCESS); 4359 } 4360 } 4361 } 4362 /* 1) See if a specialized converter is known to the current matrix and the desired class */ 4363 for (i = 0; i < 3; i++) { 4364 PetscCall(PetscStrncpy(convname, "MatConvert_", sizeof(convname))); 4365 PetscCall(PetscStrlcat(convname, ((PetscObject)mat)->type_name, sizeof(convname))); 4366 PetscCall(PetscStrlcat(convname, "_", sizeof(convname))); 4367 PetscCall(PetscStrlcat(convname, prefix[i], sizeof(convname))); 4368 PetscCall(PetscStrlcat(convname, issame ? ((PetscObject)mat)->type_name : newtype, sizeof(convname))); 4369 PetscCall(PetscStrlcat(convname, "_C", sizeof(convname))); 4370 PetscCall(PetscObjectQueryFunction((PetscObject)mat, convname, &conv)); 4371 PetscCall(PetscInfo(mat, "Check specialized (1) %s (%s) -> %d\n", convname, ((PetscObject)mat)->type_name, !!conv)); 4372 if (conv) goto foundconv; 4373 } 4374 4375 /* 2) See if a specialized converter is known to the desired matrix class. */ 4376 PetscCall(MatCreate(PetscObjectComm((PetscObject)mat), &B)); 4377 PetscCall(MatSetSizes(B, mat->rmap->n, mat->cmap->n, mat->rmap->N, mat->cmap->N)); 4378 PetscCall(MatSetType(B, newtype)); 4379 for (i = 0; i < 3; i++) { 4380 PetscCall(PetscStrncpy(convname, "MatConvert_", sizeof(convname))); 4381 PetscCall(PetscStrlcat(convname, ((PetscObject)mat)->type_name, sizeof(convname))); 4382 PetscCall(PetscStrlcat(convname, "_", sizeof(convname))); 4383 PetscCall(PetscStrlcat(convname, prefix[i], sizeof(convname))); 4384 PetscCall(PetscStrlcat(convname, newtype, sizeof(convname))); 4385 PetscCall(PetscStrlcat(convname, "_C", sizeof(convname))); 4386 PetscCall(PetscObjectQueryFunction((PetscObject)B, convname, &conv)); 4387 PetscCall(PetscInfo(mat, "Check specialized (2) %s (%s) -> %d\n", convname, ((PetscObject)B)->type_name, !!conv)); 4388 if (conv) { 4389 PetscCall(MatDestroy(&B)); 4390 goto foundconv; 4391 } 4392 } 4393 4394 /* 3) See if a good general converter is registered for the desired class */ 4395 conv = B->ops->convertfrom; 4396 PetscCall(PetscInfo(mat, "Check convertfrom (%s) -> %d\n", ((PetscObject)B)->type_name, !!conv)); 4397 PetscCall(MatDestroy(&B)); 4398 if (conv) goto foundconv; 4399 4400 /* 4) See if a good general converter is known for the current matrix */ 4401 if (mat->ops->convert) conv = mat->ops->convert; 4402 PetscCall(PetscInfo(mat, "Check general convert (%s) -> %d\n", ((PetscObject)mat)->type_name, !!conv)); 4403 if (conv) goto foundconv; 4404 4405 /* 5) Use a really basic converter. */ 4406 PetscCall(PetscInfo(mat, "Using MatConvert_Basic\n")); 4407 conv = MatConvert_Basic; 4408 4409 foundconv: 4410 PetscCall(PetscLogEventBegin(MAT_Convert, mat, 0, 0, 0)); 4411 PetscCall((*conv)(mat, newtype, reuse, M)); 4412 if (mat->rmap->mapping && mat->cmap->mapping && !(*M)->rmap->mapping && !(*M)->cmap->mapping) { 4413 /* the block sizes must be same if the mappings are copied over */ 4414 (*M)->rmap->bs = mat->rmap->bs; 4415 (*M)->cmap->bs = mat->cmap->bs; 4416 PetscCall(PetscObjectReference((PetscObject)mat->rmap->mapping)); 4417 PetscCall(PetscObjectReference((PetscObject)mat->cmap->mapping)); 4418 (*M)->rmap->mapping = mat->rmap->mapping; 4419 (*M)->cmap->mapping = mat->cmap->mapping; 4420 } 4421 (*M)->stencil.dim = mat->stencil.dim; 4422 (*M)->stencil.noc = mat->stencil.noc; 4423 for (i = 0; i <= mat->stencil.dim; i++) { 4424 (*M)->stencil.dims[i] = mat->stencil.dims[i]; 4425 (*M)->stencil.starts[i] = mat->stencil.starts[i]; 4426 } 4427 PetscCall(PetscLogEventEnd(MAT_Convert, mat, 0, 0, 0)); 4428 } 4429 PetscCall(PetscObjectStateIncrease((PetscObject)*M)); 4430 4431 /* Copy Mat options */ 4432 if (issymmetric == PETSC_BOOL3_TRUE) PetscCall(MatSetOption(*M, MAT_SYMMETRIC, PETSC_TRUE)); 4433 else if (issymmetric == PETSC_BOOL3_FALSE) PetscCall(MatSetOption(*M, MAT_SYMMETRIC, PETSC_FALSE)); 4434 if (ishermitian == PETSC_BOOL3_TRUE) PetscCall(MatSetOption(*M, MAT_HERMITIAN, PETSC_TRUE)); 4435 else if (ishermitian == PETSC_BOOL3_FALSE) PetscCall(MatSetOption(*M, MAT_HERMITIAN, PETSC_FALSE)); 4436 PetscFunctionReturn(PETSC_SUCCESS); 4437 } 4438 4439 /*@C 4440 MatFactorGetSolverType - Returns name of the package providing the factorization routines 4441 4442 Not Collective 4443 4444 Input Parameter: 4445 . mat - the matrix, must be a factored matrix 4446 4447 Output Parameter: 4448 . type - the string name of the package (do not free this string) 4449 4450 Level: intermediate 4451 4452 Fortran Note: 4453 Pass in an empty string and the package name will be copied into it. Make sure the string is long enough. 4454 4455 .seealso: `Mat`, [Matrix Factorization](sec_matfactor), `MatGetFactor()`, `MatSolverType`, `MatCopy()`, `MatDuplicate()`, `MatGetFactorAvailable()`, `MatGetFactor()` 4456 @*/ 4457 PetscErrorCode MatFactorGetSolverType(Mat mat, MatSolverType *type) 4458 { 4459 PetscErrorCode (*conv)(Mat, MatSolverType *); 4460 4461 PetscFunctionBegin; 4462 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 4463 PetscValidType(mat, 1); 4464 PetscValidPointer(type, 2); 4465 PetscCheck(mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Only for factored matrix"); 4466 PetscCall(PetscObjectQueryFunction((PetscObject)mat, "MatFactorGetSolverType_C", &conv)); 4467 if (conv) PetscCall((*conv)(mat, type)); 4468 else *type = MATSOLVERPETSC; 4469 PetscFunctionReturn(PETSC_SUCCESS); 4470 } 4471 4472 typedef struct _MatSolverTypeForSpecifcType *MatSolverTypeForSpecifcType; 4473 struct _MatSolverTypeForSpecifcType { 4474 MatType mtype; 4475 /* no entry for MAT_FACTOR_NONE */ 4476 PetscErrorCode (*createfactor[MAT_FACTOR_NUM_TYPES - 1])(Mat, MatFactorType, Mat *); 4477 MatSolverTypeForSpecifcType next; 4478 }; 4479 4480 typedef struct _MatSolverTypeHolder *MatSolverTypeHolder; 4481 struct _MatSolverTypeHolder { 4482 char *name; 4483 MatSolverTypeForSpecifcType handlers; 4484 MatSolverTypeHolder next; 4485 }; 4486 4487 static MatSolverTypeHolder MatSolverTypeHolders = NULL; 4488 4489 /*@C 4490 MatSolverTypeRegister - Registers a `MatSolverType` that works for a particular matrix type 4491 4492 Input Parameters: 4493 + package - name of the package, for example petsc or superlu 4494 . mtype - the matrix type that works with this package 4495 . ftype - the type of factorization supported by the package 4496 - createfactor - routine that will create the factored matrix ready to be used 4497 4498 Level: developer 4499 4500 .seealso: `Mat`, [Matrix Factorization](sec_matfactor), `MatFactorGetSolverType()`, `MatCopy()`, `MatDuplicate()`, `MatGetFactorAvailable()`, `MatGetFactor()` 4501 @*/ 4502 PetscErrorCode MatSolverTypeRegister(MatSolverType package, MatType mtype, MatFactorType ftype, PetscErrorCode (*createfactor)(Mat, MatFactorType, Mat *)) 4503 { 4504 MatSolverTypeHolder next = MatSolverTypeHolders, prev = NULL; 4505 PetscBool flg; 4506 MatSolverTypeForSpecifcType inext, iprev = NULL; 4507 4508 PetscFunctionBegin; 4509 PetscCall(MatInitializePackage()); 4510 if (!next) { 4511 PetscCall(PetscNew(&MatSolverTypeHolders)); 4512 PetscCall(PetscStrallocpy(package, &MatSolverTypeHolders->name)); 4513 PetscCall(PetscNew(&MatSolverTypeHolders->handlers)); 4514 PetscCall(PetscStrallocpy(mtype, (char **)&MatSolverTypeHolders->handlers->mtype)); 4515 MatSolverTypeHolders->handlers->createfactor[(int)ftype - 1] = createfactor; 4516 PetscFunctionReturn(PETSC_SUCCESS); 4517 } 4518 while (next) { 4519 PetscCall(PetscStrcasecmp(package, next->name, &flg)); 4520 if (flg) { 4521 PetscCheck(next->handlers, PETSC_COMM_SELF, PETSC_ERR_PLIB, "MatSolverTypeHolder is missing handlers"); 4522 inext = next->handlers; 4523 while (inext) { 4524 PetscCall(PetscStrcasecmp(mtype, inext->mtype, &flg)); 4525 if (flg) { 4526 inext->createfactor[(int)ftype - 1] = createfactor; 4527 PetscFunctionReturn(PETSC_SUCCESS); 4528 } 4529 iprev = inext; 4530 inext = inext->next; 4531 } 4532 PetscCall(PetscNew(&iprev->next)); 4533 PetscCall(PetscStrallocpy(mtype, (char **)&iprev->next->mtype)); 4534 iprev->next->createfactor[(int)ftype - 1] = createfactor; 4535 PetscFunctionReturn(PETSC_SUCCESS); 4536 } 4537 prev = next; 4538 next = next->next; 4539 } 4540 PetscCall(PetscNew(&prev->next)); 4541 PetscCall(PetscStrallocpy(package, &prev->next->name)); 4542 PetscCall(PetscNew(&prev->next->handlers)); 4543 PetscCall(PetscStrallocpy(mtype, (char **)&prev->next->handlers->mtype)); 4544 prev->next->handlers->createfactor[(int)ftype - 1] = createfactor; 4545 PetscFunctionReturn(PETSC_SUCCESS); 4546 } 4547 4548 /*@C 4549 MatSolverTypeGet - Gets the function that creates the factor matrix if it exist 4550 4551 Input Parameters: 4552 + type - name of the package, for example petsc or superlu 4553 . ftype - the type of factorization supported by the type 4554 - mtype - the matrix type that works with this type 4555 4556 Output Parameters: 4557 + foundtype - `PETSC_TRUE` if the type was registered 4558 . foundmtype - `PETSC_TRUE` if the type supports the requested mtype 4559 - createfactor - routine that will create the factored matrix ready to be used or `NULL` if not found 4560 4561 Level: developer 4562 4563 .seealso: `Mat`, `MatFactorType`, `MatType`, `MatCopy()`, `MatDuplicate()`, `MatGetFactorAvailable()`, `MatSolverTypeRegister()`, `MatGetFactor()` 4564 @*/ 4565 PetscErrorCode MatSolverTypeGet(MatSolverType type, MatType mtype, MatFactorType ftype, PetscBool *foundtype, PetscBool *foundmtype, PetscErrorCode (**createfactor)(Mat, MatFactorType, Mat *)) 4566 { 4567 MatSolverTypeHolder next = MatSolverTypeHolders; 4568 PetscBool flg; 4569 MatSolverTypeForSpecifcType inext; 4570 4571 PetscFunctionBegin; 4572 if (foundtype) *foundtype = PETSC_FALSE; 4573 if (foundmtype) *foundmtype = PETSC_FALSE; 4574 if (createfactor) *createfactor = NULL; 4575 4576 if (type) { 4577 while (next) { 4578 PetscCall(PetscStrcasecmp(type, next->name, &flg)); 4579 if (flg) { 4580 if (foundtype) *foundtype = PETSC_TRUE; 4581 inext = next->handlers; 4582 while (inext) { 4583 PetscCall(PetscStrbeginswith(mtype, inext->mtype, &flg)); 4584 if (flg) { 4585 if (foundmtype) *foundmtype = PETSC_TRUE; 4586 if (createfactor) *createfactor = inext->createfactor[(int)ftype - 1]; 4587 PetscFunctionReturn(PETSC_SUCCESS); 4588 } 4589 inext = inext->next; 4590 } 4591 } 4592 next = next->next; 4593 } 4594 } else { 4595 while (next) { 4596 inext = next->handlers; 4597 while (inext) { 4598 PetscCall(PetscStrcmp(mtype, inext->mtype, &flg)); 4599 if (flg && inext->createfactor[(int)ftype - 1]) { 4600 if (foundtype) *foundtype = PETSC_TRUE; 4601 if (foundmtype) *foundmtype = PETSC_TRUE; 4602 if (createfactor) *createfactor = inext->createfactor[(int)ftype - 1]; 4603 PetscFunctionReturn(PETSC_SUCCESS); 4604 } 4605 inext = inext->next; 4606 } 4607 next = next->next; 4608 } 4609 /* try with base classes inext->mtype */ 4610 next = MatSolverTypeHolders; 4611 while (next) { 4612 inext = next->handlers; 4613 while (inext) { 4614 PetscCall(PetscStrbeginswith(mtype, inext->mtype, &flg)); 4615 if (flg && inext->createfactor[(int)ftype - 1]) { 4616 if (foundtype) *foundtype = PETSC_TRUE; 4617 if (foundmtype) *foundmtype = PETSC_TRUE; 4618 if (createfactor) *createfactor = inext->createfactor[(int)ftype - 1]; 4619 PetscFunctionReturn(PETSC_SUCCESS); 4620 } 4621 inext = inext->next; 4622 } 4623 next = next->next; 4624 } 4625 } 4626 PetscFunctionReturn(PETSC_SUCCESS); 4627 } 4628 4629 PetscErrorCode MatSolverTypeDestroy(void) 4630 { 4631 MatSolverTypeHolder next = MatSolverTypeHolders, prev; 4632 MatSolverTypeForSpecifcType inext, iprev; 4633 4634 PetscFunctionBegin; 4635 while (next) { 4636 PetscCall(PetscFree(next->name)); 4637 inext = next->handlers; 4638 while (inext) { 4639 PetscCall(PetscFree(inext->mtype)); 4640 iprev = inext; 4641 inext = inext->next; 4642 PetscCall(PetscFree(iprev)); 4643 } 4644 prev = next; 4645 next = next->next; 4646 PetscCall(PetscFree(prev)); 4647 } 4648 MatSolverTypeHolders = NULL; 4649 PetscFunctionReturn(PETSC_SUCCESS); 4650 } 4651 4652 /*@C 4653 MatFactorGetCanUseOrdering - Indicates if the factorization can use the ordering provided in `MatLUFactorSymbolic()`, `MatCholeskyFactorSymbolic()` 4654 4655 Logically Collective 4656 4657 Input Parameters: 4658 . mat - the matrix 4659 4660 Output Parameters: 4661 . flg - `PETSC_TRUE` if uses the ordering 4662 4663 Level: developer 4664 4665 Note: 4666 Most internal PETSc factorizations use the ordering passed to the factorization routine but external 4667 packages do not, thus we want to skip generating the ordering when it is not needed or used. 4668 4669 .seealso: `Mat`, [Matrix Factorization](sec_matfactor), `MatCopy()`, `MatDuplicate()`, `MatGetFactorAvailable()`, `MatGetFactor()`, `MatLUFactorSymbolic()`, `MatCholeskyFactorSymbolic()` 4670 @*/ 4671 PetscErrorCode MatFactorGetCanUseOrdering(Mat mat, PetscBool *flg) 4672 { 4673 PetscFunctionBegin; 4674 *flg = mat->canuseordering; 4675 PetscFunctionReturn(PETSC_SUCCESS); 4676 } 4677 4678 /*@C 4679 MatFactorGetPreferredOrdering - The preferred ordering for a particular matrix factor object 4680 4681 Logically Collective 4682 4683 Input Parameters: 4684 . mat - the matrix obtained with `MatGetFactor()` 4685 4686 Output Parameters: 4687 . otype - the preferred type 4688 4689 Level: developer 4690 4691 .seealso: `Mat`, [Matrix Factorization](sec_matfactor), `MatFactorType`, `MatOrderingType`, `MatCopy()`, `MatDuplicate()`, `MatGetFactorAvailable()`, `MatGetFactor()`, `MatLUFactorSymbolic()`, `MatCholeskyFactorSymbolic()` 4692 @*/ 4693 PetscErrorCode MatFactorGetPreferredOrdering(Mat mat, MatFactorType ftype, MatOrderingType *otype) 4694 { 4695 PetscFunctionBegin; 4696 *otype = mat->preferredordering[ftype]; 4697 PetscCheck(*otype, PETSC_COMM_SELF, PETSC_ERR_PLIB, "MatFactor did not have a preferred ordering"); 4698 PetscFunctionReturn(PETSC_SUCCESS); 4699 } 4700 4701 /*@C 4702 MatGetFactor - Returns a matrix suitable to calls to MatXXFactorSymbolic() 4703 4704 Collective 4705 4706 Input Parameters: 4707 + mat - the matrix 4708 . type - name of solver type, for example, superlu, petsc (to use PETSc's default) 4709 - ftype - factor type, `MAT_FACTOR_LU`, `MAT_FACTOR_CHOLESKY`, `MAT_FACTOR_ICC`, `MAT_FACTOR_ILU`, `MAT_FACTOR_QR` 4710 4711 Output Parameters: 4712 . f - the factor matrix used with MatXXFactorSymbolic() calls 4713 4714 Options Database Key: 4715 . -mat_factor_bind_factorization <host, device> - Where to do matrix factorization? Default is device (might consume more device memory. 4716 One can choose host to save device memory). Currently only supported with `MATSEQAIJCUSPARSE` matrices. 4717 4718 Level: intermediate 4719 4720 Notes: 4721 Users usually access the factorization solvers via `KSP` 4722 4723 Some PETSc matrix formats have alternative solvers available that are contained in alternative packages 4724 such as pastix, superlu, mumps etc. 4725 4726 PETSc must have been ./configure to use the external solver, using the option --download-package 4727 4728 Some of the packages have options for controlling the factorization, these are in the form -prefix_mat_packagename_packageoption 4729 where prefix is normally obtained from the calling `KSP`/`PC`. If `MatGetFactor()` is called directly one can set 4730 call `MatSetOptionsPrefixFactor()` on the originating matrix or `MatSetOptionsPrefix()` on the resulting factor matrix. 4731 4732 Developer Note: 4733 This should actually be called `MatCreateFactor()` since it creates a new factor object 4734 4735 .seealso: `Mat`, [Matrix Factorization](sec_matfactor), `KSP`, `MatSolverType`, `MatFactorType`, `MatCopy()`, `MatDuplicate()`, `MatGetFactorAvailable()`, `MatFactorGetCanUseOrdering()`, `MatSolverTypeRegister()`, 4736 `MAT_FACTOR_LU`, `MAT_FACTOR_CHOLESKY`, `MAT_FACTOR_ICC`, `MAT_FACTOR_ILU`, `MAT_FACTOR_QR` 4737 @*/ 4738 PetscErrorCode MatGetFactor(Mat mat, MatSolverType type, MatFactorType ftype, Mat *f) 4739 { 4740 PetscBool foundtype, foundmtype; 4741 PetscErrorCode (*conv)(Mat, MatFactorType, Mat *); 4742 4743 PetscFunctionBegin; 4744 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 4745 PetscValidType(mat, 1); 4746 4747 PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 4748 MatCheckPreallocated(mat, 1); 4749 4750 PetscCall(MatSolverTypeGet(type, ((PetscObject)mat)->type_name, ftype, &foundtype, &foundmtype, &conv)); 4751 if (!foundtype) { 4752 if (type) { 4753 SETERRQ(PetscObjectComm((PetscObject)mat), PETSC_ERR_MISSING_FACTOR, "Could not locate solver type %s for factorization type %s and matrix type %s. Perhaps you must ./configure with --download-%s", type, MatFactorTypes[ftype], 4754 ((PetscObject)mat)->type_name, type); 4755 } else { 4756 SETERRQ(PetscObjectComm((PetscObject)mat), PETSC_ERR_MISSING_FACTOR, "Could not locate a solver type for factorization type %s and matrix type %s.", MatFactorTypes[ftype], ((PetscObject)mat)->type_name); 4757 } 4758 } 4759 PetscCheck(foundmtype, PetscObjectComm((PetscObject)mat), PETSC_ERR_MISSING_FACTOR, "MatSolverType %s does not support matrix type %s", type, ((PetscObject)mat)->type_name); 4760 PetscCheck(conv, PetscObjectComm((PetscObject)mat), PETSC_ERR_MISSING_FACTOR, "MatSolverType %s does not support factorization type %s for matrix type %s", type, MatFactorTypes[ftype], ((PetscObject)mat)->type_name); 4761 4762 PetscCall((*conv)(mat, ftype, f)); 4763 if (mat->factorprefix) PetscCall(MatSetOptionsPrefix(*f, mat->factorprefix)); 4764 PetscFunctionReturn(PETSC_SUCCESS); 4765 } 4766 4767 /*@C 4768 MatGetFactorAvailable - Returns a a flag if matrix supports particular type and factor type 4769 4770 Not Collective 4771 4772 Input Parameters: 4773 + mat - the matrix 4774 . type - name of solver type, for example, superlu, petsc (to use PETSc's default) 4775 - ftype - factor type, `MAT_FACTOR_LU`, `MAT_FACTOR_CHOLESKY`, `MAT_FACTOR_ICC`, `MAT_FACTOR_ILU`, `MAT_FACTOR_QR` 4776 4777 Output Parameter: 4778 . flg - PETSC_TRUE if the factorization is available 4779 4780 Level: intermediate 4781 4782 Notes: 4783 Some PETSc matrix formats have alternative solvers available that are contained in alternative packages 4784 such as pastix, superlu, mumps etc. 4785 4786 PETSc must have been ./configure to use the external solver, using the option --download-package 4787 4788 Developer Note: 4789 This should actually be called MatCreateFactorAvailable() since MatGetFactor() creates a new factor object 4790 4791 .seealso: `Mat`, [Matrix Factorization](sec_matfactor), `MatSolverType`, `MatFactorType`, `MatGetFactor()`, `MatCopy()`, `MatDuplicate()`, `MatGetFactor()`, `MatSolverTypeRegister()`, 4792 `MAT_FACTOR_LU`, `MAT_FACTOR_CHOLESKY`, `MAT_FACTOR_ICC`, `MAT_FACTOR_ILU`, `MAT_FACTOR_QR` 4793 @*/ 4794 PetscErrorCode MatGetFactorAvailable(Mat mat, MatSolverType type, MatFactorType ftype, PetscBool *flg) 4795 { 4796 PetscErrorCode (*gconv)(Mat, MatFactorType, Mat *); 4797 4798 PetscFunctionBegin; 4799 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 4800 PetscValidType(mat, 1); 4801 PetscValidBoolPointer(flg, 4); 4802 4803 PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 4804 MatCheckPreallocated(mat, 1); 4805 4806 PetscCall(MatSolverTypeGet(type, ((PetscObject)mat)->type_name, ftype, NULL, NULL, &gconv)); 4807 *flg = gconv ? PETSC_TRUE : PETSC_FALSE; 4808 PetscFunctionReturn(PETSC_SUCCESS); 4809 } 4810 4811 /*@ 4812 MatDuplicate - Duplicates a matrix including the non-zero structure. 4813 4814 Collective 4815 4816 Input Parameters: 4817 + mat - the matrix 4818 - op - One of `MAT_DO_NOT_COPY_VALUES`, `MAT_COPY_VALUES`, or `MAT_SHARE_NONZERO_PATTERN`. 4819 See the manual page for `MatDuplicateOption()` for an explanation of these options. 4820 4821 Output Parameter: 4822 . M - pointer to place new matrix 4823 4824 Level: intermediate 4825 4826 Notes: 4827 You cannot change the nonzero pattern for the parent or child matrix if you use `MAT_SHARE_NONZERO_PATTERN`. 4828 4829 May be called with an unassembled input `Mat` if `MAT_DO_NOT_COPY_VALUES` is used, in which case the output `Mat` is unassembled as well. 4830 4831 When original mat is a product of matrix operation, e.g., an output of `MatMatMult()` or `MatCreateSubMatrix()`, only the simple matrix data structure of mat 4832 is duplicated and the internal data structures created for the reuse of previous matrix operations are not duplicated. 4833 User should not use `MatDuplicate()` to create new matrix M if M is intended to be reused as the product of matrix operation. 4834 4835 .seealso: `Mat`, `MatCopy()`, `MatConvert()`, `MatDuplicateOption` 4836 @*/ 4837 PetscErrorCode MatDuplicate(Mat mat, MatDuplicateOption op, Mat *M) 4838 { 4839 Mat B; 4840 VecType vtype; 4841 PetscInt i; 4842 PetscObject dm; 4843 void (*viewf)(void); 4844 4845 PetscFunctionBegin; 4846 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 4847 PetscValidType(mat, 1); 4848 PetscValidPointer(M, 3); 4849 PetscCheck(op != MAT_COPY_VALUES || mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "MAT_COPY_VALUES not allowed for unassembled matrix"); 4850 PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 4851 MatCheckPreallocated(mat, 1); 4852 4853 *M = NULL; 4854 PetscCall(PetscLogEventBegin(MAT_Convert, mat, 0, 0, 0)); 4855 PetscUseTypeMethod(mat, duplicate, op, M); 4856 PetscCall(PetscLogEventEnd(MAT_Convert, mat, 0, 0, 0)); 4857 B = *M; 4858 4859 PetscCall(MatGetOperation(mat, MATOP_VIEW, &viewf)); 4860 if (viewf) PetscCall(MatSetOperation(B, MATOP_VIEW, viewf)); 4861 PetscCall(MatGetVecType(mat, &vtype)); 4862 PetscCall(MatSetVecType(B, vtype)); 4863 4864 B->stencil.dim = mat->stencil.dim; 4865 B->stencil.noc = mat->stencil.noc; 4866 for (i = 0; i <= mat->stencil.dim; i++) { 4867 B->stencil.dims[i] = mat->stencil.dims[i]; 4868 B->stencil.starts[i] = mat->stencil.starts[i]; 4869 } 4870 4871 B->nooffproczerorows = mat->nooffproczerorows; 4872 B->nooffprocentries = mat->nooffprocentries; 4873 4874 PetscCall(PetscObjectQuery((PetscObject)mat, "__PETSc_dm", &dm)); 4875 if (dm) PetscCall(PetscObjectCompose((PetscObject)B, "__PETSc_dm", dm)); 4876 PetscCall(PetscObjectStateIncrease((PetscObject)B)); 4877 PetscFunctionReturn(PETSC_SUCCESS); 4878 } 4879 4880 /*@ 4881 MatGetDiagonal - Gets the diagonal of a matrix as a `Vec` 4882 4883 Logically Collective 4884 4885 Input Parameters: 4886 + mat - the matrix 4887 - v - the vector for storing the diagonal 4888 4889 Output Parameter: 4890 . v - the diagonal of the matrix 4891 4892 Level: intermediate 4893 4894 Note: 4895 Currently only correct in parallel for square matrices. 4896 4897 .seealso: `Mat`, `Vec`, `MatGetRow()`, `MatCreateSubMatrices()`, `MatCreateSubMatrix()`, `MatGetRowMaxAbs()` 4898 @*/ 4899 PetscErrorCode MatGetDiagonal(Mat mat, Vec v) 4900 { 4901 PetscFunctionBegin; 4902 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 4903 PetscValidType(mat, 1); 4904 PetscValidHeaderSpecific(v, VEC_CLASSID, 2); 4905 PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 4906 MatCheckPreallocated(mat, 1); 4907 4908 PetscUseTypeMethod(mat, getdiagonal, v); 4909 PetscCall(PetscObjectStateIncrease((PetscObject)v)); 4910 PetscFunctionReturn(PETSC_SUCCESS); 4911 } 4912 4913 /*@C 4914 MatGetRowMin - Gets the minimum value (of the real part) of each 4915 row of the matrix 4916 4917 Logically Collective 4918 4919 Input Parameter: 4920 . mat - the matrix 4921 4922 Output Parameters: 4923 + v - the vector for storing the maximums 4924 - idx - the indices of the column found for each row (optional) 4925 4926 Level: intermediate 4927 4928 Note: 4929 The result of this call are the same as if one converted the matrix to dense format 4930 and found the minimum value in each row (i.e. the implicit zeros are counted as zeros). 4931 4932 This code is only implemented for a couple of matrix formats. 4933 4934 .seealso: `Mat`, `MatGetDiagonal()`, `MatCreateSubMatrices()`, `MatCreateSubMatrix()`, `MatGetRowMaxAbs()`, `MatGetRowMinAbs()`, 4935 `MatGetRowMax()` 4936 @*/ 4937 PetscErrorCode MatGetRowMin(Mat mat, Vec v, PetscInt idx[]) 4938 { 4939 PetscFunctionBegin; 4940 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 4941 PetscValidType(mat, 1); 4942 PetscValidHeaderSpecific(v, VEC_CLASSID, 2); 4943 PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 4944 4945 if (!mat->cmap->N) { 4946 PetscCall(VecSet(v, PETSC_MAX_REAL)); 4947 if (idx) { 4948 PetscInt i, m = mat->rmap->n; 4949 for (i = 0; i < m; i++) idx[i] = -1; 4950 } 4951 } else { 4952 MatCheckPreallocated(mat, 1); 4953 } 4954 PetscUseTypeMethod(mat, getrowmin, v, idx); 4955 PetscCall(PetscObjectStateIncrease((PetscObject)v)); 4956 PetscFunctionReturn(PETSC_SUCCESS); 4957 } 4958 4959 /*@C 4960 MatGetRowMinAbs - Gets the minimum value (in absolute value) of each 4961 row of the matrix 4962 4963 Logically Collective 4964 4965 Input Parameter: 4966 . mat - the matrix 4967 4968 Output Parameters: 4969 + v - the vector for storing the minimums 4970 - idx - the indices of the column found for each row (or `NULL` if not needed) 4971 4972 Level: intermediate 4973 4974 Notes: 4975 if a row is completely empty or has only 0.0 values then the idx[] value for that 4976 row is 0 (the first column). 4977 4978 This code is only implemented for a couple of matrix formats. 4979 4980 .seealso: `Mat`, `MatGetDiagonal()`, `MatCreateSubMatrices()`, `MatCreateSubMatrix()`, `MatGetRowMax()`, `MatGetRowMaxAbs()`, `MatGetRowMin()` 4981 @*/ 4982 PetscErrorCode MatGetRowMinAbs(Mat mat, Vec v, PetscInt idx[]) 4983 { 4984 PetscFunctionBegin; 4985 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 4986 PetscValidType(mat, 1); 4987 PetscValidHeaderSpecific(v, VEC_CLASSID, 2); 4988 PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 4989 PetscCheck(!mat->factortype, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 4990 4991 if (!mat->cmap->N) { 4992 PetscCall(VecSet(v, 0.0)); 4993 if (idx) { 4994 PetscInt i, m = mat->rmap->n; 4995 for (i = 0; i < m; i++) idx[i] = -1; 4996 } 4997 } else { 4998 MatCheckPreallocated(mat, 1); 4999 if (idx) PetscCall(PetscArrayzero(idx, mat->rmap->n)); 5000 PetscUseTypeMethod(mat, getrowminabs, v, idx); 5001 } 5002 PetscCall(PetscObjectStateIncrease((PetscObject)v)); 5003 PetscFunctionReturn(PETSC_SUCCESS); 5004 } 5005 5006 /*@C 5007 MatGetRowMax - Gets the maximum value (of the real part) of each 5008 row of the matrix 5009 5010 Logically Collective 5011 5012 Input Parameter: 5013 . mat - the matrix 5014 5015 Output Parameters: 5016 + v - the vector for storing the maximums 5017 - idx - the indices of the column found for each row (optional) 5018 5019 Level: intermediate 5020 5021 Notes: 5022 The result of this call are the same as if one converted the matrix to dense format 5023 and found the minimum value in each row (i.e. the implicit zeros are counted as zeros). 5024 5025 This code is only implemented for a couple of matrix formats. 5026 5027 .seealso: `Mat`, `MatGetDiagonal()`, `MatCreateSubMatrices()`, `MatCreateSubMatrix()`, `MatGetRowMaxAbs()`, `MatGetRowMin()`, `MatGetRowMinAbs()` 5028 @*/ 5029 PetscErrorCode MatGetRowMax(Mat mat, Vec v, PetscInt idx[]) 5030 { 5031 PetscFunctionBegin; 5032 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 5033 PetscValidType(mat, 1); 5034 PetscValidHeaderSpecific(v, VEC_CLASSID, 2); 5035 PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 5036 5037 if (!mat->cmap->N) { 5038 PetscCall(VecSet(v, PETSC_MIN_REAL)); 5039 if (idx) { 5040 PetscInt i, m = mat->rmap->n; 5041 for (i = 0; i < m; i++) idx[i] = -1; 5042 } 5043 } else { 5044 MatCheckPreallocated(mat, 1); 5045 PetscUseTypeMethod(mat, getrowmax, v, idx); 5046 } 5047 PetscCall(PetscObjectStateIncrease((PetscObject)v)); 5048 PetscFunctionReturn(PETSC_SUCCESS); 5049 } 5050 5051 /*@C 5052 MatGetRowMaxAbs - Gets the maximum value (in absolute value) of each 5053 row of the matrix 5054 5055 Logically Collective 5056 5057 Input Parameter: 5058 . mat - the matrix 5059 5060 Output Parameters: 5061 + v - the vector for storing the maximums 5062 - idx - the indices of the column found for each row (or `NULL` if not needed) 5063 5064 Level: intermediate 5065 5066 Notes: 5067 if a row is completely empty or has only 0.0 values then the idx[] value for that 5068 row is 0 (the first column). 5069 5070 This code is only implemented for a couple of matrix formats. 5071 5072 .seealso: `Mat`, `MatGetDiagonal()`, `MatCreateSubMatrices()`, `MatCreateSubMatrix()`, `MatGetRowMax()`, `MatGetRowMin()`, `MatGetRowMinAbs()` 5073 @*/ 5074 PetscErrorCode MatGetRowMaxAbs(Mat mat, Vec v, PetscInt idx[]) 5075 { 5076 PetscFunctionBegin; 5077 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 5078 PetscValidType(mat, 1); 5079 PetscValidHeaderSpecific(v, VEC_CLASSID, 2); 5080 PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 5081 5082 if (!mat->cmap->N) { 5083 PetscCall(VecSet(v, 0.0)); 5084 if (idx) { 5085 PetscInt i, m = mat->rmap->n; 5086 for (i = 0; i < m; i++) idx[i] = -1; 5087 } 5088 } else { 5089 MatCheckPreallocated(mat, 1); 5090 if (idx) PetscCall(PetscArrayzero(idx, mat->rmap->n)); 5091 PetscUseTypeMethod(mat, getrowmaxabs, v, idx); 5092 } 5093 PetscCall(PetscObjectStateIncrease((PetscObject)v)); 5094 PetscFunctionReturn(PETSC_SUCCESS); 5095 } 5096 5097 /*@ 5098 MatGetRowSum - Gets the sum of each row of the matrix 5099 5100 Logically or Neighborhood Collective 5101 5102 Input Parameters: 5103 . mat - the matrix 5104 5105 Output Parameter: 5106 . v - the vector for storing the sum of rows 5107 5108 Level: intermediate 5109 5110 Notes: 5111 This code is slow since it is not currently specialized for different formats 5112 5113 .seealso: `Mat`, `MatGetDiagonal()`, `MatCreateSubMatrices()`, `MatCreateSubMatrix()`, `MatGetRowMax()`, `MatGetRowMin()`, `MatGetRowMaxAbs()`, `MatGetRowMinAbs()` 5114 @*/ 5115 PetscErrorCode MatGetRowSum(Mat mat, Vec v) 5116 { 5117 Vec ones; 5118 5119 PetscFunctionBegin; 5120 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 5121 PetscValidType(mat, 1); 5122 PetscValidHeaderSpecific(v, VEC_CLASSID, 2); 5123 PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 5124 MatCheckPreallocated(mat, 1); 5125 PetscCall(MatCreateVecs(mat, &ones, NULL)); 5126 PetscCall(VecSet(ones, 1.)); 5127 PetscCall(MatMult(mat, ones, v)); 5128 PetscCall(VecDestroy(&ones)); 5129 PetscFunctionReturn(PETSC_SUCCESS); 5130 } 5131 5132 /*@ 5133 MatTransposeSetPrecursor - Set the matrix from which the second matrix will receive numerical transpose data with a call to `MatTranspose`(A,`MAT_REUSE_MATRIX`,&B) 5134 when B was not obtained with `MatTranspose`(A,`MAT_INITIAL_MATRIX`,&B) 5135 5136 Collective 5137 5138 Input Parameter: 5139 . mat - the matrix to provide the transpose 5140 5141 Output Parameter: 5142 . mat - the matrix to contain the transpose; it MUST have the nonzero structure of the transpose of A or the code will crash or generate incorrect results 5143 5144 Level: advanced 5145 5146 Note: 5147 Normally he use of `MatTranspose`(A,`MAT_REUSE_MATRIX`,&B) requires that B was obtained with a call to `MatTranspose`(A,`MAT_INITIAL_MATRIX`,&B). This 5148 routine allows bypassing that call. 5149 5150 .seealso: `Mat`, `MatTransposeSymbolic()`, `MatTranspose()`, `MatMultTranspose()`, `MatMultTransposeAdd()`, `MatIsTranspose()`, `MatReuse`, `MAT_INITIAL_MATRIX`, `MAT_REUSE_MATRIX`, `MAT_INPLACE_MATRIX` 5151 @*/ 5152 PetscErrorCode MatTransposeSetPrecursor(Mat mat, Mat B) 5153 { 5154 PetscContainer rB = NULL; 5155 MatParentState *rb = NULL; 5156 5157 PetscFunctionBegin; 5158 PetscCall(PetscNew(&rb)); 5159 rb->id = ((PetscObject)mat)->id; 5160 rb->state = 0; 5161 PetscCall(MatGetNonzeroState(mat, &rb->nonzerostate)); 5162 PetscCall(PetscContainerCreate(PetscObjectComm((PetscObject)B), &rB)); 5163 PetscCall(PetscContainerSetPointer(rB, rb)); 5164 PetscCall(PetscContainerSetUserDestroy(rB, PetscContainerUserDestroyDefault)); 5165 PetscCall(PetscObjectCompose((PetscObject)B, "MatTransposeParent", (PetscObject)rB)); 5166 PetscCall(PetscObjectDereference((PetscObject)rB)); 5167 PetscFunctionReturn(PETSC_SUCCESS); 5168 } 5169 5170 /*@ 5171 MatTranspose - Computes an in-place or out-of-place transpose of a matrix. 5172 5173 Collective 5174 5175 Input Parameters: 5176 + mat - the matrix to transpose 5177 - reuse - either `MAT_INITIAL_MATRIX`, `MAT_REUSE_MATRIX`, or `MAT_INPLACE_MATRIX` 5178 5179 Output Parameter: 5180 . B - the transpose 5181 5182 Level: intermediate 5183 5184 Notes: 5185 If you use `MAT_INPLACE_MATRIX` then you must pass in &mat for B 5186 5187 `MAT_REUSE_MATRIX` uses the B matrix obtained from a previous call to this function with `MAT_INITIAL_MATRIX`. If you already have a matrix to contain the 5188 transpose, call `MatTransposeSetPrecursor`(mat,B) before calling this routine. 5189 5190 If the nonzero structure of mat changed from the previous call to this function with the same matrices an error will be generated for some matrix types. 5191 5192 Consider using `MatCreateTranspose()` instead if you only need a matrix that behaves like the transpose, but don't need the storage to be changed. 5193 5194 If mat is unchanged from the last call this function returns immediately without recomputing the result 5195 5196 If you only need the symbolic transpose, and not the numerical values, use `MatTransposeSymbolic()` 5197 5198 .seealso: `Mat`, `MatTransposeSetPrecursor()`, `MatMultTranspose()`, `MatMultTransposeAdd()`, `MatIsTranspose()`, `MatReuse`, `MAT_INITIAL_MATRIX`, `MAT_REUSE_MATRIX`, `MAT_INPLACE_MATRIX`, 5199 `MatTransposeSymbolic()` 5200 @*/ 5201 PetscErrorCode MatTranspose(Mat mat, MatReuse reuse, Mat *B) 5202 { 5203 PetscContainer rB = NULL; 5204 MatParentState *rb = NULL; 5205 5206 PetscFunctionBegin; 5207 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 5208 PetscValidType(mat, 1); 5209 PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 5210 PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 5211 PetscCheck(reuse != MAT_INPLACE_MATRIX || mat == *B, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "MAT_INPLACE_MATRIX requires last matrix to match first"); 5212 PetscCheck(reuse != MAT_REUSE_MATRIX || mat != *B, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Perhaps you mean MAT_INPLACE_MATRIX"); 5213 MatCheckPreallocated(mat, 1); 5214 if (reuse == MAT_REUSE_MATRIX) { 5215 PetscCall(PetscObjectQuery((PetscObject)*B, "MatTransposeParent", (PetscObject *)&rB)); 5216 PetscCheck(rB, PetscObjectComm((PetscObject)*B), PETSC_ERR_ARG_WRONG, "Reuse matrix used was not generated from call to MatTranspose(). Suggest MatTransposeSetPrecursor()."); 5217 PetscCall(PetscContainerGetPointer(rB, (void **)&rb)); 5218 PetscCheck(rb->id == ((PetscObject)mat)->id, PetscObjectComm((PetscObject)*B), PETSC_ERR_ARG_WRONG, "Reuse matrix used was not generated from input matrix"); 5219 if (rb->state == ((PetscObject)mat)->state) PetscFunctionReturn(PETSC_SUCCESS); 5220 } 5221 5222 PetscCall(PetscLogEventBegin(MAT_Transpose, mat, 0, 0, 0)); 5223 if (reuse != MAT_INPLACE_MATRIX || mat->symmetric != PETSC_BOOL3_TRUE) { 5224 PetscUseTypeMethod(mat, transpose, reuse, B); 5225 PetscCall(PetscObjectStateIncrease((PetscObject)*B)); 5226 } 5227 PetscCall(PetscLogEventEnd(MAT_Transpose, mat, 0, 0, 0)); 5228 5229 if (reuse == MAT_INITIAL_MATRIX) PetscCall(MatTransposeSetPrecursor(mat, *B)); 5230 if (reuse != MAT_INPLACE_MATRIX) { 5231 PetscCall(PetscObjectQuery((PetscObject)*B, "MatTransposeParent", (PetscObject *)&rB)); 5232 PetscCall(PetscContainerGetPointer(rB, (void **)&rb)); 5233 rb->state = ((PetscObject)mat)->state; 5234 rb->nonzerostate = mat->nonzerostate; 5235 } 5236 PetscFunctionReturn(PETSC_SUCCESS); 5237 } 5238 5239 /*@ 5240 MatTransposeSymbolic - Computes the symbolic part of the transpose of a matrix. 5241 5242 Collective 5243 5244 Input Parameters: 5245 . A - the matrix to transpose 5246 5247 Output Parameter: 5248 . B - the transpose. This is a complete matrix but the numerical portion is invalid. One can call `MatTranspose`(A,`MAT_REUSE_MATRIX`,&B) to compute the 5249 numerical portion. 5250 5251 Level: intermediate 5252 5253 Note: 5254 This is not supported for many matrix types, use `MatTranspose()` in those cases 5255 5256 .seealso: `Mat`, `MatTransposeSetPrecursor()`, `MatTranspose()`, `MatMultTranspose()`, `MatMultTransposeAdd()`, `MatIsTranspose()`, `MatReuse`, `MAT_INITIAL_MATRIX`, `MAT_REUSE_MATRIX`, `MAT_INPLACE_MATRIX` 5257 @*/ 5258 PetscErrorCode MatTransposeSymbolic(Mat A, Mat *B) 5259 { 5260 PetscFunctionBegin; 5261 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 5262 PetscValidType(A, 1); 5263 PetscCheck(A->assembled, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 5264 PetscCheck(!A->factortype, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 5265 PetscCheck(A->ops->transposesymbolic, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Mat type %s", ((PetscObject)A)->type_name); 5266 PetscCall(PetscLogEventBegin(MAT_Transpose, A, 0, 0, 0)); 5267 PetscCall((*A->ops->transposesymbolic)(A, B)); 5268 PetscCall(PetscLogEventEnd(MAT_Transpose, A, 0, 0, 0)); 5269 5270 PetscCall(MatTransposeSetPrecursor(A, *B)); 5271 PetscFunctionReturn(PETSC_SUCCESS); 5272 } 5273 5274 PetscErrorCode MatTransposeCheckNonzeroState_Private(Mat A, Mat B) 5275 { 5276 PetscContainer rB; 5277 MatParentState *rb; 5278 5279 PetscFunctionBegin; 5280 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 5281 PetscValidType(A, 1); 5282 PetscCheck(A->assembled, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 5283 PetscCheck(!A->factortype, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 5284 PetscCall(PetscObjectQuery((PetscObject)B, "MatTransposeParent", (PetscObject *)&rB)); 5285 PetscCheck(rB, PetscObjectComm((PetscObject)B), PETSC_ERR_ARG_WRONG, "Reuse matrix used was not generated from call to MatTranspose()"); 5286 PetscCall(PetscContainerGetPointer(rB, (void **)&rb)); 5287 PetscCheck(rb->id == ((PetscObject)A)->id, PetscObjectComm((PetscObject)B), PETSC_ERR_ARG_WRONG, "Reuse matrix used was not generated from input matrix"); 5288 PetscCheck(rb->nonzerostate == A->nonzerostate, PetscObjectComm((PetscObject)B), PETSC_ERR_ARG_WRONGSTATE, "Reuse matrix has changed nonzero structure"); 5289 PetscFunctionReturn(PETSC_SUCCESS); 5290 } 5291 5292 /*@ 5293 MatIsTranspose - Test whether a matrix is another one's transpose, 5294 or its own, in which case it tests symmetry. 5295 5296 Collective 5297 5298 Input Parameters: 5299 + A - the matrix to test 5300 - B - the matrix to test against, this can equal the first parameter 5301 5302 Output Parameters: 5303 . flg - the result 5304 5305 Level: intermediate 5306 5307 Notes: 5308 Only available for `MATAIJ` matrices. 5309 5310 The sequential algorithm has a running time of the order of the number of nonzeros; the parallel 5311 test involves parallel copies of the block-offdiagonal parts of the matrix. 5312 5313 .seealso: `Mat`, `MatTranspose()`, `MatIsSymmetric()`, `MatIsHermitian()` 5314 @*/ 5315 PetscErrorCode MatIsTranspose(Mat A, Mat B, PetscReal tol, PetscBool *flg) 5316 { 5317 PetscErrorCode (*f)(Mat, Mat, PetscReal, PetscBool *), (*g)(Mat, Mat, PetscReal, PetscBool *); 5318 5319 PetscFunctionBegin; 5320 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 5321 PetscValidHeaderSpecific(B, MAT_CLASSID, 2); 5322 PetscValidBoolPointer(flg, 4); 5323 PetscCall(PetscObjectQueryFunction((PetscObject)A, "MatIsTranspose_C", &f)); 5324 PetscCall(PetscObjectQueryFunction((PetscObject)B, "MatIsTranspose_C", &g)); 5325 *flg = PETSC_FALSE; 5326 if (f && g) { 5327 PetscCheck(f == g, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_NOTSAMETYPE, "Matrices do not have the same comparator for symmetry test"); 5328 PetscCall((*f)(A, B, tol, flg)); 5329 } else { 5330 MatType mattype; 5331 5332 PetscCall(MatGetType(f ? B : A, &mattype)); 5333 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Matrix of type %s does not support checking for transpose", mattype); 5334 } 5335 PetscFunctionReturn(PETSC_SUCCESS); 5336 } 5337 5338 /*@ 5339 MatHermitianTranspose - Computes an in-place or out-of-place Hermitian transpose of a matrix in complex conjugate. 5340 5341 Collective 5342 5343 Input Parameters: 5344 + mat - the matrix to transpose and complex conjugate 5345 - reuse - either `MAT_INITIAL_MATRIX`, `MAT_REUSE_MATRIX`, or `MAT_INPLACE_MATRIX` 5346 5347 Output Parameter: 5348 . B - the Hermitian transpose 5349 5350 Level: intermediate 5351 5352 .seealso: `Mat`, `MatTranspose()`, `MatMultTranspose()`, `MatMultTransposeAdd()`, `MatIsTranspose()`, `MatReuse` 5353 @*/ 5354 PetscErrorCode MatHermitianTranspose(Mat mat, MatReuse reuse, Mat *B) 5355 { 5356 PetscFunctionBegin; 5357 PetscCall(MatTranspose(mat, reuse, B)); 5358 #if defined(PETSC_USE_COMPLEX) 5359 PetscCall(MatConjugate(*B)); 5360 #endif 5361 PetscFunctionReturn(PETSC_SUCCESS); 5362 } 5363 5364 /*@ 5365 MatIsHermitianTranspose - Test whether a matrix is another one's Hermitian transpose, 5366 5367 Collective 5368 5369 Input Parameters: 5370 + A - the matrix to test 5371 - B - the matrix to test against, this can equal the first parameter 5372 5373 Output Parameters: 5374 . flg - the result 5375 5376 Level: intermediate 5377 5378 Notes: 5379 Only available for `MATAIJ` matrices. 5380 5381 The sequential algorithm 5382 has a running time of the order of the number of nonzeros; the parallel 5383 test involves parallel copies of the block-offdiagonal parts of the matrix. 5384 5385 .seealso: `Mat`, `MatTranspose()`, `MatIsSymmetric()`, `MatIsHermitian()`, `MatIsTranspose()` 5386 @*/ 5387 PetscErrorCode MatIsHermitianTranspose(Mat A, Mat B, PetscReal tol, PetscBool *flg) 5388 { 5389 PetscErrorCode (*f)(Mat, Mat, PetscReal, PetscBool *), (*g)(Mat, Mat, PetscReal, PetscBool *); 5390 5391 PetscFunctionBegin; 5392 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 5393 PetscValidHeaderSpecific(B, MAT_CLASSID, 2); 5394 PetscValidBoolPointer(flg, 4); 5395 PetscCall(PetscObjectQueryFunction((PetscObject)A, "MatIsHermitianTranspose_C", &f)); 5396 PetscCall(PetscObjectQueryFunction((PetscObject)B, "MatIsHermitianTranspose_C", &g)); 5397 if (f && g) { 5398 PetscCheck(f != g, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_NOTSAMETYPE, "Matrices do not have the same comparator for Hermitian test"); 5399 PetscCall((*f)(A, B, tol, flg)); 5400 } 5401 PetscFunctionReturn(PETSC_SUCCESS); 5402 } 5403 5404 /*@ 5405 MatPermute - Creates a new matrix with rows and columns permuted from the 5406 original. 5407 5408 Collective 5409 5410 Input Parameters: 5411 + mat - the matrix to permute 5412 . row - row permutation, each processor supplies only the permutation for its rows 5413 - col - column permutation, each processor supplies only the permutation for its columns 5414 5415 Output Parameters: 5416 . B - the permuted matrix 5417 5418 Level: advanced 5419 5420 Note: 5421 The index sets map from row/col of permuted matrix to row/col of original matrix. 5422 The index sets should be on the same communicator as mat and have the same local sizes. 5423 5424 Developer Note: 5425 If you want to implement `MatPermute()` for a matrix type, and your approach doesn't 5426 exploit the fact that row and col are permutations, consider implementing the 5427 more general `MatCreateSubMatrix()` instead. 5428 5429 .seealso: `Mat`, `MatGetOrdering()`, `ISAllGather()`, `MatCreateSubMatrix()` 5430 @*/ 5431 PetscErrorCode MatPermute(Mat mat, IS row, IS col, Mat *B) 5432 { 5433 PetscFunctionBegin; 5434 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 5435 PetscValidType(mat, 1); 5436 PetscValidHeaderSpecific(row, IS_CLASSID, 2); 5437 PetscValidHeaderSpecific(col, IS_CLASSID, 3); 5438 PetscValidPointer(B, 4); 5439 PetscCheckSameComm(mat, 1, row, 2); 5440 if (row != col) PetscCheckSameComm(row, 2, col, 3); 5441 PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 5442 PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 5443 PetscCheck(mat->ops->permute || mat->ops->createsubmatrix, PETSC_COMM_SELF, PETSC_ERR_SUP, "MatPermute not available for Mat type %s", ((PetscObject)mat)->type_name); 5444 MatCheckPreallocated(mat, 1); 5445 5446 if (mat->ops->permute) { 5447 PetscUseTypeMethod(mat, permute, row, col, B); 5448 PetscCall(PetscObjectStateIncrease((PetscObject)*B)); 5449 } else { 5450 PetscCall(MatCreateSubMatrix(mat, row, col, MAT_INITIAL_MATRIX, B)); 5451 } 5452 PetscFunctionReturn(PETSC_SUCCESS); 5453 } 5454 5455 /*@ 5456 MatEqual - Compares two matrices. 5457 5458 Collective 5459 5460 Input Parameters: 5461 + A - the first matrix 5462 - B - the second matrix 5463 5464 Output Parameter: 5465 . flg - `PETSC_TRUE` if the matrices are equal; `PETSC_FALSE` otherwise. 5466 5467 Level: intermediate 5468 5469 .seealso: `Mat` 5470 @*/ 5471 PetscErrorCode MatEqual(Mat A, Mat B, PetscBool *flg) 5472 { 5473 PetscFunctionBegin; 5474 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 5475 PetscValidHeaderSpecific(B, MAT_CLASSID, 2); 5476 PetscValidType(A, 1); 5477 PetscValidType(B, 2); 5478 PetscValidBoolPointer(flg, 3); 5479 PetscCheckSameComm(A, 1, B, 2); 5480 MatCheckPreallocated(A, 1); 5481 MatCheckPreallocated(B, 2); 5482 PetscCheck(A->assembled, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 5483 PetscCheck(B->assembled, PetscObjectComm((PetscObject)B), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 5484 PetscCheck(A->rmap->N == B->rmap->N && A->cmap->N == B->cmap->N, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_SIZ, "Mat A,Mat B: global dim %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT, A->rmap->N, B->rmap->N, A->cmap->N, 5485 B->cmap->N); 5486 if (A->ops->equal && A->ops->equal == B->ops->equal) { 5487 PetscUseTypeMethod(A, equal, B, flg); 5488 } else { 5489 PetscCall(MatMultEqual(A, B, 10, flg)); 5490 } 5491 PetscFunctionReturn(PETSC_SUCCESS); 5492 } 5493 5494 /*@ 5495 MatDiagonalScale - Scales a matrix on the left and right by diagonal 5496 matrices that are stored as vectors. Either of the two scaling 5497 matrices can be `NULL`. 5498 5499 Collective 5500 5501 Input Parameters: 5502 + mat - the matrix to be scaled 5503 . l - the left scaling vector (or `NULL`) 5504 - r - the right scaling vector (or `NULL`) 5505 5506 Level: intermediate 5507 5508 Note: 5509 `MatDiagonalScale()` computes A = LAR, where 5510 L = a diagonal matrix (stored as a vector), R = a diagonal matrix (stored as a vector) 5511 The L scales the rows of the matrix, the R scales the columns of the matrix. 5512 5513 .seealso: `Mat`, `MatScale()`, `MatShift()`, `MatDiagonalSet()` 5514 @*/ 5515 PetscErrorCode MatDiagonalScale(Mat mat, Vec l, Vec r) 5516 { 5517 PetscFunctionBegin; 5518 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 5519 PetscValidType(mat, 1); 5520 if (l) { 5521 PetscValidHeaderSpecific(l, VEC_CLASSID, 2); 5522 PetscCheckSameComm(mat, 1, l, 2); 5523 } 5524 if (r) { 5525 PetscValidHeaderSpecific(r, VEC_CLASSID, 3); 5526 PetscCheckSameComm(mat, 1, r, 3); 5527 } 5528 PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 5529 PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 5530 MatCheckPreallocated(mat, 1); 5531 if (!l && !r) PetscFunctionReturn(PETSC_SUCCESS); 5532 5533 PetscCall(PetscLogEventBegin(MAT_Scale, mat, 0, 0, 0)); 5534 PetscUseTypeMethod(mat, diagonalscale, l, r); 5535 PetscCall(PetscLogEventEnd(MAT_Scale, mat, 0, 0, 0)); 5536 PetscCall(PetscObjectStateIncrease((PetscObject)mat)); 5537 if (l != r) mat->symmetric = PETSC_BOOL3_FALSE; 5538 PetscFunctionReturn(PETSC_SUCCESS); 5539 } 5540 5541 /*@ 5542 MatScale - Scales all elements of a matrix by a given number. 5543 5544 Logically Collective 5545 5546 Input Parameters: 5547 + mat - the matrix to be scaled 5548 - a - the scaling value 5549 5550 Output Parameter: 5551 . mat - the scaled matrix 5552 5553 Level: intermediate 5554 5555 .seealso: `Mat`, `MatDiagonalScale()` 5556 @*/ 5557 PetscErrorCode MatScale(Mat mat, PetscScalar a) 5558 { 5559 PetscFunctionBegin; 5560 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 5561 PetscValidType(mat, 1); 5562 PetscCheck(a == (PetscScalar)1.0 || mat->ops->scale, PETSC_COMM_SELF, PETSC_ERR_SUP, "Mat type %s", ((PetscObject)mat)->type_name); 5563 PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 5564 PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 5565 PetscValidLogicalCollectiveScalar(mat, a, 2); 5566 MatCheckPreallocated(mat, 1); 5567 5568 PetscCall(PetscLogEventBegin(MAT_Scale, mat, 0, 0, 0)); 5569 if (a != (PetscScalar)1.0) { 5570 PetscUseTypeMethod(mat, scale, a); 5571 PetscCall(PetscObjectStateIncrease((PetscObject)mat)); 5572 } 5573 PetscCall(PetscLogEventEnd(MAT_Scale, mat, 0, 0, 0)); 5574 PetscFunctionReturn(PETSC_SUCCESS); 5575 } 5576 5577 /*@ 5578 MatNorm - Calculates various norms of a matrix. 5579 5580 Collective 5581 5582 Input Parameters: 5583 + mat - the matrix 5584 - type - the type of norm, `NORM_1`, `NORM_FROBENIUS`, `NORM_INFINITY` 5585 5586 Output Parameter: 5587 . nrm - the resulting norm 5588 5589 Level: intermediate 5590 5591 .seealso: `Mat` 5592 @*/ 5593 PetscErrorCode MatNorm(Mat mat, NormType type, PetscReal *nrm) 5594 { 5595 PetscFunctionBegin; 5596 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 5597 PetscValidType(mat, 1); 5598 PetscValidRealPointer(nrm, 3); 5599 5600 PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 5601 PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 5602 MatCheckPreallocated(mat, 1); 5603 5604 PetscUseTypeMethod(mat, norm, type, nrm); 5605 PetscFunctionReturn(PETSC_SUCCESS); 5606 } 5607 5608 /* 5609 This variable is used to prevent counting of MatAssemblyBegin() that 5610 are called from within a MatAssemblyEnd(). 5611 */ 5612 static PetscInt MatAssemblyEnd_InUse = 0; 5613 /*@ 5614 MatAssemblyBegin - Begins assembling the matrix. This routine should 5615 be called after completing all calls to `MatSetValues()`. 5616 5617 Collective 5618 5619 Input Parameters: 5620 + mat - the matrix 5621 - type - type of assembly, either `MAT_FLUSH_ASSEMBLY` or `MAT_FINAL_ASSEMBLY` 5622 5623 Level: beginner 5624 5625 Notes: 5626 `MatSetValues()` generally caches the values that belong to other MPI ranks. The matrix is ready to 5627 use only after `MatAssemblyBegin()` and `MatAssemblyEnd()` have been called. 5628 5629 Use `MAT_FLUSH_ASSEMBLY` when switching between `ADD_VALUES` and `INSERT_VALUES` 5630 in `MatSetValues()`; use `MAT_FINAL_ASSEMBLY` for the final assembly before 5631 using the matrix. 5632 5633 ALL processes that share a matrix MUST call `MatAssemblyBegin()` and `MatAssemblyEnd()` the SAME NUMBER of times, and each time with the 5634 same flag of `MAT_FLUSH_ASSEMBLY` or `MAT_FINAL_ASSEMBLY` for all processes. Thus you CANNOT locally change from `ADD_VALUES` to `INSERT_VALUES`, that is 5635 a global collective operation requiring all processes that share the matrix. 5636 5637 Space for preallocated nonzeros that is not filled by a call to `MatSetValues()` or a related routine are compressed 5638 out by assembly. If you intend to use that extra space on a subsequent assembly, be sure to insert explicit zeros 5639 before `MAT_FINAL_ASSEMBLY` so the space is not compressed out. 5640 5641 .seealso: `Mat`, `MatAssemblyEnd()`, `MatSetValues()`, `MatAssembled()` 5642 @*/ 5643 PetscErrorCode MatAssemblyBegin(Mat mat, MatAssemblyType type) 5644 { 5645 PetscFunctionBegin; 5646 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 5647 PetscValidType(mat, 1); 5648 MatCheckPreallocated(mat, 1); 5649 PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix.\nDid you forget to call MatSetUnfactored()?"); 5650 if (mat->assembled) { 5651 mat->was_assembled = PETSC_TRUE; 5652 mat->assembled = PETSC_FALSE; 5653 } 5654 5655 if (!MatAssemblyEnd_InUse) { 5656 PetscCall(PetscLogEventBegin(MAT_AssemblyBegin, mat, 0, 0, 0)); 5657 PetscTryTypeMethod(mat, assemblybegin, type); 5658 PetscCall(PetscLogEventEnd(MAT_AssemblyBegin, mat, 0, 0, 0)); 5659 } else PetscTryTypeMethod(mat, assemblybegin, type); 5660 PetscFunctionReturn(PETSC_SUCCESS); 5661 } 5662 5663 /*@ 5664 MatAssembled - Indicates if a matrix has been assembled and is ready for 5665 use; for example, in matrix-vector product. 5666 5667 Not Collective 5668 5669 Input Parameter: 5670 . mat - the matrix 5671 5672 Output Parameter: 5673 . assembled - `PETSC_TRUE` or `PETSC_FALSE` 5674 5675 Level: advanced 5676 5677 .seealso: `Mat`, `MatAssemblyEnd()`, `MatSetValues()`, `MatAssemblyBegin()` 5678 @*/ 5679 PetscErrorCode MatAssembled(Mat mat, PetscBool *assembled) 5680 { 5681 PetscFunctionBegin; 5682 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 5683 PetscValidBoolPointer(assembled, 2); 5684 *assembled = mat->assembled; 5685 PetscFunctionReturn(PETSC_SUCCESS); 5686 } 5687 5688 /*@ 5689 MatAssemblyEnd - Completes assembling the matrix. This routine should 5690 be called after `MatAssemblyBegin()`. 5691 5692 Collective on Mat 5693 5694 Input Parameters: 5695 + mat - the matrix 5696 - type - type of assembly, either `MAT_FLUSH_ASSEMBLY` or `MAT_FINAL_ASSEMBLY` 5697 5698 Options Database Keys: 5699 + -mat_view ::ascii_info - Prints info on matrix at conclusion of `MatEndAssembly()` 5700 . -mat_view ::ascii_info_detail - Prints more detailed info 5701 . -mat_view - Prints matrix in ASCII format 5702 . -mat_view ::ascii_matlab - Prints matrix in Matlab format 5703 . -mat_view draw - draws nonzero structure of matrix, using `MatView()` and `PetscDrawOpenX()`. 5704 . -display <name> - Sets display name (default is host) 5705 . -draw_pause <sec> - Sets number of seconds to pause after display 5706 . -mat_view socket - Sends matrix to socket, can be accessed from Matlab (See [Using MATLAB with PETSc](ch_matlab)) 5707 . -viewer_socket_machine <machine> - Machine to use for socket 5708 . -viewer_socket_port <port> - Port number to use for socket 5709 - -mat_view binary:filename[:append] - Save matrix to file in binary format 5710 5711 Level: beginner 5712 5713 .seealso: `Mat`, `MatAssemblyBegin()`, `MatSetValues()`, `PetscDrawOpenX()`, `PetscDrawCreate()`, `MatView()`, `MatAssembled()`, `PetscViewerSocketOpen()` 5714 @*/ 5715 PetscErrorCode MatAssemblyEnd(Mat mat, MatAssemblyType type) 5716 { 5717 static PetscInt inassm = 0; 5718 PetscBool flg = PETSC_FALSE; 5719 5720 PetscFunctionBegin; 5721 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 5722 PetscValidType(mat, 1); 5723 5724 inassm++; 5725 MatAssemblyEnd_InUse++; 5726 if (MatAssemblyEnd_InUse == 1) { /* Do the logging only the first time through */ 5727 PetscCall(PetscLogEventBegin(MAT_AssemblyEnd, mat, 0, 0, 0)); 5728 PetscTryTypeMethod(mat, assemblyend, type); 5729 PetscCall(PetscLogEventEnd(MAT_AssemblyEnd, mat, 0, 0, 0)); 5730 } else PetscTryTypeMethod(mat, assemblyend, type); 5731 5732 /* Flush assembly is not a true assembly */ 5733 if (type != MAT_FLUSH_ASSEMBLY) { 5734 if (mat->num_ass) { 5735 if (!mat->symmetry_eternal) { 5736 mat->symmetric = PETSC_BOOL3_UNKNOWN; 5737 mat->hermitian = PETSC_BOOL3_UNKNOWN; 5738 } 5739 if (!mat->structural_symmetry_eternal && mat->ass_nonzerostate != mat->nonzerostate) mat->structurally_symmetric = PETSC_BOOL3_UNKNOWN; 5740 if (!mat->spd_eternal) mat->spd = PETSC_BOOL3_UNKNOWN; 5741 } 5742 mat->num_ass++; 5743 mat->assembled = PETSC_TRUE; 5744 mat->ass_nonzerostate = mat->nonzerostate; 5745 } 5746 5747 mat->insertmode = NOT_SET_VALUES; 5748 MatAssemblyEnd_InUse--; 5749 PetscCall(PetscObjectStateIncrease((PetscObject)mat)); 5750 if (inassm == 1 && type != MAT_FLUSH_ASSEMBLY) { 5751 PetscCall(MatViewFromOptions(mat, NULL, "-mat_view")); 5752 5753 if (mat->checksymmetryonassembly) { 5754 PetscCall(MatIsSymmetric(mat, mat->checksymmetrytol, &flg)); 5755 if (flg) { 5756 PetscCall(PetscPrintf(PetscObjectComm((PetscObject)mat), "Matrix is symmetric (tolerance %g)\n", (double)mat->checksymmetrytol)); 5757 } else { 5758 PetscCall(PetscPrintf(PetscObjectComm((PetscObject)mat), "Matrix is not symmetric (tolerance %g)\n", (double)mat->checksymmetrytol)); 5759 } 5760 } 5761 if (mat->nullsp && mat->checknullspaceonassembly) PetscCall(MatNullSpaceTest(mat->nullsp, mat, NULL)); 5762 } 5763 inassm--; 5764 PetscFunctionReturn(PETSC_SUCCESS); 5765 } 5766 5767 /*@ 5768 MatSetOption - Sets a parameter option for a matrix. Some options 5769 may be specific to certain storage formats. Some options 5770 determine how values will be inserted (or added). Sorted, 5771 row-oriented input will generally assemble the fastest. The default 5772 is row-oriented. 5773 5774 Logically Collective for certain operations, such as `MAT_SPD`, not collective for `MAT_ROW_ORIENTED`, see `MatOption` 5775 5776 Input Parameters: 5777 + mat - the matrix 5778 . option - the option, one of those listed below (and possibly others), 5779 - flg - turn the option on (`PETSC_TRUE`) or off (`PETSC_FALSE`) 5780 5781 Options Describing Matrix Structure: 5782 + `MAT_SPD` - symmetric positive definite 5783 . `MAT_SYMMETRIC` - symmetric in terms of both structure and value 5784 . `MAT_HERMITIAN` - transpose is the complex conjugation 5785 . `MAT_STRUCTURALLY_SYMMETRIC` - symmetric nonzero structure 5786 . `MAT_SYMMETRY_ETERNAL` - indicates the symmetry (or Hermitian structure) or its absence will persist through any changes to the matrix 5787 . `MAT_STRUCTURAL_SYMMETRY_ETERNAL` - indicates the structural symmetry or its absence will persist through any changes to the matrix 5788 - `MAT_SPD_ETERNAL` - indicates the value of `MAT_SPD` (true or false) will persist through any changes to the matrix 5789 5790 These are not really options of the matrix, they are knowledge about the structure of the matrix that users may provide so that they 5791 do not need to be computed (usually at a high cost) 5792 5793 Options For Use with `MatSetValues()`: 5794 Insert a logically dense subblock, which can be 5795 . `MAT_ROW_ORIENTED` - row-oriented (default) 5796 5797 Note these options reflect the data you pass in with `MatSetValues()`; it has 5798 nothing to do with how the data is stored internally in the matrix 5799 data structure. 5800 5801 When (re)assembling a matrix, we can restrict the input for 5802 efficiency/debugging purposes. These options include 5803 + `MAT_NEW_NONZERO_LOCATIONS` - additional insertions will be allowed if they generate a new nonzero (slow) 5804 . `MAT_FORCE_DIAGONAL_ENTRIES` - forces diagonal entries to be allocated 5805 . `MAT_IGNORE_OFF_PROC_ENTRIES` - drops off-processor entries 5806 . `MAT_NEW_NONZERO_LOCATION_ERR` - generates an error for new matrix entry 5807 . `MAT_USE_HASH_TABLE` - uses a hash table to speed up matrix assembly 5808 . `MAT_NO_OFF_PROC_ENTRIES` - you know each process will only set values for its own rows, will generate an error if 5809 any process sets values for another process. This avoids all reductions in the MatAssembly routines and thus improves 5810 performance for very large process counts. 5811 - `MAT_SUBSET_OFF_PROC_ENTRIES` - you know that the first assembly after setting this flag will set a superset 5812 of the off-process entries required for all subsequent assemblies. This avoids a rendezvous step in the MatAssembly 5813 functions, instead sending only neighbor messages. 5814 5815 Level: intermediate 5816 5817 Notes: 5818 Except for `MAT_UNUSED_NONZERO_LOCATION_ERR` and `MAT_ROW_ORIENTED` all processes that share the matrix must pass the same value in flg! 5819 5820 Some options are relevant only for particular matrix types and 5821 are thus ignored by others. Other options are not supported by 5822 certain matrix types and will generate an error message if set. 5823 5824 If using Fortran to compute a matrix, one may need to 5825 use the column-oriented option (or convert to the row-oriented 5826 format). 5827 5828 `MAT_NEW_NONZERO_LOCATIONS` set to `PETSC_FALSE` indicates that any add or insertion 5829 that would generate a new entry in the nonzero structure is instead 5830 ignored. Thus, if memory has not already been allocated for this particular 5831 data, then the insertion is ignored. For dense matrices, in which 5832 the entire array is allocated, no entries are ever ignored. 5833 Set after the first `MatAssemblyEnd()`. If this option is set then the MatAssemblyBegin/End() processes has one less global reduction 5834 5835 `MAT_NEW_NONZERO_LOCATION_ERR` set to PETSC_TRUE indicates that any add or insertion 5836 that would generate a new entry in the nonzero structure instead produces 5837 an error. (Currently supported for `MATAIJ` and `MATBAIJ` formats only.) If this option is set then the `MatAssemblyBegin()`/`MatAssemblyEnd()` processes has one less global reduction 5838 5839 `MAT_NEW_NONZERO_ALLOCATION_ERR` set to `PETSC_TRUE` indicates that any add or insertion 5840 that would generate a new entry that has not been preallocated will 5841 instead produce an error. (Currently supported for `MATAIJ` and `MATBAIJ` formats 5842 only.) This is a useful flag when debugging matrix memory preallocation. 5843 If this option is set then the `MatAssemblyBegin()`/`MatAssemblyEnd()` processes has one less global reduction 5844 5845 `MAT_IGNORE_OFF_PROC_ENTRIES` set to `PETSC_TRUE` indicates entries destined for 5846 other processors should be dropped, rather than stashed. 5847 This is useful if you know that the "owning" processor is also 5848 always generating the correct matrix entries, so that PETSc need 5849 not transfer duplicate entries generated on another processor. 5850 5851 `MAT_USE_HASH_TABLE` indicates that a hash table be used to improve the 5852 searches during matrix assembly. When this flag is set, the hash table 5853 is created during the first matrix assembly. This hash table is 5854 used the next time through, during `MatSetValues()`/`MatSetValuesBlocked()` 5855 to improve the searching of indices. `MAT_NEW_NONZERO_LOCATIONS` flag 5856 should be used with `MAT_USE_HASH_TABLE` flag. This option is currently 5857 supported by `MATMPIBAIJ` format only. 5858 5859 `MAT_KEEP_NONZERO_PATTERN` indicates when `MatZeroRows()` is called the zeroed entries 5860 are kept in the nonzero structure 5861 5862 `MAT_IGNORE_ZERO_ENTRIES` - for `MATAIJ` and `MATIS` matrices this will stop zero values from creating 5863 a zero location in the matrix 5864 5865 `MAT_USE_INODES` - indicates using inode version of the code - works with `MATAIJ` matrix types 5866 5867 `MAT_NO_OFF_PROC_ZERO_ROWS` - you know each process will only zero its own rows. This avoids all reductions in the 5868 zero row routines and thus improves performance for very large process counts. 5869 5870 `MAT_IGNORE_LOWER_TRIANGULAR` - For `MATSBAIJ` matrices will ignore any insertions you make in the lower triangular 5871 part of the matrix (since they should match the upper triangular part). 5872 5873 `MAT_SORTED_FULL` - each process provides exactly its local rows; all column indices for a given row are passed in a 5874 single call to `MatSetValues()`, preallocation is perfect, row oriented, `INSERT_VALUES` is used. Common 5875 with finite difference schemes with non-periodic boundary conditions. 5876 5877 Developer Note: 5878 `MAT_SYMMETRY_ETERNAL`, `MAT_STRUCTURAL_SYMMETRY_ETERNAL`, and `MAT_SPD_ETERNAL` are used by `MatAssemblyEnd()` and in other 5879 places where otherwise the value of `MAT_SYMMETRIC`, `MAT_STRUCTURAL_SYMMETRIC` or `MAT_SPD` would need to be changed back 5880 to `PETSC_BOOL3_UNKNOWN` because the matrix values had changed so the code cannot be certain that the related property had 5881 not changed. 5882 5883 .seealso: `MatOption`, `Mat`, `MatGetOption()` 5884 @*/ 5885 PetscErrorCode MatSetOption(Mat mat, MatOption op, PetscBool flg) 5886 { 5887 PetscFunctionBegin; 5888 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 5889 if (op > 0) { 5890 PetscValidLogicalCollectiveEnum(mat, op, 2); 5891 PetscValidLogicalCollectiveBool(mat, flg, 3); 5892 } 5893 5894 PetscCheck(((int)op) > MAT_OPTION_MIN && ((int)op) < MAT_OPTION_MAX, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_OUTOFRANGE, "Options %d is out of range", (int)op); 5895 5896 switch (op) { 5897 case MAT_FORCE_DIAGONAL_ENTRIES: 5898 mat->force_diagonals = flg; 5899 PetscFunctionReturn(PETSC_SUCCESS); 5900 case MAT_NO_OFF_PROC_ENTRIES: 5901 mat->nooffprocentries = flg; 5902 PetscFunctionReturn(PETSC_SUCCESS); 5903 case MAT_SUBSET_OFF_PROC_ENTRIES: 5904 mat->assembly_subset = flg; 5905 if (!mat->assembly_subset) { /* See the same logic in VecAssembly wrt VEC_SUBSET_OFF_PROC_ENTRIES */ 5906 #if !defined(PETSC_HAVE_MPIUNI) 5907 PetscCall(MatStashScatterDestroy_BTS(&mat->stash)); 5908 #endif 5909 mat->stash.first_assembly_done = PETSC_FALSE; 5910 } 5911 PetscFunctionReturn(PETSC_SUCCESS); 5912 case MAT_NO_OFF_PROC_ZERO_ROWS: 5913 mat->nooffproczerorows = flg; 5914 PetscFunctionReturn(PETSC_SUCCESS); 5915 case MAT_SPD: 5916 if (flg) { 5917 mat->spd = PETSC_BOOL3_TRUE; 5918 mat->symmetric = PETSC_BOOL3_TRUE; 5919 mat->structurally_symmetric = PETSC_BOOL3_TRUE; 5920 } else { 5921 mat->spd = PETSC_BOOL3_FALSE; 5922 } 5923 break; 5924 case MAT_SYMMETRIC: 5925 mat->symmetric = flg ? PETSC_BOOL3_TRUE : PETSC_BOOL3_FALSE; 5926 if (flg) mat->structurally_symmetric = PETSC_BOOL3_TRUE; 5927 #if !defined(PETSC_USE_COMPLEX) 5928 mat->hermitian = flg ? PETSC_BOOL3_TRUE : PETSC_BOOL3_FALSE; 5929 #endif 5930 break; 5931 case MAT_HERMITIAN: 5932 mat->hermitian = flg ? PETSC_BOOL3_TRUE : PETSC_BOOL3_FALSE; 5933 if (flg) mat->structurally_symmetric = PETSC_BOOL3_TRUE; 5934 #if !defined(PETSC_USE_COMPLEX) 5935 mat->symmetric = flg ? PETSC_BOOL3_TRUE : PETSC_BOOL3_FALSE; 5936 #endif 5937 break; 5938 case MAT_STRUCTURALLY_SYMMETRIC: 5939 mat->structurally_symmetric = flg ? PETSC_BOOL3_TRUE : PETSC_BOOL3_FALSE; 5940 break; 5941 case MAT_SYMMETRY_ETERNAL: 5942 PetscCheck(mat->symmetric != PETSC_BOOL3_UNKNOWN, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Cannot set MAT_SYMMETRY_ETERNAL without first setting MAT_SYMMETRIC to true or false"); 5943 mat->symmetry_eternal = flg; 5944 if (flg) mat->structural_symmetry_eternal = PETSC_TRUE; 5945 break; 5946 case MAT_STRUCTURAL_SYMMETRY_ETERNAL: 5947 PetscCheck(mat->structurally_symmetric != PETSC_BOOL3_UNKNOWN, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Cannot set MAT_STRUCTURAL_SYMMETRY_ETERNAL without first setting MAT_STRUCTURAL_SYMMETRIC to true or false"); 5948 mat->structural_symmetry_eternal = flg; 5949 break; 5950 case MAT_SPD_ETERNAL: 5951 PetscCheck(mat->spd != PETSC_BOOL3_UNKNOWN, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Cannot set MAT_SPD_ETERNAL without first setting MAT_SPD to true or false"); 5952 mat->spd_eternal = flg; 5953 if (flg) { 5954 mat->structural_symmetry_eternal = PETSC_TRUE; 5955 mat->symmetry_eternal = PETSC_TRUE; 5956 } 5957 break; 5958 case MAT_STRUCTURE_ONLY: 5959 mat->structure_only = flg; 5960 break; 5961 case MAT_SORTED_FULL: 5962 mat->sortedfull = flg; 5963 break; 5964 default: 5965 break; 5966 } 5967 PetscTryTypeMethod(mat, setoption, op, flg); 5968 PetscFunctionReturn(PETSC_SUCCESS); 5969 } 5970 5971 /*@ 5972 MatGetOption - Gets a parameter option that has been set for a matrix. 5973 5974 Logically Collective 5975 5976 Input Parameters: 5977 + mat - the matrix 5978 - option - the option, this only responds to certain options, check the code for which ones 5979 5980 Output Parameter: 5981 . flg - turn the option on (`PETSC_TRUE`) or off (`PETSC_FALSE`) 5982 5983 Level: intermediate 5984 5985 Notes: 5986 Can only be called after `MatSetSizes()` and `MatSetType()` have been set. 5987 5988 Certain option values may be unknown, for those use the routines `MatIsSymmetric()`, `MatIsHermitian()`, `MatIsStructurallySymmetric()`, or 5989 `MatIsSymmetricKnown()`, `MatIsHermitianKnown()`, `MatIsStructurallySymmetricKnown()` 5990 5991 .seealso: `Mat`, `MatOption`, `MatSetOption()`, `MatIsSymmetric()`, `MatIsHermitian()`, `MatIsStructurallySymmetric()`, 5992 `MatIsSymmetricKnown()`, `MatIsHermitianKnown()`, `MatIsStructurallySymmetricKnown()` 5993 @*/ 5994 PetscErrorCode MatGetOption(Mat mat, MatOption op, PetscBool *flg) 5995 { 5996 PetscFunctionBegin; 5997 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 5998 PetscValidType(mat, 1); 5999 6000 PetscCheck(((int)op) > MAT_OPTION_MIN && ((int)op) < MAT_OPTION_MAX, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_OUTOFRANGE, "Options %d is out of range", (int)op); 6001 PetscCheck(((PetscObject)mat)->type_name, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_TYPENOTSET, "Cannot get options until type and size have been set, see MatSetType() and MatSetSizes()"); 6002 6003 switch (op) { 6004 case MAT_NO_OFF_PROC_ENTRIES: 6005 *flg = mat->nooffprocentries; 6006 break; 6007 case MAT_NO_OFF_PROC_ZERO_ROWS: 6008 *flg = mat->nooffproczerorows; 6009 break; 6010 case MAT_SYMMETRIC: 6011 SETERRQ(PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Use MatIsSymmetric() or MatIsSymmetricKnown()"); 6012 break; 6013 case MAT_HERMITIAN: 6014 SETERRQ(PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Use MatIsHermitian() or MatIsHermitianKnown()"); 6015 break; 6016 case MAT_STRUCTURALLY_SYMMETRIC: 6017 SETERRQ(PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Use MatIsStructurallySymmetric() or MatIsStructurallySymmetricKnown()"); 6018 break; 6019 case MAT_SPD: 6020 SETERRQ(PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Use MatIsSPDKnown()"); 6021 break; 6022 case MAT_SYMMETRY_ETERNAL: 6023 *flg = mat->symmetry_eternal; 6024 break; 6025 case MAT_STRUCTURAL_SYMMETRY_ETERNAL: 6026 *flg = mat->symmetry_eternal; 6027 break; 6028 default: 6029 break; 6030 } 6031 PetscFunctionReturn(PETSC_SUCCESS); 6032 } 6033 6034 /*@ 6035 MatZeroEntries - Zeros all entries of a matrix. For sparse matrices 6036 this routine retains the old nonzero structure. 6037 6038 Logically Collective 6039 6040 Input Parameters: 6041 . mat - the matrix 6042 6043 Level: intermediate 6044 6045 Note: 6046 If the matrix was not preallocated then a default, likely poor preallocation will be set in the matrix, so this should be called after the preallocation phase. 6047 See the Performance chapter of the users manual for information on preallocating matrices. 6048 6049 .seealso: `Mat`, `MatZeroRows()`, `MatZeroRowsColumns()` 6050 @*/ 6051 PetscErrorCode MatZeroEntries(Mat mat) 6052 { 6053 PetscFunctionBegin; 6054 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 6055 PetscValidType(mat, 1); 6056 PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 6057 PetscCheck(mat->insertmode == NOT_SET_VALUES, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for matrices where you have set values but not yet assembled"); 6058 MatCheckPreallocated(mat, 1); 6059 6060 PetscCall(PetscLogEventBegin(MAT_ZeroEntries, mat, 0, 0, 0)); 6061 PetscUseTypeMethod(mat, zeroentries); 6062 PetscCall(PetscLogEventEnd(MAT_ZeroEntries, mat, 0, 0, 0)); 6063 PetscCall(PetscObjectStateIncrease((PetscObject)mat)); 6064 PetscFunctionReturn(PETSC_SUCCESS); 6065 } 6066 6067 /*@ 6068 MatZeroRowsColumns - Zeros all entries (except possibly the main diagonal) 6069 of a set of rows and columns of a matrix. 6070 6071 Collective 6072 6073 Input Parameters: 6074 + mat - the matrix 6075 . numRows - the number of rows/columns to zero 6076 . rows - the global row indices 6077 . diag - value put in the diagonal of the eliminated rows 6078 . x - optional vector of the solution for zeroed rows (other entries in vector are not used), these must be set before this call 6079 - b - optional vector of the right hand side, that will be adjusted by provided solution entries 6080 6081 Level: intermediate 6082 6083 Notes: 6084 This routine, along with `MatZeroRows()`, is typically used to eliminate known Dirichlet boundary conditions from a linear system. 6085 6086 For each zeroed row, the value of the corresponding `b` is set to diag times the value of the corresponding `x`. 6087 The other entries of `b` will be adjusted by the known values of `x` times the corresponding matrix entries in the columns that are being eliminated 6088 6089 If the resulting linear system is to be solved with `KSP` then one can (but does not have to) call `KSPSetInitialGuessNonzero()` to allow the 6090 Krylov method to take advantage of the known solution on the zeroed rows. 6091 6092 For the parallel case, all processes that share the matrix (i.e., 6093 those in the communicator used for matrix creation) MUST call this 6094 routine, regardless of whether any rows being zeroed are owned by 6095 them. 6096 6097 Unlike `MatZeroRows()` this does not change the nonzero structure of the matrix, it merely zeros those entries in the matrix. 6098 6099 Each processor can indicate any rows in the entire matrix to be zeroed (i.e. each process does NOT have to 6100 list only rows local to itself). 6101 6102 The option `MAT_NO_OFF_PROC_ZERO_ROWS` does not apply to this routine. 6103 6104 .seealso: `Mat`, `MatZeroRowsIS()`, `MatZeroRows()`, `MatZeroRowsLocalIS()`, `MatZeroRowsStencil()`, `MatZeroEntries()`, `MatZeroRowsLocal()`, `MatSetOption()`, 6105 `MatZeroRowsColumnsLocal()`, `MatZeroRowsColumnsLocalIS()`, `MatZeroRowsColumnsIS()`, `MatZeroRowsColumnsStencil()` 6106 @*/ 6107 PetscErrorCode MatZeroRowsColumns(Mat mat, PetscInt numRows, const PetscInt rows[], PetscScalar diag, Vec x, Vec b) 6108 { 6109 PetscFunctionBegin; 6110 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 6111 PetscValidType(mat, 1); 6112 if (numRows) PetscValidIntPointer(rows, 3); 6113 PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 6114 PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 6115 MatCheckPreallocated(mat, 1); 6116 6117 PetscUseTypeMethod(mat, zerorowscolumns, numRows, rows, diag, x, b); 6118 PetscCall(MatViewFromOptions(mat, NULL, "-mat_view")); 6119 PetscCall(PetscObjectStateIncrease((PetscObject)mat)); 6120 PetscFunctionReturn(PETSC_SUCCESS); 6121 } 6122 6123 /*@ 6124 MatZeroRowsColumnsIS - Zeros all entries (except possibly the main diagonal) 6125 of a set of rows and columns of a matrix. 6126 6127 Collective 6128 6129 Input Parameters: 6130 + mat - the matrix 6131 . is - the rows to zero 6132 . diag - value put in all diagonals of eliminated rows (0.0 will even eliminate diagonal entry) 6133 . x - optional vector of solutions for zeroed rows (other entries in vector are not used) 6134 - b - optional vector of right hand side, that will be adjusted by provided solution 6135 6136 Level: intermediate 6137 6138 Note: 6139 See `MatZeroRowsColumns()` for details on how this routine operates. 6140 6141 .seealso: `Mat`, `MatZeroRowsIS()`, `MatZeroRowsColumns()`, `MatZeroRowsLocalIS()`, `MatZeroRowsStencil()`, `MatZeroEntries()`, `MatZeroRowsLocal()`, `MatSetOption()`, 6142 `MatZeroRowsColumnsLocal()`, `MatZeroRowsColumnsLocalIS()`, `MatZeroRows()`, `MatZeroRowsColumnsStencil()` 6143 @*/ 6144 PetscErrorCode MatZeroRowsColumnsIS(Mat mat, IS is, PetscScalar diag, Vec x, Vec b) 6145 { 6146 PetscInt numRows; 6147 const PetscInt *rows; 6148 6149 PetscFunctionBegin; 6150 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 6151 PetscValidHeaderSpecific(is, IS_CLASSID, 2); 6152 PetscValidType(mat, 1); 6153 PetscValidType(is, 2); 6154 PetscCall(ISGetLocalSize(is, &numRows)); 6155 PetscCall(ISGetIndices(is, &rows)); 6156 PetscCall(MatZeroRowsColumns(mat, numRows, rows, diag, x, b)); 6157 PetscCall(ISRestoreIndices(is, &rows)); 6158 PetscFunctionReturn(PETSC_SUCCESS); 6159 } 6160 6161 /*@ 6162 MatZeroRows - Zeros all entries (except possibly the main diagonal) 6163 of a set of rows of a matrix. 6164 6165 Collective 6166 6167 Input Parameters: 6168 + mat - the matrix 6169 . numRows - the number of rows to zero 6170 . rows - the global row indices 6171 . diag - value put in the diagonal of the zeroed rows 6172 . x - optional vector of solutions for zeroed rows (other entries in vector are not used), these must be set before this call 6173 - b - optional vector of right hand side, that will be adjusted by provided solution entries 6174 6175 Level: intermediate 6176 6177 Notes: 6178 This routine, along with `MatZeroRowsColumns()`, is typically used to eliminate known Dirichlet boundary conditions from a linear system. 6179 6180 For each zeroed row, the value of the corresponding `b` is set to `diag` times the value of the corresponding `x`. 6181 6182 If the resulting linear system is to be solved with `KSP` then one can (but does not have to) call `KSPSetInitialGuessNonzero()` to allow the 6183 Krylov method to take advantage of the known solution on the zeroed rows. 6184 6185 May be followed by using a `PC` of type `PCREDISTRIBUTE` to solve the reduced problem (`PCDISTRIBUTE` completely eliminates the zeroed rows and their corresponding columns) 6186 from the matrix. 6187 6188 Unlike `MatZeroRowsColumns()` for the `MATAIJ` and `MATBAIJ` matrix formats this removes the old nonzero structure, from the eliminated rows of the matrix 6189 but does not release memory. Because of this removal matrix-vector products with the adjusted matrix will be a bit faster. For the dense and block diagonal 6190 formats this does not alter the nonzero structure. 6191 6192 If the option `MatSetOption`(mat,`MAT_KEEP_NONZERO_PATTERN`,`PETSC_TRUE`) the nonzero structure 6193 of the matrix is not changed the values are 6194 merely zeroed. 6195 6196 The user can set a value in the diagonal entry (or for the `MATAIJ` format 6197 formats can optionally remove the main diagonal entry from the 6198 nonzero structure as well, by passing 0.0 as the final argument). 6199 6200 For the parallel case, all processes that share the matrix (i.e., 6201 those in the communicator used for matrix creation) MUST call this 6202 routine, regardless of whether any rows being zeroed are owned by 6203 them. 6204 6205 Each processor can indicate any rows in the entire matrix to be zeroed (i.e. each process does NOT have to 6206 list only rows local to itself). 6207 6208 You can call `MatSetOption`(mat,`MAT_NO_OFF_PROC_ZERO_ROWS`,`PETSC_TRUE`) if each process indicates only rows it 6209 owns that are to be zeroed. This saves a global synchronization in the implementation. 6210 6211 .seealso: `Mat`, `MatZeroRowsIS()`, `MatZeroRowsColumns()`, `MatZeroRowsLocalIS()`, `MatZeroRowsStencil()`, `MatZeroEntries()`, `MatZeroRowsLocal()`, `MatSetOption()`, 6212 `MatZeroRowsColumnsLocal()`, `MatZeroRowsColumnsLocalIS()`, `MatZeroRowsColumnsIS()`, `MatZeroRowsColumnsStencil()`, `PCREDISTRIBUTE` 6213 @*/ 6214 PetscErrorCode MatZeroRows(Mat mat, PetscInt numRows, const PetscInt rows[], PetscScalar diag, Vec x, Vec b) 6215 { 6216 PetscFunctionBegin; 6217 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 6218 PetscValidType(mat, 1); 6219 if (numRows) PetscValidIntPointer(rows, 3); 6220 PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 6221 PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 6222 MatCheckPreallocated(mat, 1); 6223 6224 PetscUseTypeMethod(mat, zerorows, numRows, rows, diag, x, b); 6225 PetscCall(MatViewFromOptions(mat, NULL, "-mat_view")); 6226 PetscCall(PetscObjectStateIncrease((PetscObject)mat)); 6227 PetscFunctionReturn(PETSC_SUCCESS); 6228 } 6229 6230 /*@ 6231 MatZeroRowsIS - Zeros all entries (except possibly the main diagonal) 6232 of a set of rows of a matrix. 6233 6234 Collective on Mat 6235 6236 Input Parameters: 6237 + mat - the matrix 6238 . is - index set of rows to remove (if `NULL` then no row is removed) 6239 . diag - value put in all diagonals of eliminated rows 6240 . x - optional vector of solutions for zeroed rows (other entries in vector are not used) 6241 - b - optional vector of right hand side, that will be adjusted by provided solution 6242 6243 Level: intermediate 6244 6245 Note: 6246 See `MatZeroRows()` for details on how this routine operates. 6247 6248 .seealso: `Mat`, `MatZeroRows()`, `MatZeroRowsColumns()`, `MatZeroRowsLocalIS()`, `MatZeroRowsStencil()`, `MatZeroEntries()`, `MatZeroRowsLocal()`, `MatSetOption()`, 6249 `MatZeroRowsColumnsLocal()`, `MatZeroRowsColumnsLocalIS()`, `MatZeroRowsColumnsIS()`, `MatZeroRowsColumnsStencil()` 6250 @*/ 6251 PetscErrorCode MatZeroRowsIS(Mat mat, IS is, PetscScalar diag, Vec x, Vec b) 6252 { 6253 PetscInt numRows = 0; 6254 const PetscInt *rows = NULL; 6255 6256 PetscFunctionBegin; 6257 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 6258 PetscValidType(mat, 1); 6259 if (is) { 6260 PetscValidHeaderSpecific(is, IS_CLASSID, 2); 6261 PetscCall(ISGetLocalSize(is, &numRows)); 6262 PetscCall(ISGetIndices(is, &rows)); 6263 } 6264 PetscCall(MatZeroRows(mat, numRows, rows, diag, x, b)); 6265 if (is) PetscCall(ISRestoreIndices(is, &rows)); 6266 PetscFunctionReturn(PETSC_SUCCESS); 6267 } 6268 6269 /*@ 6270 MatZeroRowsStencil - Zeros all entries (except possibly the main diagonal) 6271 of a set of rows of a matrix. These rows must be local to the process. 6272 6273 Collective 6274 6275 Input Parameters: 6276 + mat - the matrix 6277 . numRows - the number of rows to remove 6278 . rows - the grid coordinates (and component number when dof > 1) for matrix rows 6279 . diag - value put in all diagonals of eliminated rows (0.0 will even eliminate diagonal entry) 6280 . x - optional vector of solutions for zeroed rows (other entries in vector are not used) 6281 - b - optional vector of right hand side, that will be adjusted by provided solution 6282 6283 Level: intermediate 6284 6285 Notes: 6286 See `MatZeroRows()` for details on how this routine operates. 6287 6288 The grid coordinates are across the entire grid, not just the local portion 6289 6290 For periodic boundary conditions use negative indices for values to the left (below 0; that are to be 6291 obtained by wrapping values from right edge). For values to the right of the last entry using that index plus one 6292 etc to obtain values that obtained by wrapping the values from the left edge. This does not work for anything but the 6293 `DM_BOUNDARY_PERIODIC` boundary type. 6294 6295 For indices that don't mean anything for your case (like the k index when working in 2d) or the c index when you have 6296 a single value per point) you can skip filling those indices. 6297 6298 Fortran Note: 6299 idxm and idxn should be declared as 6300 $ MatStencil idxm(4,m) 6301 and the values inserted using 6302 .vb 6303 idxm(MatStencil_i,1) = i 6304 idxm(MatStencil_j,1) = j 6305 idxm(MatStencil_k,1) = k 6306 idxm(MatStencil_c,1) = c 6307 etc 6308 .ve 6309 6310 .seealso: `Mat`, `MatZeroRowsIS()`, `MatZeroRowsColumns()`, `MatZeroRowsLocalIS()`, `MatZeroRowsl()`, `MatZeroEntries()`, `MatZeroRowsLocal()`, `MatSetOption()`, 6311 `MatZeroRowsColumnsLocal()`, `MatZeroRowsColumnsLocalIS()`, `MatZeroRowsColumnsIS()`, `MatZeroRowsColumnsStencil()` 6312 @*/ 6313 PetscErrorCode MatZeroRowsStencil(Mat mat, PetscInt numRows, const MatStencil rows[], PetscScalar diag, Vec x, Vec b) 6314 { 6315 PetscInt dim = mat->stencil.dim; 6316 PetscInt sdim = dim - (1 - (PetscInt)mat->stencil.noc); 6317 PetscInt *dims = mat->stencil.dims + 1; 6318 PetscInt *starts = mat->stencil.starts; 6319 PetscInt *dxm = (PetscInt *)rows; 6320 PetscInt *jdxm, i, j, tmp, numNewRows = 0; 6321 6322 PetscFunctionBegin; 6323 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 6324 PetscValidType(mat, 1); 6325 if (numRows) PetscValidPointer(rows, 3); 6326 6327 PetscCall(PetscMalloc1(numRows, &jdxm)); 6328 for (i = 0; i < numRows; ++i) { 6329 /* Skip unused dimensions (they are ordered k, j, i, c) */ 6330 for (j = 0; j < 3 - sdim; ++j) dxm++; 6331 /* Local index in X dir */ 6332 tmp = *dxm++ - starts[0]; 6333 /* Loop over remaining dimensions */ 6334 for (j = 0; j < dim - 1; ++j) { 6335 /* If nonlocal, set index to be negative */ 6336 if ((*dxm++ - starts[j + 1]) < 0 || tmp < 0) tmp = PETSC_MIN_INT; 6337 /* Update local index */ 6338 else tmp = tmp * dims[j] + *(dxm - 1) - starts[j + 1]; 6339 } 6340 /* Skip component slot if necessary */ 6341 if (mat->stencil.noc) dxm++; 6342 /* Local row number */ 6343 if (tmp >= 0) jdxm[numNewRows++] = tmp; 6344 } 6345 PetscCall(MatZeroRowsLocal(mat, numNewRows, jdxm, diag, x, b)); 6346 PetscCall(PetscFree(jdxm)); 6347 PetscFunctionReturn(PETSC_SUCCESS); 6348 } 6349 6350 /*@ 6351 MatZeroRowsColumnsStencil - Zeros all row and column entries (except possibly the main diagonal) 6352 of a set of rows and columns of a matrix. 6353 6354 Collective 6355 6356 Input Parameters: 6357 + mat - the matrix 6358 . numRows - the number of rows/columns to remove 6359 . rows - the grid coordinates (and component number when dof > 1) for matrix rows 6360 . diag - value put in all diagonals of eliminated rows (0.0 will even eliminate diagonal entry) 6361 . x - optional vector of solutions for zeroed rows (other entries in vector are not used) 6362 - b - optional vector of right hand side, that will be adjusted by provided solution 6363 6364 Level: intermediate 6365 6366 Notes: 6367 See `MatZeroRowsColumns()` for details on how this routine operates. 6368 6369 The grid coordinates are across the entire grid, not just the local portion 6370 6371 For periodic boundary conditions use negative indices for values to the left (below 0; that are to be 6372 obtained by wrapping values from right edge). For values to the right of the last entry using that index plus one 6373 etc to obtain values that obtained by wrapping the values from the left edge. This does not work for anything but the 6374 `DM_BOUNDARY_PERIODIC` boundary type. 6375 6376 For indices that don't mean anything for your case (like the k index when working in 2d) or the c index when you have 6377 a single value per point) you can skip filling those indices. 6378 6379 Fortran Note: 6380 In Fortran idxm and idxn should be declared as 6381 $ MatStencil idxm(4,m) 6382 and the values inserted using 6383 .vb 6384 idxm(MatStencil_i,1) = i 6385 idxm(MatStencil_j,1) = j 6386 idxm(MatStencil_k,1) = k 6387 idxm(MatStencil_c,1) = c 6388 etc 6389 .ve 6390 6391 .seealso: `Mat`, `MatZeroRowsIS()`, `MatZeroRowsColumns()`, `MatZeroRowsLocalIS()`, `MatZeroRowsStencil()`, `MatZeroEntries()`, `MatZeroRowsLocal()`, `MatSetOption()`, 6392 `MatZeroRowsColumnsLocal()`, `MatZeroRowsColumnsLocalIS()`, `MatZeroRowsColumnsIS()`, `MatZeroRows()` 6393 @*/ 6394 PetscErrorCode MatZeroRowsColumnsStencil(Mat mat, PetscInt numRows, const MatStencil rows[], PetscScalar diag, Vec x, Vec b) 6395 { 6396 PetscInt dim = mat->stencil.dim; 6397 PetscInt sdim = dim - (1 - (PetscInt)mat->stencil.noc); 6398 PetscInt *dims = mat->stencil.dims + 1; 6399 PetscInt *starts = mat->stencil.starts; 6400 PetscInt *dxm = (PetscInt *)rows; 6401 PetscInt *jdxm, i, j, tmp, numNewRows = 0; 6402 6403 PetscFunctionBegin; 6404 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 6405 PetscValidType(mat, 1); 6406 if (numRows) PetscValidPointer(rows, 3); 6407 6408 PetscCall(PetscMalloc1(numRows, &jdxm)); 6409 for (i = 0; i < numRows; ++i) { 6410 /* Skip unused dimensions (they are ordered k, j, i, c) */ 6411 for (j = 0; j < 3 - sdim; ++j) dxm++; 6412 /* Local index in X dir */ 6413 tmp = *dxm++ - starts[0]; 6414 /* Loop over remaining dimensions */ 6415 for (j = 0; j < dim - 1; ++j) { 6416 /* If nonlocal, set index to be negative */ 6417 if ((*dxm++ - starts[j + 1]) < 0 || tmp < 0) tmp = PETSC_MIN_INT; 6418 /* Update local index */ 6419 else tmp = tmp * dims[j] + *(dxm - 1) - starts[j + 1]; 6420 } 6421 /* Skip component slot if necessary */ 6422 if (mat->stencil.noc) dxm++; 6423 /* Local row number */ 6424 if (tmp >= 0) jdxm[numNewRows++] = tmp; 6425 } 6426 PetscCall(MatZeroRowsColumnsLocal(mat, numNewRows, jdxm, diag, x, b)); 6427 PetscCall(PetscFree(jdxm)); 6428 PetscFunctionReturn(PETSC_SUCCESS); 6429 } 6430 6431 /*@C 6432 MatZeroRowsLocal - Zeros all entries (except possibly the main diagonal) 6433 of a set of rows of a matrix; using local numbering of rows. 6434 6435 Collective 6436 6437 Input Parameters: 6438 + mat - the matrix 6439 . numRows - the number of rows to remove 6440 . rows - the local row indices 6441 . diag - value put in all diagonals of eliminated rows 6442 . x - optional vector of solutions for zeroed rows (other entries in vector are not used) 6443 - b - optional vector of right hand side, that will be adjusted by provided solution 6444 6445 Level: intermediate 6446 6447 Notes: 6448 Before calling `MatZeroRowsLocal()`, the user must first set the 6449 local-to-global mapping by calling MatSetLocalToGlobalMapping(), this is often already set for matrices obtained with `DMCreateMatrix()`. 6450 6451 See `MatZeroRows()` for details on how this routine operates. 6452 6453 .seealso: `Mat`, `MatZeroRowsIS()`, `MatZeroRowsColumns()`, `MatZeroRowsLocalIS()`, `MatZeroRowsStencil()`, `MatZeroEntries()`, `MatZeroRows()`, `MatSetOption()`, 6454 `MatZeroRowsColumnsLocal()`, `MatZeroRowsColumnsLocalIS()`, `MatZeroRowsColumnsIS()`, `MatZeroRowsColumnsStencil()` 6455 @*/ 6456 PetscErrorCode MatZeroRowsLocal(Mat mat, PetscInt numRows, const PetscInt rows[], PetscScalar diag, Vec x, Vec b) 6457 { 6458 PetscFunctionBegin; 6459 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 6460 PetscValidType(mat, 1); 6461 if (numRows) PetscValidIntPointer(rows, 3); 6462 PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 6463 PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 6464 MatCheckPreallocated(mat, 1); 6465 6466 if (mat->ops->zerorowslocal) { 6467 PetscUseTypeMethod(mat, zerorowslocal, numRows, rows, diag, x, b); 6468 } else { 6469 IS is, newis; 6470 const PetscInt *newRows; 6471 6472 PetscCheck(mat->rmap->mapping, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Need to provide local to global mapping to matrix first"); 6473 PetscCall(ISCreateGeneral(PETSC_COMM_SELF, numRows, rows, PETSC_COPY_VALUES, &is)); 6474 PetscCall(ISLocalToGlobalMappingApplyIS(mat->rmap->mapping, is, &newis)); 6475 PetscCall(ISGetIndices(newis, &newRows)); 6476 PetscUseTypeMethod(mat, zerorows, numRows, newRows, diag, x, b); 6477 PetscCall(ISRestoreIndices(newis, &newRows)); 6478 PetscCall(ISDestroy(&newis)); 6479 PetscCall(ISDestroy(&is)); 6480 } 6481 PetscCall(PetscObjectStateIncrease((PetscObject)mat)); 6482 PetscFunctionReturn(PETSC_SUCCESS); 6483 } 6484 6485 /*@ 6486 MatZeroRowsLocalIS - Zeros all entries (except possibly the main diagonal) 6487 of a set of rows of a matrix; using local numbering of rows. 6488 6489 Collective 6490 6491 Input Parameters: 6492 + mat - the matrix 6493 . is - index set of rows to remove 6494 . diag - value put in all diagonals of eliminated rows 6495 . x - optional vector of solutions for zeroed rows (other entries in vector are not used) 6496 - b - optional vector of right hand side, that will be adjusted by provided solution 6497 6498 Level: intermediate 6499 6500 Notes: 6501 Before calling `MatZeroRowsLocalIS()`, the user must first set the 6502 local-to-global mapping by calling `MatSetLocalToGlobalMapping()`, this is often already set for matrices obtained with `DMCreateMatrix()`. 6503 6504 See `MatZeroRows()` for details on how this routine operates. 6505 6506 .seealso: `Mat`, `MatZeroRowsIS()`, `MatZeroRowsColumns()`, `MatZeroRows()`, `MatZeroRowsStencil()`, `MatZeroEntries()`, `MatZeroRowsLocal()`, `MatSetOption()`, 6507 `MatZeroRowsColumnsLocal()`, `MatZeroRowsColumnsLocalIS()`, `MatZeroRowsColumnsIS()`, `MatZeroRowsColumnsStencil()` 6508 @*/ 6509 PetscErrorCode MatZeroRowsLocalIS(Mat mat, IS is, PetscScalar diag, Vec x, Vec b) 6510 { 6511 PetscInt numRows; 6512 const PetscInt *rows; 6513 6514 PetscFunctionBegin; 6515 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 6516 PetscValidType(mat, 1); 6517 PetscValidHeaderSpecific(is, IS_CLASSID, 2); 6518 PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 6519 PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 6520 MatCheckPreallocated(mat, 1); 6521 6522 PetscCall(ISGetLocalSize(is, &numRows)); 6523 PetscCall(ISGetIndices(is, &rows)); 6524 PetscCall(MatZeroRowsLocal(mat, numRows, rows, diag, x, b)); 6525 PetscCall(ISRestoreIndices(is, &rows)); 6526 PetscFunctionReturn(PETSC_SUCCESS); 6527 } 6528 6529 /*@ 6530 MatZeroRowsColumnsLocal - Zeros all entries (except possibly the main diagonal) 6531 of a set of rows and columns of a matrix; using local numbering of rows. 6532 6533 Collective 6534 6535 Input Parameters: 6536 + mat - the matrix 6537 . numRows - the number of rows to remove 6538 . rows - the global row indices 6539 . diag - value put in all diagonals of eliminated rows 6540 . x - optional vector of solutions for zeroed rows (other entries in vector are not used) 6541 - b - optional vector of right hand side, that will be adjusted by provided solution 6542 6543 Level: intermediate 6544 6545 Notes: 6546 Before calling `MatZeroRowsColumnsLocal()`, the user must first set the 6547 local-to-global mapping by calling `MatSetLocalToGlobalMapping()`, this is often already set for matrices obtained with `DMCreateMatrix()`. 6548 6549 See `MatZeroRowsColumns()` for details on how this routine operates. 6550 6551 .seealso: `Mat`, `MatZeroRowsIS()`, `MatZeroRowsColumns()`, `MatZeroRowsLocalIS()`, `MatZeroRowsStencil()`, `MatZeroEntries()`, `MatZeroRowsLocal()`, `MatSetOption()`, 6552 `MatZeroRows()`, `MatZeroRowsColumnsLocalIS()`, `MatZeroRowsColumnsIS()`, `MatZeroRowsColumnsStencil()` 6553 @*/ 6554 PetscErrorCode MatZeroRowsColumnsLocal(Mat mat, PetscInt numRows, const PetscInt rows[], PetscScalar diag, Vec x, Vec b) 6555 { 6556 IS is, newis; 6557 const PetscInt *newRows; 6558 6559 PetscFunctionBegin; 6560 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 6561 PetscValidType(mat, 1); 6562 if (numRows) PetscValidIntPointer(rows, 3); 6563 PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 6564 PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 6565 MatCheckPreallocated(mat, 1); 6566 6567 PetscCheck(mat->cmap->mapping, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Need to provide local to global mapping to matrix first"); 6568 PetscCall(ISCreateGeneral(PETSC_COMM_SELF, numRows, rows, PETSC_COPY_VALUES, &is)); 6569 PetscCall(ISLocalToGlobalMappingApplyIS(mat->cmap->mapping, is, &newis)); 6570 PetscCall(ISGetIndices(newis, &newRows)); 6571 PetscUseTypeMethod(mat, zerorowscolumns, numRows, newRows, diag, x, b); 6572 PetscCall(ISRestoreIndices(newis, &newRows)); 6573 PetscCall(ISDestroy(&newis)); 6574 PetscCall(ISDestroy(&is)); 6575 PetscCall(PetscObjectStateIncrease((PetscObject)mat)); 6576 PetscFunctionReturn(PETSC_SUCCESS); 6577 } 6578 6579 /*@ 6580 MatZeroRowsColumnsLocalIS - Zeros all entries (except possibly the main diagonal) 6581 of a set of rows and columns of a matrix; using local numbering of rows. 6582 6583 Collective on Mat 6584 6585 Input Parameters: 6586 + mat - the matrix 6587 . is - index set of rows to remove 6588 . diag - value put in all diagonals of eliminated rows 6589 . x - optional vector of solutions for zeroed rows (other entries in vector are not used) 6590 - b - optional vector of right hand side, that will be adjusted by provided solution 6591 6592 Level: intermediate 6593 6594 Notes: 6595 Before calling `MatZeroRowsColumnsLocalIS()`, the user must first set the 6596 local-to-global mapping by calling `MatSetLocalToGlobalMapping()`, this is often already set for matrices obtained with `DMCreateMatrix()`. 6597 6598 See `MatZeroRowsColumns()` for details on how this routine operates. 6599 6600 .seealso: `Mat`, `MatZeroRowsIS()`, `MatZeroRowsColumns()`, `MatZeroRowsLocalIS()`, `MatZeroRowsStencil()`, `MatZeroEntries()`, `MatZeroRowsLocal()`, `MatSetOption()`, 6601 `MatZeroRowsColumnsLocal()`, `MatZeroRows()`, `MatZeroRowsColumnsIS()`, `MatZeroRowsColumnsStencil()` 6602 @*/ 6603 PetscErrorCode MatZeroRowsColumnsLocalIS(Mat mat, IS is, PetscScalar diag, Vec x, Vec b) 6604 { 6605 PetscInt numRows; 6606 const PetscInt *rows; 6607 6608 PetscFunctionBegin; 6609 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 6610 PetscValidType(mat, 1); 6611 PetscValidHeaderSpecific(is, IS_CLASSID, 2); 6612 PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 6613 PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 6614 MatCheckPreallocated(mat, 1); 6615 6616 PetscCall(ISGetLocalSize(is, &numRows)); 6617 PetscCall(ISGetIndices(is, &rows)); 6618 PetscCall(MatZeroRowsColumnsLocal(mat, numRows, rows, diag, x, b)); 6619 PetscCall(ISRestoreIndices(is, &rows)); 6620 PetscFunctionReturn(PETSC_SUCCESS); 6621 } 6622 6623 /*@C 6624 MatGetSize - Returns the numbers of rows and columns in a matrix. 6625 6626 Not Collective 6627 6628 Input Parameter: 6629 . mat - the matrix 6630 6631 Level: beginner 6632 6633 Output Parameters: 6634 + m - the number of global rows 6635 - n - the number of global columns 6636 6637 Note: 6638 Both output parameters can be `NULL` on input. 6639 6640 .seealso: `Mat`, `MatSetSizes()`, `MatGetLocalSize()` 6641 @*/ 6642 PetscErrorCode MatGetSize(Mat mat, PetscInt *m, PetscInt *n) 6643 { 6644 PetscFunctionBegin; 6645 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 6646 if (m) *m = mat->rmap->N; 6647 if (n) *n = mat->cmap->N; 6648 PetscFunctionReturn(PETSC_SUCCESS); 6649 } 6650 6651 /*@C 6652 MatGetLocalSize - For most matrix formats, excluding `MATELEMENTAL` and `MATSCALAPACK`, Returns the number of local rows and local columns 6653 of a matrix. For all matrices this is the local size of the left and right vectors as returned by `MatCreateVecs()`. 6654 6655 Not Collective 6656 6657 Input Parameter: 6658 . mat - the matrix 6659 6660 Output Parameters: 6661 + m - the number of local rows, use `NULL` to not obtain this value 6662 - n - the number of local columns, use `NULL` to not obtain this value 6663 6664 Level: beginner 6665 6666 .seealso: `Mat`, `MatSetSizes()`, `MatGetSize()` 6667 @*/ 6668 PetscErrorCode MatGetLocalSize(Mat mat, PetscInt *m, PetscInt *n) 6669 { 6670 PetscFunctionBegin; 6671 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 6672 if (m) PetscValidIntPointer(m, 2); 6673 if (n) PetscValidIntPointer(n, 3); 6674 if (m) *m = mat->rmap->n; 6675 if (n) *n = mat->cmap->n; 6676 PetscFunctionReturn(PETSC_SUCCESS); 6677 } 6678 6679 /*@C 6680 MatGetOwnershipRangeColumn - Returns the range of matrix columns associated with rows of a vector one multiplies this matrix by that are owned by 6681 this processor. (The columns of the "diagonal block" for most sparse matrix formats). See :any:`<sec_matlayout>` for details on matrix layouts. 6682 6683 Not Collective, unless matrix has not been allocated, then collective 6684 6685 Input Parameter: 6686 . mat - the matrix 6687 6688 Output Parameters: 6689 + m - the global index of the first local column, use `NULL` to not obtain this value 6690 - n - one more than the global index of the last local column, use `NULL` to not obtain this value 6691 6692 Level: developer 6693 6694 .seealso: `Mat`, `MatGetOwnershipRange()`, `MatGetOwnershipRanges()`, `MatGetOwnershipRangesColumn()`, `PetscLayout` 6695 @*/ 6696 PetscErrorCode MatGetOwnershipRangeColumn(Mat mat, PetscInt *m, PetscInt *n) 6697 { 6698 PetscFunctionBegin; 6699 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 6700 PetscValidType(mat, 1); 6701 if (m) PetscValidIntPointer(m, 2); 6702 if (n) PetscValidIntPointer(n, 3); 6703 MatCheckPreallocated(mat, 1); 6704 if (m) *m = mat->cmap->rstart; 6705 if (n) *n = mat->cmap->rend; 6706 PetscFunctionReturn(PETSC_SUCCESS); 6707 } 6708 6709 /*@C 6710 MatGetOwnershipRange - For matrices that own values by row, excludes `MATELEMENTAL` and `MATSCALAPACK`, returns the range of matrix rows owned by 6711 this MPI rank. For all matrices it returns the range of matrix rows associated with rows of a vector that would contain the result of a matrix 6712 vector product with this matrix. See :any:`<sec_matlayout>` for details on matrix layouts 6713 6714 Not Collective 6715 6716 Input Parameter: 6717 . mat - the matrix 6718 6719 Output Parameters: 6720 + m - the global index of the first local row, use `NULL` to not obtain this value 6721 - n - one more than the global index of the last local row, use `NULL` to not obtain this value 6722 6723 Level: beginner 6724 6725 Note: 6726 This function requires that the matrix be preallocated. If you have not preallocated, consider using 6727 `PetscSplitOwnership`(`MPI_Comm` comm, `PetscInt` *n, `PetscInt` *N) 6728 and then `MPI_Scan()` to calculate prefix sums of the local sizes. 6729 6730 .seealso: `Mat`, `MatGetOwnershipRanges()`, `MatGetOwnershipRangeColumn()`, `MatGetOwnershipRangesColumn()`, `PetscSplitOwnership()`, `PetscSplitOwnershipBlock()`, 6731 `PetscLayout` 6732 @*/ 6733 PetscErrorCode MatGetOwnershipRange(Mat mat, PetscInt *m, PetscInt *n) 6734 { 6735 PetscFunctionBegin; 6736 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 6737 PetscValidType(mat, 1); 6738 if (m) PetscValidIntPointer(m, 2); 6739 if (n) PetscValidIntPointer(n, 3); 6740 MatCheckPreallocated(mat, 1); 6741 if (m) *m = mat->rmap->rstart; 6742 if (n) *n = mat->rmap->rend; 6743 PetscFunctionReturn(PETSC_SUCCESS); 6744 } 6745 6746 /*@C 6747 MatGetOwnershipRanges - For matrices that own values by row, excludes `MATELEMENTAL` and `MATSCALAPACK`, returns the range of matrix rows owned by 6748 each process. For all matrices it returns the ranges of matrix rows associated with rows of a vector that would contain the result of a matrix 6749 vector product with this matrix. See :any:`<sec_matlayout>` for details on matrix layouts 6750 6751 Not Collective, unless matrix has not been allocated, then collective 6752 6753 Input Parameters: 6754 . mat - the matrix 6755 6756 Output Parameters: 6757 . ranges - start of each processors portion plus one more than the total length at the end 6758 6759 Level: beginner 6760 6761 .seealso: `Mat`, `MatGetOwnershipRange()`, `MatGetOwnershipRangeColumn()`, `MatGetOwnershipRangesColumn()`, `PetscLayout` 6762 @*/ 6763 PetscErrorCode MatGetOwnershipRanges(Mat mat, const PetscInt **ranges) 6764 { 6765 PetscFunctionBegin; 6766 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 6767 PetscValidType(mat, 1); 6768 MatCheckPreallocated(mat, 1); 6769 PetscCall(PetscLayoutGetRanges(mat->rmap, ranges)); 6770 PetscFunctionReturn(PETSC_SUCCESS); 6771 } 6772 6773 /*@C 6774 MatGetOwnershipRangesColumn - Returns the ranges of matrix columns associated with rows of a vector one multiplies this vector by that are owned by 6775 each processor. (The columns of the "diagonal blocks", for most sparse matrix formats). See :any:`<sec_matlayout>` for details on matrix layouts. 6776 6777 Not Collective, unless matrix has not been allocated, then collective on Mat 6778 6779 Input Parameters: 6780 . mat - the matrix 6781 6782 Output Parameters: 6783 . ranges - start of each processors portion plus one more then the total length at the end 6784 6785 Level: beginner 6786 6787 .seealso: `Mat`, `MatGetOwnershipRange()`, `MatGetOwnershipRangeColumn()`, `MatGetOwnershipRanges()` 6788 @*/ 6789 PetscErrorCode MatGetOwnershipRangesColumn(Mat mat, const PetscInt **ranges) 6790 { 6791 PetscFunctionBegin; 6792 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 6793 PetscValidType(mat, 1); 6794 MatCheckPreallocated(mat, 1); 6795 PetscCall(PetscLayoutGetRanges(mat->cmap, ranges)); 6796 PetscFunctionReturn(PETSC_SUCCESS); 6797 } 6798 6799 /*@C 6800 MatGetOwnershipIS - Get row and column ownership of a matrices' values as index sets. For most matrices, excluding `MATELEMENTAL` and `MATSCALAPACK`, this 6801 corresponds to values returned by `MatGetOwnershipRange()`, `MatGetOwnershipRangeColumn()`. For `MATELEMENTAL` and `MATSCALAPACK` the ownership 6802 is more complicated. See :any:`<sec_matlayout>` for details on matrix layouts. 6803 6804 Not Collective 6805 6806 Input Parameter: 6807 . A - matrix 6808 6809 Output Parameters: 6810 + rows - rows in which this process owns elements, , use `NULL` to not obtain this value 6811 - cols - columns in which this process owns elements, use `NULL` to not obtain this value 6812 6813 Level: intermediate 6814 6815 .seealso: `Mat`, `MatGetOwnershipRange()`, `MatGetOwnershipRangeColumn()`, `MatSetValues()`, ``MATELEMENTAL``, ``MATSCALAPACK`` 6816 @*/ 6817 PetscErrorCode MatGetOwnershipIS(Mat A, IS *rows, IS *cols) 6818 { 6819 PetscErrorCode (*f)(Mat, IS *, IS *); 6820 6821 PetscFunctionBegin; 6822 MatCheckPreallocated(A, 1); 6823 PetscCall(PetscObjectQueryFunction((PetscObject)A, "MatGetOwnershipIS_C", &f)); 6824 if (f) { 6825 PetscCall((*f)(A, rows, cols)); 6826 } else { /* Create a standard row-based partition, each process is responsible for ALL columns in their row block */ 6827 if (rows) PetscCall(ISCreateStride(PETSC_COMM_SELF, A->rmap->n, A->rmap->rstart, 1, rows)); 6828 if (cols) PetscCall(ISCreateStride(PETSC_COMM_SELF, A->cmap->N, 0, 1, cols)); 6829 } 6830 PetscFunctionReturn(PETSC_SUCCESS); 6831 } 6832 6833 /*@C 6834 MatILUFactorSymbolic - Performs symbolic ILU factorization of a matrix obtained with `MatGetFactor()` 6835 Uses levels of fill only, not drop tolerance. Use `MatLUFactorNumeric()` 6836 to complete the factorization. 6837 6838 Collective on fact 6839 6840 Input Parameters: 6841 + fact - the factorized matrix obtained with `MatGetFactor()` 6842 . mat - the matrix 6843 . row - row permutation 6844 . column - column permutation 6845 - info - structure containing 6846 .vb 6847 levels - number of levels of fill. 6848 expected fill - as ratio of original fill. 6849 1 or 0 - indicating force fill on diagonal (improves robustness for matrices 6850 missing diagonal entries) 6851 .ve 6852 6853 Output Parameters: 6854 . fact - new matrix that has been symbolically factored 6855 6856 Level: developer 6857 6858 Notes: 6859 See [Matrix Factorization](sec_matfactor) for additional information. 6860 6861 Most users should employ the `KSP` interface for linear solvers 6862 instead of working directly with matrix algebra routines such as this. 6863 See, e.g., `KSPCreate()`. 6864 6865 Uses the definition of level of fill as in Y. Saad, 2003 6866 6867 Developer Note: 6868 The Fortran interface is not autogenerated as the 6869 interface definition cannot be generated correctly [due to `MatFactorInfo`] 6870 6871 References: 6872 . * - Y. Saad, Iterative methods for sparse linear systems Philadelphia: Society for Industrial and Applied Mathematics, 2003 6873 6874 .seealso: `Mat`, [Matrix Factorization](sec_matfactor), `MatGetFactor()`, `MatLUFactorSymbolic()`, `MatLUFactorNumeric()`, `MatCholeskyFactor()` 6875 `MatGetOrdering()`, `MatFactorInfo` 6876 @*/ 6877 PetscErrorCode MatILUFactorSymbolic(Mat fact, Mat mat, IS row, IS col, const MatFactorInfo *info) 6878 { 6879 PetscFunctionBegin; 6880 PetscValidHeaderSpecific(mat, MAT_CLASSID, 2); 6881 PetscValidType(mat, 2); 6882 if (row) PetscValidHeaderSpecific(row, IS_CLASSID, 3); 6883 if (col) PetscValidHeaderSpecific(col, IS_CLASSID, 4); 6884 PetscValidPointer(info, 5); 6885 PetscValidPointer(fact, 1); 6886 PetscCheck(info->levels >= 0, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_OUTOFRANGE, "Levels of fill negative %" PetscInt_FMT, (PetscInt)info->levels); 6887 PetscCheck(info->fill >= 1.0, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_OUTOFRANGE, "Expected fill less than 1.0 %g", (double)info->fill); 6888 if (!fact->ops->ilufactorsymbolic) { 6889 MatSolverType stype; 6890 PetscCall(MatFactorGetSolverType(fact, &stype)); 6891 SETERRQ(PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Matrix type %s symbolic ILU using solver type %s", ((PetscObject)mat)->type_name, stype); 6892 } 6893 PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 6894 PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 6895 MatCheckPreallocated(mat, 2); 6896 6897 if (!fact->trivialsymbolic) PetscCall(PetscLogEventBegin(MAT_ILUFactorSymbolic, mat, row, col, 0)); 6898 PetscCall((fact->ops->ilufactorsymbolic)(fact, mat, row, col, info)); 6899 if (!fact->trivialsymbolic) PetscCall(PetscLogEventEnd(MAT_ILUFactorSymbolic, mat, row, col, 0)); 6900 PetscFunctionReturn(PETSC_SUCCESS); 6901 } 6902 6903 /*@C 6904 MatICCFactorSymbolic - Performs symbolic incomplete 6905 Cholesky factorization for a symmetric matrix. Use 6906 `MatCholeskyFactorNumeric()` to complete the factorization. 6907 6908 Collective on fact 6909 6910 Input Parameters: 6911 + fact - the factorized matrix obtained with `MatGetFactor()` 6912 . mat - the matrix to be factored 6913 . perm - row and column permutation 6914 - info - structure containing 6915 .vb 6916 levels - number of levels of fill. 6917 expected fill - as ratio of original fill. 6918 .ve 6919 6920 Output Parameter: 6921 . fact - the factored matrix 6922 6923 Level: developer 6924 6925 Notes: 6926 Most users should employ the `KSP` interface for linear solvers 6927 instead of working directly with matrix algebra routines such as this. 6928 See, e.g., `KSPCreate()`. 6929 6930 This uses the definition of level of fill as in Y. Saad, 2003 6931 6932 Developer Note: 6933 The Fortran interface is not autogenerated as the 6934 interface definition cannot be generated correctly [due to `MatFactorInfo`] 6935 6936 References: 6937 . * - Y. Saad, Iterative methods for sparse linear systems Philadelphia: Society for Industrial and Applied Mathematics, 2003 6938 6939 .seealso: `Mat`, `MatGetFactor()`, `MatCholeskyFactorNumeric()`, `MatCholeskyFactor()`, `MatFactorInfo` 6940 @*/ 6941 PetscErrorCode MatICCFactorSymbolic(Mat fact, Mat mat, IS perm, const MatFactorInfo *info) 6942 { 6943 PetscFunctionBegin; 6944 PetscValidHeaderSpecific(mat, MAT_CLASSID, 2); 6945 PetscValidType(mat, 2); 6946 if (perm) PetscValidHeaderSpecific(perm, IS_CLASSID, 3); 6947 PetscValidPointer(info, 4); 6948 PetscValidPointer(fact, 1); 6949 PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 6950 PetscCheck(info->levels >= 0, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_OUTOFRANGE, "Levels negative %" PetscInt_FMT, (PetscInt)info->levels); 6951 PetscCheck(info->fill >= 1.0, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_OUTOFRANGE, "Expected fill less than 1.0 %g", (double)info->fill); 6952 if (!(fact)->ops->iccfactorsymbolic) { 6953 MatSolverType stype; 6954 PetscCall(MatFactorGetSolverType(fact, &stype)); 6955 SETERRQ(PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Matrix type %s symbolic ICC using solver type %s", ((PetscObject)mat)->type_name, stype); 6956 } 6957 PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 6958 MatCheckPreallocated(mat, 2); 6959 6960 if (!fact->trivialsymbolic) PetscCall(PetscLogEventBegin(MAT_ICCFactorSymbolic, mat, perm, 0, 0)); 6961 PetscCall((fact->ops->iccfactorsymbolic)(fact, mat, perm, info)); 6962 if (!fact->trivialsymbolic) PetscCall(PetscLogEventEnd(MAT_ICCFactorSymbolic, mat, perm, 0, 0)); 6963 PetscFunctionReturn(PETSC_SUCCESS); 6964 } 6965 6966 /*@C 6967 MatCreateSubMatrices - Extracts several submatrices from a matrix. If submat 6968 points to an array of valid matrices, they may be reused to store the new 6969 submatrices. 6970 6971 Collective 6972 6973 Input Parameters: 6974 + mat - the matrix 6975 . n - the number of submatrixes to be extracted (on this processor, may be zero) 6976 . irow, icol - index sets of rows and columns to extract 6977 - scall - either `MAT_INITIAL_MATRIX` or `MAT_REUSE_MATRIX` 6978 6979 Output Parameter: 6980 . submat - the array of submatrices 6981 6982 Level: advanced 6983 6984 Notes: 6985 `MatCreateSubMatrices()` can extract ONLY sequential submatrices 6986 (from both sequential and parallel matrices). Use `MatCreateSubMatrix()` 6987 to extract a parallel submatrix. 6988 6989 Some matrix types place restrictions on the row and column 6990 indices, such as that they be sorted or that they be equal to each other. 6991 6992 The index sets may not have duplicate entries. 6993 6994 When extracting submatrices from a parallel matrix, each processor can 6995 form a different submatrix by setting the rows and columns of its 6996 individual index sets according to the local submatrix desired. 6997 6998 When finished using the submatrices, the user should destroy 6999 them with `MatDestroySubMatrices()`. 7000 7001 `MAT_REUSE_MATRIX` can only be used when the nonzero structure of the 7002 original matrix has not changed from that last call to `MatCreateSubMatrices()`. 7003 7004 This routine creates the matrices in submat; you should NOT create them before 7005 calling it. It also allocates the array of matrix pointers submat. 7006 7007 For `MATBAIJ` matrices the index sets must respect the block structure, that is if they 7008 request one row/column in a block, they must request all rows/columns that are in 7009 that block. For example, if the block size is 2 you cannot request just row 0 and 7010 column 0. 7011 7012 Fortran Note: 7013 The Fortran interface is slightly different from that given below; it 7014 requires one to pass in as submat a `Mat` (integer) array of size at least n+1. 7015 7016 .seealso: `Mat`, `MatDestroySubMatrices()`, `MatCreateSubMatrix()`, `MatGetRow()`, `MatGetDiagonal()`, `MatReuse` 7017 @*/ 7018 PetscErrorCode MatCreateSubMatrices(Mat mat, PetscInt n, const IS irow[], const IS icol[], MatReuse scall, Mat *submat[]) 7019 { 7020 PetscInt i; 7021 PetscBool eq; 7022 7023 PetscFunctionBegin; 7024 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 7025 PetscValidType(mat, 1); 7026 if (n) { 7027 PetscValidPointer(irow, 3); 7028 for (i = 0; i < n; i++) PetscValidHeaderSpecific(irow[i], IS_CLASSID, 3); 7029 PetscValidPointer(icol, 4); 7030 for (i = 0; i < n; i++) PetscValidHeaderSpecific(icol[i], IS_CLASSID, 4); 7031 } 7032 PetscValidPointer(submat, 6); 7033 if (n && scall == MAT_REUSE_MATRIX) { 7034 PetscValidPointer(*submat, 6); 7035 for (i = 0; i < n; i++) PetscValidHeaderSpecific((*submat)[i], MAT_CLASSID, 6); 7036 } 7037 PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 7038 PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 7039 MatCheckPreallocated(mat, 1); 7040 PetscCall(PetscLogEventBegin(MAT_CreateSubMats, mat, 0, 0, 0)); 7041 PetscUseTypeMethod(mat, createsubmatrices, n, irow, icol, scall, submat); 7042 PetscCall(PetscLogEventEnd(MAT_CreateSubMats, mat, 0, 0, 0)); 7043 for (i = 0; i < n; i++) { 7044 (*submat)[i]->factortype = MAT_FACTOR_NONE; /* in case in place factorization was previously done on submatrix */ 7045 PetscCall(ISEqualUnsorted(irow[i], icol[i], &eq)); 7046 if (eq) PetscCall(MatPropagateSymmetryOptions(mat, (*submat)[i])); 7047 #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) || defined(PETSC_HAVE_HIP) 7048 if (mat->boundtocpu && mat->bindingpropagates) { 7049 PetscCall(MatBindToCPU((*submat)[i], PETSC_TRUE)); 7050 PetscCall(MatSetBindingPropagates((*submat)[i], PETSC_TRUE)); 7051 } 7052 #endif 7053 } 7054 PetscFunctionReturn(PETSC_SUCCESS); 7055 } 7056 7057 /*@C 7058 MatCreateSubMatricesMPI - Extracts MPI submatrices across a sub communicator of mat (by pairs of `IS` that may live on subcomms). 7059 7060 Collective 7061 7062 Input Parameters: 7063 + mat - the matrix 7064 . n - the number of submatrixes to be extracted 7065 . irow, icol - index sets of rows and columns to extract 7066 - scall - either `MAT_INITIAL_MATRIX` or `MAT_REUSE_MATRIX` 7067 7068 Output Parameter: 7069 . submat - the array of submatrices 7070 7071 Level: advanced 7072 7073 Note: 7074 This is used by `PCGASM` 7075 7076 .seealso: `Mat`, `PCGASM`, `MatCreateSubMatrices()`, `MatCreateSubMatrix()`, `MatGetRow()`, `MatGetDiagonal()`, `MatReuse` 7077 @*/ 7078 PetscErrorCode MatCreateSubMatricesMPI(Mat mat, PetscInt n, const IS irow[], const IS icol[], MatReuse scall, Mat *submat[]) 7079 { 7080 PetscInt i; 7081 PetscBool eq; 7082 7083 PetscFunctionBegin; 7084 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 7085 PetscValidType(mat, 1); 7086 if (n) { 7087 PetscValidPointer(irow, 3); 7088 PetscValidHeaderSpecific(*irow, IS_CLASSID, 3); 7089 PetscValidPointer(icol, 4); 7090 PetscValidHeaderSpecific(*icol, IS_CLASSID, 4); 7091 } 7092 PetscValidPointer(submat, 6); 7093 if (n && scall == MAT_REUSE_MATRIX) { 7094 PetscValidPointer(*submat, 6); 7095 PetscValidHeaderSpecific(**submat, MAT_CLASSID, 6); 7096 } 7097 PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 7098 PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 7099 MatCheckPreallocated(mat, 1); 7100 7101 PetscCall(PetscLogEventBegin(MAT_CreateSubMats, mat, 0, 0, 0)); 7102 PetscUseTypeMethod(mat, createsubmatricesmpi, n, irow, icol, scall, submat); 7103 PetscCall(PetscLogEventEnd(MAT_CreateSubMats, mat, 0, 0, 0)); 7104 for (i = 0; i < n; i++) { 7105 PetscCall(ISEqualUnsorted(irow[i], icol[i], &eq)); 7106 if (eq) PetscCall(MatPropagateSymmetryOptions(mat, (*submat)[i])); 7107 } 7108 PetscFunctionReturn(PETSC_SUCCESS); 7109 } 7110 7111 /*@C 7112 MatDestroyMatrices - Destroys an array of matrices. 7113 7114 Collective 7115 7116 Input Parameters: 7117 + n - the number of local matrices 7118 - mat - the matrices (note that this is a pointer to the array of matrices) 7119 7120 Level: advanced 7121 7122 Note: 7123 Frees not only the matrices, but also the array that contains the matrices 7124 7125 Fortran Note: 7126 Will not free the array. 7127 7128 .seealso: `Mat`, `MatCreateSubMatrices()` `MatDestroySubMatrices()` 7129 @*/ 7130 PetscErrorCode MatDestroyMatrices(PetscInt n, Mat *mat[]) 7131 { 7132 PetscInt i; 7133 7134 PetscFunctionBegin; 7135 if (!*mat) PetscFunctionReturn(PETSC_SUCCESS); 7136 PetscCheck(n >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Trying to destroy negative number of matrices %" PetscInt_FMT, n); 7137 PetscValidPointer(mat, 2); 7138 7139 for (i = 0; i < n; i++) PetscCall(MatDestroy(&(*mat)[i])); 7140 7141 /* memory is allocated even if n = 0 */ 7142 PetscCall(PetscFree(*mat)); 7143 PetscFunctionReturn(PETSC_SUCCESS); 7144 } 7145 7146 /*@C 7147 MatDestroySubMatrices - Destroys a set of matrices obtained with `MatCreateSubMatrices()`. 7148 7149 Collective 7150 7151 Input Parameters: 7152 + n - the number of local matrices 7153 - mat - the matrices (note that this is a pointer to the array of matrices, just to match the calling 7154 sequence of MatCreateSubMatrices()) 7155 7156 Level: advanced 7157 7158 Note: 7159 Frees not only the matrices, but also the array that contains the matrices 7160 7161 Fortran Note: 7162 Will not free the array. 7163 7164 .seealso: `Mat`, `MatCreateSubMatrices()`, `MatDestroyMatrices()` 7165 @*/ 7166 PetscErrorCode MatDestroySubMatrices(PetscInt n, Mat *mat[]) 7167 { 7168 Mat mat0; 7169 7170 PetscFunctionBegin; 7171 if (!*mat) PetscFunctionReturn(PETSC_SUCCESS); 7172 /* mat[] is an array of length n+1, see MatCreateSubMatrices_xxx() */ 7173 PetscCheck(n >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Trying to destroy negative number of matrices %" PetscInt_FMT, n); 7174 PetscValidPointer(mat, 2); 7175 7176 mat0 = (*mat)[0]; 7177 if (mat0 && mat0->ops->destroysubmatrices) { 7178 PetscCall((mat0->ops->destroysubmatrices)(n, mat)); 7179 } else { 7180 PetscCall(MatDestroyMatrices(n, mat)); 7181 } 7182 PetscFunctionReturn(PETSC_SUCCESS); 7183 } 7184 7185 /*@C 7186 MatGetSeqNonzeroStructure - Extracts the nonzero structure from a matrix and stores it, in its entirety, on each process 7187 7188 Collective 7189 7190 Input Parameters: 7191 . mat - the matrix 7192 7193 Output Parameter: 7194 . matstruct - the sequential matrix with the nonzero structure of mat 7195 7196 Level: developer 7197 7198 .seealso: `Mat`, `MatDestroySeqNonzeroStructure()`, `MatCreateSubMatrices()`, `MatDestroyMatrices()` 7199 @*/ 7200 PetscErrorCode MatGetSeqNonzeroStructure(Mat mat, Mat *matstruct) 7201 { 7202 PetscFunctionBegin; 7203 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 7204 PetscValidPointer(matstruct, 2); 7205 7206 PetscValidType(mat, 1); 7207 PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 7208 MatCheckPreallocated(mat, 1); 7209 7210 PetscCall(PetscLogEventBegin(MAT_GetSeqNonzeroStructure, mat, 0, 0, 0)); 7211 PetscUseTypeMethod(mat, getseqnonzerostructure, matstruct); 7212 PetscCall(PetscLogEventEnd(MAT_GetSeqNonzeroStructure, mat, 0, 0, 0)); 7213 PetscFunctionReturn(PETSC_SUCCESS); 7214 } 7215 7216 /*@C 7217 MatDestroySeqNonzeroStructure - Destroys matrix obtained with `MatGetSeqNonzeroStructure()`. 7218 7219 Collective 7220 7221 Input Parameters: 7222 . mat - the matrix (note that this is a pointer to the array of matrices, just to match the calling 7223 sequence of `MatGetSequentialNonzeroStructure()`) 7224 7225 Level: advanced 7226 7227 Note: 7228 Frees not only the matrices, but also the array that contains the matrices 7229 7230 .seealso: `Mat`, `MatGetSeqNonzeroStructure()` 7231 @*/ 7232 PetscErrorCode MatDestroySeqNonzeroStructure(Mat *mat) 7233 { 7234 PetscFunctionBegin; 7235 PetscValidPointer(mat, 1); 7236 PetscCall(MatDestroy(mat)); 7237 PetscFunctionReturn(PETSC_SUCCESS); 7238 } 7239 7240 /*@ 7241 MatIncreaseOverlap - Given a set of submatrices indicated by index sets, 7242 replaces the index sets by larger ones that represent submatrices with 7243 additional overlap. 7244 7245 Collective 7246 7247 Input Parameters: 7248 + mat - the matrix 7249 . n - the number of index sets 7250 . is - the array of index sets (these index sets will changed during the call) 7251 - ov - the additional overlap requested 7252 7253 Options Database Key: 7254 . -mat_increase_overlap_scalable - use a scalable algorithm to compute the overlap (supported by MPIAIJ matrix) 7255 7256 Level: developer 7257 7258 Note: 7259 The computed overlap preserves the matrix block sizes when the blocks are square. 7260 That is: if a matrix nonzero for a given block would increase the overlap all columns associated with 7261 that block are included in the overlap regardless of whether each specific column would increase the overlap. 7262 7263 .seealso: `Mat`, `PCASM`, `MatSetBlockSize()`, `MatIncreaseOverlapSplit()`, `MatCreateSubMatrices()` 7264 @*/ 7265 PetscErrorCode MatIncreaseOverlap(Mat mat, PetscInt n, IS is[], PetscInt ov) 7266 { 7267 PetscInt i, bs, cbs; 7268 7269 PetscFunctionBegin; 7270 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 7271 PetscValidType(mat, 1); 7272 PetscValidLogicalCollectiveInt(mat, n, 2); 7273 PetscCheck(n >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must have one or more domains, you have %" PetscInt_FMT, n); 7274 if (n) { 7275 PetscValidPointer(is, 3); 7276 for (i = 0; i < n; i++) PetscValidHeaderSpecific(is[i], IS_CLASSID, 3); 7277 } 7278 PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 7279 PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 7280 MatCheckPreallocated(mat, 1); 7281 7282 if (!ov || !n) PetscFunctionReturn(PETSC_SUCCESS); 7283 PetscCall(PetscLogEventBegin(MAT_IncreaseOverlap, mat, 0, 0, 0)); 7284 PetscUseTypeMethod(mat, increaseoverlap, n, is, ov); 7285 PetscCall(PetscLogEventEnd(MAT_IncreaseOverlap, mat, 0, 0, 0)); 7286 PetscCall(MatGetBlockSizes(mat, &bs, &cbs)); 7287 if (bs == cbs) { 7288 for (i = 0; i < n; i++) PetscCall(ISSetBlockSize(is[i], bs)); 7289 } 7290 PetscFunctionReturn(PETSC_SUCCESS); 7291 } 7292 7293 PetscErrorCode MatIncreaseOverlapSplit_Single(Mat, IS *, PetscInt); 7294 7295 /*@ 7296 MatIncreaseOverlapSplit - Given a set of submatrices indicated by index sets across 7297 a sub communicator, replaces the index sets by larger ones that represent submatrices with 7298 additional overlap. 7299 7300 Collective 7301 7302 Input Parameters: 7303 + mat - the matrix 7304 . n - the number of index sets 7305 . is - the array of index sets (these index sets will changed during the call) 7306 - ov - the additional overlap requested 7307 7308 ` Options Database Key: 7309 . -mat_increase_overlap_scalable - use a scalable algorithm to compute the overlap (supported by MPIAIJ matrix) 7310 7311 Level: developer 7312 7313 .seealso: `Mat`, `MatCreateSubMatrices()`, `MatIncreaseOverlap()` 7314 @*/ 7315 PetscErrorCode MatIncreaseOverlapSplit(Mat mat, PetscInt n, IS is[], PetscInt ov) 7316 { 7317 PetscInt i; 7318 7319 PetscFunctionBegin; 7320 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 7321 PetscValidType(mat, 1); 7322 PetscCheck(n >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must have one or more domains, you have %" PetscInt_FMT, n); 7323 if (n) { 7324 PetscValidPointer(is, 3); 7325 PetscValidHeaderSpecific(*is, IS_CLASSID, 3); 7326 } 7327 PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 7328 PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 7329 MatCheckPreallocated(mat, 1); 7330 if (!ov) PetscFunctionReturn(PETSC_SUCCESS); 7331 PetscCall(PetscLogEventBegin(MAT_IncreaseOverlap, mat, 0, 0, 0)); 7332 for (i = 0; i < n; i++) PetscCall(MatIncreaseOverlapSplit_Single(mat, &is[i], ov)); 7333 PetscCall(PetscLogEventEnd(MAT_IncreaseOverlap, mat, 0, 0, 0)); 7334 PetscFunctionReturn(PETSC_SUCCESS); 7335 } 7336 7337 /*@ 7338 MatGetBlockSize - Returns the matrix block size. 7339 7340 Not Collective 7341 7342 Input Parameter: 7343 . mat - the matrix 7344 7345 Output Parameter: 7346 . bs - block size 7347 7348 Level: intermediate 7349 7350 Notes: 7351 Block row formats are `MATBAIJ` and `MATSBAIJ` ALWAYS have square block storage in the matrix. 7352 7353 If the block size has not been set yet this routine returns 1. 7354 7355 .seealso: `Mat`, `MATBAIJ`, `MATSBAIJ`, `MatCreateSeqBAIJ()`, `MatCreateBAIJ()`, `MatGetBlockSizes()` 7356 @*/ 7357 PetscErrorCode MatGetBlockSize(Mat mat, PetscInt *bs) 7358 { 7359 PetscFunctionBegin; 7360 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 7361 PetscValidIntPointer(bs, 2); 7362 *bs = PetscAbs(mat->rmap->bs); 7363 PetscFunctionReturn(PETSC_SUCCESS); 7364 } 7365 7366 /*@ 7367 MatGetBlockSizes - Returns the matrix block row and column sizes. 7368 7369 Not Collective 7370 7371 Input Parameter: 7372 . mat - the matrix 7373 7374 Output Parameters: 7375 + rbs - row block size 7376 - cbs - column block size 7377 7378 Level: intermediate 7379 7380 Notes: 7381 Block row formats are `MATBAIJ` and `MATSBAIJ` ALWAYS have square block storage in the matrix. 7382 If you pass a different block size for the columns than the rows, the row block size determines the square block storage. 7383 7384 If a block size has not been set yet this routine returns 1. 7385 7386 .seealso: `Mat`, `MATBAIJ`, `MATSBAIJ`, `MatCreateSeqBAIJ()`, `MatCreateBAIJ()`, `MatGetBlockSize()`, `MatSetBlockSize()`, `MatSetBlockSizes()` 7387 @*/ 7388 PetscErrorCode MatGetBlockSizes(Mat mat, PetscInt *rbs, PetscInt *cbs) 7389 { 7390 PetscFunctionBegin; 7391 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 7392 if (rbs) PetscValidIntPointer(rbs, 2); 7393 if (cbs) PetscValidIntPointer(cbs, 3); 7394 if (rbs) *rbs = PetscAbs(mat->rmap->bs); 7395 if (cbs) *cbs = PetscAbs(mat->cmap->bs); 7396 PetscFunctionReturn(PETSC_SUCCESS); 7397 } 7398 7399 /*@ 7400 MatSetBlockSize - Sets the matrix block size. 7401 7402 Logically Collective 7403 7404 Input Parameters: 7405 + mat - the matrix 7406 - bs - block size 7407 7408 Level: intermediate 7409 7410 Notes: 7411 Block row formats are `MATBAIJ` and `MATSBAIJ` formats ALWAYS have square block storage in the matrix. 7412 This must be called before `MatSetUp()` or MatXXXSetPreallocation() (or will default to 1) and the block size cannot be changed later. 7413 7414 For `MATAIJ` matrix format, this function can be called at a later stage, provided that the specified block size 7415 is compatible with the matrix local sizes. 7416 7417 .seealso: `Mat`, `MATBAIJ`, `MATSBAIJ`, `MATAIJ`, `MatCreateSeqBAIJ()`, `MatCreateBAIJ()`, `MatGetBlockSize()`, `MatSetBlockSizes()`, `MatGetBlockSizes()` 7418 @*/ 7419 PetscErrorCode MatSetBlockSize(Mat mat, PetscInt bs) 7420 { 7421 PetscFunctionBegin; 7422 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 7423 PetscValidLogicalCollectiveInt(mat, bs, 2); 7424 PetscCall(MatSetBlockSizes(mat, bs, bs)); 7425 PetscFunctionReturn(PETSC_SUCCESS); 7426 } 7427 7428 typedef struct { 7429 PetscInt n; 7430 IS *is; 7431 Mat *mat; 7432 PetscObjectState nonzerostate; 7433 Mat C; 7434 } EnvelopeData; 7435 7436 static PetscErrorCode EnvelopeDataDestroy(EnvelopeData *edata) 7437 { 7438 for (PetscInt i = 0; i < edata->n; i++) PetscCall(ISDestroy(&edata->is[i])); 7439 PetscCall(PetscFree(edata->is)); 7440 PetscCall(PetscFree(edata)); 7441 return PETSC_SUCCESS; 7442 } 7443 7444 /* 7445 MatComputeVariableBlockEnvelope - Given a matrix whose nonzeros are in blocks along the diagonal this computes and stores 7446 the sizes of these blocks in the matrix. An individual block may lie over several processes. 7447 7448 Collective 7449 7450 Input Parameter: 7451 . mat - the matrix 7452 7453 Notes: 7454 There can be zeros within the blocks 7455 7456 The blocks can overlap between processes, including laying on more than two processes 7457 7458 .seealso: `Mat`, `MatInvertVariableBlockEnvelope()`, `MatSetVariableBlockSizes()` 7459 */ 7460 static PetscErrorCode MatComputeVariableBlockEnvelope(Mat mat) 7461 { 7462 PetscInt n, *sizes, *starts, i = 0, env = 0, tbs = 0, lblocks = 0, rstart, II, ln = 0, cnt = 0, cstart, cend; 7463 PetscInt *diag, *odiag, sc; 7464 VecScatter scatter; 7465 PetscScalar *seqv; 7466 const PetscScalar *parv; 7467 const PetscInt *ia, *ja; 7468 PetscBool set, flag, done; 7469 Mat AA = mat, A; 7470 MPI_Comm comm; 7471 PetscMPIInt rank, size, tag; 7472 MPI_Status status; 7473 PetscContainer container; 7474 EnvelopeData *edata; 7475 Vec seq, par; 7476 IS isglobal; 7477 7478 PetscFunctionBegin; 7479 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 7480 PetscCall(MatIsSymmetricKnown(mat, &set, &flag)); 7481 if (!set || !flag) { 7482 /* TOO: only needs nonzero structure of transpose */ 7483 PetscCall(MatTranspose(mat, MAT_INITIAL_MATRIX, &AA)); 7484 PetscCall(MatAXPY(AA, 1.0, mat, DIFFERENT_NONZERO_PATTERN)); 7485 } 7486 PetscCall(MatAIJGetLocalMat(AA, &A)); 7487 PetscCall(MatGetRowIJ(A, 0, PETSC_FALSE, PETSC_FALSE, &n, &ia, &ja, &done)); 7488 PetscCheck(done, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Unable to get IJ structure from matrix"); 7489 7490 PetscCall(MatGetLocalSize(mat, &n, NULL)); 7491 PetscCall(PetscObjectGetNewTag((PetscObject)mat, &tag)); 7492 PetscCall(PetscObjectGetComm((PetscObject)mat, &comm)); 7493 PetscCallMPI(MPI_Comm_size(comm, &size)); 7494 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 7495 7496 PetscCall(PetscMalloc2(n, &sizes, n, &starts)); 7497 7498 if (rank > 0) { 7499 PetscCallMPI(MPI_Recv(&env, 1, MPIU_INT, rank - 1, tag, comm, &status)); 7500 PetscCallMPI(MPI_Recv(&tbs, 1, MPIU_INT, rank - 1, tag, comm, &status)); 7501 } 7502 PetscCall(MatGetOwnershipRange(mat, &rstart, NULL)); 7503 for (i = 0; i < n; i++) { 7504 env = PetscMax(env, ja[ia[i + 1] - 1]); 7505 II = rstart + i; 7506 if (env == II) { 7507 starts[lblocks] = tbs; 7508 sizes[lblocks++] = 1 + II - tbs; 7509 tbs = 1 + II; 7510 } 7511 } 7512 if (rank < size - 1) { 7513 PetscCallMPI(MPI_Send(&env, 1, MPIU_INT, rank + 1, tag, comm)); 7514 PetscCallMPI(MPI_Send(&tbs, 1, MPIU_INT, rank + 1, tag, comm)); 7515 } 7516 7517 PetscCall(MatRestoreRowIJ(A, 0, PETSC_FALSE, PETSC_FALSE, &n, &ia, &ja, &done)); 7518 if (!set || !flag) PetscCall(MatDestroy(&AA)); 7519 PetscCall(MatDestroy(&A)); 7520 7521 PetscCall(PetscNew(&edata)); 7522 PetscCall(MatGetNonzeroState(mat, &edata->nonzerostate)); 7523 edata->n = lblocks; 7524 /* create IS needed for extracting blocks from the original matrix */ 7525 PetscCall(PetscMalloc1(lblocks, &edata->is)); 7526 for (PetscInt i = 0; i < lblocks; i++) PetscCall(ISCreateStride(PETSC_COMM_SELF, sizes[i], starts[i], 1, &edata->is[i])); 7527 7528 /* Create the resulting inverse matrix structure with preallocation information */ 7529 PetscCall(MatCreate(PetscObjectComm((PetscObject)mat), &edata->C)); 7530 PetscCall(MatSetSizes(edata->C, mat->rmap->n, mat->cmap->n, mat->rmap->N, mat->cmap->N)); 7531 PetscCall(MatSetBlockSizesFromMats(edata->C, mat, mat)); 7532 PetscCall(MatSetType(edata->C, MATAIJ)); 7533 7534 /* Communicate the start and end of each row, from each block to the correct rank */ 7535 /* TODO: Use PetscSF instead of VecScatter */ 7536 for (PetscInt i = 0; i < lblocks; i++) ln += sizes[i]; 7537 PetscCall(VecCreateSeq(PETSC_COMM_SELF, 2 * ln, &seq)); 7538 PetscCall(VecGetArrayWrite(seq, &seqv)); 7539 for (PetscInt i = 0; i < lblocks; i++) { 7540 for (PetscInt j = 0; j < sizes[i]; j++) { 7541 seqv[cnt] = starts[i]; 7542 seqv[cnt + 1] = starts[i] + sizes[i]; 7543 cnt += 2; 7544 } 7545 } 7546 PetscCall(VecRestoreArrayWrite(seq, &seqv)); 7547 PetscCallMPI(MPI_Scan(&cnt, &sc, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)mat))); 7548 sc -= cnt; 7549 PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)mat), 2 * mat->rmap->n, 2 * mat->rmap->N, &par)); 7550 PetscCall(ISCreateStride(PETSC_COMM_SELF, cnt, sc, 1, &isglobal)); 7551 PetscCall(VecScatterCreate(seq, NULL, par, isglobal, &scatter)); 7552 PetscCall(ISDestroy(&isglobal)); 7553 PetscCall(VecScatterBegin(scatter, seq, par, INSERT_VALUES, SCATTER_FORWARD)); 7554 PetscCall(VecScatterEnd(scatter, seq, par, INSERT_VALUES, SCATTER_FORWARD)); 7555 PetscCall(VecScatterDestroy(&scatter)); 7556 PetscCall(VecDestroy(&seq)); 7557 PetscCall(MatGetOwnershipRangeColumn(mat, &cstart, &cend)); 7558 PetscCall(PetscMalloc2(mat->rmap->n, &diag, mat->rmap->n, &odiag)); 7559 PetscCall(VecGetArrayRead(par, &parv)); 7560 cnt = 0; 7561 PetscCall(MatGetSize(mat, NULL, &n)); 7562 for (PetscInt i = 0; i < mat->rmap->n; i++) { 7563 PetscInt start, end, d = 0, od = 0; 7564 7565 start = (PetscInt)PetscRealPart(parv[cnt]); 7566 end = (PetscInt)PetscRealPart(parv[cnt + 1]); 7567 cnt += 2; 7568 7569 if (start < cstart) { 7570 od += cstart - start + n - cend; 7571 d += cend - cstart; 7572 } else if (start < cend) { 7573 od += n - cend; 7574 d += cend - start; 7575 } else od += n - start; 7576 if (end <= cstart) { 7577 od -= cstart - end + n - cend; 7578 d -= cend - cstart; 7579 } else if (end < cend) { 7580 od -= n - cend; 7581 d -= cend - end; 7582 } else od -= n - end; 7583 7584 odiag[i] = od; 7585 diag[i] = d; 7586 } 7587 PetscCall(VecRestoreArrayRead(par, &parv)); 7588 PetscCall(VecDestroy(&par)); 7589 PetscCall(MatXAIJSetPreallocation(edata->C, mat->rmap->bs, diag, odiag, NULL, NULL)); 7590 PetscCall(PetscFree2(diag, odiag)); 7591 PetscCall(PetscFree2(sizes, starts)); 7592 7593 PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &container)); 7594 PetscCall(PetscContainerSetPointer(container, edata)); 7595 PetscCall(PetscContainerSetUserDestroy(container, (PetscErrorCode(*)(void *))EnvelopeDataDestroy)); 7596 PetscCall(PetscObjectCompose((PetscObject)mat, "EnvelopeData", (PetscObject)container)); 7597 PetscCall(PetscObjectDereference((PetscObject)container)); 7598 PetscFunctionReturn(PETSC_SUCCESS); 7599 } 7600 7601 /*@ 7602 MatInvertVariableBlockEnvelope - set matrix C to be the inverted block diagonal of matrix A 7603 7604 Collective 7605 7606 Input Parameters: 7607 . A - the matrix 7608 7609 Output Parameters: 7610 . C - matrix with inverted block diagonal of A. This matrix should be created and may have its type set. 7611 7612 Level: advanced 7613 7614 Note: 7615 For efficiency the matrix A should have all the nonzero entries clustered in smallish blocks along the diagonal. 7616 7617 .seealso: `Mat`, `MatInvertBlockDiagonal()`, `MatComputeBlockDiagonal()` 7618 @*/ 7619 PetscErrorCode MatInvertVariableBlockEnvelope(Mat A, MatReuse reuse, Mat *C) 7620 { 7621 PetscContainer container; 7622 EnvelopeData *edata; 7623 PetscObjectState nonzerostate; 7624 7625 PetscFunctionBegin; 7626 PetscCall(PetscObjectQuery((PetscObject)A, "EnvelopeData", (PetscObject *)&container)); 7627 if (!container) { 7628 PetscCall(MatComputeVariableBlockEnvelope(A)); 7629 PetscCall(PetscObjectQuery((PetscObject)A, "EnvelopeData", (PetscObject *)&container)); 7630 } 7631 PetscCall(PetscContainerGetPointer(container, (void **)&edata)); 7632 PetscCall(MatGetNonzeroState(A, &nonzerostate)); 7633 PetscCheck(nonzerostate <= edata->nonzerostate, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Cannot handle changes to matrix nonzero structure"); 7634 PetscCheck(reuse != MAT_REUSE_MATRIX || *C == edata->C, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "C matrix must be the same as previously output"); 7635 7636 PetscCall(MatCreateSubMatrices(A, edata->n, edata->is, edata->is, MAT_INITIAL_MATRIX, &edata->mat)); 7637 *C = edata->C; 7638 7639 for (PetscInt i = 0; i < edata->n; i++) { 7640 Mat D; 7641 PetscScalar *dvalues; 7642 7643 PetscCall(MatConvert(edata->mat[i], MATSEQDENSE, MAT_INITIAL_MATRIX, &D)); 7644 PetscCall(MatSetOption(*C, MAT_ROW_ORIENTED, PETSC_FALSE)); 7645 PetscCall(MatSeqDenseInvert(D)); 7646 PetscCall(MatDenseGetArray(D, &dvalues)); 7647 PetscCall(MatSetValuesIS(*C, edata->is[i], edata->is[i], dvalues, INSERT_VALUES)); 7648 PetscCall(MatDestroy(&D)); 7649 } 7650 PetscCall(MatDestroySubMatrices(edata->n, &edata->mat)); 7651 PetscCall(MatAssemblyBegin(*C, MAT_FINAL_ASSEMBLY)); 7652 PetscCall(MatAssemblyEnd(*C, MAT_FINAL_ASSEMBLY)); 7653 PetscFunctionReturn(PETSC_SUCCESS); 7654 } 7655 7656 /*@ 7657 MatSetVariableBlockSizes - Sets diagonal point-blocks of the matrix that need not be of the same size 7658 7659 Logically Collective 7660 7661 Input Parameters: 7662 + mat - the matrix 7663 . nblocks - the number of blocks on this process, each block can only exist on a single process 7664 - bsizes - the block sizes 7665 7666 Level: intermediate 7667 7668 Notes: 7669 Currently used by `PCVPBJACOBI` for `MATAIJ` matrices 7670 7671 Each variable point-block set of degrees of freedom must live on a single MPI rank. That is a point block cannot straddle two MPI ranks. 7672 7673 .seealso: `Mat`, `MatCreateSeqBAIJ()`, `MatCreateBAIJ()`, `MatGetBlockSize()`, `MatSetBlockSizes()`, `MatGetBlockSizes()`, `MatGetVariableBlockSizes()`, 7674 `MatComputeVariableBlockEnvelope()`, `PCVPBJACOBI` 7675 @*/ 7676 PetscErrorCode MatSetVariableBlockSizes(Mat mat, PetscInt nblocks, PetscInt *bsizes) 7677 { 7678 PetscInt i, ncnt = 0, nlocal; 7679 7680 PetscFunctionBegin; 7681 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 7682 PetscCheck(nblocks >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Number of local blocks must be great than or equal to zero"); 7683 PetscCall(MatGetLocalSize(mat, &nlocal, NULL)); 7684 for (i = 0; i < nblocks; i++) ncnt += bsizes[i]; 7685 PetscCheck(ncnt == nlocal, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Sum of local block sizes %" PetscInt_FMT " does not equal local size of matrix %" PetscInt_FMT, ncnt, nlocal); 7686 PetscCall(PetscFree(mat->bsizes)); 7687 mat->nblocks = nblocks; 7688 PetscCall(PetscMalloc1(nblocks, &mat->bsizes)); 7689 PetscCall(PetscArraycpy(mat->bsizes, bsizes, nblocks)); 7690 PetscFunctionReturn(PETSC_SUCCESS); 7691 } 7692 7693 /*@C 7694 MatGetVariableBlockSizes - Gets a diagonal blocks of the matrix that need not be of the same size 7695 7696 Logically Collective; No Fortran Support 7697 7698 Input Parameter: 7699 . mat - the matrix 7700 7701 Output Parameters: 7702 + nblocks - the number of blocks on this process 7703 - bsizes - the block sizes 7704 7705 Level: intermediate 7706 7707 .seealso: `Mat`, `MatCreateSeqBAIJ()`, `MatCreateBAIJ()`, `MatGetBlockSize()`, `MatSetBlockSizes()`, `MatGetBlockSizes()`, `MatSetVariableBlockSizes()`, `MatComputeVariableBlockEnvelope()` 7708 @*/ 7709 PetscErrorCode MatGetVariableBlockSizes(Mat mat, PetscInt *nblocks, const PetscInt **bsizes) 7710 { 7711 PetscFunctionBegin; 7712 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 7713 *nblocks = mat->nblocks; 7714 *bsizes = mat->bsizes; 7715 PetscFunctionReturn(PETSC_SUCCESS); 7716 } 7717 7718 /*@ 7719 MatSetBlockSizes - Sets the matrix block row and column sizes. 7720 7721 Logically Collective 7722 7723 Input Parameters: 7724 + mat - the matrix 7725 . rbs - row block size 7726 - cbs - column block size 7727 7728 Level: intermediate 7729 7730 Notes: 7731 Block row formats are `MATBAIJ` and `MATSBAIJ`. These formats ALWAYS have square block storage in the matrix. 7732 If you pass a different block size for the columns than the rows, the row block size determines the square block storage. 7733 This must be called before `MatSetUp()` or MatXXXSetPreallocation() (or will default to 1) and the block size cannot be changed later. 7734 7735 For `MATAIJ` matrix this function can be called at a later stage, provided that the specified block sizes 7736 are compatible with the matrix local sizes. 7737 7738 The row and column block size determine the blocksize of the "row" and "column" vectors returned by `MatCreateVecs()`. 7739 7740 .seealso: `Mat`, `MatCreateSeqBAIJ()`, `MatCreateBAIJ()`, `MatGetBlockSize()`, `MatSetBlockSize()`, `MatGetBlockSizes()` 7741 @*/ 7742 PetscErrorCode MatSetBlockSizes(Mat mat, PetscInt rbs, PetscInt cbs) 7743 { 7744 PetscFunctionBegin; 7745 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 7746 PetscValidLogicalCollectiveInt(mat, rbs, 2); 7747 PetscValidLogicalCollectiveInt(mat, cbs, 3); 7748 PetscTryTypeMethod(mat, setblocksizes, rbs, cbs); 7749 if (mat->rmap->refcnt) { 7750 ISLocalToGlobalMapping l2g = NULL; 7751 PetscLayout nmap = NULL; 7752 7753 PetscCall(PetscLayoutDuplicate(mat->rmap, &nmap)); 7754 if (mat->rmap->mapping) PetscCall(ISLocalToGlobalMappingDuplicate(mat->rmap->mapping, &l2g)); 7755 PetscCall(PetscLayoutDestroy(&mat->rmap)); 7756 mat->rmap = nmap; 7757 mat->rmap->mapping = l2g; 7758 } 7759 if (mat->cmap->refcnt) { 7760 ISLocalToGlobalMapping l2g = NULL; 7761 PetscLayout nmap = NULL; 7762 7763 PetscCall(PetscLayoutDuplicate(mat->cmap, &nmap)); 7764 if (mat->cmap->mapping) PetscCall(ISLocalToGlobalMappingDuplicate(mat->cmap->mapping, &l2g)); 7765 PetscCall(PetscLayoutDestroy(&mat->cmap)); 7766 mat->cmap = nmap; 7767 mat->cmap->mapping = l2g; 7768 } 7769 PetscCall(PetscLayoutSetBlockSize(mat->rmap, rbs)); 7770 PetscCall(PetscLayoutSetBlockSize(mat->cmap, cbs)); 7771 PetscFunctionReturn(PETSC_SUCCESS); 7772 } 7773 7774 /*@ 7775 MatSetBlockSizesFromMats - Sets the matrix block row and column sizes to match a pair of matrices 7776 7777 Logically Collective 7778 7779 Input Parameters: 7780 + mat - the matrix 7781 . fromRow - matrix from which to copy row block size 7782 - fromCol - matrix from which to copy column block size (can be same as fromRow) 7783 7784 Level: developer 7785 7786 .seealso: `Mat`, `MatCreateSeqBAIJ()`, `MatCreateBAIJ()`, `MatGetBlockSize()`, `MatSetBlockSizes()` 7787 @*/ 7788 PetscErrorCode MatSetBlockSizesFromMats(Mat mat, Mat fromRow, Mat fromCol) 7789 { 7790 PetscFunctionBegin; 7791 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 7792 PetscValidHeaderSpecific(fromRow, MAT_CLASSID, 2); 7793 PetscValidHeaderSpecific(fromCol, MAT_CLASSID, 3); 7794 if (fromRow->rmap->bs > 0) PetscCall(PetscLayoutSetBlockSize(mat->rmap, fromRow->rmap->bs)); 7795 if (fromCol->cmap->bs > 0) PetscCall(PetscLayoutSetBlockSize(mat->cmap, fromCol->cmap->bs)); 7796 PetscFunctionReturn(PETSC_SUCCESS); 7797 } 7798 7799 /*@ 7800 MatResidual - Default routine to calculate the residual r = b - Ax 7801 7802 Collective 7803 7804 Input Parameters: 7805 + mat - the matrix 7806 . b - the right-hand-side 7807 - x - the approximate solution 7808 7809 Output Parameter: 7810 . r - location to store the residual 7811 7812 Level: developer 7813 7814 .seealso: `Mat`, `MatMult()`, `MatMultAdd()`, `PCMGSetResidual()` 7815 @*/ 7816 PetscErrorCode MatResidual(Mat mat, Vec b, Vec x, Vec r) 7817 { 7818 PetscFunctionBegin; 7819 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 7820 PetscValidHeaderSpecific(b, VEC_CLASSID, 2); 7821 PetscValidHeaderSpecific(x, VEC_CLASSID, 3); 7822 PetscValidHeaderSpecific(r, VEC_CLASSID, 4); 7823 PetscValidType(mat, 1); 7824 MatCheckPreallocated(mat, 1); 7825 PetscCall(PetscLogEventBegin(MAT_Residual, mat, 0, 0, 0)); 7826 if (!mat->ops->residual) { 7827 PetscCall(MatMult(mat, x, r)); 7828 PetscCall(VecAYPX(r, -1.0, b)); 7829 } else { 7830 PetscUseTypeMethod(mat, residual, b, x, r); 7831 } 7832 PetscCall(PetscLogEventEnd(MAT_Residual, mat, 0, 0, 0)); 7833 PetscFunctionReturn(PETSC_SUCCESS); 7834 } 7835 7836 /*MC 7837 MatGetRowIJF90 - Obtains the compressed row storage i and j indices for the local rows of a sparse matrix 7838 7839 Synopsis: 7840 MatGetRowIJF90(Mat A, PetscInt shift, PetscBool symmetric, PetscBool inodecompressed, PetscInt n, {PetscInt, pointer :: ia(:)}, {PetscInt, pointer :: ja(:)}, PetscBool done,integer ierr) 7841 7842 Not Collective 7843 7844 Input Parameters: 7845 + A - the matrix 7846 . shift - 0 or 1 indicating we want the indices starting at 0 or 1 7847 . symmetric - `PETSC_TRUE` or `PETSC_FALSE` indicating the matrix data structure should be symmetrized 7848 - inodecompressed - `PETSC_TRUE` or `PETSC_FALSE` indicats if the nonzero structure of the 7849 inodes or the nonzero elements is wanted. For `MATBAIJ` matrices the compressed version is 7850 always used. 7851 7852 Output Parameters: 7853 + n - number of local rows in the (possibly compressed) matrix 7854 . ia - the row pointers; that is ia[0] = 0, ia[row] = ia[row-1] + number of elements in that row of the matrix 7855 . ja - the column indices 7856 - done - indicates if the routine actually worked and returned appropriate ia[] and ja[] arrays; callers 7857 are responsible for handling the case when done == `PETSC_FALSE` and ia and ja are not set 7858 7859 Level: developer 7860 7861 Note: 7862 Use `MatRestoreRowIJF90()` when you no longer need access to the data 7863 7864 .seealso: [](sec_fortranarrays), `Mat`, `MATMPIAIJ`, `MatGetRowIJ()`, `MatRestoreRowIJ()`, `MatRestoreRowIJF90()` 7865 M*/ 7866 7867 /*MC 7868 MatRestoreRowIJF90 - restores the compressed row storage i and j indices for the local rows of a sparse matrix obtained with `MatGetRowIJF90()` 7869 7870 Synopsis: 7871 MatRestoreRowIJF90(Mat A, PetscInt shift, PetscBool symmetric, PetscBool inodecompressed, PetscInt n, {PetscInt, pointer :: ia(:)}, {PetscInt, pointer :: ja(:)}, PetscBool done,integer ierr) 7872 7873 Not Collective 7874 7875 Input Parameters: 7876 + A - the matrix 7877 . shift - 0 or 1 indicating we want the indices starting at 0 or 1 7878 . symmetric - `PETSC_TRUE` or `PETSC_FALSE` indicating the matrix data structure should be symmetrized 7879 inodecompressed - `PETSC_TRUE` or `PETSC_FALSE` indicats if the nonzero structure of the 7880 inodes or the nonzero elements is wanted. For `MATBAIJ` matrices the compressed version is 7881 always used. 7882 . n - number of local rows in the (possibly compressed) matrix 7883 . ia - the row pointers; that is ia[0] = 0, ia[row] = ia[row-1] + number of elements in that row of the matrix 7884 . ja - the column indices 7885 - done - indicates if the routine actually worked and returned appropriate ia[] and ja[] arrays; callers 7886 are responsible for handling the case when done == `PETSC_FALSE` and ia and ja are not set 7887 7888 Level: developer 7889 7890 .seealso: [](sec_fortranarrays), `Mat`, `MATMPIAIJ`, `MatGetRowIJ()`, `MatRestoreRowIJ()`, `MatGetRowIJF90()` 7891 M*/ 7892 7893 /*@C 7894 MatGetRowIJ - Returns the compressed row storage i and j indices for the local rows of a sparse matrix 7895 7896 Collective 7897 7898 Input Parameters: 7899 + mat - the matrix 7900 . shift - 0 or 1 indicating we want the indices starting at 0 or 1 7901 . symmetric - `PETSC_TRUE` or `PETSC_FALSE` indicating the matrix data structure should be symmetrized 7902 - inodecompressed - `PETSC_TRUE` or `PETSC_FALSE` indicats if the nonzero structure of the 7903 inodes or the nonzero elements is wanted. For `MATBAIJ` matrices the compressed version is 7904 always used. 7905 7906 Output Parameters: 7907 + n - number of local rows in the (possibly compressed) matrix, use `NULL` if not needed 7908 . ia - the row pointers; that is ia[0] = 0, ia[row] = ia[row-1] + number of elements in that row of the matrix, use `NULL` if not needed 7909 . ja - the column indices, use `NULL` if not needed 7910 - done - indicates if the routine actually worked and returned appropriate ia[] and ja[] arrays; callers 7911 are responsible for handling the case when done == `PETSC_FALSE` and ia and ja are not set 7912 7913 Level: developer 7914 7915 Notes: 7916 You CANNOT change any of the ia[] or ja[] values. 7917 7918 Use `MatRestoreRowIJ()` when you are finished accessing the ia[] and ja[] values. 7919 7920 Fortran Notes: 7921 In Fortran use 7922 .vb 7923 PetscInt, pointer :: ia(:),ja(:) 7924 call MatGetRowIJF90(mat,shift,symmetric,inodecompressed,n,ia,ja,done,ierr) 7925 ! Access the ith and jth entries via ia(i) and ja(j) 7926 .ve 7927 `MatGetRowIJ()` Fortran binding is deprecated (since PETSc 3.19), use `MatGetRowIJF90()` 7928 7929 .seealso: `Mat`, `MATAIJ`, `MatGetRowIJF90()`, `MatGetColumnIJ()`, `MatRestoreRowIJ()`, `MatSeqAIJGetArray()` 7930 @*/ 7931 PetscErrorCode MatGetRowIJ(Mat mat, PetscInt shift, PetscBool symmetric, PetscBool inodecompressed, PetscInt *n, const PetscInt *ia[], const PetscInt *ja[], PetscBool *done) 7932 { 7933 PetscFunctionBegin; 7934 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 7935 PetscValidType(mat, 1); 7936 if (n) PetscValidIntPointer(n, 5); 7937 if (ia) PetscValidPointer(ia, 6); 7938 if (ja) PetscValidPointer(ja, 7); 7939 if (done) PetscValidBoolPointer(done, 8); 7940 MatCheckPreallocated(mat, 1); 7941 if (!mat->ops->getrowij && done) *done = PETSC_FALSE; 7942 else { 7943 if (done) *done = PETSC_TRUE; 7944 PetscCall(PetscLogEventBegin(MAT_GetRowIJ, mat, 0, 0, 0)); 7945 PetscUseTypeMethod(mat, getrowij, shift, symmetric, inodecompressed, n, ia, ja, done); 7946 PetscCall(PetscLogEventEnd(MAT_GetRowIJ, mat, 0, 0, 0)); 7947 } 7948 PetscFunctionReturn(PETSC_SUCCESS); 7949 } 7950 7951 /*@C 7952 MatGetColumnIJ - Returns the compressed column storage i and j indices for sequential matrices. 7953 7954 Collective 7955 7956 Input Parameters: 7957 + mat - the matrix 7958 . shift - 1 or zero indicating we want the indices starting at 0 or 1 7959 . symmetric - `PETSC_TRUE` or `PETSC_FALSE` indicating the matrix data structure should be 7960 symmetrized 7961 . inodecompressed - `PETSC_TRUE` or `PETSC_FALSE` indicating if the nonzero structure of the 7962 inodes or the nonzero elements is wanted. For `MATBAIJ` matrices the compressed version is 7963 always used. 7964 . n - number of columns in the (possibly compressed) matrix 7965 . ia - the column pointers; that is ia[0] = 0, ia[col] = i[col-1] + number of elements in that col of the matrix 7966 - ja - the row indices 7967 7968 Output Parameters: 7969 . done - `PETSC_TRUE` or `PETSC_FALSE`, indicating whether the values have been returned 7970 7971 Level: developer 7972 7973 .seealso: `Mat`, `MatGetRowIJ()`, `MatRestoreColumnIJ()` 7974 @*/ 7975 PetscErrorCode MatGetColumnIJ(Mat mat, PetscInt shift, PetscBool symmetric, PetscBool inodecompressed, PetscInt *n, const PetscInt *ia[], const PetscInt *ja[], PetscBool *done) 7976 { 7977 PetscFunctionBegin; 7978 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 7979 PetscValidType(mat, 1); 7980 PetscValidIntPointer(n, 5); 7981 if (ia) PetscValidPointer(ia, 6); 7982 if (ja) PetscValidPointer(ja, 7); 7983 PetscValidBoolPointer(done, 8); 7984 MatCheckPreallocated(mat, 1); 7985 if (!mat->ops->getcolumnij) *done = PETSC_FALSE; 7986 else { 7987 *done = PETSC_TRUE; 7988 PetscUseTypeMethod(mat, getcolumnij, shift, symmetric, inodecompressed, n, ia, ja, done); 7989 } 7990 PetscFunctionReturn(PETSC_SUCCESS); 7991 } 7992 7993 /*@C 7994 MatRestoreRowIJ - Call after you are completed with the ia,ja indices obtained with `MatGetRowIJ()`. 7995 7996 Collective 7997 7998 Input Parameters: 7999 + mat - the matrix 8000 . shift - 1 or zero indicating we want the indices starting at 0 or 1 8001 . symmetric - `PETSC_TRUE` or `PETSC_FALSE` indicating the matrix data structure should be symmetrized 8002 . inodecompressed - `PETSC_TRUE` or `PETSC_FALSE` indicating if the nonzero structure of the 8003 inodes or the nonzero elements is wanted. For `MATBAIJ` matrices the compressed version is 8004 always used. 8005 . n - size of (possibly compressed) matrix 8006 . ia - the row pointers 8007 - ja - the column indices 8008 8009 Output Parameters: 8010 . done - `PETSC_TRUE` or `PETSC_FALSE` indicated that the values have been returned 8011 8012 Level: developer 8013 8014 Note: 8015 This routine zeros out `n`, `ia`, and `ja`. This is to prevent accidental 8016 us of the array after it has been restored. If you pass `NULL`, it will 8017 not zero the pointers. Use of ia or ja after `MatRestoreRowIJ()` is invalid. 8018 8019 Fortran Note: 8020 `MatRestoreRowIJ()` Fortran binding is deprecated (since PETSc 3.19), use `MatRestoreRowIJF90()` 8021 .seealso: `Mat`, `MatGetRowIJ()`, `MatRestoreRowIJF90()`, `MatRestoreColumnIJ()` 8022 @*/ 8023 PetscErrorCode MatRestoreRowIJ(Mat mat, PetscInt shift, PetscBool symmetric, PetscBool inodecompressed, PetscInt *n, const PetscInt *ia[], const PetscInt *ja[], PetscBool *done) 8024 { 8025 PetscFunctionBegin; 8026 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 8027 PetscValidType(mat, 1); 8028 if (ia) PetscValidPointer(ia, 6); 8029 if (ja) PetscValidPointer(ja, 7); 8030 if (done) PetscValidBoolPointer(done, 8); 8031 MatCheckPreallocated(mat, 1); 8032 8033 if (!mat->ops->restorerowij && done) *done = PETSC_FALSE; 8034 else { 8035 if (done) *done = PETSC_TRUE; 8036 PetscUseTypeMethod(mat, restorerowij, shift, symmetric, inodecompressed, n, ia, ja, done); 8037 if (n) *n = 0; 8038 if (ia) *ia = NULL; 8039 if (ja) *ja = NULL; 8040 } 8041 PetscFunctionReturn(PETSC_SUCCESS); 8042 } 8043 8044 /*@C 8045 MatRestoreColumnIJ - Call after you are completed with the ia,ja indices obtained with `MatGetColumnIJ()`. 8046 8047 Collective on Mat 8048 8049 Input Parameters: 8050 + mat - the matrix 8051 . shift - 1 or zero indicating we want the indices starting at 0 or 1 8052 . symmetric - `PETSC_TRUE` or `PETSC_FALSE` indicating the matrix data structure should be symmetrized 8053 - inodecompressed - `PETSC_TRUE` or `PETSC_FALSE` indicating if the nonzero structure of the 8054 inodes or the nonzero elements is wanted. For `MATBAIJ` matrices the compressed version is 8055 always used. 8056 8057 Output Parameters: 8058 + n - size of (possibly compressed) matrix 8059 . ia - the column pointers 8060 . ja - the row indices 8061 - done - `PETSC_TRUE` or `PETSC_FALSE` indicated that the values have been returned 8062 8063 Level: developer 8064 8065 .seealso: `Mat`, `MatGetColumnIJ()`, `MatRestoreRowIJ()` 8066 @*/ 8067 PetscErrorCode MatRestoreColumnIJ(Mat mat, PetscInt shift, PetscBool symmetric, PetscBool inodecompressed, PetscInt *n, const PetscInt *ia[], const PetscInt *ja[], PetscBool *done) 8068 { 8069 PetscFunctionBegin; 8070 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 8071 PetscValidType(mat, 1); 8072 if (ia) PetscValidPointer(ia, 6); 8073 if (ja) PetscValidPointer(ja, 7); 8074 PetscValidBoolPointer(done, 8); 8075 MatCheckPreallocated(mat, 1); 8076 8077 if (!mat->ops->restorecolumnij) *done = PETSC_FALSE; 8078 else { 8079 *done = PETSC_TRUE; 8080 PetscUseTypeMethod(mat, restorecolumnij, shift, symmetric, inodecompressed, n, ia, ja, done); 8081 if (n) *n = 0; 8082 if (ia) *ia = NULL; 8083 if (ja) *ja = NULL; 8084 } 8085 PetscFunctionReturn(PETSC_SUCCESS); 8086 } 8087 8088 /*@C 8089 MatColoringPatch -Used inside matrix coloring routines that use `MatGetRowIJ()` and/or `MatGetColumnIJ()`. 8090 8091 Collective 8092 8093 Input Parameters: 8094 + mat - the matrix 8095 . ncolors - max color value 8096 . n - number of entries in colorarray 8097 - colorarray - array indicating color for each column 8098 8099 Output Parameters: 8100 . iscoloring - coloring generated using colorarray information 8101 8102 Level: developer 8103 8104 .seealso: `Mat`, `MatGetRowIJ()`, `MatGetColumnIJ()` 8105 @*/ 8106 PetscErrorCode MatColoringPatch(Mat mat, PetscInt ncolors, PetscInt n, ISColoringValue colorarray[], ISColoring *iscoloring) 8107 { 8108 PetscFunctionBegin; 8109 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 8110 PetscValidType(mat, 1); 8111 PetscValidIntPointer(colorarray, 4); 8112 PetscValidPointer(iscoloring, 5); 8113 MatCheckPreallocated(mat, 1); 8114 8115 if (!mat->ops->coloringpatch) { 8116 PetscCall(ISColoringCreate(PetscObjectComm((PetscObject)mat), ncolors, n, colorarray, PETSC_OWN_POINTER, iscoloring)); 8117 } else { 8118 PetscUseTypeMethod(mat, coloringpatch, ncolors, n, colorarray, iscoloring); 8119 } 8120 PetscFunctionReturn(PETSC_SUCCESS); 8121 } 8122 8123 /*@ 8124 MatSetUnfactored - Resets a factored matrix to be treated as unfactored. 8125 8126 Logically Collective 8127 8128 Input Parameter: 8129 . mat - the factored matrix to be reset 8130 8131 Level: developer 8132 8133 Notes: 8134 This routine should be used only with factored matrices formed by in-place 8135 factorization via ILU(0) (or by in-place LU factorization for the `MATSEQDENSE` 8136 format). This option can save memory, for example, when solving nonlinear 8137 systems with a matrix-free Newton-Krylov method and a matrix-based, in-place 8138 ILU(0) preconditioner. 8139 8140 Note that one can specify in-place ILU(0) factorization by calling 8141 .vb 8142 PCType(pc,PCILU); 8143 PCFactorSeUseInPlace(pc); 8144 .ve 8145 or by using the options -pc_type ilu -pc_factor_in_place 8146 8147 In-place factorization ILU(0) can also be used as a local 8148 solver for the blocks within the block Jacobi or additive Schwarz 8149 methods (runtime option: -sub_pc_factor_in_place). See Users-Manual: ch_pc 8150 for details on setting local solver options. 8151 8152 Most users should employ the `KSP` interface for linear solvers 8153 instead of working directly with matrix algebra routines such as this. 8154 See, e.g., `KSPCreate()`. 8155 8156 .seealso: `Mat`, `PCFactorSetUseInPlace()`, `PCFactorGetUseInPlace()` 8157 @*/ 8158 PetscErrorCode MatSetUnfactored(Mat mat) 8159 { 8160 PetscFunctionBegin; 8161 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 8162 PetscValidType(mat, 1); 8163 MatCheckPreallocated(mat, 1); 8164 mat->factortype = MAT_FACTOR_NONE; 8165 if (!mat->ops->setunfactored) PetscFunctionReturn(PETSC_SUCCESS); 8166 PetscUseTypeMethod(mat, setunfactored); 8167 PetscFunctionReturn(PETSC_SUCCESS); 8168 } 8169 8170 /*MC 8171 MatDenseGetArrayF90 - Accesses a matrix array from Fortran 8172 8173 Synopsis: 8174 MatDenseGetArrayF90(Mat x,{Scalar, pointer :: xx_v(:,:)},integer ierr) 8175 8176 Not collective 8177 8178 Input Parameter: 8179 . x - matrix 8180 8181 Output Parameters: 8182 + xx_v - the Fortran pointer to the array 8183 - ierr - error code 8184 8185 Example of Usage: 8186 .vb 8187 PetscScalar, pointer xx_v(:,:) 8188 .... 8189 call MatDenseGetArrayF90(x,xx_v,ierr) 8190 a = xx_v(3) 8191 call MatDenseRestoreArrayF90(x,xx_v,ierr) 8192 .ve 8193 8194 Level: advanced 8195 8196 .seealso: `Mat`, `MatDenseRestoreArrayF90()`, `MatDenseGetArray()`, `MatDenseRestoreArray()`, `MatSeqAIJGetArrayF90()` 8197 M*/ 8198 8199 /*MC 8200 MatDenseRestoreArrayF90 - Restores a matrix array that has been 8201 accessed with `MatDenseGetArrayF90()`. 8202 8203 Synopsis: 8204 MatDenseRestoreArrayF90(Mat x,{Scalar, pointer :: xx_v(:,:)},integer ierr) 8205 8206 Not collective 8207 8208 Input Parameters: 8209 + x - matrix 8210 - xx_v - the Fortran90 pointer to the array 8211 8212 Output Parameter: 8213 . ierr - error code 8214 8215 Example of Usage: 8216 .vb 8217 PetscScalar, pointer xx_v(:,:) 8218 .... 8219 call MatDenseGetArrayF90(x,xx_v,ierr) 8220 a = xx_v(3) 8221 call MatDenseRestoreArrayF90(x,xx_v,ierr) 8222 .ve 8223 8224 Level: advanced 8225 8226 .seealso: `Mat`, `MatDenseGetArrayF90()`, `MatDenseGetArray()`, `MatDenseRestoreArray()`, `MatSeqAIJRestoreArrayF90()` 8227 M*/ 8228 8229 /*MC 8230 MatSeqAIJGetArrayF90 - Accesses a matrix array from Fortran. 8231 8232 Synopsis: 8233 MatSeqAIJGetArrayF90(Mat x,{Scalar, pointer :: xx_v(:)},integer ierr) 8234 8235 Not collective 8236 8237 Input Parameter: 8238 . x - matrix 8239 8240 Output Parameters: 8241 + xx_v - the Fortran pointer to the array 8242 - ierr - error code 8243 8244 Example of Usage: 8245 .vb 8246 PetscScalar, pointer xx_v(:) 8247 .... 8248 call MatSeqAIJGetArrayF90(x,xx_v,ierr) 8249 a = xx_v(3) 8250 call MatSeqAIJRestoreArrayF90(x,xx_v,ierr) 8251 .ve 8252 8253 Level: advanced 8254 8255 .seealso: `Mat`, `MatSeqAIJRestoreArrayF90()`, `MatSeqAIJGetArray()`, `MatSeqAIJRestoreArray()`, `MatDenseGetArrayF90()` 8256 M*/ 8257 8258 /*MC 8259 MatSeqAIJRestoreArrayF90 - Restores a matrix array that has been 8260 accessed with `MatSeqAIJGetArrayF90()`. 8261 8262 Synopsis: 8263 MatSeqAIJRestoreArrayF90(Mat x,{Scalar, pointer :: xx_v(:)},integer ierr) 8264 8265 Not collective 8266 8267 Input Parameters: 8268 + x - matrix 8269 - xx_v - the Fortran90 pointer to the array 8270 8271 Output Parameter: 8272 . ierr - error code 8273 8274 Example of Usage: 8275 .vb 8276 PetscScalar, pointer xx_v(:) 8277 .... 8278 call MatSeqAIJGetArrayF90(x,xx_v,ierr) 8279 a = xx_v(3) 8280 call MatSeqAIJRestoreArrayF90(x,xx_v,ierr) 8281 .ve 8282 8283 Level: advanced 8284 8285 .seealso: `Mat`, `MatSeqAIJGetArrayF90()`, `MatSeqAIJGetArray()`, `MatSeqAIJRestoreArray()`, `MatDenseRestoreArrayF90()` 8286 M*/ 8287 8288 /*@ 8289 MatCreateSubMatrix - Gets a single submatrix on the same number of processors 8290 as the original matrix. 8291 8292 Collective 8293 8294 Input Parameters: 8295 + mat - the original matrix 8296 . isrow - parallel IS containing the rows this processor should obtain 8297 . iscol - parallel IS containing all columns you wish to keep. Each process should list the columns that will be in IT's "diagonal part" in the new matrix. 8298 - cll - either `MAT_INITIAL_MATRIX` or `MAT_REUSE_MATRIX` 8299 8300 Output Parameter: 8301 . newmat - the new submatrix, of the same type as the old 8302 8303 Level: advanced 8304 8305 Notes: 8306 The submatrix will be able to be multiplied with vectors using the same layout as iscol. 8307 8308 Some matrix types place restrictions on the row and column indices, such 8309 as that they be sorted or that they be equal to each other. For `MATBAIJ` and `MATSBAIJ` matrices the indices must include all rows/columns of a block; 8310 for example, if the block size is 3 one cannot select the 0 and 2 rows without selecting the 1 row. 8311 8312 The index sets may not have duplicate entries. 8313 8314 The first time this is called you should use a cll of `MAT_INITIAL_MATRIX`, 8315 the `MatCreateSubMatrix()` routine will create the newmat for you. Any additional calls 8316 to this routine with a mat of the same nonzero structure and with a call of `MAT_REUSE_MATRIX` 8317 will reuse the matrix generated the first time. You should call `MatDestroy()` on newmat when 8318 you are finished using it. 8319 8320 The communicator of the newly obtained matrix is ALWAYS the same as the communicator of 8321 the input matrix. 8322 8323 If iscol is `NULL` then all columns are obtained (not supported in Fortran). 8324 8325 Example usage: 8326 Consider the following 8x8 matrix with 34 non-zero values, that is 8327 assembled across 3 processors. Let's assume that proc0 owns 3 rows, 8328 proc1 owns 3 rows, proc2 owns 2 rows. This division can be shown 8329 as follows: 8330 8331 .vb 8332 1 2 0 | 0 3 0 | 0 4 8333 Proc0 0 5 6 | 7 0 0 | 8 0 8334 9 0 10 | 11 0 0 | 12 0 8335 ------------------------------------- 8336 13 0 14 | 15 16 17 | 0 0 8337 Proc1 0 18 0 | 19 20 21 | 0 0 8338 0 0 0 | 22 23 0 | 24 0 8339 ------------------------------------- 8340 Proc2 25 26 27 | 0 0 28 | 29 0 8341 30 0 0 | 31 32 33 | 0 34 8342 .ve 8343 8344 Suppose isrow = [0 1 | 4 | 6 7] and iscol = [1 2 | 3 4 5 | 6]. The resulting submatrix is 8345 8346 .vb 8347 2 0 | 0 3 0 | 0 8348 Proc0 5 6 | 7 0 0 | 8 8349 ------------------------------- 8350 Proc1 18 0 | 19 20 21 | 0 8351 ------------------------------- 8352 Proc2 26 27 | 0 0 28 | 29 8353 0 0 | 31 32 33 | 0 8354 .ve 8355 8356 .seealso: `Mat`, `MatCreateSubMatrices()`, `MatCreateSubMatricesMPI()`, `MatCreateSubMatrixVirtual()`, `MatSubMatrixVirtualUpdate()` 8357 @*/ 8358 PetscErrorCode MatCreateSubMatrix(Mat mat, IS isrow, IS iscol, MatReuse cll, Mat *newmat) 8359 { 8360 PetscMPIInt size; 8361 Mat *local; 8362 IS iscoltmp; 8363 PetscBool flg; 8364 8365 PetscFunctionBegin; 8366 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 8367 PetscValidHeaderSpecific(isrow, IS_CLASSID, 2); 8368 if (iscol) PetscValidHeaderSpecific(iscol, IS_CLASSID, 3); 8369 PetscValidPointer(newmat, 5); 8370 if (cll == MAT_REUSE_MATRIX) PetscValidHeaderSpecific(*newmat, MAT_CLASSID, 5); 8371 PetscValidType(mat, 1); 8372 PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 8373 PetscCheck(cll != MAT_IGNORE_MATRIX, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Cannot use MAT_IGNORE_MATRIX"); 8374 8375 MatCheckPreallocated(mat, 1); 8376 PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)mat), &size)); 8377 8378 if (!iscol || isrow == iscol) { 8379 PetscBool stride; 8380 PetscMPIInt grabentirematrix = 0, grab; 8381 PetscCall(PetscObjectTypeCompare((PetscObject)isrow, ISSTRIDE, &stride)); 8382 if (stride) { 8383 PetscInt first, step, n, rstart, rend; 8384 PetscCall(ISStrideGetInfo(isrow, &first, &step)); 8385 if (step == 1) { 8386 PetscCall(MatGetOwnershipRange(mat, &rstart, &rend)); 8387 if (rstart == first) { 8388 PetscCall(ISGetLocalSize(isrow, &n)); 8389 if (n == rend - rstart) grabentirematrix = 1; 8390 } 8391 } 8392 } 8393 PetscCall(MPIU_Allreduce(&grabentirematrix, &grab, 1, MPI_INT, MPI_MIN, PetscObjectComm((PetscObject)mat))); 8394 if (grab) { 8395 PetscCall(PetscInfo(mat, "Getting entire matrix as submatrix\n")); 8396 if (cll == MAT_INITIAL_MATRIX) { 8397 *newmat = mat; 8398 PetscCall(PetscObjectReference((PetscObject)mat)); 8399 } 8400 PetscFunctionReturn(PETSC_SUCCESS); 8401 } 8402 } 8403 8404 if (!iscol) { 8405 PetscCall(ISCreateStride(PetscObjectComm((PetscObject)mat), mat->cmap->n, mat->cmap->rstart, 1, &iscoltmp)); 8406 } else { 8407 iscoltmp = iscol; 8408 } 8409 8410 /* if original matrix is on just one processor then use submatrix generated */ 8411 if (mat->ops->createsubmatrices && !mat->ops->createsubmatrix && size == 1 && cll == MAT_REUSE_MATRIX) { 8412 PetscCall(MatCreateSubMatrices(mat, 1, &isrow, &iscoltmp, MAT_REUSE_MATRIX, &newmat)); 8413 goto setproperties; 8414 } else if (mat->ops->createsubmatrices && !mat->ops->createsubmatrix && size == 1) { 8415 PetscCall(MatCreateSubMatrices(mat, 1, &isrow, &iscoltmp, MAT_INITIAL_MATRIX, &local)); 8416 *newmat = *local; 8417 PetscCall(PetscFree(local)); 8418 goto setproperties; 8419 } else if (!mat->ops->createsubmatrix) { 8420 /* Create a new matrix type that implements the operation using the full matrix */ 8421 PetscCall(PetscLogEventBegin(MAT_CreateSubMat, mat, 0, 0, 0)); 8422 switch (cll) { 8423 case MAT_INITIAL_MATRIX: 8424 PetscCall(MatCreateSubMatrixVirtual(mat, isrow, iscoltmp, newmat)); 8425 break; 8426 case MAT_REUSE_MATRIX: 8427 PetscCall(MatSubMatrixVirtualUpdate(*newmat, mat, isrow, iscoltmp)); 8428 break; 8429 default: 8430 SETERRQ(PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_OUTOFRANGE, "Invalid MatReuse, must be either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX"); 8431 } 8432 PetscCall(PetscLogEventEnd(MAT_CreateSubMat, mat, 0, 0, 0)); 8433 goto setproperties; 8434 } 8435 8436 PetscCall(PetscLogEventBegin(MAT_CreateSubMat, mat, 0, 0, 0)); 8437 PetscUseTypeMethod(mat, createsubmatrix, isrow, iscoltmp, cll, newmat); 8438 PetscCall(PetscLogEventEnd(MAT_CreateSubMat, mat, 0, 0, 0)); 8439 8440 setproperties: 8441 PetscCall(ISEqualUnsorted(isrow, iscoltmp, &flg)); 8442 if (flg) PetscCall(MatPropagateSymmetryOptions(mat, *newmat)); 8443 if (!iscol) PetscCall(ISDestroy(&iscoltmp)); 8444 if (*newmat && cll == MAT_INITIAL_MATRIX) PetscCall(PetscObjectStateIncrease((PetscObject)*newmat)); 8445 PetscFunctionReturn(PETSC_SUCCESS); 8446 } 8447 8448 /*@ 8449 MatPropagateSymmetryOptions - Propagates symmetry options set on a matrix to another matrix 8450 8451 Not Collective 8452 8453 Input Parameters: 8454 + A - the matrix we wish to propagate options from 8455 - B - the matrix we wish to propagate options to 8456 8457 Level: beginner 8458 8459 Note: 8460 Propagates the options associated to `MAT_SYMMETRY_ETERNAL`, `MAT_STRUCTURALLY_SYMMETRIC`, `MAT_HERMITIAN`, `MAT_SPD`, `MAT_SYMMETRIC`, and `MAT_STRUCTURAL_SYMMETRY_ETERNAL` 8461 8462 .seealso: `Mat`, `MatSetOption()`, `MatIsSymmetricKnown()`, `MatIsSPDKnown()`, `MatIsHermitianKnown()`, MatIsStructurallySymmetricKnown()` 8463 @*/ 8464 PetscErrorCode MatPropagateSymmetryOptions(Mat A, Mat B) 8465 { 8466 PetscFunctionBegin; 8467 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 8468 PetscValidHeaderSpecific(B, MAT_CLASSID, 2); 8469 B->symmetry_eternal = A->symmetry_eternal; 8470 B->structural_symmetry_eternal = A->structural_symmetry_eternal; 8471 B->symmetric = A->symmetric; 8472 B->structurally_symmetric = A->structurally_symmetric; 8473 B->spd = A->spd; 8474 B->hermitian = A->hermitian; 8475 PetscFunctionReturn(PETSC_SUCCESS); 8476 } 8477 8478 /*@ 8479 MatStashSetInitialSize - sets the sizes of the matrix stash, that is 8480 used during the assembly process to store values that belong to 8481 other processors. 8482 8483 Not Collective 8484 8485 Input Parameters: 8486 + mat - the matrix 8487 . size - the initial size of the stash. 8488 - bsize - the initial size of the block-stash(if used). 8489 8490 Options Database Keys: 8491 + -matstash_initial_size <size> or <size0,size1,...sizep-1> 8492 - -matstash_block_initial_size <bsize> or <bsize0,bsize1,...bsizep-1> 8493 8494 Level: intermediate 8495 8496 Notes: 8497 The block-stash is used for values set with `MatSetValuesBlocked()` while 8498 the stash is used for values set with `MatSetValues()` 8499 8500 Run with the option -info and look for output of the form 8501 MatAssemblyBegin_MPIXXX:Stash has MM entries, uses nn mallocs. 8502 to determine the appropriate value, MM, to use for size and 8503 MatAssemblyBegin_MPIXXX:Block-Stash has BMM entries, uses nn mallocs. 8504 to determine the value, BMM to use for bsize 8505 8506 .seealso: `MatAssemblyBegin()`, `MatAssemblyEnd()`, `Mat`, `MatStashGetInfo()` 8507 @*/ 8508 PetscErrorCode MatStashSetInitialSize(Mat mat, PetscInt size, PetscInt bsize) 8509 { 8510 PetscFunctionBegin; 8511 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 8512 PetscValidType(mat, 1); 8513 PetscCall(MatStashSetInitialSize_Private(&mat->stash, size)); 8514 PetscCall(MatStashSetInitialSize_Private(&mat->bstash, bsize)); 8515 PetscFunctionReturn(PETSC_SUCCESS); 8516 } 8517 8518 /*@ 8519 MatInterpolateAdd - w = y + A*x or A'*x depending on the shape of 8520 the matrix 8521 8522 Neighbor-wise Collective 8523 8524 Input Parameters: 8525 + mat - the matrix 8526 . x,y - the vectors 8527 - w - where the result is stored 8528 8529 Level: intermediate 8530 8531 Notes: 8532 `w` may be the same vector as `y`. 8533 8534 This allows one to use either the restriction or interpolation (its transpose) 8535 matrix to do the interpolation 8536 8537 .seealso: `Mat`, `MatMultAdd()`, `MatMultTransposeAdd()`, `MatRestrict()`, `PCMG` 8538 @*/ 8539 PetscErrorCode MatInterpolateAdd(Mat A, Vec x, Vec y, Vec w) 8540 { 8541 PetscInt M, N, Ny; 8542 8543 PetscFunctionBegin; 8544 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 8545 PetscValidHeaderSpecific(x, VEC_CLASSID, 2); 8546 PetscValidHeaderSpecific(y, VEC_CLASSID, 3); 8547 PetscValidHeaderSpecific(w, VEC_CLASSID, 4); 8548 PetscCall(MatGetSize(A, &M, &N)); 8549 PetscCall(VecGetSize(y, &Ny)); 8550 if (M == Ny) { 8551 PetscCall(MatMultAdd(A, x, y, w)); 8552 } else { 8553 PetscCall(MatMultTransposeAdd(A, x, y, w)); 8554 } 8555 PetscFunctionReturn(PETSC_SUCCESS); 8556 } 8557 8558 /*@ 8559 MatInterpolate - y = A*x or A'*x depending on the shape of 8560 the matrix 8561 8562 Neighbor-wise Collective 8563 8564 Input Parameters: 8565 + mat - the matrix 8566 - x,y - the vectors 8567 8568 Level: intermediate 8569 8570 Note: 8571 This allows one to use either the restriction or interpolation (its transpose) 8572 matrix to do the interpolation 8573 8574 .seealso: `Mat`, `MatMultAdd()`, `MatMultTransposeAdd()`, `MatRestrict()`, `PCMG` 8575 @*/ 8576 PetscErrorCode MatInterpolate(Mat A, Vec x, Vec y) 8577 { 8578 PetscInt M, N, Ny; 8579 8580 PetscFunctionBegin; 8581 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 8582 PetscValidHeaderSpecific(x, VEC_CLASSID, 2); 8583 PetscValidHeaderSpecific(y, VEC_CLASSID, 3); 8584 PetscCall(MatGetSize(A, &M, &N)); 8585 PetscCall(VecGetSize(y, &Ny)); 8586 if (M == Ny) { 8587 PetscCall(MatMult(A, x, y)); 8588 } else { 8589 PetscCall(MatMultTranspose(A, x, y)); 8590 } 8591 PetscFunctionReturn(PETSC_SUCCESS); 8592 } 8593 8594 /*@ 8595 MatRestrict - y = A*x or A'*x 8596 8597 Neighbor-wise Collective on Mat 8598 8599 Input Parameters: 8600 + mat - the matrix 8601 - x,y - the vectors 8602 8603 Level: intermediate 8604 8605 Note: 8606 This allows one to use either the restriction or interpolation (its transpose) 8607 matrix to do the restriction 8608 8609 .seealso: `Mat`, `MatMultAdd()`, `MatMultTransposeAdd()`, `MatInterpolate()`, `PCMG` 8610 @*/ 8611 PetscErrorCode MatRestrict(Mat A, Vec x, Vec y) 8612 { 8613 PetscInt M, N, Ny; 8614 8615 PetscFunctionBegin; 8616 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 8617 PetscValidHeaderSpecific(x, VEC_CLASSID, 2); 8618 PetscValidHeaderSpecific(y, VEC_CLASSID, 3); 8619 PetscCall(MatGetSize(A, &M, &N)); 8620 PetscCall(VecGetSize(y, &Ny)); 8621 if (M == Ny) { 8622 PetscCall(MatMult(A, x, y)); 8623 } else { 8624 PetscCall(MatMultTranspose(A, x, y)); 8625 } 8626 PetscFunctionReturn(PETSC_SUCCESS); 8627 } 8628 8629 /*@ 8630 MatMatInterpolateAdd - Y = W + A*X or W + A'*X 8631 8632 Neighbor-wise Collective on Mat 8633 8634 Input Parameters: 8635 + mat - the matrix 8636 - w, x - the input dense matrices 8637 8638 Output Parameters: 8639 . y - the output dense matrix 8640 8641 Level: intermediate 8642 8643 Note: 8644 This allows one to use either the restriction or interpolation (its transpose) 8645 matrix to do the interpolation. y matrix can be reused if already created with the proper sizes, 8646 otherwise it will be recreated. y must be initialized to `NULL` if not supplied. 8647 8648 .seealso: `Mat`, `MatInterpolateAdd()`, `MatMatInterpolate()`, `MatMatRestrict()`, `PCMG` 8649 @*/ 8650 PetscErrorCode MatMatInterpolateAdd(Mat A, Mat x, Mat w, Mat *y) 8651 { 8652 PetscInt M, N, Mx, Nx, Mo, My = 0, Ny = 0; 8653 PetscBool trans = PETSC_TRUE; 8654 MatReuse reuse = MAT_INITIAL_MATRIX; 8655 8656 PetscFunctionBegin; 8657 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 8658 PetscValidHeaderSpecific(x, MAT_CLASSID, 2); 8659 PetscValidType(x, 2); 8660 if (w) PetscValidHeaderSpecific(w, MAT_CLASSID, 3); 8661 if (*y) PetscValidHeaderSpecific(*y, MAT_CLASSID, 4); 8662 PetscCall(MatGetSize(A, &M, &N)); 8663 PetscCall(MatGetSize(x, &Mx, &Nx)); 8664 if (N == Mx) trans = PETSC_FALSE; 8665 else PetscCheck(M == Mx, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Size mismatch: A %" PetscInt_FMT "x%" PetscInt_FMT ", X %" PetscInt_FMT "x%" PetscInt_FMT, M, N, Mx, Nx); 8666 Mo = trans ? N : M; 8667 if (*y) { 8668 PetscCall(MatGetSize(*y, &My, &Ny)); 8669 if (Mo == My && Nx == Ny) { 8670 reuse = MAT_REUSE_MATRIX; 8671 } else { 8672 PetscCheck(w || *y != w, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Cannot reuse y and w, size mismatch: A %" PetscInt_FMT "x%" PetscInt_FMT ", X %" PetscInt_FMT "x%" PetscInt_FMT ", Y %" PetscInt_FMT "x%" PetscInt_FMT, M, N, Mx, Nx, My, Ny); 8673 PetscCall(MatDestroy(y)); 8674 } 8675 } 8676 8677 if (w && *y == w) { /* this is to minimize changes in PCMG */ 8678 PetscBool flg; 8679 8680 PetscCall(PetscObjectQuery((PetscObject)*y, "__MatMatIntAdd_w", (PetscObject *)&w)); 8681 if (w) { 8682 PetscInt My, Ny, Mw, Nw; 8683 8684 PetscCall(PetscObjectTypeCompare((PetscObject)*y, ((PetscObject)w)->type_name, &flg)); 8685 PetscCall(MatGetSize(*y, &My, &Ny)); 8686 PetscCall(MatGetSize(w, &Mw, &Nw)); 8687 if (!flg || My != Mw || Ny != Nw) w = NULL; 8688 } 8689 if (!w) { 8690 PetscCall(MatDuplicate(*y, MAT_COPY_VALUES, &w)); 8691 PetscCall(PetscObjectCompose((PetscObject)*y, "__MatMatIntAdd_w", (PetscObject)w)); 8692 PetscCall(PetscObjectDereference((PetscObject)w)); 8693 } else { 8694 PetscCall(MatCopy(*y, w, UNKNOWN_NONZERO_PATTERN)); 8695 } 8696 } 8697 if (!trans) { 8698 PetscCall(MatMatMult(A, x, reuse, PETSC_DEFAULT, y)); 8699 } else { 8700 PetscCall(MatTransposeMatMult(A, x, reuse, PETSC_DEFAULT, y)); 8701 } 8702 if (w) PetscCall(MatAXPY(*y, 1.0, w, UNKNOWN_NONZERO_PATTERN)); 8703 PetscFunctionReturn(PETSC_SUCCESS); 8704 } 8705 8706 /*@ 8707 MatMatInterpolate - Y = A*X or A'*X 8708 8709 Neighbor-wise Collective on Mat 8710 8711 Input Parameters: 8712 + mat - the matrix 8713 - x - the input dense matrix 8714 8715 Output Parameters: 8716 . y - the output dense matrix 8717 8718 Level: intermediate 8719 8720 Note: 8721 This allows one to use either the restriction or interpolation (its transpose) 8722 matrix to do the interpolation. y matrix can be reused if already created with the proper sizes, 8723 otherwise it will be recreated. y must be initialized to `NULL` if not supplied. 8724 8725 .seealso: `Mat`, `MatInterpolate()`, `MatRestrict()`, `MatMatRestrict()`, `PCMG` 8726 @*/ 8727 PetscErrorCode MatMatInterpolate(Mat A, Mat x, Mat *y) 8728 { 8729 PetscFunctionBegin; 8730 PetscCall(MatMatInterpolateAdd(A, x, NULL, y)); 8731 PetscFunctionReturn(PETSC_SUCCESS); 8732 } 8733 8734 /*@ 8735 MatMatRestrict - Y = A*X or A'*X 8736 8737 Neighbor-wise Collective on Mat 8738 8739 Input Parameters: 8740 + mat - the matrix 8741 - x - the input dense matrix 8742 8743 Output Parameters: 8744 . y - the output dense matrix 8745 8746 Level: intermediate 8747 8748 Note: 8749 This allows one to use either the restriction or interpolation (its transpose) 8750 matrix to do the restriction. y matrix can be reused if already created with the proper sizes, 8751 otherwise it will be recreated. y must be initialized to `NULL` if not supplied. 8752 8753 .seealso: `Mat`, `MatRestrict()`, `MatInterpolate()`, `MatMatInterpolate()`, `PCMG` 8754 @*/ 8755 PetscErrorCode MatMatRestrict(Mat A, Mat x, Mat *y) 8756 { 8757 PetscFunctionBegin; 8758 PetscCall(MatMatInterpolateAdd(A, x, NULL, y)); 8759 PetscFunctionReturn(PETSC_SUCCESS); 8760 } 8761 8762 /*@ 8763 MatGetNullSpace - retrieves the null space of a matrix. 8764 8765 Logically Collective 8766 8767 Input Parameters: 8768 + mat - the matrix 8769 - nullsp - the null space object 8770 8771 Level: developer 8772 8773 .seealso: `Mat`, `MatCreate()`, `MatNullSpaceCreate()`, `MatSetNearNullSpace()`, `MatSetNullSpace()`, `MatNullSpace` 8774 @*/ 8775 PetscErrorCode MatGetNullSpace(Mat mat, MatNullSpace *nullsp) 8776 { 8777 PetscFunctionBegin; 8778 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 8779 PetscValidPointer(nullsp, 2); 8780 *nullsp = (mat->symmetric == PETSC_BOOL3_TRUE && !mat->nullsp) ? mat->transnullsp : mat->nullsp; 8781 PetscFunctionReturn(PETSC_SUCCESS); 8782 } 8783 8784 /*@ 8785 MatSetNullSpace - attaches a null space to a matrix. 8786 8787 Logically Collective 8788 8789 Input Parameters: 8790 + mat - the matrix 8791 - nullsp - the null space object 8792 8793 Level: advanced 8794 8795 Notes: 8796 This null space is used by the `KSP` linear solvers to solve singular systems. 8797 8798 Overwrites any previous null space that may have been attached. You can remove the null space from the matrix object by calling this routine with an nullsp of `NULL` 8799 8800 For inconsistent singular systems (linear systems where the right hand side is not in the range of the operator) the `KSP` residuals will not converge to 8801 to zero but the linear system will still be solved in a least squares sense. 8802 8803 The fundamental theorem of linear algebra (Gilbert Strang, Introduction to Applied Mathematics, page 72) states that 8804 the domain of a matrix A (from R^n to R^m (m rows, n columns) R^n = the direct sum of the null space of A, n(A), + the range of A^T, R(A^T). 8805 Similarly R^m = direct sum n(A^T) + R(A). Hence the linear system A x = b has a solution only if b in R(A) (or correspondingly b is orthogonal to 8806 n(A^T)) and if x is a solution then x + alpha n(A) is a solution for any alpha. The minimum norm solution is orthogonal to n(A). For problems without a solution 8807 the solution that minimizes the norm of the residual (the least squares solution) can be obtained by solving A x = \hat{b} where \hat{b} is b orthogonalized to the n(A^T). 8808 This \hat{b} can be obtained by calling MatNullSpaceRemove() with the null space of the transpose of the matrix. 8809 8810 If the matrix is known to be symmetric because it is an `MATSBAIJ` matrix or one as called 8811 `MatSetOption`(mat,`MAT_SYMMETRIC` or possibly `MAT_SYMMETRY_ETERNAL`,`PETSC_TRUE`); this 8812 routine also automatically calls `MatSetTransposeNullSpace()`. 8813 8814 The user should call `MatNullSpaceDestroy()`. 8815 8816 .seealso: `Mat`, `MatCreate()`, `MatNullSpaceCreate()`, `MatSetNearNullSpace()`, `MatGetNullSpace()`, `MatSetTransposeNullSpace()`, `MatGetTransposeNullSpace()`, `MatNullSpaceRemove()`, 8817 `KSPSetPCSide()` 8818 @*/ 8819 PetscErrorCode MatSetNullSpace(Mat mat, MatNullSpace nullsp) 8820 { 8821 PetscFunctionBegin; 8822 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 8823 if (nullsp) PetscValidHeaderSpecific(nullsp, MAT_NULLSPACE_CLASSID, 2); 8824 if (nullsp) PetscCall(PetscObjectReference((PetscObject)nullsp)); 8825 PetscCall(MatNullSpaceDestroy(&mat->nullsp)); 8826 mat->nullsp = nullsp; 8827 if (mat->symmetric == PETSC_BOOL3_TRUE) PetscCall(MatSetTransposeNullSpace(mat, nullsp)); 8828 PetscFunctionReturn(PETSC_SUCCESS); 8829 } 8830 8831 /*@ 8832 MatGetTransposeNullSpace - retrieves the null space of the transpose of a matrix. 8833 8834 Logically Collective 8835 8836 Input Parameters: 8837 + mat - the matrix 8838 - nullsp - the null space object 8839 8840 Level: developer 8841 8842 .seealso: `Mat`, `MatNullSpace`, `MatCreate()`, `MatNullSpaceCreate()`, `MatSetNearNullSpace()`, `MatSetTransposeNullSpace()`, `MatSetNullSpace()`, `MatGetNullSpace()` 8843 @*/ 8844 PetscErrorCode MatGetTransposeNullSpace(Mat mat, MatNullSpace *nullsp) 8845 { 8846 PetscFunctionBegin; 8847 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 8848 PetscValidType(mat, 1); 8849 PetscValidPointer(nullsp, 2); 8850 *nullsp = (mat->symmetric == PETSC_BOOL3_TRUE && !mat->transnullsp) ? mat->nullsp : mat->transnullsp; 8851 PetscFunctionReturn(PETSC_SUCCESS); 8852 } 8853 8854 /*@ 8855 MatSetTransposeNullSpace - attaches the null space of a transpose of a matrix to the matrix 8856 8857 Logically Collective 8858 8859 Input Parameters: 8860 + mat - the matrix 8861 - nullsp - the null space object 8862 8863 Level: advanced 8864 8865 Notes: 8866 This allows solving singular linear systems defined by the transpose of the matrix using `KSP` solvers with left preconditioning. 8867 8868 See `MatSetNullSpace()` 8869 8870 .seealso: `Mat`, `MatNullSpace`, `MatCreate()`, `MatNullSpaceCreate()`, `MatSetNearNullSpace()`, `MatGetNullSpace()`, `MatSetNullSpace()`, `MatGetTransposeNullSpace()`, `MatNullSpaceRemove()`, `KSPSetPCSide()` 8871 @*/ 8872 PetscErrorCode MatSetTransposeNullSpace(Mat mat, MatNullSpace nullsp) 8873 { 8874 PetscFunctionBegin; 8875 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 8876 if (nullsp) PetscValidHeaderSpecific(nullsp, MAT_NULLSPACE_CLASSID, 2); 8877 if (nullsp) PetscCall(PetscObjectReference((PetscObject)nullsp)); 8878 PetscCall(MatNullSpaceDestroy(&mat->transnullsp)); 8879 mat->transnullsp = nullsp; 8880 PetscFunctionReturn(PETSC_SUCCESS); 8881 } 8882 8883 /*@ 8884 MatSetNearNullSpace - attaches a null space to a matrix, which is often the null space (rigid body modes) of the operator without boundary conditions 8885 This null space will be used to provide near null space vectors to a multigrid preconditioner built from this matrix. 8886 8887 Logically Collective 8888 8889 Input Parameters: 8890 + mat - the matrix 8891 - nullsp - the null space object 8892 8893 Level: advanced 8894 8895 Notes: 8896 Overwrites any previous near null space that may have been attached 8897 8898 You can remove the null space by calling this routine with an nullsp of `NULL` 8899 8900 .seealso: `Mat`, `MatNullSpace`, `MatCreate()`, `MatNullSpaceCreate()`, `MatSetNullSpace()`, `MatNullSpaceCreateRigidBody()`, `MatGetNearNullSpace()` 8901 @*/ 8902 PetscErrorCode MatSetNearNullSpace(Mat mat, MatNullSpace nullsp) 8903 { 8904 PetscFunctionBegin; 8905 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 8906 PetscValidType(mat, 1); 8907 if (nullsp) PetscValidHeaderSpecific(nullsp, MAT_NULLSPACE_CLASSID, 2); 8908 MatCheckPreallocated(mat, 1); 8909 if (nullsp) PetscCall(PetscObjectReference((PetscObject)nullsp)); 8910 PetscCall(MatNullSpaceDestroy(&mat->nearnullsp)); 8911 mat->nearnullsp = nullsp; 8912 PetscFunctionReturn(PETSC_SUCCESS); 8913 } 8914 8915 /*@ 8916 MatGetNearNullSpace - Get null space attached with `MatSetNearNullSpace()` 8917 8918 Not Collective 8919 8920 Input Parameter: 8921 . mat - the matrix 8922 8923 Output Parameter: 8924 . nullsp - the null space object, `NULL` if not set 8925 8926 Level: advanced 8927 8928 .seealso: `Mat`, `MatNullSpace`, `MatSetNearNullSpace()`, `MatGetNullSpace()`, `MatNullSpaceCreate()` 8929 @*/ 8930 PetscErrorCode MatGetNearNullSpace(Mat mat, MatNullSpace *nullsp) 8931 { 8932 PetscFunctionBegin; 8933 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 8934 PetscValidType(mat, 1); 8935 PetscValidPointer(nullsp, 2); 8936 MatCheckPreallocated(mat, 1); 8937 *nullsp = mat->nearnullsp; 8938 PetscFunctionReturn(PETSC_SUCCESS); 8939 } 8940 8941 /*@C 8942 MatICCFactor - Performs in-place incomplete Cholesky factorization of matrix. 8943 8944 Collective 8945 8946 Input Parameters: 8947 + mat - the matrix 8948 . row - row/column permutation 8949 . fill - expected fill factor >= 1.0 8950 - level - level of fill, for ICC(k) 8951 8952 Notes: 8953 Probably really in-place only when level of fill is zero, otherwise allocates 8954 new space to store factored matrix and deletes previous memory. 8955 8956 Most users should employ the `KSP` interface for linear solvers 8957 instead of working directly with matrix algebra routines such as this. 8958 See, e.g., `KSPCreate()`. 8959 8960 Level: developer 8961 8962 Developer Note: 8963 The Fortran interface is not autogenerated as the 8964 interface definition cannot be generated correctly [due to `MatFactorInfo`] 8965 8966 .seealso: `Mat`, `MatGetFactor()`, `MatICCFactorSymbolic()`, `MatLUFactorNumeric()`, `MatCholeskyFactor()` 8967 @*/ 8968 PetscErrorCode MatICCFactor(Mat mat, IS row, const MatFactorInfo *info) 8969 { 8970 PetscFunctionBegin; 8971 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 8972 PetscValidType(mat, 1); 8973 if (row) PetscValidHeaderSpecific(row, IS_CLASSID, 2); 8974 PetscValidPointer(info, 3); 8975 PetscCheck(mat->rmap->N == mat->cmap->N, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONG, "matrix must be square"); 8976 PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 8977 PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 8978 MatCheckPreallocated(mat, 1); 8979 PetscUseTypeMethod(mat, iccfactor, row, info); 8980 PetscCall(PetscObjectStateIncrease((PetscObject)mat)); 8981 PetscFunctionReturn(PETSC_SUCCESS); 8982 } 8983 8984 /*@ 8985 MatDiagonalScaleLocal - Scales columns of a matrix given the scaling values including the 8986 ghosted ones. 8987 8988 Not Collective 8989 8990 Input Parameters: 8991 + mat - the matrix 8992 - diag - the diagonal values, including ghost ones 8993 8994 Level: developer 8995 8996 Notes: 8997 Works only for `MATMPIAIJ` and `MATMPIBAIJ` matrices 8998 8999 This allows one to avoid during communication to perform the scaling that must be done with `MatDiagonalScale()` 9000 9001 .seealso: `Mat`, `MatDiagonalScale()` 9002 @*/ 9003 PetscErrorCode MatDiagonalScaleLocal(Mat mat, Vec diag) 9004 { 9005 PetscMPIInt size; 9006 9007 PetscFunctionBegin; 9008 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 9009 PetscValidHeaderSpecific(diag, VEC_CLASSID, 2); 9010 PetscValidType(mat, 1); 9011 9012 PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Matrix must be already assembled"); 9013 PetscCall(PetscLogEventBegin(MAT_Scale, mat, 0, 0, 0)); 9014 PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)mat), &size)); 9015 if (size == 1) { 9016 PetscInt n, m; 9017 PetscCall(VecGetSize(diag, &n)); 9018 PetscCall(MatGetSize(mat, NULL, &m)); 9019 PetscCheck(m == n, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only supported for sequential matrices when no ghost points/periodic conditions"); 9020 PetscCall(MatDiagonalScale(mat, NULL, diag)); 9021 } else { 9022 PetscUseMethod(mat, "MatDiagonalScaleLocal_C", (Mat, Vec), (mat, diag)); 9023 } 9024 PetscCall(PetscLogEventEnd(MAT_Scale, mat, 0, 0, 0)); 9025 PetscCall(PetscObjectStateIncrease((PetscObject)mat)); 9026 PetscFunctionReturn(PETSC_SUCCESS); 9027 } 9028 9029 /*@ 9030 MatGetInertia - Gets the inertia from a factored matrix 9031 9032 Collective 9033 9034 Input Parameter: 9035 . mat - the matrix 9036 9037 Output Parameters: 9038 + nneg - number of negative eigenvalues 9039 . nzero - number of zero eigenvalues 9040 - npos - number of positive eigenvalues 9041 9042 Level: advanced 9043 9044 Note: 9045 Matrix must have been factored by `MatCholeskyFactor()` 9046 9047 .seealso: `Mat`, `MatGetFactor()`, `MatCholeskyFactor()` 9048 @*/ 9049 PetscErrorCode MatGetInertia(Mat mat, PetscInt *nneg, PetscInt *nzero, PetscInt *npos) 9050 { 9051 PetscFunctionBegin; 9052 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 9053 PetscValidType(mat, 1); 9054 PetscCheck(mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Unfactored matrix"); 9055 PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Numeric factor mat is not assembled"); 9056 PetscUseTypeMethod(mat, getinertia, nneg, nzero, npos); 9057 PetscFunctionReturn(PETSC_SUCCESS); 9058 } 9059 9060 /*@C 9061 MatSolves - Solves A x = b, given a factored matrix, for a collection of vectors 9062 9063 Neighbor-wise Collective 9064 9065 Input Parameters: 9066 + mat - the factored matrix obtained with `MatGetFactor()` 9067 - b - the right-hand-side vectors 9068 9069 Output Parameter: 9070 . x - the result vectors 9071 9072 Level: developer 9073 9074 Note: 9075 The vectors `b` and `x` cannot be the same. I.e., one cannot 9076 call `MatSolves`(A,x,x). 9077 9078 .seealso: `Mat`, `Vecs`, `MatSolveAdd()`, `MatSolveTranspose()`, `MatSolveTransposeAdd()`, `MatSolve()` 9079 @*/ 9080 PetscErrorCode MatSolves(Mat mat, Vecs b, Vecs x) 9081 { 9082 PetscFunctionBegin; 9083 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 9084 PetscValidType(mat, 1); 9085 PetscCheck(x != b, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_IDN, "x and b must be different vectors"); 9086 PetscCheck(mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Unfactored matrix"); 9087 if (!mat->rmap->N && !mat->cmap->N) PetscFunctionReturn(PETSC_SUCCESS); 9088 9089 MatCheckPreallocated(mat, 1); 9090 PetscCall(PetscLogEventBegin(MAT_Solves, mat, 0, 0, 0)); 9091 PetscUseTypeMethod(mat, solves, b, x); 9092 PetscCall(PetscLogEventEnd(MAT_Solves, mat, 0, 0, 0)); 9093 PetscFunctionReturn(PETSC_SUCCESS); 9094 } 9095 9096 /*@ 9097 MatIsSymmetric - Test whether a matrix is symmetric 9098 9099 Collective 9100 9101 Input Parameters: 9102 + A - the matrix to test 9103 - tol - difference between value and its transpose less than this amount counts as equal (use 0.0 for exact transpose) 9104 9105 Output Parameters: 9106 . flg - the result 9107 9108 Level: intermediate 9109 9110 Notes: 9111 For real numbers `MatIsSymmetric()` and `MatIsHermitian()` return identical results 9112 9113 If the matrix does not yet know if it is symmetric or not this can be an expensive operation, also available `MatIsSymmetricKnown()` 9114 9115 One can declare that a matrix is symmetric with `MatSetOption`(mat,`MAT_SYMMETRIC`,`PETSC_TRUE`) and if it is known to remain symmetric 9116 after changes to the matrices values one can call `MatSetOption`(mat,`MAT_SYMMETRY_ETERNAL`,`PETSC_TRUE`) 9117 9118 .seealso: `Mat`, `MatTranspose()`, `MatIsTranspose()`, `MatIsHermitian()`, `MatIsStructurallySymmetric()`, `MatSetOption()`, `MatIsSymmetricKnown()`, 9119 `MAT_SYMMETRIC`, `MAT_SYMMETRY_ETERNAL`, `MatSetOption()` 9120 @*/ 9121 PetscErrorCode MatIsSymmetric(Mat A, PetscReal tol, PetscBool *flg) 9122 { 9123 PetscFunctionBegin; 9124 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 9125 PetscValidBoolPointer(flg, 3); 9126 9127 if (A->symmetric == PETSC_BOOL3_TRUE) *flg = PETSC_TRUE; 9128 else if (A->symmetric == PETSC_BOOL3_FALSE) *flg = PETSC_FALSE; 9129 else { 9130 if (!A->ops->issymmetric) { 9131 MatType mattype; 9132 PetscCall(MatGetType(A, &mattype)); 9133 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Matrix of type %s does not support checking for symmetric", mattype); 9134 } 9135 PetscUseTypeMethod(A, issymmetric, tol, flg); 9136 if (!tol) PetscCall(MatSetOption(A, MAT_SYMMETRIC, *flg)); 9137 } 9138 PetscFunctionReturn(PETSC_SUCCESS); 9139 } 9140 9141 /*@ 9142 MatIsHermitian - Test whether a matrix is Hermitian 9143 9144 Collective on Mat 9145 9146 Input Parameters: 9147 + A - the matrix to test 9148 - tol - difference between value and its transpose less than this amount counts as equal (use 0.0 for exact Hermitian) 9149 9150 Output Parameters: 9151 . flg - the result 9152 9153 Level: intermediate 9154 9155 Notes: 9156 For real numbers `MatIsSymmetric()` and `MatIsHermitian()` return identical results 9157 9158 If the matrix does not yet know if it is Hermitian or not this can be an expensive operation, also available `MatIsHermitianKnown()` 9159 9160 One can declare that a matrix is Hermitian with `MatSetOption`(mat,`MAT_HERMITIAN`,`PETSC_TRUE`) and if it is known to remain Hermitian 9161 after changes to the matrices values one can call `MatSetOption`(mat,`MAT_SYMEMTRY_ETERNAL`,`PETSC_TRUE`) 9162 9163 .seealso: `Mat`, `MatTranspose()`, `MatIsTranspose()`, `MatIsHermitianKnown()`, `MatIsStructurallySymmetric()`, `MatSetOption()`, 9164 `MatIsSymmetricKnown()`, `MatIsSymmetric()`, `MAT_HERMITIAN`, `MAT_SYMMETRY_ETERNAL`, `MatSetOption()` 9165 @*/ 9166 PetscErrorCode MatIsHermitian(Mat A, PetscReal tol, PetscBool *flg) 9167 { 9168 PetscFunctionBegin; 9169 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 9170 PetscValidBoolPointer(flg, 3); 9171 9172 if (A->hermitian == PETSC_BOOL3_TRUE) *flg = PETSC_TRUE; 9173 else if (A->hermitian == PETSC_BOOL3_FALSE) *flg = PETSC_FALSE; 9174 else { 9175 if (!A->ops->ishermitian) { 9176 MatType mattype; 9177 PetscCall(MatGetType(A, &mattype)); 9178 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Matrix of type %s does not support checking for hermitian", mattype); 9179 } 9180 PetscUseTypeMethod(A, ishermitian, tol, flg); 9181 if (!tol) PetscCall(MatSetOption(A, MAT_HERMITIAN, *flg)); 9182 } 9183 PetscFunctionReturn(PETSC_SUCCESS); 9184 } 9185 9186 /*@ 9187 MatIsSymmetricKnown - Checks if a matrix knows if it is symmetric or not and its symmetric state 9188 9189 Not Collective 9190 9191 Input Parameter: 9192 . A - the matrix to check 9193 9194 Output Parameters: 9195 + set - `PETSC_TRUE` if the matrix knows its symmetry state (this tells you if the next flag is valid) 9196 - flg - the result (only valid if set is `PETSC_TRUE`) 9197 9198 Level: advanced 9199 9200 Notes: 9201 Does not check the matrix values directly, so this may return unknown (set = `PETSC_FALSE`). Use `MatIsSymmetric()` 9202 if you want it explicitly checked 9203 9204 One can declare that a matrix is symmetric with `MatSetOption`(mat,`MAT_SYMMETRIC`,`PETSC_TRUE`) and if it is known to remain symmetric 9205 after changes to the matrices values one can call `MatSetOption`(mat,`MAT_SYMMETRY_ETERNAL`,`PETSC_TRUE`) 9206 9207 .seealso: `Mat`, `MAT_SYMMETRY_ETERNAL`, `MatTranspose()`, `MatIsTranspose()`, `MatIsHermitian()`, `MatIsStructurallySymmetric()`, `MatSetOption()`, `MatIsSymmetric()`, `MatIsHermitianKnown()` 9208 @*/ 9209 PetscErrorCode MatIsSymmetricKnown(Mat A, PetscBool *set, PetscBool *flg) 9210 { 9211 PetscFunctionBegin; 9212 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 9213 PetscValidBoolPointer(set, 2); 9214 PetscValidBoolPointer(flg, 3); 9215 if (A->symmetric != PETSC_BOOL3_UNKNOWN) { 9216 *set = PETSC_TRUE; 9217 *flg = PetscBool3ToBool(A->symmetric); 9218 } else { 9219 *set = PETSC_FALSE; 9220 } 9221 PetscFunctionReturn(PETSC_SUCCESS); 9222 } 9223 9224 /*@ 9225 MatIsSPDKnown - Checks if a matrix knows if it is symmetric positive definite or not and its symmetric positive definite state 9226 9227 Not Collective 9228 9229 Input Parameter: 9230 . A - the matrix to check 9231 9232 Output Parameters: 9233 + set - `PETSC_TRUE` if the matrix knows its symmetric positive definite state (this tells you if the next flag is valid) 9234 - flg - the result (only valid if set is `PETSC_TRUE`) 9235 9236 Level: advanced 9237 9238 Notes: 9239 Does not check the matrix values directly, so this may return unknown (set = `PETSC_FALSE`). 9240 9241 One can declare that a matrix is SPD with `MatSetOption`(mat,`MAT_SPD`,`PETSC_TRUE`) and if it is known to remain SPD 9242 after changes to the matrices values one can call `MatSetOption`(mat,`MAT_SPD_ETERNAL`,`PETSC_TRUE`) 9243 9244 .seealso: `Mat`, `MAT_SPD_ETERNAL`, `MAT_SPD`, `MatTranspose()`, `MatIsTranspose()`, `MatIsHermitian()`, `MatIsStructurallySymmetric()`, `MatSetOption()`, `MatIsSymmetric()`, `MatIsHermitianKnown()` 9245 @*/ 9246 PetscErrorCode MatIsSPDKnown(Mat A, PetscBool *set, PetscBool *flg) 9247 { 9248 PetscFunctionBegin; 9249 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 9250 PetscValidBoolPointer(set, 2); 9251 PetscValidBoolPointer(flg, 3); 9252 if (A->spd != PETSC_BOOL3_UNKNOWN) { 9253 *set = PETSC_TRUE; 9254 *flg = PetscBool3ToBool(A->spd); 9255 } else { 9256 *set = PETSC_FALSE; 9257 } 9258 PetscFunctionReturn(PETSC_SUCCESS); 9259 } 9260 9261 /*@ 9262 MatIsHermitianKnown - Checks if a matrix knows if it is Hermitian or not and its Hermitian state 9263 9264 Not Collective 9265 9266 Input Parameter: 9267 . A - the matrix to check 9268 9269 Output Parameters: 9270 + set - `PETSC_TRUE` if the matrix knows its Hermitian state (this tells you if the next flag is valid) 9271 - flg - the result (only valid if set is `PETSC_TRUE`) 9272 9273 Level: advanced 9274 9275 Notes: 9276 Does not check the matrix values directly, so this may return unknown (set = `PETSC_FALSE`). Use `MatIsHermitian()` 9277 if you want it explicitly checked 9278 9279 One can declare that a matrix is Hermitian with `MatSetOption`(mat,`MAT_HERMITIAN`,`PETSC_TRUE`) and if it is known to remain Hermitian 9280 after changes to the matrices values one can call `MatSetOption`(mat,`MAT_SYMMETRY_ETERNAL`,`PETSC_TRUE`) 9281 9282 .seealso: `Mat`, `MAT_SYMMETRY_ETERNAL`, `MAT_HERMITIAN`, `MatTranspose()`, `MatIsTranspose()`, `MatIsHermitian()`, `MatIsStructurallySymmetric()`, `MatSetOption()`, `MatIsSymmetric()` 9283 @*/ 9284 PetscErrorCode MatIsHermitianKnown(Mat A, PetscBool *set, PetscBool *flg) 9285 { 9286 PetscFunctionBegin; 9287 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 9288 PetscValidBoolPointer(set, 2); 9289 PetscValidBoolPointer(flg, 3); 9290 if (A->hermitian != PETSC_BOOL3_UNKNOWN) { 9291 *set = PETSC_TRUE; 9292 *flg = PetscBool3ToBool(A->hermitian); 9293 } else { 9294 *set = PETSC_FALSE; 9295 } 9296 PetscFunctionReturn(PETSC_SUCCESS); 9297 } 9298 9299 /*@ 9300 MatIsStructurallySymmetric - Test whether a matrix is structurally symmetric 9301 9302 Collective on Mat 9303 9304 Input Parameter: 9305 . A - the matrix to test 9306 9307 Output Parameters: 9308 . flg - the result 9309 9310 Level: intermediate 9311 9312 Notes: 9313 If the matrix does yet know it is structurally symmetric this can be an expensive operation, also available `MatIsStructurallySymmetricKnown()` 9314 9315 One can declare that a matrix is structurally symmetric with `MatSetOption`(mat,`MAT_STRUCTURALLY_SYMMETRIC`,`PETSC_TRUE`) and if it is known to remain structurally 9316 symmetric after changes to the matrices values one can call `MatSetOption`(mat,`MAT_STRUCTURAL_SYMMETRY_ETERNAL`,`PETSC_TRUE`) 9317 9318 .seealso: `Mat`, `MAT_STRUCTURALLY_SYMMETRIC`, `MAT_STRUCTURAL_SYMMETRY_ETERNAL`, `MatTranspose()`, `MatIsTranspose()`, `MatIsHermitian()`, `MatIsSymmetric()`, `MatSetOption()`, `MatIsStructurallySymmetricKnown()` 9319 @*/ 9320 PetscErrorCode MatIsStructurallySymmetric(Mat A, PetscBool *flg) 9321 { 9322 PetscFunctionBegin; 9323 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 9324 PetscValidBoolPointer(flg, 2); 9325 if (A->structurally_symmetric != PETSC_BOOL3_UNKNOWN) { 9326 *flg = PetscBool3ToBool(A->structurally_symmetric); 9327 } else { 9328 PetscUseTypeMethod(A, isstructurallysymmetric, flg); 9329 PetscCall(MatSetOption(A, MAT_STRUCTURALLY_SYMMETRIC, *flg)); 9330 } 9331 PetscFunctionReturn(PETSC_SUCCESS); 9332 } 9333 9334 /*@ 9335 MatIsStructurallySymmetricKnown - Checks if a matrix knows if it is structurally symmetric or not and its structurally symmetric state 9336 9337 Not Collective 9338 9339 Input Parameter: 9340 . A - the matrix to check 9341 9342 Output Parameters: 9343 + set - PETSC_TRUE if the matrix knows its structurally symmetric state (this tells you if the next flag is valid) 9344 - flg - the result (only valid if set is PETSC_TRUE) 9345 9346 Level: advanced 9347 9348 Notes: 9349 One can declare that a matrix is structurally symmetric with `MatSetOption`(mat,`MAT_STRUCTURALLY_SYMMETRIC`,`PETSC_TRUE`) and if it is known to remain structurally 9350 symmetric after changes to the matrices values one can call `MatSetOption`(mat,`MAT_STRUCTURAL_SYMMETRY_ETERNAL`,`PETSC_TRUE`) 9351 9352 Use `MatIsStructurallySymmetric()` to explicitly check if a matrix is structurally symmetric (this is an expensive operation) 9353 9354 .seealso: `Mat`, `MAT_STRUCTURALLY_SYMMETRIC`, `MatTranspose()`, `MatIsTranspose()`, `MatIsHermitian()`, `MatIsStructurallySymmetric()`, `MatSetOption()`, `MatIsSymmetric()`, `MatIsHermitianKnown()` 9355 @*/ 9356 PetscErrorCode MatIsStructurallySymmetricKnown(Mat A, PetscBool *set, PetscBool *flg) 9357 { 9358 PetscFunctionBegin; 9359 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 9360 PetscValidBoolPointer(set, 2); 9361 PetscValidBoolPointer(flg, 3); 9362 if (A->structurally_symmetric != PETSC_BOOL3_UNKNOWN) { 9363 *set = PETSC_TRUE; 9364 *flg = PetscBool3ToBool(A->structurally_symmetric); 9365 } else { 9366 *set = PETSC_FALSE; 9367 } 9368 PetscFunctionReturn(PETSC_SUCCESS); 9369 } 9370 9371 /*@ 9372 MatStashGetInfo - Gets how many values are currently in the matrix stash, i.e. need 9373 to be communicated to other processors during the `MatAssemblyBegin()`/`MatAssemblyEnd()` process 9374 9375 Not collective 9376 9377 Input Parameter: 9378 . mat - the matrix 9379 9380 Output Parameters: 9381 + nstash - the size of the stash 9382 . reallocs - the number of additional mallocs incurred. 9383 . bnstash - the size of the block stash 9384 - breallocs - the number of additional mallocs incurred.in the block stash 9385 9386 Level: advanced 9387 9388 .seealso: `MatAssemblyBegin()`, `MatAssemblyEnd()`, `Mat`, `MatStashSetInitialSize()` 9389 @*/ 9390 PetscErrorCode MatStashGetInfo(Mat mat, PetscInt *nstash, PetscInt *reallocs, PetscInt *bnstash, PetscInt *breallocs) 9391 { 9392 PetscFunctionBegin; 9393 PetscCall(MatStashGetInfo_Private(&mat->stash, nstash, reallocs)); 9394 PetscCall(MatStashGetInfo_Private(&mat->bstash, bnstash, breallocs)); 9395 PetscFunctionReturn(PETSC_SUCCESS); 9396 } 9397 9398 /*@C 9399 MatCreateVecs - Get vector(s) compatible with the matrix, i.e. with the same 9400 parallel layout, `PetscLayout` for rows and columns 9401 9402 Collective 9403 9404 Input Parameter: 9405 . mat - the matrix 9406 9407 Output Parameters: 9408 + right - (optional) vector that the matrix can be multiplied against 9409 - left - (optional) vector that the matrix vector product can be stored in 9410 9411 Level: advanced 9412 9413 Notes: 9414 The blocksize of the returned vectors is determined by the row and column block sizes set with `MatSetBlockSizes()` or the single blocksize (same for both) set by `MatSetBlockSize()`. 9415 9416 These are new vectors which are not owned by the mat, they should be destroyed in `VecDestroy()` when no longer needed 9417 9418 .seealso: `Mat`, `Vec`, `VecCreate()`, `VecDestroy()`, `DMCreateGlobalVector()` 9419 @*/ 9420 PetscErrorCode MatCreateVecs(Mat mat, Vec *right, Vec *left) 9421 { 9422 PetscFunctionBegin; 9423 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 9424 PetscValidType(mat, 1); 9425 if (mat->ops->getvecs) { 9426 PetscUseTypeMethod(mat, getvecs, right, left); 9427 } else { 9428 PetscInt rbs, cbs; 9429 PetscCall(MatGetBlockSizes(mat, &rbs, &cbs)); 9430 if (right) { 9431 PetscCheck(mat->cmap->n >= 0, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "PetscLayout for columns not yet setup"); 9432 PetscCall(VecCreate(PetscObjectComm((PetscObject)mat), right)); 9433 PetscCall(VecSetSizes(*right, mat->cmap->n, PETSC_DETERMINE)); 9434 PetscCall(VecSetBlockSize(*right, cbs)); 9435 PetscCall(VecSetType(*right, mat->defaultvectype)); 9436 #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) || defined(PETSC_HAVE_HIP) 9437 if (mat->boundtocpu && mat->bindingpropagates) { 9438 PetscCall(VecSetBindingPropagates(*right, PETSC_TRUE)); 9439 PetscCall(VecBindToCPU(*right, PETSC_TRUE)); 9440 } 9441 #endif 9442 PetscCall(PetscLayoutReference(mat->cmap, &(*right)->map)); 9443 } 9444 if (left) { 9445 PetscCheck(mat->rmap->n >= 0, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "PetscLayout for rows not yet setup"); 9446 PetscCall(VecCreate(PetscObjectComm((PetscObject)mat), left)); 9447 PetscCall(VecSetSizes(*left, mat->rmap->n, PETSC_DETERMINE)); 9448 PetscCall(VecSetBlockSize(*left, rbs)); 9449 PetscCall(VecSetType(*left, mat->defaultvectype)); 9450 #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) || defined(PETSC_HAVE_HIP) 9451 if (mat->boundtocpu && mat->bindingpropagates) { 9452 PetscCall(VecSetBindingPropagates(*left, PETSC_TRUE)); 9453 PetscCall(VecBindToCPU(*left, PETSC_TRUE)); 9454 } 9455 #endif 9456 PetscCall(PetscLayoutReference(mat->rmap, &(*left)->map)); 9457 } 9458 } 9459 PetscFunctionReturn(PETSC_SUCCESS); 9460 } 9461 9462 /*@C 9463 MatFactorInfoInitialize - Initializes a `MatFactorInfo` data structure 9464 with default values. 9465 9466 Not Collective 9467 9468 Input Parameters: 9469 . info - the `MatFactorInfo` data structure 9470 9471 Level: developer 9472 9473 Notes: 9474 The solvers are generally used through the `KSP` and `PC` objects, for example 9475 `PCLU`, `PCILU`, `PCCHOLESKY`, `PCICC` 9476 9477 Once the data structure is initialized one may change certain entries as desired for the particular factorization to be performed 9478 9479 Developer Note: 9480 The Fortran interface is not autogenerated as the 9481 interface definition cannot be generated correctly [due to `MatFactorInfo`] 9482 9483 .seealso: `Mat`, `MatGetFactor()`, `MatFactorInfo` 9484 @*/ 9485 PetscErrorCode MatFactorInfoInitialize(MatFactorInfo *info) 9486 { 9487 PetscFunctionBegin; 9488 PetscCall(PetscMemzero(info, sizeof(MatFactorInfo))); 9489 PetscFunctionReturn(PETSC_SUCCESS); 9490 } 9491 9492 /*@ 9493 MatFactorSetSchurIS - Set indices corresponding to the Schur complement you wish to have computed 9494 9495 Collective 9496 9497 Input Parameters: 9498 + mat - the factored matrix 9499 - is - the index set defining the Schur indices (0-based) 9500 9501 Level: advanced 9502 9503 Notes: 9504 Call `MatFactorSolveSchurComplement()` or `MatFactorSolveSchurComplementTranspose()` after this call to solve a Schur complement system. 9505 9506 You can call `MatFactorGetSchurComplement()` or `MatFactorCreateSchurComplement()` after this call. 9507 9508 This functionality is only supported for `MATSOLVERMUMPS` and `MATSOLVERMKL_PARDISO` 9509 9510 .seealso: `Mat`, `MatGetFactor()`, `MatFactorGetSchurComplement()`, `MatFactorRestoreSchurComplement()`, `MatFactorCreateSchurComplement()`, `MatFactorSolveSchurComplement()`, 9511 `MatFactorSolveSchurComplementTranspose()`, `MatFactorSolveSchurComplement()`, `MATSOLVERMUMPS`, `MATSOLVERMKL_PARDISO` 9512 @*/ 9513 PetscErrorCode MatFactorSetSchurIS(Mat mat, IS is) 9514 { 9515 PetscErrorCode (*f)(Mat, IS); 9516 9517 PetscFunctionBegin; 9518 PetscValidType(mat, 1); 9519 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 9520 PetscValidType(is, 2); 9521 PetscValidHeaderSpecific(is, IS_CLASSID, 2); 9522 PetscCheckSameComm(mat, 1, is, 2); 9523 PetscCheck(mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Only for factored matrix"); 9524 PetscCall(PetscObjectQueryFunction((PetscObject)mat, "MatFactorSetSchurIS_C", &f)); 9525 PetscCheck(f, PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "The selected MatSolverType does not support Schur complement computation. You should use MATSOLVERMUMPS or MATSOLVERMKL_PARDISO"); 9526 PetscCall(MatDestroy(&mat->schur)); 9527 PetscCall((*f)(mat, is)); 9528 PetscCheck(mat->schur, PetscObjectComm((PetscObject)mat), PETSC_ERR_PLIB, "Schur complement has not been created"); 9529 PetscFunctionReturn(PETSC_SUCCESS); 9530 } 9531 9532 /*@ 9533 MatFactorCreateSchurComplement - Create a Schur complement matrix object using Schur data computed during the factorization step 9534 9535 Logically Collective 9536 9537 Input Parameters: 9538 + F - the factored matrix obtained by calling `MatGetFactor()` 9539 . S - location where to return the Schur complement, can be `NULL` 9540 - status - the status of the Schur complement matrix, can be `NULL` 9541 9542 Level: advanced 9543 9544 Notes: 9545 You must call `MatFactorSetSchurIS()` before calling this routine. 9546 9547 This functionality is only supported for `MATSOLVERMUMPS` and `MATSOLVERMKL_PARDISO` 9548 9549 The routine provides a copy of the Schur matrix stored within the solver data structures. 9550 The caller must destroy the object when it is no longer needed. 9551 If `MatFactorInvertSchurComplement()` has been called, the routine gets back the inverse. 9552 9553 Use `MatFactorGetSchurComplement()` to get access to the Schur complement matrix inside the factored matrix instead of making a copy of it (which this function does) 9554 9555 See `MatCreateSchurComplement()` or `MatGetSchurComplement()` for ways to create virtual or approximate Schur complements. 9556 9557 Developer Note: 9558 The reason this routine exists is because the representation of the Schur complement within the factor matrix may be different than a standard PETSc 9559 matrix representation and we normally do not want to use the time or memory to make a copy as a regular PETSc matrix. 9560 9561 .seealso: `Mat`, `MatGetFactor()`, `MatFactorSetSchurIS()`, `MatFactorGetSchurComplement()`, `MatFactorSchurStatus`, `MATSOLVERMUMPS`, `MATSOLVERMKL_PARDISO` 9562 @*/ 9563 PetscErrorCode MatFactorCreateSchurComplement(Mat F, Mat *S, MatFactorSchurStatus *status) 9564 { 9565 PetscFunctionBegin; 9566 PetscValidHeaderSpecific(F, MAT_CLASSID, 1); 9567 if (S) PetscValidPointer(S, 2); 9568 if (status) PetscValidPointer(status, 3); 9569 if (S) { 9570 PetscErrorCode (*f)(Mat, Mat *); 9571 9572 PetscCall(PetscObjectQueryFunction((PetscObject)F, "MatFactorCreateSchurComplement_C", &f)); 9573 if (f) { 9574 PetscCall((*f)(F, S)); 9575 } else { 9576 PetscCall(MatDuplicate(F->schur, MAT_COPY_VALUES, S)); 9577 } 9578 } 9579 if (status) *status = F->schur_status; 9580 PetscFunctionReturn(PETSC_SUCCESS); 9581 } 9582 9583 /*@ 9584 MatFactorGetSchurComplement - Gets access to a Schur complement matrix using the current Schur data within a factored matrix 9585 9586 Logically Collective 9587 9588 Input Parameters: 9589 + F - the factored matrix obtained by calling `MatGetFactor()` 9590 . *S - location where to return the Schur complement, can be `NULL` 9591 - status - the status of the Schur complement matrix, can be `NULL` 9592 9593 Level: advanced 9594 9595 Notes: 9596 You must call `MatFactorSetSchurIS()` before calling this routine. 9597 9598 Schur complement mode is currently implemented for sequential matrices with factor type of `MATSOLVERMUMPS` 9599 9600 The routine returns a the Schur Complement stored within the data structures of the solver. 9601 9602 If `MatFactorInvertSchurComplement()` has previously been called, the returned matrix is actually the inverse of the Schur complement. 9603 9604 The returned matrix should not be destroyed; the caller should call `MatFactorRestoreSchurComplement()` when the object is no longer needed. 9605 9606 Use `MatFactorCreateSchurComplement()` to create a copy of the Schur complement matrix that is within a factored matrix 9607 9608 See `MatCreateSchurComplement()` or `MatGetSchurComplement()` for ways to create virtual or approximate Schur complements. 9609 9610 .seealso: `Mat`, `MatGetFactor()`, `MatFactorSetSchurIS()`, `MatFactorRestoreSchurComplement()`, `MatFactorCreateSchurComplement()`, `MatFactorSchurStatus` 9611 @*/ 9612 PetscErrorCode MatFactorGetSchurComplement(Mat F, Mat *S, MatFactorSchurStatus *status) 9613 { 9614 PetscFunctionBegin; 9615 PetscValidHeaderSpecific(F, MAT_CLASSID, 1); 9616 if (S) PetscValidPointer(S, 2); 9617 if (status) PetscValidPointer(status, 3); 9618 if (S) *S = F->schur; 9619 if (status) *status = F->schur_status; 9620 PetscFunctionReturn(PETSC_SUCCESS); 9621 } 9622 9623 static PetscErrorCode MatFactorUpdateSchurStatus_Private(Mat F) 9624 { 9625 Mat S = F->schur; 9626 9627 PetscFunctionBegin; 9628 switch (F->schur_status) { 9629 case MAT_FACTOR_SCHUR_UNFACTORED: // fall-through 9630 case MAT_FACTOR_SCHUR_INVERTED: 9631 if (S) { 9632 S->ops->solve = NULL; 9633 S->ops->matsolve = NULL; 9634 S->ops->solvetranspose = NULL; 9635 S->ops->matsolvetranspose = NULL; 9636 S->ops->solveadd = NULL; 9637 S->ops->solvetransposeadd = NULL; 9638 S->factortype = MAT_FACTOR_NONE; 9639 PetscCall(PetscFree(S->solvertype)); 9640 } 9641 case MAT_FACTOR_SCHUR_FACTORED: // fall-through 9642 break; 9643 default: 9644 SETERRQ(PetscObjectComm((PetscObject)F), PETSC_ERR_SUP, "Unhandled MatFactorSchurStatus %d", F->schur_status); 9645 } 9646 PetscFunctionReturn(PETSC_SUCCESS); 9647 } 9648 9649 /*@ 9650 MatFactorRestoreSchurComplement - Restore the Schur complement matrix object obtained from a call to `MatFactorGetSchurComplement()` 9651 9652 Logically Collective 9653 9654 Input Parameters: 9655 + F - the factored matrix obtained by calling `MatGetFactor()` 9656 . *S - location where the Schur complement is stored 9657 - status - the status of the Schur complement matrix (see `MatFactorSchurStatus`) 9658 9659 Level: advanced 9660 9661 .seealso: `Mat`, `MatGetFactor()`, `MatFactorSetSchurIS()`, `MatFactorRestoreSchurComplement()`, `MatFactorCreateSchurComplement()`, `MatFactorSchurStatus` 9662 @*/ 9663 PetscErrorCode MatFactorRestoreSchurComplement(Mat F, Mat *S, MatFactorSchurStatus status) 9664 { 9665 PetscFunctionBegin; 9666 PetscValidHeaderSpecific(F, MAT_CLASSID, 1); 9667 if (S) { 9668 PetscValidHeaderSpecific(*S, MAT_CLASSID, 2); 9669 *S = NULL; 9670 } 9671 F->schur_status = status; 9672 PetscCall(MatFactorUpdateSchurStatus_Private(F)); 9673 PetscFunctionReturn(PETSC_SUCCESS); 9674 } 9675 9676 /*@ 9677 MatFactorSolveSchurComplementTranspose - Solve the transpose of the Schur complement system computed during the factorization step 9678 9679 Logically Collective 9680 9681 Input Parameters: 9682 + F - the factored matrix obtained by calling `MatGetFactor()` 9683 . rhs - location where the right hand side of the Schur complement system is stored 9684 - sol - location where the solution of the Schur complement system has to be returned 9685 9686 Level: advanced 9687 9688 Notes: 9689 The sizes of the vectors should match the size of the Schur complement 9690 9691 Must be called after `MatFactorSetSchurIS()` 9692 9693 .seealso: `Mat`, `MatGetFactor()`, `MatFactorSetSchurIS()`, `MatFactorSolveSchurComplement()` 9694 @*/ 9695 PetscErrorCode MatFactorSolveSchurComplementTranspose(Mat F, Vec rhs, Vec sol) 9696 { 9697 PetscFunctionBegin; 9698 PetscValidType(F, 1); 9699 PetscValidType(rhs, 2); 9700 PetscValidType(sol, 3); 9701 PetscValidHeaderSpecific(F, MAT_CLASSID, 1); 9702 PetscValidHeaderSpecific(rhs, VEC_CLASSID, 2); 9703 PetscValidHeaderSpecific(sol, VEC_CLASSID, 3); 9704 PetscCheckSameComm(F, 1, rhs, 2); 9705 PetscCheckSameComm(F, 1, sol, 3); 9706 PetscCall(MatFactorFactorizeSchurComplement(F)); 9707 switch (F->schur_status) { 9708 case MAT_FACTOR_SCHUR_FACTORED: 9709 PetscCall(MatSolveTranspose(F->schur, rhs, sol)); 9710 break; 9711 case MAT_FACTOR_SCHUR_INVERTED: 9712 PetscCall(MatMultTranspose(F->schur, rhs, sol)); 9713 break; 9714 default: 9715 SETERRQ(PetscObjectComm((PetscObject)F), PETSC_ERR_SUP, "Unhandled MatFactorSchurStatus %d", F->schur_status); 9716 } 9717 PetscFunctionReturn(PETSC_SUCCESS); 9718 } 9719 9720 /*@ 9721 MatFactorSolveSchurComplement - Solve the Schur complement system computed during the factorization step 9722 9723 Logically Collective 9724 9725 Input Parameters: 9726 + F - the factored matrix obtained by calling `MatGetFactor()` 9727 . rhs - location where the right hand side of the Schur complement system is stored 9728 - sol - location where the solution of the Schur complement system has to be returned 9729 9730 Level: advanced 9731 9732 Notes: 9733 The sizes of the vectors should match the size of the Schur complement 9734 9735 Must be called after `MatFactorSetSchurIS()` 9736 9737 .seealso: `Mat`, `MatGetFactor()`, `MatFactorSetSchurIS()`, `MatFactorSolveSchurComplementTranspose()` 9738 @*/ 9739 PetscErrorCode MatFactorSolveSchurComplement(Mat F, Vec rhs, Vec sol) 9740 { 9741 PetscFunctionBegin; 9742 PetscValidType(F, 1); 9743 PetscValidType(rhs, 2); 9744 PetscValidType(sol, 3); 9745 PetscValidHeaderSpecific(F, MAT_CLASSID, 1); 9746 PetscValidHeaderSpecific(rhs, VEC_CLASSID, 2); 9747 PetscValidHeaderSpecific(sol, VEC_CLASSID, 3); 9748 PetscCheckSameComm(F, 1, rhs, 2); 9749 PetscCheckSameComm(F, 1, sol, 3); 9750 PetscCall(MatFactorFactorizeSchurComplement(F)); 9751 switch (F->schur_status) { 9752 case MAT_FACTOR_SCHUR_FACTORED: 9753 PetscCall(MatSolve(F->schur, rhs, sol)); 9754 break; 9755 case MAT_FACTOR_SCHUR_INVERTED: 9756 PetscCall(MatMult(F->schur, rhs, sol)); 9757 break; 9758 default: 9759 SETERRQ(PetscObjectComm((PetscObject)F), PETSC_ERR_SUP, "Unhandled MatFactorSchurStatus %d", F->schur_status); 9760 } 9761 PetscFunctionReturn(PETSC_SUCCESS); 9762 } 9763 9764 PETSC_EXTERN PetscErrorCode MatSeqDenseInvertFactors_Private(Mat); 9765 #if PetscDefined(HAVE_CUDA) 9766 PETSC_EXTERN PetscErrorCode MatSeqDenseCUDAInvertFactors_Private(Mat); 9767 #endif 9768 9769 /* Schur status updated in the interface */ 9770 static PetscErrorCode MatFactorInvertSchurComplement_Private(Mat F) 9771 { 9772 Mat S = F->schur; 9773 9774 PetscFunctionBegin; 9775 if (S) { 9776 PetscMPIInt size; 9777 PetscBool isdense, isdensecuda; 9778 9779 PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)S), &size)); 9780 PetscCheck(size <= 1, PetscObjectComm((PetscObject)S), PETSC_ERR_SUP, "Not yet implemented"); 9781 PetscCall(PetscObjectTypeCompare((PetscObject)S, MATSEQDENSE, &isdense)); 9782 PetscCall(PetscObjectTypeCompare((PetscObject)S, MATSEQDENSECUDA, &isdensecuda)); 9783 PetscCheck(isdense || isdensecuda, PetscObjectComm((PetscObject)S), PETSC_ERR_SUP, "Not implemented for type %s", ((PetscObject)S)->type_name); 9784 PetscCall(PetscLogEventBegin(MAT_FactorInvS, F, 0, 0, 0)); 9785 if (isdense) { 9786 PetscCall(MatSeqDenseInvertFactors_Private(S)); 9787 } else if (isdensecuda) { 9788 #if defined(PETSC_HAVE_CUDA) 9789 PetscCall(MatSeqDenseCUDAInvertFactors_Private(S)); 9790 #endif 9791 } 9792 PetscCall(PetscLogEventEnd(MAT_FactorInvS, F, 0, 0, 0)); 9793 } 9794 PetscFunctionReturn(PETSC_SUCCESS); 9795 } 9796 9797 /*@ 9798 MatFactorInvertSchurComplement - Invert the Schur complement matrix computed during the factorization step 9799 9800 Logically Collective on F 9801 9802 Input Parameters: 9803 . F - the factored matrix obtained by calling `MatGetFactor()` 9804 9805 Level: advanced 9806 9807 Notes: 9808 Must be called after `MatFactorSetSchurIS()`. 9809 9810 Call `MatFactorGetSchurComplement()` or `MatFactorCreateSchurComplement()` AFTER this call to actually compute the inverse and get access to it. 9811 9812 .seealso: `Mat`, `MatGetFactor()`, `MatFactorSetSchurIS()`, `MatFactorGetSchurComplement()`, `MatFactorCreateSchurComplement()` 9813 @*/ 9814 PetscErrorCode MatFactorInvertSchurComplement(Mat F) 9815 { 9816 PetscFunctionBegin; 9817 PetscValidType(F, 1); 9818 PetscValidHeaderSpecific(F, MAT_CLASSID, 1); 9819 if (F->schur_status == MAT_FACTOR_SCHUR_INVERTED) PetscFunctionReturn(PETSC_SUCCESS); 9820 PetscCall(MatFactorFactorizeSchurComplement(F)); 9821 PetscCall(MatFactorInvertSchurComplement_Private(F)); 9822 F->schur_status = MAT_FACTOR_SCHUR_INVERTED; 9823 PetscFunctionReturn(PETSC_SUCCESS); 9824 } 9825 9826 /*@ 9827 MatFactorFactorizeSchurComplement - Factorize the Schur complement matrix computed during the factorization step 9828 9829 Logically Collective 9830 9831 Input Parameters: 9832 . F - the factored matrix obtained by calling `MatGetFactor()` 9833 9834 Level: advanced 9835 9836 Note: 9837 Must be called after `MatFactorSetSchurIS()` 9838 9839 .seealso: `Mat`, `MatGetFactor()`, `MatFactorSetSchurIS()`, `MatFactorInvertSchurComplement()` 9840 @*/ 9841 PetscErrorCode MatFactorFactorizeSchurComplement(Mat F) 9842 { 9843 MatFactorInfo info; 9844 9845 PetscFunctionBegin; 9846 PetscValidType(F, 1); 9847 PetscValidHeaderSpecific(F, MAT_CLASSID, 1); 9848 if (F->schur_status == MAT_FACTOR_SCHUR_INVERTED || F->schur_status == MAT_FACTOR_SCHUR_FACTORED) PetscFunctionReturn(PETSC_SUCCESS); 9849 PetscCall(PetscLogEventBegin(MAT_FactorFactS, F, 0, 0, 0)); 9850 if (F->factortype == MAT_FACTOR_CHOLESKY) { /* LDL^t regarded as Cholesky */ 9851 PetscCall(MatCholeskyFactor(F->schur, NULL, &info)); 9852 } else { 9853 PetscCall(MatLUFactor(F->schur, NULL, NULL, &info)); 9854 } 9855 PetscCall(PetscLogEventEnd(MAT_FactorFactS, F, 0, 0, 0)); 9856 F->schur_status = MAT_FACTOR_SCHUR_FACTORED; 9857 PetscFunctionReturn(PETSC_SUCCESS); 9858 } 9859 9860 /*@ 9861 MatPtAP - Creates the matrix product C = P^T * A * P 9862 9863 Neighbor-wise Collective 9864 9865 Input Parameters: 9866 + A - the matrix 9867 . P - the projection matrix 9868 . scall - either `MAT_INITIAL_MATRIX` or `MAT_REUSE_MATRIX` 9869 - fill - expected fill as ratio of nnz(C)/(nnz(A) + nnz(P)), use `PETSC_DEFAULT` if you do not have a good estimate 9870 if the result is a dense matrix this is irrelevant 9871 9872 Output Parameters: 9873 . C - the product matrix 9874 9875 Level: intermediate 9876 9877 Notes: 9878 C will be created and must be destroyed by the user with `MatDestroy()`. 9879 9880 An alternative approach to this function is to use `MatProductCreate()` and set the desired options before the computation is done 9881 9882 Developer Note: 9883 For matrix types without special implementation the function fallbacks to `MatMatMult()` followed by `MatTransposeMatMult()`. 9884 9885 .seealso: `Mat`, `MatProductCreate()`, `MatMatMult()`, `MatRARt()` 9886 @*/ 9887 PetscErrorCode MatPtAP(Mat A, Mat P, MatReuse scall, PetscReal fill, Mat *C) 9888 { 9889 PetscFunctionBegin; 9890 if (scall == MAT_REUSE_MATRIX) MatCheckProduct(*C, 5); 9891 PetscCheck(scall != MAT_INPLACE_MATRIX, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Inplace product not supported"); 9892 9893 if (scall == MAT_INITIAL_MATRIX) { 9894 PetscCall(MatProductCreate(A, P, NULL, C)); 9895 PetscCall(MatProductSetType(*C, MATPRODUCT_PtAP)); 9896 PetscCall(MatProductSetAlgorithm(*C, "default")); 9897 PetscCall(MatProductSetFill(*C, fill)); 9898 9899 (*C)->product->api_user = PETSC_TRUE; 9900 PetscCall(MatProductSetFromOptions(*C)); 9901 PetscCheck((*C)->ops->productsymbolic, PetscObjectComm((PetscObject)(*C)), PETSC_ERR_SUP, "MatProduct %s not supported for A %s and P %s", MatProductTypes[MATPRODUCT_PtAP], ((PetscObject)A)->type_name, ((PetscObject)P)->type_name); 9902 PetscCall(MatProductSymbolic(*C)); 9903 } else { /* scall == MAT_REUSE_MATRIX */ 9904 PetscCall(MatProductReplaceMats(A, P, NULL, *C)); 9905 } 9906 9907 PetscCall(MatProductNumeric(*C)); 9908 (*C)->symmetric = A->symmetric; 9909 (*C)->spd = A->spd; 9910 PetscFunctionReturn(PETSC_SUCCESS); 9911 } 9912 9913 /*@ 9914 MatRARt - Creates the matrix product C = R * A * R^T 9915 9916 Neighbor-wise Collective 9917 9918 Input Parameters: 9919 + A - the matrix 9920 . R - the projection matrix 9921 . scall - either `MAT_INITIAL_MATRIX` or `MAT_REUSE_MATRIX` 9922 - fill - expected fill as ratio of nnz(C)/nnz(A), use `PETSC_DEFAULT` if you do not have a good estimate 9923 if the result is a dense matrix this is irrelevant 9924 9925 Output Parameters: 9926 . C - the product matrix 9927 9928 Level: intermediate 9929 9930 Notes: 9931 C will be created and must be destroyed by the user with `MatDestroy()`. 9932 9933 An alternative approach to this function is to use `MatProductCreate()` and set the desired options before the computation is done 9934 9935 This routine is currently only implemented for pairs of `MATAIJ` matrices and classes 9936 which inherit from `MATAIJ`. Due to PETSc sparse matrix block row distribution among processes, 9937 parallel MatRARt is implemented via explicit transpose of R, which could be very expensive. 9938 We recommend using MatPtAP(). 9939 9940 .seealso: `Mat`, `MatProductCreate()`, `MatMatMult()`, `MatPtAP()` 9941 @*/ 9942 PetscErrorCode MatRARt(Mat A, Mat R, MatReuse scall, PetscReal fill, Mat *C) 9943 { 9944 PetscFunctionBegin; 9945 if (scall == MAT_REUSE_MATRIX) MatCheckProduct(*C, 5); 9946 PetscCheck(scall != MAT_INPLACE_MATRIX, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Inplace product not supported"); 9947 9948 if (scall == MAT_INITIAL_MATRIX) { 9949 PetscCall(MatProductCreate(A, R, NULL, C)); 9950 PetscCall(MatProductSetType(*C, MATPRODUCT_RARt)); 9951 PetscCall(MatProductSetAlgorithm(*C, "default")); 9952 PetscCall(MatProductSetFill(*C, fill)); 9953 9954 (*C)->product->api_user = PETSC_TRUE; 9955 PetscCall(MatProductSetFromOptions(*C)); 9956 PetscCheck((*C)->ops->productsymbolic, PetscObjectComm((PetscObject)(*C)), PETSC_ERR_SUP, "MatProduct %s not supported for A %s and R %s", MatProductTypes[MATPRODUCT_RARt], ((PetscObject)A)->type_name, ((PetscObject)R)->type_name); 9957 PetscCall(MatProductSymbolic(*C)); 9958 } else { /* scall == MAT_REUSE_MATRIX */ 9959 PetscCall(MatProductReplaceMats(A, R, NULL, *C)); 9960 } 9961 9962 PetscCall(MatProductNumeric(*C)); 9963 if (A->symmetric == PETSC_BOOL3_TRUE) PetscCall(MatSetOption(*C, MAT_SYMMETRIC, PETSC_TRUE)); 9964 PetscFunctionReturn(PETSC_SUCCESS); 9965 } 9966 9967 static PetscErrorCode MatProduct_Private(Mat A, Mat B, MatReuse scall, PetscReal fill, MatProductType ptype, Mat *C) 9968 { 9969 PetscFunctionBegin; 9970 PetscCheck(scall != MAT_INPLACE_MATRIX, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Inplace product not supported"); 9971 9972 if (scall == MAT_INITIAL_MATRIX) { 9973 PetscCall(PetscInfo(A, "Calling MatProduct API with MAT_INITIAL_MATRIX and product type %s\n", MatProductTypes[ptype])); 9974 PetscCall(MatProductCreate(A, B, NULL, C)); 9975 PetscCall(MatProductSetType(*C, ptype)); 9976 PetscCall(MatProductSetAlgorithm(*C, MATPRODUCTALGORITHMDEFAULT)); 9977 PetscCall(MatProductSetFill(*C, fill)); 9978 9979 (*C)->product->api_user = PETSC_TRUE; 9980 PetscCall(MatProductSetFromOptions(*C)); 9981 PetscCall(MatProductSymbolic(*C)); 9982 } else { /* scall == MAT_REUSE_MATRIX */ 9983 Mat_Product *product = (*C)->product; 9984 PetscBool isdense; 9985 9986 PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)(*C), &isdense, MATSEQDENSE, MATMPIDENSE, "")); 9987 if (isdense && product && product->type != ptype) { 9988 PetscCall(MatProductClear(*C)); 9989 product = NULL; 9990 } 9991 PetscCall(PetscInfo(A, "Calling MatProduct API with MAT_REUSE_MATRIX %s product present and product type %s\n", product ? "with" : "without", MatProductTypes[ptype])); 9992 if (!product) { /* user provide the dense matrix *C without calling MatProductCreate() or reusing it from previous calls */ 9993 PetscCheck(isdense, PetscObjectComm((PetscObject)(*C)), PETSC_ERR_SUP, "Call MatProductCreate() first"); 9994 PetscCall(MatProductCreate_Private(A, B, NULL, *C)); 9995 product = (*C)->product; 9996 product->fill = fill; 9997 product->api_user = PETSC_TRUE; 9998 product->clear = PETSC_TRUE; 9999 10000 PetscCall(MatProductSetType(*C, ptype)); 10001 PetscCall(MatProductSetFromOptions(*C)); 10002 PetscCheck((*C)->ops->productsymbolic, PetscObjectComm((PetscObject)(*C)), PETSC_ERR_SUP, "MatProduct %s not supported for %s and %s", MatProductTypes[ptype], ((PetscObject)A)->type_name, ((PetscObject)B)->type_name); 10003 PetscCall(MatProductSymbolic(*C)); 10004 } else { /* user may change input matrices A or B when REUSE */ 10005 PetscCall(MatProductReplaceMats(A, B, NULL, *C)); 10006 } 10007 } 10008 PetscCall(MatProductNumeric(*C)); 10009 PetscFunctionReturn(PETSC_SUCCESS); 10010 } 10011 10012 /*@ 10013 MatMatMult - Performs matrix-matrix multiplication C=A*B. 10014 10015 Neighbor-wise Collective 10016 10017 Input Parameters: 10018 + A - the left matrix 10019 . B - the right matrix 10020 . scall - either `MAT_INITIAL_MATRIX` or `MAT_REUSE_MATRIX` 10021 - fill - expected fill as ratio of nnz(C)/(nnz(A) + nnz(B)), use `PETSC_DEFAULT` if you do not have a good estimate 10022 if the result is a dense matrix this is irrelevant 10023 10024 Output Parameters: 10025 . C - the product matrix 10026 10027 Notes: 10028 Unless scall is `MAT_REUSE_MATRIX` C will be created. 10029 10030 `MAT_REUSE_MATRIX` can only be used if the matrices A and B have the same nonzero pattern as in the previous call and C was obtained from a previous 10031 call to this function with `MAT_INITIAL_MATRIX`. 10032 10033 To determine the correct fill value, run with -info and search for the string "Fill ratio" to see the value actually needed. 10034 10035 In the special case where matrix B (and hence C) are dense you can create the correctly sized matrix C yourself and then call this routine with `MAT_REUSE_MATRIX`, 10036 rather than first having `MatMatMult()` create it for you. You can NEVER do this if the matrix C is sparse. 10037 10038 Example of Usage: 10039 .vb 10040 MatProductCreate(A,B,NULL,&C); 10041 MatProductSetType(C,MATPRODUCT_AB); 10042 MatProductSymbolic(C); 10043 MatProductNumeric(C); // compute C=A * B 10044 MatProductReplaceMats(A1,B1,NULL,C); // compute C=A1 * B1 10045 MatProductNumeric(C); 10046 MatProductReplaceMats(A2,NULL,NULL,C); // compute C=A2 * B1 10047 MatProductNumeric(C); 10048 .ve 10049 10050 Level: intermediate 10051 10052 .seealso: `Mat`, `MatProductType`, `MATPRODUCT_AB`, `MatTransposeMatMult()`, `MatMatTransposeMult()`, `MatPtAP()`, `MatProductCreate()`, `MatProductSymbolic()`, `MatProductReplaceMats()`, `MatProductNumeric()` 10053 @*/ 10054 PetscErrorCode MatMatMult(Mat A, Mat B, MatReuse scall, PetscReal fill, Mat *C) 10055 { 10056 PetscFunctionBegin; 10057 PetscCall(MatProduct_Private(A, B, scall, fill, MATPRODUCT_AB, C)); 10058 PetscFunctionReturn(PETSC_SUCCESS); 10059 } 10060 10061 /*@ 10062 MatMatTransposeMult - Performs matrix-matrix multiplication C=A*B^T. 10063 10064 Neighbor-wise Collective 10065 10066 Input Parameters: 10067 + A - the left matrix 10068 . B - the right matrix 10069 . scall - either `MAT_INITIAL_MATRIX` or `MAT_REUSE_MATRIX` 10070 - fill - expected fill as ratio of nnz(C)/(nnz(A) + nnz(B)), use `PETSC_DEFAULT` if not known 10071 10072 Output Parameters: 10073 . C - the product matrix 10074 10075 Level: intermediate 10076 10077 Notes: 10078 C will be created if `MAT_INITIAL_MATRIX` and must be destroyed by the user with `MatDestroy()`. 10079 10080 `MAT_REUSE_MATRIX` can only be used if the matrices A and B have the same nonzero pattern as in the previous call 10081 10082 To determine the correct fill value, run with -info and search for the string "Fill ratio" to see the value 10083 actually needed. 10084 10085 This routine is currently only implemented for pairs of `MATSEQAIJ` matrices, for the `MATSEQDENSE` class, 10086 and for pairs of `MATMPIDENSE` matrices. 10087 10088 This routine is shorthand for using `MatProductCreate()` with the `MatProductType` of `MATPRODUCT_ABt` 10089 10090 Options Database Keys: 10091 . -matmattransmult_mpidense_mpidense_via {allgatherv,cyclic} - Choose between algorithms for `MATMPIDENSE` matrices: the 10092 first redundantly copies the transposed B matrix on each process and requiers O(log P) communication complexity; 10093 the second never stores more than one portion of the B matrix at a time by requires O(P) communication complexity. 10094 10095 .seealso: `Mat`, `MatProductCreate()`, `MATPRODUCT_ABt`, `MatMatMult()`, `MatTransposeMatMult()` `MatPtAP()`, `MatProductCreate()`, `MatProductAlgorithm`, `MatProductType`, `MATPRODUCT_ABt` 10096 @*/ 10097 PetscErrorCode MatMatTransposeMult(Mat A, Mat B, MatReuse scall, PetscReal fill, Mat *C) 10098 { 10099 PetscFunctionBegin; 10100 PetscCall(MatProduct_Private(A, B, scall, fill, MATPRODUCT_ABt, C)); 10101 if (A == B) PetscCall(MatSetOption(*C, MAT_SYMMETRIC, PETSC_TRUE)); 10102 PetscFunctionReturn(PETSC_SUCCESS); 10103 } 10104 10105 /*@ 10106 MatTransposeMatMult - Performs matrix-matrix multiplication C=A^T*B. 10107 10108 Neighbor-wise Collective 10109 10110 Input Parameters: 10111 + A - the left matrix 10112 . B - the right matrix 10113 . scall - either `MAT_INITIAL_MATRIX` or `MAT_REUSE_MATRIX` 10114 - fill - expected fill as ratio of nnz(C)/(nnz(A) + nnz(B)), use `PETSC_DEFAULT` if not known 10115 10116 Output Parameters: 10117 . C - the product matrix 10118 10119 Level: intermediate 10120 10121 Notes: 10122 C will be created if `MAT_INITIAL_MATRIX` and must be destroyed by the user with `MatDestroy()`. 10123 10124 `MAT_REUSE_MATRIX` can only be used if the matrices A and B have the same nonzero pattern as in the previous call. 10125 10126 This routine is shorthand for using `MatProductCreate()` with the `MatProductType` of `MATPRODUCT_AtB` 10127 10128 To determine the correct fill value, run with -info and search for the string "Fill ratio" to see the value 10129 actually needed. 10130 10131 This routine is currently implemented for pairs of `MATAIJ` matrices and pairs of `MATSEQDENSE` matrices and classes 10132 which inherit from `MATSEQAIJ`. C will be of the same type as the input matrices. 10133 10134 .seealso: `Mat`, `MatProductCreate()`, `MATPRODUCT_AtB`, `MatMatMult()`, `MatMatTransposeMult()`, `MatPtAP()` 10135 @*/ 10136 PetscErrorCode MatTransposeMatMult(Mat A, Mat B, MatReuse scall, PetscReal fill, Mat *C) 10137 { 10138 PetscFunctionBegin; 10139 PetscCall(MatProduct_Private(A, B, scall, fill, MATPRODUCT_AtB, C)); 10140 PetscFunctionReturn(PETSC_SUCCESS); 10141 } 10142 10143 /*@ 10144 MatMatMatMult - Performs matrix-matrix-matrix multiplication D=A*B*C. 10145 10146 Neighbor-wise Collective 10147 10148 Input Parameters: 10149 + A - the left matrix 10150 . B - the middle matrix 10151 . C - the right matrix 10152 . scall - either `MAT_INITIAL_MATRIX` or `MAT_REUSE_MATRIX` 10153 - fill - expected fill as ratio of nnz(D)/(nnz(A) + nnz(B)+nnz(C)), use `PETSC_DEFAULT` if you do not have a good estimate 10154 if the result is a dense matrix this is irrelevant 10155 10156 Output Parameters: 10157 . D - the product matrix 10158 10159 Level: intermediate 10160 10161 Notes: 10162 Unless scall is `MAT_REUSE_MATRIX` D will be created. 10163 10164 `MAT_REUSE_MATRIX` can only be used if the matrices A, B and C have the same nonzero pattern as in the previous call 10165 10166 This routine is shorthand for using `MatProductCreate()` with the `MatProductType` of `MATPRODUCT_ABC` 10167 10168 To determine the correct fill value, run with -info and search for the string "Fill ratio" to see the value 10169 actually needed. 10170 10171 If you have many matrices with the same non-zero structure to multiply, you 10172 should use `MAT_REUSE_MATRIX` in all calls but the first 10173 10174 .seealso: `Mat`, `MatProductCreate()`, `MATPRODUCT_ABC`, `MatMatMult`, `MatPtAP()`, `MatMatTransposeMult()`, `MatTransposeMatMult()` 10175 @*/ 10176 PetscErrorCode MatMatMatMult(Mat A, Mat B, Mat C, MatReuse scall, PetscReal fill, Mat *D) 10177 { 10178 PetscFunctionBegin; 10179 if (scall == MAT_REUSE_MATRIX) MatCheckProduct(*D, 6); 10180 PetscCheck(scall != MAT_INPLACE_MATRIX, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Inplace product not supported"); 10181 10182 if (scall == MAT_INITIAL_MATRIX) { 10183 PetscCall(MatProductCreate(A, B, C, D)); 10184 PetscCall(MatProductSetType(*D, MATPRODUCT_ABC)); 10185 PetscCall(MatProductSetAlgorithm(*D, "default")); 10186 PetscCall(MatProductSetFill(*D, fill)); 10187 10188 (*D)->product->api_user = PETSC_TRUE; 10189 PetscCall(MatProductSetFromOptions(*D)); 10190 PetscCheck((*D)->ops->productsymbolic, PetscObjectComm((PetscObject)(*D)), PETSC_ERR_SUP, "MatProduct %s not supported for A %s, B %s and C %s", MatProductTypes[MATPRODUCT_ABC], ((PetscObject)A)->type_name, ((PetscObject)B)->type_name, 10191 ((PetscObject)C)->type_name); 10192 PetscCall(MatProductSymbolic(*D)); 10193 } else { /* user may change input matrices when REUSE */ 10194 PetscCall(MatProductReplaceMats(A, B, C, *D)); 10195 } 10196 PetscCall(MatProductNumeric(*D)); 10197 PetscFunctionReturn(PETSC_SUCCESS); 10198 } 10199 10200 /*@ 10201 MatCreateRedundantMatrix - Create redundant matrices and put them into processors of subcommunicators. 10202 10203 Collective 10204 10205 Input Parameters: 10206 + mat - the matrix 10207 . nsubcomm - the number of subcommunicators (= number of redundant parallel or sequential matrices) 10208 . subcomm - MPI communicator split from the communicator where mat resides in (or `MPI_COMM_NULL` if nsubcomm is used) 10209 - reuse - either `MAT_INITIAL_MATRIX` or `MAT_REUSE_MATRIX` 10210 10211 Output Parameter: 10212 . matredundant - redundant matrix 10213 10214 Level: advanced 10215 10216 Notes: 10217 `MAT_REUSE_MATRIX` can only be used when the nonzero structure of the 10218 original matrix has not changed from that last call to MatCreateRedundantMatrix(). 10219 10220 This routine creates the duplicated matrices in the subcommunicators; you should NOT create them before 10221 calling it. 10222 10223 `PetscSubcommCreate()` can be used to manage the creation of the subcomm but need not be. 10224 10225 .seealso: `Mat`, `MatDestroy()`, `PetscSubcommCreate()`, `PetscSubComm` 10226 @*/ 10227 PetscErrorCode MatCreateRedundantMatrix(Mat mat, PetscInt nsubcomm, MPI_Comm subcomm, MatReuse reuse, Mat *matredundant) 10228 { 10229 MPI_Comm comm; 10230 PetscMPIInt size; 10231 PetscInt mloc_sub, nloc_sub, rstart, rend, M = mat->rmap->N, N = mat->cmap->N, bs = mat->rmap->bs; 10232 Mat_Redundant *redund = NULL; 10233 PetscSubcomm psubcomm = NULL; 10234 MPI_Comm subcomm_in = subcomm; 10235 Mat *matseq; 10236 IS isrow, iscol; 10237 PetscBool newsubcomm = PETSC_FALSE; 10238 10239 PetscFunctionBegin; 10240 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 10241 if (nsubcomm && reuse == MAT_REUSE_MATRIX) { 10242 PetscValidPointer(*matredundant, 5); 10243 PetscValidHeaderSpecific(*matredundant, MAT_CLASSID, 5); 10244 } 10245 10246 PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)mat), &size)); 10247 if (size == 1 || nsubcomm == 1) { 10248 if (reuse == MAT_INITIAL_MATRIX) { 10249 PetscCall(MatDuplicate(mat, MAT_COPY_VALUES, matredundant)); 10250 } else { 10251 PetscCheck(*matredundant != mat, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "MAT_REUSE_MATRIX means reuse the matrix passed in as the final argument, not the original matrix"); 10252 PetscCall(MatCopy(mat, *matredundant, SAME_NONZERO_PATTERN)); 10253 } 10254 PetscFunctionReturn(PETSC_SUCCESS); 10255 } 10256 10257 PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 10258 PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 10259 MatCheckPreallocated(mat, 1); 10260 10261 PetscCall(PetscLogEventBegin(MAT_RedundantMat, mat, 0, 0, 0)); 10262 if (subcomm_in == MPI_COMM_NULL && reuse == MAT_INITIAL_MATRIX) { /* get subcomm if user does not provide subcomm */ 10263 /* create psubcomm, then get subcomm */ 10264 PetscCall(PetscObjectGetComm((PetscObject)mat, &comm)); 10265 PetscCallMPI(MPI_Comm_size(comm, &size)); 10266 PetscCheck(nsubcomm >= 1 && nsubcomm <= size, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "nsubcomm must between 1 and %d", size); 10267 10268 PetscCall(PetscSubcommCreate(comm, &psubcomm)); 10269 PetscCall(PetscSubcommSetNumber(psubcomm, nsubcomm)); 10270 PetscCall(PetscSubcommSetType(psubcomm, PETSC_SUBCOMM_CONTIGUOUS)); 10271 PetscCall(PetscSubcommSetFromOptions(psubcomm)); 10272 PetscCall(PetscCommDuplicate(PetscSubcommChild(psubcomm), &subcomm, NULL)); 10273 newsubcomm = PETSC_TRUE; 10274 PetscCall(PetscSubcommDestroy(&psubcomm)); 10275 } 10276 10277 /* get isrow, iscol and a local sequential matrix matseq[0] */ 10278 if (reuse == MAT_INITIAL_MATRIX) { 10279 mloc_sub = PETSC_DECIDE; 10280 nloc_sub = PETSC_DECIDE; 10281 if (bs < 1) { 10282 PetscCall(PetscSplitOwnership(subcomm, &mloc_sub, &M)); 10283 PetscCall(PetscSplitOwnership(subcomm, &nloc_sub, &N)); 10284 } else { 10285 PetscCall(PetscSplitOwnershipBlock(subcomm, bs, &mloc_sub, &M)); 10286 PetscCall(PetscSplitOwnershipBlock(subcomm, bs, &nloc_sub, &N)); 10287 } 10288 PetscCallMPI(MPI_Scan(&mloc_sub, &rend, 1, MPIU_INT, MPI_SUM, subcomm)); 10289 rstart = rend - mloc_sub; 10290 PetscCall(ISCreateStride(PETSC_COMM_SELF, mloc_sub, rstart, 1, &isrow)); 10291 PetscCall(ISCreateStride(PETSC_COMM_SELF, N, 0, 1, &iscol)); 10292 } else { /* reuse == MAT_REUSE_MATRIX */ 10293 PetscCheck(*matredundant != mat, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "MAT_REUSE_MATRIX means reuse the matrix passed in as the final argument, not the original matrix"); 10294 /* retrieve subcomm */ 10295 PetscCall(PetscObjectGetComm((PetscObject)(*matredundant), &subcomm)); 10296 redund = (*matredundant)->redundant; 10297 isrow = redund->isrow; 10298 iscol = redund->iscol; 10299 matseq = redund->matseq; 10300 } 10301 PetscCall(MatCreateSubMatrices(mat, 1, &isrow, &iscol, reuse, &matseq)); 10302 10303 /* get matredundant over subcomm */ 10304 if (reuse == MAT_INITIAL_MATRIX) { 10305 PetscCall(MatCreateMPIMatConcatenateSeqMat(subcomm, matseq[0], nloc_sub, reuse, matredundant)); 10306 10307 /* create a supporting struct and attach it to C for reuse */ 10308 PetscCall(PetscNew(&redund)); 10309 (*matredundant)->redundant = redund; 10310 redund->isrow = isrow; 10311 redund->iscol = iscol; 10312 redund->matseq = matseq; 10313 if (newsubcomm) { 10314 redund->subcomm = subcomm; 10315 } else { 10316 redund->subcomm = MPI_COMM_NULL; 10317 } 10318 } else { 10319 PetscCall(MatCreateMPIMatConcatenateSeqMat(subcomm, matseq[0], PETSC_DECIDE, reuse, matredundant)); 10320 } 10321 #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA) || defined(PETSC_HAVE_HIP) 10322 if (matseq[0]->boundtocpu && matseq[0]->bindingpropagates) { 10323 PetscCall(MatBindToCPU(*matredundant, PETSC_TRUE)); 10324 PetscCall(MatSetBindingPropagates(*matredundant, PETSC_TRUE)); 10325 } 10326 #endif 10327 PetscCall(PetscLogEventEnd(MAT_RedundantMat, mat, 0, 0, 0)); 10328 PetscFunctionReturn(PETSC_SUCCESS); 10329 } 10330 10331 /*@C 10332 MatGetMultiProcBlock - Create multiple 'parallel submatrices' from 10333 a given `Mat`. Each submatrix can span multiple procs. 10334 10335 Collective 10336 10337 Input Parameters: 10338 + mat - the matrix 10339 . subcomm - the sub communicator obtained as if by `MPI_Comm_split(PetscObjectComm((PetscObject)mat))` 10340 - scall - either `MAT_INITIAL_MATRIX` or `MAT_REUSE_MATRIX` 10341 10342 Output Parameter: 10343 . subMat - parallel sub-matrices each spanning a given `subcomm` 10344 10345 Level: advanced 10346 10347 Notes: 10348 The submatrix partition across processors is dictated by `subComm` a 10349 communicator obtained by `MPI_comm_split()` or via `PetscSubcommCreate()`. The `subComm` 10350 is not restricted to be grouped with consecutive original ranks. 10351 10352 Due the `MPI_Comm_split()` usage, the parallel layout of the submatrices 10353 map directly to the layout of the original matrix [wrt the local 10354 row,col partitioning]. So the original 'DiagonalMat' naturally maps 10355 into the 'DiagonalMat' of the `subMat`, hence it is used directly from 10356 the `subMat`. However the offDiagMat looses some columns - and this is 10357 reconstructed with `MatSetValues()` 10358 10359 This is used by `PCBJACOBI` when a single block spans multiple MPI ranks 10360 10361 .seealso: `Mat`, `MatCreateRedundantMatrix()`, `MatCreateSubMatrices()`, `PCBJACOBI` 10362 @*/ 10363 PetscErrorCode MatGetMultiProcBlock(Mat mat, MPI_Comm subComm, MatReuse scall, Mat *subMat) 10364 { 10365 PetscMPIInt commsize, subCommSize; 10366 10367 PetscFunctionBegin; 10368 PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)mat), &commsize)); 10369 PetscCallMPI(MPI_Comm_size(subComm, &subCommSize)); 10370 PetscCheck(subCommSize <= commsize, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_OUTOFRANGE, "CommSize %d < SubCommZize %d", commsize, subCommSize); 10371 10372 PetscCheck(scall != MAT_REUSE_MATRIX || *subMat != mat, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "MAT_REUSE_MATRIX means reuse the matrix passed in as the final argument, not the original matrix"); 10373 PetscCall(PetscLogEventBegin(MAT_GetMultiProcBlock, mat, 0, 0, 0)); 10374 PetscUseTypeMethod(mat, getmultiprocblock, subComm, scall, subMat); 10375 PetscCall(PetscLogEventEnd(MAT_GetMultiProcBlock, mat, 0, 0, 0)); 10376 PetscFunctionReturn(PETSC_SUCCESS); 10377 } 10378 10379 /*@ 10380 MatGetLocalSubMatrix - Gets a reference to a submatrix specified in local numbering 10381 10382 Not Collective 10383 10384 Input Parameters: 10385 + mat - matrix to extract local submatrix from 10386 . isrow - local row indices for submatrix 10387 - iscol - local column indices for submatrix 10388 10389 Output Parameter: 10390 . submat - the submatrix 10391 10392 Level: intermediate 10393 10394 Notes: 10395 `submat` should be disposed of with `MatRestoreLocalSubMatrix()`. 10396 10397 Depending on the format of `mat`, the returned submat may not implement `MatMult()`. Its communicator may be 10398 the same as mat, it may be `PETSC_COMM_SELF`, or some other subcomm of `mat`'s. 10399 10400 `submat` always implements `MatSetValuesLocal()`. If `isrow` and `iscol` have the same block size, then 10401 `MatSetValuesBlockedLocal()` will also be implemented. 10402 10403 `mat` must have had a `ISLocalToGlobalMapping` provided to it with `MatSetLocalToGlobalMapping()`. 10404 Matrices obtained with `DMCreateMatrix()` generally already have the local to global mapping provided. 10405 10406 .seealso: `Mat`, `MatRestoreLocalSubMatrix()`, `MatCreateLocalRef()`, `MatSetLocalToGlobalMapping()` 10407 @*/ 10408 PetscErrorCode MatGetLocalSubMatrix(Mat mat, IS isrow, IS iscol, Mat *submat) 10409 { 10410 PetscFunctionBegin; 10411 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 10412 PetscValidHeaderSpecific(isrow, IS_CLASSID, 2); 10413 PetscValidHeaderSpecific(iscol, IS_CLASSID, 3); 10414 PetscCheckSameComm(isrow, 2, iscol, 3); 10415 PetscValidPointer(submat, 4); 10416 PetscCheck(mat->rmap->mapping, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Matrix must have local to global mapping provided before this call"); 10417 10418 if (mat->ops->getlocalsubmatrix) { 10419 PetscUseTypeMethod(mat, getlocalsubmatrix, isrow, iscol, submat); 10420 } else { 10421 PetscCall(MatCreateLocalRef(mat, isrow, iscol, submat)); 10422 } 10423 PetscFunctionReturn(PETSC_SUCCESS); 10424 } 10425 10426 /*@ 10427 MatRestoreLocalSubMatrix - Restores a reference to a submatrix specified in local numbering obtained with `MatGetLocalSubMatrix()` 10428 10429 Not Collective 10430 10431 Input Parameters: 10432 + mat - matrix to extract local submatrix from 10433 . isrow - local row indices for submatrix 10434 . iscol - local column indices for submatrix 10435 - submat - the submatrix 10436 10437 Level: intermediate 10438 10439 .seealso: `Mat`, `MatGetLocalSubMatrix()` 10440 @*/ 10441 PetscErrorCode MatRestoreLocalSubMatrix(Mat mat, IS isrow, IS iscol, Mat *submat) 10442 { 10443 PetscFunctionBegin; 10444 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 10445 PetscValidHeaderSpecific(isrow, IS_CLASSID, 2); 10446 PetscValidHeaderSpecific(iscol, IS_CLASSID, 3); 10447 PetscCheckSameComm(isrow, 2, iscol, 3); 10448 PetscValidPointer(submat, 4); 10449 if (*submat) PetscValidHeaderSpecific(*submat, MAT_CLASSID, 4); 10450 10451 if (mat->ops->restorelocalsubmatrix) { 10452 PetscUseTypeMethod(mat, restorelocalsubmatrix, isrow, iscol, submat); 10453 } else { 10454 PetscCall(MatDestroy(submat)); 10455 } 10456 *submat = NULL; 10457 PetscFunctionReturn(PETSC_SUCCESS); 10458 } 10459 10460 /*@ 10461 MatFindZeroDiagonals - Finds all the rows of a matrix that have zero or no diagonal entry in the matrix 10462 10463 Collective 10464 10465 Input Parameter: 10466 . mat - the matrix 10467 10468 Output Parameter: 10469 . is - if any rows have zero diagonals this contains the list of them 10470 10471 Level: developer 10472 10473 .seealso: `Mat`, `MatMultTranspose()`, `MatMultAdd()`, `MatMultTransposeAdd()` 10474 @*/ 10475 PetscErrorCode MatFindZeroDiagonals(Mat mat, IS *is) 10476 { 10477 PetscFunctionBegin; 10478 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 10479 PetscValidType(mat, 1); 10480 PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 10481 PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 10482 10483 if (!mat->ops->findzerodiagonals) { 10484 Vec diag; 10485 const PetscScalar *a; 10486 PetscInt *rows; 10487 PetscInt rStart, rEnd, r, nrow = 0; 10488 10489 PetscCall(MatCreateVecs(mat, &diag, NULL)); 10490 PetscCall(MatGetDiagonal(mat, diag)); 10491 PetscCall(MatGetOwnershipRange(mat, &rStart, &rEnd)); 10492 PetscCall(VecGetArrayRead(diag, &a)); 10493 for (r = 0; r < rEnd - rStart; ++r) 10494 if (a[r] == 0.0) ++nrow; 10495 PetscCall(PetscMalloc1(nrow, &rows)); 10496 nrow = 0; 10497 for (r = 0; r < rEnd - rStart; ++r) 10498 if (a[r] == 0.0) rows[nrow++] = r + rStart; 10499 PetscCall(VecRestoreArrayRead(diag, &a)); 10500 PetscCall(VecDestroy(&diag)); 10501 PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)mat), nrow, rows, PETSC_OWN_POINTER, is)); 10502 } else { 10503 PetscUseTypeMethod(mat, findzerodiagonals, is); 10504 } 10505 PetscFunctionReturn(PETSC_SUCCESS); 10506 } 10507 10508 /*@ 10509 MatFindOffBlockDiagonalEntries - Finds all the rows of a matrix that have entries outside of the main diagonal block (defined by the matrix block size) 10510 10511 Collective 10512 10513 Input Parameter: 10514 . mat - the matrix 10515 10516 Output Parameter: 10517 . is - contains the list of rows with off block diagonal entries 10518 10519 Level: developer 10520 10521 .seealso: `Mat`, `MatMultTranspose()`, `MatMultAdd()`, `MatMultTransposeAdd()` 10522 @*/ 10523 PetscErrorCode MatFindOffBlockDiagonalEntries(Mat mat, IS *is) 10524 { 10525 PetscFunctionBegin; 10526 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 10527 PetscValidType(mat, 1); 10528 PetscCheck(mat->assembled, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 10529 PetscCheck(!mat->factortype, PetscObjectComm((PetscObject)mat), PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 10530 10531 PetscUseTypeMethod(mat, findoffblockdiagonalentries, is); 10532 PetscFunctionReturn(PETSC_SUCCESS); 10533 } 10534 10535 /*@C 10536 MatInvertBlockDiagonal - Inverts the block diagonal entries. 10537 10538 Collective; No Fortran Support 10539 10540 Input Parameters: 10541 . mat - the matrix 10542 10543 Output Parameters: 10544 . values - the block inverses in column major order (FORTRAN-like) 10545 10546 Level: advanced 10547 10548 Notes: 10549 The size of the blocks is determined by the block size of the matrix. 10550 10551 The blocks never overlap between two MPI ranks, use `MatInvertVariableBlockEnvelope()` for that case 10552 10553 The blocks all have the same size, use `MatInvertVariableBlockDiagonal()` for variable block size 10554 10555 .seealso: `Mat`, `MatInvertVariableBlockEnvelope()`, `MatInvertBlockDiagonalMat()` 10556 @*/ 10557 PetscErrorCode MatInvertBlockDiagonal(Mat mat, const PetscScalar **values) 10558 { 10559 PetscFunctionBegin; 10560 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 10561 PetscCheck(mat->assembled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 10562 PetscCheck(!mat->factortype, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 10563 PetscUseTypeMethod(mat, invertblockdiagonal, values); 10564 PetscFunctionReturn(PETSC_SUCCESS); 10565 } 10566 10567 /*@C 10568 MatInvertVariableBlockDiagonal - Inverts the point block diagonal entries. 10569 10570 Collective; No Fortran Support 10571 10572 Input Parameters: 10573 + mat - the matrix 10574 . nblocks - the number of blocks on the process, set with `MatSetVariableBlockSizes()` 10575 - bsizes - the size of each block on the process, set with `MatSetVariableBlockSizes()` 10576 10577 Output Parameters: 10578 . values - the block inverses in column major order (FORTRAN-like) 10579 10580 Level: advanced 10581 10582 Notes: 10583 Use `MatInvertBlockDiagonal()` if all blocks have the same size 10584 10585 The blocks never overlap between two MPI ranks, use `MatInvertVariableBlockEnvelope()` for that case 10586 10587 .seealso: `Mat`, `MatInvertBlockDiagonal()`, `MatSetVariableBlockSizes()`, `MatInvertVariableBlockEnvelope()` 10588 @*/ 10589 PetscErrorCode MatInvertVariableBlockDiagonal(Mat mat, PetscInt nblocks, const PetscInt *bsizes, PetscScalar *values) 10590 { 10591 PetscFunctionBegin; 10592 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 10593 PetscCheck(mat->assembled, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix"); 10594 PetscCheck(!mat->factortype, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for factored matrix"); 10595 PetscUseTypeMethod(mat, invertvariableblockdiagonal, nblocks, bsizes, values); 10596 PetscFunctionReturn(PETSC_SUCCESS); 10597 } 10598 10599 /*@ 10600 MatInvertBlockDiagonalMat - set the values of matrix C to be the inverted block diagonal of matrix A 10601 10602 Collective on Mat 10603 10604 Input Parameters: 10605 + A - the matrix 10606 - C - matrix with inverted block diagonal of `A`. This matrix should be created and may have its type set. 10607 10608 Level: advanced 10609 10610 Note: 10611 The blocksize of the matrix is used to determine the blocks on the diagonal of `C` 10612 10613 .seealso: `Mat`, `MatInvertBlockDiagonal()` 10614 @*/ 10615 PetscErrorCode MatInvertBlockDiagonalMat(Mat A, Mat C) 10616 { 10617 const PetscScalar *vals; 10618 PetscInt *dnnz; 10619 PetscInt m, rstart, rend, bs, i, j; 10620 10621 PetscFunctionBegin; 10622 PetscCall(MatInvertBlockDiagonal(A, &vals)); 10623 PetscCall(MatGetBlockSize(A, &bs)); 10624 PetscCall(MatGetLocalSize(A, &m, NULL)); 10625 PetscCall(MatSetLayouts(C, A->rmap, A->cmap)); 10626 PetscCall(PetscMalloc1(m / bs, &dnnz)); 10627 for (j = 0; j < m / bs; j++) dnnz[j] = 1; 10628 PetscCall(MatXAIJSetPreallocation(C, bs, dnnz, NULL, NULL, NULL)); 10629 PetscCall(PetscFree(dnnz)); 10630 PetscCall(MatGetOwnershipRange(C, &rstart, &rend)); 10631 PetscCall(MatSetOption(C, MAT_ROW_ORIENTED, PETSC_FALSE)); 10632 for (i = rstart / bs; i < rend / bs; i++) PetscCall(MatSetValuesBlocked(C, 1, &i, 1, &i, &vals[(i - rstart / bs) * bs * bs], INSERT_VALUES)); 10633 PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY)); 10634 PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY)); 10635 PetscCall(MatSetOption(C, MAT_ROW_ORIENTED, PETSC_TRUE)); 10636 PetscFunctionReturn(PETSC_SUCCESS); 10637 } 10638 10639 /*@C 10640 MatTransposeColoringDestroy - Destroys a coloring context for matrix product C=A*B^T that was created 10641 via `MatTransposeColoringCreate()`. 10642 10643 Collective on c 10644 10645 Input Parameter: 10646 . c - coloring context 10647 10648 Level: intermediate 10649 10650 .seealso: `Mat`, `MatTransposeColoringCreate()` 10651 @*/ 10652 PetscErrorCode MatTransposeColoringDestroy(MatTransposeColoring *c) 10653 { 10654 MatTransposeColoring matcolor = *c; 10655 10656 PetscFunctionBegin; 10657 if (!matcolor) PetscFunctionReturn(PETSC_SUCCESS); 10658 if (--((PetscObject)matcolor)->refct > 0) { 10659 matcolor = NULL; 10660 PetscFunctionReturn(PETSC_SUCCESS); 10661 } 10662 10663 PetscCall(PetscFree3(matcolor->ncolumns, matcolor->nrows, matcolor->colorforrow)); 10664 PetscCall(PetscFree(matcolor->rows)); 10665 PetscCall(PetscFree(matcolor->den2sp)); 10666 PetscCall(PetscFree(matcolor->colorforcol)); 10667 PetscCall(PetscFree(matcolor->columns)); 10668 if (matcolor->brows > 0) PetscCall(PetscFree(matcolor->lstart)); 10669 PetscCall(PetscHeaderDestroy(c)); 10670 PetscFunctionReturn(PETSC_SUCCESS); 10671 } 10672 10673 /*@C 10674 MatTransColoringApplySpToDen - Given a symbolic matrix product C=A*B^T for which 10675 a `MatTransposeColoring` context has been created, computes a dense B^T by applying 10676 `MatTransposeColoring` to sparse B. 10677 10678 Collective on coloring 10679 10680 Input Parameters: 10681 + B - sparse matrix B 10682 . Btdense - symbolic dense matrix B^T 10683 - coloring - coloring context created with `MatTransposeColoringCreate()` 10684 10685 Output Parameter: 10686 . Btdense - dense matrix B^T 10687 10688 Level: developer 10689 10690 Note: 10691 These are used internally for some implementations of `MatRARt()` 10692 10693 .seealso: `Mat`, `MatTransposeColoringCreate()`, `MatTransposeColoringDestroy()`, `MatTransColoringApplyDenToSp()` 10694 @*/ 10695 PetscErrorCode MatTransColoringApplySpToDen(MatTransposeColoring coloring, Mat B, Mat Btdense) 10696 { 10697 PetscFunctionBegin; 10698 PetscValidHeaderSpecific(B, MAT_CLASSID, 2); 10699 PetscValidHeaderSpecific(Btdense, MAT_CLASSID, 3); 10700 PetscValidHeaderSpecific(coloring, MAT_TRANSPOSECOLORING_CLASSID, 1); 10701 10702 PetscCall((B->ops->transcoloringapplysptoden)(coloring, B, Btdense)); 10703 PetscFunctionReturn(PETSC_SUCCESS); 10704 } 10705 10706 /*@C 10707 MatTransColoringApplyDenToSp - Given a symbolic matrix product Csp=A*B^T for which 10708 a `MatTransposeColoring` context has been created and a dense matrix Cden=A*Btdense 10709 in which Btdens is obtained from `MatTransColoringApplySpToDen()`, recover sparse matrix 10710 Csp from Cden. 10711 10712 Collective 10713 10714 Input Parameters: 10715 + coloring - coloring context created with `MatTransposeColoringCreate()` 10716 - Cden - matrix product of a sparse matrix and a dense matrix Btdense 10717 10718 Output Parameter: 10719 . Csp - sparse matrix 10720 10721 Level: developer 10722 10723 Note: 10724 These are used internally for some implementations of `MatRARt()` 10725 10726 .seealso: `Mat`, `MatTransposeColoringCreate()`, `MatTransposeColoringDestroy()`, `MatTransColoringApplySpToDen()` 10727 @*/ 10728 PetscErrorCode MatTransColoringApplyDenToSp(MatTransposeColoring matcoloring, Mat Cden, Mat Csp) 10729 { 10730 PetscFunctionBegin; 10731 PetscValidHeaderSpecific(matcoloring, MAT_TRANSPOSECOLORING_CLASSID, 1); 10732 PetscValidHeaderSpecific(Cden, MAT_CLASSID, 2); 10733 PetscValidHeaderSpecific(Csp, MAT_CLASSID, 3); 10734 10735 PetscCall((Csp->ops->transcoloringapplydentosp)(matcoloring, Cden, Csp)); 10736 PetscCall(MatAssemblyBegin(Csp, MAT_FINAL_ASSEMBLY)); 10737 PetscCall(MatAssemblyEnd(Csp, MAT_FINAL_ASSEMBLY)); 10738 PetscFunctionReturn(PETSC_SUCCESS); 10739 } 10740 10741 /*@C 10742 MatTransposeColoringCreate - Creates a matrix coloring context for the matrix product C=A*B^T. 10743 10744 Collective 10745 10746 Input Parameters: 10747 + mat - the matrix product C 10748 - iscoloring - the coloring of the matrix; usually obtained with `MatColoringCreate()` or `DMCreateColoring()` 10749 10750 Output Parameter: 10751 . color - the new coloring context 10752 10753 Level: intermediate 10754 10755 .seealso: `Mat`, `MatTransposeColoringDestroy()`, `MatTransColoringApplySpToDen()`, 10756 `MatTransColoringApplyDenToSp()` 10757 @*/ 10758 PetscErrorCode MatTransposeColoringCreate(Mat mat, ISColoring iscoloring, MatTransposeColoring *color) 10759 { 10760 MatTransposeColoring c; 10761 MPI_Comm comm; 10762 10763 PetscFunctionBegin; 10764 PetscCall(PetscLogEventBegin(MAT_TransposeColoringCreate, mat, 0, 0, 0)); 10765 PetscCall(PetscObjectGetComm((PetscObject)mat, &comm)); 10766 PetscCall(PetscHeaderCreate(c, MAT_TRANSPOSECOLORING_CLASSID, "MatTransposeColoring", "Matrix product C=A*B^T via coloring", "Mat", comm, MatTransposeColoringDestroy, NULL)); 10767 10768 c->ctype = iscoloring->ctype; 10769 PetscUseTypeMethod(mat, transposecoloringcreate, iscoloring, c); 10770 10771 *color = c; 10772 PetscCall(PetscLogEventEnd(MAT_TransposeColoringCreate, mat, 0, 0, 0)); 10773 PetscFunctionReturn(PETSC_SUCCESS); 10774 } 10775 10776 /*@ 10777 MatGetNonzeroState - Returns a 64 bit integer representing the current state of nonzeros in the matrix. If the 10778 matrix has had no new nonzero locations added to (or removed from) the matrix since the previous call then the value will be the 10779 same, otherwise it will be larger 10780 10781 Not Collective 10782 10783 Input Parameter: 10784 . A - the matrix 10785 10786 Output Parameter: 10787 . state - the current state 10788 10789 Level: intermediate 10790 10791 Notes: 10792 You can only compare states from two different calls to the SAME matrix, you cannot compare calls between 10793 different matrices 10794 10795 Use `PetscObjectStateGet()` to check for changes to the numerical values in a matrix 10796 10797 Use the result of `PetscObjectGetId()` to compare if a previously checked matrix is the same as the current matrix, do not compare object pointers. 10798 10799 .seealso: `Mat`, `PetscObjectStateGet()`, `PetscObjectGetId()` 10800 @*/ 10801 PetscErrorCode MatGetNonzeroState(Mat mat, PetscObjectState *state) 10802 { 10803 PetscFunctionBegin; 10804 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 10805 *state = mat->nonzerostate; 10806 PetscFunctionReturn(PETSC_SUCCESS); 10807 } 10808 10809 /*@ 10810 MatCreateMPIMatConcatenateSeqMat - Creates a single large PETSc matrix by concatenating sequential 10811 matrices from each processor 10812 10813 Collective 10814 10815 Input Parameters: 10816 + comm - the communicators the parallel matrix will live on 10817 . seqmat - the input sequential matrices 10818 . n - number of local columns (or `PETSC_DECIDE`) 10819 - reuse - either `MAT_INITIAL_MATRIX` or `MAT_REUSE_MATRIX` 10820 10821 Output Parameter: 10822 . mpimat - the parallel matrix generated 10823 10824 Level: developer 10825 10826 Note: 10827 The number of columns of the matrix in EACH processor MUST be the same. 10828 10829 .seealso: `Mat` 10830 @*/ 10831 PetscErrorCode MatCreateMPIMatConcatenateSeqMat(MPI_Comm comm, Mat seqmat, PetscInt n, MatReuse reuse, Mat *mpimat) 10832 { 10833 PetscMPIInt size; 10834 10835 PetscFunctionBegin; 10836 PetscCallMPI(MPI_Comm_size(comm, &size)); 10837 if (size == 1) { 10838 if (reuse == MAT_INITIAL_MATRIX) { 10839 PetscCall(MatDuplicate(seqmat, MAT_COPY_VALUES, mpimat)); 10840 } else { 10841 PetscCall(MatCopy(seqmat, *mpimat, SAME_NONZERO_PATTERN)); 10842 } 10843 PetscFunctionReturn(PETSC_SUCCESS); 10844 } 10845 10846 PetscCheck(reuse != MAT_REUSE_MATRIX || seqmat != *mpimat, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "MAT_REUSE_MATRIX means reuse the matrix passed in as the final argument, not the original matrix"); 10847 10848 PetscCall(PetscLogEventBegin(MAT_Merge, seqmat, 0, 0, 0)); 10849 PetscCall((*seqmat->ops->creatempimatconcatenateseqmat)(comm, seqmat, n, reuse, mpimat)); 10850 PetscCall(PetscLogEventEnd(MAT_Merge, seqmat, 0, 0, 0)); 10851 PetscFunctionReturn(PETSC_SUCCESS); 10852 } 10853 10854 /*@ 10855 MatSubdomainsCreateCoalesce - Creates index subdomains by coalescing adjacent ranks' ownership ranges. 10856 10857 Collective 10858 10859 Input Parameters: 10860 + A - the matrix to create subdomains from 10861 - N - requested number of subdomains 10862 10863 Output Parameters: 10864 + n - number of subdomains resulting on this rank 10865 - iss - `IS` list with indices of subdomains on this rank 10866 10867 Level: advanced 10868 10869 Note: 10870 The number of subdomains must be smaller than the communicator size 10871 10872 .seealso: `Mat`, `IS` 10873 @*/ 10874 PetscErrorCode MatSubdomainsCreateCoalesce(Mat A, PetscInt N, PetscInt *n, IS *iss[]) 10875 { 10876 MPI_Comm comm, subcomm; 10877 PetscMPIInt size, rank, color; 10878 PetscInt rstart, rend, k; 10879 10880 PetscFunctionBegin; 10881 PetscCall(PetscObjectGetComm((PetscObject)A, &comm)); 10882 PetscCallMPI(MPI_Comm_size(comm, &size)); 10883 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 10884 PetscCheck(N >= 1 && N < (PetscInt)size, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "number of subdomains must be > 0 and < %d, got N = %" PetscInt_FMT, size, N); 10885 *n = 1; 10886 k = ((PetscInt)size) / N + ((PetscInt)size % N > 0); /* There are up to k ranks to a color */ 10887 color = rank / k; 10888 PetscCallMPI(MPI_Comm_split(comm, color, rank, &subcomm)); 10889 PetscCall(PetscMalloc1(1, iss)); 10890 PetscCall(MatGetOwnershipRange(A, &rstart, &rend)); 10891 PetscCall(ISCreateStride(subcomm, rend - rstart, rstart, 1, iss[0])); 10892 PetscCallMPI(MPI_Comm_free(&subcomm)); 10893 PetscFunctionReturn(PETSC_SUCCESS); 10894 } 10895 10896 /*@ 10897 MatGalerkin - Constructs the coarse grid problem matrix via Galerkin projection. 10898 10899 If the interpolation and restriction operators are the same, uses `MatPtAP()`. 10900 If they are not the same, uses `MatMatMatMult()`. 10901 10902 Once the coarse grid problem is constructed, correct for interpolation operators 10903 that are not of full rank, which can legitimately happen in the case of non-nested 10904 geometric multigrid. 10905 10906 Input Parameters: 10907 + restrct - restriction operator 10908 . dA - fine grid matrix 10909 . interpolate - interpolation operator 10910 . reuse - either `MAT_INITIAL_MATRIX` or `MAT_REUSE_MATRIX` 10911 - fill - expected fill, use `PETSC_DEFAULT` if you do not have a good estimate 10912 10913 Output Parameters: 10914 . A - the Galerkin coarse matrix 10915 10916 Options Database Key: 10917 . -pc_mg_galerkin <both,pmat,mat,none> - for what matrices the Galerkin process should be used 10918 10919 Level: developer 10920 10921 .seealso: `Mat`, `MatPtAP()`, `MatMatMatMult()` 10922 @*/ 10923 PetscErrorCode MatGalerkin(Mat restrct, Mat dA, Mat interpolate, MatReuse reuse, PetscReal fill, Mat *A) 10924 { 10925 IS zerorows; 10926 Vec diag; 10927 10928 PetscFunctionBegin; 10929 PetscCheck(reuse != MAT_INPLACE_MATRIX, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Inplace product not supported"); 10930 /* Construct the coarse grid matrix */ 10931 if (interpolate == restrct) { 10932 PetscCall(MatPtAP(dA, interpolate, reuse, fill, A)); 10933 } else { 10934 PetscCall(MatMatMatMult(restrct, dA, interpolate, reuse, fill, A)); 10935 } 10936 10937 /* If the interpolation matrix is not of full rank, A will have zero rows. 10938 This can legitimately happen in the case of non-nested geometric multigrid. 10939 In that event, we set the rows of the matrix to the rows of the identity, 10940 ignoring the equations (as the RHS will also be zero). */ 10941 10942 PetscCall(MatFindZeroRows(*A, &zerorows)); 10943 10944 if (zerorows != NULL) { /* if there are any zero rows */ 10945 PetscCall(MatCreateVecs(*A, &diag, NULL)); 10946 PetscCall(MatGetDiagonal(*A, diag)); 10947 PetscCall(VecISSet(diag, zerorows, 1.0)); 10948 PetscCall(MatDiagonalSet(*A, diag, INSERT_VALUES)); 10949 PetscCall(VecDestroy(&diag)); 10950 PetscCall(ISDestroy(&zerorows)); 10951 } 10952 PetscFunctionReturn(PETSC_SUCCESS); 10953 } 10954 10955 /*@C 10956 MatSetOperation - Allows user to set a matrix operation for any matrix type 10957 10958 Logically Collective 10959 10960 Input Parameters: 10961 + mat - the matrix 10962 . op - the name of the operation 10963 - f - the function that provides the operation 10964 10965 Level: developer 10966 10967 Usage: 10968 .vb 10969 extern PetscErrorCode usermult(Mat, Vec, Vec); 10970 10971 PetscCall(MatCreateXXX(comm, ..., &A)); 10972 PetscCall(MatSetOperation(A, MATOP_MULT, (PetscVoidFunction)usermult)); 10973 .ve 10974 10975 Notes: 10976 See the file `include/petscmat.h` for a complete list of matrix 10977 operations, which all have the form MATOP_<OPERATION>, where 10978 <OPERATION> is the name (in all capital letters) of the 10979 user interface routine (e.g., `MatMult()` -> `MATOP_MULT`). 10980 10981 All user-provided functions (except for `MATOP_DESTROY`) should have the same calling 10982 sequence as the usual matrix interface routines, since they 10983 are intended to be accessed via the usual matrix interface 10984 routines, e.g., 10985 .vb 10986 MatMult(Mat, Vec, Vec) -> usermult(Mat, Vec, Vec) 10987 .ve 10988 10989 In particular each function MUST return `PETSC_SUCCESS` on success and 10990 nonzero on failure. 10991 10992 This routine is distinct from `MatShellSetOperation()` in that it can be called on any matrix type. 10993 10994 .seealso: `Mat`, `MatGetOperation()`, `MatCreateShell()`, `MatShellSetContext()`, `MatShellSetOperation()` 10995 @*/ 10996 PetscErrorCode MatSetOperation(Mat mat, MatOperation op, void (*f)(void)) 10997 { 10998 PetscFunctionBegin; 10999 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 11000 if (op == MATOP_VIEW && !mat->ops->viewnative && f != (void (*)(void))(mat->ops->view)) mat->ops->viewnative = mat->ops->view; 11001 (((void (**)(void))mat->ops)[op]) = f; 11002 PetscFunctionReturn(PETSC_SUCCESS); 11003 } 11004 11005 /*@C 11006 MatGetOperation - Gets a matrix operation for any matrix type. 11007 11008 Not Collective 11009 11010 Input Parameters: 11011 + mat - the matrix 11012 - op - the name of the operation 11013 11014 Output Parameter: 11015 . f - the function that provides the operation 11016 11017 Level: developer 11018 11019 Usage: 11020 $ PetscErrorCode (*usermult)(Mat,Vec,Vec); 11021 $ MatGetOperation(A,MATOP_MULT,(void(**)(void))&usermult); 11022 11023 Notes: 11024 See the file include/petscmat.h for a complete list of matrix 11025 operations, which all have the form MATOP_<OPERATION>, where 11026 <OPERATION> is the name (in all capital letters) of the 11027 user interface routine (e.g., `MatMult()` -> `MATOP_MULT`). 11028 11029 This routine is distinct from `MatShellGetOperation()` in that it can be called on any matrix type. 11030 11031 .seealso: `Mat`, `MatSetOperation()`, `MatCreateShell()`, `MatShellGetContext()`, `MatShellGetOperation()` 11032 @*/ 11033 PetscErrorCode MatGetOperation(Mat mat, MatOperation op, void (**f)(void)) 11034 { 11035 PetscFunctionBegin; 11036 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 11037 *f = (((void (**)(void))mat->ops)[op]); 11038 PetscFunctionReturn(PETSC_SUCCESS); 11039 } 11040 11041 /*@ 11042 MatHasOperation - Determines whether the given matrix supports the particular operation. 11043 11044 Not Collective 11045 11046 Input Parameters: 11047 + mat - the matrix 11048 - op - the operation, for example, `MATOP_GET_DIAGONAL` 11049 11050 Output Parameter: 11051 . has - either `PETSC_TRUE` or `PETSC_FALSE` 11052 11053 Level: advanced 11054 11055 Note: 11056 See `MatSetOperation()` for additional discussion on naming convention and usage of `op`. 11057 11058 .seealso: `Mat`, `MatCreateShell()`, `MatGetOperation()`, `MatSetOperation()` 11059 @*/ 11060 PetscErrorCode MatHasOperation(Mat mat, MatOperation op, PetscBool *has) 11061 { 11062 PetscFunctionBegin; 11063 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 11064 PetscValidBoolPointer(has, 3); 11065 if (mat->ops->hasoperation) { 11066 PetscUseTypeMethod(mat, hasoperation, op, has); 11067 } else { 11068 if (((void **)mat->ops)[op]) *has = PETSC_TRUE; 11069 else { 11070 *has = PETSC_FALSE; 11071 if (op == MATOP_CREATE_SUBMATRIX) { 11072 PetscMPIInt size; 11073 11074 PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)mat), &size)); 11075 if (size == 1) PetscCall(MatHasOperation(mat, MATOP_CREATE_SUBMATRICES, has)); 11076 } 11077 } 11078 } 11079 PetscFunctionReturn(PETSC_SUCCESS); 11080 } 11081 11082 /*@ 11083 MatHasCongruentLayouts - Determines whether the rows and columns layouts of the matrix are congruent 11084 11085 Collective 11086 11087 Input Parameters: 11088 . mat - the matrix 11089 11090 Output Parameter: 11091 . cong - either `PETSC_TRUE` or `PETSC_FALSE` 11092 11093 Level: beginner 11094 11095 .seealso: `Mat`, `MatCreate()`, `MatSetSizes()`, `PetscLayout` 11096 @*/ 11097 PetscErrorCode MatHasCongruentLayouts(Mat mat, PetscBool *cong) 11098 { 11099 PetscFunctionBegin; 11100 PetscValidHeaderSpecific(mat, MAT_CLASSID, 1); 11101 PetscValidType(mat, 1); 11102 PetscValidBoolPointer(cong, 2); 11103 if (!mat->rmap || !mat->cmap) { 11104 *cong = mat->rmap == mat->cmap ? PETSC_TRUE : PETSC_FALSE; 11105 PetscFunctionReturn(PETSC_SUCCESS); 11106 } 11107 if (mat->congruentlayouts == PETSC_DECIDE) { /* first time we compare rows and cols layouts */ 11108 PetscCall(PetscLayoutSetUp(mat->rmap)); 11109 PetscCall(PetscLayoutSetUp(mat->cmap)); 11110 PetscCall(PetscLayoutCompare(mat->rmap, mat->cmap, cong)); 11111 if (*cong) mat->congruentlayouts = 1; 11112 else mat->congruentlayouts = 0; 11113 } else *cong = mat->congruentlayouts ? PETSC_TRUE : PETSC_FALSE; 11114 PetscFunctionReturn(PETSC_SUCCESS); 11115 } 11116 11117 PetscErrorCode MatSetInf(Mat A) 11118 { 11119 PetscFunctionBegin; 11120 PetscUseTypeMethod(A, setinf); 11121 PetscFunctionReturn(PETSC_SUCCESS); 11122 } 11123 11124 /*@C 11125 MatCreateGraph - create a scalar matrix (that is a matrix with one vertex for each block vertex in the original matrix), for use in graph algorithms 11126 and possibly removes small values from the graph structure. 11127 11128 Collective 11129 11130 Input Parameters: 11131 + A - the matrix 11132 . sym - `PETSC_TRUE` indicates that the graph should be symmetrized 11133 . scale - `PETSC_TRUE` indicates that the graph edge weights should be symmetrically scaled with the diagonal entry 11134 - filter - filter value - < 0: does nothing; == 0: removes only 0.0 entries; otherwise: removes entries with abs(entries) <= value 11135 11136 Output Parameter: 11137 . graph - the resulting graph 11138 11139 Level: advanced 11140 11141 .seealso: `Mat`, `MatCreate()`, `PCGAMG` 11142 @*/ 11143 PetscErrorCode MatCreateGraph(Mat A, PetscBool sym, PetscBool scale, PetscReal filter, Mat *graph) 11144 { 11145 PetscFunctionBegin; 11146 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 11147 PetscValidType(A, 1); 11148 PetscValidLogicalCollectiveBool(A, scale, 3); 11149 PetscValidPointer(graph, 5); 11150 PetscUseTypeMethod(A, creategraph, sym, scale, filter, graph); 11151 PetscFunctionReturn(PETSC_SUCCESS); 11152 } 11153 11154 /*@ 11155 MatEliminateZeros - eliminate the nondiagonal zero entries in place from the nonzero structure of a sparse `Mat` in place, 11156 meaning the same memory is used for the matrix, and no new memory is allocated. 11157 11158 Collective on mat 11159 11160 Input Parameter: 11161 . A - the matrix 11162 11163 Output Parameter: 11164 . A - the matrix 11165 11166 Level: intermediate 11167 11168 .seealso: `Mat`, `MatCreate()`, `MatCreateGraph()`, `MatChop()` 11169 @*/ 11170 PetscErrorCode MatEliminateZeros(Mat A) 11171 { 11172 PetscFunctionBegin; 11173 PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 11174 PetscUseTypeMethod(A, eliminatezeros); 11175 PetscFunctionReturn(PETSC_SUCCESS); 11176 } 11177