1 #ifndef lint 2 static char vcid[] = "$Id: matrix.c,v 1.126 1996/01/12 22:06:54 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() 536 .seealso: MatCholeskyFactorNumeric() 537 @*/ 538 int MatCholeskyFactor(Mat mat,IS perm,double f) 539 { 540 int ierr; 541 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 542 if (!mat->ops.choleskyfactor) SETERRQ(PETSC_ERR_SUP,"MatCholeskyFactor"); 543 PLogEventBegin(MAT_CholeskyFactor,mat,perm,0,0); 544 ierr = (*mat->ops.choleskyfactor)(mat,perm,f); CHKERRQ(ierr); 545 PLogEventEnd(MAT_CholeskyFactor,mat,perm,0,0); 546 return 0; 547 } 548 /*@ 549 MatCholeskyFactorSymbolic - Performs symbolic Cholesky factorization 550 of a symmetric matrix. 551 552 Input Parameters: 553 . mat - the matrix 554 . perm - row and column permutations 555 . f - expected fill as ratio of original 556 557 Output Parameter: 558 . fact - the factored matrix 559 560 Notes: 561 See MatLUFactorSymbolic() for the nonsymmetric case. See also 562 MatCholeskyFactor() and MatCholeskyFactorNumeric(). 563 564 .keywords: matrix, factor, factorization, symbolic, Cholesky 565 566 .seealso: MatLUFactorSymbolic() 567 .seealso: MatCholeskyFactor(), MatCholeskyFactorNumeric() 568 @*/ 569 int MatCholeskyFactorSymbolic(Mat mat,IS perm,double f,Mat *fact) 570 { 571 int ierr; 572 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 573 if (!fact) 574 SETERRQ(1,"MatCholeskyFactorSymbolic:Missing factor matrix argument"); 575 if (!mat->ops.choleskyfactorsymbolic)SETERRQ(PETSC_ERR_SUP,"MatCholeskyFactorSymbolic"); 576 PLogEventBegin(MAT_CholeskyFactorSymbolic,mat,perm,0,0); 577 ierr = (*mat->ops.choleskyfactorsymbolic)(mat,perm,f,fact); CHKERRQ(ierr); 578 PLogEventEnd(MAT_CholeskyFactorSymbolic,mat,perm,0,0); 579 return 0; 580 } 581 /*@ 582 MatCholeskyFactorNumeric - Performs numeric Cholesky factorization 583 of a symmetric matrix. Call this routine after first calling 584 MatCholeskyFactorSymbolic(). 585 586 Input Parameter: 587 . mat - the initial matrix 588 589 Output Parameter: 590 . fact - the factored matrix 591 592 .keywords: matrix, factor, numeric, Cholesky 593 594 .seealso: MatCholeskyFactorSymbolic(), MatCholeskyFactor() 595 .seealso: MatLUFactorNumeric() 596 @*/ 597 int MatCholeskyFactorNumeric(Mat mat,Mat *fact) 598 { 599 int ierr; 600 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 601 if (!fact) 602 SETERRQ(1,"MatCholeskyFactorNumeric:Missing factor matrix argument"); 603 if (!mat->ops.choleskyfactornumeric) SETERRQ(PETSC_ERR_SUP,"MatCholeskyFactorNumeric"); 604 PLogEventBegin(MAT_CholeskyFactorNumeric,mat,*fact,0,0); 605 ierr = (*mat->ops.choleskyfactornumeric)(mat,fact); CHKERRQ(ierr); 606 PLogEventEnd(MAT_CholeskyFactorNumeric,mat,*fact,0,0); 607 return 0; 608 } 609 /* ----------------------------------------------------------------*/ 610 /*@ 611 MatSolve - Solves A x = b, given a factored matrix. 612 613 Input Parameters: 614 . mat - the factored matrix 615 . b - the right-hand-side vector 616 617 Output Parameter: 618 . x - the result vector 619 620 .keywords: matrix, linear system, solve, LU, Cholesky, triangular solve 621 622 .seealso: MatSolveAdd(), MatSolveTrans(), MatSolveTransAdd() 623 @*/ 624 int MatSolve(Mat mat,Vec b,Vec x) 625 { 626 int ierr; 627 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 628 PETSCVALIDHEADERSPECIFIC(b,VEC_COOKIE); PETSCVALIDHEADERSPECIFIC(x,VEC_COOKIE); 629 if (!mat->factor) SETERRQ(1,"MatSolve:Unfactored matrix"); 630 if (!mat->ops.solve) SETERRQ(PETSC_ERR_SUP,"MatSolve"); 631 PLogEventBegin(MAT_Solve,mat,b,x,0); 632 ierr = (*mat->ops.solve)(mat,b,x); CHKERRQ(ierr); 633 PLogEventEnd(MAT_Solve,mat,b,x,0); 634 return 0; 635 } 636 637 /* @ 638 MatForwardSolve - Solves L x = b, given a factored matrix, A = LU. 639 640 Input Parameters: 641 . mat - the factored matrix 642 . b - the right-hand-side vector 643 644 Output Parameter: 645 . x - the result vector 646 647 Notes: 648 MatSolve() should be used for most applications, as it performs 649 a forward solve followed by a backward solve. 650 651 .keywords: matrix, forward, LU, Cholesky, triangular solve 652 653 .seealso: MatSolve(), MatBackwardSolve() 654 @ */ 655 int MatForwardSolve(Mat mat,Vec b,Vec x) 656 { 657 int ierr; 658 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 659 PETSCVALIDHEADERSPECIFIC(b,VEC_COOKIE); PETSCVALIDHEADERSPECIFIC(x,VEC_COOKIE); 660 if (!mat->factor) SETERRQ(1,"MatForwardSolve:Unfactored matrix"); 661 if (!mat->ops.forwardsolve) SETERRQ(PETSC_ERR_SUP,"MatForwardSolve"); 662 PLogEventBegin(MAT_ForwardSolve,mat,b,x,0); 663 ierr = (*mat->ops.forwardsolve)(mat,b,x); CHKERRQ(ierr); 664 PLogEventEnd(MAT_ForwardSolve,mat,b,x,0); 665 return 0; 666 } 667 668 /* @ 669 MatBackwardSolve - Solves U x = b, given a factored matrix, A = LU. 670 671 Input Parameters: 672 . mat - the factored matrix 673 . b - the right-hand-side vector 674 675 Output Parameter: 676 . x - the result vector 677 678 Notes: 679 MatSolve() should be used for most applications, as it performs 680 a forward solve followed by a backward solve. 681 682 .keywords: matrix, backward, LU, Cholesky, triangular solve 683 684 .seealso: MatSolve(), MatForwardSolve() 685 @ */ 686 int MatBackwardSolve(Mat mat,Vec b,Vec x) 687 { 688 int ierr; 689 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 690 PETSCVALIDHEADERSPECIFIC(b,VEC_COOKIE); PETSCVALIDHEADERSPECIFIC(x,VEC_COOKIE); 691 if (!mat->factor) SETERRQ(1,"MatBackwardSolve:Unfactored matrix"); 692 if (!mat->ops.backwardsolve) SETERRQ(PETSC_ERR_SUP,"MatBackwardSolve"); 693 PLogEventBegin(MAT_BackwardSolve,mat,b,x,0); 694 ierr = (*mat->ops.backwardsolve)(mat,b,x); CHKERRQ(ierr); 695 PLogEventEnd(MAT_BackwardSolve,mat,b,x,0); 696 return 0; 697 } 698 699 /*@ 700 MatSolveAdd - Computes x = y + inv(A)*b, given a factored matrix. 701 702 Input Parameters: 703 . mat - the factored matrix 704 . b - the right-hand-side vector 705 . y - the vector to be added to 706 707 Output Parameter: 708 . x - the result vector 709 710 .keywords: matrix, linear system, solve, LU, Cholesky, add 711 712 .seealso: MatSolve(), MatSolveTrans(), MatSolveTransAdd() 713 @*/ 714 int MatSolveAdd(Mat mat,Vec b,Vec y,Vec x) 715 { 716 Scalar one = 1.0; 717 Vec tmp; 718 int ierr; 719 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);PETSCVALIDHEADERSPECIFIC(y,VEC_COOKIE); 720 PETSCVALIDHEADERSPECIFIC(b,VEC_COOKIE); PETSCVALIDHEADERSPECIFIC(x,VEC_COOKIE); 721 if (!mat->factor) SETERRQ(1,"MatSolveAdd:Unfactored matrix"); 722 PLogEventBegin(MAT_SolveAdd,mat,b,x,y); 723 if (mat->ops.solveadd) { 724 ierr = (*mat->ops.solveadd)(mat,b,y,x); CHKERRQ(ierr); 725 } 726 else { 727 /* do the solve then the add manually */ 728 if (x != y) { 729 ierr = MatSolve(mat,b,x); CHKERRQ(ierr); 730 ierr = VecAXPY(&one,y,x); CHKERRQ(ierr); 731 } 732 else { 733 ierr = VecDuplicate(x,&tmp); CHKERRQ(ierr); 734 PLogObjectParent(mat,tmp); 735 ierr = VecCopy(x,tmp); CHKERRQ(ierr); 736 ierr = MatSolve(mat,b,x); CHKERRQ(ierr); 737 ierr = VecAXPY(&one,tmp,x); CHKERRQ(ierr); 738 ierr = VecDestroy(tmp); CHKERRQ(ierr); 739 } 740 } 741 PLogEventEnd(MAT_SolveAdd,mat,b,x,y); 742 return 0; 743 } 744 /*@ 745 MatSolveTrans - Solves A' x = b, given a factored matrix. 746 747 Input Parameters: 748 . mat - the factored matrix 749 . b - the right-hand-side vector 750 751 Output Parameter: 752 . x - the result vector 753 754 .keywords: matrix, linear system, solve, LU, Cholesky, transpose 755 756 .seealso: MatSolve(), MatSolveAdd(), MatSolveTransAdd() 757 @*/ 758 int MatSolveTrans(Mat mat,Vec b,Vec x) 759 { 760 int ierr; 761 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 762 PETSCVALIDHEADERSPECIFIC(b,VEC_COOKIE); PETSCVALIDHEADERSPECIFIC(x,VEC_COOKIE); 763 if (!mat->factor) SETERRQ(1,"MatSolveTrans:Unfactored matrix"); 764 if (!mat->ops.solvetrans) SETERRQ(PETSC_ERR_SUP,"MatSolveTrans"); 765 PLogEventBegin(MAT_SolveTrans,mat,b,x,0); 766 ierr = (*mat->ops.solvetrans)(mat,b,x); CHKERRQ(ierr); 767 PLogEventEnd(MAT_SolveTrans,mat,b,x,0); 768 return 0; 769 } 770 /*@ 771 MatSolveTransAdd - Computes x = y + inv(trans(A)) b, given a 772 factored matrix. 773 774 Input Parameters: 775 . mat - the factored matrix 776 . b - the right-hand-side vector 777 . y - the vector to be added to 778 779 Output Parameter: 780 . x - the result vector 781 782 .keywords: matrix, linear system, solve, LU, Cholesky, transpose, add 783 784 .seealso: MatSolve(), MatSolveAdd(), MatSolveTrans() 785 @*/ 786 int MatSolveTransAdd(Mat mat,Vec b,Vec y,Vec x) 787 { 788 Scalar one = 1.0; 789 int ierr; 790 Vec tmp; 791 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE);PETSCVALIDHEADERSPECIFIC(y,VEC_COOKIE); 792 PETSCVALIDHEADERSPECIFIC(b,VEC_COOKIE); PETSCVALIDHEADERSPECIFIC(x,VEC_COOKIE); 793 if (!mat->factor) SETERRQ(1,"MatSolveTransAdd:Unfactored matrix"); 794 PLogEventBegin(MAT_SolveTransAdd,mat,b,x,y); 795 if (mat->ops.solvetransadd) { 796 ierr = (*mat->ops.solvetransadd)(mat,b,y,x); CHKERRQ(ierr); 797 } 798 else { 799 /* do the solve then the add manually */ 800 if (x != y) { 801 ierr = MatSolveTrans(mat,b,x); CHKERRQ(ierr); 802 ierr = VecAXPY(&one,y,x); CHKERRQ(ierr); 803 } 804 else { 805 ierr = VecDuplicate(x,&tmp); CHKERRQ(ierr); 806 PLogObjectParent(mat,tmp); 807 ierr = VecCopy(x,tmp); CHKERRQ(ierr); 808 ierr = MatSolveTrans(mat,b,x); CHKERRQ(ierr); 809 ierr = VecAXPY(&one,tmp,x); CHKERRQ(ierr); 810 ierr = VecDestroy(tmp); CHKERRQ(ierr); 811 } 812 } 813 PLogEventEnd(MAT_SolveTransAdd,mat,b,x,y); 814 return 0; 815 } 816 /* ----------------------------------------------------------------*/ 817 818 /*@ 819 MatRelax - Computes one relaxation sweep. 820 821 Input Parameters: 822 . mat - the matrix 823 . b - the right hand side 824 . omega - the relaxation factor 825 . flag - flag indicating the type of SOR, one of 826 $ SOR_FORWARD_SWEEP 827 $ SOR_BACKWARD_SWEEP 828 $ SOR_SYMMETRIC_SWEEP (SSOR method) 829 $ SOR_LOCAL_FORWARD_SWEEP 830 $ SOR_LOCAL_BACKWARD_SWEEP 831 $ SOR_LOCAL_SYMMETRIC_SWEEP (local SSOR) 832 $ SOR_APPLY_UPPER, SOR_APPLY_LOWER - applies 833 $ upper/lower triangular part of matrix to 834 $ vector (with omega) 835 $ SOR_ZERO_INITIAL_GUESS - zero initial guess 836 . shift - diagonal shift 837 . its - the number of iterations 838 839 Output Parameters: 840 . x - the solution (can contain an initial guess) 841 842 Notes: 843 SOR_LOCAL_FORWARD_SWEEP, SOR_LOCAL_BACKWARD_SWEEP, and 844 SOR_LOCAL_SYMMETRIC_SWEEP perform seperate independent smoothings 845 on each processor. 846 847 Application programmers will not generally use MatRelax() directly, 848 but instead will employ the SLES/PC interface. 849 850 Notes for Advanced Users: 851 The flags are implemented as bitwise inclusive or operations. 852 For example, use (SOR_ZERO_INITIAL_GUESS | SOR_SYMMETRIC_SWEEP) 853 to specify a zero initial guess for SSOR. 854 855 .keywords: matrix, relax, relaxation, sweep 856 @*/ 857 int MatRelax(Mat mat,Vec b,double omega,MatSORType flag,double shift, 858 int its,Vec x) 859 { 860 int ierr; 861 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 862 PETSCVALIDHEADERSPECIFIC(b,VEC_COOKIE); PETSCVALIDHEADERSPECIFIC(x,VEC_COOKIE); 863 if (!mat->ops.relax) SETERRQ(PETSC_ERR_SUP,"MatRelax"); 864 PLogEventBegin(MAT_Relax,mat,b,x,0); 865 ierr =(*mat->ops.relax)(mat,b,omega,flag,shift,its,x); CHKERRQ(ierr); 866 PLogEventEnd(MAT_Relax,mat,b,x,0); 867 return 0; 868 } 869 870 /* 871 Default matrix copy routine. 872 */ 873 int MatCopy_Basic(Mat A,Mat B) 874 { 875 int ierr,i,rstart,rend,nz,*cwork; 876 Scalar *vwork; 877 878 ierr = MatZeroEntries(B); CHKERRQ(ierr); 879 ierr = MatGetOwnershipRange(A,&rstart,&rend); CHKERRQ(ierr); 880 for (i=rstart; i<rend; i++) { 881 ierr = MatGetRow(A,i,&nz,&cwork,&vwork); CHKERRQ(ierr); 882 ierr = MatSetValues(B,1,&i,nz,cwork,vwork,INSERT_VALUES); CHKERRQ(ierr); 883 ierr = MatRestoreRow(A,i,&nz,&cwork,&vwork); CHKERRQ(ierr); 884 } 885 ierr = MatAssemblyBegin(B,FINAL_ASSEMBLY); CHKERRQ(ierr); 886 ierr = MatAssemblyEnd(B,FINAL_ASSEMBLY); CHKERRQ(ierr); 887 return 0; 888 } 889 890 /*@C 891 MatCopy - Copys a matrix to another matrix. 892 893 Input Parameters: 894 . A - the matrix 895 896 Output Parameter: 897 . B - where the copy is put 898 899 Notes: 900 MatCopy() copies the matrix entries of a matrix to another existing 901 matrix (after first zeroing the second matrix). A related routine is 902 MatConvert(), which first creates a new matrix and then copies the data. 903 904 .keywords: matrix, copy, convert 905 906 .seealso: MatConvert() 907 @*/ 908 int MatCopy(Mat A,Mat B) 909 { 910 int ierr; 911 PETSCVALIDHEADERSPECIFIC(A,MAT_COOKIE); 912 PETSCVALIDHEADERSPECIFIC(B,MAT_COOKIE); 913 914 PLogEventBegin(MAT_Copy,A,B,0,0); 915 if (A->ops.copy) { 916 ierr = (*A->ops.copy)(A,B); CHKERRQ(ierr); 917 } 918 else { /* generic conversion */ 919 ierr = MatCopy_Basic(A,B); CHKERRQ(ierr); 920 } 921 PLogEventEnd(MAT_Copy,A,B,0,0); 922 return 0; 923 } 924 925 /*@C 926 MatConvert - Converts a matrix to another matrix, either of the same 927 or different type. 928 929 Input Parameters: 930 . mat - the matrix 931 . newtype - new matrix type. Use MATSAME to create a new matrix of the 932 same type as the original matrix. 933 934 Output Parameter: 935 . M - pointer to place new matrix 936 937 Notes: 938 MatConvert() first creates a new matrix and then copies the data from 939 the first matrix. A related routine is MatCopy(), which copies the matrix 940 entries of one matrix to another already existing matrix context. 941 942 .keywords: matrix, copy, convert 943 944 .seealso: MatCopy() 945 @*/ 946 int MatConvert(Mat mat,MatType newtype,Mat *M) 947 { 948 int ierr; 949 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 950 if (!M) SETERRQ(1,"MatConvert:Bad new matrix address"); 951 PLogEventBegin(MAT_Convert,mat,0,0,0); 952 if (newtype == mat->type || newtype == MATSAME) { 953 if (mat->ops.convertsametype) { /* customized copy */ 954 ierr = (*mat->ops.convertsametype)(mat,M,COPY_VALUES); CHKERRQ(ierr); 955 } 956 } 957 else if (mat->ops.convert) { /* customized conversion */ 958 ierr = (*mat->ops.convert)(mat,newtype,M); CHKERRQ(ierr); 959 } 960 else { /* generic conversion */ 961 ierr = MatConvert_Basic(mat,newtype,M); CHKERRQ(ierr); 962 } 963 PLogEventEnd(MAT_Convert,mat,0,0,0); 964 return 0; 965 } 966 967 /*@ 968 MatGetDiagonal - Gets the diagonal of a matrix. 969 970 Input Parameters: 971 . mat - the matrix 972 973 Output Parameters: 974 . v - the vector for storing the diagonal 975 976 .keywords: matrix, get, diagonal 977 @*/ 978 int MatGetDiagonal(Mat mat,Vec v) 979 { 980 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 981 PETSCVALIDHEADERSPECIFIC(v,VEC_COOKIE); 982 if (mat->ops.getdiagonal) return (*mat->ops.getdiagonal)(mat,v); 983 SETERRQ(PETSC_ERR_SUP,"MatGetDiagonal"); 984 } 985 986 /*@C 987 MatTranspose - Computes an in-place or out-of-place transpose of a matrix. 988 989 Input Parameters: 990 . mat - the matrix to transpose 991 992 Output Parameters: 993 . B - the transpose - pass in zero for an in-place transpose 994 995 .keywords: matrix, transpose 996 @*/ 997 int MatTranspose(Mat mat,Mat *B) 998 { 999 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 1000 if (mat->ops.transpose) return (*mat->ops.transpose)(mat,B); 1001 SETERRQ(PETSC_ERR_SUP,"MatTranspose"); 1002 } 1003 1004 /*@ 1005 MatEqual - Compares two matrices. Returns 1 if two matrices are equal. 1006 1007 Input Parameters: 1008 . mat1 - the first matrix 1009 . mat2 - the second matrix 1010 1011 Returns: 1012 Returns 1 if the matrices are equal; returns 0 otherwise. 1013 1014 .keywords: matrix, equal, equivalent 1015 @*/ 1016 int MatEqual(Mat mat1,Mat mat2) 1017 { 1018 PETSCVALIDHEADERSPECIFIC(mat1,MAT_COOKIE); PETSCVALIDHEADERSPECIFIC(mat2,MAT_COOKIE); 1019 if (mat1->ops.equal) return (*mat1->ops.equal)(mat1,mat2); 1020 SETERRQ(PETSC_ERR_SUP,"MatEqual"); 1021 } 1022 1023 /*@ 1024 MatScale - Scales a matrix on the left and right by diagonal 1025 matrices that are stored as vectors. Either of the two scaling 1026 matrices can be null. 1027 1028 Input Parameters: 1029 . mat - the matrix to be scaled 1030 . l - the left scaling vector 1031 . r - the right scaling vector 1032 1033 .keywords: matrix, scale 1034 @*/ 1035 int MatScale(Mat mat,Vec l,Vec r) 1036 { 1037 int ierr; 1038 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 1039 if (!mat->ops.scale) SETERRQ(PETSC_ERR_SUP,"MatScale"); 1040 if (l) PETSCVALIDHEADERSPECIFIC(l,VEC_COOKIE); 1041 if (r) PETSCVALIDHEADERSPECIFIC(r,VEC_COOKIE); 1042 PLogEventBegin(MAT_Scale,mat,0,0,0); 1043 ierr = (*mat->ops.scale)(mat,l,r); CHKERRQ(ierr); 1044 PLogEventEnd(MAT_Scale,mat,0,0,0); 1045 return 0; 1046 } 1047 1048 /*@ 1049 MatNorm - Calculates various norms of a matrix. 1050 1051 Input Parameters: 1052 . mat - the matrix 1053 . type - the type of norm, NORM_1, NORM_2, NORM_FROBENIUS, NORM_INFINITY 1054 1055 Output Parameters: 1056 . norm - the resulting norm 1057 1058 .keywords: matrix, norm, Frobenius 1059 @*/ 1060 int MatNorm(Mat mat,NormType type,double *norm) 1061 { 1062 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 1063 if (!norm) SETERRQ(1,"MatNorm:bad addess for value"); 1064 if (mat->ops.norm) return (*mat->ops.norm)(mat,type,norm); 1065 SETERRQ(PETSC_ERR_SUP,"MatNorm:Not for this matrix type"); 1066 } 1067 1068 /*@ 1069 MatAssemblyBegin - Begins assembling the matrix. This routine should 1070 be called after completing all calls to MatSetValues(). 1071 1072 Input Parameters: 1073 . mat - the matrix 1074 . type - type of assembly, either FLUSH_ASSEMBLY or FINAL_ASSEMBLY 1075 1076 Notes: 1077 MatSetValues() generally caches the values. The matrix is ready to 1078 use only after MatAssemblyBegin() and MatAssemblyEnd() have been called. 1079 Use FLUSH_ASSEMBLY when switching between ADD_VALUES and SetValues; use 1080 FINAL_ASSEMBLY for the final assembly before the matrix is used. 1081 1082 .keywords: matrix, assembly, assemble, begin 1083 1084 .seealso: MatAssemblyEnd(), MatSetValues() 1085 @*/ 1086 int MatAssemblyBegin(Mat mat,MatAssemblyType type) 1087 { 1088 int ierr; 1089 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 1090 PLogEventBegin(MAT_AssemblyBegin,mat,0,0,0); 1091 if (mat->ops.assemblybegin) {ierr = (*mat->ops.assemblybegin)(mat,type); CHKERRQ(ierr);} 1092 PLogEventEnd(MAT_AssemblyBegin,mat,0,0,0); 1093 return 0; 1094 } 1095 1096 /*@ 1097 MatAssemblyEnd - Completes assembling the matrix. This routine should 1098 be called after all calls to MatSetValues() and after MatAssemblyBegin(). 1099 1100 Input Parameters: 1101 . mat - the matrix 1102 . type - type of assembly, either FLUSH_ASSEMBLY or FINAL_ASSEMBLY 1103 1104 Options Database Keys: 1105 $ -mat_view_draw : Draw nonzero structure of matrix at conclusion of MatEndAssembly(), 1106 using MatView() and DrawOpenX(). 1107 $ -mat_view_info : Prints info on matrix. 1108 $ -mat_view_info_detailed: More detailed information. 1109 $ -mat_view_ascii : Prints matrix out in ascii. 1110 $ -display <name> : Set display name (default is host) 1111 $ -draw_pause <sec> : Set number of seconds to pause after display 1112 1113 Note: 1114 MatSetValues() generally caches the values. The matrix is ready to 1115 use only after MatAssemblyBegin() and MatAssemblyEnd() have been called. 1116 Use FLUSH_ASSEMBLY when switching between ADD_VALUES and SetValues; use 1117 FINAL_ASSEMBLY for the final assembly before the matrix is used. 1118 1119 .keywords: matrix, assembly, assemble, end 1120 1121 .seealso: MatAssemblyBegin(), MatSetValues() 1122 @*/ 1123 int MatAssemblyEnd(Mat mat,MatAssemblyType type) 1124 { 1125 int ierr,flg; 1126 static int inassm = 0; 1127 1128 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 1129 inassm++; 1130 PLogEventBegin(MAT_AssemblyEnd,mat,0,0,0); 1131 if (mat->ops.assemblyend) {ierr = (*mat->ops.assemblyend)(mat,type); CHKERRQ(ierr);} 1132 PLogEventEnd(MAT_AssemblyEnd,mat,0,0,0); 1133 if (inassm == 1) { 1134 ierr = OptionsHasName(PETSC_NULL,"-mat_view_info",&flg); CHKERRQ(ierr); 1135 if (flg) { 1136 Viewer viewer; 1137 ierr = ViewerFileOpenASCII(mat->comm,"stdout",&viewer);CHKERRQ(ierr); 1138 ierr = ViewerFileSetFormat(viewer,FILE_FORMAT_INFO,0);CHKERRQ(ierr); 1139 ierr = MatView(mat,viewer); CHKERRQ(ierr); 1140 ierr = ViewerDestroy(viewer); CHKERRQ(ierr); 1141 } 1142 ierr = OptionsHasName(PETSC_NULL,"-mat_view_info_detailed",&flg); CHKERRQ(ierr); 1143 if (flg) { 1144 Viewer viewer; 1145 ierr = ViewerFileOpenASCII(mat->comm,"stdout",&viewer);CHKERRQ(ierr); 1146 ierr = ViewerFileSetFormat(viewer,FILE_FORMAT_INFO_DETAILED,0);CHKERRQ(ierr); 1147 ierr = MatView(mat,viewer); CHKERRQ(ierr); 1148 ierr = ViewerDestroy(viewer); CHKERRQ(ierr); 1149 } 1150 ierr = OptionsHasName(PETSC_NULL,"-mat_view_draw",&flg); CHKERRQ(ierr); 1151 if (flg) { 1152 Draw win; 1153 ierr = DrawOpenX(mat->comm,0,0,0,0,300,300,&win); CHKERRQ(ierr); 1154 ierr = MatView(mat,(Viewer)win); CHKERRQ(ierr); 1155 ierr = DrawSyncFlush(win); CHKERRQ(ierr); 1156 ierr = DrawDestroy(win); CHKERRQ(ierr); 1157 } 1158 } 1159 inassm--; 1160 return 0; 1161 } 1162 1163 /*@ 1164 MatCompress - Tries to store the matrix in as little space as 1165 possible. May fail if memory is already fully used, since it 1166 tries to allocate new space. 1167 1168 Input Parameters: 1169 . mat - the matrix 1170 1171 .keywords: matrix, compress 1172 @*/ 1173 int MatCompress(Mat mat) 1174 { 1175 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 1176 if (mat->ops.compress) return (*mat->ops.compress)(mat); 1177 return 0; 1178 } 1179 /*@ 1180 MatSetOption - Sets a parameter option for a matrix. Some options 1181 may be specific to certain storage formats. Some options 1182 determine how values will be inserted (or added). Sorted, 1183 row-oriented input will generally assemble the fastest. The default 1184 is row-oriented, nonsorted input. 1185 1186 Input Parameters: 1187 . mat - the matrix 1188 . option - the option, one of the following: 1189 $ ROW_ORIENTED 1190 $ COLUMN_ORIENTED, 1191 $ ROWS_SORTED, 1192 $ COLUMNS_SORTED, 1193 $ NO_NEW_NONZERO_LOCATIONS, 1194 $ YES_NEW_NONZERO_LOCATIONS, 1195 $ SYMMETRIC_MATRIX, 1196 $ STRUCTURALLY_SYMMETRIC_MATRIX, 1197 $ NO_NEW_DIAGONALS, 1198 $ YES_NEW_DIAGONALS, 1199 $ and possibly others. 1200 1201 Notes: 1202 Some options are relevant only for particular matrix types and 1203 are thus ignored by others. Other options are not supported by 1204 certain matrix types and will generate an error message if set. 1205 1206 If using a Fortran 77 module to compute a matrix, one may need to 1207 use the column-oriented option (or convert to the row-oriented 1208 format). 1209 1210 NO_NEW_NONZERO_LOCATIONS indicates that any add or insertion 1211 that will generate a new entry in the nonzero structure is ignored. 1212 What this means is if memory is not allocated for this particular 1213 lot, then the insertion is ignored. For dense matrices, where 1214 the entire array is allocated, no entries are ever ignored. 1215 1216 .keywords: matrix, option, row-oriented, column-oriented, sorted, nonzero 1217 @*/ 1218 int MatSetOption(Mat mat,MatOption op) 1219 { 1220 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 1221 if (mat->ops.setoption) return (*mat->ops.setoption)(mat,op); 1222 return 0; 1223 } 1224 1225 /*@ 1226 MatZeroEntries - Zeros all entries of a matrix. For sparse matrices 1227 this routine retains the old nonzero structure. 1228 1229 Input Parameters: 1230 . mat - the matrix 1231 1232 .keywords: matrix, zero, entries 1233 1234 .seealso: MatZeroRows() 1235 @*/ 1236 int MatZeroEntries(Mat mat) 1237 { 1238 int ierr; 1239 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 1240 if (!mat->ops.zeroentries) SETERRQ(PETSC_ERR_SUP,"MatZeroEntries"); 1241 PLogEventBegin(MAT_ZeroEntries,mat,0,0,0); 1242 ierr = (*mat->ops.zeroentries)(mat); CHKERRQ(ierr); 1243 PLogEventEnd(MAT_ZeroEntries,mat,0,0,0); 1244 return 0; 1245 } 1246 1247 /*@ 1248 MatZeroRows - Zeros all entries (except possibly the main diagonal) 1249 of a set of rows of a matrix. 1250 1251 Input Parameters: 1252 . mat - the matrix 1253 . is - index set of rows to remove 1254 . diag - pointer to value put in all diagonals of eliminated rows. 1255 Note that diag is not a pointer to an array, but merely a 1256 pointer to a single value. 1257 1258 Notes: 1259 For the AIJ and row matrix formats this removes the old nonzero 1260 structure, but does not release memory. For the dense and block 1261 diagonal formats this does not alter the nonzero structure. 1262 1263 The user can set a value in the diagonal entry (or for the AIJ and 1264 row formats can optionally remove the main diagonal entry from the 1265 nonzero structure as well, by passing a null pointer as the final 1266 argument). 1267 1268 .keywords: matrix, zero, rows, boundary conditions 1269 1270 .seealso: MatZeroEntries(), MatGetSubMatrix(), MatGetSubMatrixInPlace() 1271 @*/ 1272 int MatZeroRows(Mat mat,IS is, Scalar *diag) 1273 { 1274 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 1275 if (mat->ops.zerorows) return (*mat->ops.zerorows)(mat,is,diag); 1276 SETERRQ(PETSC_ERR_SUP,"MatZeroRows"); 1277 } 1278 1279 /*@ 1280 MatGetSize - Returns the numbers of rows and columns in a matrix. 1281 1282 Input Parameter: 1283 . mat - the matrix 1284 1285 Output Parameters: 1286 . m - the number of global rows 1287 . n - the number of global columns 1288 1289 .keywords: matrix, dimension, size, rows, columns, global, get 1290 1291 .seealso: MatGetLocalSize() 1292 @*/ 1293 int MatGetSize(Mat mat,int *m,int* n) 1294 { 1295 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 1296 if (!m || !n) SETERRQ(1,"MatGetSize:Bad address for result"); 1297 return (*mat->ops.getsize)(mat,m,n); 1298 } 1299 1300 /*@ 1301 MatGetLocalSize - Returns the number of rows and columns in a matrix 1302 stored locally. This information may be implementation dependent, so 1303 use with care. 1304 1305 Input Parameters: 1306 . mat - the matrix 1307 1308 Output Parameters: 1309 . m - the number of local rows 1310 . n - the number of local columns 1311 1312 .keywords: matrix, dimension, size, local, rows, columns, get 1313 1314 .seealso: MatGetSize() 1315 @*/ 1316 int MatGetLocalSize(Mat mat,int *m,int* n) 1317 { 1318 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 1319 if (!m || !n) SETERRQ(1,"MatGetLocalSize:Bad address for result"); 1320 return (*mat->ops.getlocalsize)(mat,m,n); 1321 } 1322 1323 /*@ 1324 MatGetOwnershipRange - Returns the range of matrix rows owned by 1325 this processor, assuming that the matrix is laid out with the first 1326 n1 rows on the first processor, the next n2 rows on the second, etc. 1327 For certain parallel layouts this range may not be well-defined. 1328 1329 Input Parameters: 1330 . mat - the matrix 1331 1332 Output Parameters: 1333 . m - the first local row 1334 . n - one more then the last local row 1335 1336 .keywords: matrix, get, range, ownership 1337 @*/ 1338 int MatGetOwnershipRange(Mat mat,int *m,int* n) 1339 { 1340 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 1341 if (!m || !n) SETERRQ(1,"MatGetOwnershipRange:Bad address for result"); 1342 if (mat->ops.getownershiprange) return (*mat->ops.getownershiprange)(mat,m,n); 1343 SETERRQ(PETSC_ERR_SUP,"MatGetOwnershipRange"); 1344 } 1345 1346 /*@ 1347 MatILUFactorSymbolic - Performs symbolic ILU factorization of a matrix. 1348 Uses levels of fill only, not drop tolerance. Use MatLUFactorNumeric() 1349 to complete the factorization. 1350 1351 Input Parameters: 1352 . mat - the matrix 1353 . row - row permutation 1354 . column - column permutation 1355 . fill - number of levels of fill 1356 . f - expected fill as ratio of original fill 1357 1358 Output Parameters: 1359 . fact - puts factor 1360 1361 .keywords: matrix, factor, incomplete, ILU, symbolic, fill 1362 1363 .seealso: MatLUFactorSymbolic(), MatLUFactorNumeric() 1364 @*/ 1365 int MatILUFactorSymbolic(Mat mat,IS row,IS col,double f,int fill,Mat *fact) 1366 { 1367 int ierr,flg; 1368 1369 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 1370 if (fill < 0) SETERRQ(1,"MatILUFactorSymbolic:Levels of fill negative"); 1371 if (!fact) SETERRQ(1,"MatILUFactorSymbolic:Fact argument is missing"); 1372 if (!mat->ops.ilufactorsymbolic) SETERRQ(PETSC_ERR_SUP,"MatILUFactorSymbolic"); 1373 ierr = OptionsGetDouble(PETSC_NULL,"-mat_ilu_fill",&f,&flg); CHKERRQ(ierr); 1374 PLogEventBegin(MAT_ILUFactorSymbolic,mat,row,col,0); 1375 ierr = (*mat->ops.ilufactorsymbolic)(mat,row,col,f,fill,fact); CHKERRQ(ierr); 1376 PLogEventEnd(MAT_ILUFactorSymbolic,mat,row,col,0); 1377 return 0; 1378 } 1379 1380 /*@ 1381 MatIncompleteCholeskyFactorSymbolic - Performs symbolic incomplete 1382 Cholesky factorization for a symmetric matrix. Use 1383 MatCholeskyFactorNumeric() to complete the factorization. 1384 1385 Input Parameters: 1386 . mat - the matrix 1387 . perm - row and column permutation 1388 . fill - levels of fill 1389 . f - expected fill as ratio of original fill 1390 1391 Output Parameter: 1392 . fact - the factored matrix 1393 1394 Note: Currently only no-fill factorization is supported. 1395 1396 .keywords: matrix, factor, incomplete, ICC, Cholesky, symbolic, fill 1397 1398 .seealso: MatCholeskyFactorNumeric(), MatCholeskyFactor() 1399 @*/ 1400 int MatIncompleteCholeskyFactorSymbolic(Mat mat,IS perm,double f,int fill, 1401 Mat *fact) 1402 { 1403 int ierr; 1404 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 1405 if (fill < 0) SETERRQ(1,"MatIncompleteCholeskyFactorSymbolic:Fill negative"); 1406 if (!fact) SETERRQ(1,"MatIncompleteCholeskyFactorSymbolic:Missing fact argument"); 1407 if (!mat->ops.incompletecholeskyfactorsymbolic) 1408 SETERRQ(PETSC_ERR_SUP,"MatIncompleteCholeskyFactorSymbolic"); 1409 PLogEventBegin(MAT_IncompleteCholeskyFactorSymbolic,mat,perm,0,0); 1410 ierr = (*mat->ops.incompletecholeskyfactorsymbolic)(mat,perm,f,fill,fact); 1411 CHKERRQ(ierr); 1412 PLogEventEnd(MAT_IncompleteCholeskyFactorSymbolic,mat,perm,0,0); 1413 return 0; 1414 } 1415 1416 /*@C 1417 MatGetArray - Returns a pointer to the element values in the matrix. 1418 This routine is implementation dependent, and may not even work for 1419 certain matrix types. 1420 1421 Input Parameter: 1422 . mat - the matrix 1423 1424 Output Parameter: 1425 . v - the location of the values 1426 1427 Fortran Note: 1428 The Fortran interface is slightly different from that given below. 1429 See the users manual and petsc/src/mat/examples for details. 1430 1431 .keywords: matrix, array, elements, values 1432 @*/ 1433 int MatGetArray(Mat mat,Scalar **v) 1434 { 1435 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 1436 if (!v) SETERRQ(1,"MatGetArray:Bad input, array pointer location"); 1437 if (!mat->ops.getarray) SETERRQ(PETSC_ERR_SUP,"MatGetArraye"); 1438 return (*mat->ops.getarray)(mat,v); 1439 } 1440 1441 /*@C 1442 MatGetSubMatrix - Extracts a submatrix from a matrix. If submat points 1443 to a valid matrix, it may be reused. 1444 1445 Input Parameters: 1446 . mat - the matrix 1447 . irow, icol - index sets of rows and columns to extract 1448 . scall - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX 1449 1450 Output Parameter: 1451 . submat - the submatrix 1452 1453 Notes: 1454 MatGetSubMatrix() can be useful in setting boundary conditions. 1455 1456 Use MatGetSubMatrices() to extract multiple submatrices. 1457 1458 .keywords: matrix, get, submatrix, boundary conditions 1459 1460 .seealso: MatZeroRows(), MatGetSubMatrixInPlace(), MatGetSubMatrices() 1461 @*/ 1462 int MatGetSubMatrix(Mat mat,IS irow,IS icol,MatGetSubMatrixCall scall,Mat *submat) 1463 { 1464 int ierr; 1465 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 1466 if (scall == MAT_REUSE_MATRIX) { 1467 PETSCVALIDHEADERSPECIFIC(*submat,MAT_COOKIE); 1468 } 1469 if (!mat->ops.getsubmatrix) SETERRQ(PETSC_ERR_SUP,"MatGetSubMatrix"); 1470 PLogEventBegin(MAT_GetSubMatrix,mat,irow,icol,0); 1471 ierr = (*mat->ops.getsubmatrix)(mat,irow,icol,scall,submat); CHKERRQ(ierr); 1472 PLogEventEnd(MAT_GetSubMatrix,mat,irow,icol,0); 1473 return 0; 1474 } 1475 1476 /*@C 1477 MatGetSubMatrices - Extracts several submatrices from a matrix. If submat 1478 points to an array of valid matrices, it may be reused. 1479 1480 Input Parameters: 1481 . mat - the matrix 1482 . irow, icol - index sets of rows and columns to extract 1483 1484 Output Parameter: 1485 . submat - the submatrices 1486 1487 Note: 1488 Use MatGetSubMatrix() for extracting a sinble submatrix. 1489 1490 .keywords: matrix, get, submatrix, submatrices 1491 1492 .seealso: MatGetSubMatrix() 1493 @*/ 1494 int MatGetSubMatrices(Mat mat,int n, IS *irow,IS *icol,MatGetSubMatrixCall scall, 1495 Mat **submat) 1496 { 1497 int ierr; 1498 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 1499 if (!mat->ops.getsubmatrices) SETERRQ(PETSC_ERR_SUP,"MatGetSubMatrices"); 1500 PLogEventBegin(MAT_GetSubMatrices,mat,0,0,0); 1501 ierr = (*mat->ops.getsubmatrices)(mat,n,irow,icol,scall,submat); CHKERRQ(ierr); 1502 PLogEventEnd(MAT_GetSubMatrices,mat,0,0,0); 1503 return 0; 1504 } 1505 1506 /*@ 1507 MatGetSubMatrixInPlace - Extracts a submatrix from a matrix, returning 1508 the submatrix in place of the original matrix. 1509 1510 Input Parameters: 1511 . mat - the matrix 1512 . irow, icol - index sets of rows and columns to extract 1513 1514 .keywords: matrix, get, submatrix, boundary conditions, in-place 1515 1516 .seealso: MatZeroRows(), MatGetSubMatrix() 1517 @*/ 1518 int MatGetSubMatrixInPlace(Mat mat,IS irow,IS icol) 1519 { 1520 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 1521 if (!mat->ops.getsubmatrixinplace) SETERRQ(PETSC_ERR_SUP,"MatGetSubmatrixInPlace"); 1522 return (*mat->ops.getsubmatrixinplace)(mat,irow,icol); 1523 } 1524 1525 /*@ 1526 MatIncreaseOverlap - Given a set of submatrices indicated by index sets, 1527 replaces the index by larger ones that represent submatrices with more 1528 overlap. 1529 1530 Input Parameters: 1531 . mat - the matrix 1532 . n - the number of index sets 1533 . is - the array of pointers to index sets 1534 . ov - the additional overlap requested 1535 1536 .keywords: matrix, overlap, Schwarz 1537 1538 .seealso: MatGetSubMatrices() 1539 @*/ 1540 int MatIncreaseOverlap(Mat mat,int n, IS *is, int ov) 1541 { 1542 int ierr; 1543 PETSCVALIDHEADERSPECIFIC(mat,MAT_COOKIE); 1544 if (!mat->ops.increaseoverlap) SETERRQ(PETSC_ERR_SUP,"MatIncreaseOverlap"); 1545 ierr = (*mat->ops.increaseoverlap)(mat,n,is,ov); CHKERRQ(ierr); 1546 return 0; 1547 } 1548 1549 /*@ 1550 MatPrintHelp - Prints all the options for the matrix. 1551 1552 Input Parameter: 1553 . mat - the matrix 1554 1555 Options Database Keys: 1556 $ -help, -h 1557 1558 .keywords: mat, help 1559 1560 .seealso: MatCreate(), MatCreateXXX() 1561 @*/ 1562 int MatPrintHelp(Mat mat) 1563 { 1564 static int called = 0; 1565 MPI_Comm comm = mat->comm; 1566 1567 if (!called) { 1568 MPIU_printf(comm,"General matrix options:\n"); 1569 MPIU_printf(comm," -mat_view_info : view basic matrix info during MatAssemblyEnd()\n"); 1570 MPIU_printf(comm," -mat_view_info_detailed : view detailed matrix info during MatAssemblyEnd()\n"); 1571 MPIU_printf(comm," -mat_view_draw : draw nonzero matrix structure during MatAssemblyEnd()\n"); 1572 MPIU_printf(comm," -draw_pause <sec> : set seconds of display pause\n"); 1573 MPIU_printf(comm," -display <name> : set alternate display\n"); 1574 called = 1; 1575 } 1576 if (mat->ops.printhelp) (*mat->ops.printhelp)(mat); 1577 return 0; 1578 } 1579