1be1d678aSKris Buschelman 28016bdd1SSatish Balay /* 3d9653453SSatish Balay Support for the parallel BAIJ matrix vector multiply 48016bdd1SSatish Balay */ 5c6db04a5SJed Brown #include <../src/mat/impls/baij/mpi/mpibaij.h> 6af0996ceSBarry Smith #include <petsc/private/isimpl.h> /* needed because accesses data structure of ISLocalToGlobalMapping directly */ 7bba1ac68SSatish Balay 8dfbe8321SBarry Smith PetscErrorCode MatSetUpMultiply_MPIBAIJ(Mat mat) 98016bdd1SSatish Balay { 10d9653453SSatish Balay Mat_MPIBAIJ *baij = (Mat_MPIBAIJ*)mat->data; 11d9653453SSatish Balay Mat_SeqBAIJ *B = (Mat_SeqBAIJ*)(baij->B->data); 126849ba73SBarry Smith PetscErrorCode ierr; 13b24ad042SBarry Smith PetscInt i,j,*aj = B->j,ec = 0,*garray; 14d0f46423SBarry Smith PetscInt bs = mat->rmap->bs,*stmp; 158016bdd1SSatish Balay IS from,to; 168016bdd1SSatish Balay Vec gvec; 17aa482453SBarry Smith #if defined(PETSC_USE_CTABLE) 180f5bd95cSBarry Smith PetscTable gid1_lid1; 190f5bd95cSBarry Smith PetscTablePosition tpos; 20b24ad042SBarry Smith PetscInt gid,lid; 216f531f54SSatish Balay #else 22b24ad042SBarry Smith PetscInt Nbs = baij->Nbs,*indices; 2373a2e727SSatish Balay #endif 248016bdd1SSatish Balay 253a40ed3dSBarry Smith PetscFunctionBegin; 26aa482453SBarry Smith #if defined(PETSC_USE_CTABLE) 2773a2e727SSatish Balay /* use a table - Mark Adams */ 28e23dfa41SBarry Smith ierr = PetscTableCreate(B->mbs,baij->Nbs+1,&gid1_lid1);CHKERRQ(ierr); 2973a2e727SSatish Balay for (i=0; i<B->mbs; i++) { 3073a2e727SSatish Balay for (j=0; j<B->ilen[i]; j++) { 31b24ad042SBarry Smith PetscInt data,gid1 = aj[B->i[i]+j] + 1; 320f5bd95cSBarry Smith ierr = PetscTableFind(gid1_lid1,gid1,&data);CHKERRQ(ierr); 33fa46199cSSatish Balay if (!data) { 3473a2e727SSatish Balay /* one based table */ 353861aac3SJed Brown ierr = PetscTableAdd(gid1_lid1,gid1,++ec,INSERT_VALUES);CHKERRQ(ierr); 3673a2e727SSatish Balay } 3773a2e727SSatish Balay } 3873a2e727SSatish Balay } 3973a2e727SSatish Balay /* form array of columns we need */ 40854ce69bSBarry Smith ierr = PetscMalloc1(ec+1,&garray);CHKERRQ(ierr); 410f5bd95cSBarry Smith ierr = PetscTableGetHeadPosition(gid1_lid1,&tpos);CHKERRQ(ierr); 4273a2e727SSatish Balay while (tpos) { 430f5bd95cSBarry Smith ierr = PetscTableGetNext(gid1_lid1,&tpos,&gid,&lid);CHKERRQ(ierr); 4473a2e727SSatish Balay gid--; lid--; 4573a2e727SSatish Balay garray[lid] = gid; 4673a2e727SSatish Balay } 470064e2bbSSatish Balay ierr = PetscSortInt(ec,garray);CHKERRQ(ierr); 480f5bd95cSBarry Smith ierr = PetscTableRemoveAll(gid1_lid1);CHKERRQ(ierr); 4973a2e727SSatish Balay for (i=0; i<ec; i++) { 503861aac3SJed Brown ierr = PetscTableAdd(gid1_lid1,garray[i]+1,i+1,INSERT_VALUES);CHKERRQ(ierr); 5173a2e727SSatish Balay } 5273a2e727SSatish Balay /* compact out the extra columns in B */ 5373a2e727SSatish Balay for (i=0; i<B->mbs; i++) { 5473a2e727SSatish Balay for (j=0; j<B->ilen[i]; j++) { 55b24ad042SBarry Smith PetscInt gid1 = aj[B->i[i] + j] + 1; 560f5bd95cSBarry Smith ierr = PetscTableFind(gid1_lid1,gid1,&lid);CHKERRQ(ierr); 57fa46199cSSatish Balay lid--; 5873a2e727SSatish Balay aj[B->i[i]+j] = lid; 5973a2e727SSatish Balay } 6073a2e727SSatish Balay } 6173a2e727SSatish Balay B->nbs = ec; 62*ca5434daSLawrence Mitchell ierr = PetscLayoutDestroy(&baij->B->cmap);CHKERRQ(ierr); 63*ca5434daSLawrence Mitchell ierr = PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)baij->B),ec*mat->rmap->bs,ec*mat->rmap->bs,mat->rmap->bs,&baij->B->cmap);CHKERRQ(ierr); 646bc0bbbfSBarry Smith ierr = PetscTableDestroy(&gid1_lid1);CHKERRQ(ierr); 6573a2e727SSatish Balay #else 66435da068SBarry Smith /* Make an array as long as the number of columns */ 67d9653453SSatish Balay /* mark those columns that are in baij->B */ 681795a4d1SJed Brown ierr = PetscCalloc1(Nbs+1,&indices);CHKERRQ(ierr); 69d9653453SSatish Balay for (i=0; i<B->mbs; i++) { 708016bdd1SSatish Balay for (j=0; j<B->ilen[i]; j++) { 71d9653453SSatish Balay if (!indices[aj[B->i[i] + j]]) ec++; 72d9653453SSatish Balay indices[aj[B->i[i] + j]] = 1; 738016bdd1SSatish Balay } 748016bdd1SSatish Balay } 758016bdd1SSatish Balay 768016bdd1SSatish Balay /* form array of columns we need */ 77854ce69bSBarry Smith ierr = PetscMalloc1(ec+1,&garray);CHKERRQ(ierr); 788016bdd1SSatish Balay ec = 0; 79d9653453SSatish Balay for (i=0; i<Nbs; i++) { 800bdbc534SSatish Balay if (indices[i]) { 810bdbc534SSatish Balay garray[ec++] = i; 820bdbc534SSatish Balay } 838016bdd1SSatish Balay } 848016bdd1SSatish Balay 858016bdd1SSatish Balay /* make indices now point into garray */ 868016bdd1SSatish Balay for (i=0; i<ec; i++) { 87d9653453SSatish Balay indices[garray[i]] = i; 888016bdd1SSatish Balay } 898016bdd1SSatish Balay 908016bdd1SSatish Balay /* compact out the extra columns in B */ 91d9653453SSatish Balay for (i=0; i<B->mbs; i++) { 928016bdd1SSatish Balay for (j=0; j<B->ilen[i]; j++) { 93d9653453SSatish Balay aj[B->i[i] + j] = indices[aj[B->i[i] + j]]; 948016bdd1SSatish Balay } 958016bdd1SSatish Balay } 96d9653453SSatish Balay B->nbs = ec; 97*ca5434daSLawrence Mitchell ierr = PetscLayoutDestroy(&baij->B->cmap);CHKERRQ(ierr); 98*ca5434daSLawrence Mitchell ierr = PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)baij->B),ec*mat->rmap->bs,ec*mat->rmap->bs,mat->rmap->bs,&baij->B->cmap);CHKERRQ(ierr); 99606d414cSSatish Balay ierr = PetscFree(indices);CHKERRQ(ierr); 10073a2e727SSatish Balay #endif 1018016bdd1SSatish Balay 1028016bdd1SSatish Balay /* create local vector that is used to scatter into */ 103029af93fSBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,ec*bs,&baij->lvec);CHKERRQ(ierr); 1048016bdd1SSatish Balay 105c16cb8f2SBarry Smith /* create two temporary index sets for building scatter-gather */ 106deff0451SBarry Smith ierr = ISCreateBlock(PETSC_COMM_SELF,bs,ec,garray,PETSC_COPY_VALUES,&from);CHKERRQ(ierr); 107c16cb8f2SBarry Smith 108854ce69bSBarry Smith ierr = PetscMalloc1(ec+1,&stmp);CHKERRQ(ierr); 10926fbe8dcSKarl Rupp for (i=0; i<ec; i++) stmp[i] = i; 110deff0451SBarry Smith ierr = ISCreateBlock(PETSC_COMM_SELF,bs,ec,stmp,PETSC_OWN_POINTER,&to);CHKERRQ(ierr); 1118016bdd1SSatish Balay 1128016bdd1SSatish Balay /* create temporary global vector to generate scatter context */ 113ce94432eSBarry Smith ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)mat),1,mat->cmap->n,mat->cmap->N,NULL,&gvec);CHKERRQ(ierr); 1148016bdd1SSatish Balay 1159448b7f1SJunchao Zhang ierr = VecScatterCreate(gvec,from,baij->lvec,to,&baij->Mvctx);CHKERRQ(ierr); 11690f02eecSBarry Smith 1173bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)baij->Mvctx);CHKERRQ(ierr); 1183bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)baij->lvec);CHKERRQ(ierr); 1193bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)from);CHKERRQ(ierr); 1203bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)to);CHKERRQ(ierr); 12126fbe8dcSKarl Rupp 122d9653453SSatish Balay baij->garray = garray; 12326fbe8dcSKarl Rupp 1243bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)mat,(ec+1)*sizeof(PetscInt));CHKERRQ(ierr); 1256bf464f9SBarry Smith ierr = ISDestroy(&from);CHKERRQ(ierr); 1266bf464f9SBarry Smith ierr = ISDestroy(&to);CHKERRQ(ierr); 1276bf464f9SBarry Smith ierr = VecDestroy(&gvec);CHKERRQ(ierr); 1283a40ed3dSBarry Smith PetscFunctionReturn(0); 1298016bdd1SSatish Balay } 1308016bdd1SSatish Balay 1318016bdd1SSatish Balay /* 132d9653453SSatish Balay Takes the local part of an already assembled MPIBAIJ matrix 1338016bdd1SSatish Balay and disassembles it. This is to allow new nonzeros into the matrix 1348016bdd1SSatish Balay that require more communication in the matrix vector multiply. 1358016bdd1SSatish Balay Thus certain data-structures must be rebuilt. 1368016bdd1SSatish Balay 1378016bdd1SSatish Balay Kind of slow! But that's what application programmers get when 1388016bdd1SSatish Balay they are sloppy. 1398016bdd1SSatish Balay */ 140ab9863d7SBarry Smith PetscErrorCode MatDisAssemble_MPIBAIJ(Mat A) 1418016bdd1SSatish Balay { 142d9653453SSatish Balay Mat_MPIBAIJ *baij = (Mat_MPIBAIJ*)A->data; 143d9653453SSatish Balay Mat B = baij->B,Bnew; 144d9653453SSatish Balay Mat_SeqBAIJ *Bbaij = (Mat_SeqBAIJ*)B->data; 145dfbe8321SBarry Smith PetscErrorCode ierr; 146d0f46423SBarry Smith PetscInt i,j,mbs=Bbaij->mbs,n = A->cmap->N,col,*garray=baij->garray; 147d0f46423SBarry Smith PetscInt bs2 = baij->bs2,*nz,ec,m = A->rmap->n; 1483eda8832SBarry Smith MatScalar *a = Bbaij->a; 149dd6ea824SBarry Smith MatScalar *atmp; 15097e5c40aSBarry Smith 1518016bdd1SSatish Balay 1523a40ed3dSBarry Smith PetscFunctionBegin; 1538016bdd1SSatish Balay /* free stuff related to matrix-vec multiply */ 154b0a32e0cSBarry Smith ierr = VecGetSize(baij->lvec,&ec);CHKERRQ(ierr); /* needed for PetscLogObjectMemory below */ 1556bf464f9SBarry Smith ierr = VecDestroy(&baij->lvec);CHKERRQ(ierr); baij->lvec = 0; 1566bf464f9SBarry Smith ierr = VecScatterDestroy(&baij->Mvctx);CHKERRQ(ierr); baij->Mvctx = 0; 157d9653453SSatish Balay if (baij->colmap) { 158aa482453SBarry Smith #if defined(PETSC_USE_CTABLE) 1596bc0bbbfSBarry Smith ierr = PetscTableDestroy(&baij->colmap);CHKERRQ(ierr); 16048e59246SSatish Balay #else 161606d414cSSatish Balay ierr = PetscFree(baij->colmap);CHKERRQ(ierr); 1623bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)A,-Bbaij->nbs*sizeof(PetscInt));CHKERRQ(ierr); 16348e59246SSatish Balay #endif 1648016bdd1SSatish Balay } 1658016bdd1SSatish Balay 1668016bdd1SSatish Balay /* make sure that B is assembled so we can access its values */ 1676d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1683eda8832SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1698016bdd1SSatish Balay 1708016bdd1SSatish Balay /* invent new B and copy stuff over */ 171785e854fSJed Brown ierr = PetscMalloc1(mbs,&nz);CHKERRQ(ierr); 172d9653453SSatish Balay for (i=0; i<mbs; i++) { 173d9653453SSatish Balay nz[i] = Bbaij->i[i+1]-Bbaij->i[i]; 1748016bdd1SSatish Balay } 175ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)B),&Bnew);CHKERRQ(ierr); 176f69a0ea3SMatthew Knepley ierr = MatSetSizes(Bnew,m,n,m,n);CHKERRQ(ierr); 1777adad957SLisandro Dalcin ierr = MatSetType(Bnew,((PetscObject)B)->type_name);CHKERRQ(ierr); 178d0f46423SBarry Smith ierr = MatSeqBAIJSetPreallocation(Bnew,B->rmap->bs,0,nz);CHKERRQ(ierr); 179b38c15b3SStefano Zampini if (Bbaij->nonew >= 0) { /* Inherit insertion error options (if positive). */ 180b38c15b3SStefano Zampini ((Mat_SeqBAIJ*)Bnew->data)->nonew = Bbaij->nonew; 181b38c15b3SStefano Zampini } 18226fbe8dcSKarl Rupp 1834e0d8c25SBarry Smith ierr = MatSetOption(Bnew,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr); 18477341eacSDmitry Karpeev /* 18577341eacSDmitry Karpeev Ensure that B's nonzerostate is monotonically increasing. 18677341eacSDmitry Karpeev Or should this follow the MatSetValuesBlocked() loop to preserve B's nonzerstate across a MatDisAssemble() call? 18777341eacSDmitry Karpeev */ 18877341eacSDmitry Karpeev Bnew->nonzerostate = B->nonzerostate; 189d9653453SSatish Balay 190bba1ac68SSatish Balay for (i=0; i<mbs; i++) { 191bba1ac68SSatish Balay for (j=Bbaij->i[i]; j<Bbaij->i[i+1]; j++) { 192bba1ac68SSatish Balay col = garray[Bbaij->j[j]]; 1933eda8832SBarry Smith atmp = a + j*bs2; 194bba1ac68SSatish Balay ierr = MatSetValuesBlocked_SeqBAIJ(Bnew,1,&i,1,&col,atmp,B->insertmode);CHKERRQ(ierr); 1958016bdd1SSatish Balay } 1968016bdd1SSatish Balay } 1974e0d8c25SBarry Smith ierr = MatSetOption(Bnew,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr); 198bba1ac68SSatish Balay 199bba1ac68SSatish Balay ierr = PetscFree(nz);CHKERRQ(ierr); 200606d414cSSatish Balay ierr = PetscFree(baij->garray);CHKERRQ(ierr); 2013bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)A,-ec*sizeof(PetscInt));CHKERRQ(ierr); 2026bf464f9SBarry Smith ierr = MatDestroy(&B);CHKERRQ(ierr); 2033bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)A,(PetscObject)Bnew);CHKERRQ(ierr); 20426fbe8dcSKarl Rupp 205d9653453SSatish Balay baij->B = Bnew; 2068016bdd1SSatish Balay A->was_assembled = PETSC_FALSE; 2076a719282SBarry Smith A->assembled = PETSC_FALSE; 2083a40ed3dSBarry Smith PetscFunctionReturn(0); 2098016bdd1SSatish Balay } 2108016bdd1SSatish Balay 21104f1ad80SBarry Smith /* ugly stuff added for Glenn someday we should fix this up */ 21204f1ad80SBarry Smith 21326fbe8dcSKarl Rupp static PetscInt *uglyrmapd = 0,*uglyrmapo = 0; /* mapping from the local ordering to the "diagonal" and "off-diagonal" parts of the local matrix */ 2142cd6534aSBarry Smith static Vec uglydd = 0,uglyoo = 0; /* work vectors used to scale the two parts of the local matrix */ 21504f1ad80SBarry Smith 21604f1ad80SBarry Smith 217dfbe8321SBarry Smith PetscErrorCode MatMPIBAIJDiagonalScaleLocalSetUp(Mat inA,Vec scale) 21804f1ad80SBarry Smith { 21904f1ad80SBarry Smith Mat_MPIBAIJ *ina = (Mat_MPIBAIJ*) inA->data; /*access private part of matrix */ 22004f1ad80SBarry Smith Mat_SeqBAIJ *B = (Mat_SeqBAIJ*)ina->B->data; 221dfbe8321SBarry Smith PetscErrorCode ierr; 222d0f46423SBarry Smith PetscInt bs = inA->rmap->bs,i,n,nt,j,cstart,cend,no,*garray = ina->garray,*lindices; 223b24ad042SBarry Smith PetscInt *r_rmapd,*r_rmapo; 22404f1ad80SBarry Smith 22504f1ad80SBarry Smith PetscFunctionBegin; 22604f1ad80SBarry Smith ierr = MatGetOwnershipRange(inA,&cstart,&cend);CHKERRQ(ierr); 2270298fd71SBarry Smith ierr = MatGetSize(ina->A,NULL,&n);CHKERRQ(ierr); 228854ce69bSBarry Smith ierr = PetscCalloc1(inA->rmap->mapping->n+1,&r_rmapd);CHKERRQ(ierr); 22904f1ad80SBarry Smith nt = 0; 23045b6f7e9SBarry Smith for (i=0; i<inA->rmap->mapping->n; i++) { 23145b6f7e9SBarry Smith if (inA->rmap->mapping->indices[i]*bs >= cstart && inA->rmap->mapping->indices[i]*bs < cend) { 23204f1ad80SBarry Smith nt++; 23345b6f7e9SBarry Smith r_rmapd[i] = inA->rmap->mapping->indices[i] + 1; 23404f1ad80SBarry Smith } 23504f1ad80SBarry Smith } 236e32f2f54SBarry Smith if (nt*bs != n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hmm nt*bs %D n %D",nt*bs,n); 237854ce69bSBarry Smith ierr = PetscMalloc1(n+1,&uglyrmapd);CHKERRQ(ierr); 23845b6f7e9SBarry Smith for (i=0; i<inA->rmap->mapping->n; i++) { 23904f1ad80SBarry Smith if (r_rmapd[i]) { 24004f1ad80SBarry Smith for (j=0; j<bs; j++) { 24104f1ad80SBarry Smith uglyrmapd[(r_rmapd[i]-1)*bs+j-cstart] = i*bs + j; 24204f1ad80SBarry Smith } 24304f1ad80SBarry Smith } 24404f1ad80SBarry Smith } 24504f1ad80SBarry Smith ierr = PetscFree(r_rmapd);CHKERRQ(ierr); 24604f1ad80SBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,n,&uglydd);CHKERRQ(ierr); 24704f1ad80SBarry Smith 248854ce69bSBarry Smith ierr = PetscCalloc1(ina->Nbs+1,&lindices);CHKERRQ(ierr); 24904f1ad80SBarry Smith for (i=0; i<B->nbs; i++) { 25004f1ad80SBarry Smith lindices[garray[i]] = i+1; 25104f1ad80SBarry Smith } 25245b6f7e9SBarry Smith no = inA->rmap->mapping->n - nt; 253854ce69bSBarry Smith ierr = PetscCalloc1(inA->rmap->mapping->n+1,&r_rmapo);CHKERRQ(ierr); 25404f1ad80SBarry Smith nt = 0; 25545b6f7e9SBarry Smith for (i=0; i<inA->rmap->mapping->n; i++) { 25645b6f7e9SBarry Smith if (lindices[inA->rmap->mapping->indices[i]]) { 25704f1ad80SBarry Smith nt++; 25845b6f7e9SBarry Smith r_rmapo[i] = lindices[inA->rmap->mapping->indices[i]]; 25904f1ad80SBarry Smith } 26004f1ad80SBarry Smith } 261e32f2f54SBarry Smith if (nt > no) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hmm nt %D no %D",nt,n); 26204f1ad80SBarry Smith ierr = PetscFree(lindices);CHKERRQ(ierr); 263854ce69bSBarry Smith ierr = PetscMalloc1(nt*bs+1,&uglyrmapo);CHKERRQ(ierr); 26445b6f7e9SBarry Smith for (i=0; i<inA->rmap->mapping->n; i++) { 26504f1ad80SBarry Smith if (r_rmapo[i]) { 26604f1ad80SBarry Smith for (j=0; j<bs; j++) { 26704f1ad80SBarry Smith uglyrmapo[(r_rmapo[i]-1)*bs+j] = i*bs + j; 26804f1ad80SBarry Smith } 26904f1ad80SBarry Smith } 27004f1ad80SBarry Smith } 27104f1ad80SBarry Smith ierr = PetscFree(r_rmapo);CHKERRQ(ierr); 27204f1ad80SBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,nt*bs,&uglyoo);CHKERRQ(ierr); 27304f1ad80SBarry Smith PetscFunctionReturn(0); 27404f1ad80SBarry Smith } 27504f1ad80SBarry Smith 2767087cfbeSBarry Smith PetscErrorCode MatMPIBAIJDiagonalScaleLocal(Mat A,Vec scale) 27704f1ad80SBarry Smith { 27892b32695SKris Buschelman /* This routine should really be abandoned as it duplicates MatDiagonalScaleLocal */ 2794ac538c5SBarry Smith PetscErrorCode ierr; 28092b32695SKris Buschelman 28192b32695SKris Buschelman PetscFunctionBegin; 2824ac538c5SBarry Smith ierr = PetscTryMethod(A,"MatDiagonalScaleLocal_C",(Mat,Vec),(A,scale));CHKERRQ(ierr); 28392b32695SKris Buschelman PetscFunctionReturn(0); 28492b32695SKris Buschelman } 28592b32695SKris Buschelman 2867087cfbeSBarry Smith PetscErrorCode MatDiagonalScaleLocal_MPIBAIJ(Mat A,Vec scale) 28792b32695SKris Buschelman { 28804f1ad80SBarry Smith Mat_MPIBAIJ *a = (Mat_MPIBAIJ*) A->data; /*access private part of matrix */ 289dfbe8321SBarry Smith PetscErrorCode ierr; 290b24ad042SBarry Smith PetscInt n,i; 291bca11509SBarry Smith PetscScalar *d,*o; 292bca11509SBarry Smith const PetscScalar *s; 29304f1ad80SBarry Smith 29404f1ad80SBarry Smith PetscFunctionBegin; 2952cd6534aSBarry Smith if (!uglyrmapd) { 2962cd6534aSBarry Smith ierr = MatMPIBAIJDiagonalScaleLocalSetUp(A,scale);CHKERRQ(ierr); 2972cd6534aSBarry Smith } 2982cd6534aSBarry Smith 299bca11509SBarry Smith ierr = VecGetArrayRead(scale,&s);CHKERRQ(ierr); 30004f1ad80SBarry Smith 30104f1ad80SBarry Smith ierr = VecGetLocalSize(uglydd,&n);CHKERRQ(ierr); 3021ebc52fbSHong Zhang ierr = VecGetArray(uglydd,&d);CHKERRQ(ierr); 30304f1ad80SBarry Smith for (i=0; i<n; i++) { 30404f1ad80SBarry Smith d[i] = s[uglyrmapd[i]]; /* copy "diagonal" (true local) portion of scale into dd vector */ 30504f1ad80SBarry Smith } 3061ebc52fbSHong Zhang ierr = VecRestoreArray(uglydd,&d);CHKERRQ(ierr); 30704f1ad80SBarry Smith /* column scale "diagonal" portion of local matrix */ 3080298fd71SBarry Smith ierr = MatDiagonalScale(a->A,NULL,uglydd);CHKERRQ(ierr); 30904f1ad80SBarry Smith 31004f1ad80SBarry Smith ierr = VecGetLocalSize(uglyoo,&n);CHKERRQ(ierr); 3111ebc52fbSHong Zhang ierr = VecGetArray(uglyoo,&o);CHKERRQ(ierr); 31204f1ad80SBarry Smith for (i=0; i<n; i++) { 31304f1ad80SBarry Smith o[i] = s[uglyrmapo[i]]; /* copy "off-diagonal" portion of scale into oo vector */ 31404f1ad80SBarry Smith } 315bca11509SBarry Smith ierr = VecRestoreArrayRead(scale,&s);CHKERRQ(ierr); 3161ebc52fbSHong Zhang ierr = VecRestoreArray(uglyoo,&o);CHKERRQ(ierr); 31704f1ad80SBarry Smith /* column scale "off-diagonal" portion of local matrix */ 3180298fd71SBarry Smith ierr = MatDiagonalScale(a->B,NULL,uglyoo);CHKERRQ(ierr); 31904f1ad80SBarry Smith PetscFunctionReturn(0); 32004f1ad80SBarry Smith } 321