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