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