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