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