1ee1cef2cSJed Brown 2af0996ceSBarry Smith #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 lwork,rwork; /* work vectors inside the scatters */ 7*11c5f74dSStefano Zampini Vec lwork2,rwork2; /* work vectors inside the scatters */ 8ee1cef2cSJed Brown VecScatter lrestrict,rprolong; 9ee1cef2cSJed Brown Mat A; 1054e05e6cSHong Zhang } Mat_SubVirtual; 11ee1cef2cSJed Brown 12bddd1ffdSAlp Dener static PetscErrorCode MatScale_SubMatrix(Mat N,PetscScalar a) 13bddd1ffdSAlp Dener { 14bddd1ffdSAlp Dener Mat_SubVirtual *Na = (Mat_SubVirtual*)N->data; 15bddd1ffdSAlp Dener PetscErrorCode ierr; 16bddd1ffdSAlp Dener 17bddd1ffdSAlp Dener PetscFunctionBegin; 18*11c5f74dSStefano Zampini ierr = MatScale(Na->A,a);CHKERRQ(ierr); 19bddd1ffdSAlp Dener PetscFunctionReturn(0); 20bddd1ffdSAlp Dener } 21bddd1ffdSAlp Dener 22bddd1ffdSAlp Dener static PetscErrorCode MatShift_SubMatrix(Mat N,PetscScalar a) 23bddd1ffdSAlp Dener { 24bddd1ffdSAlp Dener Mat_SubVirtual *Na = (Mat_SubVirtual*)N->data; 25bddd1ffdSAlp Dener PetscErrorCode ierr; 26bddd1ffdSAlp Dener 27bddd1ffdSAlp Dener PetscFunctionBegin; 28*11c5f74dSStefano Zampini ierr = MatShift(Na->A,a);CHKERRQ(ierr); 29ee1cef2cSJed Brown PetscFunctionReturn(0); 30ee1cef2cSJed Brown } 31ee1cef2cSJed Brown 32ee1cef2cSJed Brown static PetscErrorCode MatDiagonalScale_SubMatrix(Mat N,Vec left,Vec right) 33ee1cef2cSJed Brown { 3454e05e6cSHong Zhang Mat_SubVirtual *Na = (Mat_SubVirtual*)N->data; 35ee1cef2cSJed Brown PetscErrorCode ierr; 36ee1cef2cSJed Brown 37ee1cef2cSJed Brown PetscFunctionBegin; 38ee1cef2cSJed Brown if (right) { 39*11c5f74dSStefano Zampini ierr = VecZeroEntries(Na->rwork);CHKERRQ(ierr); 40*11c5f74dSStefano Zampini ierr = VecScatterBegin(Na->rprolong,right,Na->rwork,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 41*11c5f74dSStefano Zampini ierr = VecScatterEnd(Na->rprolong,right,Na->rwork,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 42ee1cef2cSJed Brown } 43*11c5f74dSStefano Zampini if (left) { 44*11c5f74dSStefano Zampini ierr = VecZeroEntries(Na->lwork);CHKERRQ(ierr); 45*11c5f74dSStefano Zampini ierr = VecScatterBegin(Na->lrestrict,left,Na->lwork,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 46*11c5f74dSStefano Zampini ierr = VecScatterEnd(Na->lrestrict,left,Na->lwork,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 47ee1cef2cSJed Brown } 48*11c5f74dSStefano Zampini ierr = MatDiagonalScale(Na->A,left ? Na->lwork : NULL,right ? Na->rwork : NULL);CHKERRQ(ierr); 49*11c5f74dSStefano Zampini PetscFunctionReturn(0); 50*11c5f74dSStefano Zampini } 51*11c5f74dSStefano Zampini 52*11c5f74dSStefano Zampini static PetscErrorCode MatGetDiagonal_SubMatrix(Mat N,Vec d) 53*11c5f74dSStefano Zampini { 54*11c5f74dSStefano Zampini Mat_SubVirtual *Na = (Mat_SubVirtual*)N->data; 55*11c5f74dSStefano Zampini PetscErrorCode ierr; 56*11c5f74dSStefano Zampini 57*11c5f74dSStefano Zampini PetscFunctionBegin; 58*11c5f74dSStefano Zampini ierr = MatGetDiagonal(Na->A,Na->rwork);CHKERRQ(ierr); 59*11c5f74dSStefano Zampini ierr = VecScatterBegin(Na->rprolong,Na->rwork,d,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 60*11c5f74dSStefano Zampini ierr = VecScatterEnd(Na->rprolong,Na->rwork,d,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 61ee1cef2cSJed Brown PetscFunctionReturn(0); 62ee1cef2cSJed Brown } 63ee1cef2cSJed Brown 64ee1cef2cSJed Brown static PetscErrorCode MatMult_SubMatrix(Mat N,Vec x,Vec y) 65ee1cef2cSJed Brown { 6654e05e6cSHong Zhang Mat_SubVirtual *Na = (Mat_SubVirtual*)N->data; 67ee1cef2cSJed Brown PetscErrorCode ierr; 68ee1cef2cSJed Brown 69ee1cef2cSJed Brown PetscFunctionBegin; 70ee1cef2cSJed Brown ierr = VecZeroEntries(Na->rwork);CHKERRQ(ierr); 71*11c5f74dSStefano Zampini ierr = VecScatterBegin(Na->rprolong,x,Na->rwork,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 72*11c5f74dSStefano Zampini ierr = VecScatterEnd(Na->rprolong,x,Na->rwork,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 73ee1cef2cSJed Brown ierr = MatMult(Na->A,Na->rwork,Na->lwork);CHKERRQ(ierr); 74ee1cef2cSJed Brown ierr = VecScatterBegin(Na->lrestrict,Na->lwork,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 75ee1cef2cSJed Brown ierr = VecScatterEnd(Na->lrestrict,Na->lwork,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 76ee1cef2cSJed Brown PetscFunctionReturn(0); 77ee1cef2cSJed Brown } 78ee1cef2cSJed Brown 79ee1cef2cSJed Brown static PetscErrorCode MatMultAdd_SubMatrix(Mat N,Vec v1,Vec v2,Vec v3) 80ee1cef2cSJed Brown { 8154e05e6cSHong Zhang Mat_SubVirtual *Na = (Mat_SubVirtual*)N->data; 82ee1cef2cSJed Brown PetscErrorCode ierr; 83ee1cef2cSJed Brown 84ee1cef2cSJed Brown PetscFunctionBegin; 85ee1cef2cSJed Brown ierr = VecZeroEntries(Na->rwork);CHKERRQ(ierr); 86*11c5f74dSStefano Zampini ierr = VecScatterBegin(Na->rprolong,v1,Na->rwork,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 87*11c5f74dSStefano Zampini ierr = VecScatterEnd(Na->rprolong,v1,Na->rwork,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 88*11c5f74dSStefano Zampini if (v1 == v2) { 89*11c5f74dSStefano Zampini ierr = MatMultAdd(Na->A,Na->rwork,Na->rwork,Na->lwork);CHKERRQ(ierr); 90*11c5f74dSStefano Zampini } else if (v2 == v3) { 91*11c5f74dSStefano Zampini ierr = VecZeroEntries(Na->lwork);CHKERRQ(ierr); 92*11c5f74dSStefano Zampini ierr = VecScatterBegin(Na->lrestrict,v2,Na->lwork,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 93*11c5f74dSStefano Zampini ierr = VecScatterEnd(Na->lrestrict,v2,Na->lwork,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 94*11c5f74dSStefano Zampini ierr = MatMultAdd(Na->A,Na->rwork,Na->lwork,Na->lwork);CHKERRQ(ierr); 951e7f0bbdSJed Brown } else { 96*11c5f74dSStefano Zampini if (!Na->lwork2) { 97*11c5f74dSStefano Zampini ierr = VecDuplicate(Na->lwork,&Na->lwork2);CHKERRQ(ierr); 98*11c5f74dSStefano Zampini } else { 99*11c5f74dSStefano Zampini ierr = VecZeroEntries(Na->lwork2);CHKERRQ(ierr); 100*11c5f74dSStefano Zampini } 101*11c5f74dSStefano Zampini ierr = VecScatterBegin(Na->lrestrict,v2,Na->lwork2,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 102*11c5f74dSStefano Zampini ierr = VecScatterEnd(Na->lrestrict,v2,Na->lwork2,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 103*11c5f74dSStefano Zampini ierr = MatMultAdd(Na->A,Na->rwork,Na->lwork2,Na->lwork);CHKERRQ(ierr); 104*11c5f74dSStefano Zampini } 105ee1cef2cSJed Brown ierr = VecScatterBegin(Na->lrestrict,Na->lwork,v3,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 106ee1cef2cSJed Brown ierr = VecScatterEnd(Na->lrestrict,Na->lwork,v3,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 107ee1cef2cSJed Brown PetscFunctionReturn(0); 108ee1cef2cSJed Brown } 109ee1cef2cSJed Brown 110ee1cef2cSJed Brown static PetscErrorCode MatMultTranspose_SubMatrix(Mat N,Vec x,Vec y) 111ee1cef2cSJed Brown { 11254e05e6cSHong Zhang Mat_SubVirtual *Na = (Mat_SubVirtual*)N->data; 113ee1cef2cSJed Brown PetscErrorCode ierr; 114ee1cef2cSJed Brown 115ee1cef2cSJed Brown PetscFunctionBegin; 116ee1cef2cSJed Brown ierr = VecZeroEntries(Na->lwork);CHKERRQ(ierr); 117*11c5f74dSStefano Zampini ierr = VecScatterBegin(Na->lrestrict,x,Na->lwork,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 118*11c5f74dSStefano Zampini ierr = VecScatterEnd(Na->lrestrict,x,Na->lwork,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 119ee1cef2cSJed Brown ierr = MatMultTranspose(Na->A,Na->lwork,Na->rwork);CHKERRQ(ierr); 1201e7f0bbdSJed Brown ierr = VecScatterBegin(Na->rprolong,Na->rwork,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1211e7f0bbdSJed Brown ierr = VecScatterEnd(Na->rprolong,Na->rwork,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 122ee1cef2cSJed Brown PetscFunctionReturn(0); 123ee1cef2cSJed Brown } 124ee1cef2cSJed Brown 125ee1cef2cSJed Brown static PetscErrorCode MatMultTransposeAdd_SubMatrix(Mat N,Vec v1,Vec v2,Vec v3) 126ee1cef2cSJed Brown { 12754e05e6cSHong Zhang Mat_SubVirtual *Na = (Mat_SubVirtual*)N->data; 128ee1cef2cSJed Brown PetscErrorCode ierr; 129ee1cef2cSJed Brown 130ee1cef2cSJed Brown PetscFunctionBegin; 131ee1cef2cSJed Brown ierr = VecZeroEntries(Na->lwork);CHKERRQ(ierr); 132*11c5f74dSStefano Zampini ierr = VecScatterBegin(Na->lrestrict,v1,Na->lwork,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 133*11c5f74dSStefano Zampini ierr = VecScatterEnd(Na->lrestrict,v1,Na->lwork,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 134*11c5f74dSStefano Zampini if (v1 == v2) { 135*11c5f74dSStefano Zampini ierr = MatMultTransposeAdd(Na->A,Na->lwork,Na->lwork,Na->rwork);CHKERRQ(ierr); 136*11c5f74dSStefano Zampini } else if (v2 == v3) { 137*11c5f74dSStefano Zampini ierr = VecZeroEntries(Na->rwork);CHKERRQ(ierr); 138*11c5f74dSStefano Zampini ierr = VecScatterBegin(Na->rprolong,v2,Na->rwork,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 139*11c5f74dSStefano Zampini ierr = VecScatterEnd(Na->rprolong,v2,Na->rwork,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 140*11c5f74dSStefano Zampini ierr = MatMultTransposeAdd(Na->A,Na->lwork,Na->rwork,Na->rwork);CHKERRQ(ierr); 1411e7f0bbdSJed Brown } else { 142*11c5f74dSStefano Zampini if (!Na->rwork2) { 143*11c5f74dSStefano Zampini ierr = VecDuplicate(Na->rwork,&Na->rwork2);CHKERRQ(ierr); 144*11c5f74dSStefano Zampini } else { 145*11c5f74dSStefano Zampini ierr = VecZeroEntries(Na->rwork2);CHKERRQ(ierr); 146*11c5f74dSStefano Zampini } 147*11c5f74dSStefano Zampini ierr = VecScatterBegin(Na->rprolong,v2,Na->rwork2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 148*11c5f74dSStefano Zampini ierr = VecScatterEnd(Na->rprolong,v2,Na->rwork2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 149*11c5f74dSStefano Zampini ierr = MatMultTransposeAdd(Na->A,Na->lwork,Na->rwork2,Na->rwork);CHKERRQ(ierr); 150*11c5f74dSStefano Zampini } 1511e7f0bbdSJed Brown ierr = VecScatterBegin(Na->rprolong,Na->rwork,v3,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 1521e7f0bbdSJed Brown ierr = VecScatterEnd(Na->rprolong,Na->rwork,v3,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 153ee1cef2cSJed Brown PetscFunctionReturn(0); 154ee1cef2cSJed Brown } 155ee1cef2cSJed Brown 156ee1cef2cSJed Brown static PetscErrorCode MatDestroy_SubMatrix(Mat N) 157ee1cef2cSJed Brown { 15854e05e6cSHong Zhang Mat_SubVirtual *Na = (Mat_SubVirtual*)N->data; 159ee1cef2cSJed Brown PetscErrorCode ierr; 160ee1cef2cSJed Brown 161ee1cef2cSJed Brown PetscFunctionBegin; 1626bf464f9SBarry Smith ierr = ISDestroy(&Na->isrow);CHKERRQ(ierr); 1636bf464f9SBarry Smith ierr = ISDestroy(&Na->iscol);CHKERRQ(ierr); 1646bf464f9SBarry Smith ierr = VecDestroy(&Na->lwork);CHKERRQ(ierr); 1656bf464f9SBarry Smith ierr = VecDestroy(&Na->rwork);CHKERRQ(ierr); 166*11c5f74dSStefano Zampini ierr = VecDestroy(&Na->lwork2);CHKERRQ(ierr); 167*11c5f74dSStefano Zampini ierr = VecDestroy(&Na->rwork2);CHKERRQ(ierr); 1686bf464f9SBarry Smith ierr = VecScatterDestroy(&Na->lrestrict);CHKERRQ(ierr); 1696bf464f9SBarry Smith ierr = VecScatterDestroy(&Na->rprolong);CHKERRQ(ierr); 1706bf464f9SBarry Smith ierr = MatDestroy(&Na->A);CHKERRQ(ierr); 1716bf464f9SBarry Smith ierr = PetscFree(N->data);CHKERRQ(ierr); 172ee1cef2cSJed Brown PetscFunctionReturn(0); 173ee1cef2cSJed Brown } 174ee1cef2cSJed Brown 175ee1cef2cSJed Brown /*@ 17654e05e6cSHong Zhang MatCreateSubMatrixVirtual - Creates a virtual matrix that acts as a submatrix 177ee1cef2cSJed Brown 178ee1cef2cSJed Brown Collective on Mat 179ee1cef2cSJed Brown 180ee1cef2cSJed Brown Input Parameters: 181ee1cef2cSJed Brown + A - matrix that we will extract a submatrix of 182ee1cef2cSJed Brown . isrow - rows to be present in the submatrix 183ee1cef2cSJed Brown - iscol - columns to be present in the submatrix 184ee1cef2cSJed Brown 185ee1cef2cSJed Brown Output Parameters: 186ee1cef2cSJed Brown . newmat - new matrix 187ee1cef2cSJed Brown 188ee1cef2cSJed Brown Level: developer 189ee1cef2cSJed Brown 190ee1cef2cSJed Brown Notes: 1917dae84e0SHong Zhang Most will use MatCreateSubMatrix which provides a more efficient representation if it is available. 192ee1cef2cSJed Brown 19354e05e6cSHong Zhang .seealso: MatCreateSubMatrix(), MatSubMatrixVirtualUpdate() 194ee1cef2cSJed Brown @*/ 19554e05e6cSHong Zhang PetscErrorCode MatCreateSubMatrixVirtual(Mat A,IS isrow,IS iscol,Mat *newmat) 196ee1cef2cSJed Brown { 197ee1cef2cSJed Brown Vec left,right; 198ee1cef2cSJed Brown PetscInt m,n; 199ee1cef2cSJed Brown Mat N; 20054e05e6cSHong Zhang Mat_SubVirtual *Na; 201ee1cef2cSJed Brown PetscErrorCode ierr; 202ee1cef2cSJed Brown 203ee1cef2cSJed Brown PetscFunctionBegin; 2040700a824SBarry Smith PetscValidHeaderSpecific(A,MAT_CLASSID,1); 2050700a824SBarry Smith PetscValidHeaderSpecific(isrow,IS_CLASSID,2); 2060700a824SBarry Smith PetscValidHeaderSpecific(iscol,IS_CLASSID,3); 207ee1cef2cSJed Brown PetscValidPointer(newmat,4); 208ee1cef2cSJed Brown *newmat = 0; 209ee1cef2cSJed Brown 210ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)A),&N);CHKERRQ(ierr); 211ee1cef2cSJed Brown ierr = ISGetLocalSize(isrow,&m);CHKERRQ(ierr); 212ee1cef2cSJed Brown ierr = ISGetLocalSize(iscol,&n);CHKERRQ(ierr); 213ee1cef2cSJed Brown ierr = MatSetSizes(N,m,n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 214ee1cef2cSJed Brown ierr = PetscObjectChangeTypeName((PetscObject)N,MATSUBMATRIX);CHKERRQ(ierr); 215ee1cef2cSJed Brown 216b00a9115SJed Brown ierr = PetscNewLog(N,&Na);CHKERRQ(ierr); 217ee1cef2cSJed Brown N->data = (void*)Na; 218*11c5f74dSStefano Zampini 219ee1cef2cSJed Brown ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr); 220ee1cef2cSJed Brown ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr); 221ee1cef2cSJed Brown Na->isrow = isrow; 222ee1cef2cSJed Brown Na->iscol = iscol; 223*11c5f74dSStefano Zampini 224*11c5f74dSStefano Zampini ierr = PetscFree(N->defaultvectype);CHKERRQ(ierr); 225*11c5f74dSStefano Zampini ierr = PetscStrallocpy(A->defaultvectype,&N->defaultvectype);CHKERRQ(ierr); 226*11c5f74dSStefano Zampini /* Do not use MatConvert directly since MatShell has a duplicate operation which does not increase 227*11c5f74dSStefano Zampini the reference count of the context. This is a problem if A is already of type MATSHELL */ 228*11c5f74dSStefano Zampini ierr = MatConvertFrom_Shell(A,MATSHELL,MAT_INITIAL_MATRIX,&Na->A);CHKERRQ(ierr); 229ee1cef2cSJed Brown 230ee1cef2cSJed Brown N->ops->destroy = MatDestroy_SubMatrix; 231ee1cef2cSJed Brown N->ops->mult = MatMult_SubMatrix; 232ee1cef2cSJed Brown N->ops->multadd = MatMultAdd_SubMatrix; 233ee1cef2cSJed Brown N->ops->multtranspose = MatMultTranspose_SubMatrix; 234ee1cef2cSJed Brown N->ops->multtransposeadd = MatMultTransposeAdd_SubMatrix; 235ee1cef2cSJed Brown N->ops->scale = MatScale_SubMatrix; 236ee1cef2cSJed Brown N->ops->diagonalscale = MatDiagonalScale_SubMatrix; 237bddd1ffdSAlp Dener N->ops->shift = MatShift_SubMatrix; 238*11c5f74dSStefano Zampini N->ops->convert = MatConvert_Shell; 239*11c5f74dSStefano Zampini N->ops->getdiagonal = MatGetDiagonal_SubMatrix; 240ee1cef2cSJed Brown 24133d57670SJed Brown ierr = MatSetBlockSizesFromMats(N,A,A);CHKERRQ(ierr); 24226283091SBarry Smith ierr = PetscLayoutSetUp(N->rmap);CHKERRQ(ierr); 24326283091SBarry Smith ierr = PetscLayoutSetUp(N->cmap);CHKERRQ(ierr); 244ee1cef2cSJed Brown 2452a7a6963SBarry Smith ierr = MatCreateVecs(A,&Na->rwork,&Na->lwork);CHKERRQ(ierr); 246*11c5f74dSStefano Zampini ierr = MatCreateVecs(N,&right,&left);CHKERRQ(ierr); 2479448b7f1SJunchao Zhang ierr = VecScatterCreate(Na->lwork,isrow,left,NULL,&Na->lrestrict);CHKERRQ(ierr); 2489448b7f1SJunchao Zhang ierr = VecScatterCreate(right,NULL,Na->rwork,iscol,&Na->rprolong);CHKERRQ(ierr); 2496bf464f9SBarry Smith ierr = VecDestroy(&left);CHKERRQ(ierr); 2506bf464f9SBarry Smith ierr = VecDestroy(&right);CHKERRQ(ierr); 251be7c243fSJed Brown ierr = MatSetUp(N);CHKERRQ(ierr); 25226fbe8dcSKarl Rupp 253*11c5f74dSStefano Zampini N->assembled = PETSC_TRUE; 254ee1cef2cSJed Brown *newmat = N; 255ee1cef2cSJed Brown PetscFunctionReturn(0); 256ee1cef2cSJed Brown } 257ee1cef2cSJed Brown 258ee1cef2cSJed Brown 259ee1cef2cSJed Brown /*@ 26054e05e6cSHong Zhang MatSubMatrixVirtualUpdate - Updates a submatrix 261ee1cef2cSJed Brown 262ee1cef2cSJed Brown Collective on Mat 263ee1cef2cSJed Brown 264ee1cef2cSJed Brown Input Parameters: 265ee1cef2cSJed Brown + N - submatrix to update 266ee1cef2cSJed Brown . A - full matrix in the submatrix 267ee1cef2cSJed Brown . isrow - rows in the update (same as the first time the submatrix was created) 268ee1cef2cSJed Brown - iscol - columns in the update (same as the first time the submatrix was created) 269ee1cef2cSJed Brown 270ee1cef2cSJed Brown Level: developer 271ee1cef2cSJed Brown 272ee1cef2cSJed Brown Notes: 2737dae84e0SHong Zhang Most will use MatCreateSubMatrix which provides a more efficient representation if it is available. 274ee1cef2cSJed Brown 27554e05e6cSHong Zhang .seealso: MatCreateSubMatrixVirtual() 276ee1cef2cSJed Brown @*/ 27754e05e6cSHong Zhang PetscErrorCode MatSubMatrixVirtualUpdate(Mat N,Mat A,IS isrow,IS iscol) 278ee1cef2cSJed Brown { 279ee1cef2cSJed Brown PetscErrorCode ierr; 280ace3abfcSBarry Smith PetscBool flg; 28154e05e6cSHong Zhang Mat_SubVirtual *Na; 282ee1cef2cSJed Brown 283ee1cef2cSJed Brown PetscFunctionBegin; 2840700a824SBarry Smith PetscValidHeaderSpecific(N,MAT_CLASSID,1); 2850700a824SBarry Smith PetscValidHeaderSpecific(A,MAT_CLASSID,2); 2860700a824SBarry Smith PetscValidHeaderSpecific(isrow,IS_CLASSID,3); 2870700a824SBarry Smith PetscValidHeaderSpecific(iscol,IS_CLASSID,4); 288251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)N,MATSUBMATRIX,&flg);CHKERRQ(ierr); 289ce94432eSBarry Smith if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix has wrong type"); 290ee1cef2cSJed Brown 29154e05e6cSHong Zhang Na = (Mat_SubVirtual*)N->data; 292ee1cef2cSJed Brown ierr = ISEqual(isrow,Na->isrow,&flg);CHKERRQ(ierr); 293e32f2f54SBarry Smith if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Cannot update submatrix with different row indices"); 294ee1cef2cSJed Brown ierr = ISEqual(iscol,Na->iscol,&flg);CHKERRQ(ierr); 295e32f2f54SBarry Smith if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Cannot update submatrix with different column indices"); 296ee1cef2cSJed Brown 297*11c5f74dSStefano Zampini ierr = PetscFree(N->defaultvectype);CHKERRQ(ierr); 298*11c5f74dSStefano Zampini ierr = PetscStrallocpy(A->defaultvectype,&N->defaultvectype);CHKERRQ(ierr); 2996bf464f9SBarry Smith ierr = MatDestroy(&Na->A);CHKERRQ(ierr); 300*11c5f74dSStefano Zampini /* Do not use MatConvert directly since MatShell has a duplicate operation which does not increase 301*11c5f74dSStefano Zampini the reference count of the context. This is a problem if A is already of type MATSHELL */ 302*11c5f74dSStefano Zampini ierr = MatConvertFrom_Shell(A,MATSHELL,MAT_INITIAL_MATRIX,&Na->A);CHKERRQ(ierr); 303ee1cef2cSJed Brown PetscFunctionReturn(0); 304ee1cef2cSJed Brown } 305