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