1 2 /* 3 This is where the abstract matrix operations are defined 4 */ 5 6 #include "petsc.h" 7 #include "matimpl.h" /*I "mat.h" I*/ 8 #include "vec/vecimpl.h" 9 10 /*@ 11 MatGetReordering - 12 13 Input Parameters: 14 . mat - the matrix 15 . type - type of reordering; one of:ORDER_ND,ORDER_1WD,ORDER_RCM,ORDER_QMD 16 17 OutPut Parameters: 18 . rperm - row permutation indices 19 . cperm - column permutation indices 20 21 If the column permutations and row permutations are the same, then 22 returns 0 in cperm. 23 @*/ 24 int MatGetReordering(Mat mat,int type,IS *rperm,IS *cperm) 25 { 26 VALIDHEADER(mat,MAT_COOKIE); 27 if (!mat->ops->order) SETERR(1,"Cannot reorder this matrix"); 28 return (*mat->ops->order)(mat,type,rperm,cperm); 29 } 30 31 /*@ 32 MatGetRow - gets a row of a matrix. This routine is provided 33 for people who need to have direct address to the 34 structure of a matrix, we hope that we provide 35 enough high level matrix routines so that few users 36 need it. You MUST call MatRestoreRow() for each row 37 that you get to insure that your application does 38 not bleed memory. 39 40 Input Parameters: 41 . mat - the matrix 42 . row - the row to get 43 44 OutPut Parameters: 45 . ncols, cols - the number of nonzeros and their columns 46 . vals - if nonzero the column values 47 @*/ 48 int MatGetRow(Mat mat,int row,int *ncols,int **cols,Scalar **vals) 49 { 50 VALIDHEADER(mat,MAT_COOKIE); 51 return (*mat->ops->getrow)(mat,row,ncols,cols,vals); 52 } 53 54 /*@ 55 MatRestoreRow - frees any temporary space malloced by 56 MatGetRow(). 57 58 Input Parameters: 59 . mat - the matrix 60 . row - the row to get 61 . ncols, cols - the number of nonzeros and their columns 62 . vals - if nonzero the column values 63 @*/ 64 int MatRestoreRow(Mat mat,int row,int *ncols,int **cols,Scalar **vals) 65 { 66 VALIDHEADER(mat,MAT_COOKIE); 67 if (!mat->ops->restorerow) return 0; 68 return (*mat->ops->restorerow)(mat,row,ncols,cols,vals); 69 } 70 /*@ 71 MatView - visualize a matrix object. 72 73 Input Parameters: 74 . mat - the matrix 75 . ptr - a visualization context 76 @*/ 77 int MatView(Mat mat,Viewer ptr) 78 { 79 VALIDHEADER(mat,MAT_COOKIE); 80 return (*mat->view)((PetscObject)mat,ptr); 81 } 82 /*@ 83 MatDestroy - frees space taken by matrix 84 85 Input Parameters: 86 . mat the matrix 87 @*/ 88 int MatDestroy(Mat mat) 89 { 90 VALIDHEADER(mat,MAT_COOKIE); 91 return (*mat->destroy)((PetscObject)mat); 92 } 93 /*@ 94 MatValidMatrix - returns 1 if a valid matrix else 0 95 96 Input Parameter: 97 . m - the matrix to check 98 @*/ 99 int MatValidMatrix(Mat m) 100 { 101 if (!m) return 0; 102 if (m->cookie != MAT_COOKIE) return 0; 103 return 1; 104 } 105 106 /*@ 107 MatInsertValues - inserts a block of values into a matrix 108 109 Input Parameters: 110 . mat - the matrix 111 . v - a logically two dimensional array of values, Fortran ordering 112 . m, indexm - the number of rows and their indices 113 . n, indexn - the number of columns and their indices 114 @*/ 115 int MatInsertValues(Mat mat,Scalar *v,int m,int *indexm,int n,int *indexn) 116 { 117 VALIDHEADER(mat,MAT_COOKIE); 118 return (*mat->ops->insert)(mat,v,m,indexm,n,indexn); 119 } 120 /*@ 121 MatAddValues - adds a block of values into a matrix 122 123 Input Parameters: 124 . mat - the matrix 125 . v - a logically two dimensional array of values, Fortran ordering 126 . m, indexm - the number of rows and their indices 127 . n, indexn - the number of columns and their indices 128 @*/ 129 int MatAddValues(Mat mat,Scalar *v,int m,int *indexm,int n,int *indexn) 130 { 131 VALIDHEADER(mat,MAT_COOKIE); 132 return (*mat->ops->insertadd)(mat,v,m,indexm,n,indexn); 133 } 134 /* --------------------------------------------------------*/ 135 /*@ 136 MatMult - Matrix vector multiply 137 138 Input Parameters: 139 . mat -the matrix 140 . x - the vector to be multilplied 141 Output Parameters: 142 . y - the result 143 @*/ 144 int MatMult(Mat mat,Vec x,Vec y) 145 { 146 VALIDHEADER(mat,MAT_COOKIE); 147 VALIDHEADER(x,VEC_COOKIE);VALIDHEADER(y,VEC_COOKIE); 148 return (*mat->ops->mult)(mat,x,y); 149 } 150 /*@ 151 MatMultTrans - Matrix transpose times a vector 152 153 Input Parameters: 154 . mat -the matrix 155 . x - the vector to be multilplied 156 Output Parameters: 157 . y - the result 158 @*/ 159 int MatMultTrans(Mat mat,Vec x,Vec y) 160 { 161 VALIDHEADER(mat,MAT_COOKIE); 162 VALIDHEADER(x,VEC_COOKIE); VALIDHEADER(y,VEC_COOKIE); 163 return (*mat->ops->multtrans)(mat,x,y); 164 } 165 /*@ 166 MatMultAdd - v3 = v2 + A v1 167 168 Input Parameters: 169 . mat -the matrix 170 . v1, v2 - the vectors 171 Output Parameters: 172 . v3 - the result 173 @*/ 174 int MatMultAdd(Mat mat,Vec v1,Vec v2,Vec v3) 175 { 176 VALIDHEADER(mat,MAT_COOKIE);VALIDHEADER(v1,VEC_COOKIE); 177 VALIDHEADER(v2,VEC_COOKIE); VALIDHEADER(v3,VEC_COOKIE); 178 return (*mat->ops->multadd)(mat,v1,v2,v3); 179 } 180 /*@ 181 MatMultTransAdd - v3 = v2 + A' v1 182 183 Input Parameters: 184 . mat -the matrix 185 . v1, v2 - the vectors 186 Output Parameters: 187 . v3 - the result 188 @*/ 189 int MatMultTransAdd(Mat mat,Vec v1,Vec v2,Vec v3) 190 { 191 VALIDHEADER(mat,MAT_COOKIE); VALIDHEADER(v1,VEC_COOKIE); 192 VALIDHEADER(v2,VEC_COOKIE); VALIDHEADER(v3,VEC_COOKIE); 193 return (*mat->ops->multtransadd)(mat,v1,v2,v3); 194 } 195 /* ------------------------------------------------------------*/ 196 /*@ 197 MatNonZeros - returns the number of nonzeros in a matrix 198 199 Input Parameters: 200 . mat - the matrix 201 202 Output Parameters: 203 . returns the number of nonzeros 204 @*/ 205 int MatNonZeros(Mat mat,int *nz) 206 { 207 VALIDHEADER(mat,MAT_COOKIE); 208 return (*mat->ops->NZ)(mat,nz); 209 } 210 /*@ 211 MatMemoryUsed - returns the amount of memory used to store matrix. 212 213 Input Parameters: 214 . mat - the matrix 215 216 Output Parameters: 217 . mem - memory used 218 @*/ 219 int MatMemoryUsed(Mat mat,int *mem) 220 { 221 VALIDHEADER(mat,MAT_COOKIE); 222 return (*mat->ops->memory)(mat,mem); 223 } 224 /* ----------------------------------------------------------*/ 225 /*@ 226 MatLUFactor - performs an inplace LU factorization of matrix. 227 See MatLUFactorSymbolic() and MatLUFactorNumeric() for 228 out-of-place factorization. See MatCholeskyFactor() 229 for symmetric, positive definite case. 230 231 Input Parameters: 232 . mat - the matrix 233 . row, col - row and column permutations 234 235 @*/ 236 int MatLUFactor(Mat mat,IS row,IS col) 237 { 238 VALIDHEADER(mat,MAT_COOKIE); 239 if (mat->ops->lufactor) return (*mat->ops->lufactor)(mat,row,col); 240 SETERR(1,"No MatLUFactor for implementation"); 241 } 242 /*@ 243 MatLUFactorSymbolic - performs a symbolic LU factorization of matrix. 244 See MatLUFactor() for in-place factorization. See 245 MatCholeskyFactor() for symmetric, positive definite case. 246 247 Input Parameters: 248 . mat - the matrix 249 . row, col - row and column permutations 250 251 Output Parameters: 252 . fact - puts factor else does inplace factorization 253 @*/ 254 int MatLUFactorSymbolic(Mat mat,IS row,IS col,Mat *fact) 255 { 256 VALIDHEADER(mat,MAT_COOKIE); 257 if (!fact) SETERR(1,"Missing slot for symbolic factorization"); 258 if (mat->ops->lufactorsymbolic) return (*mat->ops->lufactorsymbolic)(mat,row, 259 col,fact); 260 SETERR(1,"No MatLUFactorSymbolic for implementation"); 261 } 262 /*@ 263 MatLUFactorNumeric - performs a numeric LU factorization of matrix. 264 See MatLUFactor() for in-place factorization. See 265 MatCholeskyFactor() for symmetric, positive definite case. 266 267 Input Parameters: 268 . mat - the matrix 269 . row, col - row and column permutations 270 271 Output Parameters: 272 . fact - must have been obtained with MatLUFactorSymbolic() 273 @*/ 274 int MatLUFactorNumeric(Mat mat,Mat fact) 275 { 276 VALIDHEADER(mat,MAT_COOKIE); 277 if (!fact) SETERR(1,"Missing symbolic factorization"); 278 if (mat->ops->lufactornumeric) return (*mat->ops->lufactornumeric)(mat,fact); 279 SETERR(1,"No MatLUFactorNumeric for implementation"); 280 } 281 /*@ 282 MatCholeskyFactor - performs an inplace CC' factorization of matrix. 283 See also MatLUFactor(), MatCholeskyFactorSymbolic() and 284 MatCholeskyFactorNumeric(). 285 286 Input Parameters: 287 . mat - the matrix 288 . perm - row and column permutations 289 290 @*/ 291 int MatCholeskyFactor(Mat mat,IS perm) 292 { 293 VALIDHEADER(mat,MAT_COOKIE); 294 if (mat->ops->chfactor) return (*mat->ops->chfactor)(mat,perm); 295 SETERR(1,"No MatCholeskyFactor for implementation"); 296 } 297 /*@ 298 MatCholeskyFactorSymbolic - performs a symbolic Cholesky factorization. 299 See also MatLUFactor(), MatCholeskyFactorSymbolic() and 300 MatCholeskyFactorNumeric(). 301 302 Input Parameters: 303 . mat - the matrix 304 . perm - row and column permutations 305 306 @*/ 307 int MatCholeskyFactorSymbolic(Mat mat,IS perm,Mat *fact) 308 { 309 VALIDHEADER(mat,MAT_COOKIE); 310 if (!fact) SETERR(1,"Missing slot for symbolic factorization"); 311 if (mat->ops->chfactorsymbolic) return (*mat->ops->chfactorsymbolic)(mat, 312 perm,fact); 313 SETERR(1,"No MatCholeskyFactorSymbolic for implementation"); 314 } 315 /*@ 316 MatCholeskyFactorNumeric - performs a numeric Cholesky factorization. 317 See also MatLUFactor(), MatCholeskyFactor() and 318 MatCholeskyFactorSymbolic(). 319 320 Input Parameters: 321 . mat - the matrix 322 323 324 @*/ 325 int MatCholeskyFactorNumeric(Mat mat,Mat fact) 326 { 327 VALIDHEADER(mat,MAT_COOKIE); 328 if (!fact) SETERR(1,"Missing symbolic factorization"); 329 if (mat->ops->chfactornumeric) return (*mat->ops->chfactornumeric)(mat, 330 fact); 331 SETERR(1,"No MatCholeskyFactorNumeric for implementation"); 332 } 333 /* ----------------------------------------------------------------*/ 334 /*@ 335 MatSolve - solve A x = b 336 337 Input Parameters: 338 . mat - the factored matrix 339 . b - the right hand side 340 341 Output Parameter: 342 . x- the result 343 @*/ 344 int MatSolve(Mat mat,Vec b,Vec x) 345 { 346 VALIDHEADER(mat,MAT_COOKIE); 347 VALIDHEADER(b,VEC_COOKIE); VALIDHEADER(x,VEC_COOKIE); 348 if (!mat->factor) SETERR(1,"Attempt solve on nonfactored matrix"); 349 if (mat->ops->solve) return (*mat->ops->solve)(mat,b,x); 350 SETERR(1,"No MatSolve for implementation"); 351 } 352 /*@ 353 MatSolveAdd - x = y + A^-1 b 354 355 Input Parameters: 356 . mat - the factored matrix 357 . b - the right hand side 358 . y - te vector to be added to 359 360 Output Parameter: 361 . x- the result 362 @*/ 363 int MatSolveAdd(Mat mat,Vec b,Vec y,Vec x) 364 { 365 Scalar one = 1.0; 366 Vec tmp; 367 int ierr; 368 VALIDHEADER(mat,MAT_COOKIE);VALIDHEADER(y,VEC_COOKIE); 369 VALIDHEADER(b,VEC_COOKIE); VALIDHEADER(x,VEC_COOKIE); 370 if (!mat->factor) SETERR(1,"Attempt solve on nonfactored matrix"); 371 if (mat->ops->solveadd) return (*mat->ops->solveadd)(mat,b,y,x); 372 /* do the solve then the add manually */ 373 if (x != y) { 374 ierr = MatSolve(mat,b,x); CHKERR(ierr); 375 ierr = VecAXPY(&one,y,x); CHKERR(ierr); 376 return 0; 377 } 378 else { 379 ierr = VecCreate(x,&tmp); CHKERR(ierr); 380 ierr = VecCopy(x,tmp); CHKERR(ierr); 381 ierr = MatSolve(mat,b,x); CHKERR(ierr); 382 ierr = VecAXPY(&one,tmp,x); CHKERR(ierr); 383 ierr = VecDestroy(tmp); CHKERR(ierr); 384 return 0; 385 } 386 } 387 /*@ 388 MatSolveTrans - solve A' x = b 389 390 Input Parameters: 391 . mat - the factored matrix 392 . b - the right hand side 393 394 Output Parameter: 395 . x- the result 396 @*/ 397 int MatSolveTrans(Mat mat,Vec b,Vec x) 398 { 399 VALIDHEADER(mat,MAT_COOKIE); 400 VALIDHEADER(b,VEC_COOKIE); VALIDHEADER(x,VEC_COOKIE); 401 if (!mat->factor) SETERR(1,"Attempt solve on nonfactored matrix"); 402 if (mat->ops->solvetrans) return (*mat->ops->solvetrans)(mat,b,x); 403 SETERR(1,"No MatSolveTrans for implementation"); 404 } 405 /*@ 406 MatSolveTransAdd - x = y + A^-T b 407 408 Input Parameters: 409 . mat - the factored matrix 410 . b - the right hand side 411 . y - te vector to be added to 412 413 Output Parameter: 414 . x- the result 415 @*/ 416 int MatSolveTransAdd(Mat mat,Vec b,Vec y,Vec x) 417 { 418 Scalar one = 1.0; 419 int ierr; 420 Vec tmp; 421 VALIDHEADER(mat,MAT_COOKIE);VALIDHEADER(y,VEC_COOKIE); 422 VALIDHEADER(b,VEC_COOKIE); VALIDHEADER(x,VEC_COOKIE); 423 if (!mat->factor) SETERR(1,"Attempt solve on nonfactored matrix"); 424 if (mat->ops->solvetransadd) return (*mat->ops->solvetransadd)(mat,b,y,x); 425 /* do the solve then the add manually */ 426 if (x != y) { 427 ierr = MatSolveTrans(mat,b,x); CHKERR(ierr); 428 ierr = VecAXPY(&one,y,x); CHKERR(ierr); 429 return 0; 430 } 431 else { 432 ierr = VecCreate(x,&tmp); CHKERR(ierr); 433 ierr = VecCopy(x,tmp); CHKERR(ierr); 434 ierr = MatSolveTrans(mat,b,x); CHKERR(ierr); 435 ierr = VecAXPY(&one,tmp,x); CHKERR(ierr); 436 ierr = VecDestroy(tmp); CHKERR(ierr); 437 return 0; 438 } 439 } 440 /* ----------------------------------------------------------------*/ 441 442 /*@ 443 MatRelax - one relaxation sweep 444 445 Input Parameters: 446 . mat - the matrix 447 . b - the right hand side 448 . omega - the relaxation factor 449 . flag - or together from SOR_FORWARD_SWEEP, SOR_BACKWARD_SWEEP, 450 . SOR_ZERO_INITIAL_GUESS. (SOR_SYMMETRIC_SWEEP is 451 . equivalent to SOR_FORWARD_SWEEP | SOR_BACKWARD_SWEEP) 452 . is - optional index set indicating ordering 453 . its - the number of iterations 454 455 Output Parameters: 456 . x - the solution - can contain initial guess 457 @*/ 458 int MatRelax(Mat mat,Vec b,double omega,int flag,IS is,int its,Vec x) 459 { 460 VALIDHEADER(mat,MAT_COOKIE); 461 VALIDHEADER(b,VEC_COOKIE); VALIDHEADER(x,VEC_COOKIE); 462 if (flag < 1 || flag > 7) SETERR(1,"Bad relaxation type"); 463 if (mat->ops->relax) return (*mat->ops->relax)(mat,b,omega,flag,is,its,x); 464 SETERR(1,"No MatRelax for implementation"); 465 } 466 467 /* ----------------------------------------------------------------*/ 468 /*@ 469 MatCopy - copies a matrix 470 471 Input Parameters: 472 . mat - the matrix 473 474 Output Parameters: 475 . M - pointer to place new matrix 476 @*/ 477 int MatCopy(Mat mat,Mat *M) 478 { 479 VALIDHEADER(mat,MAT_COOKIE); 480 if (mat->ops->copy) return (*mat->ops->copy)(mat,M); 481 SETERR(1,"No MatCopy for implementation"); 482 } 483 484 /*@ 485 MatGetDiagonal - get diagonal of matrix 486 487 Input Parameters: 488 . mat - the matrix 489 490 Output Parameters: 491 . v - the vector to store the diagonal in 492 @*/ 493 int MatGetDiagonal(Mat mat,Vec v) 494 { 495 VALIDHEADER(mat,MAT_COOKIE); 496 VALIDHEADER(v,VEC_COOKIE); 497 if (mat->ops->getdiag) return (*mat->ops->getdiag)(mat,v); 498 SETERR(1,"No MatGetDiagonal for implementaion"); 499 } 500 501 /*@ 502 MatTranspose - in place transpose of matrix 503 504 Input Parameters: 505 . mat - the matrix to transpose 506 @*/ 507 int MatTranspose(Mat mat) 508 { 509 VALIDHEADER(mat,MAT_COOKIE); 510 if (mat->ops->trans) return (*mat->ops->trans)(mat); 511 SETERR(1,"No MatTranspose for implementaion"); 512 } 513 514 /*@ 515 MatEqual - returns 1 if two matrices are equal 516 517 Input Parameters: 518 . mat1, mat2 - the matrices 519 520 OutPut Parameter: 521 . returns 1 if matrices are equal 522 @*/ 523 int MatEqual(Mat mat1,Mat mat2) 524 { 525 VALIDHEADER(mat1,MAT_COOKIE); VALIDHEADER(mat2,MAT_COOKIE); 526 if (mat1->ops->equal) return (*mat1->ops->equal)(mat1,mat2); 527 SETERR(1,"No MatEqual for implementaion"); 528 } 529 530 /*@ 531 MatScale - scale a matrix on the left and right by diagonal 532 matrices stored as vectors. Either of the two may be 533 null. 534 535 Input Parameters: 536 . mat - the matrix to be scaled 537 . l,r - the left and right scaling vectors 538 @*/ 539 int MatScale(Mat mat,Vec l,Vec r) 540 { 541 VALIDHEADER(mat,MAT_COOKIE); 542 VALIDHEADER(l,VEC_COOKIE); VALIDHEADER(r,VEC_COOKIE); 543 if (mat->ops->scale) return (*mat->ops->scale)(mat,l,r); 544 SETERR(1,"No MatScale for implementaion"); 545 } 546 547 /*@ 548 MatNorm - calculate various norms of a matrix 549 550 Input Parameters: 551 . mat - the matrix 552 . type - the type of norm, NORM_1, NORM_2, NORM_FROBENIUS, NORM_INFINITY 553 554 Output Parameters: 555 . norm - the resulting norm 556 @*/ 557 int MatNorm(Mat mat,int type,double *norm) 558 { 559 VALIDHEADER(mat,MAT_COOKIE); 560 if (mat->ops->norm) return (*mat->ops->norm)(mat,type,norm); 561 SETERR(1,"No MatNorm for implementaion"); 562 } 563 564 /*@ 565 MatBeginAssembly - begin assemblying the matrix. This should 566 be called after all the calls to MatInsertValues() 567 and MatAddValues(). 568 569 Input Parameters: 570 . mat - the matrix 571 572 Note: when you call MatInsertValues() it generally caches the values 573 only after you have called MatBeginAssembly() and 574 MatEndAssembly() is the matrix ready to be used. 575 @*/ 576 int MatBeginAssembly(Mat mat) 577 { 578 VALIDHEADER(mat,MAT_COOKIE); 579 if (mat->ops->bassembly) return (*mat->ops->bassembly)(mat); 580 return 0; 581 } 582 /*@ 583 MatEndAssembly - begin assemblying the matrix. This should 584 be called after all the calls to MatInsertValues() 585 and MatAddValues() and after MatBeginAssembly(). 586 587 Input Parameters: 588 . mat - the matrix 589 590 Note: when you call MatInsertValues() it generally caches the values 591 only after you have called MatBeginAssembly() and 592 MatEndAssembly() is the matrix ready to be used. 593 @*/ 594 int MatEndAssembly(Mat mat) 595 { 596 VALIDHEADER(mat,MAT_COOKIE); 597 if (mat->ops->eassembly) return (*mat->ops->eassembly)(mat); 598 return 0; 599 } 600 /*@ 601 MatCompress - trys to store matrix in as little space as 602 possible. May fail if memory is already fully used as 603 it tries to allocate new space. 604 605 Input Parameters: 606 . mat - the matrix 607 608 @*/ 609 int MatCompress(Mat mat) 610 { 611 VALIDHEADER(mat,MAT_COOKIE); 612 if (mat->ops->compress) return (*mat->ops->compress)(mat); 613 return 0; 614 } 615 /*@ 616 MatSetInsertOption - set parameter option on how you will insert 617 (or add) values. In general sorted, row oriented input will assemble 618 fastest. Defaults to row oriented, nonsorted input. Note that 619 if you are using a Fortran 77 module to compute an element 620 stiffness you may need to modify the module or use column 621 oriented input. 622 623 Input Parameters: 624 . mat - the matrix 625 . option - one of ROW_ORIENTED, COLUMN_ORIENTED, ROWS_SORTED, 626 COLUMNS_SORTED, NO_NEW_NONZERO_LOCATIONS 627 628 NO_NEW_NONZERO_LOCATIONS indicates that any add or insertion 629 that will generate a new entry in the nonzero structure is ignored. 630 What this means is if memory is not allocated for this particular 631 lot then the insertion is ignored. For dense matrices where 632 the entire array is allocated no entries are every ignored. This 633 may not be a good idea?? 634 635 @*/ 636 int MatSetInsertOption(Mat mat,int op) 637 { 638 VALIDHEADER(mat,MAT_COOKIE); 639 if (op < 0 || op > 7) SETERR(1,"Invalid option to MatSetInsertOption"); 640 if (mat->ops->insopt) return (*mat->ops->insopt)(mat,op); 641 return 0; 642 } 643 644 /*@ 645 MatZeroEntries - zeros all entries in a matrix. For sparse matrices 646 it keeps the old nonzero structure. 647 648 Input Parameters: 649 . mat - the matrix 650 651 @*/ 652 int MatZeroEntries(Mat mat) 653 { 654 VALIDHEADER(mat,MAT_COOKIE); 655 if (mat->ops->zeroentries) return (*mat->ops->zeroentries)(mat); 656 SETERR(1,"No MatZeroEntries for implementation"); 657 } 658 659 /*@ 660 MatZeroRow - zeros all entries in a row of a matrix. For sparse matrices 661 it removes the old nonzero structure. 662 663 Input Parameters: 664 . mat - the matrix 665 . row - the row to zap 666 667 @*/ 668 int MatZeroRow(int row,Mat mat) 669 { 670 VALIDHEADER(mat,MAT_COOKIE); 671 if (mat->ops->zerorow) return (*mat->ops->zerorow)(row,mat); 672 SETERR(1,"No MatZeroRow for implementation"); 673 } 674 675 /*@ 676 MatScatterBegin - Begins a scatter of one matrix into another. 677 This is much more general than a standard scatter, it include 678 scatter, gather and cmobinations. In a parallel enviroment 679 it can be used to simultaneous gather from a global matrix 680 into local and vice-versa. 681 682 Input Parameters: 683 . mat - the matrix 684 . is1,is2 - the rows and columns to snag 685 . is3,is4 - rows and columns to drop snagged in 686 687 Output Paramters: 688 . out - matrix to receive 689 . ctx - optional pointer to MatScatterCtx - allows reuse of communcation 690 pattern if same scatter used more than once. 691 692 Keywords: matrix, scatter, gather 693 @*/ 694 int MatScatterBegin(Mat mat,IS is1,IS is2,Mat out,IS is3, IS is4, 695 MatScatterCtx *ctx) 696 { 697 VALIDHEADER(mat,MAT_COOKIE); VALIDHEADER(out,MAT_COOKIE); 698 if (mat->ops->scatbegin) 699 return (*mat->ops->scatbegin)(mat,is1,is2,out,is3,is4,ctx); 700 SETERR(1,"No MatScatterBegin for implementation"); 701 } 702 703