1 2 #ifndef lint 3 static char vcid[] = "$Id: matrix.c,v 1.93 1995/10/07 20:47:53 curfman Exp curfman $"; 4 #endif 5 6 /* 7 This is where the abstract matrix operations are defined 8 */ 9 10 #include "petsc.h" 11 #include "matimpl.h" /*I "mat.h" I*/ 12 #include "vec/vecimpl.h" 13 #include "pinclude/pviewer.h" 14 15 /*@C 16 MatGetReordering - Gets a reordering for a matrix to reduce fill or to 17 improve numerical stability of LU factorization. 18 19 Input Parameters: 20 . mat - the matrix 21 . type - type of reordering, one of the following: 22 $ ORDER_NATURAL - Natural 23 $ ORDER_ND - Nested Dissection 24 $ ORDER_1WD - One-way Dissection 25 $ ORDER_RCM - Reverse Cuthill-McGee 26 $ ORDER_QMD - Quotient Minimum Degree 27 28 Output Parameters: 29 . rperm - row permutation indices 30 . cperm - column permutation indices 31 32 Options Database Keys: 33 To specify the ordering through the options database, use one of 34 the following 35 $ -mat_order natural, -mat_order nd, -mat_order 1wd, 36 $ -mat_order rcm, -mat_order qmd 37 38 Notes: 39 If the column permutations and row permutations are the same, 40 then MatGetReordering() returns 0 in cperm. 41 42 The user can define additional orderings; see MatReorderingRegister(). 43 44 .keywords: matrix, set, ordering, factorization, direct, ILU, LU, 45 fill, reordering, natural, Nested Dissection, 46 One-way Dissection, Cholesky, Reverse Cuthill-McGee, 47 Quotient Minimum Degree 48 49 .seealso: MatGetReorderingTypeFromOptions(), MatReorderingRegister() 50 @*/ 51 int MatGetReordering(Mat mat,MatOrdering type,IS *rperm,IS *cperm) 52 { 53 int ierr; 54 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 55 if (!mat->ops.getreordering) {*rperm = 0; *cperm = 0; return 0;} 56 PLogEventBegin(MAT_GetReordering,mat,0,0,0); 57 ierr = MatGetReorderingTypeFromOptions(0,&type); CHKERRQ(ierr); 58 ierr = (*mat->ops.getreordering)(mat,type,rperm,cperm); CHKERRQ(ierr); 59 PLogEventEnd(MAT_GetReordering,mat,0,0,0); 60 return 0; 61 } 62 63 /*@C 64 MatGetRow - Gets a row of a matrix. You MUST call MatRestoreRow() 65 for each row that you get to ensure that your application does 66 not bleed memory. 67 68 Input Parameters: 69 . mat - the matrix 70 . row - the row to get 71 72 Output Parameters: 73 . ncols - the number of nonzeros in the row 74 . cols - if nonzero, the column numbers 75 . vals - if nonzero, the values 76 77 Notes: 78 This routine is provided for people who need to have direct access 79 to the structure of a matrix. We hope that we provide enough 80 high-level matrix routines that few users will need it. 81 82 For better efficiency, set cols and/or vals to zero if you do not 83 wish to extract these quantities. 84 85 .keywords: matrix, row, get, extract 86 87 .seealso: MatRestoreRow() 88 @*/ 89 int MatGetRow(Mat mat,int row,int *ncols,int **cols,Scalar **vals) 90 { 91 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 92 return (*mat->ops.getrow)(mat,row,ncols,cols,vals); 93 } 94 95 /*@C 96 MatRestoreRow - Frees any temporary space allocated by MatGetRow(). 97 98 Input Parameters: 99 . mat - the matrix 100 . row - the row to get 101 . ncols, cols - the number of nonzeros and their columns 102 . vals - if nonzero the column values 103 104 .keywords: matrix, row, restore 105 106 .seealso: MatGetRow() 107 @*/ 108 int MatRestoreRow(Mat mat,int row,int *ncols,int **cols,Scalar **vals) 109 { 110 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 111 if (!mat->ops.restorerow) return 0; 112 return (*mat->ops.restorerow)(mat,row,ncols,cols,vals); 113 } 114 /*@ 115 MatView - Visualizes a matrix object. 116 117 Input Parameters: 118 . mat - the matrix 119 . ptr - visualization context 120 121 Notes: 122 The available visualization contexts include 123 $ STDOUT_VIEWER_SELF - standard output (default) 124 $ STDOUT_VIEWER_WORLD - synchronized standard 125 $ output where only the first processor opens 126 $ the file. All other processors send their 127 $ data to the first processor to print. 128 129 The user can open alternative vistualization contexts with 130 $ ViewerFileOpenASCII() - output to a specified file 131 $ DrawOpenX() - output nonzero matrix structure to 132 $ an X window display 133 $ ViewerMatlabOpen() - output matrix to Matlab viewer. 134 $ Currently only the sequential dense and AIJ 135 $ matrix types support the Matlab viewer. 136 137 The user can call ViewerFileSetFormat() to specify the output 138 format. Available formats include 139 $ FILE_FORMAT_DEFAULT - default, prints matrix contents 140 $ FILE_FORMAT_IMPL - implementation-specific format 141 $ (which is in many cases the same as the default) 142 $ FILE_FORMAT_INFO - basic information about the matrix 143 $ size and structure (not the matrix entries) 144 145 .keywords: matrix, view, visualize 146 147 .seealso: ViewerFileSetFormat(), ViewerFileOpenASCII(), DrawOpenX(), 148 ViewerMatlabOpen() 149 @*/ 150 int MatView(Mat mat,Viewer ptr) 151 { 152 int format, ierr, rows, cols,nz, nzalloc, mem; 153 FILE *fd; 154 char *cstring; 155 PetscObject vobj = (PetscObject) ptr; 156 157 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 158 if (!ptr) { /* so that viewers may be used from debuggers */ 159 ptr = STDOUT_VIEWER_SELF; vobj = (PetscObject) ptr; 160 } 161 ierr = ViewerFileGetFormat_Private(ptr,&format); CHKERRQ(ierr); 162 ierr = ViewerFileGetPointer_Private(ptr,&fd); CHKERRQ(ierr); 163 if (vobj->cookie == VIEWER_COOKIE && format == FILE_FORMAT_INFO && 164 (vobj->type == ASCII_FILE_VIEWER || vobj->type == ASCII_FILES_VIEWER)) { 165 MPIU_fprintf(mat->comm,fd,"Matrix Object:\n"); 166 ierr = MatGetName(mat,&cstring); CHKERRQ(ierr); 167 ierr = MatGetSize(mat,&rows,&cols); CHKERRQ(ierr); 168 MPIU_fprintf(mat->comm,fd, 169 " type=%s, rows=%d, cols=%d\n",cstring,rows,cols); 170 if (mat->ops.getinfo) { 171 ierr = MatGetInfo(mat,MAT_LOCAL,&nz,&nzalloc,&mem); CHKERRQ(ierr); 172 MPIU_fprintf(mat->comm,fd, 173 " nonzeros=%d, allocated nonzeros=%d\n",nz,nzalloc); 174 } 175 } 176 if (mat->view) {ierr = (*mat->view)((PetscObject)mat,ptr); CHKERRQ(ierr);} 177 return 0; 178 } 179 /*@C 180 MatDestroy - Frees space taken by a matrix. 181 182 Input Parameter: 183 . mat - the matrix 184 185 .keywords: matrix, destroy 186 @*/ 187 int MatDestroy(Mat mat) 188 { 189 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 190 return (*mat->destroy)((PetscObject)mat); 191 } 192 /*@ 193 MatValidMatrix - Returns 1 if a valid matrix else 0. 194 195 Input Parameter: 196 . m - the matrix to check 197 198 .keywords: matrix, valid 199 @*/ 200 int MatValidMatrix(Mat m) 201 { 202 if (!m) return 0; 203 if (m->cookie != MAT_COOKIE) return 0; 204 return 1; 205 } 206 207 /*@ 208 MatSetValues - Inserts or adds a block of values into a matrix. 209 These values may be cached, so MatAssemblyBegin() and MatAssemblyEnd() 210 MUST be called after all calls to MatSetValues() have been completed. 211 212 Input Parameters: 213 . mat - the matrix 214 . v - a logically two-dimensional array of values 215 . m, indexm - the number of rows and their global indices 216 . n, indexn - the number of columns and their global indices 217 . addv - either ADD_VALUES or INSERT_VALUES, where 218 $ ADD_VALUES - adds values to any existing entries 219 $ INSERT_VALUES - replaces existing entries with new values 220 221 Notes: 222 By default the values, v, are row-oriented and unsorted. 223 See MatSetOptions() for other options. 224 225 Calls to MatSetValues() with the INSERT_VALUES and ADD_VALUES 226 options cannot be mixed without intervening calls to the assembly 227 routines. 228 229 .keywords: matrix, insert, add, set, values 230 231 .seealso: MatSetOptions(), MatAssemblyBegin(), MatAssemblyEnd() 232 @*/ 233 int MatSetValues(Mat mat,int m,int *idxm,int n,int *idxn,Scalar *v, 234 InsertMode addv) 235 { 236 int ierr; 237 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 238 PLogEventBegin(MAT_SetValues,mat,0,0,0); 239 ierr = (*mat->ops.setvalues)(mat,m,idxm,n,idxn,v,addv);CHKERRQ(ierr); 240 PLogEventEnd(MAT_SetValues,mat,0,0,0); 241 return 0; 242 } 243 244 /* --------------------------------------------------------*/ 245 /*@ 246 MatMult - Computes matrix-vector product. 247 248 Input Parameters: 249 . mat - the matrix 250 . x - the vector to be multilplied 251 252 Output Parameters: 253 . y - the result 254 255 .keywords: matrix, multiply, matrix-vector product 256 257 .seealso: MatMultTrans(), MatMultAdd(), MatMultTransAdd() 258 @*/ 259 int MatMult(Mat mat,Vec x,Vec y) 260 { 261 int ierr; 262 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 263 PETSCVALIDHEADERSPECIFIC(x,VEC_COOKIE);PETSCVALIDHEADERSPECIFIC(y,VEC_COOKIE); 264 PLogEventBegin(MAT_Mult,mat,x,y,0); 265 ierr = (*mat->ops.mult)(mat,x,y); CHKERRQ(ierr); 266 PLogEventEnd(MAT_Mult,mat,x,y,0); 267 return 0; 268 } 269 /*@ 270 MatMultTrans - Computes matrix transpose times a vector. 271 272 Input Parameters: 273 . mat - the matrix 274 . x - the vector to be multilplied 275 276 Output Parameters: 277 . y - the result 278 279 .keywords: matrix, multiply, matrix-vector product, transpose 280 281 .seealso: MatMult(), MatMultAdd(), MatMultTransAdd() 282 @*/ 283 int MatMultTrans(Mat mat,Vec x,Vec y) 284 { 285 int ierr; 286 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 287 PETSCVALIDHEADERSPECIFIC(x,VEC_COOKIE); PETSCVALIDHEADERSPECIFIC(y,VEC_COOKIE); 288 PLogEventBegin(MAT_MultTrans,mat,x,y,0); 289 ierr = (*mat->ops.multtrans)(mat,x,y); CHKERRQ(ierr); 290 PLogEventEnd(MAT_MultTrans,mat,x,y,0); 291 return 0; 292 } 293 /*@ 294 MatMultAdd - Computes v3 = v2 + A * v1. 295 296 Input Parameters: 297 . mat - the matrix 298 . v1, v2 - the vectors 299 300 Output Parameters: 301 . v3 - the result 302 303 .keywords: matrix, multiply, matrix-vector product, add 304 305 .seealso: MatMultTrans(), MatMult(), MatMultTransAdd() 306 @*/ 307 int MatMultAdd(Mat mat,Vec v1,Vec v2,Vec v3) 308 { 309 int ierr; 310 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);PETSCVALIDHEADERSPECIFIC(v1,VEC_COOKIE); 311 PETSCVALIDHEADERSPECIFIC(v2,VEC_COOKIE); PETSCVALIDHEADERSPECIFIC(v3,VEC_COOKIE); 312 PLogEventBegin(MAT_MultAdd,mat,v1,v2,v3); 313 ierr = (*mat->ops.multadd)(mat,v1,v2,v3); CHKERRQ(ierr); 314 PLogEventEnd(MAT_MultAdd,mat,v1,v2,v3); 315 return 0; 316 } 317 /*@ 318 MatMultTransAdd - Computes v3 = v2 + A' * v1. 319 320 Input Parameters: 321 . mat - the matrix 322 . v1, v2 - the vectors 323 324 Output Parameters: 325 . v3 - the result 326 327 .keywords: matrix, multiply, matrix-vector product, transpose, add 328 329 .seealso: MatMultTrans(), MatMultAdd(), MatMult() 330 @*/ 331 int MatMultTransAdd(Mat mat,Vec v1,Vec v2,Vec v3) 332 { 333 int ierr; 334 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); PETSCVALIDHEADERSPECIFIC(v1,VEC_COOKIE); 335 PETSCVALIDHEADERSPECIFIC(v2,VEC_COOKIE); PETSCVALIDHEADERSPECIFIC(v3,VEC_COOKIE); 336 if (!mat->ops.multtransadd) SETERRQ(PETSC_ERR_SUP,"MatMultTransAdd"); 337 PLogEventBegin(MAT_MultTransAdd,mat,v1,v2,v3); 338 ierr = (*mat->ops.multtransadd)(mat,v1,v2,v3); CHKERRQ(ierr); 339 PLogEventEnd(MAT_MultTransAdd,mat,v1,v2,v3); 340 return 0; 341 } 342 /* ------------------------------------------------------------*/ 343 /*@ 344 MatGetInfo - Returns information about matrix storage (number of 345 nonzeros, memory). 346 347 Input Parameters: 348 . mat - the matrix 349 350 Output Parameters: 351 . flag - flag indicating the type of parameters to be returned 352 $ flag = MAT_LOCAL: local matrix 353 $ flag = MAT_GLOBAL_MAX: maximum over all processors 354 $ flag = MAT_GLOBAL_SUM: sum over all processors 355 . nz - the number of nonzeros 356 . nzalloc - the number of allocated nonzeros 357 . mem - the memory used (in bytes) 358 359 .keywords: matrix, get, info, storage, nonzeros, memory 360 @*/ 361 int MatGetInfo(Mat mat,MatInfoType flag,int *nz,int *nzalloc,int *mem) 362 { 363 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 364 if (!mat->ops.getinfo) SETERRQ(PETSC_ERR_SUP,"MatGetInfo"); 365 return (*mat->ops.getinfo)(mat,flag,nz,nzalloc,mem); 366 } 367 /* ----------------------------------------------------------*/ 368 /*@ 369 MatLUFactor - Performs in-place LU factorization of matrix. 370 371 Input Parameters: 372 . mat - the matrix 373 . row - row permutation 374 . col - column permutation 375 . f - expected fill as ratio of original fill. 376 377 .keywords: matrix, factor, LU, in-place 378 379 .seealso: MatLUFactorSymbolic(), MatLUFactorNumeric(), MatCholeskyFactor() 380 @*/ 381 int MatLUFactor(Mat mat,IS row,IS col,double f) 382 { 383 int ierr; 384 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 385 if (!mat->ops.lufactor) SETERRQ(PETSC_ERR_SUP,"MatLUFactor"); 386 PLogEventBegin(MAT_LUFactor,mat,row,col,0); 387 ierr = (*mat->ops.lufactor)(mat,row,col,f); CHKERRQ(ierr); 388 PLogEventEnd(MAT_LUFactor,mat,row,col,0); 389 return 0; 390 } 391 /*@ 392 MatILUFactor - Performs in-place ILU factorization of matrix. 393 394 Input Parameters: 395 . mat - the matrix 396 . row - row permutation 397 . col - column permutation 398 . f - expected fill as ratio of original fill. 399 . level - number of levels of fill. 400 401 Note: probably really only in-place when level is zero. 402 .keywords: matrix, factor, ILU, in-place 403 404 .seealso: MatILUFactorSymbolic(), MatLUFactorNumeric(), MatCholeskyFactor() 405 @*/ 406 int MatILUFactor(Mat mat,IS row,IS col,double f,int level) 407 { 408 int ierr; 409 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 410 if (!mat->ops.ilufactor) SETERRQ(PETSC_ERR_SUP,"MatILUFactor"); 411 PLogEventBegin(MAT_ILUFactor,mat,row,col,0); 412 ierr = (*mat->ops.ilufactor)(mat,row,col,f,level); CHKERRQ(ierr); 413 PLogEventEnd(MAT_ILUFactor,mat,row,col,0); 414 return 0; 415 } 416 417 /*@ 418 MatLUFactorSymbolic - Performs symbolic LU factorization of matrix. 419 Call this routine before calling MatLUFactorNumeric(). 420 421 Input Parameters: 422 . mat - the matrix 423 . row, col - row and column permutations 424 . f - expected fill as ratio of the original number of nonzeros, 425 for example 3.0; choosing this parameter well can result in 426 more efficient use of time and space. 427 428 Output Parameters: 429 . fact - new matrix that has been symbolically factored 430 431 .keywords: matrix, factor, LU, symbol CHKERRQ(ierr);ic 432 433 .seealso: MatLUFactor(), MatLUFactorNumeric(), MatCholeskyFactor() 434 @*/ 435 int MatLUFactorSymbolic(Mat mat,IS row,IS col,double f,Mat *fact) 436 { 437 int ierr; 438 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 439 if (!fact) SETERRQ(1,"MatLUFactorSymbolic:Missing factor matrix argument"); 440 if (!mat->ops.lufactorsymbolic) SETERRQ(PETSC_ERR_SUP,"MatLUFactorSymbolic"); 441 OptionsGetDouble(0,"-mat_lu_fill",&f); 442 PLogEventBegin(MAT_LUFactorSymbolic,mat,row,col,0); 443 ierr = (*mat->ops.lufactorsymbolic)(mat,row,col,f,fact); CHKERRQ(ierr); 444 PLogEventEnd(MAT_LUFactorSymbolic,mat,row,col,0); 445 return 0; 446 } 447 /*@ 448 MatLUFactorNumeric - Performs numeric LU factorization of a matrix. 449 Call this routine after first calling MatLUFactorSymbolic(). 450 451 Input Parameters: 452 . mat - the matrix 453 . row, col - row and column permutations 454 455 Output Parameters: 456 . fact - symbolically factored matrix that must have been generated 457 by MatLUFactorSymbolic() 458 459 Notes: 460 See MatLUFactor() for in-place factorization. See 461 MatCholeskyFactorNumeric() for the symmetric, positive definite case. 462 463 .keywords: matrix, factor, LU, numeric 464 465 .seealso: MatLUFactorSymbolic(), MatLUFactor(), MatCholeskyFactor() 466 @*/ 467 int MatLUFactorNumeric(Mat mat,Mat *fact) 468 { 469 int ierr; 470 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 471 if (!fact) SETERRQ(1,"MatLUFactorNumeric:Missing factor matrix argument"); 472 if (!mat->ops.lufactornumeric) SETERRQ(PETSC_ERR_SUP,"MatLUFactorNumeric"); 473 PLogEventBegin(MAT_LUFactorNumeric,mat,*fact,0,0); 474 ierr = (*mat->ops.lufactornumeric)(mat,fact); CHKERRQ(ierr); 475 PLogEventEnd(MAT_LUFactorNumeric,mat,*fact,0,0); 476 return 0; 477 } 478 /*@ 479 MatCholeskyFactor - Performs in-place Cholesky factorization of a 480 symmetric matrix. 481 482 Input Parameters: 483 . mat - the matrix 484 . perm - row and column permutations 485 . f - expected fill as ratio of original fill 486 487 Notes: 488 See MatLUFactor() for the nonsymmetric case. See also 489 MatCholeskyFactorSymbolic(), and MatCholeskyFactorNumeric(). 490 491 .keywords: matrix, factor, in-place, Cholesky 492 493 .seealso: MatLUFactor(), MatCholeskyFactorSymbolic() 494 .seealso: MatCholeskyFactorNumeric() 495 @*/ 496 int MatCholeskyFactor(Mat mat,IS perm,double f) 497 { 498 int ierr; 499 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 500 if (!mat->ops.choleskyfactor) SETERRQ(PETSC_ERR_SUP,"MatCholeskyFactor"); 501 PLogEventBegin(MAT_CholeskyFactor,mat,perm,0,0); 502 ierr = (*mat->ops.choleskyfactor)(mat,perm,f); CHKERRQ(ierr); 503 PLogEventEnd(MAT_CholeskyFactor,mat,perm,0,0); 504 return 0; 505 } 506 /*@ 507 MatCholeskyFactorSymbolic - Performs symbolic Cholesky factorization 508 of a symmetric matrix. 509 510 Input Parameters: 511 . mat - the matrix 512 . perm - row and column permutations 513 . f - expected fill as ratio of original 514 515 Output Parameter: 516 . fact - the factored matrix 517 518 Notes: 519 See MatLUFactorSymbolic() for the nonsymmetric case. See also 520 MatCholeskyFactor() and MatCholeskyFactorNumeric(). 521 522 .keywords: matrix, factor, factorization, symbolic, Cholesky 523 524 .seealso: MatLUFactorSymbolic() 525 .seealso: MatCholeskyFactor(), MatCholeskyFactorNumeric() 526 @*/ 527 int MatCholeskyFactorSymbolic(Mat mat,IS perm,double f,Mat *fact) 528 { 529 int ierr; 530 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 531 if (!fact) 532 SETERRQ(1,"MatCholeskyFactorSymbolic:Missing factor matrix argument"); 533 if (!mat->ops.choleskyfactorsymbolic)SETERRQ(PETSC_ERR_SUP,"MatCholeskyFactorSymbolic"); 534 PLogEventBegin(MAT_CholeskyFactorSymbolic,mat,perm,0,0); 535 ierr = (*mat->ops.choleskyfactorsymbolic)(mat,perm,f,fact); CHKERRQ(ierr); 536 PLogEventEnd(MAT_CholeskyFactorSymbolic,mat,perm,0,0); 537 return 0; 538 } 539 /*@ 540 MatCholeskyFactorNumeric - Performs numeric Cholesky factorization 541 of a symmetric matrix. Call this routine after first calling 542 MatCholeskyFactorSymbolic(). 543 544 Input Parameter: 545 . mat - the initial matrix 546 547 Output Parameter: 548 . fact - the factored matrix 549 550 .keywords: matrix, factor, numeric, Cholesky 551 552 .seealso: MatCholeskyFactorSymbolic(), MatCholeskyFactor() 553 .seealso: MatLUFactorNumeric() 554 @*/ 555 int MatCholeskyFactorNumeric(Mat mat,Mat *fact) 556 { 557 int ierr; 558 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 559 if (!fact) 560 SETERRQ(1,"MatCholeskyFactorNumeric:Missing factor matrix argument"); 561 if (!mat->ops.choleskyfactornumeric) SETERRQ(PETSC_ERR_SUP,"MatCholeskyFactorNumeric"); 562 PLogEventBegin(MAT_CholeskyFactorNumeric,mat,*fact,0,0); 563 ierr = (*mat->ops.choleskyfactornumeric)(mat,fact); CHKERRQ(ierr); 564 PLogEventEnd(MAT_CholeskyFactorNumeric,mat,*fact,0,0); 565 return 0; 566 } 567 /* ----------------------------------------------------------------*/ 568 /*@ 569 MatSolve - Solves A x = b, given a factored matrix. 570 571 Input Parameters: 572 . mat - the factored matrix 573 . b - the right-hand-side vector 574 575 Output Parameter: 576 . x - the result vector 577 578 .keywords: matrix, linear system, solve, LU, Cholesky, triangular solve 579 580 .seealso: MatSolveAdd(), MatSolveTrans(), MatSolveTransAdd() 581 @*/ 582 int MatSolve(Mat mat,Vec b,Vec x) 583 { 584 int ierr; 585 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 586 PETSCVALIDHEADERSPECIFIC(b,VEC_COOKIE); PETSCVALIDHEADERSPECIFIC(x,VEC_COOKIE); 587 if (!mat->factor) SETERRQ(1,"MatSolve:Unfactored matrix"); 588 if (!mat->ops.solve) SETERRQ(PETSC_ERR_SUP,"MatSolve"); 589 PLogEventBegin(MAT_Solve,mat,b,x,0); 590 ierr = (*mat->ops.solve)(mat,b,x); CHKERRQ(ierr); 591 PLogEventEnd(MAT_Solve,mat,b,x,0); 592 return 0; 593 } 594 595 /* @ 596 MatForwardSolve - Solves L x = b, given a factored matrix, A = LU. 597 598 Input Parameters: 599 . mat - the factored matrix 600 . b - the right-hand-side vector 601 602 Output Parameter: 603 . x - the result vector 604 605 Notes: 606 MatSolve() should be used for most applications, as it performs 607 a forward solve followed by a backward solve. 608 609 .keywords: matrix, forward, LU, Cholesky, triangular solve 610 611 .seealso: MatSolve(), MatBackwardSolve() 612 @ */ 613 int MatForwardSolve(Mat mat,Vec b,Vec x) 614 { 615 int ierr; 616 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 617 PETSCVALIDHEADERSPECIFIC(b,VEC_COOKIE); PETSCVALIDHEADERSPECIFIC(x,VEC_COOKIE); 618 if (!mat->factor) SETERRQ(1,"MatForwardSolve:Unfactored matrix"); 619 if (!mat->ops.forwardsolve) SETERRQ(PETSC_ERR_SUP,"MatForwardSolve"); 620 PLogEventBegin(MAT_ForwardSolve,mat,b,x,0); 621 ierr = (*mat->ops.forwardsolve)(mat,b,x); CHKERRQ(ierr); 622 PLogEventEnd(MAT_ForwardSolve,mat,b,x,0); 623 return 0; 624 } 625 626 /* @ 627 MatBackwardSolve - Solves U x = b, given a factored matrix, A = LU. 628 629 Input Parameters: 630 . mat - the factored matrix 631 . b - the right-hand-side vector 632 633 Output Parameter: 634 . x - the result vector 635 636 Notes: 637 MatSolve() should be used for most applications, as it performs 638 a forward solve followed by a backward solve. 639 640 .keywords: matrix, backward, LU, Cholesky, triangular solve 641 642 .seealso: MatSolve(), MatForwardSolve() 643 @ */ 644 int MatBackwardSolve(Mat mat,Vec b,Vec x) 645 { 646 int ierr; 647 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 648 PETSCVALIDHEADERSPECIFIC(b,VEC_COOKIE); PETSCVALIDHEADERSPECIFIC(x,VEC_COOKIE); 649 if (!mat->factor) SETERRQ(1,"MatBackwardSolve:Unfactored matrix"); 650 if (!mat->ops.backwardsolve) SETERRQ(PETSC_ERR_SUP,"MatBackwardSolve"); 651 PLogEventBegin(MAT_BackwardSolve,mat,b,x,0); 652 ierr = (*mat->ops.backwardsolve)(mat,b,x); CHKERRQ(ierr); 653 PLogEventEnd(MAT_BackwardSolve,mat,b,x,0); 654 return 0; 655 } 656 657 /*@ 658 MatSolveAdd - Computes x = y + inv(A)*b, given a factored matrix. 659 660 Input Parameters: 661 . mat - the factored matrix 662 . b - the right-hand-side vector 663 . y - the vector to be added to 664 665 Output Parameter: 666 . x - the result vector 667 668 .keywords: matrix, linear system, solve, LU, Cholesky, add 669 670 .seealso: MatSolve(), MatSolveTrans(), MatSolveTransAdd() 671 @*/ 672 int MatSolveAdd(Mat mat,Vec b,Vec y,Vec x) 673 { 674 Scalar one = 1.0; 675 Vec tmp; 676 int ierr; 677 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);PETSCVALIDHEADERSPECIFIC(y,VEC_COOKIE); 678 PETSCVALIDHEADERSPECIFIC(b,VEC_COOKIE); PETSCVALIDHEADERSPECIFIC(x,VEC_COOKIE); 679 if (!mat->factor) SETERRQ(1,"MatSolveAdd:Unfactored matrix"); 680 PLogEventBegin(MAT_SolveAdd,mat,b,x,y); 681 if (mat->ops.solveadd) { 682 ierr = (*mat->ops.solveadd)(mat,b,y,x); CHKERRQ(ierr); 683 } 684 else { 685 /* do the solve then the add manually */ 686 if (x != y) { 687 ierr = MatSolve(mat,b,x); CHKERRQ(ierr); 688 ierr = VecAXPY(&one,y,x); CHKERRQ(ierr); 689 } 690 else { 691 ierr = VecDuplicate(x,&tmp); CHKERRQ(ierr); 692 PLogObjectParent(mat,tmp); 693 ierr = VecCopy(x,tmp); CHKERRQ(ierr); 694 ierr = MatSolve(mat,b,x); CHKERRQ(ierr); 695 ierr = VecAXPY(&one,tmp,x); CHKERRQ(ierr); 696 ierr = VecDestroy(tmp); CHKERRQ(ierr); 697 } 698 } 699 PLogEventEnd(MAT_SolveAdd,mat,b,x,y); 700 return 0; 701 } 702 /*@ 703 MatSolveTrans - Solves A' x = b, given a factored matrix. 704 705 Input Parameters: 706 . mat - the factored matrix 707 . b - the right-hand-side vector 708 709 Output Parameter: 710 . x - the result vector 711 712 .keywords: matrix, linear system, solve, LU, Cholesky, transpose 713 714 .seealso: MatSolve(), MatSolveAdd(), MatSolveTransAdd() 715 @*/ 716 int MatSolveTrans(Mat mat,Vec b,Vec x) 717 { 718 int ierr; 719 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 720 PETSCVALIDHEADERSPECIFIC(b,VEC_COOKIE); PETSCVALIDHEADERSPECIFIC(x,VEC_COOKIE); 721 if (!mat->factor) SETERRQ(1,"MatSolveTrans:Unfactored matrix"); 722 if (!mat->ops.solvetrans) SETERRQ(PETSC_ERR_SUP,"MatSolveTrans"); 723 PLogEventBegin(MAT_SolveTrans,mat,b,x,0); 724 ierr = (*mat->ops.solvetrans)(mat,b,x); CHKERRQ(ierr); 725 PLogEventEnd(MAT_SolveTrans,mat,b,x,0); 726 return 0; 727 } 728 /*@ 729 MatSolveTransAdd - Computes x = y + inv(trans(A)) b, given a 730 factored matrix. 731 732 Input Parameters: 733 . mat - the factored matrix 734 . b - the right-hand-side vector 735 . y - the vector to be added to 736 737 Output Parameter: 738 . x - the result vector 739 740 .keywords: matrix, linear system, solve, LU, Cholesky, transpose, add 741 742 .seealso: MatSolve(), MatSolveAdd(), MatSolveTrans() 743 @*/ 744 int MatSolveTransAdd(Mat mat,Vec b,Vec y,Vec x) 745 { 746 Scalar one = 1.0; 747 int ierr; 748 Vec tmp; 749 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);PETSCVALIDHEADERSPECIFIC(y,VEC_COOKIE); 750 PETSCVALIDHEADERSPECIFIC(b,VEC_COOKIE); PETSCVALIDHEADERSPECIFIC(x,VEC_COOKIE); 751 if (!mat->factor) SETERRQ(1,"MatSolveTransAdd:Unfactored matrix"); 752 PLogEventBegin(MAT_SolveTransAdd,mat,b,x,y); 753 if (mat->ops.solvetransadd) { 754 ierr = (*mat->ops.solvetransadd)(mat,b,y,x); CHKERRQ(ierr); 755 } 756 else { 757 /* do the solve then the add manually */ 758 if (x != y) { 759 ierr = MatSolveTrans(mat,b,x); CHKERRQ(ierr); 760 ierr = VecAXPY(&one,y,x); CHKERRQ(ierr); 761 } 762 else { 763 ierr = VecDuplicate(x,&tmp); CHKERRQ(ierr); 764 PLogObjectParent(mat,tmp); 765 ierr = VecCopy(x,tmp); CHKERRQ(ierr); 766 ierr = MatSolveTrans(mat,b,x); CHKERRQ(ierr); 767 ierr = VecAXPY(&one,tmp,x); CHKERRQ(ierr); 768 ierr = VecDestroy(tmp); CHKERRQ(ierr); 769 } 770 } 771 PLogEventEnd(MAT_SolveTransAdd,mat,b,x,y); 772 return 0; 773 } 774 /* ----------------------------------------------------------------*/ 775 776 /*@ 777 MatRelax - Computes one relaxation sweep. 778 779 Input Parameters: 780 . mat - the matrix 781 . b - the right hand side 782 . omega - the relaxation factor 783 . flag - flag indicating the type of SOR, one of 784 $ SOR_FORWARD_SWEEP 785 $ SOR_BACKWARD_SWEEP 786 $ SOR_SYMMETRIC_SWEEP (SSOR method) 787 $ SOR_LOCAL_FORWARD_SWEEP 788 $ SOR_LOCAL_BACKWARD_SWEEP 789 $ SOR_LOCAL_SYMMETRIC_SWEEP (local SSOR) 790 $ SOR_APPLY_UPPER, SOR_APPLY_LOWER - applies 791 $ upper/lower triangular part of matrix to 792 $ vector (with omega) 793 $ SOR_ZERO_INITIAL_GUESS - zero initial guess 794 . shift - diagonal shift 795 . its - the number of iterations 796 797 Output Parameters: 798 . x - the solution (can contain an initial guess) 799 800 Notes: 801 SOR_LOCAL_FORWARD_SWEEP, SOR_LOCAL_BACKWARD_SWEEP, and 802 SOR_LOCAL_SYMMETRIC_SWEEP perform seperate independent smoothings 803 on each processor. 804 805 Application programmers will not generally use MatRelax() directly, 806 but instead will employ the SLES/PC interface. 807 808 Notes for Advanced Users: 809 The flags are implemented as bitwise inclusive or operations. 810 For example, use (SOR_ZERO_INITIAL_GUESS | SOR_SYMMETRIC_SWEEP) 811 to specify a zero initial guess for SSOR. 812 813 .keywords: matrix, relax, relaxation, sweep 814 @*/ 815 int MatRelax(Mat mat,Vec b,double omega,MatSORType flag,double shift, 816 int its,Vec x) 817 { 818 int ierr; 819 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 820 PETSCVALIDHEADERSPECIFIC(b,VEC_COOKIE); PETSCVALIDHEADERSPECIFIC(x,VEC_COOKIE); 821 if (!mat->ops.relax) SETERRQ(PETSC_ERR_SUP,"MatRelax"); 822 PLogEventBegin(MAT_Relax,mat,b,x,0); 823 ierr =(*mat->ops.relax)(mat,b,omega,flag,shift,its,x); CHKERRQ(ierr); 824 PLogEventEnd(MAT_Relax,mat,b,x,0); 825 return 0; 826 } 827 828 /*@C 829 MatConvert - Converts a matrix to another matrix, either of the same 830 or different type. 831 832 Input Parameters: 833 . mat - the matrix 834 . newtype - new matrix type. Use MATSAME to create a new matrix of the 835 same type as the original matrix. 836 837 Output Parameter: 838 . M - pointer to place new matrix 839 840 .keywords: matrix, copy, convert 841 @*/ 842 int MatConvert(Mat mat,MatType newtype,Mat *M) 843 { 844 int ierr; 845 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 846 if (!M) SETERRQ(1,"MatConvert:Bad address"); 847 if (newtype == mat->type || newtype == MATSAME) { 848 if (mat->ops.copyprivate) { 849 PLogEventBegin(MAT_Convert,mat,0,0,0); 850 ierr = (*mat->ops.copyprivate)(mat,M); CHKERRQ(ierr); 851 PLogEventEnd(MAT_Convert,mat,0,0,0); 852 return 0; 853 } 854 } 855 if (!mat->ops.convert) SETERRQ(PETSC_ERR_SUP,"MatConvert"); 856 PLogEventBegin(MAT_Convert,mat,0,0,0); 857 if (mat->ops.convert) { 858 ierr = (*mat->ops.convert)(mat,newtype,M); CHKERRQ(ierr); 859 } 860 PLogEventEnd(MAT_Convert,mat,0,0,0); 861 return 0; 862 } 863 864 /*@ 865 MatGetDiagonal - Gets the diagonal of a matrix. 866 867 Input Parameters: 868 . mat - the matrix 869 870 Output Parameters: 871 . v - the vector for storing the diagonal 872 873 .keywords: matrix, get, diagonal 874 @*/ 875 int MatGetDiagonal(Mat mat,Vec v) 876 { 877 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 878 PETSCVALIDHEADERSPECIFIC(v,VEC_COOKIE); 879 if (mat->ops.getdiagonal) return (*mat->ops.getdiagonal)(mat,v); 880 SETERRQ(PETSC_ERR_SUP,"MatGetDiagonal"); 881 } 882 883 /*@ 884 MatTranspose - Computes an in-place or out-of-place transpose of a matrix. 885 886 Input Parameters: 887 . mat - the matrix to transpose 888 889 Output Parameters: 890 . B - the transpose - pass in zero for an in-place transpose 891 892 .keywords: matrix, transpose 893 @*/ 894 int MatTranspose(Mat mat,Mat *B) 895 { 896 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 897 if (mat->ops.transpose) return (*mat->ops.transpose)(mat,B); 898 SETERRQ(PETSC_ERR_SUP,"MatTranspose"); 899 } 900 901 /*@ 902 MatEqual - Compares two matrices. Returns 1 if two matrices are equal. 903 904 Input Parameters: 905 . mat1 - the first matrix 906 . mat2 - the second matrix 907 908 Returns: 909 Returns 1 if the matrices are equal; returns 0 otherwise. 910 911 .keywords: matrix, equal, equivalent 912 @*/ 913 int MatEqual(Mat mat1,Mat mat2) 914 { 915 PETSCVALIDHEADERSPECIFIC(mat1,MAT_COOKIE); PETSCVALIDHEADERSPECIFIC(mat2,MAT_COOKIE); 916 if (mat1->ops.equal) return (*mat1->ops.equal)(mat1,mat2); 917 SETERRQ(PETSC_ERR_SUP,"MatEqual"); 918 } 919 920 /*@ 921 MatScale - Scales a matrix on the left and right by diagonal 922 matrices that are stored as vectors. Either of the two scaling 923 matrices can be null. 924 925 Input Parameters: 926 . mat - the matrix to be scaled 927 . l - the left scaling vector 928 . r - the right scaling vector 929 930 .keywords: matrix, scale 931 @*/ 932 int MatScale(Mat mat,Vec l,Vec r) 933 { 934 int ierr; 935 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 936 if (!mat->ops.scale) SETERRQ(PETSC_ERR_SUP,"MatScale"); 937 if (l) PETSCVALIDHEADERSPECIFIC(l,VEC_COOKIE); 938 if (r) PETSCVALIDHEADERSPECIFIC(r,VEC_COOKIE); 939 PLogEventBegin(MAT_Scale,mat,0,0,0); 940 ierr = (*mat->ops.scale)(mat,l,r); CHKERRQ(ierr); 941 PLogEventEnd(MAT_Scale,mat,0,0,0); 942 return 0; 943 } 944 945 /*@ 946 MatNorm - Calculates various norms of a matrix. 947 948 Input Parameters: 949 . mat - the matrix 950 . type - the type of norm, NORM_1, NORM_2, NORM_FROBENIUS, NORM_INFINITY 951 952 Output Parameters: 953 . norm - the resulting norm 954 955 .keywords: matrix, norm, Frobenius 956 @*/ 957 int MatNorm(Mat mat,MatNormType type,double *norm) 958 { 959 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 960 if (!norm) SETERRQ(1,"MatNorm:bad addess for value"); 961 if (mat->ops.norm) return (*mat->ops.norm)(mat,type,norm); 962 SETERRQ(PETSC_ERR_SUP,"MatNorm:Not for this matrix type"); 963 } 964 965 /*@ 966 MatAssemblyBegin - Begins assembling the matrix. This routine should 967 be called after completing all calls to MatSetValues(). 968 969 Input Parameters: 970 . mat - the matrix 971 . type - type of assembly, either FLUSH_ASSEMBLY or FINAL_ASSEMBLY 972 973 Notes: 974 MatSetValues() generally caches the values. The matrix is ready to 975 use only after MatAssemblyBegin() and MatAssemblyEnd() have been called. 976 Use FLUSH_ASSEMBLY when switching between ADD_VALUES and SetValues; use 977 FINAL_ASSEMBLY for the final assembly before the matrix is used. 978 979 .keywords: matrix, assembly, assemble, begin 980 981 .seealso: MatAssemblyEnd(), MatSetValues() 982 @*/ 983 int MatAssemblyBegin(Mat mat,MatAssemblyType type) 984 { 985 int ierr; 986 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 987 PLogEventBegin(MAT_AssemblyBegin,mat,0,0,0); 988 if (mat->ops.assemblybegin) {ierr = (*mat->ops.assemblybegin)(mat,type); CHKERRQ(ierr);} 989 PLogEventEnd(MAT_AssemblyBegin,mat,0,0,0); 990 return 0; 991 } 992 /*@ 993 MatAssemblyEnd - Completes assembling the matrix. This routine should 994 be called after all calls to MatSetValues() and after MatAssemblyBegin(). 995 996 Input Parameters: 997 . mat - the matrix 998 . type - type of assembly, either FLUSH_ASSEMBLY or FINAL_ASSEMBLY 999 1000 Note: 1001 MatSetValues() generally caches the values. The matrix is ready to 1002 use only after MatAssemblyBegin() and MatAssemblyEnd() have been called. 1003 Use FLUSH_ASSEMBLY when switching between ADD_VALUES and SetValues; use 1004 FINAL_ASSEMBLY for the final assembly before the matrix is used. 1005 1006 .keywords: matrix, assembly, assemble, end 1007 1008 .seealso: MatAssemblyBegin(), MatSetValues() 1009 @*/ 1010 int MatAssemblyEnd(Mat mat,MatAssemblyType type) 1011 { 1012 int ierr; 1013 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 1014 PLogEventBegin(MAT_AssemblyEnd,mat,0,0,0); 1015 if (mat->ops.assemblyend) {ierr = (*mat->ops.assemblyend)(mat,type); CHKERRQ(ierr);} 1016 PLogEventEnd(MAT_AssemblyEnd,mat,0,0,0); 1017 return 0; 1018 } 1019 /*@ 1020 MatCompress - Tries to store the matrix in as little space as 1021 possible. May fail if memory is already fully used, since it 1022 tries to allocate new space. 1023 1024 Input Parameters: 1025 . mat - the matrix 1026 1027 .keywords: matrix, compress 1028 @*/ 1029 int MatCompress(Mat mat) 1030 { 1031 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 1032 if (mat->ops.compress) return (*mat->ops.compress)(mat); 1033 return 0; 1034 } 1035 /*@ 1036 MatSetOption - Sets a parameter option for a matrix. Some options 1037 may be specific to certain storage formats. Some options 1038 determine how values will be inserted (or added). Sorted, 1039 row-oriented input will generally assemble the fastest. The default 1040 is row-oriented, nonsorted input. 1041 1042 Input Parameters: 1043 . mat - the matrix 1044 . option - the option, one of the following: 1045 $ ROW_ORIENTED, COLUMN_ORIENTED, ROWS_SORTED, 1046 $ COLUMNS_SORTED, NO_NEW_NONZERO_LOCATIONS, 1047 $ SYMMETRIC_MATRIX, 1048 $ YES_NEW_NONZERO_LOCATIONS, and possibly others 1049 1050 Notes: 1051 If using a Fortran 77 module to compute a matrix, one may need to 1052 use the column-oriented option (or convert to the row-oriented 1053 format). 1054 1055 NO_NEW_NONZERO_LOCATIONS indicates that any add or insertion 1056 that will generate a new entry in the nonzero structure is ignored. 1057 What this means is if memory is not allocated for this particular 1058 lot, then the insertion is ignored. For dense matrices where 1059 the entire array is allocated no entries are ever ignored. This 1060 may not be a good idea?? 1061 1062 .keywords: matrix, option, row-oriented, column-oriented, sorted, nonzero 1063 @*/ 1064 int MatSetOption(Mat mat,MatOption op) 1065 { 1066 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 1067 if (mat->ops.setoption) return (*mat->ops.setoption)(mat,op); 1068 return 0; 1069 } 1070 1071 /*@ 1072 MatZeroEntries - Zeros all entries of a matrix. For sparse matrices 1073 this routine retains the old nonzero structure. 1074 1075 Input Parameters: 1076 . mat - the matrix 1077 1078 .keywords: matrix, zero, entries 1079 1080 .seealso: MatZeroRows() 1081 @*/ 1082 int MatZeroEntries(Mat mat) 1083 { 1084 int ierr; 1085 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 1086 if (!mat->ops.zeroentries) SETERRQ(PETSC_ERR_SUP,"MatZeroEntries"); 1087 PLogEventBegin(MAT_ZeroEntries,mat,0,0,0); 1088 ierr = (*mat->ops.zeroentries)(mat); CHKERRQ(ierr); 1089 PLogEventEnd(MAT_ZeroEntries,mat,0,0,0); 1090 return 0; 1091 } 1092 1093 /*@ 1094 MatZeroRows - Zeros all entries (except possibly the main diagonal) 1095 of a set of rows of a matrix. 1096 1097 Input Parameters: 1098 . mat - the matrix 1099 . is - index set of rows to remove 1100 . diag - pointer to value put in all diagonals of eliminated rows. 1101 Note that diag is not a pointer to an array, but merely a 1102 pointer to a single value. 1103 1104 Notes: 1105 For the AIJ and row matrix formats this removes the old nonzero 1106 structure, but does not release memory. For the dense and block 1107 diagonal formats this does not alter the nonzero structure. 1108 1109 The user can set a value in the diagonal entry (or for the AIJ and 1110 row formats can optionally remove the main diagonal entry from the 1111 nonzero structure as well, by passing a null pointer as the final 1112 argument). 1113 1114 .keywords: matrix, zero, rows, boundary conditions 1115 1116 .seealso: MatZeroEntries(), MatGetSubMatrix(), MatGetSubMatrixInPlace() 1117 @*/ 1118 int MatZeroRows(Mat mat,IS is, Scalar *diag) 1119 { 1120 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 1121 if (mat->ops.zerorows) return (*mat->ops.zerorows)(mat,is,diag); 1122 SETERRQ(PETSC_ERR_SUP,"MatZeroRows"); 1123 } 1124 1125 /*@ 1126 MatGetSize - Returns the numbers of rows and columns in a matrix. 1127 1128 Input Parameter: 1129 . mat - the matrix 1130 1131 Output Parameters: 1132 . m - the number of global rows 1133 . n - the number of global columns 1134 1135 .keywords: matrix, dimension, size, rows, columns, global, get 1136 1137 .seealso: MatGetLocalSize() 1138 @*/ 1139 int MatGetSize(Mat mat,int *m,int* n) 1140 { 1141 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 1142 if (!m || !n) SETERRQ(1,"MatGetSize:Bad address for result"); 1143 return (*mat->ops.getsize)(mat,m,n); 1144 } 1145 1146 /*@ 1147 MatGetLocalSize - Returns the number of rows and columns in a matrix 1148 stored locally. This information may be implementation dependent, so 1149 use with care. 1150 1151 Input Parameters: 1152 . mat - the matrix 1153 1154 Output Parameters: 1155 . m - the number of local rows 1156 . n - the number of local columns 1157 1158 .keywords: matrix, dimension, size, local, rows, columns, get 1159 1160 .seealso: MatGetSize() 1161 @*/ 1162 int MatGetLocalSize(Mat mat,int *m,int* n) 1163 { 1164 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 1165 if (!m || !n) SETERRQ(1,"MatGetLocalSize:Bad address for result"); 1166 return (*mat->ops.getlocalsize)(mat,m,n); 1167 } 1168 1169 /*@ 1170 MatGetOwnershipRange - Returns the range of matrix rows owned by 1171 this processor, assuming that the matrix is laid out with the first 1172 n1 rows on the first processor, the next n2 rows on the second, etc. 1173 For certain parallel layouts this range may not be well-defined. 1174 1175 Input Parameters: 1176 . mat - the matrix 1177 1178 Output Parameters: 1179 . m - the first local row 1180 . n - one more then the last local row 1181 1182 .keywords: matrix, get, range, ownership 1183 @*/ 1184 int MatGetOwnershipRange(Mat mat,int *m,int* n) 1185 { 1186 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 1187 if (!m || !n) SETERRQ(1,"MatGetOwnershipRange:Bad address for result"); 1188 if (mat->ops.getownershiprange) return (*mat->ops.getownershiprange)(mat,m,n); 1189 SETERRQ(PETSC_ERR_SUP,"MatGetOwnershipRange"); 1190 } 1191 1192 /*@ 1193 MatILUFactorSymbolic - Performs symbolic ILU factorization of a matrix. 1194 Uses levels of fill only, not drop tolerance. Use MatLUFactorNumeric() 1195 to complete the factorization. 1196 1197 Input Parameters: 1198 . mat - the matrix 1199 . row - row permutation 1200 . column - column permutation 1201 . fill - number of levels of fill 1202 . f - expected fill as ratio of original fill 1203 1204 Output Parameters: 1205 . fact - puts factor 1206 1207 .keywords: matrix, factor, incomplete, ILU, symbolic, fill 1208 1209 .seealso: MatLUFactorSymbolic(), MatLUFactorNumeric() 1210 @*/ 1211 int MatILUFactorSymbolic(Mat mat,IS row,IS col,double f,int fill,Mat *fact) 1212 { 1213 int ierr; 1214 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 1215 if (fill < 0) SETERRQ(1,"MatILUFactorSymbolic:Levels of fill negative"); 1216 if (!fact) SETERRQ(1,"MatILUFactorSymbolic:Fact argument is missing"); 1217 if (!mat->ops.ilufactorsymbolic) SETERRQ(PETSC_ERR_SUP,"MatILUFactorSymbolic"); 1218 PLogEventBegin(MAT_ILUFactorSymbolic,mat,row,col,0); 1219 ierr = (*mat->ops.ilufactorsymbolic)(mat,row,col,f,fill,fact); 1220 CHKERRQ(ierr); 1221 PLogEventEnd(MAT_ILUFactorSymbolic,mat,row,col,0); 1222 return 0; 1223 } 1224 1225 /*@ 1226 MatIncompleteCholeskyFactorSymbolic - Performs symbolic incomplete 1227 Cholesky factorization for a symmetric matrix. Use 1228 MatCholeskyFactorNumeric() to complete the factorization. 1229 1230 Input Parameters: 1231 . mat - the matrix 1232 . perm - row and column permutation 1233 . fill - levels of fill 1234 . f - expected fill as ratio of original fill 1235 1236 Output Parameter: 1237 . fact - the factored matrix 1238 1239 Note: Currently only no-fill factorization is supported. 1240 1241 .keywords: matrix, factor, incomplete, ICC, Cholesky, symbolic, fill 1242 1243 .seealso: MatCholeskyFactorNumeric(), MatCholeskyFactor() 1244 @*/ 1245 int MatIncompleteCholeskyFactorSymbolic(Mat mat,IS perm,double f,int fill, 1246 Mat *fact) 1247 { 1248 int ierr; 1249 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 1250 if (fill < 0) SETERRQ(1,"MatIncompleteCholeskyFactorSymbolic:Fill negative"); 1251 if (!fact) SETERRQ(1,"MatIncompleteCholeskyFactorSymbolic:Missing fact argument"); 1252 if (!mat->ops.incompletecholeskyfactorsymbolic) 1253 SETERRQ(PETSC_ERR_SUP,"MatIncompleteCholeskyFactorSymbolic"); 1254 PLogEventBegin(MAT_IncompleteCholeskyFactorSymbolic,mat,perm,0,0); 1255 ierr = (*mat->ops.incompletecholeskyfactorsymbolic)(mat,perm,f,fill,fact); 1256 CHKERRQ(ierr); 1257 PLogEventEnd(MAT_IncompleteCholeskyFactorSymbolic,mat,perm,0,0); 1258 return 0; 1259 } 1260 1261 /*@C 1262 MatGetArray - Returns a pointer to the element values in the matrix. 1263 This routine is implementation dependent, and may not even work for 1264 certain matrix types. 1265 1266 Input Parameter: 1267 . mat - the matrix 1268 1269 Output Parameter: 1270 . v - the location of the values 1271 1272 .keywords: matrix, array, elements, values 1273 @*/ 1274 int MatGetArray(Mat mat,Scalar **v) 1275 { 1276 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 1277 if (!v) SETERRQ(1,"MatGetArray:Bad input, array pointer location"); 1278 if (!mat->ops.getarray) SETERRQ(PETSC_ERR_SUP,"MatGetArraye"); 1279 return (*mat->ops.getarray)(mat,v); 1280 } 1281 1282 /*@ 1283 MatGetSubMatrix - Extracts a submatrix from a matrix. 1284 1285 Input Parameters: 1286 . mat - the matrix 1287 . irow, icol - index sets of rows and columns to extract 1288 1289 Output Parameter: 1290 . submat - the submatrix 1291 1292 Notes: 1293 MatGetSubMatrix() can be useful in setting boundary conditions. 1294 1295 .keywords: matrix, get, submatrix, boundary conditions 1296 1297 .seealso: MatZeroRows(), MatGetSubMatrixInPlace() 1298 @*/ 1299 int MatGetSubMatrix(Mat mat,IS irow,IS icol,Mat *submat) 1300 { 1301 int ierr; 1302 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 1303 if (!mat->ops.getsubmatrix) SETERRQ(PETSC_ERR_SUP,"MatGetSubMatrix"); 1304 PLogEventBegin(MAT_GetSubMatrix,mat,irow,icol,0); 1305 ierr = (*mat->ops.getsubmatrix)(mat,irow,icol,submat); CHKERRQ(ierr); 1306 PLogEventEnd(MAT_GetSubMatrix,mat,irow,icol,0); 1307 return 0; 1308 } 1309 1310 /*@ 1311 MatGetSubMatrixInPlace - Extracts a submatrix from a matrix, returning 1312 the submatrix in place of the original matrix. 1313 1314 Input Parameters: 1315 . mat - the matrix 1316 . irow, icol - index sets of rows and columns to extract 1317 1318 .keywords: matrix, get, submatrix, boundary conditions, in-place 1319 1320 .seealso: MatZeroRows(), MatGetSubMatrix() 1321 @*/ 1322 int MatGetSubMatrixInPlace(Mat mat,IS irow,IS icol) 1323 { 1324 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 1325 if (!mat->ops.getsubmatrixinplace) SETERRQ(PETSC_ERR_SUP,"MatGetSubmatrixInPlace"); 1326 return (*mat->ops.getsubmatrixinplace)(mat,irow,icol); 1327 } 1328 1329 /*@ 1330 MatGetType - Returns the type of the matrix, one of MATSEQDENSE, MATSEQAIJ, etc. 1331 1332 Input Parameter: 1333 . mat - the matrix 1334 1335 Ouput Parameter: 1336 . type - the matrix type 1337 1338 Notes: 1339 The matrix types are given in petsc/include/mat.h 1340 1341 .keywords: matrix, get, type 1342 1343 .seealso: MatGetName() 1344 @*/ 1345 int MatGetType(Mat mat,MatType *type) 1346 { 1347 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 1348 *type = (MatType) mat->type; 1349 return 0; 1350 } 1351