1 #ifdef PETSC_RCS_HEADER 2 static char vcid[] = "$Id: matrix.c,v 1.345 1999/09/28 21:28:34 balay Exp bsmith $"; 3 #endif 4 5 /* 6 This is where the abstract matrix operations are defined 7 */ 8 9 #include "src/mat/matimpl.h" /*I "mat.h" I*/ 10 #include "src/vec/vecimpl.h" 11 12 #undef __FUNC__ 13 #define __FUNC__ "MatGetRow" 14 /*@C 15 MatGetRow - Gets a row of a matrix. You MUST call MatRestoreRow() 16 for each row that you get to ensure that your application does 17 not bleed memory. 18 19 Not Collective 20 21 Input Parameters: 22 + mat - the matrix 23 - row - the row to get 24 25 Output Parameters: 26 + ncols - the number of nonzeros in the row 27 . cols - if nonzero, the column numbers 28 - vals - if nonzero, the values 29 30 Notes: 31 This routine is provided for people who need to have direct access 32 to the structure of a matrix. We hope that we provide enough 33 high-level matrix routines that few users will need it. 34 35 MatGetRow() always returns 0-based column indices, regardless of 36 whether the internal representation is 0-based (default) or 1-based. 37 38 For better efficiency, set cols and/or vals to PETSC_NULL if you do 39 not wish to extract these quantities. 40 41 The user can only examine the values extracted with MatGetRow(); 42 the values cannot be altered. To change the matrix entries, one 43 must use MatSetValues(). 44 45 You can only have one call to MatGetRow() outstanding for a particular 46 matrix at a time. 47 48 Fortran Notes: 49 The calling sequence from Fortran is 50 .vb 51 MatGetRow(matrix,row,ncols,cols,values,ierr) 52 Mat matrix (input) 53 integer row (input) 54 integer ncols (output) 55 integer cols(maxcols) (output) 56 double precision (or double complex) values(maxcols) output 57 .ve 58 where maxcols >= maximum nonzeros in any row of the matrix. 59 60 Caution: 61 Do not try to change the contents of the output arrays (cols and vals). 62 In some cases, this may corrupt the matrix. 63 64 Level: advanced 65 66 .keywords: matrix, row, get, extract 67 68 .seealso: MatRestoreRow(), MatSetValues(), MatGetValues(), MatGetSubmatrices(), MatGetDiagonal() 69 @*/ 70 int MatGetRow(Mat mat,int row,int *ncols,int **cols,Scalar **vals) 71 { 72 int ierr; 73 74 PetscFunctionBegin; 75 PetscValidHeaderSpecific(mat,MAT_COOKIE); 76 PetscValidIntPointer(ncols); 77 if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix"); 78 if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix"); 79 if (!mat->ops->getrow) SETERRQ(PETSC_ERR_SUP,0,""); 80 PLogEventBegin(MAT_GetRow,mat,0,0,0); 81 ierr = (*mat->ops->getrow)(mat,row,ncols,cols,vals);CHKERRQ(ierr); 82 PLogEventEnd(MAT_GetRow,mat,0,0,0); 83 PetscFunctionReturn(0); 84 } 85 86 #undef __FUNC__ 87 #define __FUNC__ "MatRestoreRow" 88 /*@C 89 MatRestoreRow - Frees any temporary space allocated by MatGetRow(). 90 91 Not Collective 92 93 Input Parameters: 94 + mat - the matrix 95 . row - the row to get 96 . ncols, cols - the number of nonzeros and their columns 97 - vals - if nonzero the column values 98 99 Notes: 100 This routine should be called after you have finished examining the entries. 101 102 Fortran Notes: 103 The calling sequence from Fortran is 104 .vb 105 MatRestoreRow(matrix,row,ncols,cols,values,ierr) 106 Mat matrix (input) 107 integer row (input) 108 integer ncols (output) 109 integer cols(maxcols) (output) 110 double precision (or double complex) values(maxcols) output 111 .ve 112 Where maxcols >= maximum nonzeros in any row of the matrix. 113 114 In Fortran MatRestoreRow() MUST be called after MatGetRow() 115 before another call to MatGetRow() can be made. 116 117 Level: advanced 118 119 .keywords: matrix, row, restore 120 121 .seealso: MatGetRow() 122 @*/ 123 int MatRestoreRow(Mat mat,int row,int *ncols,int **cols,Scalar **vals) 124 { 125 int ierr; 126 127 PetscFunctionBegin; 128 PetscValidHeaderSpecific(mat,MAT_COOKIE); 129 PetscValidIntPointer(ncols); 130 if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix"); 131 if (!mat->ops->restorerow) PetscFunctionReturn(0); 132 ierr = (*mat->ops->restorerow)(mat,row,ncols,cols,vals);CHKERRQ(ierr); 133 PetscFunctionReturn(0); 134 } 135 136 #undef __FUNC__ 137 #define __FUNC__ "MatView" 138 /*@C 139 MatView - Visualizes a matrix object. 140 141 Collective on Mat unless Viewer is VIEWER_STDOUT_SELF 142 143 Input Parameters: 144 + mat - the matrix 145 - ptr - visualization context 146 147 Notes: 148 The available visualization contexts include 149 + VIEWER_STDOUT_SELF - standard output (default) 150 . VIEWER_STDOUT_WORLD - synchronized standard 151 output where only the first processor opens 152 the file. All other processors send their 153 data to the first processor to print. 154 - VIEWER_DRAW_WORLD - graphical display of nonzero structure 155 156 The user can open alternative visualization contexts with 157 + ViewerASCIIOpen() - Outputs matrix to a specified file 158 . ViewerBinaryOpen() - Outputs matrix in binary to a 159 specified file; corresponding input uses MatLoad() 160 . ViewerDrawOpen() - Outputs nonzero matrix structure to 161 an X window display 162 - ViewerSocketOpen() - Outputs matrix to Socket viewer. 163 Currently only the sequential dense and AIJ 164 matrix types support the Socket viewer. 165 166 The user can call ViewerSetFormat() to specify the output 167 format of ASCII printed objects (when using VIEWER_STDOUT_SELF, 168 VIEWER_STDOUT_WORLD and ViewerASCIIOpen). Available formats include 169 + VIEWER_FORMAT_ASCII_DEFAULT - default, prints matrix contents 170 . VIEWER_FORMAT_ASCII_MATLAB - prints matrix contents in Matlab format 171 . VIEWER_FORMAT_ASCII_DENSE - prints entire matrix including zeros 172 . VIEWER_FORMAT_ASCII_COMMON - prints matrix contents, using a sparse 173 format common among all matrix types 174 . VIEWER_FORMAT_ASCII_IMPL - prints matrix contents, using an implementation-specific 175 format (which is in many cases the same as the default) 176 . VIEWER_FORMAT_ASCII_INFO - prints basic information about the matrix 177 size and structure (not the matrix entries) 178 - VIEWER_FORMAT_ASCII_INFO_LONG - prints more detailed information about 179 the matrix structure 180 181 Level: beginner 182 183 .keywords: matrix, view, visualize, output, print, write, draw 184 185 .seealso: ViewerSetFormat(), ViewerASCIIOpen(), ViewerDrawOpen(), 186 ViewerSocketOpen(), ViewerBinaryOpen(), MatLoad() 187 @*/ 188 int MatView(Mat mat,Viewer viewer) 189 { 190 int format, ierr, rows, cols; 191 char *cstr; 192 ViewerType vtype; 193 194 PetscFunctionBegin; 195 PetscValidHeaderSpecific(mat,MAT_COOKIE); 196 if (viewer) PetscValidHeaderSpecific(viewer,VIEWER_COOKIE); 197 198 if (!viewer) { 199 viewer = VIEWER_STDOUT_SELF; 200 } 201 202 ierr = ViewerGetType(viewer,&vtype);CHKERRQ(ierr); 203 if (PetscTypeCompare(vtype,ASCII_VIEWER)) { 204 ierr = ViewerGetFormat(viewer,&format);CHKERRQ(ierr); 205 if (format == VIEWER_FORMAT_ASCII_INFO || format == VIEWER_FORMAT_ASCII_INFO_LONG) { 206 ierr = ViewerASCIIPrintf(viewer,"Matrix Object:\n");CHKERRQ(ierr); 207 ierr = MatGetType(mat,PETSC_NULL,&cstr);CHKERRQ(ierr); 208 ierr = MatGetSize(mat,&rows,&cols);CHKERRQ(ierr); 209 ierr = ViewerASCIIPrintf(viewer," type=%s, rows=%d, cols=%d\n",cstr,rows,cols);CHKERRQ(ierr); 210 if (mat->ops->getinfo) { 211 MatInfo info; 212 ierr = MatGetInfo(mat,MAT_GLOBAL_SUM,&info);CHKERRQ(ierr); 213 ierr = ViewerASCIIPrintf(viewer," total: nonzeros=%d, allocated nonzeros=%d\n", 214 (int)info.nz_used,(int)info.nz_allocated);CHKERRQ(ierr); 215 } 216 } 217 } 218 if (mat->ops->view) { 219 ierr = ViewerASCIIPushTab(viewer);CHKERRQ(ierr); 220 ierr = (*mat->ops->view)(mat,viewer);CHKERRQ(ierr); 221 ierr = ViewerASCIIPopTab(viewer);CHKERRQ(ierr); 222 } 223 PetscFunctionReturn(0); 224 } 225 226 #undef __FUNC__ 227 #define __FUNC__ "MatScaleSystem" 228 /*@C 229 MatScaleSystem - Scale a vector solution and right hand side to 230 match the scaling of a scaled matrix. 231 232 Collective on Mat 233 234 Input Parameter: 235 + mat - the matrix 236 . x - solution vector (or PETSC_NULL) 237 + b - right hand side vector (or PETSC_NULL) 238 239 240 Notes: 241 For AIJ, BAIJ, and BDiag matrix formats, the matrices are not 242 internally scaled, so this does nothing. For MPIROWBS it 243 permutes and diagonally scales. 244 245 The SLES methods automatically call this routine when required 246 (via PCPreSolve()) so it is rarely used directly. 247 248 Level: Developer 249 250 .keywords: matrix, scale 251 252 .seealso: MatUseScaledForm(), MatUnScaleSystem() 253 @*/ 254 int MatScaleSystem(Mat mat,Vec x,Vec b) 255 { 256 int ierr; 257 258 PetscFunctionBegin; 259 PetscValidHeaderSpecific(mat,MAT_COOKIE); 260 if (x) {PetscValidHeaderSpecific(x,VEC_COOKIE);PetscCheckSameComm(mat,x);} 261 if (b) {PetscValidHeaderSpecific(b,VEC_COOKIE);PetscCheckSameComm(mat,b);} 262 263 if (mat->ops->scalesystem) { 264 ierr = (*mat->ops->scalesystem)(mat,x,b);CHKERRQ(ierr); 265 } 266 PetscFunctionReturn(0); 267 } 268 269 #undef __FUNC__ 270 #define __FUNC__ "MatUnScaleSystem" 271 /*@C 272 MatUnScaleSystem - Unscales a vector solution and right hand side to 273 match the original scaling of a scaled matrix. 274 275 Collective on Mat 276 277 Input Parameter: 278 + mat - the matrix 279 . x - solution vector (or PETSC_NULL) 280 + b - right hand side vector (or PETSC_NULL) 281 282 283 Notes: 284 For AIJ, BAIJ, and BDiag matrix formats, the matrices are not 285 internally scaled, so this does nothing. For MPIROWBS it 286 permutes and diagonally scales. 287 288 The SLES methods automatically call this routine when required 289 (via PCPreSolve()) so it is rarely used directly. 290 291 Level: Developer 292 293 .keywords: matrix, scale 294 295 .seealso: MatUseScaledForm(), MatScaleSystem() 296 @*/ 297 int MatUnScaleSystem(Mat mat,Vec x,Vec b) 298 { 299 int ierr; 300 301 PetscFunctionBegin; 302 PetscValidHeaderSpecific(mat,MAT_COOKIE); 303 if (x) {PetscValidHeaderSpecific(x,VEC_COOKIE);PetscCheckSameComm(mat,x);} 304 if (b) {PetscValidHeaderSpecific(b,VEC_COOKIE);PetscCheckSameComm(mat,x);} 305 if (mat->ops->unscalesystem) { 306 ierr = (*mat->ops->unscalesystem)(mat,x,b);CHKERRQ(ierr); 307 } 308 PetscFunctionReturn(0); 309 } 310 311 #undef __FUNC__ 312 #define __FUNC__ "MatUseScaledForm" 313 /*@C 314 MatUseScaledForm - For matrix storage formats that scale the 315 matrix (for example MPIRowBS matrices are diagonally scaled on 316 assembly) indicates matrix operations (MatMult() etc) are 317 applied using the scaled matrix. 318 319 Collective on Mat 320 321 Input Parameter: 322 + mat - the matrix 323 - scaled - PETSC_TRUE for applying the scaled, PETSC_FALSE for 324 applying the original matrix 325 326 Notes: 327 For scaled matrix formats, applying the original, unscaled matrix 328 will be slightly more expensive 329 330 Level: Developer 331 332 .keywords: matrix, scale 333 334 .seealso: MatScaleSystem(), MatUnScaleSystem() 335 @*/ 336 int MatUseScaledForm(Mat mat,PetscTruth scaled) 337 { 338 int ierr; 339 340 PetscFunctionBegin; 341 PetscValidHeaderSpecific(mat,MAT_COOKIE); 342 if (mat->ops->usescaledform) { 343 ierr = (*mat->ops->usescaledform)(mat,scaled);CHKERRQ(ierr); 344 } 345 PetscFunctionReturn(0); 346 } 347 348 #undef __FUNC__ 349 #define __FUNC__ "MatDestroy" 350 /*@C 351 MatDestroy - Frees space taken by a matrix. 352 353 Collective on Mat 354 355 Input Parameter: 356 . mat - the matrix 357 358 Level: beginner 359 360 .keywords: matrix, destroy 361 @*/ 362 int MatDestroy(Mat mat) 363 { 364 int ierr; 365 366 PetscFunctionBegin; 367 PetscValidHeaderSpecific(mat,MAT_COOKIE); 368 if (--mat->refct > 0) PetscFunctionReturn(0); 369 370 /* if memory was published with AMS then destroy it */ 371 ierr = PetscAMSDestroy(mat);CHKERRQ(ierr); 372 373 ierr = (*mat->ops->destroy)(mat);CHKERRQ(ierr); 374 PetscFunctionReturn(0); 375 } 376 377 #undef __FUNC__ 378 #define __FUNC__ "MatValid" 379 /*@ 380 MatValid - Checks whether a matrix object is valid. 381 382 Collective on Mat 383 384 Input Parameter: 385 . m - the matrix to check 386 387 Output Parameter: 388 flg - flag indicating matrix status, either 389 PETSC_TRUE if matrix is valid, or PETSC_FALSE otherwise. 390 391 Level: developer 392 393 .keywords: matrix, valid 394 @*/ 395 int MatValid(Mat m,PetscTruth *flg) 396 { 397 PetscFunctionBegin; 398 PetscValidIntPointer(flg); 399 if (!m) *flg = PETSC_FALSE; 400 else if (m->cookie != MAT_COOKIE) *flg = PETSC_FALSE; 401 else *flg = PETSC_TRUE; 402 PetscFunctionReturn(0); 403 } 404 405 #undef __FUNC__ 406 #define __FUNC__ "MatSetValues" 407 /*@ 408 MatSetValues - Inserts or adds a block of values into a matrix. 409 These values may be cached, so MatAssemblyBegin() and MatAssemblyEnd() 410 MUST be called after all calls to MatSetValues() have been completed. 411 412 Not Collective 413 414 Input Parameters: 415 + mat - the matrix 416 . v - a logically two-dimensional array of values 417 . m, idxm - the number of rows and their global indices 418 . n, idxn - the number of columns and their global indices 419 - addv - either ADD_VALUES or INSERT_VALUES, where 420 ADD_VALUES adds values to any existing entries, and 421 INSERT_VALUES replaces existing entries with new values 422 423 Notes: 424 By default the values, v, are row-oriented and unsorted. 425 See MatSetOption() for other options. 426 427 Calls to MatSetValues() with the INSERT_VALUES and ADD_VALUES 428 options cannot be mixed without intervening calls to the assembly 429 routines. 430 431 MatSetValues() uses 0-based row and column numbers in Fortran 432 as well as in C. 433 434 Negative indices may be passed in idxm and idxn, these rows and columns are 435 simply ignored. This allows easily inserting element stiffness matrices 436 with homogeneous Dirchlet boundary conditions that you don't want represented 437 in the matrix. 438 439 Efficiency Alert: 440 The routine MatSetValuesBlocked() may offer much better efficiency 441 for users of block sparse formats (MATSEQBAIJ and MATMPIBAIJ). 442 443 Level: beginner 444 445 .keywords: matrix, insert, add, set, values 446 447 .seealso: MatSetOption(), MatAssemblyBegin(), MatAssemblyEnd(), MatSetValuesBlocked(), MatSetValuesLocal() 448 @*/ 449 int MatSetValues(Mat mat,int m,int *idxm,int n,int *idxn,Scalar *v,InsertMode addv) 450 { 451 int ierr; 452 453 PetscFunctionBegin; 454 if (!m || !n) PetscFunctionReturn(0); /* no values to insert */ 455 PetscValidHeaderSpecific(mat,MAT_COOKIE); 456 PetscValidIntPointer(idxm); 457 PetscValidIntPointer(idxn); 458 PetscValidScalarPointer(v); 459 if (mat->insertmode == NOT_SET_VALUES) { 460 mat->insertmode = addv; 461 } 462 #if defined(PETSC_USE_BOPT_g) 463 else if (mat->insertmode != addv) { 464 SETERRQ(PETSC_ERR_ARG_WRONGSTATE,1,"Cannot mix add values and insert values"); 465 } 466 if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix"); 467 #endif 468 469 if (mat->assembled) { 470 mat->was_assembled = PETSC_TRUE; 471 mat->assembled = PETSC_FALSE; 472 } 473 PLogEventBegin(MAT_SetValues,mat,0,0,0); 474 if (!mat->ops->setvalues) SETERRQ(PETSC_ERR_SUP,1,"Not supported for this matrix type"); 475 ierr = (*mat->ops->setvalues)(mat,m,idxm,n,idxn,v,addv);CHKERRQ(ierr); 476 PLogEventEnd(MAT_SetValues,mat,0,0,0); 477 PetscFunctionReturn(0); 478 } 479 480 #undef __FUNC__ 481 #define __FUNC__ "MatSetValuesBlocked" 482 /*@ 483 MatSetValuesBlocked - Inserts or adds a block of values into a matrix. 484 485 Not Collective 486 487 Input Parameters: 488 + mat - the matrix 489 . v - a logically two-dimensional array of values 490 . m, idxm - the number of block rows and their global block indices 491 . n, idxn - the number of block columns and their global block indices 492 - addv - either ADD_VALUES or INSERT_VALUES, where 493 ADD_VALUES adds values to any existing entries, and 494 INSERT_VALUES replaces existing entries with new values 495 496 Notes: 497 By default the values, v, are row-oriented and unsorted. So the layout of 498 v is the same as for MatSetValues(). See MatSetOption() for other options. 499 500 Calls to MatSetValuesBlocked() with the INSERT_VALUES and ADD_VALUES 501 options cannot be mixed without intervening calls to the assembly 502 routines. 503 504 MatSetValuesBlocked() uses 0-based row and column numbers in Fortran 505 as well as in C. 506 507 Negative indices may be passed in idxm and idxn, these rows and columns are 508 simply ignored. This allows easily inserting element stiffness matrices 509 with homogeneous Dirchlet boundary conditions that you don't want represented 510 in the matrix. 511 512 Each time an entry is set within a sparse matrix via MatSetValues(), 513 internal searching must be done to determine where to place the the 514 data in the matrix storage space. By instead inserting blocks of 515 entries via MatSetValuesBlocked(), the overhead of matrix assembly is 516 reduced. 517 518 Restrictions: 519 MatSetValuesBlocked() is currently supported only for the block AIJ 520 matrix format (MATSEQBAIJ and MATMPIBAIJ, which are created via 521 MatCreateSeqBAIJ() and MatCreateMPIBAIJ()). 522 523 Level: intermediate 524 525 .keywords: matrix, insert, add, set, values 526 527 .seealso: MatSetOption(), MatAssemblyBegin(), MatAssemblyEnd(), MatSetValues(), MatSetValuesBlockedLocal() 528 @*/ 529 int MatSetValuesBlocked(Mat mat,int m,int *idxm,int n,int *idxn,Scalar *v,InsertMode addv) 530 { 531 int ierr; 532 533 PetscFunctionBegin; 534 if (!m || !n) PetscFunctionReturn(0); /* no values to insert */ 535 PetscValidHeaderSpecific(mat,MAT_COOKIE); 536 PetscValidIntPointer(idxm); 537 PetscValidIntPointer(idxn); 538 PetscValidScalarPointer(v); 539 if (mat->insertmode == NOT_SET_VALUES) { 540 mat->insertmode = addv; 541 } 542 #if defined(PETSC_USE_BOPT_g) 543 else if (mat->insertmode != addv) { 544 SETERRQ(PETSC_ERR_ARG_WRONGSTATE,1,"Cannot mix add values and insert values"); 545 } 546 if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix"); 547 #endif 548 549 if (mat->assembled) { 550 mat->was_assembled = PETSC_TRUE; 551 mat->assembled = PETSC_FALSE; 552 } 553 PLogEventBegin(MAT_SetValues,mat,0,0,0); 554 if (!mat->ops->setvaluesblocked) SETERRQ(PETSC_ERR_SUP,1,"Not supported for this matrix type"); 555 ierr = (*mat->ops->setvaluesblocked)(mat,m,idxm,n,idxn,v,addv);CHKERRQ(ierr); 556 PLogEventEnd(MAT_SetValues,mat,0,0,0); 557 PetscFunctionReturn(0); 558 } 559 560 /*MC 561 MatSetValue - Set a single entry into a matrix. 562 563 Synopsis: 564 void MatSetValue(Mat m,int row,int col,Scalar value,InsertMode mode); 565 566 Not collective 567 568 Input Parameters: 569 + m - the matrix 570 . row - the row location of the entry 571 . col - the column location of the entry 572 . value - the value to insert 573 - mode - either INSERT_VALUES or ADD_VALUES 574 575 Notes: 576 For efficiency one should use MatSetValues() and set several or many 577 values simultaneously if possible. 578 579 Note that VecSetValue() does NOT return an error code (since this 580 is checked internally). 581 582 Level: beginner 583 584 .seealso: MatSetValues() 585 M*/ 586 587 #undef __FUNC__ 588 #define __FUNC__ "MatGetValues" 589 /*@ 590 MatGetValues - Gets a block of values from a matrix. 591 592 Not Collective; currently only returns a local block 593 594 Input Parameters: 595 + mat - the matrix 596 . v - a logically two-dimensional array for storing the values 597 . m, idxm - the number of rows and their global indices 598 - n, idxn - the number of columns and their global indices 599 600 Notes: 601 The user must allocate space (m*n Scalars) for the values, v. 602 The values, v, are then returned in a row-oriented format, 603 analogous to that used by default in MatSetValues(). 604 605 MatGetValues() uses 0-based row and column numbers in 606 Fortran as well as in C. 607 608 MatGetValues() requires that the matrix has been assembled 609 with MatAssemblyBegin()/MatAssemblyEnd(). Thus, calls to 610 MatSetValues() and MatGetValues() CANNOT be made in succession 611 without intermediate matrix assembly. 612 613 Level: advanced 614 615 .keywords: matrix, get, values 616 617 .seealso: MatGetRow(), MatGetSubmatrices(), MatSetValues() 618 @*/ 619 int MatGetValues(Mat mat,int m,int *idxm,int n,int *idxn,Scalar *v) 620 { 621 int ierr; 622 623 PetscFunctionBegin; 624 PetscValidHeaderSpecific(mat,MAT_COOKIE); 625 PetscValidIntPointer(idxm); 626 PetscValidIntPointer(idxn); 627 PetscValidScalarPointer(v); 628 if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix"); 629 if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix"); 630 if (!mat->ops->getvalues) SETERRQ(PETSC_ERR_SUP,0,""); 631 632 PLogEventBegin(MAT_GetValues,mat,0,0,0); 633 ierr = (*mat->ops->getvalues)(mat,m,idxm,n,idxn,v);CHKERRQ(ierr); 634 PLogEventEnd(MAT_GetValues,mat,0,0,0); 635 PetscFunctionReturn(0); 636 } 637 638 #undef __FUNC__ 639 #define __FUNC__ "MatSetLocalToGlobalMapping" 640 /*@ 641 MatSetLocalToGlobalMapping - Sets a local-to-global numbering for use by 642 the routine MatSetValuesLocal() to allow users to insert matrix entries 643 using a local (per-processor) numbering. 644 645 Not Collective 646 647 Input Parameters: 648 + x - the matrix 649 - mapping - mapping created with ISLocalToGlobalMappingCreate() 650 or ISLocalToGlobalMappingCreateIS() 651 652 Level: intermediate 653 654 .keywords: matrix, set, values, local, global, mapping 655 656 .seealso: MatAssemblyBegin(), MatAssemblyEnd(), MatSetValues(), MatSetValuesLocal() 657 @*/ 658 int MatSetLocalToGlobalMapping(Mat x,ISLocalToGlobalMapping mapping) 659 { 660 PetscFunctionBegin; 661 PetscValidHeaderSpecific(x,MAT_COOKIE); 662 PetscValidHeaderSpecific(mapping,IS_LTOGM_COOKIE); 663 if (x->mapping) { 664 SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Mapping already set for matrix"); 665 } 666 667 x->mapping = mapping; 668 PetscObjectReference((PetscObject)mapping); 669 PetscFunctionReturn(0); 670 } 671 672 #undef __FUNC__ 673 #define __FUNC__ "MatSetLocalToGlobalMappingBlock" 674 /*@ 675 MatSetLocalToGlobalMappingBlock - Sets a local-to-global numbering for use 676 by the routine MatSetValuesBlockedLocal() to allow users to insert matrix 677 entries using a local (per-processor) numbering. 678 679 Not Collective 680 681 Input Parameters: 682 + x - the matrix 683 - mapping - mapping created with ISLocalToGlobalMappingCreate() or 684 ISLocalToGlobalMappingCreateIS() 685 686 Level: intermediate 687 688 .keywords: matrix, set, values, local ordering 689 690 .seealso: MatAssemblyBegin(), MatAssemblyEnd(), MatSetValues(), MatSetValuesBlockedLocal(), 691 MatSetValuesBlocked(), MatSetValuesLocal() 692 @*/ 693 int MatSetLocalToGlobalMappingBlock(Mat x,ISLocalToGlobalMapping mapping) 694 { 695 PetscFunctionBegin; 696 PetscValidHeaderSpecific(x,MAT_COOKIE); 697 PetscValidHeaderSpecific(mapping,IS_LTOGM_COOKIE); 698 if (x->bmapping) { 699 SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Mapping already set for matrix"); 700 } 701 702 x->bmapping = mapping; 703 PetscObjectReference((PetscObject)mapping); 704 PetscFunctionReturn(0); 705 } 706 707 #undef __FUNC__ 708 #define __FUNC__ "MatSetValuesLocal" 709 /*@ 710 MatSetValuesLocal - Inserts or adds values into certain locations of a matrix, 711 using a local ordering of the nodes. 712 713 Not Collective 714 715 Input Parameters: 716 + x - the matrix 717 . nrow, irow - number of rows and their local indices 718 . ncol, icol - number of columns and their local indices 719 . y - a logically two-dimensional array of values 720 - addv - either INSERT_VALUES or ADD_VALUES, where 721 ADD_VALUES adds values to any existing entries, and 722 INSERT_VALUES replaces existing entries with new values 723 724 Notes: 725 Before calling MatSetValuesLocal(), the user must first set the 726 local-to-global mapping by calling MatSetLocalToGlobalMapping(). 727 728 Calls to MatSetValuesLocal() with the INSERT_VALUES and ADD_VALUES 729 options cannot be mixed without intervening calls to the assembly 730 routines. 731 732 These values may be cached, so MatAssemblyBegin() and MatAssemblyEnd() 733 MUST be called after all calls to MatSetValuesLocal() have been completed. 734 735 Level: intermediate 736 737 .keywords: matrix, set, values, local ordering 738 739 .seealso: MatAssemblyBegin(), MatAssemblyEnd(), MatSetValues(), MatSetLocalToGlobalMapping() 740 @*/ 741 int MatSetValuesLocal(Mat mat,int nrow,int *irow,int ncol, int *icol,Scalar *y,InsertMode addv) 742 { 743 int ierr,irowm[2048],icolm[2048]; 744 745 PetscFunctionBegin; 746 PetscValidHeaderSpecific(mat,MAT_COOKIE); 747 PetscValidIntPointer(irow); 748 PetscValidIntPointer(icol); 749 PetscValidScalarPointer(y); 750 751 if (mat->insertmode == NOT_SET_VALUES) { 752 mat->insertmode = addv; 753 } 754 #if defined(PETSC_USE_BOPT_g) 755 else if (mat->insertmode != addv) { 756 SETERRQ(PETSC_ERR_ARG_WRONGSTATE,1,"Cannot mix add values and insert values"); 757 } 758 if (!mat->mapping) { 759 SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Local to global never set with MatSetLocalToGlobalMapping()"); 760 } 761 if (nrow > 2048 || ncol > 2048) { 762 SETERRQ2(PETSC_ERR_SUP,0,"Number column/row indices must be <= 2048: are %d %d",nrow,ncol); 763 } 764 if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix"); 765 #endif 766 767 if (mat->assembled) { 768 mat->was_assembled = PETSC_TRUE; 769 mat->assembled = PETSC_FALSE; 770 } 771 PLogEventBegin(MAT_SetValues,mat,0,0,0); 772 ierr = ISLocalToGlobalMappingApply(mat->mapping,nrow,irow,irowm);CHKERRQ(ierr); 773 ierr = ISLocalToGlobalMappingApply(mat->mapping,ncol,icol,icolm);CHKERRQ(ierr); 774 ierr = (*mat->ops->setvalues)(mat,nrow,irowm,ncol,icolm,y,addv);CHKERRQ(ierr); 775 PLogEventEnd(MAT_SetValues,mat,0,0,0); 776 PetscFunctionReturn(0); 777 } 778 779 #undef __FUNC__ 780 #define __FUNC__ "MatSetValuesBlockedLocal" 781 /*@ 782 MatSetValuesBlockedLocal - Inserts or adds values into certain locations of a matrix, 783 using a local ordering of the nodes a block at a time. 784 785 Not Collective 786 787 Input Parameters: 788 + x - the matrix 789 . nrow, irow - number of rows and their local indices 790 . ncol, icol - number of columns and their local indices 791 . y - a logically two-dimensional array of values 792 - addv - either INSERT_VALUES or ADD_VALUES, where 793 ADD_VALUES adds values to any existing entries, and 794 INSERT_VALUES replaces existing entries with new values 795 796 Notes: 797 Before calling MatSetValuesBlockedLocal(), the user must first set the 798 local-to-global mapping by calling MatSetLocalToGlobalMappingBlock(), 799 where the mapping MUST be set for matrix blocks, not for matrix elements. 800 801 Calls to MatSetValuesBlockedLocal() with the INSERT_VALUES and ADD_VALUES 802 options cannot be mixed without intervening calls to the assembly 803 routines. 804 805 These values may be cached, so MatAssemblyBegin() and MatAssemblyEnd() 806 MUST be called after all calls to MatSetValuesBlockedLocal() have been completed. 807 808 Level: intermediate 809 810 .keywords: matrix, set, values, blocked, local 811 812 .seealso: MatAssemblyBegin(), MatAssemblyEnd(), MatSetValuesLocal(), MatSetLocalToGlobalMappingBlock(), MatSetValuesBlocked() 813 @*/ 814 int MatSetValuesBlockedLocal(Mat mat,int nrow,int *irow,int ncol,int *icol,Scalar *y,InsertMode addv) 815 { 816 int ierr,irowm[2048],icolm[2048]; 817 818 PetscFunctionBegin; 819 PetscValidHeaderSpecific(mat,MAT_COOKIE); 820 PetscValidIntPointer(irow); 821 PetscValidIntPointer(icol); 822 PetscValidScalarPointer(y); 823 if (mat->insertmode == NOT_SET_VALUES) { 824 mat->insertmode = addv; 825 } 826 #if defined(PETSC_USE_BOPT_g) 827 else if (mat->insertmode != addv) { 828 SETERRQ(PETSC_ERR_ARG_WRONGSTATE,1,"Cannot mix add values and insert values"); 829 } 830 if (!mat->bmapping) { 831 SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Local to global never set with MatSetLocalToGlobalMappingBlock()"); 832 } 833 if (nrow > 2048 || ncol > 2048) { 834 SETERRQ2(PETSC_ERR_SUP,0,"Number column/row indices must be <= 2048: are %d %d",nrow,ncol); 835 } 836 if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix"); 837 #endif 838 839 if (mat->assembled) { 840 mat->was_assembled = PETSC_TRUE; 841 mat->assembled = PETSC_FALSE; 842 } 843 PLogEventBegin(MAT_SetValues,mat,0,0,0); 844 ierr = ISLocalToGlobalMappingApply(mat->bmapping,nrow,irow,irowm);CHKERRQ(ierr); 845 ierr = ISLocalToGlobalMappingApply(mat->bmapping,ncol,icol,icolm);CHKERRQ(ierr); 846 ierr = (*mat->ops->setvaluesblocked)(mat,nrow,irowm,ncol,icolm,y,addv);CHKERRQ(ierr); 847 PLogEventEnd(MAT_SetValues,mat,0,0,0); 848 PetscFunctionReturn(0); 849 } 850 851 /* --------------------------------------------------------*/ 852 #undef __FUNC__ 853 #define __FUNC__ "MatMult" 854 /*@ 855 MatMult - Computes the matrix-vector product, y = Ax. 856 857 Collective on Mat and Vec 858 859 Input Parameters: 860 + mat - the matrix 861 - x - the vector to be multilplied 862 863 Output Parameters: 864 . y - the result 865 866 Notes: 867 The vectors x and y cannot be the same. I.e., one cannot 868 call MatMult(A,y,y). 869 870 Level: beginner 871 872 .keywords: matrix, multiply, matrix-vector product 873 874 .seealso: MatMultTrans(), MatMultAdd(), MatMultTransAdd() 875 @*/ 876 int MatMult(Mat mat,Vec x,Vec y) 877 { 878 int ierr; 879 880 PetscFunctionBegin; 881 PetscValidHeaderSpecific(mat,MAT_COOKIE); 882 PetscValidHeaderSpecific(x,VEC_COOKIE); 883 PetscValidHeaderSpecific(y,VEC_COOKIE); 884 PetscCheckSameComm(mat,x); 885 PetscCheckSameComm(mat,y); 886 if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix"); 887 if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix"); 888 if (x == y) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"x and y must be different vectors"); 889 if (mat->N != x->N) SETERRQ2(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec x: global dim %d %d",mat->N,x->N); 890 if (mat->M != y->N) SETERRQ2(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec y: global dim %d %d",mat->M,y->N); 891 if (mat->m != y->n) SETERRQ2(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec y: local dim %d %d",mat->m,y->n); 892 893 PLogEventBegin(MAT_Mult,mat,x,y,0); 894 ierr = (*mat->ops->mult)(mat,x,y);CHKERRQ(ierr); 895 PLogEventEnd(MAT_Mult,mat,x,y,0); 896 897 PetscFunctionReturn(0); 898 } 899 900 #undef __FUNC__ 901 #define __FUNC__ "MatMultTrans" 902 /*@ 903 MatMultTrans - Computes matrix transpose times a vector. 904 905 Collective on Mat and Vec 906 907 Input Parameters: 908 + mat - the matrix 909 - x - the vector to be multilplied 910 911 Output Parameters: 912 . y - the result 913 914 Notes: 915 The vectors x and y cannot be the same. I.e., one cannot 916 call MatMultTrans(A,y,y). 917 918 Level: beginner 919 920 .keywords: matrix, multiply, matrix-vector product, transpose 921 922 .seealso: MatMult(), MatMultAdd(), MatMultTransAdd() 923 @*/ 924 int MatMultTrans(Mat mat,Vec x,Vec y) 925 { 926 int ierr; 927 928 PetscFunctionBegin; 929 PetscValidHeaderSpecific(mat,MAT_COOKIE); 930 PetscValidHeaderSpecific(x,VEC_COOKIE); 931 PetscValidHeaderSpecific(y,VEC_COOKIE); 932 PetscCheckSameComm(mat,x); 933 PetscCheckSameComm(mat,y); 934 if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix"); 935 if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix"); 936 if (x == y) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"x and y must be different vectors"); 937 if (mat->M != x->N) SETERRQ2(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec x: global dim %d %d",mat->M,x->N); 938 if (mat->N != y->N) SETERRQ2(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec y: global dim %d %d",mat->N,y->N); 939 940 PLogEventBegin(MAT_MultTrans,mat,x,y,0); 941 ierr = (*mat->ops->multtrans)(mat,x,y);CHKERRQ(ierr); 942 PLogEventEnd(MAT_MultTrans,mat,x,y,0); 943 PetscFunctionReturn(0); 944 } 945 946 #undef __FUNC__ 947 #define __FUNC__ "MatMultAdd" 948 /*@ 949 MatMultAdd - Computes v3 = v2 + A * v1. 950 951 Collective on Mat and Vec 952 953 Input Parameters: 954 + mat - the matrix 955 - v1, v2 - the vectors 956 957 Output Parameters: 958 . v3 - the result 959 960 Notes: 961 The vectors v1 and v3 cannot be the same. I.e., one cannot 962 call MatMultAdd(A,v1,v2,v1). 963 964 Level: beginner 965 966 .keywords: matrix, multiply, matrix-vector product, add 967 968 .seealso: MatMultTrans(), MatMult(), MatMultTransAdd() 969 @*/ 970 int MatMultAdd(Mat mat,Vec v1,Vec v2,Vec v3) 971 { 972 int ierr; 973 974 PetscFunctionBegin; 975 PetscValidHeaderSpecific(mat,MAT_COOKIE); 976 PetscValidHeaderSpecific(v1,VEC_COOKIE); 977 PetscValidHeaderSpecific(v2,VEC_COOKIE); 978 PetscValidHeaderSpecific(v3,VEC_COOKIE); 979 PetscCheckSameComm(mat,v1); 980 PetscCheckSameComm(mat,v2); 981 PetscCheckSameComm(mat,v3); 982 if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix"); 983 if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix"); 984 if (mat->N != v1->N) SETERRQ2(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec v1: global dim %d %d",mat->N,v1->N); 985 if (mat->M != v2->N) SETERRQ2(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec v2: global dim %d %d",mat->M,v2->N); 986 if (mat->M != v3->N) SETERRQ2(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec v3: global dim %d %d",mat->M,v3->N); 987 if (mat->m != v3->n) SETERRQ2(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec v3: local dim %d %d",mat->m,v3->n); 988 if (mat->m != v2->n) SETERRQ2(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec v2: local dim %d %d",mat->m,v2->n); 989 if (v1 == v3) SETERRQ(PETSC_ERR_ARG_IDN,0,"v1 and v3 must be different vectors"); 990 991 PLogEventBegin(MAT_MultAdd,mat,v1,v2,v3); 992 ierr = (*mat->ops->multadd)(mat,v1,v2,v3);CHKERRQ(ierr); 993 PLogEventEnd(MAT_MultAdd,mat,v1,v2,v3); 994 PetscFunctionReturn(0); 995 } 996 997 #undef __FUNC__ 998 #define __FUNC__ "MatMultTransAdd" 999 /*@ 1000 MatMultTransAdd - Computes v3 = v2 + A' * v1. 1001 1002 Collective on Mat and Vec 1003 1004 Input Parameters: 1005 + mat - the matrix 1006 - v1, v2 - the vectors 1007 1008 Output Parameters: 1009 . v3 - the result 1010 1011 Notes: 1012 The vectors v1 and v3 cannot be the same. I.e., one cannot 1013 call MatMultTransAdd(A,v1,v2,v1). 1014 1015 Level: beginner 1016 1017 .keywords: matrix, multiply, matrix-vector product, transpose, add 1018 1019 .seealso: MatMultTrans(), MatMultAdd(), MatMult() 1020 @*/ 1021 int MatMultTransAdd(Mat mat,Vec v1,Vec v2,Vec v3) 1022 { 1023 int ierr; 1024 1025 PetscFunctionBegin; 1026 PetscValidHeaderSpecific(mat,MAT_COOKIE); 1027 PetscValidHeaderSpecific(v1,VEC_COOKIE); 1028 PetscValidHeaderSpecific(v2,VEC_COOKIE); 1029 PetscValidHeaderSpecific(v3,VEC_COOKIE); 1030 PetscCheckSameComm(mat,v1); 1031 PetscCheckSameComm(mat,v2); 1032 PetscCheckSameComm(mat,v3); 1033 if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix"); 1034 if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix"); 1035 if (!mat->ops->multtransadd) SETERRQ(PETSC_ERR_SUP,0,""); 1036 if (v1 == v3) SETERRQ(PETSC_ERR_ARG_IDN,0,"v1 and v3 must be different vectors"); 1037 if (mat->M != v1->N) SETERRQ2(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec v1: global dim %d %d",mat->M,v1->N); 1038 if (mat->N != v2->N) SETERRQ2(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec v2: global dim %d %d",mat->N,v2->N); 1039 if (mat->N != v3->N) SETERRQ2(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec v3: global dim %d %d",mat->N,v3->N); 1040 1041 PLogEventBegin(MAT_MultTransAdd,mat,v1,v2,v3); 1042 ierr = (*mat->ops->multtransadd)(mat,v1,v2,v3);CHKERRQ(ierr); 1043 PLogEventEnd(MAT_MultTransAdd,mat,v1,v2,v3); 1044 PetscFunctionReturn(0); 1045 } 1046 /* ------------------------------------------------------------*/ 1047 #undef __FUNC__ 1048 #define __FUNC__ "MatGetInfo" 1049 /*@C 1050 MatGetInfo - Returns information about matrix storage (number of 1051 nonzeros, memory, etc.). 1052 1053 Collective on Mat if MAT_GLOBAL_MAX or MAT_GLOBAL_SUM is used 1054 as the flag 1055 1056 Input Parameters: 1057 . mat - the matrix 1058 1059 Output Parameters: 1060 + flag - flag indicating the type of parameters to be returned 1061 (MAT_LOCAL - local matrix, MAT_GLOBAL_MAX - maximum over all processors, 1062 MAT_GLOBAL_SUM - sum over all processors) 1063 - info - matrix information context 1064 1065 Notes: 1066 The MatInfo context contains a variety of matrix data, including 1067 number of nonzeros allocated and used, number of mallocs during 1068 matrix assembly, etc. Additional information for factored matrices 1069 is provided (such as the fill ratio, number of mallocs during 1070 factorization, etc.). Much of this info is printed to STDOUT 1071 when using the runtime options 1072 $ -log_info -mat_view_info 1073 1074 Example for C/C++ Users: 1075 See the file ${PETSC_DIR}/include/mat.h for a complete list of 1076 data within the MatInfo context. For example, 1077 .vb 1078 MatInfo info; 1079 Mat A; 1080 double mal, nz_a, nz_u; 1081 1082 MatGetInfo(A,MAT_LOCAL,&info); 1083 mal = info.mallocs; 1084 nz_a = info.nz_allocated; 1085 .ve 1086 1087 Example for Fortran Users: 1088 Fortran users should declare info as a double precision 1089 array of dimension MAT_INFO_SIZE, and then extract the parameters 1090 of interest. See the file ${PETSC_DIR}/include/finclude/mat.h 1091 a complete list of parameter names. 1092 .vb 1093 double precision info(MAT_INFO_SIZE) 1094 double precision mal, nz_a 1095 Mat A 1096 integer ierr 1097 1098 call MatGetInfo(A,MAT_LOCAL,info,ierr) 1099 mal = info(MAT_INFO_MALLOCS) 1100 nz_a = info(MAT_INFO_NZ_ALLOCATED) 1101 .ve 1102 1103 Level: intermediate 1104 1105 .keywords: matrix, get, info, storage, nonzeros, memory, fill 1106 @*/ 1107 int MatGetInfo(Mat mat,MatInfoType flag,MatInfo *info) 1108 { 1109 int ierr; 1110 1111 PetscFunctionBegin; 1112 PetscValidHeaderSpecific(mat,MAT_COOKIE); 1113 PetscValidPointer(info); 1114 if (!mat->ops->getinfo) SETERRQ(PETSC_ERR_SUP,0,""); 1115 ierr = (*mat->ops->getinfo)(mat,flag,info);CHKERRQ(ierr); 1116 PetscFunctionReturn(0); 1117 } 1118 1119 /* ----------------------------------------------------------*/ 1120 #undef __FUNC__ 1121 #define __FUNC__ "MatILUDTFactor" 1122 /*@ 1123 MatILUDTFactor - Performs a drop tolerance ILU factorization. 1124 1125 Collective on Mat 1126 1127 Input Parameters: 1128 + mat - the matrix 1129 . dt - the drop tolerance 1130 . maxnz - the maximum number of nonzeros per row allowed 1131 . row - row permutation 1132 - col - column permutation 1133 1134 Output Parameters: 1135 . fact - the factored matrix 1136 1137 Notes: 1138 Most users should employ the simplified SLES interface for linear solvers 1139 instead of working directly with matrix algebra routines such as this. 1140 See, e.g., SLESCreate(). 1141 1142 Level: developer 1143 1144 .keywords: matrix, factor, LU, in-place 1145 1146 .seealso: MatLUFactorSymbolic(), MatLUFactorNumeric(), MatCholeskyFactor() 1147 @*/ 1148 int MatILUDTFactor(Mat mat,double dt,int maxnz,IS row,IS col,Mat *fact) 1149 { 1150 int ierr; 1151 1152 PetscFunctionBegin; 1153 PetscValidHeaderSpecific(mat,MAT_COOKIE); 1154 PetscValidPointer(fact); 1155 if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix"); 1156 if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix"); 1157 if (!mat->ops->iludtfactor) SETERRQ(PETSC_ERR_SUP,0,""); 1158 1159 PLogEventBegin(MAT_ILUFactor,mat,row,col,0); 1160 ierr = (*mat->ops->iludtfactor)(mat,dt,maxnz,row,col,fact);CHKERRQ(ierr); 1161 PLogEventEnd(MAT_ILUFactor,mat,row,col,0); 1162 1163 PetscFunctionReturn(0); 1164 } 1165 1166 #undef __FUNC__ 1167 #define __FUNC__ "MatLUFactor" 1168 /*@ 1169 MatLUFactor - Performs in-place LU factorization of matrix. 1170 1171 Collective on Mat 1172 1173 Input Parameters: 1174 + mat - the matrix 1175 . row - row permutation 1176 . col - column permutation 1177 - f - expected fill as ratio of original fill. 1178 1179 Notes: 1180 Most users should employ the simplified SLES interface for linear solvers 1181 instead of working directly with matrix algebra routines such as this. 1182 See, e.g., SLESCreate(). 1183 1184 Level: developer 1185 1186 .keywords: matrix, factor, LU, in-place 1187 1188 .seealso: MatLUFactorSymbolic(), MatLUFactorNumeric(), MatCholeskyFactor(), 1189 MatGetOrdering() 1190 1191 @*/ 1192 int MatLUFactor(Mat mat,IS row,IS col,double f) 1193 { 1194 int ierr; 1195 1196 PetscFunctionBegin; 1197 PetscValidHeaderSpecific(mat,MAT_COOKIE); 1198 if (mat->M != mat->N) SETERRQ(PETSC_ERR_ARG_WRONG,0,"matrix must be square"); 1199 if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix"); 1200 if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix"); 1201 if (!mat->ops->lufactor) SETERRQ(PETSC_ERR_SUP,0,""); 1202 1203 PLogEventBegin(MAT_LUFactor,mat,row,col,0); 1204 ierr = (*mat->ops->lufactor)(mat,row,col,f);CHKERRQ(ierr); 1205 PLogEventEnd(MAT_LUFactor,mat,row,col,0); 1206 PetscFunctionReturn(0); 1207 } 1208 1209 #undef __FUNC__ 1210 #define __FUNC__ "MatILUFactor" 1211 /*@ 1212 MatILUFactor - Performs in-place ILU factorization of matrix. 1213 1214 Collective on Mat 1215 1216 Input Parameters: 1217 + mat - the matrix 1218 . row - row permutation 1219 . col - column permutation 1220 - info - structure containing 1221 $ levels - number of levels of fill. 1222 $ expected fill - as ratio of original fill. 1223 $ 1 or 0 - indicating force fill on diagonal (improves robustness for matrices 1224 missing diagonal entries) 1225 1226 Notes: 1227 Probably really in-place only when level of fill is zero, otherwise allocates 1228 new space to store factored matrix and deletes previous memory. 1229 1230 Most users should employ the simplified SLES interface for linear solvers 1231 instead of working directly with matrix algebra routines such as this. 1232 See, e.g., SLESCreate(). 1233 1234 Level: developer 1235 1236 .keywords: matrix, factor, ILU, in-place 1237 1238 .seealso: MatILUFactorSymbolic(), MatLUFactorNumeric(), MatCholeskyFactor() 1239 @*/ 1240 int MatILUFactor(Mat mat,IS row,IS col,MatILUInfo *info) 1241 { 1242 int ierr; 1243 1244 PetscFunctionBegin; 1245 PetscValidHeaderSpecific(mat,MAT_COOKIE); 1246 if (mat->M != mat->N) SETERRQ(PETSC_ERR_ARG_WRONG,0,"matrix must be square"); 1247 if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix"); 1248 if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix"); 1249 if (!mat->ops->ilufactor) SETERRQ(PETSC_ERR_SUP,0,""); 1250 1251 PLogEventBegin(MAT_ILUFactor,mat,row,col,0); 1252 ierr = (*mat->ops->ilufactor)(mat,row,col,info);CHKERRQ(ierr); 1253 PLogEventEnd(MAT_ILUFactor,mat,row,col,0); 1254 PetscFunctionReturn(0); 1255 } 1256 1257 #undef __FUNC__ 1258 #define __FUNC__ "MatLUFactorSymbolic" 1259 /*@ 1260 MatLUFactorSymbolic - Performs symbolic LU factorization of matrix. 1261 Call this routine before calling MatLUFactorNumeric(). 1262 1263 Collective on Mat 1264 1265 Input Parameters: 1266 + mat - the matrix 1267 . row, col - row and column permutations 1268 - f - expected fill as ratio of the original number of nonzeros, 1269 for example 3.0; choosing this parameter well can result in 1270 more efficient use of time and space. Run with the option -log_info 1271 to determine an optimal value to use 1272 1273 Output Parameter: 1274 . fact - new matrix that has been symbolically factored 1275 1276 Notes: 1277 See the users manual for additional information about 1278 choosing the fill factor for better efficiency. 1279 1280 Most users should employ the simplified SLES interface for linear solvers 1281 instead of working directly with matrix algebra routines such as this. 1282 See, e.g., SLESCreate(). 1283 1284 Level: developer 1285 1286 .keywords: matrix, factor, LU, symbolic, fill 1287 1288 .seealso: MatLUFactor(), MatLUFactorNumeric(), MatCholeskyFactor() 1289 @*/ 1290 int MatLUFactorSymbolic(Mat mat,IS row,IS col,double f,Mat *fact) 1291 { 1292 int ierr; 1293 1294 PetscFunctionBegin; 1295 PetscValidHeaderSpecific(mat,MAT_COOKIE); 1296 if (mat->M != mat->N) SETERRQ(PETSC_ERR_ARG_WRONG,0,"matrix must be square"); 1297 PetscValidPointer(fact); 1298 if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix"); 1299 if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix"); 1300 if (!mat->ops->lufactorsymbolic) SETERRQ(PETSC_ERR_SUP,0,""); 1301 1302 PLogEventBegin(MAT_LUFactorSymbolic,mat,row,col,0); 1303 ierr = (*mat->ops->lufactorsymbolic)(mat,row,col,f,fact);CHKERRQ(ierr); 1304 PLogEventEnd(MAT_LUFactorSymbolic,mat,row,col,0); 1305 PetscFunctionReturn(0); 1306 } 1307 1308 #undef __FUNC__ 1309 #define __FUNC__ "MatLUFactorNumeric" 1310 /*@ 1311 MatLUFactorNumeric - Performs numeric LU factorization of a matrix. 1312 Call this routine after first calling MatLUFactorSymbolic(). 1313 1314 Collective on Mat 1315 1316 Input Parameters: 1317 + mat - the matrix 1318 - row, col - row and column permutations 1319 1320 Output Parameters: 1321 . fact - symbolically factored matrix that must have been generated 1322 by MatLUFactorSymbolic() 1323 1324 Notes: 1325 See MatLUFactor() for in-place factorization. See 1326 MatCholeskyFactorNumeric() for the symmetric, positive definite case. 1327 1328 Most users should employ the simplified SLES interface for linear solvers 1329 instead of working directly with matrix algebra routines such as this. 1330 See, e.g., SLESCreate(). 1331 1332 Level: developer 1333 1334 .keywords: matrix, factor, LU, numeric 1335 1336 .seealso: MatLUFactorSymbolic(), MatLUFactor(), MatCholeskyFactor() 1337 @*/ 1338 int MatLUFactorNumeric(Mat mat,Mat *fact) 1339 { 1340 int ierr,flg; 1341 1342 PetscFunctionBegin; 1343 PetscValidHeaderSpecific(mat,MAT_COOKIE); 1344 PetscValidPointer(fact); 1345 PetscValidHeaderSpecific(*fact,MAT_COOKIE); 1346 if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix"); 1347 if (mat->M != (*fact)->M || mat->N != (*fact)->N) { 1348 SETERRQ4(PETSC_ERR_ARG_SIZ,0,"Mat mat,Mat *fact: global dimensions are different %d should = %d %d should = %d", 1349 mat->M,(*fact)->M,mat->N,(*fact)->N); 1350 } 1351 if (!(*fact)->ops->lufactornumeric) SETERRQ(PETSC_ERR_SUP,0,""); 1352 1353 PLogEventBegin(MAT_LUFactorNumeric,mat,*fact,0,0); 1354 ierr = (*(*fact)->ops->lufactornumeric)(mat,fact);CHKERRQ(ierr); 1355 PLogEventEnd(MAT_LUFactorNumeric,mat,*fact,0,0); 1356 ierr = OptionsHasName(PETSC_NULL,"-mat_view_draw",&flg);CHKERRQ(ierr); 1357 if (flg) { 1358 ierr = OptionsHasName(0,"-mat_view_contour",&flg);CHKERRQ(ierr); 1359 if (flg) { 1360 ViewerPushFormat(VIEWER_DRAW_(mat->comm),VIEWER_FORMAT_DRAW_CONTOUR,0);CHKERRQ(ierr); 1361 } 1362 ierr = MatView(*fact,VIEWER_DRAW_(mat->comm));CHKERRQ(ierr); 1363 ierr = ViewerFlush(VIEWER_DRAW_(mat->comm));CHKERRQ(ierr); 1364 if (flg) { 1365 ViewerPopFormat(VIEWER_DRAW_(mat->comm));CHKERRQ(ierr); 1366 } 1367 } 1368 PetscFunctionReturn(0); 1369 } 1370 1371 #undef __FUNC__ 1372 #define __FUNC__ "MatCholeskyFactor" 1373 /*@ 1374 MatCholeskyFactor - Performs in-place Cholesky factorization of a 1375 symmetric matrix. 1376 1377 Collective on Mat 1378 1379 Input Parameters: 1380 + mat - the matrix 1381 . perm - row and column permutations 1382 - f - expected fill as ratio of original fill 1383 1384 Notes: 1385 See MatLUFactor() for the nonsymmetric case. See also 1386 MatCholeskyFactorSymbolic(), and MatCholeskyFactorNumeric(). 1387 1388 Most users should employ the simplified SLES interface for linear solvers 1389 instead of working directly with matrix algebra routines such as this. 1390 See, e.g., SLESCreate(). 1391 1392 Level: developer 1393 1394 .keywords: matrix, factor, in-place, Cholesky 1395 1396 .seealso: MatLUFactor(), MatCholeskyFactorSymbolic(), MatCholeskyFactorNumeric() 1397 MatGetOrdering() 1398 1399 @*/ 1400 int MatCholeskyFactor(Mat mat,IS perm,double f) 1401 { 1402 int ierr; 1403 1404 PetscFunctionBegin; 1405 PetscValidHeaderSpecific(mat,MAT_COOKIE); 1406 if (mat->M != mat->N) SETERRQ(PETSC_ERR_ARG_WRONG,0,"Matrix must be square"); 1407 if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix"); 1408 if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix"); 1409 if (!mat->ops->choleskyfactor) SETERRQ(PETSC_ERR_SUP,0,""); 1410 1411 PLogEventBegin(MAT_CholeskyFactor,mat,perm,0,0); 1412 ierr = (*mat->ops->choleskyfactor)(mat,perm,f);CHKERRQ(ierr); 1413 PLogEventEnd(MAT_CholeskyFactor,mat,perm,0,0); 1414 PetscFunctionReturn(0); 1415 } 1416 1417 #undef __FUNC__ 1418 #define __FUNC__ "MatCholeskyFactorSymbolic" 1419 /*@ 1420 MatCholeskyFactorSymbolic - Performs symbolic Cholesky factorization 1421 of a symmetric matrix. 1422 1423 Collective on Mat 1424 1425 Input Parameters: 1426 + mat - the matrix 1427 . perm - row and column permutations 1428 - f - expected fill as ratio of original 1429 1430 Output Parameter: 1431 . fact - the factored matrix 1432 1433 Notes: 1434 See MatLUFactorSymbolic() for the nonsymmetric case. See also 1435 MatCholeskyFactor() and MatCholeskyFactorNumeric(). 1436 1437 Most users should employ the simplified SLES interface for linear solvers 1438 instead of working directly with matrix algebra routines such as this. 1439 See, e.g., SLESCreate(). 1440 1441 Level: developer 1442 1443 .keywords: matrix, factor, factorization, symbolic, Cholesky 1444 1445 .seealso: MatLUFactorSymbolic(), MatCholeskyFactor(), MatCholeskyFactorNumeric() 1446 MatGetOrdering() 1447 1448 @*/ 1449 int MatCholeskyFactorSymbolic(Mat mat,IS perm,double f,Mat *fact) 1450 { 1451 int ierr; 1452 1453 PetscFunctionBegin; 1454 PetscValidHeaderSpecific(mat,MAT_COOKIE); 1455 PetscValidPointer(fact); 1456 if (mat->M != mat->N) SETERRQ(PETSC_ERR_ARG_WRONG,0,"Matrix must be square"); 1457 if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix"); 1458 if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix"); 1459 if (!mat->ops->choleskyfactorsymbolic) SETERRQ(PETSC_ERR_SUP,0,""); 1460 1461 PLogEventBegin(MAT_CholeskyFactorSymbolic,mat,perm,0,0); 1462 ierr = (*mat->ops->choleskyfactorsymbolic)(mat,perm,f,fact);CHKERRQ(ierr); 1463 PLogEventEnd(MAT_CholeskyFactorSymbolic,mat,perm,0,0); 1464 PetscFunctionReturn(0); 1465 } 1466 1467 #undef __FUNC__ 1468 #define __FUNC__ "MatCholeskyFactorNumeric" 1469 /*@ 1470 MatCholeskyFactorNumeric - Performs numeric Cholesky factorization 1471 of a symmetric matrix. Call this routine after first calling 1472 MatCholeskyFactorSymbolic(). 1473 1474 Collective on Mat 1475 1476 Input Parameter: 1477 . mat - the initial matrix 1478 1479 Output Parameter: 1480 . fact - the factored matrix 1481 1482 Notes: 1483 Most users should employ the simplified SLES interface for linear solvers 1484 instead of working directly with matrix algebra routines such as this. 1485 See, e.g., SLESCreate(). 1486 1487 Level: developer 1488 1489 .keywords: matrix, factor, numeric, Cholesky 1490 1491 .seealso: MatCholeskyFactorSymbolic(), MatCholeskyFactor(), MatLUFactorNumeric() 1492 @*/ 1493 int MatCholeskyFactorNumeric(Mat mat,Mat *fact) 1494 { 1495 int ierr; 1496 1497 PetscFunctionBegin; 1498 PetscValidHeaderSpecific(mat,MAT_COOKIE); 1499 PetscValidPointer(fact); 1500 if (!mat->ops->choleskyfactornumeric) SETERRQ(PETSC_ERR_SUP,0,""); 1501 if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix"); 1502 if (mat->M != (*fact)->M || mat->N != (*fact)->N) { 1503 SETERRQ4(PETSC_ERR_ARG_SIZ,0,"Mat mat,Mat *fact: global dim %d should = %d %d should = %d", 1504 mat->M,(*fact)->M,mat->N,(*fact)->N); 1505 } 1506 1507 PLogEventBegin(MAT_CholeskyFactorNumeric,mat,*fact,0,0); 1508 ierr = (*mat->ops->choleskyfactornumeric)(mat,fact);CHKERRQ(ierr); 1509 PLogEventEnd(MAT_CholeskyFactorNumeric,mat,*fact,0,0); 1510 PetscFunctionReturn(0); 1511 } 1512 1513 /* ----------------------------------------------------------------*/ 1514 #undef __FUNC__ 1515 #define __FUNC__ "MatSolve" 1516 /*@ 1517 MatSolve - Solves A x = b, given a factored matrix. 1518 1519 Collective on Mat and Vec 1520 1521 Input Parameters: 1522 + mat - the factored matrix 1523 - b - the right-hand-side vector 1524 1525 Output Parameter: 1526 . x - the result vector 1527 1528 Notes: 1529 The vectors b and x cannot be the same. I.e., one cannot 1530 call MatSolve(A,x,x). 1531 1532 Notes: 1533 Most users should employ the simplified SLES interface for linear solvers 1534 instead of working directly with matrix algebra routines such as this. 1535 See, e.g., SLESCreate(). 1536 1537 Level: developer 1538 1539 .keywords: matrix, linear system, solve, LU, Cholesky, triangular solve 1540 1541 .seealso: MatSolveAdd(), MatSolveTrans(), MatSolveTransAdd() 1542 @*/ 1543 int MatSolve(Mat mat,Vec b,Vec x) 1544 { 1545 int ierr; 1546 1547 PetscFunctionBegin; 1548 PetscValidHeaderSpecific(mat,MAT_COOKIE); 1549 PetscValidHeaderSpecific(b,VEC_COOKIE); 1550 PetscValidHeaderSpecific(x,VEC_COOKIE); 1551 PetscCheckSameComm(mat,b); 1552 PetscCheckSameComm(mat,x); 1553 if (x == b) SETERRQ(PETSC_ERR_ARG_IDN,0,"x and b must be different vectors"); 1554 if (!mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Unfactored matrix"); 1555 if (mat->N != x->N) SETERRQ2(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec x: global dim %d %d",mat->N,x->N); 1556 if (mat->M != b->N) SETERRQ2(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec b: global dim %d %d",mat->M,b->N); 1557 if (mat->m != b->n) SETERRQ2(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec b: local dim %d %d",mat->m,b->n); 1558 if (mat->M == 0 && mat->N == 0) PetscFunctionReturn(0); 1559 1560 if (!mat->ops->solve) SETERRQ(PETSC_ERR_SUP,0,""); 1561 PLogEventBegin(MAT_Solve,mat,b,x,0); 1562 ierr = (*mat->ops->solve)(mat,b,x);CHKERRQ(ierr); 1563 PLogEventEnd(MAT_Solve,mat,b,x,0); 1564 PetscFunctionReturn(0); 1565 } 1566 1567 #undef __FUNC__ 1568 #define __FUNC__ "MatForwardSolve" 1569 /* @ 1570 MatForwardSolve - Solves L x = b, given a factored matrix, A = LU. 1571 1572 Collective on Mat and Vec 1573 1574 Input Parameters: 1575 + mat - the factored matrix 1576 - b - the right-hand-side vector 1577 1578 Output Parameter: 1579 . x - the result vector 1580 1581 Notes: 1582 MatSolve() should be used for most applications, as it performs 1583 a forward solve followed by a backward solve. 1584 1585 The vectors b and x cannot be the same. I.e., one cannot 1586 call MatForwardSolve(A,x,x). 1587 1588 Most users should employ the simplified SLES interface for linear solvers 1589 instead of working directly with matrix algebra routines such as this. 1590 See, e.g., SLESCreate(). 1591 1592 Level: developer 1593 1594 .keywords: matrix, forward, LU, Cholesky, triangular solve 1595 1596 .seealso: MatSolve(), MatBackwardSolve() 1597 @ */ 1598 int MatForwardSolve(Mat mat,Vec b,Vec x) 1599 { 1600 int ierr; 1601 1602 PetscFunctionBegin; 1603 PetscValidHeaderSpecific(mat,MAT_COOKIE); 1604 PetscValidHeaderSpecific(b,VEC_COOKIE); 1605 PetscValidHeaderSpecific(x,VEC_COOKIE); 1606 PetscCheckSameComm(mat,b); 1607 PetscCheckSameComm(mat,x); 1608 if (x == b) SETERRQ(PETSC_ERR_ARG_IDN,0,"x and b must be different vectors"); 1609 if (!mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Unfactored matrix"); 1610 if (!mat->ops->forwardsolve) SETERRQ(PETSC_ERR_SUP,0,""); 1611 if (mat->N != x->N) SETERRQ2(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec x: global dim %d %d",mat->N,x->N); 1612 if (mat->M != b->N) SETERRQ2(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec b: global dim %d %d",mat->M,b->N); 1613 if (mat->m != b->n) SETERRQ2(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec b: local dim %d %d",mat->m,b->n); 1614 1615 PLogEventBegin(MAT_ForwardSolve,mat,b,x,0); 1616 ierr = (*mat->ops->forwardsolve)(mat,b,x);CHKERRQ(ierr); 1617 PLogEventEnd(MAT_ForwardSolve,mat,b,x,0); 1618 PetscFunctionReturn(0); 1619 } 1620 1621 #undef __FUNC__ 1622 #define __FUNC__ "MatBackwardSolve" 1623 /* @ 1624 MatBackwardSolve - Solves U x = b, given a factored matrix, A = LU. 1625 1626 Collective on Mat and Vec 1627 1628 Input Parameters: 1629 + mat - the factored matrix 1630 - b - the right-hand-side vector 1631 1632 Output Parameter: 1633 . x - the result vector 1634 1635 Notes: 1636 MatSolve() should be used for most applications, as it performs 1637 a forward solve followed by a backward solve. 1638 1639 The vectors b and x cannot be the same. I.e., one cannot 1640 call MatBackwardSolve(A,x,x). 1641 1642 Most users should employ the simplified SLES interface for linear solvers 1643 instead of working directly with matrix algebra routines such as this. 1644 See, e.g., SLESCreate(). 1645 1646 Level: developer 1647 1648 .keywords: matrix, backward, LU, Cholesky, triangular solve 1649 1650 .seealso: MatSolve(), MatForwardSolve() 1651 @ */ 1652 int MatBackwardSolve(Mat mat,Vec b,Vec x) 1653 { 1654 int ierr; 1655 1656 PetscFunctionBegin; 1657 PetscValidHeaderSpecific(mat,MAT_COOKIE); 1658 PetscValidHeaderSpecific(b,VEC_COOKIE); 1659 PetscValidHeaderSpecific(x,VEC_COOKIE); 1660 PetscCheckSameComm(mat,b); 1661 PetscCheckSameComm(mat,x); 1662 if (x == b) SETERRQ(PETSC_ERR_ARG_IDN,0,"x and b must be different vectors"); 1663 if (!mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Unfactored matrix"); 1664 if (!mat->ops->backwardsolve) SETERRQ(PETSC_ERR_SUP,0,""); 1665 if (mat->N != x->N) SETERRQ2(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec x: global dim %d %d",mat->N,x->N); 1666 if (mat->M != b->N) SETERRQ2(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec b: global dim %d %d",mat->M,b->N); 1667 if (mat->m != b->n) SETERRQ2(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec b: local dim %d %d",mat->m,b->n); 1668 1669 PLogEventBegin(MAT_BackwardSolve,mat,b,x,0); 1670 ierr = (*mat->ops->backwardsolve)(mat,b,x);CHKERRQ(ierr); 1671 PLogEventEnd(MAT_BackwardSolve,mat,b,x,0); 1672 PetscFunctionReturn(0); 1673 } 1674 1675 #undef __FUNC__ 1676 #define __FUNC__ "MatSolveAdd" 1677 /*@ 1678 MatSolveAdd - Computes x = y + inv(A)*b, given a factored matrix. 1679 1680 Collective on Mat and Vec 1681 1682 Input Parameters: 1683 + mat - the factored matrix 1684 . b - the right-hand-side vector 1685 - y - the vector to be added to 1686 1687 Output Parameter: 1688 . x - the result vector 1689 1690 Notes: 1691 The vectors b and x cannot be the same. I.e., one cannot 1692 call MatSolveAdd(A,x,y,x). 1693 1694 Most users should employ the simplified SLES interface for linear solvers 1695 instead of working directly with matrix algebra routines such as this. 1696 See, e.g., SLESCreate(). 1697 1698 Level: developer 1699 1700 .keywords: matrix, linear system, solve, LU, Cholesky, add 1701 1702 .seealso: MatSolve(), MatSolveTrans(), MatSolveTransAdd() 1703 @*/ 1704 int MatSolveAdd(Mat mat,Vec b,Vec y,Vec x) 1705 { 1706 Scalar one = 1.0; 1707 Vec tmp; 1708 int ierr; 1709 1710 PetscFunctionBegin; 1711 PetscValidHeaderSpecific(mat,MAT_COOKIE); 1712 PetscValidHeaderSpecific(y,VEC_COOKIE); 1713 PetscValidHeaderSpecific(b,VEC_COOKIE); 1714 PetscValidHeaderSpecific(x,VEC_COOKIE); 1715 PetscCheckSameComm(mat,b); 1716 PetscCheckSameComm(mat,y); 1717 PetscCheckSameComm(mat,x); 1718 if (x == b) SETERRQ(PETSC_ERR_ARG_IDN,0,"x and b must be different vectors"); 1719 if (!mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Unfactored matrix"); 1720 if (mat->N != x->N) SETERRQ2(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec x: global dim %d %d",mat->N,x->N); 1721 if (mat->M != b->N) SETERRQ2(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec b: global dim %d %d",mat->M,b->N); 1722 if (mat->M != y->N) SETERRQ2(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec y: global dim %d %d",mat->M,y->N); 1723 if (mat->m != b->n) SETERRQ2(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec b: local dim %d %d",mat->m,b->n); 1724 if (x->n != y->n) SETERRQ2(PETSC_ERR_ARG_SIZ,0,"Vec x,Vec y: local dim %d %d",x->n,y->n); 1725 1726 PLogEventBegin(MAT_SolveAdd,mat,b,x,y); 1727 if (mat->ops->solveadd) { 1728 ierr = (*mat->ops->solveadd)(mat,b,y,x);CHKERRQ(ierr); 1729 } else { 1730 /* do the solve then the add manually */ 1731 if (x != y) { 1732 ierr = MatSolve(mat,b,x);CHKERRQ(ierr); 1733 ierr = VecAXPY(&one,y,x);CHKERRQ(ierr); 1734 } else { 1735 ierr = VecDuplicate(x,&tmp);CHKERRQ(ierr); 1736 PLogObjectParent(mat,tmp); 1737 ierr = VecCopy(x,tmp);CHKERRQ(ierr); 1738 ierr = MatSolve(mat,b,x);CHKERRQ(ierr); 1739 ierr = VecAXPY(&one,tmp,x);CHKERRQ(ierr); 1740 ierr = VecDestroy(tmp);CHKERRQ(ierr); 1741 } 1742 } 1743 PLogEventEnd(MAT_SolveAdd,mat,b,x,y); 1744 PetscFunctionReturn(0); 1745 } 1746 1747 #undef __FUNC__ 1748 #define __FUNC__ "MatSolveTrans" 1749 /*@ 1750 MatSolveTrans - Solves A' x = b, given a factored matrix. 1751 1752 Collective on Mat and Vec 1753 1754 Input Parameters: 1755 + mat - the factored matrix 1756 - b - the right-hand-side vector 1757 1758 Output Parameter: 1759 . x - the result vector 1760 1761 Notes: 1762 The vectors b and x cannot be the same. I.e., one cannot 1763 call MatSolveTrans(A,x,x). 1764 1765 Most users should employ the simplified SLES interface for linear solvers 1766 instead of working directly with matrix algebra routines such as this. 1767 See, e.g., SLESCreate(). 1768 1769 Level: developer 1770 1771 .keywords: matrix, linear system, solve, LU, Cholesky, transpose 1772 1773 .seealso: MatSolve(), MatSolveAdd(), MatSolveTransAdd() 1774 @*/ 1775 int MatSolveTrans(Mat mat,Vec b,Vec x) 1776 { 1777 int ierr; 1778 1779 PetscFunctionBegin; 1780 PetscValidHeaderSpecific(mat,MAT_COOKIE); 1781 PetscValidHeaderSpecific(b,VEC_COOKIE); 1782 PetscValidHeaderSpecific(x,VEC_COOKIE); 1783 PetscCheckSameComm(mat,b); 1784 PetscCheckSameComm(mat,x); 1785 if (!mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Unfactored matrix"); 1786 if (x == b) SETERRQ(PETSC_ERR_ARG_IDN,0,"x and b must be different vectors"); 1787 if (!mat->ops->solvetrans) SETERRQ(PETSC_ERR_SUP,0,""); 1788 if (mat->M != x->N) SETERRQ2(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec x: global dim %d %d",mat->M,x->N); 1789 if (mat->N != b->N) SETERRQ2(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec b: global dim %d %d",mat->N,b->N); 1790 1791 PLogEventBegin(MAT_SolveTrans,mat,b,x,0); 1792 ierr = (*mat->ops->solvetrans)(mat,b,x);CHKERRQ(ierr); 1793 PLogEventEnd(MAT_SolveTrans,mat,b,x,0); 1794 PetscFunctionReturn(0); 1795 } 1796 1797 #undef __FUNC__ 1798 #define __FUNC__ "MatSolveTransAdd" 1799 /*@ 1800 MatSolveTransAdd - Computes x = y + inv(trans(A)) b, given a 1801 factored matrix. 1802 1803 Collective on Mat and Vec 1804 1805 Input Parameters: 1806 + mat - the factored matrix 1807 . b - the right-hand-side vector 1808 - y - the vector to be added to 1809 1810 Output Parameter: 1811 . x - the result vector 1812 1813 Notes: 1814 The vectors b and x cannot be the same. I.e., one cannot 1815 call MatSolveTransAdd(A,x,y,x). 1816 1817 Most users should employ the simplified SLES interface for linear solvers 1818 instead of working directly with matrix algebra routines such as this. 1819 See, e.g., SLESCreate(). 1820 1821 Level: developer 1822 1823 .keywords: matrix, linear system, solve, LU, Cholesky, transpose, add 1824 1825 .seealso: MatSolve(), MatSolveAdd(), MatSolveTrans() 1826 @*/ 1827 int MatSolveTransAdd(Mat mat,Vec b,Vec y,Vec x) 1828 { 1829 Scalar one = 1.0; 1830 int ierr; 1831 Vec tmp; 1832 1833 PetscFunctionBegin; 1834 PetscValidHeaderSpecific(mat,MAT_COOKIE); 1835 PetscValidHeaderSpecific(y,VEC_COOKIE); 1836 PetscValidHeaderSpecific(b,VEC_COOKIE); 1837 PetscValidHeaderSpecific(x,VEC_COOKIE); 1838 PetscCheckSameComm(mat,b); 1839 PetscCheckSameComm(mat,y); 1840 PetscCheckSameComm(mat,x); 1841 if (x == b) SETERRQ(PETSC_ERR_ARG_IDN,0,"x and b must be different vectors"); 1842 if (!mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Unfactored matrix"); 1843 if (mat->M != x->N) SETERRQ2(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec x: global dim %d %d",mat->M,x->N); 1844 if (mat->N != b->N) SETERRQ2(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec b: global dim %d %d",mat->N,b->N); 1845 if (mat->N != y->N) SETERRQ2(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec y: global dim %d %d",mat->N,y->N); 1846 if (x->n != y->n) SETERRQ2(PETSC_ERR_ARG_SIZ,0,"Vec x,Vec y: local dim %d %d",x->n,y->n); 1847 1848 PLogEventBegin(MAT_SolveTransAdd,mat,b,x,y); 1849 if (mat->ops->solvetransadd) { 1850 ierr = (*mat->ops->solvetransadd)(mat,b,y,x);CHKERRQ(ierr); 1851 } else { 1852 /* do the solve then the add manually */ 1853 if (x != y) { 1854 ierr = MatSolveTrans(mat,b,x);CHKERRQ(ierr); 1855 ierr = VecAXPY(&one,y,x);CHKERRQ(ierr); 1856 } else { 1857 ierr = VecDuplicate(x,&tmp);CHKERRQ(ierr); 1858 PLogObjectParent(mat,tmp); 1859 ierr = VecCopy(x,tmp);CHKERRQ(ierr); 1860 ierr = MatSolveTrans(mat,b,x);CHKERRQ(ierr); 1861 ierr = VecAXPY(&one,tmp,x);CHKERRQ(ierr); 1862 ierr = VecDestroy(tmp);CHKERRQ(ierr); 1863 } 1864 } 1865 PLogEventEnd(MAT_SolveTransAdd,mat,b,x,y); 1866 PetscFunctionReturn(0); 1867 } 1868 /* ----------------------------------------------------------------*/ 1869 1870 #undef __FUNC__ 1871 #define __FUNC__ "MatRelax" 1872 /*@ 1873 MatRelax - Computes one relaxation sweep. 1874 1875 Collective on Mat and Vec 1876 1877 Input Parameters: 1878 + mat - the matrix 1879 . b - the right hand side 1880 . omega - the relaxation factor 1881 . flag - flag indicating the type of SOR (see below) 1882 . shift - diagonal shift 1883 - its - the number of iterations 1884 1885 Output Parameters: 1886 . x - the solution (can contain an initial guess) 1887 1888 SOR Flags: 1889 . SOR_FORWARD_SWEEP - forward SOR 1890 . SOR_BACKWARD_SWEEP - backward SOR 1891 . SOR_SYMMETRIC_SWEEP - SSOR (symmetric SOR) 1892 . SOR_LOCAL_FORWARD_SWEEP - local forward SOR 1893 . SOR_LOCAL_BACKWARD_SWEEP - local forward SOR 1894 . SOR_LOCAL_SYMMETRIC_SWEEP - local SSOR 1895 . SOR_APPLY_UPPER, SOR_APPLY_LOWER - applies 1896 upper/lower triangular part of matrix to 1897 vector (with omega) 1898 . SOR_ZERO_INITIAL_GUESS - zero initial guess 1899 1900 Notes: 1901 SOR_LOCAL_FORWARD_SWEEP, SOR_LOCAL_BACKWARD_SWEEP, and 1902 SOR_LOCAL_SYMMETRIC_SWEEP perform seperate independent smoothings 1903 on each processor. 1904 1905 Application programmers will not generally use MatRelax() directly, 1906 but instead will employ the SLES/PC interface. 1907 1908 Notes for Advanced Users: 1909 The flags are implemented as bitwise inclusive or operations. 1910 For example, use (SOR_ZERO_INITIAL_GUESS | SOR_SYMMETRIC_SWEEP) 1911 to specify a zero initial guess for SSOR. 1912 1913 Most users should employ the simplified SLES interface for linear solvers 1914 instead of working directly with matrix algebra routines such as this. 1915 See, e.g., SLESCreate(). 1916 1917 Level: developer 1918 1919 .keywords: matrix, relax, relaxation, sweep 1920 @*/ 1921 int MatRelax(Mat mat,Vec b,double omega,MatSORType flag,double shift,int its,Vec x) 1922 { 1923 int ierr; 1924 1925 PetscFunctionBegin; 1926 PetscValidHeaderSpecific(mat,MAT_COOKIE); 1927 PetscValidHeaderSpecific(b,VEC_COOKIE); 1928 PetscValidHeaderSpecific(x,VEC_COOKIE); 1929 PetscCheckSameComm(mat,b); 1930 PetscCheckSameComm(mat,x); 1931 if (!mat->ops->relax) SETERRQ(PETSC_ERR_SUP,0,""); 1932 if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix"); 1933 if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix"); 1934 if (mat->N != x->N) SETERRQ2(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec x: global dim %d %d",mat->N,x->N); 1935 if (mat->M != b->N) SETERRQ2(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec b: global dim %d %d",mat->M,b->N); 1936 if (mat->m != b->n) SETERRQ2(PETSC_ERR_ARG_SIZ,0,"Mat mat,Vec b: local dim %d %d",mat->m,b->n); 1937 1938 PLogEventBegin(MAT_Relax,mat,b,x,0); 1939 ierr =(*mat->ops->relax)(mat,b,omega,flag,shift,its,x);CHKERRQ(ierr); 1940 PLogEventEnd(MAT_Relax,mat,b,x,0); 1941 PetscFunctionReturn(0); 1942 } 1943 1944 #undef __FUNC__ 1945 #define __FUNC__ "MatCopy_Basic" 1946 /* 1947 Default matrix copy routine. 1948 */ 1949 int MatCopy_Basic(Mat A,Mat B,MatStructure str) 1950 { 1951 int ierr,i,rstart,rend,nz,*cwork; 1952 Scalar *vwork; 1953 1954 PetscFunctionBegin; 1955 ierr = MatZeroEntries(B);CHKERRQ(ierr); 1956 ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 1957 for (i=rstart; i<rend; i++) { 1958 ierr = MatGetRow(A,i,&nz,&cwork,&vwork);CHKERRQ(ierr); 1959 ierr = MatSetValues(B,1,&i,nz,cwork,vwork,INSERT_VALUES);CHKERRQ(ierr); 1960 ierr = MatRestoreRow(A,i,&nz,&cwork,&vwork);CHKERRQ(ierr); 1961 } 1962 ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1963 ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1964 PetscFunctionReturn(0); 1965 } 1966 1967 #undef __FUNC__ 1968 #define __FUNC__ "MatCopy" 1969 /*@C 1970 MatCopy - Copys a matrix to another matrix. 1971 1972 Collective on Mat 1973 1974 Input Parameters: 1975 + A - the matrix 1976 - str - SAME_NONZERO_PATTERN or DIFFERENT_NONZERO_PATTERN 1977 1978 Output Parameter: 1979 . B - where the copy is put 1980 1981 Notes: 1982 If you use SAME_NONZERO_PATTERN then the zero matrices had better have the 1983 same nonzero pattern or the routine will crash. 1984 1985 MatCopy() copies the matrix entries of a matrix to another existing 1986 matrix (after first zeroing the second matrix). A related routine is 1987 MatConvert(), which first creates a new matrix and then copies the data. 1988 1989 Level: intermediate 1990 1991 .keywords: matrix, copy, convert 1992 1993 .seealso: MatConvert() 1994 @*/ 1995 int MatCopy(Mat A,Mat B,MatStructure str) 1996 { 1997 int ierr; 1998 1999 PetscFunctionBegin; 2000 PetscValidHeaderSpecific(A,MAT_COOKIE); 2001 PetscValidHeaderSpecific(B,MAT_COOKIE); 2002 PetscCheckSameComm(A,B); 2003 if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix"); 2004 if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix"); 2005 if (A->M != B->M || A->N != B->N) SETERRQ4(PETSC_ERR_ARG_SIZ,0,"Mat A,Mat B: global dim %d %d",A->M,B->M, 2006 A->N,B->N); 2007 2008 PLogEventBegin(MAT_Copy,A,B,0,0); 2009 if (A->ops->copy) { 2010 ierr = (*A->ops->copy)(A,B,str);CHKERRQ(ierr); 2011 } else { /* generic conversion */ 2012 ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr); 2013 } 2014 PLogEventEnd(MAT_Copy,A,B,0,0); 2015 PetscFunctionReturn(0); 2016 } 2017 2018 static int MatConvertersSet = 0; 2019 static int (*MatConverters[MAX_MATRIX_TYPES][MAX_MATRIX_TYPES])(Mat,MatType,Mat*) = 2020 {{0,0,0,0,0,0,0,0,0,0,0,0},{0,0,0,0,0,0,0,0,0,0,0,0}, 2021 {0,0,0,0,0,0,0,0,0,0,0,0},{0,0,0,0,0,0,0,0,0,0,0,0}, 2022 {0,0,0,0,0,0,0,0,0,0,0,0},{0,0,0,0,0,0,0,0,0,0,0,0}, 2023 {0,0,0,0,0,0,0,0,0,0,0,0},{0,0,0,0,0,0,0,0,0,0,0,0}, 2024 {0,0,0,0,0,0,0,0,0,0,0,0},{0,0,0,0,0,0,0,0,0,0,0,0}, 2025 {0,0,0,0,0,0,0,0,0,0,0,0},{0,0,0,0,0,0,0,0,0,0,0,0}}; 2026 2027 #undef __FUNC__ 2028 #define __FUNC__ "MatConvertRegister" 2029 /*@C 2030 MatConvertRegister - Allows one to register a routine that converts between 2031 two matrix types. 2032 2033 Not Collective 2034 2035 Input Parameters: 2036 + intype - the type of matrix (defined in include/mat.h), for example, MATSEQAIJ. 2037 - outtype - new matrix type, or MATSAME 2038 2039 Level: advanced 2040 2041 .seealso: MatConvertRegisterAll() 2042 @*/ 2043 int MatConvertRegister(MatType intype,MatType outtype,int (*converter)(Mat,MatType,Mat*)) 2044 { 2045 PetscFunctionBegin; 2046 MatConverters[intype][outtype] = converter; 2047 MatConvertersSet = 1; 2048 PetscFunctionReturn(0); 2049 } 2050 2051 #undef __FUNC__ 2052 #define __FUNC__ "MatConvert" 2053 /*@C 2054 MatConvert - Converts a matrix to another matrix, either of the same 2055 or different type. 2056 2057 Collective on Mat 2058 2059 Input Parameters: 2060 + mat - the matrix 2061 - newtype - new matrix type. Use MATSAME to create a new matrix of the 2062 same type as the original matrix. 2063 2064 Output Parameter: 2065 . M - pointer to place new matrix 2066 2067 Notes: 2068 MatConvert() first creates a new matrix and then copies the data from 2069 the first matrix. A related routine is MatCopy(), which copies the matrix 2070 entries of one matrix to another already existing matrix context. 2071 2072 Level: intermediate 2073 2074 .keywords: matrix, copy, convert 2075 2076 .seealso: MatCopy(), MatDuplicate() 2077 @*/ 2078 int MatConvert(Mat mat,MatType newtype,Mat *M) 2079 { 2080 int ierr; 2081 2082 PetscFunctionBegin; 2083 PetscValidHeaderSpecific(mat,MAT_COOKIE); 2084 PetscValidPointer(M); 2085 if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix"); 2086 if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix"); 2087 2088 if (newtype > MAX_MATRIX_TYPES || newtype < -1) { 2089 SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,1,"Not a valid matrix type"); 2090 } 2091 *M = 0; 2092 2093 if (!MatConvertersSet) { 2094 ierr = MatLoadRegisterAll();CHKERRQ(ierr); 2095 } 2096 2097 PLogEventBegin(MAT_Convert,mat,0,0,0); 2098 if ((newtype == mat->type || newtype == MATSAME) && mat->ops->duplicate) { 2099 ierr = (*mat->ops->duplicate)(mat,MAT_COPY_VALUES,M);CHKERRQ(ierr); 2100 } else { 2101 if (!MatConvertersSet) { 2102 ierr = MatConvertRegisterAll();CHKERRQ(ierr); 2103 } 2104 if (!MatConverters[mat->type][newtype]) { 2105 SETERRQ(PETSC_ERR_ARG_WRONG,1,"Invalid matrix type, or matrix converter not registered"); 2106 } 2107 ierr = (*MatConverters[mat->type][newtype])(mat,newtype,M);CHKERRQ(ierr); 2108 } 2109 PLogEventEnd(MAT_Convert,mat,0,0,0); 2110 PetscFunctionReturn(0); 2111 } 2112 2113 #undef __FUNC__ 2114 #define __FUNC__ "MatDuplicate" 2115 /*@C 2116 MatDuplicate - Duplicates a matrix including the non-zero structure. 2117 2118 Collective on Mat 2119 2120 Input Parameters: 2121 + mat - the matrix 2122 - op - either MAT_DO_NOT_COPY_VALUES or MAT_COPY_VALUES, cause it to copy nonzero 2123 values as well or not 2124 2125 Output Parameter: 2126 . M - pointer to place new matrix 2127 2128 Level: intermediate 2129 2130 .keywords: matrix, copy, convert, duplicate 2131 2132 .seealso: MatCopy(), MatConvert() 2133 @*/ 2134 int MatDuplicate(Mat mat,MatDuplicateOption op,Mat *M) 2135 { 2136 int ierr; 2137 2138 PetscFunctionBegin; 2139 PetscValidHeaderSpecific(mat,MAT_COOKIE); 2140 PetscValidPointer(M); 2141 if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix"); 2142 if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix"); 2143 2144 *M = 0; 2145 PLogEventBegin(MAT_Convert,mat,0,0,0); 2146 if (!mat->ops->duplicate) { 2147 SETERRQ(PETSC_ERR_SUP,1,"Not written for this matrix type"); 2148 } 2149 ierr = (*mat->ops->duplicate)(mat,op,M);CHKERRQ(ierr); 2150 PLogEventEnd(MAT_Convert,mat,0,0,0); 2151 PetscFunctionReturn(0); 2152 } 2153 2154 #undef __FUNC__ 2155 #define __FUNC__ "MatGetDiagonal" 2156 /*@ 2157 MatGetDiagonal - Gets the diagonal of a matrix. 2158 2159 Collective on Mat and Vec 2160 2161 Input Parameters: 2162 + mat - the matrix 2163 - v - the vector for storing the diagonal 2164 2165 Output Parameter: 2166 . v - the diagonal of the matrix 2167 2168 Notes: 2169 For the SeqAIJ matrix format, this routine may also be called 2170 on a LU factored matrix; in that case it routines the reciprocal of 2171 the diagonal entries in U. It returns the entries permuted by the 2172 row and column permutation used during the symbolic factorization. 2173 2174 Level: intermediate 2175 2176 .keywords: matrix, get, diagonal 2177 2178 .seealso: MatGetRow(), MatGetSubmatrices(), MatGetSubmatrix() 2179 @*/ 2180 int MatGetDiagonal(Mat mat,Vec v) 2181 { 2182 int ierr; 2183 2184 PetscFunctionBegin; 2185 PetscValidHeaderSpecific(mat,MAT_COOKIE); 2186 PetscValidHeaderSpecific(v,VEC_COOKIE); 2187 /* PetscCheckSameComm(mat,v); Could be MPI vector but Seq matrix cause of two submatrix storage */ 2188 if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix"); 2189 if (!mat->ops->getdiagonal) SETERRQ(PETSC_ERR_SUP,0,""); 2190 2191 ierr = (*mat->ops->getdiagonal)(mat,v);CHKERRQ(ierr); 2192 PetscFunctionReturn(0); 2193 } 2194 2195 #undef __FUNC__ 2196 #define __FUNC__ "MatTranspose" 2197 /*@C 2198 MatTranspose - Computes an in-place or out-of-place transpose of a matrix. 2199 2200 Collective on Mat 2201 2202 Input Parameter: 2203 . mat - the matrix to transpose 2204 2205 Output Parameters: 2206 . B - the transpose (or pass in PETSC_NULL for an in-place transpose) 2207 2208 Level: intermediate 2209 2210 .keywords: matrix, transpose 2211 2212 .seealso: MatMultTrans(), MatMultTransAdd() 2213 @*/ 2214 int MatTranspose(Mat mat,Mat *B) 2215 { 2216 int ierr; 2217 2218 PetscFunctionBegin; 2219 PetscValidHeaderSpecific(mat,MAT_COOKIE); 2220 if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix"); 2221 if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix"); 2222 if (!mat->ops->transpose) SETERRQ(PETSC_ERR_SUP,0,""); 2223 ierr = (*mat->ops->transpose)(mat,B);CHKERRQ(ierr); 2224 PetscFunctionReturn(0); 2225 } 2226 2227 #undef __FUNC__ 2228 #define __FUNC__ "MatPermute" 2229 /*@C 2230 MatPermute - Creates a new matrix with rows and columns permuted from the 2231 original. 2232 2233 Collective on Mat 2234 2235 Input Parameters: 2236 + mat - the matrix to permute 2237 . row - row permutation 2238 - col - column permutation 2239 2240 Output Parameters: 2241 . B - the permuted matrix 2242 2243 Level: advanced 2244 2245 .keywords: matrix, transpose 2246 2247 .seealso: MatGetOrdering() 2248 @*/ 2249 int MatPermute(Mat mat,IS row,IS col,Mat *B) 2250 { 2251 int ierr; 2252 2253 PetscFunctionBegin; 2254 PetscValidHeaderSpecific(mat,MAT_COOKIE); 2255 PetscValidHeaderSpecific(row,IS_COOKIE); 2256 PetscValidHeaderSpecific(col,IS_COOKIE); 2257 if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix"); 2258 if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix"); 2259 if (!mat->ops->permute) SETERRQ(PETSC_ERR_SUP,0,""); 2260 ierr = (*mat->ops->permute)(mat,row,col,B);CHKERRQ(ierr); 2261 PetscFunctionReturn(0); 2262 } 2263 2264 #undef __FUNC__ 2265 #define __FUNC__ "MatEqual" 2266 /*@ 2267 MatEqual - Compares two matrices. 2268 2269 Collective on Mat 2270 2271 Input Parameters: 2272 + A - the first matrix 2273 - B - the second matrix 2274 2275 Output Parameter: 2276 . flg - PETSC_TRUE if the matrices are equal; PETSC_FALSE otherwise. 2277 2278 Level: intermediate 2279 2280 .keywords: matrix, equal, equivalent 2281 @*/ 2282 int MatEqual(Mat A,Mat B,PetscTruth *flg) 2283 { 2284 int ierr; 2285 2286 PetscFunctionBegin; 2287 PetscValidHeaderSpecific(A,MAT_COOKIE); 2288 PetscValidHeaderSpecific(B,MAT_COOKIE); 2289 PetscValidIntPointer(flg); 2290 PetscCheckSameComm(A,B); 2291 if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix"); 2292 if (!B->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix"); 2293 if (A->M != B->M || A->N != B->N) SETERRQ4(PETSC_ERR_ARG_SIZ,0,"Mat A,Mat B: global dim %d %d %d %d", 2294 A->M,B->M,A->N,B->N); 2295 if (!A->ops->equal) SETERRQ(PETSC_ERR_SUP,0,""); 2296 ierr = (*A->ops->equal)(A,B,flg);CHKERRQ(ierr); 2297 PetscFunctionReturn(0); 2298 } 2299 2300 #undef __FUNC__ 2301 #define __FUNC__ "MatDiagonalScale" 2302 /*@ 2303 MatDiagonalScale - Scales a matrix on the left and right by diagonal 2304 matrices that are stored as vectors. Either of the two scaling 2305 matrices can be PETSC_NULL. 2306 2307 Collective on Mat 2308 2309 Input Parameters: 2310 + mat - the matrix to be scaled 2311 . l - the left scaling vector (or PETSC_NULL) 2312 - r - the right scaling vector (or PETSC_NULL) 2313 2314 Notes: 2315 MatDiagonalScale() computes A = LAR, where 2316 L = a diagonal matrix, R = a diagonal matrix 2317 2318 Level: intermediate 2319 2320 .keywords: matrix, diagonal, scale 2321 2322 .seealso: MatScale() 2323 @*/ 2324 int MatDiagonalScale(Mat mat,Vec l,Vec r) 2325 { 2326 int ierr; 2327 2328 PetscFunctionBegin; 2329 PetscValidHeaderSpecific(mat,MAT_COOKIE); 2330 if (!mat->ops->diagonalscale) SETERRQ(PETSC_ERR_SUP,0,""); 2331 if (l) {PetscValidHeaderSpecific(l,VEC_COOKIE);PetscCheckSameComm(mat,l);} 2332 if (r) {PetscValidHeaderSpecific(r,VEC_COOKIE);PetscCheckSameComm(mat,r);} 2333 if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix"); 2334 if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix"); 2335 2336 PLogEventBegin(MAT_Scale,mat,0,0,0); 2337 ierr = (*mat->ops->diagonalscale)(mat,l,r);CHKERRQ(ierr); 2338 PLogEventEnd(MAT_Scale,mat,0,0,0); 2339 PetscFunctionReturn(0); 2340 } 2341 2342 #undef __FUNC__ 2343 #define __FUNC__ "MatScale" 2344 /*@ 2345 MatScale - Scales all elements of a matrix by a given number. 2346 2347 Collective on Mat 2348 2349 Input Parameters: 2350 + mat - the matrix to be scaled 2351 - a - the scaling value 2352 2353 Output Parameter: 2354 . mat - the scaled matrix 2355 2356 Level: intermediate 2357 2358 .keywords: matrix, scale 2359 2360 .seealso: MatDiagonalScale() 2361 @*/ 2362 int MatScale(Scalar *a,Mat mat) 2363 { 2364 int ierr; 2365 2366 PetscFunctionBegin; 2367 PetscValidHeaderSpecific(mat,MAT_COOKIE); 2368 PetscValidScalarPointer(a); 2369 if (!mat->ops->scale) SETERRQ(PETSC_ERR_SUP,0,""); 2370 if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix"); 2371 if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix"); 2372 2373 PLogEventBegin(MAT_Scale,mat,0,0,0); 2374 ierr = (*mat->ops->scale)(a,mat);CHKERRQ(ierr); 2375 PLogEventEnd(MAT_Scale,mat,0,0,0); 2376 PetscFunctionReturn(0); 2377 } 2378 2379 #undef __FUNC__ 2380 #define __FUNC__ "MatNorm" 2381 /*@ 2382 MatNorm - Calculates various norms of a matrix. 2383 2384 Collective on Mat 2385 2386 Input Parameters: 2387 + mat - the matrix 2388 - type - the type of norm, NORM_1, NORM_2, NORM_FROBENIUS, NORM_INFINITY 2389 2390 Output Parameters: 2391 . norm - the resulting norm 2392 2393 Level: intermediate 2394 2395 .keywords: matrix, norm, Frobenius 2396 @*/ 2397 int MatNorm(Mat mat,NormType type,double *norm) 2398 { 2399 int ierr; 2400 2401 PetscFunctionBegin; 2402 PetscValidHeaderSpecific(mat,MAT_COOKIE); 2403 PetscValidScalarPointer(norm); 2404 2405 if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix"); 2406 if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix"); 2407 if (!mat->ops->norm) SETERRQ(PETSC_ERR_SUP,0,"Not for this matrix type"); 2408 ierr = (*mat->ops->norm)(mat,type,norm);CHKERRQ(ierr); 2409 PetscFunctionReturn(0); 2410 } 2411 2412 /* 2413 This variable is used to prevent counting of MatAssemblyBegin() that 2414 are called from within a MatAssemblyEnd(). 2415 */ 2416 static int MatAssemblyEnd_InUse = 0; 2417 #undef __FUNC__ 2418 #define __FUNC__ "MatAssemblyBegin" 2419 /*@ 2420 MatAssemblyBegin - Begins assembling the matrix. This routine should 2421 be called after completing all calls to MatSetValues(). 2422 2423 Collective on Mat 2424 2425 Input Parameters: 2426 + mat - the matrix 2427 - type - type of assembly, either MAT_FLUSH_ASSEMBLY or MAT_FINAL_ASSEMBLY 2428 2429 Notes: 2430 MatSetValues() generally caches the values. The matrix is ready to 2431 use only after MatAssemblyBegin() and MatAssemblyEnd() have been called. 2432 Use MAT_FLUSH_ASSEMBLY when switching between ADD_VALUES and INSERT_VALUES 2433 in MatSetValues(); use MAT_FINAL_ASSEMBLY for the final assembly before 2434 using the matrix. 2435 2436 Level: beginner 2437 2438 .keywords: matrix, assembly, assemble, begin 2439 2440 .seealso: MatAssemblyEnd(), MatSetValues(), MatAssembled() 2441 @*/ 2442 int MatAssemblyBegin(Mat mat,MatAssemblyType type) 2443 { 2444 int ierr; 2445 2446 PetscFunctionBegin; 2447 PetscValidHeaderSpecific(mat,MAT_COOKIE); 2448 if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix.\n did you forget to call MatSetUnfactored()?"); 2449 if (mat->assembled) { 2450 mat->was_assembled = PETSC_TRUE; 2451 mat->assembled = PETSC_FALSE; 2452 } 2453 if (!MatAssemblyEnd_InUse) { 2454 PLogEventBegin(MAT_AssemblyBegin,mat,0,0,0); 2455 if (mat->ops->assemblybegin){ierr = (*mat->ops->assemblybegin)(mat,type);CHKERRQ(ierr);} 2456 PLogEventEnd(MAT_AssemblyBegin,mat,0,0,0); 2457 } else { 2458 if (mat->ops->assemblybegin){ierr = (*mat->ops->assemblybegin)(mat,type);CHKERRQ(ierr);} 2459 } 2460 PetscFunctionReturn(0); 2461 } 2462 2463 #undef __FUNC__ 2464 #define __FUNC__ "MatAssembed" 2465 /*@ 2466 MatAssembled - Indicates if a matrix has been assembled and is ready for 2467 use; for example, in matrix-vector product. 2468 2469 Collective on Mat 2470 2471 Input Parameter: 2472 . mat - the matrix 2473 2474 Output Parameter: 2475 . assembled - PETSC_TRUE or PETSC_FALSE 2476 2477 Level: advanced 2478 2479 .keywords: matrix, assembly, assemble, begin 2480 2481 .seealso: MatAssemblyEnd(), MatSetValues(), MatAssemblyBegin() 2482 @*/ 2483 int MatAssembled(Mat mat,PetscTruth *assembled) 2484 { 2485 PetscFunctionBegin; 2486 PetscValidHeaderSpecific(mat,MAT_COOKIE); 2487 *assembled = mat->assembled; 2488 PetscFunctionReturn(0); 2489 } 2490 2491 #undef __FUNC__ 2492 #define __FUNC__ "MatView_Private" 2493 /* 2494 Processes command line options to determine if/how a matrix 2495 is to be viewed. Called by MatAssemblyEnd() and MatLoad(). 2496 */ 2497 int MatView_Private(Mat mat) 2498 { 2499 int ierr,flg; 2500 2501 PetscFunctionBegin; 2502 ierr = OptionsHasName(PETSC_NULL,"-mat_view_info",&flg);CHKERRQ(ierr); 2503 if (flg) { 2504 ierr = ViewerPushFormat(VIEWER_STDOUT_(mat->comm),VIEWER_FORMAT_ASCII_INFO,0);CHKERRQ(ierr); 2505 ierr = MatView(mat,VIEWER_STDOUT_(mat->comm));CHKERRQ(ierr); 2506 ierr = ViewerPopFormat(VIEWER_STDOUT_(mat->comm));CHKERRQ(ierr); 2507 } 2508 ierr = OptionsHasName(PETSC_NULL,"-mat_view_info_detailed",&flg);CHKERRQ(ierr); 2509 if (flg) { 2510 ierr = ViewerPushFormat(VIEWER_STDOUT_(mat->comm),VIEWER_FORMAT_ASCII_INFO_LONG,0);CHKERRQ(ierr); 2511 ierr = MatView(mat,VIEWER_STDOUT_(mat->comm));CHKERRQ(ierr); 2512 ierr = ViewerPopFormat(VIEWER_STDOUT_(mat->comm));CHKERRQ(ierr); 2513 } 2514 ierr = OptionsHasName(PETSC_NULL,"-mat_view",&flg);CHKERRQ(ierr); 2515 if (flg) { 2516 ierr = MatView(mat,VIEWER_STDOUT_(mat->comm));CHKERRQ(ierr); 2517 } 2518 ierr = OptionsHasName(PETSC_NULL,"-mat_view_matlab",&flg);CHKERRQ(ierr); 2519 if (flg) { 2520 ierr = ViewerPushFormat(VIEWER_STDOUT_(mat->comm),VIEWER_FORMAT_ASCII_MATLAB,"M");CHKERRQ(ierr); 2521 ierr = MatView(mat,VIEWER_STDOUT_(mat->comm));CHKERRQ(ierr); 2522 ierr = ViewerPopFormat(VIEWER_STDOUT_(mat->comm));CHKERRQ(ierr); 2523 } 2524 ierr = OptionsHasName(PETSC_NULL,"-mat_view_draw",&flg);CHKERRQ(ierr); 2525 if (flg) { 2526 ierr = OptionsHasName(0,"-mat_view_contour",&flg);CHKERRQ(ierr); 2527 if (flg) { 2528 ViewerPushFormat(VIEWER_DRAW_(mat->comm),VIEWER_FORMAT_DRAW_CONTOUR,0);CHKERRQ(ierr); 2529 } 2530 ierr = MatView(mat,VIEWER_DRAW_(mat->comm));CHKERRQ(ierr); 2531 ierr = ViewerFlush(VIEWER_DRAW_(mat->comm));CHKERRQ(ierr); 2532 if (flg) { 2533 ViewerPopFormat(VIEWER_DRAW_(mat->comm));CHKERRQ(ierr); 2534 } 2535 } 2536 ierr = OptionsHasName(PETSC_NULL,"-mat_view_socket",&flg);CHKERRQ(ierr); 2537 if (flg) { 2538 ierr = MatView(mat,VIEWER_SOCKET_(mat->comm));CHKERRQ(ierr); 2539 ierr = ViewerFlush(VIEWER_SOCKET_(mat->comm));CHKERRQ(ierr); 2540 } 2541 PetscFunctionReturn(0); 2542 } 2543 2544 #undef __FUNC__ 2545 #define __FUNC__ "MatAssemblyEnd" 2546 /*@ 2547 MatAssemblyEnd - Completes assembling the matrix. This routine should 2548 be called after MatAssemblyBegin(). 2549 2550 Collective on Mat 2551 2552 Input Parameters: 2553 + mat - the matrix 2554 - type - type of assembly, either MAT_FLUSH_ASSEMBLY or MAT_FINAL_ASSEMBLY 2555 2556 Options Database Keys: 2557 + -mat_view_info - Prints info on matrix at conclusion of MatEndAssembly() 2558 . -mat_view_info_detailed - Prints more detailed info 2559 . -mat_view - Prints matrix in ASCII format 2560 . -mat_view_matlab - Prints matrix in Matlab format 2561 . -mat_view_draw - Draws nonzero structure of matrix, using MatView() and DrawOpenX(). 2562 . -display <name> - Sets display name (default is host) 2563 - -draw_pause <sec> - Sets number of seconds to pause after display 2564 2565 Notes: 2566 MatSetValues() generally caches the values. The matrix is ready to 2567 use only after MatAssemblyBegin() and MatAssemblyEnd() have been called. 2568 Use MAT_FLUSH_ASSEMBLY when switching between ADD_VALUES and INSERT_VALUES 2569 in MatSetValues(); use MAT_FINAL_ASSEMBLY for the final assembly before 2570 using the matrix. 2571 2572 Level: beginner 2573 2574 .keywords: matrix, assembly, assemble, end 2575 2576 .seealso: MatAssemblyBegin(), MatSetValues(), DrawOpenX(), MatView(), MatAssembled() 2577 @*/ 2578 int MatAssemblyEnd(Mat mat,MatAssemblyType type) 2579 { 2580 int ierr; 2581 static int inassm = 0; 2582 2583 PetscFunctionBegin; 2584 PetscValidHeaderSpecific(mat,MAT_COOKIE); 2585 2586 inassm++; 2587 MatAssemblyEnd_InUse++; 2588 if (MatAssemblyEnd_InUse == 1) { /* Do the logging only the first time through */ 2589 PLogEventBegin(MAT_AssemblyEnd,mat,0,0,0); 2590 if (mat->ops->assemblyend) { 2591 ierr = (*mat->ops->assemblyend)(mat,type);CHKERRQ(ierr); 2592 } 2593 PLogEventEnd(MAT_AssemblyEnd,mat,0,0,0); 2594 } else { 2595 if (mat->ops->assemblyend) { 2596 ierr = (*mat->ops->assemblyend)(mat,type);CHKERRQ(ierr); 2597 } 2598 } 2599 2600 /* Flush assembly is not a true assembly */ 2601 if (type != MAT_FLUSH_ASSEMBLY) { 2602 mat->assembled = PETSC_TRUE; mat->num_ass++; 2603 } 2604 mat->insertmode = NOT_SET_VALUES; 2605 MatAssemblyEnd_InUse--; 2606 2607 if (inassm == 1 && type != MAT_FLUSH_ASSEMBLY) { 2608 ierr = MatView_Private(mat);CHKERRQ(ierr); 2609 } 2610 inassm--; 2611 PetscFunctionReturn(0); 2612 } 2613 2614 2615 #undef __FUNC__ 2616 #define __FUNC__ "MatCompress" 2617 /*@ 2618 MatCompress - Tries to store the matrix in as little space as 2619 possible. May fail if memory is already fully used, since it 2620 tries to allocate new space. 2621 2622 Collective on Mat 2623 2624 Input Parameters: 2625 . mat - the matrix 2626 2627 Level: advanced 2628 2629 .keywords: matrix, compress 2630 @*/ 2631 int MatCompress(Mat mat) 2632 { 2633 int ierr; 2634 2635 PetscFunctionBegin; 2636 PetscValidHeaderSpecific(mat,MAT_COOKIE); 2637 if (mat->ops->compress) {ierr = (*mat->ops->compress)(mat);CHKERRQ(ierr);} 2638 PetscFunctionReturn(0); 2639 } 2640 2641 #undef __FUNC__ 2642 #define __FUNC__ "MatSetOption" 2643 /*@ 2644 MatSetOption - Sets a parameter option for a matrix. Some options 2645 may be specific to certain storage formats. Some options 2646 determine how values will be inserted (or added). Sorted, 2647 row-oriented input will generally assemble the fastest. The default 2648 is row-oriented, nonsorted input. 2649 2650 Collective on Mat 2651 2652 Input Parameters: 2653 + mat - the matrix 2654 - option - the option, one of those listed below (and possibly others), 2655 e.g., MAT_ROWS_SORTED, MAT_NEW_NONZERO_LOCATION_ERR 2656 2657 Options Describing Matrix Structure: 2658 + MAT_SYMMETRIC - symmetric in terms of both structure and value 2659 - MAT_STRUCTURALLY_SYMMETRIC - symmetric nonzero structure 2660 2661 Options For Use with MatSetValues(): 2662 Insert a logically dense subblock, which can be 2663 + MAT_ROW_ORIENTED - row-oriented 2664 . MAT_COLUMN_ORIENTED - column-oriented 2665 . MAT_ROWS_SORTED - sorted by row 2666 . MAT_ROWS_UNSORTED - not sorted by row 2667 . MAT_COLUMNS_SORTED - sorted by column 2668 - MAT_COLUMNS_UNSORTED - not sorted by column 2669 2670 Not these options reflect the data you pass in with MatSetValues(); it has 2671 nothing to do with how the data is stored internally in the matrix 2672 data structure. 2673 2674 When (re)assembling a matrix, we can restrict the input for 2675 efficiency/debugging purposes. These options include 2676 + MAT_NO_NEW_NONZERO_LOCATIONS - additional insertions will not be 2677 allowed if they generate a new nonzero 2678 . MAT_YES_NEW_NONZERO_LOCATIONS - additional insertions will be allowed 2679 . MAT_NO_NEW_DIAGONALS - additional insertions will not be allowed if 2680 they generate a nonzero in a new diagonal (for block diagonal format only) 2681 . MAT_YES_NEW_DIAGONALS - new diagonals will be allowed (for block diagonal format only) 2682 . MAT_IGNORE_OFF_PROC_ENTRIES - drops off-processor entries 2683 . MAT_NEW_NONZERO_LOCATION_ERR - generates an error for new matrix entry 2684 - MAT_USE_HASH_TABLE - uses a hash table to speed up matrix assembly 2685 2686 Notes: 2687 Some options are relevant only for particular matrix types and 2688 are thus ignored by others. Other options are not supported by 2689 certain matrix types and will generate an error message if set. 2690 2691 If using a Fortran 77 module to compute a matrix, one may need to 2692 use the column-oriented option (or convert to the row-oriented 2693 format). 2694 2695 MAT_NO_NEW_NONZERO_LOCATIONS indicates that any add or insertion 2696 that would generate a new entry in the nonzero structure is instead 2697 ignored. Thus, if memory has not alredy been allocated for this particular 2698 data, then the insertion is ignored. For dense matrices, in which 2699 the entire array is allocated, no entries are ever ignored. 2700 2701 MAT_NEW_NONZERO_LOCATION_ERR indicates that any add or insertion 2702 that would generate a new entry in the nonzero structure instead produces 2703 an error. (Currently supported for AIJ and BAIJ formats only.) 2704 This is a useful flag when using SAME_NONZERO_PATTERN in calling 2705 SLESSetOperators() to ensure that the nonzero pattern truely does 2706 remain unchanged. 2707 2708 MAT_NEW_NONZERO_ALLOCATION_ERR indicates that any add or insertion 2709 that would generate a new entry that has not been preallocated will 2710 instead produce an error. (Currently supported for AIJ and BAIJ formats 2711 only.) This is a useful flag when debugging matrix memory preallocation. 2712 2713 MAT_IGNORE_OFF_PROC_ENTRIES indicates entries destined for 2714 other processors should be dropped, rather than stashed. 2715 This is useful if you know that the "owning" processor is also 2716 always generating the correct matrix entries, so that PETSc need 2717 not transfer duplicate entries generated on another processor. 2718 2719 MAT_USE_HASH_TABLE indicates that a hash table be used to improve the 2720 searches during matrix assembly. When this flag is set, the hash table 2721 is created during the first Matrix Assembly. This hash table is 2722 used the next time through, during MatSetVaules()/MatSetVaulesBlocked() 2723 to improve the searching of indices. MAT_NO_NEW_NONZERO_LOCATIONS flag 2724 should be used with MAT_USE_HASH_TABLE flag. This option is currently 2725 supported by MATMPIBAIJ format only. 2726 2727 Level: intermediate 2728 2729 .keywords: matrix, option, row-oriented, column-oriented, sorted, nonzero 2730 @*/ 2731 int MatSetOption(Mat mat,MatOption op) 2732 { 2733 int ierr; 2734 2735 PetscFunctionBegin; 2736 PetscValidHeaderSpecific(mat,MAT_COOKIE); 2737 if (mat->ops->setoption) {ierr = (*mat->ops->setoption)(mat,op);CHKERRQ(ierr);} 2738 PetscFunctionReturn(0); 2739 } 2740 2741 #undef __FUNC__ 2742 #define __FUNC__ "MatZeroEntries" 2743 /*@ 2744 MatZeroEntries - Zeros all entries of a matrix. For sparse matrices 2745 this routine retains the old nonzero structure. 2746 2747 Collective on Mat 2748 2749 Input Parameters: 2750 . mat - the matrix 2751 2752 Level: intermediate 2753 2754 .keywords: matrix, zero, entries 2755 2756 .seealso: MatZeroRows() 2757 @*/ 2758 int MatZeroEntries(Mat mat) 2759 { 2760 int ierr; 2761 2762 PetscFunctionBegin; 2763 PetscValidHeaderSpecific(mat,MAT_COOKIE); 2764 if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix"); 2765 if (!mat->ops->zeroentries) SETERRQ(PETSC_ERR_SUP,0,""); 2766 2767 PLogEventBegin(MAT_ZeroEntries,mat,0,0,0); 2768 ierr = (*mat->ops->zeroentries)(mat);CHKERRQ(ierr); 2769 PLogEventEnd(MAT_ZeroEntries,mat,0,0,0); 2770 PetscFunctionReturn(0); 2771 } 2772 2773 #undef __FUNC__ 2774 #define __FUNC__ "MatZeroRows" 2775 /*@C 2776 MatZeroRows - Zeros all entries (except possibly the main diagonal) 2777 of a set of rows of a matrix. 2778 2779 Collective on Mat 2780 2781 Input Parameters: 2782 + mat - the matrix 2783 . is - index set of rows to remove 2784 - diag - pointer to value put in all diagonals of eliminated rows. 2785 Note that diag is not a pointer to an array, but merely a 2786 pointer to a single value. 2787 2788 Notes: 2789 For the AIJ matrix formats this removes the old nonzero structure, 2790 but does not release memory. For the dense and block diagonal 2791 formats this does not alter the nonzero structure. 2792 2793 The user can set a value in the diagonal entry (or for the AIJ and 2794 row formats can optionally remove the main diagonal entry from the 2795 nonzero structure as well, by passing a null pointer (PETSC_NULL 2796 in C or PETSC_NULL_SCALAR in Fortran) as the final argument). 2797 2798 For the parallel case, all processes that share the matrix (i.e., 2799 those in the communicator used for matrix creation) MUST call this 2800 routine, regardless of whether any rows being zeroed are owned by 2801 them. 2802 2803 Level: intermediate 2804 2805 .keywords: matrix, zero, rows, boundary conditions 2806 2807 .seealso: MatZeroEntries(), 2808 @*/ 2809 int MatZeroRows(Mat mat,IS is, Scalar *diag) 2810 { 2811 int ierr; 2812 2813 PetscFunctionBegin; 2814 PetscValidHeaderSpecific(mat,MAT_COOKIE); 2815 PetscValidHeaderSpecific(is,IS_COOKIE); 2816 if (diag) PetscValidScalarPointer(diag); 2817 if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix"); 2818 if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix"); 2819 if (!mat->ops->zerorows) SETERRQ(PETSC_ERR_SUP,0,""); 2820 2821 ierr = (*mat->ops->zerorows)(mat,is,diag);CHKERRQ(ierr); 2822 ierr = MatView_Private(mat);CHKERRQ(ierr); 2823 PetscFunctionReturn(0); 2824 } 2825 2826 #undef __FUNC__ 2827 #define __FUNC__ "MatZeroRowsLocal" 2828 /*@C 2829 MatZeroRowsLocal - Zeros all entries (except possibly the main diagonal) 2830 of a set of rows of a matrix; using local numbering of rows. 2831 2832 Collective on Mat 2833 2834 Input Parameters: 2835 + mat - the matrix 2836 . is - index set of rows to remove 2837 - diag - pointer to value put in all diagonals of eliminated rows. 2838 Note that diag is not a pointer to an array, but merely a 2839 pointer to a single value. 2840 2841 Notes: 2842 For the AIJ matrix formats this removes the old nonzero structure, 2843 but does not release memory. For the dense and block diagonal 2844 formats this does not alter the nonzero structure. 2845 2846 The user can set a value in the diagonal entry (or for the AIJ and 2847 row formats can optionally remove the main diagonal entry from the 2848 nonzero structure as well, by passing a null pointer as the final 2849 argument). 2850 2851 Level: intermediate 2852 2853 .keywords: matrix, zero, rows, boundary conditions 2854 2855 .seealso: MatZeroEntries(), 2856 @*/ 2857 int MatZeroRowsLocal(Mat mat,IS is, Scalar *diag) 2858 { 2859 int ierr; 2860 IS newis; 2861 2862 PetscFunctionBegin; 2863 PetscValidHeaderSpecific(mat,MAT_COOKIE); 2864 PetscValidHeaderSpecific(is,IS_COOKIE); 2865 if (diag) PetscValidScalarPointer(diag); 2866 if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix"); 2867 if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix"); 2868 if (!mat->ops->zerorows) SETERRQ(PETSC_ERR_SUP,0,""); 2869 2870 ierr = ISLocalToGlobalMappingApplyIS(mat->mapping,is,&newis);CHKERRQ(ierr); 2871 ierr = (*mat->ops->zerorows)(mat,newis,diag);CHKERRQ(ierr); 2872 ierr = ISDestroy(newis);CHKERRQ(ierr); 2873 PetscFunctionReturn(0); 2874 } 2875 2876 #undef __FUNC__ 2877 #define __FUNC__ "MatGetSize" 2878 /*@ 2879 MatGetSize - Returns the numbers of rows and columns in a matrix. 2880 2881 Not Collective 2882 2883 Input Parameter: 2884 . mat - the matrix 2885 2886 Output Parameters: 2887 + m - the number of global rows 2888 - n - the number of global columns 2889 2890 Level: beginner 2891 2892 .keywords: matrix, dimension, size, rows, columns, global, get 2893 2894 .seealso: MatGetLocalSize() 2895 @*/ 2896 int MatGetSize(Mat mat,int *m,int* n) 2897 { 2898 int ierr; 2899 2900 PetscFunctionBegin; 2901 PetscValidHeaderSpecific(mat,MAT_COOKIE); 2902 ierr = (*mat->ops->getsize)(mat,m,n);CHKERRQ(ierr); 2903 PetscFunctionReturn(0); 2904 } 2905 2906 #undef __FUNC__ 2907 #define __FUNC__ "MatGetLocalSize" 2908 /*@ 2909 MatGetLocalSize - Returns the number of rows and columns in a matrix 2910 stored locally. This information may be implementation dependent, so 2911 use with care. 2912 2913 Not Collective 2914 2915 Input Parameters: 2916 . mat - the matrix 2917 2918 Output Parameters: 2919 + m - the number of local rows 2920 - n - the number of local columns 2921 2922 Level: beginner 2923 2924 .keywords: matrix, dimension, size, local, rows, columns, get 2925 2926 .seealso: MatGetSize() 2927 @*/ 2928 int MatGetLocalSize(Mat mat,int *m,int* n) 2929 { 2930 int ierr; 2931 2932 PetscFunctionBegin; 2933 PetscValidHeaderSpecific(mat,MAT_COOKIE); 2934 ierr = (*mat->ops->getlocalsize)(mat,m,n);CHKERRQ(ierr); 2935 PetscFunctionReturn(0); 2936 } 2937 2938 #undef __FUNC__ 2939 #define __FUNC__ "MatGetOwnershipRange" 2940 /*@ 2941 MatGetOwnershipRange - Returns the range of matrix rows owned by 2942 this processor, assuming that the matrix is laid out with the first 2943 n1 rows on the first processor, the next n2 rows on the second, etc. 2944 For certain parallel layouts this range may not be well defined. 2945 2946 Not Collective 2947 2948 Input Parameters: 2949 . mat - the matrix 2950 2951 Output Parameters: 2952 + m - the global index of the first local row 2953 - n - one more than the global index of the last local row 2954 2955 Level: beginner 2956 2957 .keywords: matrix, get, range, ownership 2958 @*/ 2959 int MatGetOwnershipRange(Mat mat,int *m,int* n) 2960 { 2961 int ierr; 2962 2963 PetscFunctionBegin; 2964 PetscValidHeaderSpecific(mat,MAT_COOKIE); 2965 PetscValidIntPointer(m); 2966 PetscValidIntPointer(n); 2967 if (!mat->ops->getownershiprange) SETERRQ(PETSC_ERR_SUP,0,""); 2968 ierr = (*mat->ops->getownershiprange)(mat,m,n);CHKERRQ(ierr); 2969 PetscFunctionReturn(0); 2970 } 2971 2972 #undef __FUNC__ 2973 #define __FUNC__ "MatILUFactorSymbolic" 2974 /*@ 2975 MatILUFactorSymbolic - Performs symbolic ILU factorization of a matrix. 2976 Uses levels of fill only, not drop tolerance. Use MatLUFactorNumeric() 2977 to complete the factorization. 2978 2979 Collective on Mat 2980 2981 Input Parameters: 2982 + mat - the matrix 2983 . row - row permutation 2984 . column - column permutation 2985 - info - structure containing 2986 $ levels - number of levels of fill. 2987 $ expected fill - as ratio of original fill. 2988 $ 1 or 0 - indicating force fill on diagonal (improves robustness for matrices 2989 missing diagonal entries) 2990 2991 Output Parameters: 2992 . fact - new matrix that has been symbolically factored 2993 2994 Notes: 2995 See the users manual for additional information about 2996 choosing the fill factor for better efficiency. 2997 2998 Most users should employ the simplified SLES interface for linear solvers 2999 instead of working directly with matrix algebra routines such as this. 3000 See, e.g., SLESCreate(). 3001 3002 Level: developer 3003 3004 .keywords: matrix, factor, incomplete, ILU, symbolic, fill 3005 3006 .seealso: MatLUFactorSymbolic(), MatLUFactorNumeric(), MatCholeskyFactor() 3007 MatGetOrdering() 3008 3009 @*/ 3010 int MatILUFactorSymbolic(Mat mat,IS row,IS col,MatILUInfo *info,Mat *fact) 3011 { 3012 int ierr; 3013 3014 PetscFunctionBegin; 3015 PetscValidHeaderSpecific(mat,MAT_COOKIE); 3016 PetscValidPointer(fact); 3017 if (info && info->levels < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,0,"Levels of fill negative %d",info->levels); 3018 if (info && info->fill < 1.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,0,"Expected fill less than 1.0 %g",info->fill); 3019 if (!mat->ops->ilufactorsymbolic) SETERRQ(PETSC_ERR_SUP,0,"Only MatCreateMPIRowbs() matrices support parallel ILU"); 3020 if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix"); 3021 if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix"); 3022 3023 PLogEventBegin(MAT_ILUFactorSymbolic,mat,row,col,0); 3024 ierr = (*mat->ops->ilufactorsymbolic)(mat,row,col,info,fact);CHKERRQ(ierr); 3025 PLogEventEnd(MAT_ILUFactorSymbolic,mat,row,col,0); 3026 PetscFunctionReturn(0); 3027 } 3028 3029 #undef __FUNC__ 3030 #define __FUNC__ "MatIncompleteCholeskyFactorSymbolic" 3031 /*@ 3032 MatIncompleteCholeskyFactorSymbolic - Performs symbolic incomplete 3033 Cholesky factorization for a symmetric matrix. Use 3034 MatCholeskyFactorNumeric() to complete the factorization. 3035 3036 Collective on Mat 3037 3038 Input Parameters: 3039 + mat - the matrix 3040 . perm - row and column permutation 3041 . fill - levels of fill 3042 - f - expected fill as ratio of original fill 3043 3044 Output Parameter: 3045 . fact - the factored matrix 3046 3047 Notes: 3048 Currently only no-fill factorization is supported. 3049 3050 Most users should employ the simplified SLES interface for linear solvers 3051 instead of working directly with matrix algebra routines such as this. 3052 See, e.g., SLESCreate(). 3053 3054 Level: developer 3055 3056 .keywords: matrix, factor, incomplete, ICC, Cholesky, symbolic, fill 3057 3058 .seealso: MatCholeskyFactorNumeric(), MatCholeskyFactor() 3059 @*/ 3060 int MatIncompleteCholeskyFactorSymbolic(Mat mat,IS perm,double f,int fill,Mat *fact) 3061 { 3062 int ierr; 3063 3064 PetscFunctionBegin; 3065 PetscValidHeaderSpecific(mat,MAT_COOKIE); 3066 PetscValidPointer(fact); 3067 if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix"); 3068 if (fill < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,0,"Fill negative %d",fill); 3069 if (f < 1.0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,0,"Expected fill less than 1.0 %g",f); 3070 if (!mat->ops->incompletecholeskyfactorsymbolic) SETERRQ(PETSC_ERR_SUP,0,"Currently only MatCreateMPIRowbs() matrices support ICC in parallel"); 3071 if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix"); 3072 3073 PLogEventBegin(MAT_IncompleteCholeskyFactorSymbolic,mat,perm,0,0); 3074 ierr = (*mat->ops->incompletecholeskyfactorsymbolic)(mat,perm,f,fill,fact);CHKERRQ(ierr); 3075 PLogEventEnd(MAT_IncompleteCholeskyFactorSymbolic,mat,perm,0,0); 3076 PetscFunctionReturn(0); 3077 } 3078 3079 #undef __FUNC__ 3080 #define __FUNC__ "MatGetArray" 3081 /*@C 3082 MatGetArray - Returns a pointer to the element values in the matrix. 3083 The result of this routine is dependent on the underlying matrix data 3084 structure, and may not even work for certain matrix types. You MUST 3085 call MatRestoreArray() when you no longer need to access the array. 3086 3087 Not Collective 3088 3089 Input Parameter: 3090 . mat - the matrix 3091 3092 Output Parameter: 3093 . v - the location of the values 3094 3095 Currently returns an array only for the dense formats, giving access to 3096 the local portion of the matrix in the usual Fortran column-oriented format. 3097 3098 Fortran Note: 3099 This routine is used differently from Fortran, e.g., 3100 .vb 3101 Mat mat 3102 Scalar mat_array(1) 3103 PetscOffset i_mat 3104 int ierr 3105 call MatGetArray(mat,mat_array,i_mat,ierr) 3106 3107 C Access first local entry in matrix; note that array is 3108 C treated as one dimensional 3109 value = mat_array(i_mat + 1) 3110 3111 [... other code ...] 3112 call MatRestoreArray(mat,mat_array,i_mat,ierr) 3113 .ve 3114 3115 See the Fortran chapter of the users manual and 3116 petsc/src/mat/examples/tests for details. 3117 3118 Level: advanced 3119 3120 .keywords: matrix, array, elements, values 3121 3122 .seealso: MatRestoreArray(), MatGetArrayF90() 3123 @*/ 3124 int MatGetArray(Mat mat,Scalar **v) 3125 { 3126 int ierr; 3127 3128 PetscFunctionBegin; 3129 PetscValidHeaderSpecific(mat,MAT_COOKIE); 3130 PetscValidPointer(v); 3131 if (!mat->ops->getarray) SETERRQ(PETSC_ERR_SUP,0,""); 3132 ierr = (*mat->ops->getarray)(mat,v);CHKERRQ(ierr); 3133 PetscFunctionReturn(0); 3134 } 3135 3136 #undef __FUNC__ 3137 #define __FUNC__ "MatRestoreArray" 3138 /*@C 3139 MatRestoreArray - Restores the matrix after MatGetArray() has been called. 3140 3141 Not Collective 3142 3143 Input Parameter: 3144 + mat - the matrix 3145 - v - the location of the values 3146 3147 Fortran Note: 3148 This routine is used differently from Fortran, e.g., 3149 .vb 3150 Mat mat 3151 Scalar mat_array(1) 3152 PetscOffset i_mat 3153 int ierr 3154 call MatGetArray(mat,mat_array,i_mat,ierr) 3155 3156 C Access first local entry in matrix; note that array is 3157 C treated as one dimensional 3158 value = mat_array(i_mat + 1) 3159 3160 [... other code ...] 3161 call MatRestoreArray(mat,mat_array,i_mat,ierr) 3162 .ve 3163 3164 See the Fortran chapter of the users manual and 3165 petsc/src/mat/examples/tests for details 3166 3167 Level: advanced 3168 3169 .keywords: matrix, array, elements, values, restore 3170 3171 .seealso: MatGetArray(), MatRestoreArrayF90() 3172 @*/ 3173 int MatRestoreArray(Mat mat,Scalar **v) 3174 { 3175 int ierr; 3176 3177 PetscFunctionBegin; 3178 PetscValidHeaderSpecific(mat,MAT_COOKIE); 3179 PetscValidPointer(v); 3180 if (!mat->ops->restorearray) SETERRQ(PETSC_ERR_SUP,0,""); 3181 ierr = (*mat->ops->restorearray)(mat,v);CHKERRQ(ierr); 3182 PetscFunctionReturn(0); 3183 } 3184 3185 #undef __FUNC__ 3186 #define __FUNC__ "MatGetSubMatrices" 3187 /*@C 3188 MatGetSubMatrices - Extracts several submatrices from a matrix. If submat 3189 points to an array of valid matrices, they may be reused to store the new 3190 submatrices. 3191 3192 Collective on Mat 3193 3194 Input Parameters: 3195 + mat - the matrix 3196 . n - the number of submatrixes to be extracted 3197 . irow, icol - index sets of rows and columns to extract 3198 - scall - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX 3199 3200 Output Parameter: 3201 . submat - the array of submatrices 3202 3203 Notes: 3204 MatGetSubMatrices() can extract only sequential submatrices 3205 (from both sequential and parallel matrices). Use MatGetSubMatrix() 3206 to extract a parallel submatrix. 3207 3208 When extracting submatrices from a parallel matrix, each processor can 3209 form a different submatrix by setting the rows and columns of its 3210 individual index sets according to the local submatrix desired. 3211 3212 When finished using the submatrices, the user should destroy 3213 them with MatDestroySubMatrices(). 3214 3215 MAT_REUSE_MATRIX can only be used when the nonzero structure of the 3216 original matrix has not changed from that last call to MatGetSubMatrices(). 3217 3218 Fortran Note: 3219 The Fortran interface is slightly different from that given below; it 3220 requires one to pass in as submat a Mat (integer) array of size at least m. 3221 3222 Level: advanced 3223 3224 .keywords: matrix, get, submatrix, submatrices 3225 3226 .seealso: MatDestroyMatrices(), MatGetSubMatrix(), MatGetRow(), MatGetDiagonal() 3227 @*/ 3228 int MatGetSubMatrices(Mat mat,int n,IS *irow,IS *icol,MatReuse scall,Mat **submat) 3229 { 3230 int ierr; 3231 3232 PetscFunctionBegin; 3233 PetscValidHeaderSpecific(mat,MAT_COOKIE); 3234 if (!mat->ops->getsubmatrices) SETERRQ(PETSC_ERR_SUP,0,""); 3235 if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix"); 3236 3237 PLogEventBegin(MAT_GetSubMatrices,mat,0,0,0); 3238 ierr = (*mat->ops->getsubmatrices)(mat,n,irow,icol,scall,submat);CHKERRQ(ierr); 3239 PLogEventEnd(MAT_GetSubMatrices,mat,0,0,0); 3240 3241 PetscFunctionReturn(0); 3242 } 3243 3244 #undef __FUNC__ 3245 #define __FUNC__ "MatDestroyMatrices" 3246 /*@C 3247 MatDestroyMatrices - Destroys a set of matrices obtained with MatGetSubMatrices(). 3248 3249 Collective on Mat 3250 3251 Input Parameters: 3252 + n - the number of local matrices 3253 - mat - the matrices 3254 3255 Level: advanced 3256 3257 .keywords: matrix, destroy, submatrix, submatrices 3258 3259 .seealso: MatGetSubMatrices() 3260 @*/ 3261 int MatDestroyMatrices(int n,Mat **mat) 3262 { 3263 int ierr,i; 3264 3265 PetscFunctionBegin; 3266 if (n < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,1,"Trying to destroy negative number of matrices %d",n); 3267 PetscValidPointer(mat); 3268 for ( i=0; i<n; i++ ) { 3269 ierr = MatDestroy((*mat)[i]);CHKERRQ(ierr); 3270 } 3271 if (n) {ierr = PetscFree(*mat);CHKERRQ(ierr);} 3272 PetscFunctionReturn(0); 3273 } 3274 3275 #undef __FUNC__ 3276 #define __FUNC__ "MatIncreaseOverlap" 3277 /*@ 3278 MatIncreaseOverlap - Given a set of submatrices indicated by index sets, 3279 replaces the index sets by larger ones that represent submatrices with 3280 additional overlap. 3281 3282 Collective on Mat 3283 3284 Input Parameters: 3285 + mat - the matrix 3286 . n - the number of index sets 3287 . is - the array of pointers to index sets 3288 - ov - the additional overlap requested 3289 3290 Level: developer 3291 3292 .keywords: matrix, overlap, Schwarz 3293 3294 .seealso: MatGetSubMatrices() 3295 @*/ 3296 int MatIncreaseOverlap(Mat mat,int n, IS *is,int ov) 3297 { 3298 int ierr; 3299 3300 PetscFunctionBegin; 3301 PetscValidHeaderSpecific(mat,MAT_COOKIE); 3302 if (!mat->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for unassembled matrix"); 3303 if (mat->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix"); 3304 3305 if (ov == 0) PetscFunctionReturn(0); 3306 if (!mat->ops->increaseoverlap) SETERRQ(PETSC_ERR_SUP,0,""); 3307 PLogEventBegin(MAT_IncreaseOverlap,mat,0,0,0); 3308 ierr = (*mat->ops->increaseoverlap)(mat,n,is,ov);CHKERRQ(ierr); 3309 PLogEventEnd(MAT_IncreaseOverlap,mat,0,0,0); 3310 PetscFunctionReturn(0); 3311 } 3312 3313 #undef __FUNC__ 3314 #define __FUNC__ "MatPrintHelp" 3315 /*@ 3316 MatPrintHelp - Prints all the options for the matrix. 3317 3318 Collective on Mat 3319 3320 Input Parameter: 3321 . mat - the matrix 3322 3323 Options Database Keys: 3324 + -help - Prints matrix options 3325 - -h - Prints matrix options 3326 3327 Level: developer 3328 3329 .keywords: mat, help 3330 3331 .seealso: MatCreate(), MatCreateXXX() 3332 @*/ 3333 int MatPrintHelp(Mat mat) 3334 { 3335 static int called = 0; 3336 int ierr; 3337 MPI_Comm comm; 3338 3339 PetscFunctionBegin; 3340 PetscValidHeaderSpecific(mat,MAT_COOKIE); 3341 3342 comm = mat->comm; 3343 if (!called) { 3344 ierr = (*PetscHelpPrintf)(comm,"General matrix options:\n");CHKERRQ(ierr); 3345 ierr = (*PetscHelpPrintf)(comm," -mat_view_info: view basic matrix info during MatAssemblyEnd()\n");CHKERRQ(ierr); 3346 ierr = (*PetscHelpPrintf)(comm," -mat_view_info_detailed: view detailed matrix info during MatAssemblyEnd()\n");CHKERRQ(ierr); 3347 ierr = (*PetscHelpPrintf)(comm," -mat_view_draw: draw nonzero matrix structure during MatAssemblyEnd()\n");CHKERRQ(ierr); 3348 ierr = (*PetscHelpPrintf)(comm," -draw_pause <sec>: set seconds of display pause\n");CHKERRQ(ierr); 3349 ierr = (*PetscHelpPrintf)(comm," -display <name>: set alternate display\n");CHKERRQ(ierr); 3350 called = 1; 3351 } 3352 if (mat->ops->printhelp) { 3353 ierr = (*mat->ops->printhelp)(mat);CHKERRQ(ierr); 3354 } 3355 PetscFunctionReturn(0); 3356 } 3357 3358 #undef __FUNC__ 3359 #define __FUNC__ "MatGetBlockSize" 3360 /*@ 3361 MatGetBlockSize - Returns the matrix block size; useful especially for the 3362 block row and block diagonal formats. 3363 3364 Not Collective 3365 3366 Input Parameter: 3367 . mat - the matrix 3368 3369 Output Parameter: 3370 . bs - block size 3371 3372 Notes: 3373 Block diagonal formats are MATSEQBDIAG, MATMPIBDIAG. 3374 Block row formats are MATSEQBAIJ, MATMPIBAIJ 3375 3376 Level: intermediate 3377 3378 .keywords: matrix, get, block, size 3379 3380 .seealso: MatCreateSeqBAIJ(), MatCreateMPIBAIJ(), MatCreateSeqBDiag(), MatCreateMPIBDiag() 3381 @*/ 3382 int MatGetBlockSize(Mat mat,int *bs) 3383 { 3384 int ierr; 3385 3386 PetscFunctionBegin; 3387 PetscValidHeaderSpecific(mat,MAT_COOKIE); 3388 PetscValidIntPointer(bs); 3389 if (!mat->ops->getblocksize) SETERRQ(PETSC_ERR_SUP,0,""); 3390 ierr = (*mat->ops->getblocksize)(mat,bs);CHKERRQ(ierr); 3391 PetscFunctionReturn(0); 3392 } 3393 3394 #undef __FUNC__ 3395 #define __FUNC__ "MatGetRowIJ" 3396 /*@C 3397 MatGetRowIJ - Returns the compressed row storage i and j indices for sequential matrices. 3398 3399 Collective on Mat 3400 3401 Input Parameters: 3402 + mat - the matrix 3403 . shift - 0 or 1 indicating we want the indices starting at 0 or 1 3404 - symmetric - PETSC_TRUE or PETSC_FALSE indicating the matrix data structure should be 3405 symmetrized 3406 3407 Output Parameters: 3408 + n - number of rows in the (possibly compressed) matrix 3409 . ia - the row pointers 3410 . ja - the column indices 3411 - done - PETSC_TRUE or PETSC_FALSE, indicating whether the values have been returned 3412 3413 Level: developer 3414 3415 .seealso: MatGetColumnIJ(), MatRestoreRowIJ() 3416 @*/ 3417 int MatGetRowIJ(Mat mat,int shift,PetscTruth symmetric,int *n,int **ia,int** ja,PetscTruth *done) 3418 { 3419 int ierr; 3420 3421 PetscFunctionBegin; 3422 PetscValidHeaderSpecific(mat,MAT_COOKIE); 3423 if (ia) PetscValidIntPointer(ia); 3424 if (ja) PetscValidIntPointer(ja); 3425 PetscValidIntPointer(done); 3426 if (!mat->ops->getrowij) *done = PETSC_FALSE; 3427 else { 3428 *done = PETSC_TRUE; 3429 ierr = (*mat->ops->getrowij)(mat,shift,symmetric,n,ia,ja,done);CHKERRQ(ierr); 3430 } 3431 PetscFunctionReturn(0); 3432 } 3433 3434 #undef __FUNC__ 3435 #define __FUNC__ "MatGetColumnIJ" 3436 /*@C 3437 MatGetColumnIJ - Returns the compressed column storage i and j indices for sequential matrices. 3438 3439 Collective on Mat 3440 3441 Input Parameters: 3442 + mat - the matrix 3443 . shift - 1 or zero indicating we want the indices starting at 0 or 1 3444 - symmetric - PETSC_TRUE or PETSC_FALSE indicating the matrix data structure should be 3445 symmetrized 3446 3447 Output Parameters: 3448 + n - number of columns in the (possibly compressed) matrix 3449 . ia - the column pointers 3450 . ja - the row indices 3451 - done - PETSC_TRUE or PETSC_FALSE, indicating whether the values have been returned 3452 3453 Level: developer 3454 3455 .seealso: MatGetRowIJ(), MatRestoreColumnIJ() 3456 @*/ 3457 int MatGetColumnIJ(Mat mat,int shift,PetscTruth symmetric,int *n,int **ia,int** ja,PetscTruth *done) 3458 { 3459 int ierr; 3460 3461 PetscFunctionBegin; 3462 PetscValidHeaderSpecific(mat,MAT_COOKIE); 3463 if (ia) PetscValidIntPointer(ia); 3464 if (ja) PetscValidIntPointer(ja); 3465 PetscValidIntPointer(done); 3466 3467 if (!mat->ops->getcolumnij) *done = PETSC_FALSE; 3468 else { 3469 *done = PETSC_TRUE; 3470 ierr = (*mat->ops->getcolumnij)(mat,shift,symmetric,n,ia,ja,done);CHKERRQ(ierr); 3471 } 3472 PetscFunctionReturn(0); 3473 } 3474 3475 #undef __FUNC__ 3476 #define __FUNC__ "MatRestoreRowIJ" 3477 /*@C 3478 MatRestoreRowIJ - Call after you are completed with the ia,ja indices obtained with 3479 MatGetRowIJ(). 3480 3481 Collective on Mat 3482 3483 Input Parameters: 3484 + mat - the matrix 3485 . shift - 1 or zero indicating we want the indices starting at 0 or 1 3486 - symmetric - PETSC_TRUE or PETSC_FALSE indicating the matrix data structure should be 3487 symmetrized 3488 3489 Output Parameters: 3490 + n - size of (possibly compressed) matrix 3491 . ia - the row pointers 3492 . ja - the column indices 3493 - done - PETSC_TRUE or PETSC_FALSE indicated that the values have been returned 3494 3495 Level: developer 3496 3497 .seealso: MatGetRowIJ(), MatRestoreColumnIJ() 3498 @*/ 3499 int MatRestoreRowIJ(Mat mat,int shift,PetscTruth symmetric,int *n,int **ia,int** ja,PetscTruth *done) 3500 { 3501 int ierr; 3502 3503 PetscFunctionBegin; 3504 PetscValidHeaderSpecific(mat,MAT_COOKIE); 3505 if (ia) PetscValidIntPointer(ia); 3506 if (ja) PetscValidIntPointer(ja); 3507 PetscValidIntPointer(done); 3508 3509 if (!mat->ops->restorerowij) *done = PETSC_FALSE; 3510 else { 3511 *done = PETSC_TRUE; 3512 ierr = (*mat->ops->restorerowij)(mat,shift,symmetric,n,ia,ja,done);CHKERRQ(ierr); 3513 } 3514 PetscFunctionReturn(0); 3515 } 3516 3517 #undef __FUNC__ 3518 #define __FUNC__ "MatRestoreColumnIJ" 3519 /*@C 3520 MatRestoreColumnIJ - Call after you are completed with the ia,ja indices obtained with 3521 MatGetColumnIJ(). 3522 3523 Collective on Mat 3524 3525 Input Parameters: 3526 + mat - the matrix 3527 . shift - 1 or zero indicating we want the indices starting at 0 or 1 3528 - symmetric - PETSC_TRUE or PETSC_FALSE indicating the matrix data structure should be 3529 symmetrized 3530 3531 Output Parameters: 3532 + n - size of (possibly compressed) matrix 3533 . ia - the column pointers 3534 . ja - the row indices 3535 - done - PETSC_TRUE or PETSC_FALSE indicated that the values have been returned 3536 3537 Level: developer 3538 3539 .seealso: MatGetColumnIJ(), MatRestoreRowIJ() 3540 @*/ 3541 int MatRestoreColumnIJ(Mat mat,int shift,PetscTruth symmetric,int *n,int **ia,int** ja,PetscTruth *done) 3542 { 3543 int ierr; 3544 3545 PetscFunctionBegin; 3546 PetscValidHeaderSpecific(mat,MAT_COOKIE); 3547 if (ia) PetscValidIntPointer(ia); 3548 if (ja) PetscValidIntPointer(ja); 3549 PetscValidIntPointer(done); 3550 3551 if (!mat->ops->restorecolumnij) *done = PETSC_FALSE; 3552 else { 3553 *done = PETSC_TRUE; 3554 ierr = (*mat->ops->restorecolumnij)(mat,shift,symmetric,n,ia,ja,done);CHKERRQ(ierr); 3555 } 3556 PetscFunctionReturn(0); 3557 } 3558 3559 #undef __FUNC__ 3560 #define __FUNC__ "MatColoringPatch" 3561 /*@C 3562 MatColoringPatch -Used inside matrix coloring routines that 3563 use MatGetRowIJ() and/or MatGetColumnIJ(). 3564 3565 Collective on Mat 3566 3567 Input Parameters: 3568 + mat - the matrix 3569 . n - number of colors 3570 - colorarray - array indicating color for each column 3571 3572 Output Parameters: 3573 . iscoloring - coloring generated using colorarray information 3574 3575 Level: developer 3576 3577 .seealso: MatGetRowIJ(), MatGetColumnIJ() 3578 3579 .keywords: mat, coloring, patch 3580 @*/ 3581 int MatColoringPatch(Mat mat,int n,int *colorarray,ISColoring *iscoloring) 3582 { 3583 int ierr; 3584 3585 PetscFunctionBegin; 3586 PetscValidHeaderSpecific(mat,MAT_COOKIE); 3587 PetscValidIntPointer(colorarray); 3588 3589 if (!mat->ops->coloringpatch) {SETERRQ(PETSC_ERR_SUP,0,"");} 3590 else { 3591 ierr = (*mat->ops->coloringpatch)(mat,n,colorarray,iscoloring);CHKERRQ(ierr); 3592 } 3593 PetscFunctionReturn(0); 3594 } 3595 3596 3597 #undef __FUNC__ 3598 #define __FUNC__ "MatSetUnfactored" 3599 /*@ 3600 MatSetUnfactored - Resets a factored matrix to be treated as unfactored. 3601 3602 Collective on Mat 3603 3604 Input Parameter: 3605 . mat - the factored matrix to be reset 3606 3607 Notes: 3608 This routine should be used only with factored matrices formed by in-place 3609 factorization via ILU(0) (or by in-place LU factorization for the MATSEQDENSE 3610 format). This option can save memory, for example, when solving nonlinear 3611 systems with a matrix-free Newton-Krylov method and a matrix-based, in-place 3612 ILU(0) preconditioner. 3613 3614 Note that one can specify in-place ILU(0) factorization by calling 3615 .vb 3616 PCType(pc,PCILU); 3617 PCILUSeUseInPlace(pc); 3618 .ve 3619 or by using the options -pc_type ilu -pc_ilu_in_place 3620 3621 In-place factorization ILU(0) can also be used as a local 3622 solver for the blocks within the block Jacobi or additive Schwarz 3623 methods (runtime option: -sub_pc_ilu_in_place). See the discussion 3624 of these preconditioners in the users manual for details on setting 3625 local solver options. 3626 3627 Most users should employ the simplified SLES interface for linear solvers 3628 instead of working directly with matrix algebra routines such as this. 3629 See, e.g., SLESCreate(). 3630 3631 Level: developer 3632 3633 .seealso: PCILUSetUseInPlace(), PCLUSetUseInPlace() 3634 3635 .keywords: matrix-free, in-place ILU, in-place LU 3636 @*/ 3637 int MatSetUnfactored(Mat mat) 3638 { 3639 int ierr; 3640 3641 PetscFunctionBegin; 3642 PetscValidHeaderSpecific(mat,MAT_COOKIE); 3643 mat->factor = 0; 3644 if (!mat->ops->setunfactored) PetscFunctionReturn(0); 3645 ierr = (*mat->ops->setunfactored)(mat);CHKERRQ(ierr); 3646 PetscFunctionReturn(0); 3647 } 3648 3649 #undef __FUNC__ 3650 #define __FUNC__ "MatGetType" 3651 /*@C 3652 MatGetType - Gets the matrix type and name (as a string) from the matrix. 3653 3654 Not Collective 3655 3656 Input Parameter: 3657 . mat - the matrix 3658 3659 Output Parameter: 3660 + type - the matrix type (or use PETSC_NULL) 3661 - name - name of matrix type (or use PETSC_NULL) 3662 3663 Level: intermediate 3664 3665 .keywords: matrix, get, type, name 3666 @*/ 3667 int MatGetType(Mat mat,MatType *type,char **name) 3668 { 3669 int itype = (int)mat->type; 3670 char *matname[10]; 3671 3672 PetscFunctionBegin; 3673 PetscValidHeaderSpecific(mat,MAT_COOKIE); 3674 3675 if (type) *type = (MatType) mat->type; 3676 if (name) { 3677 /* Note: Be sure that this list corresponds to the enum in mat.h */ 3678 matname[0] = "MATSEQDENSE"; 3679 matname[1] = "MATSEQAIJ"; 3680 matname[2] = "MATMPIAIJ"; 3681 matname[3] = "MATSHELL"; 3682 matname[4] = "MATMPIROWBS"; 3683 matname[5] = "MATSEQBDIAG"; 3684 matname[6] = "MATMPIBDIAG"; 3685 matname[7] = "MATMPIDENSE"; 3686 matname[8] = "MATSEQBAIJ"; 3687 matname[9] = "MATMPIBAIJ"; 3688 3689 if (itype < 0 || itype > 9) *name = "Unknown matrix type"; 3690 else *name = matname[itype]; 3691 } 3692 PetscFunctionReturn(0); 3693 } 3694 3695 /*MC 3696 MatGetArrayF90 - Accesses a matrix array from Fortran90. 3697 3698 Synopsis: 3699 MatGetArrayF90(Mat x,{Scalar, pointer :: xx_v(:)},integer ierr) 3700 3701 Not collective 3702 3703 Input Parameter: 3704 . x - matrix 3705 3706 Output Parameters: 3707 + xx_v - the Fortran90 pointer to the array 3708 - ierr - error code 3709 3710 Example of Usage: 3711 .vb 3712 Scalar, pointer xx_v(:) 3713 .... 3714 call MatGetArrayF90(x,xx_v,ierr) 3715 a = xx_v(3) 3716 call MatRestoreArrayF90(x,xx_v,ierr) 3717 .ve 3718 3719 Notes: 3720 Not yet supported for all F90 compilers 3721 3722 Level: advanced 3723 3724 .seealso: MatRestoreArrayF90(), MatGetArray(), MatRestoreArray() 3725 3726 .keywords: matrix, array, f90 3727 M*/ 3728 3729 /*MC 3730 MatRestoreArrayF90 - Restores a matrix array that has been 3731 accessed with MatGetArrayF90(). 3732 3733 Synopsis: 3734 MatRestoreArrayF90(Mat x,{Scalar, pointer :: xx_v(:)},integer ierr) 3735 3736 Not collective 3737 3738 Input Parameters: 3739 + x - matrix 3740 - xx_v - the Fortran90 pointer to the array 3741 3742 Output Parameter: 3743 . ierr - error code 3744 3745 Example of Usage: 3746 .vb 3747 Scalar, pointer xx_v(:) 3748 .... 3749 call MatGetArrayF90(x,xx_v,ierr) 3750 a = xx_v(3) 3751 call MatRestoreArrayF90(x,xx_v,ierr) 3752 .ve 3753 3754 Notes: 3755 Not yet supported for all F90 compilers 3756 3757 Level: advanced 3758 3759 .seealso: MatGetArrayF90(), MatGetArray(), MatRestoreArray() 3760 3761 .keywords: matrix, array, f90 3762 M*/ 3763 3764 3765 #undef __FUNC__ 3766 #define __FUNC__ "MatGetSubMatrix" 3767 /*@ 3768 MatGetSubMatrix - Gets a single submatrix on the same number of processors 3769 as the original matrix. 3770 3771 Collective on Mat 3772 3773 Input Parameters: 3774 + mat - the original matrix 3775 . isrow - rows this processor should obtain 3776 . iscol - columns for all processors you wish to keep 3777 . csize - number of columns "local" to this processor (does nothing for sequential 3778 matrices). This should match the result from VecGetLocalSize(x,...) if you 3779 plan to use the matrix in a A*x; alternatively, you can use PETSC_DECIDE 3780 - cll - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX 3781 3782 Output Parameter: 3783 . newmat - the new submatrix, of the same type as the old 3784 3785 Level: advanced 3786 3787 .keywords: matrix, get, submatrix, submatrices 3788 3789 .seealso: MatGetSubMatrices() 3790 @*/ 3791 int MatGetSubMatrix(Mat mat,IS isrow,IS iscol,int csize,MatReuse cll,Mat *newmat) 3792 { 3793 int ierr, size; 3794 Mat *local; 3795 3796 PetscFunctionBegin; 3797 ierr = MPI_Comm_size(mat->comm,&size);CHKERRQ(ierr); 3798 3799 /* if original matrix is on just one processor then use submatrix generated */ 3800 if (!mat->ops->getsubmatrix && size == 1 && cll == MAT_REUSE_MATRIX) { 3801 ierr = MatGetSubMatrices(mat,1,&isrow,&iscol,MAT_REUSE_MATRIX,&newmat);CHKERRQ(ierr); 3802 PetscFunctionReturn(0); 3803 } else if (!mat->ops->getsubmatrix && size == 1) { 3804 ierr = MatGetSubMatrices(mat,1,&isrow,&iscol,MAT_INITIAL_MATRIX,&local);CHKERRQ(ierr); 3805 *newmat = *local; 3806 ierr = PetscFree(local);CHKERRQ(ierr); 3807 PetscFunctionReturn(0); 3808 } 3809 3810 if (!mat->ops->getsubmatrix) SETERRQ(PETSC_ERR_SUP,0,"Not currently implemented"); 3811 ierr = (*mat->ops->getsubmatrix)(mat,isrow,iscol,csize,cll,newmat);CHKERRQ(ierr); 3812 PetscFunctionReturn(0); 3813 } 3814 3815 #undef __FUNC__ 3816 #define __FUNC__ "MatGetMaps" 3817 /*@C 3818 MatGetMaps - Returns the maps associated with the matrix. 3819 3820 Not Collective 3821 3822 Input Parameter: 3823 . mat - the matrix 3824 3825 Output Parameters: 3826 + rmap - the row (right) map 3827 - cmap - the column (left) map 3828 3829 Level: developer 3830 3831 .keywords: matrix, get, map 3832 @*/ 3833 int MatGetMaps(Mat mat,Map *rmap,Map *cmap) 3834 { 3835 int ierr; 3836 3837 PetscFunctionBegin; 3838 PetscValidHeaderSpecific(mat,MAT_COOKIE); 3839 ierr = (*mat->ops->getmaps)(mat,rmap,cmap);CHKERRQ(ierr); 3840 PetscFunctionReturn(0); 3841 } 3842 3843 /* 3844 Version that works for all PETSc matrices 3845 */ 3846 #undef __FUNC__ 3847 #define __FUNC__ "MatGetMaps_Petsc" 3848 int MatGetMaps_Petsc(Mat mat,Map *rmap,Map *cmap) 3849 { 3850 PetscFunctionBegin; 3851 if (rmap) *rmap = mat->rmap; 3852 if (cmap) *cmap = mat->cmap; 3853 PetscFunctionReturn(0); 3854 } 3855 3856 #undef __FUNC__ 3857 #define __FUNC__ "MatSetStashInitialSize" 3858 /*@ 3859 MatSetStashInitialSize - sets the sizes of the matrix stash, that is 3860 used during the assembly process to store values that belong to 3861 other processors. 3862 3863 Not Collective 3864 3865 Input Parameters: 3866 + mat - the matrix 3867 . size - the initial size of the stash. 3868 - bsize - the initial size of the block-stash(if used). 3869 3870 Options Database Keys: 3871 + -matstash_initial_size <size> or <size0,size1,...sizep-1> 3872 - -matstash_block_initial_size <bsize> or <bsize0,bsize1,...bsizep-1> 3873 3874 Level: intermediate 3875 3876 Notes: 3877 The block-stash is used for values set with VecSetValuesBlocked() while 3878 the stash is used for values set with VecSetValues() 3879 3880 Run with the option -log_info and look for output of the form 3881 MatAssemblyBegin_MPIXXX:Stash has MM entries, uses nn mallocs. 3882 to determine the appropriate value, MM, to use for size and 3883 MatAssemblyBegin_MPIXXX:Block-Stash has BMM entries, uses nn mallocs. 3884 to determine the value, BMM to use for bsize 3885 3886 .keywords: matrix, stash, assembly 3887 @*/ 3888 int MatSetStashInitialSize(Mat mat,int size, int bsize) 3889 { 3890 int ierr; 3891 3892 PetscFunctionBegin; 3893 PetscValidHeaderSpecific(mat,MAT_COOKIE); 3894 ierr = MatStashSetInitialSize_Private(&mat->stash,size);CHKERRQ(ierr); 3895 ierr = MatStashSetInitialSize_Private(&mat->bstash,bsize);CHKERRQ(ierr); 3896 PetscFunctionReturn(0); 3897 } 3898