1ee1cef2cSJed Brown 2b45d2f2cSJed Brown #include <petsc-private/matimpl.h> /*I "petscmat.h" I*/ 3ee1cef2cSJed Brown 4ee1cef2cSJed Brown typedef struct { 5ee1cef2cSJed Brown IS isrow,iscol; /* rows and columns in submatrix, only used to check consistency */ 6ee1cef2cSJed Brown Vec left,right; /* optional scaling */ 7ee1cef2cSJed Brown Vec olwork,orwork; /* work vectors outside the scatters, only touched by PreScale and only created if needed*/ 8ee1cef2cSJed Brown Vec lwork,rwork; /* work vectors inside the scatters */ 9ee1cef2cSJed Brown VecScatter lrestrict,rprolong; 10ee1cef2cSJed Brown Mat A; 11ee1cef2cSJed Brown PetscScalar scale; 12ee1cef2cSJed Brown } Mat_SubMatrix; 13ee1cef2cSJed Brown 14ee1cef2cSJed Brown #undef __FUNCT__ 15ee1cef2cSJed Brown #define __FUNCT__ "PreScaleLeft" 16ee1cef2cSJed Brown static PetscErrorCode PreScaleLeft(Mat N,Vec x,Vec *xx) 17ee1cef2cSJed Brown { 18ee1cef2cSJed Brown Mat_SubMatrix *Na = (Mat_SubMatrix*)N->data; 19ee1cef2cSJed Brown PetscErrorCode ierr; 20ee1cef2cSJed Brown 21ee1cef2cSJed Brown PetscFunctionBegin; 22ee1cef2cSJed Brown if (!Na->left) { 23ee1cef2cSJed Brown *xx = x; 24ee1cef2cSJed Brown } else { 25ee1cef2cSJed Brown if (!Na->olwork) { 26ee1cef2cSJed Brown ierr = VecDuplicate(Na->left,&Na->olwork);CHKERRQ(ierr); 27ee1cef2cSJed Brown } 281e7f0bbdSJed Brown ierr = VecPointwiseMult(Na->olwork,x,Na->left);CHKERRQ(ierr); 29ee1cef2cSJed Brown *xx = Na->olwork; 30ee1cef2cSJed Brown } 31ee1cef2cSJed Brown PetscFunctionReturn(0); 32ee1cef2cSJed Brown } 33ee1cef2cSJed Brown 34ee1cef2cSJed Brown #undef __FUNCT__ 35ee1cef2cSJed Brown #define __FUNCT__ "PreScaleRight" 36ee1cef2cSJed Brown static PetscErrorCode PreScaleRight(Mat N,Vec x,Vec *xx) 37ee1cef2cSJed Brown { 38ee1cef2cSJed Brown Mat_SubMatrix *Na = (Mat_SubMatrix*)N->data; 39ee1cef2cSJed Brown PetscErrorCode ierr; 40ee1cef2cSJed Brown 41ee1cef2cSJed Brown PetscFunctionBegin; 42ee1cef2cSJed Brown if (!Na->right) { 43ee1cef2cSJed Brown *xx = x; 44ee1cef2cSJed Brown } else { 45ee1cef2cSJed Brown if (!Na->orwork) { 46ee1cef2cSJed Brown ierr = VecDuplicate(Na->right,&Na->orwork);CHKERRQ(ierr); 47ee1cef2cSJed Brown } 481e7f0bbdSJed Brown ierr = VecPointwiseMult(Na->orwork,x,Na->right);CHKERRQ(ierr); 49ee1cef2cSJed Brown *xx = Na->orwork; 50ee1cef2cSJed Brown } 51ee1cef2cSJed Brown PetscFunctionReturn(0); 52ee1cef2cSJed Brown } 53ee1cef2cSJed Brown 54ee1cef2cSJed Brown #undef __FUNCT__ 55ee1cef2cSJed Brown #define __FUNCT__ "PostScaleLeft" 56ee1cef2cSJed Brown static PetscErrorCode PostScaleLeft(Mat N,Vec x) 57ee1cef2cSJed Brown { 58ee1cef2cSJed Brown Mat_SubMatrix *Na = (Mat_SubMatrix*)N->data; 59ee1cef2cSJed Brown PetscErrorCode ierr; 60ee1cef2cSJed Brown 61ee1cef2cSJed Brown PetscFunctionBegin; 62ee1cef2cSJed Brown if (Na->left) { 63ee1cef2cSJed Brown ierr = VecPointwiseMult(x,x,Na->left);CHKERRQ(ierr); 64ee1cef2cSJed Brown } 65ee1cef2cSJed Brown PetscFunctionReturn(0); 66ee1cef2cSJed Brown } 67ee1cef2cSJed Brown 68ee1cef2cSJed Brown #undef __FUNCT__ 69ee1cef2cSJed Brown #define __FUNCT__ "PostScaleRight" 70ee1cef2cSJed Brown static PetscErrorCode PostScaleRight(Mat N,Vec x) 71ee1cef2cSJed Brown { 72ee1cef2cSJed Brown Mat_SubMatrix *Na = (Mat_SubMatrix*)N->data; 73ee1cef2cSJed Brown PetscErrorCode ierr; 74ee1cef2cSJed Brown 75ee1cef2cSJed Brown PetscFunctionBegin; 76ee1cef2cSJed Brown if (Na->right) { 77ee1cef2cSJed Brown ierr = VecPointwiseMult(x,x,Na->right);CHKERRQ(ierr); 78ee1cef2cSJed Brown } 79ee1cef2cSJed Brown PetscFunctionReturn(0); 80ee1cef2cSJed Brown } 81ee1cef2cSJed Brown 82ee1cef2cSJed Brown #undef __FUNCT__ 83ee1cef2cSJed Brown #define __FUNCT__ "MatScale_SubMatrix" 84ee1cef2cSJed Brown static PetscErrorCode MatScale_SubMatrix(Mat N,PetscScalar scale) 85ee1cef2cSJed Brown { 86ee1cef2cSJed Brown Mat_SubMatrix *Na = (Mat_SubMatrix*)N->data; 87ee1cef2cSJed Brown 88ee1cef2cSJed Brown PetscFunctionBegin; 89ee1cef2cSJed Brown Na->scale *= scale; 90ee1cef2cSJed Brown PetscFunctionReturn(0); 91ee1cef2cSJed Brown } 92ee1cef2cSJed Brown 93ee1cef2cSJed Brown #undef __FUNCT__ 94ee1cef2cSJed Brown #define __FUNCT__ "MatDiagonalScale_SubMatrix" 95ee1cef2cSJed Brown static PetscErrorCode MatDiagonalScale_SubMatrix(Mat N,Vec left,Vec right) 96ee1cef2cSJed Brown { 97ee1cef2cSJed Brown Mat_SubMatrix *Na = (Mat_SubMatrix*)N->data; 98ee1cef2cSJed Brown PetscErrorCode ierr; 99ee1cef2cSJed Brown 100ee1cef2cSJed Brown PetscFunctionBegin; 101ee1cef2cSJed Brown if (left) { 102ee1cef2cSJed Brown if (!Na->left) { 103ee1cef2cSJed Brown ierr = VecDuplicate(left,&Na->left);CHKERRQ(ierr); 104ee1cef2cSJed Brown ierr = VecCopy(left,Na->left);CHKERRQ(ierr); 105ee1cef2cSJed Brown } else { 106ee1cef2cSJed Brown ierr = VecPointwiseMult(Na->left,left,Na->left);CHKERRQ(ierr); 107ee1cef2cSJed Brown } 108ee1cef2cSJed Brown } 109ee1cef2cSJed Brown if (right) { 110ee1cef2cSJed Brown if (!Na->right) { 111ee1cef2cSJed Brown ierr = VecDuplicate(right,&Na->right);CHKERRQ(ierr); 112ee1cef2cSJed Brown ierr = VecCopy(right,Na->right);CHKERRQ(ierr); 113ee1cef2cSJed Brown } else { 114ee1cef2cSJed Brown ierr = VecPointwiseMult(Na->right,right,Na->right);CHKERRQ(ierr); 115ee1cef2cSJed Brown } 116ee1cef2cSJed Brown } 117ee1cef2cSJed Brown PetscFunctionReturn(0); 118ee1cef2cSJed Brown } 119ee1cef2cSJed Brown 120ee1cef2cSJed Brown #undef __FUNCT__ 121ee1cef2cSJed Brown #define __FUNCT__ "MatMult_SubMatrix" 122ee1cef2cSJed Brown static PetscErrorCode MatMult_SubMatrix(Mat N,Vec x,Vec y) 123ee1cef2cSJed Brown { 124ee1cef2cSJed Brown Mat_SubMatrix *Na = (Mat_SubMatrix*)N->data; 125b7784de6SSatish Balay Vec xx = 0; 126ee1cef2cSJed Brown PetscErrorCode ierr; 127ee1cef2cSJed Brown 128ee1cef2cSJed Brown PetscFunctionBegin; 129ee1cef2cSJed Brown ierr = PreScaleRight(N,x,&xx);CHKERRQ(ierr); 130ee1cef2cSJed Brown ierr = VecZeroEntries(Na->rwork);CHKERRQ(ierr); 131ee1cef2cSJed Brown ierr = VecScatterBegin(Na->rprolong,xx,Na->rwork,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 132ee1cef2cSJed Brown ierr = VecScatterEnd (Na->rprolong,xx,Na->rwork,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 133ee1cef2cSJed Brown ierr = MatMult(Na->A,Na->rwork,Na->lwork);CHKERRQ(ierr); 134ee1cef2cSJed Brown ierr = VecScatterBegin(Na->lrestrict,Na->lwork,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 135ee1cef2cSJed Brown ierr = VecScatterEnd (Na->lrestrict,Na->lwork,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 136ee1cef2cSJed Brown ierr = PostScaleLeft(N,y);CHKERRQ(ierr); 137ee1cef2cSJed Brown ierr = VecScale(y,Na->scale);CHKERRQ(ierr); 138ee1cef2cSJed Brown PetscFunctionReturn(0); 139ee1cef2cSJed Brown } 140ee1cef2cSJed Brown 141ee1cef2cSJed Brown #undef __FUNCT__ 142ee1cef2cSJed Brown #define __FUNCT__ "MatMultAdd_SubMatrix" 143ee1cef2cSJed Brown static PetscErrorCode MatMultAdd_SubMatrix(Mat N,Vec v1,Vec v2,Vec v3) 144ee1cef2cSJed Brown { 145ee1cef2cSJed Brown Mat_SubMatrix *Na = (Mat_SubMatrix*)N->data; 146b7784de6SSatish Balay Vec xx = 0; 147ee1cef2cSJed Brown PetscErrorCode ierr; 148ee1cef2cSJed Brown 149ee1cef2cSJed Brown PetscFunctionBegin; 150ee1cef2cSJed Brown ierr = PreScaleRight(N,v1,&xx);CHKERRQ(ierr); 151ee1cef2cSJed Brown ierr = VecZeroEntries(Na->rwork);CHKERRQ(ierr); 152ee1cef2cSJed Brown ierr = VecScatterBegin(Na->rprolong,xx,Na->rwork,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 153ee1cef2cSJed Brown ierr = VecScatterEnd (Na->rprolong,xx,Na->rwork,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 154ee1cef2cSJed Brown ierr = MatMult(Na->A,Na->rwork,Na->lwork);CHKERRQ(ierr); 1551e7f0bbdSJed Brown if (v2 == v3) { 156d4a378daSJed Brown if (Na->scale == (PetscScalar)1.0 && !Na->left) { 1571e7f0bbdSJed Brown ierr = VecScatterBegin(Na->lrestrict,Na->lwork,v3,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1581e7f0bbdSJed Brown ierr = VecScatterEnd (Na->lrestrict,Na->lwork,v3,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1591e7f0bbdSJed Brown } else { 1601e7f0bbdSJed Brown if (!Na->olwork) {ierr = VecDuplicate(v3,&Na->olwork);CHKERRQ(ierr);} 1611e7f0bbdSJed Brown ierr = VecScatterBegin(Na->lrestrict,Na->lwork,Na->olwork,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1621e7f0bbdSJed Brown ierr = VecScatterEnd (Na->lrestrict,Na->lwork,Na->olwork,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1631e7f0bbdSJed Brown ierr = PostScaleLeft(N,Na->olwork);CHKERRQ(ierr); 1641e7f0bbdSJed Brown ierr = VecAXPY(v3,Na->scale,Na->olwork);CHKERRQ(ierr); 1651e7f0bbdSJed Brown } 1661e7f0bbdSJed Brown } else { 167ee1cef2cSJed Brown ierr = VecScatterBegin(Na->lrestrict,Na->lwork,v3,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 168ee1cef2cSJed Brown ierr = VecScatterEnd (Na->lrestrict,Na->lwork,v3,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 169ee1cef2cSJed Brown ierr = PostScaleLeft(N,v3);CHKERRQ(ierr); 170ee1cef2cSJed Brown ierr = VecAYPX(v3,Na->scale,v2);CHKERRQ(ierr); 1711e7f0bbdSJed Brown } 172ee1cef2cSJed Brown PetscFunctionReturn(0); 173ee1cef2cSJed Brown } 174ee1cef2cSJed Brown 175ee1cef2cSJed Brown #undef __FUNCT__ 176ee1cef2cSJed Brown #define __FUNCT__ "MatMultTranspose_SubMatrix" 177ee1cef2cSJed Brown static PetscErrorCode MatMultTranspose_SubMatrix(Mat N,Vec x,Vec y) 178ee1cef2cSJed Brown { 179ee1cef2cSJed Brown Mat_SubMatrix *Na = (Mat_SubMatrix*)N->data; 180b7784de6SSatish Balay Vec xx = 0; 181ee1cef2cSJed Brown PetscErrorCode ierr; 182ee1cef2cSJed Brown 183ee1cef2cSJed Brown PetscFunctionBegin; 184ee1cef2cSJed Brown ierr = PreScaleLeft(N,x,&xx);CHKERRQ(ierr); 185ee1cef2cSJed Brown ierr = VecZeroEntries(Na->lwork);CHKERRQ(ierr); 186ee1cef2cSJed Brown ierr = VecScatterBegin(Na->lrestrict,xx,Na->lwork,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 187ee1cef2cSJed Brown ierr = VecScatterEnd (Na->lrestrict,xx,Na->lwork,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 188ee1cef2cSJed Brown ierr = MatMultTranspose(Na->A,Na->lwork,Na->rwork);CHKERRQ(ierr); 1891e7f0bbdSJed Brown ierr = VecScatterBegin(Na->rprolong,Na->rwork,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1901e7f0bbdSJed Brown ierr = VecScatterEnd (Na->rprolong,Na->rwork,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 191ee1cef2cSJed Brown ierr = PostScaleRight(N,y);CHKERRQ(ierr); 192ee1cef2cSJed Brown ierr = VecScale(y,Na->scale);CHKERRQ(ierr); 193ee1cef2cSJed Brown PetscFunctionReturn(0); 194ee1cef2cSJed Brown } 195ee1cef2cSJed Brown 196ee1cef2cSJed Brown #undef __FUNCT__ 197ee1cef2cSJed Brown #define __FUNCT__ "MatMultTransposeAdd_SubMatrix" 198ee1cef2cSJed Brown static PetscErrorCode MatMultTransposeAdd_SubMatrix(Mat N,Vec v1,Vec v2,Vec v3) 199ee1cef2cSJed Brown { 200ee1cef2cSJed Brown Mat_SubMatrix *Na = (Mat_SubMatrix*)N->data; 201b7784de6SSatish Balay Vec xx = 0; 202ee1cef2cSJed Brown PetscErrorCode ierr; 203ee1cef2cSJed Brown 204ee1cef2cSJed Brown PetscFunctionBegin; 205ee1cef2cSJed Brown ierr = PreScaleLeft(N,v1,&xx);CHKERRQ(ierr); 206ee1cef2cSJed Brown ierr = VecZeroEntries(Na->lwork);CHKERRQ(ierr); 207ee1cef2cSJed Brown ierr = VecScatterBegin(Na->lrestrict,xx,Na->lwork,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 208ee1cef2cSJed Brown ierr = VecScatterEnd (Na->lrestrict,xx,Na->lwork,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 209ee1cef2cSJed Brown ierr = MatMultTranspose(Na->A,Na->lwork,Na->rwork);CHKERRQ(ierr); 2101e7f0bbdSJed Brown if (v2 == v3) { 211d4a378daSJed Brown if (Na->scale == (PetscScalar)1.0 && !Na->right) { 2121e7f0bbdSJed Brown ierr = VecScatterBegin(Na->rprolong,Na->rwork,v3,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 2131e7f0bbdSJed Brown ierr = VecScatterEnd (Na->rprolong,Na->rwork,v3,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 2141e7f0bbdSJed Brown } else { 2151e7f0bbdSJed Brown if (!Na->orwork) {ierr = VecDuplicate(v3,&Na->orwork);CHKERRQ(ierr);} 2161e7f0bbdSJed Brown ierr = VecScatterBegin(Na->rprolong,Na->rwork,Na->orwork,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 2171e7f0bbdSJed Brown ierr = VecScatterEnd (Na->rprolong,Na->rwork,Na->orwork,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 2181e7f0bbdSJed Brown ierr = PostScaleRight(N,Na->orwork);CHKERRQ(ierr); 2191e7f0bbdSJed Brown ierr = VecAXPY(v3,Na->scale,Na->orwork);CHKERRQ(ierr); 2201e7f0bbdSJed Brown } 2211e7f0bbdSJed Brown } else { 2221e7f0bbdSJed Brown ierr = VecScatterBegin(Na->rprolong,Na->rwork,v3,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 2231e7f0bbdSJed Brown ierr = VecScatterEnd (Na->rprolong,Na->rwork,v3,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 224ee1cef2cSJed Brown ierr = PostScaleRight(N,v3);CHKERRQ(ierr); 225ee1cef2cSJed Brown ierr = VecAYPX(v3,Na->scale,v2);CHKERRQ(ierr); 2261e7f0bbdSJed Brown } 227ee1cef2cSJed Brown PetscFunctionReturn(0); 228ee1cef2cSJed Brown } 229ee1cef2cSJed Brown 230ee1cef2cSJed Brown #undef __FUNCT__ 231ee1cef2cSJed Brown #define __FUNCT__ "MatDestroy_SubMatrix" 232ee1cef2cSJed Brown static PetscErrorCode MatDestroy_SubMatrix(Mat N) 233ee1cef2cSJed Brown { 234ee1cef2cSJed Brown Mat_SubMatrix *Na = (Mat_SubMatrix*)N->data; 235ee1cef2cSJed Brown PetscErrorCode ierr; 236ee1cef2cSJed Brown 237ee1cef2cSJed Brown PetscFunctionBegin; 2386bf464f9SBarry Smith ierr = ISDestroy(&Na->isrow);CHKERRQ(ierr); 2396bf464f9SBarry Smith ierr = ISDestroy(&Na->iscol);CHKERRQ(ierr); 2406bf464f9SBarry Smith ierr = VecDestroy(&Na->left);CHKERRQ(ierr); 2416bf464f9SBarry Smith ierr = VecDestroy(&Na->right);CHKERRQ(ierr); 2426bf464f9SBarry Smith ierr = VecDestroy(&Na->olwork);CHKERRQ(ierr); 2436bf464f9SBarry Smith ierr = VecDestroy(&Na->orwork);CHKERRQ(ierr); 2446bf464f9SBarry Smith ierr = VecDestroy(&Na->lwork);CHKERRQ(ierr); 2456bf464f9SBarry Smith ierr = VecDestroy(&Na->rwork);CHKERRQ(ierr); 2466bf464f9SBarry Smith ierr = VecScatterDestroy(&Na->lrestrict);CHKERRQ(ierr); 2476bf464f9SBarry Smith ierr = VecScatterDestroy(&Na->rprolong);CHKERRQ(ierr); 2486bf464f9SBarry Smith ierr = MatDestroy(&Na->A);CHKERRQ(ierr); 2496bf464f9SBarry Smith ierr = PetscFree(N->data);CHKERRQ(ierr); 250ee1cef2cSJed Brown PetscFunctionReturn(0); 251ee1cef2cSJed Brown } 252ee1cef2cSJed Brown 253ee1cef2cSJed Brown #undef __FUNCT__ 254ee1cef2cSJed Brown #define __FUNCT__ "MatCreateSubMatrix" 255ee1cef2cSJed Brown /*@ 256ee1cef2cSJed Brown MatCreateSubMatrix - Creates a composite matrix that acts as a submatrix 257ee1cef2cSJed Brown 258ee1cef2cSJed Brown Collective on Mat 259ee1cef2cSJed Brown 260ee1cef2cSJed Brown Input Parameters: 261ee1cef2cSJed Brown + A - matrix that we will extract a submatrix of 262ee1cef2cSJed Brown . isrow - rows to be present in the submatrix 263ee1cef2cSJed Brown - iscol - columns to be present in the submatrix 264ee1cef2cSJed Brown 265ee1cef2cSJed Brown Output Parameters: 266ee1cef2cSJed Brown . newmat - new matrix 267ee1cef2cSJed Brown 268ee1cef2cSJed Brown Level: developer 269ee1cef2cSJed Brown 270ee1cef2cSJed Brown Notes: 271ee1cef2cSJed Brown Most will use MatGetSubMatrix which provides a more efficient representation if it is available. 272ee1cef2cSJed Brown 273ee1cef2cSJed Brown .seealso: MatGetSubMatrix(), MatSubMatrixUpdate() 274ee1cef2cSJed Brown @*/ 2757087cfbeSBarry Smith PetscErrorCode MatCreateSubMatrix(Mat A,IS isrow,IS iscol,Mat *newmat) 276ee1cef2cSJed Brown { 277ee1cef2cSJed Brown Vec left,right; 278ee1cef2cSJed Brown PetscInt m,n; 279ee1cef2cSJed Brown Mat N; 280ee1cef2cSJed Brown Mat_SubMatrix *Na; 281ee1cef2cSJed Brown PetscErrorCode ierr; 282ee1cef2cSJed Brown 283ee1cef2cSJed Brown PetscFunctionBegin; 2840700a824SBarry Smith PetscValidHeaderSpecific(A,MAT_CLASSID,1); 2850700a824SBarry Smith PetscValidHeaderSpecific(isrow,IS_CLASSID,2); 2860700a824SBarry Smith PetscValidHeaderSpecific(iscol,IS_CLASSID,3); 287ee1cef2cSJed Brown PetscValidPointer(newmat,4); 288ee1cef2cSJed Brown *newmat = 0; 289ee1cef2cSJed Brown 290*ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)A),&N);CHKERRQ(ierr); 291ee1cef2cSJed Brown ierr = ISGetLocalSize(isrow,&m);CHKERRQ(ierr); 292ee1cef2cSJed Brown ierr = ISGetLocalSize(iscol,&n);CHKERRQ(ierr); 293ee1cef2cSJed Brown ierr = MatSetSizes(N,m,n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 294ee1cef2cSJed Brown ierr = PetscObjectChangeTypeName((PetscObject)N,MATSUBMATRIX);CHKERRQ(ierr); 295ee1cef2cSJed Brown 296ee1cef2cSJed Brown ierr = PetscNewLog(N,Mat_SubMatrix,&Na);CHKERRQ(ierr); 297ee1cef2cSJed Brown N->data = (void*)Na; 298ee1cef2cSJed Brown ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 299ee1cef2cSJed Brown ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr); 300ee1cef2cSJed Brown ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr); 301ee1cef2cSJed Brown Na->A = A; 302ee1cef2cSJed Brown Na->isrow = isrow; 303ee1cef2cSJed Brown Na->iscol = iscol; 304ee1cef2cSJed Brown Na->scale = 1.0; 305ee1cef2cSJed Brown 306ee1cef2cSJed Brown N->ops->destroy = MatDestroy_SubMatrix; 307ee1cef2cSJed Brown N->ops->mult = MatMult_SubMatrix; 308ee1cef2cSJed Brown N->ops->multadd = MatMultAdd_SubMatrix; 309ee1cef2cSJed Brown N->ops->multtranspose = MatMultTranspose_SubMatrix; 310ee1cef2cSJed Brown N->ops->multtransposeadd = MatMultTransposeAdd_SubMatrix; 311ee1cef2cSJed Brown N->ops->scale = MatScale_SubMatrix; 312ee1cef2cSJed Brown N->ops->diagonalscale = MatDiagonalScale_SubMatrix; 313ee1cef2cSJed Brown 31426283091SBarry Smith ierr = PetscLayoutSetBlockSize(N->rmap,A->rmap->bs);CHKERRQ(ierr); 31526283091SBarry Smith ierr = PetscLayoutSetBlockSize(N->cmap,A->cmap->bs);CHKERRQ(ierr); 31626283091SBarry Smith ierr = PetscLayoutSetUp(N->rmap);CHKERRQ(ierr); 31726283091SBarry Smith ierr = PetscLayoutSetUp(N->cmap);CHKERRQ(ierr); 318ee1cef2cSJed Brown 319ee1cef2cSJed Brown ierr = MatGetVecs(A,&Na->rwork,&Na->lwork);CHKERRQ(ierr); 320*ce94432eSBarry Smith ierr = VecCreate(PetscObjectComm((PetscObject)isrow),&left);CHKERRQ(ierr); 321*ce94432eSBarry Smith ierr = VecCreate(PetscObjectComm((PetscObject)iscol),&right);CHKERRQ(ierr); 322ee1cef2cSJed Brown ierr = VecSetSizes(left,m,PETSC_DETERMINE);CHKERRQ(ierr); 323ee1cef2cSJed Brown ierr = VecSetSizes(right,n,PETSC_DETERMINE);CHKERRQ(ierr); 324ee1cef2cSJed Brown ierr = VecSetUp(left);CHKERRQ(ierr); 325ee1cef2cSJed Brown ierr = VecSetUp(right);CHKERRQ(ierr); 3260298fd71SBarry Smith ierr = VecScatterCreate(Na->lwork,isrow,left,NULL,&Na->lrestrict);CHKERRQ(ierr); 3270298fd71SBarry Smith ierr = VecScatterCreate(right,NULL,Na->rwork,iscol,&Na->rprolong);CHKERRQ(ierr); 3286bf464f9SBarry Smith ierr = VecDestroy(&left);CHKERRQ(ierr); 3296bf464f9SBarry Smith ierr = VecDestroy(&right);CHKERRQ(ierr); 330ee1cef2cSJed Brown 331be7c243fSJed Brown N->assembled = PETSC_TRUE; 33226fbe8dcSKarl Rupp 333be7c243fSJed Brown ierr = MatSetUp(N);CHKERRQ(ierr); 33426fbe8dcSKarl Rupp 335ee1cef2cSJed Brown *newmat = N; 336ee1cef2cSJed Brown PetscFunctionReturn(0); 337ee1cef2cSJed Brown } 338ee1cef2cSJed Brown 339ee1cef2cSJed Brown 340ee1cef2cSJed Brown #undef __FUNCT__ 341ee1cef2cSJed Brown #define __FUNCT__ "MatSubMatrixUpdate" 342ee1cef2cSJed Brown /*@ 343ee1cef2cSJed Brown MatSubMatrixUpdate - Updates a submatrix 344ee1cef2cSJed Brown 345ee1cef2cSJed Brown Collective on Mat 346ee1cef2cSJed Brown 347ee1cef2cSJed Brown Input Parameters: 348ee1cef2cSJed Brown + N - submatrix to update 349ee1cef2cSJed Brown . A - full matrix in the submatrix 350ee1cef2cSJed Brown . isrow - rows in the update (same as the first time the submatrix was created) 351ee1cef2cSJed Brown - iscol - columns in the update (same as the first time the submatrix was created) 352ee1cef2cSJed Brown 353ee1cef2cSJed Brown Level: developer 354ee1cef2cSJed Brown 355ee1cef2cSJed Brown Notes: 356ee1cef2cSJed Brown Most will use MatGetSubMatrix which provides a more efficient representation if it is available. 357ee1cef2cSJed Brown 358ee1cef2cSJed Brown .seealso: MatGetSubMatrix(), MatCreateSubMatrix() 359ee1cef2cSJed Brown @*/ 3607087cfbeSBarry Smith PetscErrorCode MatSubMatrixUpdate(Mat N,Mat A,IS isrow,IS iscol) 361ee1cef2cSJed Brown { 362ee1cef2cSJed Brown PetscErrorCode ierr; 363ace3abfcSBarry Smith PetscBool flg; 364ee1cef2cSJed Brown Mat_SubMatrix *Na; 365ee1cef2cSJed Brown 366ee1cef2cSJed Brown PetscFunctionBegin; 3670700a824SBarry Smith PetscValidHeaderSpecific(N,MAT_CLASSID,1); 3680700a824SBarry Smith PetscValidHeaderSpecific(A,MAT_CLASSID,2); 3690700a824SBarry Smith PetscValidHeaderSpecific(isrow,IS_CLASSID,3); 3700700a824SBarry Smith PetscValidHeaderSpecific(iscol,IS_CLASSID,4); 371251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)N,MATSUBMATRIX,&flg);CHKERRQ(ierr); 372*ce94432eSBarry Smith if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix has wrong type"); 373ee1cef2cSJed Brown 374ee1cef2cSJed Brown Na = (Mat_SubMatrix*)N->data; 375ee1cef2cSJed Brown ierr = ISEqual(isrow,Na->isrow,&flg);CHKERRQ(ierr); 376e32f2f54SBarry Smith if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Cannot update submatrix with different row indices"); 377ee1cef2cSJed Brown ierr = ISEqual(iscol,Na->iscol,&flg);CHKERRQ(ierr); 378e32f2f54SBarry Smith if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Cannot update submatrix with different column indices"); 379ee1cef2cSJed Brown 380ee1cef2cSJed Brown ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 3816bf464f9SBarry Smith ierr = MatDestroy(&Na->A);CHKERRQ(ierr); 382ee1cef2cSJed Brown Na->A = A; 383ee1cef2cSJed Brown 384ee1cef2cSJed Brown Na->scale = 1.0; 3856bf464f9SBarry Smith ierr = VecDestroy(&Na->left);CHKERRQ(ierr); 3866bf464f9SBarry Smith ierr = VecDestroy(&Na->right);CHKERRQ(ierr); 387ee1cef2cSJed Brown PetscFunctionReturn(0); 388ee1cef2cSJed Brown } 389