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> 6186d4ecdSBarry Smith #include <petsc-private/isimpl.h> /* needed because accesses data structure of ISLocalToGlobalMapping directly */ 7bba1ac68SSatish Balay 809573ac7SBarry Smith extern PetscErrorCode MatSetValuesBlocked_SeqBAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 98016bdd1SSatish Balay 104a2ae208SSatish Balay #undef __FUNCT__ 114a2ae208SSatish Balay #define __FUNCT__ "MatSetUpMultiply_MPIBAIJ" 12dfbe8321SBarry Smith PetscErrorCode MatSetUpMultiply_MPIBAIJ(Mat mat) 138016bdd1SSatish Balay { 14d9653453SSatish Balay Mat_MPIBAIJ *baij = (Mat_MPIBAIJ*)mat->data; 15d9653453SSatish Balay Mat_SeqBAIJ *B = (Mat_SeqBAIJ*)(baij->B->data); 166849ba73SBarry Smith PetscErrorCode ierr; 17b24ad042SBarry Smith PetscInt i,j,*aj = B->j,ec = 0,*garray; 18d0f46423SBarry Smith PetscInt bs = mat->rmap->bs,*stmp; 198016bdd1SSatish Balay IS from,to; 208016bdd1SSatish Balay Vec gvec; 21aa482453SBarry Smith #if defined(PETSC_USE_CTABLE) 220f5bd95cSBarry Smith PetscTable gid1_lid1; 230f5bd95cSBarry Smith PetscTablePosition tpos; 24b24ad042SBarry Smith PetscInt gid,lid; 256f531f54SSatish Balay #else 26b24ad042SBarry Smith PetscInt Nbs = baij->Nbs,*indices; 2773a2e727SSatish Balay #endif 288016bdd1SSatish Balay 293a40ed3dSBarry Smith PetscFunctionBegin; 30aa482453SBarry Smith #if defined(PETSC_USE_CTABLE) 3173a2e727SSatish Balay /* use a table - Mark Adams */ 32e23dfa41SBarry Smith ierr = PetscTableCreate(B->mbs,baij->Nbs+1,&gid1_lid1);CHKERRQ(ierr); 3373a2e727SSatish Balay for (i=0; i<B->mbs; i++) { 3473a2e727SSatish Balay for (j=0; j<B->ilen[i]; j++) { 35b24ad042SBarry Smith PetscInt data,gid1 = aj[B->i[i]+j] + 1; 360f5bd95cSBarry Smith ierr = PetscTableFind(gid1_lid1,gid1,&data);CHKERRQ(ierr); 37fa46199cSSatish Balay if (!data) { 3873a2e727SSatish Balay /* one based table */ 393861aac3SJed Brown ierr = PetscTableAdd(gid1_lid1,gid1,++ec,INSERT_VALUES);CHKERRQ(ierr); 4073a2e727SSatish Balay } 4173a2e727SSatish Balay } 4273a2e727SSatish Balay } 4373a2e727SSatish Balay /* form array of columns we need */ 44785e854fSJed Brown ierr = PetscMalloc1((ec+1),&garray);CHKERRQ(ierr); 450f5bd95cSBarry Smith ierr = PetscTableGetHeadPosition(gid1_lid1,&tpos);CHKERRQ(ierr); 4673a2e727SSatish Balay while (tpos) { 470f5bd95cSBarry Smith ierr = PetscTableGetNext(gid1_lid1,&tpos,&gid,&lid);CHKERRQ(ierr); 4873a2e727SSatish Balay gid--; lid--; 4973a2e727SSatish Balay garray[lid] = gid; 5073a2e727SSatish Balay } 510064e2bbSSatish Balay ierr = PetscSortInt(ec,garray);CHKERRQ(ierr); 520f5bd95cSBarry Smith ierr = PetscTableRemoveAll(gid1_lid1);CHKERRQ(ierr); 5373a2e727SSatish Balay for (i=0; i<ec; i++) { 543861aac3SJed Brown ierr = PetscTableAdd(gid1_lid1,garray[i]+1,i+1,INSERT_VALUES);CHKERRQ(ierr); 5573a2e727SSatish Balay } 5673a2e727SSatish Balay /* compact out the extra columns in B */ 5773a2e727SSatish Balay for (i=0; i<B->mbs; i++) { 5873a2e727SSatish Balay for (j=0; j<B->ilen[i]; j++) { 59b24ad042SBarry Smith PetscInt gid1 = aj[B->i[i] + j] + 1; 600f5bd95cSBarry Smith ierr = PetscTableFind(gid1_lid1,gid1,&lid);CHKERRQ(ierr); 61fa46199cSSatish Balay lid--; 6273a2e727SSatish Balay aj[B->i[i]+j] = lid; 6373a2e727SSatish Balay } 6473a2e727SSatish Balay } 6573a2e727SSatish Balay B->nbs = ec; 66d0f46423SBarry Smith baij->B->cmap->n = baij->B->cmap->N = ec*mat->rmap->bs; 6726fbe8dcSKarl Rupp 6826283091SBarry Smith ierr = PetscLayoutSetUp((baij->B->cmap));CHKERRQ(ierr); 696bc0bbbfSBarry Smith ierr = PetscTableDestroy(&gid1_lid1);CHKERRQ(ierr); 7073a2e727SSatish Balay #else 71435da068SBarry Smith /* Make an array as long as the number of columns */ 72d9653453SSatish Balay /* mark those columns that are in baij->B */ 731795a4d1SJed Brown ierr = PetscCalloc1(Nbs+1,&indices);CHKERRQ(ierr); 74d9653453SSatish Balay for (i=0; i<B->mbs; i++) { 758016bdd1SSatish Balay for (j=0; j<B->ilen[i]; j++) { 76d9653453SSatish Balay if (!indices[aj[B->i[i] + j]]) ec++; 77d9653453SSatish Balay indices[aj[B->i[i] + j]] = 1; 788016bdd1SSatish Balay } 798016bdd1SSatish Balay } 808016bdd1SSatish Balay 818016bdd1SSatish Balay /* form array of columns we need */ 82785e854fSJed Brown ierr = PetscMalloc1((ec+1),&garray);CHKERRQ(ierr); 838016bdd1SSatish Balay ec = 0; 84d9653453SSatish Balay for (i=0; i<Nbs; i++) { 850bdbc534SSatish Balay if (indices[i]) { 860bdbc534SSatish Balay garray[ec++] = i; 870bdbc534SSatish Balay } 888016bdd1SSatish Balay } 898016bdd1SSatish Balay 908016bdd1SSatish Balay /* make indices now point into garray */ 918016bdd1SSatish Balay for (i=0; i<ec; i++) { 92d9653453SSatish Balay indices[garray[i]] = i; 938016bdd1SSatish Balay } 948016bdd1SSatish Balay 958016bdd1SSatish Balay /* compact out the extra columns in B */ 96d9653453SSatish Balay for (i=0; i<B->mbs; i++) { 978016bdd1SSatish Balay for (j=0; j<B->ilen[i]; j++) { 98d9653453SSatish Balay aj[B->i[i] + j] = indices[aj[B->i[i] + j]]; 998016bdd1SSatish Balay } 1008016bdd1SSatish Balay } 101d9653453SSatish Balay B->nbs = ec; 102d0f46423SBarry Smith baij->B->cmap->n = baij->B->cmap->N = ec*mat->rmap->bs; 10326fbe8dcSKarl Rupp 10426283091SBarry Smith ierr = PetscLayoutSetUp((baij->B->cmap));CHKERRQ(ierr); 105606d414cSSatish Balay ierr = PetscFree(indices);CHKERRQ(ierr); 10673a2e727SSatish Balay #endif 1078016bdd1SSatish Balay 1088016bdd1SSatish Balay /* create local vector that is used to scatter into */ 109029af93fSBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,ec*bs,&baij->lvec);CHKERRQ(ierr); 1108016bdd1SSatish Balay 111c16cb8f2SBarry Smith /* create two temporary index sets for building scatter-gather */ 112deff0451SBarry Smith ierr = ISCreateBlock(PETSC_COMM_SELF,bs,ec,garray,PETSC_COPY_VALUES,&from);CHKERRQ(ierr); 113c16cb8f2SBarry Smith 114785e854fSJed Brown ierr = PetscMalloc1((ec+1),&stmp);CHKERRQ(ierr); 11526fbe8dcSKarl Rupp for (i=0; i<ec; i++) stmp[i] = i; 116deff0451SBarry Smith ierr = ISCreateBlock(PETSC_COMM_SELF,bs,ec,stmp,PETSC_OWN_POINTER,&to);CHKERRQ(ierr); 1178016bdd1SSatish Balay 1188016bdd1SSatish Balay /* create temporary global vector to generate scatter context */ 119ce94432eSBarry Smith ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)mat),1,mat->cmap->n,mat->cmap->N,NULL,&gvec);CHKERRQ(ierr); 1208016bdd1SSatish Balay 121d9653453SSatish Balay ierr = VecScatterCreate(gvec,from,baij->lvec,to,&baij->Mvctx);CHKERRQ(ierr); 12290f02eecSBarry Smith 1233bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)baij->Mvctx);CHKERRQ(ierr); 1243bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)baij->lvec);CHKERRQ(ierr); 1253bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)from);CHKERRQ(ierr); 1263bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)to);CHKERRQ(ierr); 12726fbe8dcSKarl Rupp 128d9653453SSatish Balay baij->garray = garray; 12926fbe8dcSKarl Rupp 1303bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)mat,(ec+1)*sizeof(PetscInt));CHKERRQ(ierr); 1316bf464f9SBarry Smith ierr = ISDestroy(&from);CHKERRQ(ierr); 1326bf464f9SBarry Smith ierr = ISDestroy(&to);CHKERRQ(ierr); 1336bf464f9SBarry Smith ierr = VecDestroy(&gvec);CHKERRQ(ierr); 1343a40ed3dSBarry Smith PetscFunctionReturn(0); 1358016bdd1SSatish Balay } 1368016bdd1SSatish Balay 1378016bdd1SSatish Balay /* 138d9653453SSatish Balay Takes the local part of an already assembled MPIBAIJ matrix 1398016bdd1SSatish Balay and disassembles it. This is to allow new nonzeros into the matrix 1408016bdd1SSatish Balay that require more communication in the matrix vector multiply. 1418016bdd1SSatish Balay Thus certain data-structures must be rebuilt. 1428016bdd1SSatish Balay 1438016bdd1SSatish Balay Kind of slow! But that's what application programmers get when 1448016bdd1SSatish Balay they are sloppy. 1458016bdd1SSatish Balay */ 1464a2ae208SSatish Balay #undef __FUNCT__ 147ab9863d7SBarry Smith #define __FUNCT__ "MatDisAssemble_MPIBAIJ" 148ab9863d7SBarry Smith PetscErrorCode MatDisAssemble_MPIBAIJ(Mat A) 1498016bdd1SSatish Balay { 150d9653453SSatish Balay Mat_MPIBAIJ *baij = (Mat_MPIBAIJ*)A->data; 151d9653453SSatish Balay Mat B = baij->B,Bnew; 152d9653453SSatish Balay Mat_SeqBAIJ *Bbaij = (Mat_SeqBAIJ*)B->data; 153dfbe8321SBarry Smith PetscErrorCode ierr; 154d0f46423SBarry Smith PetscInt i,j,mbs=Bbaij->mbs,n = A->cmap->N,col,*garray=baij->garray; 155d0f46423SBarry Smith PetscInt bs2 = baij->bs2,*nz,ec,m = A->rmap->n; 1563eda8832SBarry Smith MatScalar *a = Bbaij->a; 157dd6ea824SBarry Smith MatScalar *atmp; 15897e5c40aSBarry Smith 1598016bdd1SSatish Balay 1603a40ed3dSBarry Smith PetscFunctionBegin; 1618016bdd1SSatish Balay /* free stuff related to matrix-vec multiply */ 162b0a32e0cSBarry Smith ierr = VecGetSize(baij->lvec,&ec);CHKERRQ(ierr); /* needed for PetscLogObjectMemory below */ 1636bf464f9SBarry Smith ierr = VecDestroy(&baij->lvec);CHKERRQ(ierr); baij->lvec = 0; 1646bf464f9SBarry Smith ierr = VecScatterDestroy(&baij->Mvctx);CHKERRQ(ierr); baij->Mvctx = 0; 165d9653453SSatish Balay if (baij->colmap) { 166aa482453SBarry Smith #if defined(PETSC_USE_CTABLE) 1676bc0bbbfSBarry Smith ierr = PetscTableDestroy(&baij->colmap);CHKERRQ(ierr); 16848e59246SSatish Balay #else 169606d414cSSatish Balay ierr = PetscFree(baij->colmap);CHKERRQ(ierr); 1703bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)A,-Bbaij->nbs*sizeof(PetscInt));CHKERRQ(ierr); 17148e59246SSatish Balay #endif 1728016bdd1SSatish Balay } 1738016bdd1SSatish Balay 1748016bdd1SSatish Balay /* make sure that B is assembled so we can access its values */ 1756d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1763eda8832SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1778016bdd1SSatish Balay 1788016bdd1SSatish Balay /* invent new B and copy stuff over */ 179785e854fSJed Brown ierr = PetscMalloc1(mbs,&nz);CHKERRQ(ierr); 180d9653453SSatish Balay for (i=0; i<mbs; i++) { 181d9653453SSatish Balay nz[i] = Bbaij->i[i+1]-Bbaij->i[i]; 1828016bdd1SSatish Balay } 183ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)B),&Bnew);CHKERRQ(ierr); 184f69a0ea3SMatthew Knepley ierr = MatSetSizes(Bnew,m,n,m,n);CHKERRQ(ierr); 1857adad957SLisandro Dalcin ierr = MatSetType(Bnew,((PetscObject)B)->type_name);CHKERRQ(ierr); 186d0f46423SBarry Smith ierr = MatSeqBAIJSetPreallocation(Bnew,B->rmap->bs,0,nz);CHKERRQ(ierr); 18726fbe8dcSKarl Rupp 1882576faa2SJed Brown ((Mat_SeqBAIJ*)Bnew->data)->nonew = Bbaij->nonew; /* Inherit insertion error options. */ 18926fbe8dcSKarl Rupp 1904e0d8c25SBarry Smith ierr = MatSetOption(Bnew,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr); 191*77341eacSDmitry Karpeev /* 192*77341eacSDmitry Karpeev Ensure that B's nonzerostate is monotonically increasing. 193*77341eacSDmitry Karpeev Or should this follow the MatSetValuesBlocked() loop to preserve B's nonzerstate across a MatDisAssemble() call? 194*77341eacSDmitry Karpeev */ 195*77341eacSDmitry Karpeev Bnew->nonzerostate = B->nonzerostate; 196d9653453SSatish Balay 197bba1ac68SSatish Balay for (i=0; i<mbs; i++) { 198bba1ac68SSatish Balay for (j=Bbaij->i[i]; j<Bbaij->i[i+1]; j++) { 199bba1ac68SSatish Balay col = garray[Bbaij->j[j]]; 2003eda8832SBarry Smith atmp = a + j*bs2; 201bba1ac68SSatish Balay ierr = MatSetValuesBlocked_SeqBAIJ(Bnew,1,&i,1,&col,atmp,B->insertmode);CHKERRQ(ierr); 2028016bdd1SSatish Balay } 2038016bdd1SSatish Balay } 2044e0d8c25SBarry Smith ierr = MatSetOption(Bnew,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr); 205bba1ac68SSatish Balay 206bba1ac68SSatish Balay ierr = PetscFree(nz);CHKERRQ(ierr); 207606d414cSSatish Balay ierr = PetscFree(baij->garray);CHKERRQ(ierr); 2083bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)A,-ec*sizeof(PetscInt));CHKERRQ(ierr); 2096bf464f9SBarry Smith ierr = MatDestroy(&B);CHKERRQ(ierr); 2103bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)A,(PetscObject)Bnew);CHKERRQ(ierr); 21126fbe8dcSKarl Rupp 212d9653453SSatish Balay baij->B = Bnew; 2138016bdd1SSatish Balay A->was_assembled = PETSC_FALSE; 2146a719282SBarry Smith A->assembled = PETSC_FALSE; 2153a40ed3dSBarry Smith PetscFunctionReturn(0); 2168016bdd1SSatish Balay } 2178016bdd1SSatish Balay 21804f1ad80SBarry Smith /* ugly stuff added for Glenn someday we should fix this up */ 21904f1ad80SBarry Smith 22026fbe8dcSKarl Rupp static PetscInt *uglyrmapd = 0,*uglyrmapo = 0; /* mapping from the local ordering to the "diagonal" and "off-diagonal" parts of the local matrix */ 2212cd6534aSBarry Smith static Vec uglydd = 0,uglyoo = 0; /* work vectors used to scale the two parts of the local matrix */ 22204f1ad80SBarry Smith 22304f1ad80SBarry Smith 224348d1a9dSSatish Balay #undef __FUNCT__ 225348d1a9dSSatish Balay #define __FUNCT__ "MatMPIBAIJDiagonalScaleLocalSetUp" 226dfbe8321SBarry Smith PetscErrorCode MatMPIBAIJDiagonalScaleLocalSetUp(Mat inA,Vec scale) 22704f1ad80SBarry Smith { 22804f1ad80SBarry Smith Mat_MPIBAIJ *ina = (Mat_MPIBAIJ*) inA->data; /*access private part of matrix */ 22904f1ad80SBarry Smith Mat_SeqBAIJ *B = (Mat_SeqBAIJ*)ina->B->data; 230dfbe8321SBarry Smith PetscErrorCode ierr; 231d0f46423SBarry Smith PetscInt bs = inA->rmap->bs,i,n,nt,j,cstart,cend,no,*garray = ina->garray,*lindices; 232b24ad042SBarry Smith PetscInt *r_rmapd,*r_rmapo; 23304f1ad80SBarry Smith 23404f1ad80SBarry Smith PetscFunctionBegin; 23504f1ad80SBarry Smith ierr = MatGetOwnershipRange(inA,&cstart,&cend);CHKERRQ(ierr); 2360298fd71SBarry Smith ierr = MatGetSize(ina->A,NULL,&n);CHKERRQ(ierr); 23745b6f7e9SBarry Smith ierr = PetscCalloc1((inA->rmap->mapping->n+1),&r_rmapd);CHKERRQ(ierr); 23804f1ad80SBarry Smith nt = 0; 23945b6f7e9SBarry Smith for (i=0; i<inA->rmap->mapping->n; i++) { 24045b6f7e9SBarry Smith if (inA->rmap->mapping->indices[i]*bs >= cstart && inA->rmap->mapping->indices[i]*bs < cend) { 24104f1ad80SBarry Smith nt++; 24245b6f7e9SBarry Smith r_rmapd[i] = inA->rmap->mapping->indices[i] + 1; 24304f1ad80SBarry Smith } 24404f1ad80SBarry Smith } 245e32f2f54SBarry Smith if (nt*bs != n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hmm nt*bs %D n %D",nt*bs,n); 246785e854fSJed Brown ierr = PetscMalloc1((n+1),&uglyrmapd);CHKERRQ(ierr); 24745b6f7e9SBarry Smith for (i=0; i<inA->rmap->mapping->n; i++) { 24804f1ad80SBarry Smith if (r_rmapd[i]) { 24904f1ad80SBarry Smith for (j=0; j<bs; j++) { 25004f1ad80SBarry Smith uglyrmapd[(r_rmapd[i]-1)*bs+j-cstart] = i*bs + j; 25104f1ad80SBarry Smith } 25204f1ad80SBarry Smith } 25304f1ad80SBarry Smith } 25404f1ad80SBarry Smith ierr = PetscFree(r_rmapd);CHKERRQ(ierr); 25504f1ad80SBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,n,&uglydd);CHKERRQ(ierr); 25604f1ad80SBarry Smith 2571795a4d1SJed Brown ierr = PetscCalloc1((ina->Nbs+1),&lindices);CHKERRQ(ierr); 25804f1ad80SBarry Smith for (i=0; i<B->nbs; i++) { 25904f1ad80SBarry Smith lindices[garray[i]] = i+1; 26004f1ad80SBarry Smith } 26145b6f7e9SBarry Smith no = inA->rmap->mapping->n - nt; 26245b6f7e9SBarry Smith ierr = PetscCalloc1((inA->rmap->mapping->n+1),&r_rmapo);CHKERRQ(ierr); 26304f1ad80SBarry Smith nt = 0; 26445b6f7e9SBarry Smith for (i=0; i<inA->rmap->mapping->n; i++) { 26545b6f7e9SBarry Smith if (lindices[inA->rmap->mapping->indices[i]]) { 26604f1ad80SBarry Smith nt++; 26745b6f7e9SBarry Smith r_rmapo[i] = lindices[inA->rmap->mapping->indices[i]]; 26804f1ad80SBarry Smith } 26904f1ad80SBarry Smith } 270e32f2f54SBarry Smith if (nt > no) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hmm nt %D no %D",nt,n); 27104f1ad80SBarry Smith ierr = PetscFree(lindices);CHKERRQ(ierr); 272785e854fSJed Brown ierr = PetscMalloc1((nt*bs+1),&uglyrmapo);CHKERRQ(ierr); 27345b6f7e9SBarry Smith for (i=0; i<inA->rmap->mapping->n; i++) { 27404f1ad80SBarry Smith if (r_rmapo[i]) { 27504f1ad80SBarry Smith for (j=0; j<bs; j++) { 27604f1ad80SBarry Smith uglyrmapo[(r_rmapo[i]-1)*bs+j] = i*bs + j; 27704f1ad80SBarry Smith } 27804f1ad80SBarry Smith } 27904f1ad80SBarry Smith } 28004f1ad80SBarry Smith ierr = PetscFree(r_rmapo);CHKERRQ(ierr); 28104f1ad80SBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,nt*bs,&uglyoo);CHKERRQ(ierr); 28204f1ad80SBarry Smith PetscFunctionReturn(0); 28304f1ad80SBarry Smith } 28404f1ad80SBarry Smith 285348d1a9dSSatish Balay #undef __FUNCT__ 286348d1a9dSSatish Balay #define __FUNCT__ "MatMPIBAIJDiagonalScaleLocal" 2877087cfbeSBarry Smith PetscErrorCode MatMPIBAIJDiagonalScaleLocal(Mat A,Vec scale) 28804f1ad80SBarry Smith { 28992b32695SKris Buschelman /* This routine should really be abandoned as it duplicates MatDiagonalScaleLocal */ 2904ac538c5SBarry Smith PetscErrorCode ierr; 29192b32695SKris Buschelman 29292b32695SKris Buschelman PetscFunctionBegin; 2934ac538c5SBarry Smith ierr = PetscTryMethod(A,"MatDiagonalScaleLocal_C",(Mat,Vec),(A,scale));CHKERRQ(ierr); 29492b32695SKris Buschelman PetscFunctionReturn(0); 29592b32695SKris Buschelman } 29692b32695SKris Buschelman 29792b32695SKris Buschelman #undef __FUNCT__ 29892b32695SKris Buschelman #define __FUNCT__ "MatDiagonalScaleLocal_MPIBAIJ" 2997087cfbeSBarry Smith PetscErrorCode MatDiagonalScaleLocal_MPIBAIJ(Mat A,Vec scale) 30092b32695SKris Buschelman { 30104f1ad80SBarry Smith Mat_MPIBAIJ *a = (Mat_MPIBAIJ*) A->data; /*access private part of matrix */ 302dfbe8321SBarry Smith PetscErrorCode ierr; 303b24ad042SBarry Smith PetscInt n,i; 30404f1ad80SBarry Smith PetscScalar *d,*o,*s; 30504f1ad80SBarry Smith 30604f1ad80SBarry Smith PetscFunctionBegin; 3072cd6534aSBarry Smith if (!uglyrmapd) { 3082cd6534aSBarry Smith ierr = MatMPIBAIJDiagonalScaleLocalSetUp(A,scale);CHKERRQ(ierr); 3092cd6534aSBarry Smith } 3102cd6534aSBarry Smith 3111ebc52fbSHong Zhang ierr = VecGetArray(scale,&s);CHKERRQ(ierr); 31204f1ad80SBarry Smith 31304f1ad80SBarry Smith ierr = VecGetLocalSize(uglydd,&n);CHKERRQ(ierr); 3141ebc52fbSHong Zhang ierr = VecGetArray(uglydd,&d);CHKERRQ(ierr); 31504f1ad80SBarry Smith for (i=0; i<n; i++) { 31604f1ad80SBarry Smith d[i] = s[uglyrmapd[i]]; /* copy "diagonal" (true local) portion of scale into dd vector */ 31704f1ad80SBarry Smith } 3181ebc52fbSHong Zhang ierr = VecRestoreArray(uglydd,&d);CHKERRQ(ierr); 31904f1ad80SBarry Smith /* column scale "diagonal" portion of local matrix */ 3200298fd71SBarry Smith ierr = MatDiagonalScale(a->A,NULL,uglydd);CHKERRQ(ierr); 32104f1ad80SBarry Smith 32204f1ad80SBarry Smith ierr = VecGetLocalSize(uglyoo,&n);CHKERRQ(ierr); 3231ebc52fbSHong Zhang ierr = VecGetArray(uglyoo,&o);CHKERRQ(ierr); 32404f1ad80SBarry Smith for (i=0; i<n; i++) { 32504f1ad80SBarry Smith o[i] = s[uglyrmapo[i]]; /* copy "off-diagonal" portion of scale into oo vector */ 32604f1ad80SBarry Smith } 3271ebc52fbSHong Zhang ierr = VecRestoreArray(scale,&s);CHKERRQ(ierr); 3281ebc52fbSHong Zhang ierr = VecRestoreArray(uglyoo,&o);CHKERRQ(ierr); 32904f1ad80SBarry Smith /* column scale "off-diagonal" portion of local matrix */ 3300298fd71SBarry Smith ierr = MatDiagonalScale(a->B,NULL,uglyoo);CHKERRQ(ierr); 33104f1ad80SBarry Smith PetscFunctionReturn(0); 33204f1ad80SBarry Smith } 333