1ee1cef2cSJed Brown #define PETSCMAT_DLL 2ee1cef2cSJed Brown 3ee1cef2cSJed Brown #include "private/matimpl.h" /*I "petscmat.h" I*/ 4ee1cef2cSJed Brown 5ee1cef2cSJed Brown typedef struct { 6ee1cef2cSJed Brown IS isrow,iscol; /* rows and columns in submatrix, only used to check consistency */ 7ee1cef2cSJed Brown Vec left,right; /* optional scaling */ 8ee1cef2cSJed Brown Vec olwork,orwork; /* work vectors outside the scatters, only touched by PreScale and only created if needed*/ 9ee1cef2cSJed Brown Vec lwork,rwork; /* work vectors inside the scatters */ 10ee1cef2cSJed Brown VecScatter lrestrict,rprolong; 11ee1cef2cSJed Brown Mat A; 12ee1cef2cSJed Brown PetscScalar scale; 13ee1cef2cSJed Brown } Mat_SubMatrix; 14ee1cef2cSJed Brown 15ee1cef2cSJed Brown #undef __FUNCT__ 16ee1cef2cSJed Brown #define __FUNCT__ "PreScaleLeft" 17ee1cef2cSJed Brown static PetscErrorCode PreScaleLeft(Mat N,Vec x,Vec *xx) 18ee1cef2cSJed Brown { 19ee1cef2cSJed Brown Mat_SubMatrix *Na = (Mat_SubMatrix*)N->data; 20ee1cef2cSJed Brown PetscErrorCode ierr; 21ee1cef2cSJed Brown 22ee1cef2cSJed Brown PetscFunctionBegin; 23ee1cef2cSJed Brown if (!Na->left) { 24ee1cef2cSJed Brown *xx = x; 25ee1cef2cSJed Brown } else { 26ee1cef2cSJed Brown if (!Na->olwork) { 27ee1cef2cSJed Brown ierr = VecDuplicate(Na->left,&Na->olwork);CHKERRQ(ierr); 28ee1cef2cSJed Brown } 29*1e7f0bbdSJed Brown ierr = VecPointwiseMult(Na->olwork,x,Na->left);CHKERRQ(ierr); 30ee1cef2cSJed Brown *xx = Na->olwork; 31ee1cef2cSJed Brown } 32ee1cef2cSJed Brown PetscFunctionReturn(0); 33ee1cef2cSJed Brown } 34ee1cef2cSJed Brown 35ee1cef2cSJed Brown #undef __FUNCT__ 36ee1cef2cSJed Brown #define __FUNCT__ "PreScaleRight" 37ee1cef2cSJed Brown static PetscErrorCode PreScaleRight(Mat N,Vec x,Vec *xx) 38ee1cef2cSJed Brown { 39ee1cef2cSJed Brown Mat_SubMatrix *Na = (Mat_SubMatrix*)N->data; 40ee1cef2cSJed Brown PetscErrorCode ierr; 41ee1cef2cSJed Brown 42ee1cef2cSJed Brown PetscFunctionBegin; 43ee1cef2cSJed Brown if (!Na->right) { 44ee1cef2cSJed Brown *xx = x; 45ee1cef2cSJed Brown } else { 46ee1cef2cSJed Brown if (!Na->orwork) { 47ee1cef2cSJed Brown ierr = VecDuplicate(Na->right,&Na->orwork);CHKERRQ(ierr); 48ee1cef2cSJed Brown } 49*1e7f0bbdSJed Brown ierr = VecPointwiseMult(Na->orwork,x,Na->right);CHKERRQ(ierr); 50ee1cef2cSJed Brown *xx = Na->orwork; 51ee1cef2cSJed Brown } 52ee1cef2cSJed Brown PetscFunctionReturn(0); 53ee1cef2cSJed Brown } 54ee1cef2cSJed Brown 55ee1cef2cSJed Brown #undef __FUNCT__ 56ee1cef2cSJed Brown #define __FUNCT__ "PostScaleLeft" 57ee1cef2cSJed Brown static PetscErrorCode PostScaleLeft(Mat N,Vec x) 58ee1cef2cSJed Brown { 59ee1cef2cSJed Brown Mat_SubMatrix *Na = (Mat_SubMatrix*)N->data; 60ee1cef2cSJed Brown PetscErrorCode ierr; 61ee1cef2cSJed Brown 62ee1cef2cSJed Brown PetscFunctionBegin; 63ee1cef2cSJed Brown if (Na->left) { 64ee1cef2cSJed Brown ierr = VecPointwiseMult(x,x,Na->left);CHKERRQ(ierr); 65ee1cef2cSJed Brown } 66ee1cef2cSJed Brown PetscFunctionReturn(0); 67ee1cef2cSJed Brown } 68ee1cef2cSJed Brown 69ee1cef2cSJed Brown #undef __FUNCT__ 70ee1cef2cSJed Brown #define __FUNCT__ "PostScaleRight" 71ee1cef2cSJed Brown static PetscErrorCode PostScaleRight(Mat N,Vec x) 72ee1cef2cSJed Brown { 73ee1cef2cSJed Brown Mat_SubMatrix *Na = (Mat_SubMatrix*)N->data; 74ee1cef2cSJed Brown PetscErrorCode ierr; 75ee1cef2cSJed Brown 76ee1cef2cSJed Brown PetscFunctionBegin; 77ee1cef2cSJed Brown if (Na->right) { 78ee1cef2cSJed Brown ierr = VecPointwiseMult(x,x,Na->right);CHKERRQ(ierr); 79ee1cef2cSJed Brown } 80ee1cef2cSJed Brown PetscFunctionReturn(0); 81ee1cef2cSJed Brown } 82ee1cef2cSJed Brown 83ee1cef2cSJed Brown #undef __FUNCT__ 84ee1cef2cSJed Brown #define __FUNCT__ "MatScale_SubMatrix" 85ee1cef2cSJed Brown static PetscErrorCode MatScale_SubMatrix(Mat N,PetscScalar scale) 86ee1cef2cSJed Brown { 87ee1cef2cSJed Brown Mat_SubMatrix *Na = (Mat_SubMatrix*)N->data; 88ee1cef2cSJed Brown 89ee1cef2cSJed Brown PetscFunctionBegin; 90ee1cef2cSJed Brown Na->scale *= scale; 91ee1cef2cSJed Brown PetscFunctionReturn(0); 92ee1cef2cSJed Brown } 93ee1cef2cSJed Brown 94ee1cef2cSJed Brown #undef __FUNCT__ 95ee1cef2cSJed Brown #define __FUNCT__ "MatDiagonalScale_SubMatrix" 96ee1cef2cSJed Brown static PetscErrorCode MatDiagonalScale_SubMatrix(Mat N,Vec left,Vec right) 97ee1cef2cSJed Brown { 98ee1cef2cSJed Brown Mat_SubMatrix *Na = (Mat_SubMatrix*)N->data; 99ee1cef2cSJed Brown PetscErrorCode ierr; 100ee1cef2cSJed Brown 101ee1cef2cSJed Brown PetscFunctionBegin; 102ee1cef2cSJed Brown if (left) { 103ee1cef2cSJed Brown if (!Na->left) { 104ee1cef2cSJed Brown ierr = VecDuplicate(left,&Na->left);CHKERRQ(ierr); 105ee1cef2cSJed Brown ierr = VecCopy(left,Na->left);CHKERRQ(ierr); 106ee1cef2cSJed Brown } else { 107ee1cef2cSJed Brown ierr = VecPointwiseMult(Na->left,left,Na->left);CHKERRQ(ierr); 108ee1cef2cSJed Brown } 109ee1cef2cSJed Brown } 110ee1cef2cSJed Brown if (right) { 111ee1cef2cSJed Brown if (!Na->right) { 112ee1cef2cSJed Brown ierr = VecDuplicate(right,&Na->right);CHKERRQ(ierr); 113ee1cef2cSJed Brown ierr = VecCopy(right,Na->right);CHKERRQ(ierr); 114ee1cef2cSJed Brown } else { 115ee1cef2cSJed Brown ierr = VecPointwiseMult(Na->right,right,Na->right);CHKERRQ(ierr); 116ee1cef2cSJed Brown } 117ee1cef2cSJed Brown } 118ee1cef2cSJed Brown PetscFunctionReturn(0); 119ee1cef2cSJed Brown } 120ee1cef2cSJed Brown 121ee1cef2cSJed Brown #undef __FUNCT__ 122ee1cef2cSJed Brown #define __FUNCT__ "MatMult_SubMatrix" 123ee1cef2cSJed Brown static PetscErrorCode MatMult_SubMatrix(Mat N,Vec x,Vec y) 124ee1cef2cSJed Brown { 125ee1cef2cSJed Brown Mat_SubMatrix *Na = (Mat_SubMatrix*)N->data; 126b7784de6SSatish Balay Vec xx=0; 127ee1cef2cSJed Brown PetscErrorCode ierr; 128ee1cef2cSJed Brown 129ee1cef2cSJed Brown PetscFunctionBegin; 130ee1cef2cSJed Brown ierr = PreScaleRight(N,x,&xx);CHKERRQ(ierr); 131ee1cef2cSJed Brown ierr = VecZeroEntries(Na->rwork);CHKERRQ(ierr); 132ee1cef2cSJed Brown ierr = VecScatterBegin(Na->rprolong,xx,Na->rwork,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 133ee1cef2cSJed Brown ierr = VecScatterEnd (Na->rprolong,xx,Na->rwork,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 134ee1cef2cSJed Brown ierr = MatMult(Na->A,Na->rwork,Na->lwork);CHKERRQ(ierr); 135ee1cef2cSJed Brown ierr = VecScatterBegin(Na->lrestrict,Na->lwork,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 136ee1cef2cSJed Brown ierr = VecScatterEnd (Na->lrestrict,Na->lwork,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 137ee1cef2cSJed Brown ierr = PostScaleLeft(N,y);CHKERRQ(ierr); 138ee1cef2cSJed Brown ierr = VecScale(y,Na->scale);CHKERRQ(ierr); 139ee1cef2cSJed Brown PetscFunctionReturn(0); 140ee1cef2cSJed Brown } 141ee1cef2cSJed Brown 142ee1cef2cSJed Brown #undef __FUNCT__ 143ee1cef2cSJed Brown #define __FUNCT__ "MatMultAdd_SubMatrix" 144ee1cef2cSJed Brown static PetscErrorCode MatMultAdd_SubMatrix(Mat N,Vec v1,Vec v2,Vec v3) 145ee1cef2cSJed Brown { 146ee1cef2cSJed Brown Mat_SubMatrix *Na = (Mat_SubMatrix*)N->data; 147b7784de6SSatish Balay Vec xx=0; 148ee1cef2cSJed Brown PetscErrorCode ierr; 149ee1cef2cSJed Brown 150ee1cef2cSJed Brown PetscFunctionBegin; 151ee1cef2cSJed Brown ierr = PreScaleRight(N,v1,&xx);CHKERRQ(ierr); 152ee1cef2cSJed Brown ierr = VecZeroEntries(Na->rwork);CHKERRQ(ierr); 153ee1cef2cSJed Brown ierr = VecScatterBegin(Na->rprolong,xx,Na->rwork,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 154ee1cef2cSJed Brown ierr = VecScatterEnd (Na->rprolong,xx,Na->rwork,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 155ee1cef2cSJed Brown ierr = MatMult(Na->A,Na->rwork,Na->lwork);CHKERRQ(ierr); 156*1e7f0bbdSJed Brown if (v2 == v3) { 157*1e7f0bbdSJed Brown if (Na->scale == 1.0 && !Na->left) { 158*1e7f0bbdSJed Brown ierr = VecScatterBegin(Na->lrestrict,Na->lwork,v3,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 159*1e7f0bbdSJed Brown ierr = VecScatterEnd (Na->lrestrict,Na->lwork,v3,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 160*1e7f0bbdSJed Brown } else { 161*1e7f0bbdSJed Brown if (!Na->olwork) {ierr = VecDuplicate(v3,&Na->olwork);CHKERRQ(ierr);} 162*1e7f0bbdSJed Brown ierr = VecScatterBegin(Na->lrestrict,Na->lwork,Na->olwork,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 163*1e7f0bbdSJed Brown ierr = VecScatterEnd (Na->lrestrict,Na->lwork,Na->olwork,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 164*1e7f0bbdSJed Brown ierr = PostScaleLeft(N,Na->olwork);CHKERRQ(ierr); 165*1e7f0bbdSJed Brown ierr = VecAXPY(v3,Na->scale,Na->olwork);CHKERRQ(ierr); 166*1e7f0bbdSJed Brown } 167*1e7f0bbdSJed Brown } else { 168ee1cef2cSJed Brown ierr = VecScatterBegin(Na->lrestrict,Na->lwork,v3,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 169ee1cef2cSJed Brown ierr = VecScatterEnd (Na->lrestrict,Na->lwork,v3,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 170ee1cef2cSJed Brown ierr = PostScaleLeft(N,v3);CHKERRQ(ierr); 171ee1cef2cSJed Brown ierr = VecAYPX(v3,Na->scale,v2);CHKERRQ(ierr); 172*1e7f0bbdSJed Brown } 173ee1cef2cSJed Brown PetscFunctionReturn(0); 174ee1cef2cSJed Brown } 175ee1cef2cSJed Brown 176ee1cef2cSJed Brown #undef __FUNCT__ 177ee1cef2cSJed Brown #define __FUNCT__ "MatMultTranspose_SubMatrix" 178ee1cef2cSJed Brown static PetscErrorCode MatMultTranspose_SubMatrix(Mat N,Vec x,Vec y) 179ee1cef2cSJed Brown { 180ee1cef2cSJed Brown Mat_SubMatrix *Na = (Mat_SubMatrix*)N->data; 181b7784de6SSatish Balay Vec xx=0; 182ee1cef2cSJed Brown PetscErrorCode ierr; 183ee1cef2cSJed Brown 184ee1cef2cSJed Brown PetscFunctionBegin; 185ee1cef2cSJed Brown ierr = PreScaleLeft(N,x,&xx);CHKERRQ(ierr); 186ee1cef2cSJed Brown ierr = VecZeroEntries(Na->lwork);CHKERRQ(ierr); 187ee1cef2cSJed Brown ierr = VecScatterBegin(Na->lrestrict,xx,Na->lwork,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 188ee1cef2cSJed Brown ierr = VecScatterEnd (Na->lrestrict,xx,Na->lwork,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 189ee1cef2cSJed Brown ierr = MatMultTranspose(Na->A,Na->lwork,Na->rwork);CHKERRQ(ierr); 190*1e7f0bbdSJed Brown ierr = VecScatterBegin(Na->rprolong,Na->rwork,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 191*1e7f0bbdSJed Brown ierr = VecScatterEnd (Na->rprolong,Na->rwork,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 192ee1cef2cSJed Brown ierr = PostScaleRight(N,y);CHKERRQ(ierr); 193ee1cef2cSJed Brown ierr = VecScale(y,Na->scale);CHKERRQ(ierr); 194ee1cef2cSJed Brown PetscFunctionReturn(0); 195ee1cef2cSJed Brown } 196ee1cef2cSJed Brown 197ee1cef2cSJed Brown #undef __FUNCT__ 198ee1cef2cSJed Brown #define __FUNCT__ "MatMultTransposeAdd_SubMatrix" 199ee1cef2cSJed Brown static PetscErrorCode MatMultTransposeAdd_SubMatrix(Mat N,Vec v1,Vec v2,Vec v3) 200ee1cef2cSJed Brown { 201ee1cef2cSJed Brown Mat_SubMatrix *Na = (Mat_SubMatrix*)N->data; 202b7784de6SSatish Balay Vec xx =0; 203ee1cef2cSJed Brown PetscErrorCode ierr; 204ee1cef2cSJed Brown 205ee1cef2cSJed Brown PetscFunctionBegin; 206ee1cef2cSJed Brown ierr = PreScaleLeft(N,v1,&xx);CHKERRQ(ierr); 207ee1cef2cSJed Brown ierr = VecZeroEntries(Na->lwork);CHKERRQ(ierr); 208ee1cef2cSJed Brown ierr = VecScatterBegin(Na->lrestrict,xx,Na->lwork,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 209ee1cef2cSJed Brown ierr = VecScatterEnd (Na->lrestrict,xx,Na->lwork,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 210ee1cef2cSJed Brown ierr = MatMultTranspose(Na->A,Na->lwork,Na->rwork);CHKERRQ(ierr); 211*1e7f0bbdSJed Brown if (v2 == v3) { 212*1e7f0bbdSJed Brown if (Na->scale == 1.0 && !Na->right) { 213*1e7f0bbdSJed Brown ierr = VecScatterBegin(Na->rprolong,Na->rwork,v3,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 214*1e7f0bbdSJed Brown ierr = VecScatterEnd (Na->rprolong,Na->rwork,v3,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 215*1e7f0bbdSJed Brown } else { 216*1e7f0bbdSJed Brown if (!Na->orwork) {ierr = VecDuplicate(v3,&Na->orwork);CHKERRQ(ierr);} 217*1e7f0bbdSJed Brown ierr = VecScatterBegin(Na->rprolong,Na->rwork,Na->orwork,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 218*1e7f0bbdSJed Brown ierr = VecScatterEnd (Na->rprolong,Na->rwork,Na->orwork,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 219*1e7f0bbdSJed Brown ierr = PostScaleRight(N,Na->orwork);CHKERRQ(ierr); 220*1e7f0bbdSJed Brown ierr = VecAXPY(v3,Na->scale,Na->orwork);CHKERRQ(ierr); 221*1e7f0bbdSJed Brown } 222*1e7f0bbdSJed Brown } else { 223*1e7f0bbdSJed Brown ierr = VecScatterBegin(Na->rprolong,Na->rwork,v3,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 224*1e7f0bbdSJed Brown ierr = VecScatterEnd (Na->rprolong,Na->rwork,v3,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 225ee1cef2cSJed Brown ierr = PostScaleRight(N,v3);CHKERRQ(ierr); 226ee1cef2cSJed Brown ierr = VecAYPX(v3,Na->scale,v2);CHKERRQ(ierr); 227*1e7f0bbdSJed Brown } 228ee1cef2cSJed Brown PetscFunctionReturn(0); 229ee1cef2cSJed Brown } 230ee1cef2cSJed Brown 231ee1cef2cSJed Brown #undef __FUNCT__ 232ee1cef2cSJed Brown #define __FUNCT__ "MatDestroy_SubMatrix" 233ee1cef2cSJed Brown static PetscErrorCode MatDestroy_SubMatrix(Mat N) 234ee1cef2cSJed Brown { 235ee1cef2cSJed Brown Mat_SubMatrix *Na = (Mat_SubMatrix*)N->data; 236ee1cef2cSJed Brown PetscErrorCode ierr; 237ee1cef2cSJed Brown 238ee1cef2cSJed Brown PetscFunctionBegin; 239ee1cef2cSJed Brown ierr = ISDestroy(Na->isrow);CHKERRQ(ierr); 240ee1cef2cSJed Brown ierr = ISDestroy(Na->iscol);CHKERRQ(ierr); 241ee1cef2cSJed Brown if (Na->left) {ierr = VecDestroy(Na->left);CHKERRQ(ierr);} 242ee1cef2cSJed Brown if (Na->right) {ierr = VecDestroy(Na->right);CHKERRQ(ierr);} 243ee1cef2cSJed Brown if (Na->olwork) {ierr = VecDestroy(Na->olwork);CHKERRQ(ierr);} 244ee1cef2cSJed Brown if (Na->orwork) {ierr = VecDestroy(Na->orwork);CHKERRQ(ierr);} 245ee1cef2cSJed Brown ierr = VecDestroy(Na->lwork);CHKERRQ(ierr); 246ee1cef2cSJed Brown ierr = VecDestroy(Na->rwork);CHKERRQ(ierr); 247ee1cef2cSJed Brown ierr = VecScatterDestroy(Na->lrestrict);CHKERRQ(ierr); 248ee1cef2cSJed Brown ierr = VecScatterDestroy(Na->rprolong);CHKERRQ(ierr); 249ee1cef2cSJed Brown ierr = MatDestroy(Na->A);CHKERRQ(ierr); 250ee1cef2cSJed Brown ierr = PetscFree(Na);CHKERRQ(ierr); 251ee1cef2cSJed Brown PetscFunctionReturn(0); 252ee1cef2cSJed Brown } 253ee1cef2cSJed Brown 254ee1cef2cSJed Brown #undef __FUNCT__ 255ee1cef2cSJed Brown #define __FUNCT__ "MatCreateSubMatrix" 256ee1cef2cSJed Brown /*@ 257ee1cef2cSJed Brown MatCreateSubMatrix - Creates a composite matrix that acts as a submatrix 258ee1cef2cSJed Brown 259ee1cef2cSJed Brown Collective on Mat 260ee1cef2cSJed Brown 261ee1cef2cSJed Brown Input Parameters: 262ee1cef2cSJed Brown + A - matrix that we will extract a submatrix of 263ee1cef2cSJed Brown . isrow - rows to be present in the submatrix 264ee1cef2cSJed Brown - iscol - columns to be present in the submatrix 265ee1cef2cSJed Brown 266ee1cef2cSJed Brown Output Parameters: 267ee1cef2cSJed Brown . newmat - new matrix 268ee1cef2cSJed Brown 269ee1cef2cSJed Brown Level: developer 270ee1cef2cSJed Brown 271ee1cef2cSJed Brown Notes: 272ee1cef2cSJed Brown Most will use MatGetSubMatrix which provides a more efficient representation if it is available. 273ee1cef2cSJed Brown 274ee1cef2cSJed Brown .seealso: MatGetSubMatrix(), MatSubMatrixUpdate() 275ee1cef2cSJed Brown @*/ 276ee1cef2cSJed Brown PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSubMatrix(Mat A,IS isrow,IS iscol,Mat *newmat) 277ee1cef2cSJed Brown { 278ee1cef2cSJed Brown Vec left,right; 279ee1cef2cSJed Brown PetscInt m,n; 280ee1cef2cSJed Brown Mat N; 281ee1cef2cSJed Brown Mat_SubMatrix *Na; 282ee1cef2cSJed Brown PetscErrorCode ierr; 283ee1cef2cSJed Brown 284ee1cef2cSJed Brown PetscFunctionBegin; 285ee1cef2cSJed Brown PetscValidHeaderSpecific(A,MAT_COOKIE,1); 286ee1cef2cSJed Brown PetscValidHeaderSpecific(isrow,IS_COOKIE,2); 287ee1cef2cSJed Brown PetscValidHeaderSpecific(iscol,IS_COOKIE,3); 288ee1cef2cSJed Brown PetscValidPointer(newmat,4); 289ee1cef2cSJed Brown *newmat = 0; 290ee1cef2cSJed Brown 291ee1cef2cSJed Brown ierr = MatCreate(((PetscObject)A)->comm,&N);CHKERRQ(ierr); 292ee1cef2cSJed Brown ierr = ISGetLocalSize(isrow,&m);CHKERRQ(ierr); 293ee1cef2cSJed Brown ierr = ISGetLocalSize(iscol,&n);CHKERRQ(ierr); 294ee1cef2cSJed Brown ierr = MatSetSizes(N,m,n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 295ee1cef2cSJed Brown ierr = PetscObjectChangeTypeName((PetscObject)N,MATSUBMATRIX);CHKERRQ(ierr); 296ee1cef2cSJed Brown 297ee1cef2cSJed Brown ierr = PetscNewLog(N,Mat_SubMatrix,&Na);CHKERRQ(ierr); 298ee1cef2cSJed Brown N->data = (void*)Na; 299ee1cef2cSJed Brown ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 300ee1cef2cSJed Brown ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr); 301ee1cef2cSJed Brown ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr); 302ee1cef2cSJed Brown Na->A = A; 303ee1cef2cSJed Brown Na->isrow = isrow; 304ee1cef2cSJed Brown Na->iscol = iscol; 305ee1cef2cSJed Brown Na->scale = 1.0; 306ee1cef2cSJed Brown 307ee1cef2cSJed Brown N->ops->destroy = MatDestroy_SubMatrix; 308ee1cef2cSJed Brown N->ops->mult = MatMult_SubMatrix; 309ee1cef2cSJed Brown N->ops->multadd = MatMultAdd_SubMatrix; 310ee1cef2cSJed Brown N->ops->multtranspose = MatMultTranspose_SubMatrix; 311ee1cef2cSJed Brown N->ops->multtransposeadd = MatMultTransposeAdd_SubMatrix; 312ee1cef2cSJed Brown N->ops->scale = MatScale_SubMatrix; 313ee1cef2cSJed Brown N->ops->diagonalscale = MatDiagonalScale_SubMatrix; 314ee1cef2cSJed Brown 315ee1cef2cSJed Brown N->assembled = PETSC_TRUE; 316ee1cef2cSJed Brown 31726283091SBarry Smith ierr = PetscLayoutSetBlockSize(N->rmap,A->rmap->bs);CHKERRQ(ierr); 31826283091SBarry Smith ierr = PetscLayoutSetBlockSize(N->cmap,A->cmap->bs);CHKERRQ(ierr); 31926283091SBarry Smith ierr = PetscLayoutSetUp(N->rmap);CHKERRQ(ierr); 32026283091SBarry Smith ierr = PetscLayoutSetUp(N->cmap);CHKERRQ(ierr); 321ee1cef2cSJed Brown 322ee1cef2cSJed Brown ierr = MatGetVecs(A,&Na->rwork,&Na->lwork);CHKERRQ(ierr); 323ee1cef2cSJed Brown ierr = VecCreate(((PetscObject)isrow)->comm,&left);CHKERRQ(ierr); 324ee1cef2cSJed Brown ierr = VecCreate(((PetscObject)iscol)->comm,&right);CHKERRQ(ierr); 325ee1cef2cSJed Brown ierr = VecSetSizes(left,m,PETSC_DETERMINE);CHKERRQ(ierr); 326ee1cef2cSJed Brown ierr = VecSetSizes(right,n,PETSC_DETERMINE);CHKERRQ(ierr); 327ee1cef2cSJed Brown ierr = VecSetUp(left);CHKERRQ(ierr); 328ee1cef2cSJed Brown ierr = VecSetUp(right);CHKERRQ(ierr); 329ee1cef2cSJed Brown ierr = VecScatterCreate(Na->lwork,isrow,left,PETSC_NULL,&Na->lrestrict);CHKERRQ(ierr); 330ee1cef2cSJed Brown ierr = VecScatterCreate(right,PETSC_NULL,Na->rwork,iscol,&Na->rprolong);CHKERRQ(ierr); 331ee1cef2cSJed Brown ierr = VecDestroy(left);CHKERRQ(ierr); 332ee1cef2cSJed Brown ierr = VecDestroy(right);CHKERRQ(ierr); 333ee1cef2cSJed Brown 334ee1cef2cSJed Brown *newmat = N; 335ee1cef2cSJed Brown PetscFunctionReturn(0); 336ee1cef2cSJed Brown } 337ee1cef2cSJed Brown 338ee1cef2cSJed Brown 339ee1cef2cSJed Brown #undef __FUNCT__ 340ee1cef2cSJed Brown #define __FUNCT__ "MatSubMatrixUpdate" 341ee1cef2cSJed Brown /*@ 342ee1cef2cSJed Brown MatSubMatrixUpdate - Updates a submatrix 343ee1cef2cSJed Brown 344ee1cef2cSJed Brown Collective on Mat 345ee1cef2cSJed Brown 346ee1cef2cSJed Brown Input Parameters: 347ee1cef2cSJed Brown + N - submatrix to update 348ee1cef2cSJed Brown . A - full matrix in the submatrix 349ee1cef2cSJed Brown . isrow - rows in the update (same as the first time the submatrix was created) 350ee1cef2cSJed Brown - iscol - columns in the update (same as the first time the submatrix was created) 351ee1cef2cSJed Brown 352ee1cef2cSJed Brown Level: developer 353ee1cef2cSJed Brown 354ee1cef2cSJed Brown Notes: 355ee1cef2cSJed Brown Most will use MatGetSubMatrix which provides a more efficient representation if it is available. 356ee1cef2cSJed Brown 357ee1cef2cSJed Brown .seealso: MatGetSubMatrix(), MatCreateSubMatrix() 358ee1cef2cSJed Brown @*/ 359ee1cef2cSJed Brown PetscErrorCode PETSCMAT_DLLEXPORT MatSubMatrixUpdate(Mat N,Mat A,IS isrow,IS iscol) 360ee1cef2cSJed Brown { 361ee1cef2cSJed Brown PetscErrorCode ierr; 362ee1cef2cSJed Brown PetscTruth flg; 363ee1cef2cSJed Brown Mat_SubMatrix *Na; 364ee1cef2cSJed Brown 365ee1cef2cSJed Brown PetscFunctionBegin; 366ee1cef2cSJed Brown PetscValidHeaderSpecific(N,MAT_COOKIE,1); 367ee1cef2cSJed Brown PetscValidHeaderSpecific(A,MAT_COOKIE,2); 368ee1cef2cSJed Brown PetscValidHeaderSpecific(isrow,IS_COOKIE,3); 369ee1cef2cSJed Brown PetscValidHeaderSpecific(iscol,IS_COOKIE,4); 370ee1cef2cSJed Brown ierr = PetscTypeCompare((PetscObject)N,MATSUBMATRIX,&flg);CHKERRQ(ierr); 371ee1cef2cSJed Brown if (!flg) SETERRQ(PETSC_ERR_ARG_WRONG,"Matrix has wrong type"); 372ee1cef2cSJed Brown 373ee1cef2cSJed Brown Na = (Mat_SubMatrix*)N->data; 374ee1cef2cSJed Brown ierr = ISEqual(isrow,Na->isrow,&flg);CHKERRQ(ierr); 375ee1cef2cSJed Brown if (!flg) SETERRQ(PETSC_ERR_ARG_INCOMP,"Cannot update submatrix with different row indices"); 376ee1cef2cSJed Brown ierr = ISEqual(iscol,Na->iscol,&flg);CHKERRQ(ierr); 377ee1cef2cSJed Brown if (!flg) SETERRQ(PETSC_ERR_ARG_INCOMP,"Cannot update submatrix with different column indices"); 378ee1cef2cSJed Brown 379ee1cef2cSJed Brown ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 380ee1cef2cSJed Brown ierr = MatDestroy(Na->A);CHKERRQ(ierr); 381ee1cef2cSJed Brown Na->A = A; 382ee1cef2cSJed Brown 383ee1cef2cSJed Brown Na->scale = 1.0; 384ee1cef2cSJed Brown if (Na->left) {ierr = VecDestroy(Na->left);CHKERRQ(ierr);} 385ee1cef2cSJed Brown if (Na->right) {ierr = VecDestroy(Na->right);CHKERRQ(ierr);} 386ee1cef2cSJed Brown PetscFunctionReturn(0); 387ee1cef2cSJed Brown } 388