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