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