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> 6bba1ac68SSatish Balay 709573ac7SBarry Smith extern PetscErrorCode MatSetValuesBlocked_SeqBAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode); 88016bdd1SSatish Balay 94a2ae208SSatish Balay #undef __FUNCT__ 104a2ae208SSatish Balay #define __FUNCT__ "MatSetUpMultiply_MPIBAIJ" 11dfbe8321SBarry Smith PetscErrorCode MatSetUpMultiply_MPIBAIJ(Mat mat) 128016bdd1SSatish Balay { 13d9653453SSatish Balay Mat_MPIBAIJ *baij = (Mat_MPIBAIJ*)mat->data; 14d9653453SSatish Balay Mat_SeqBAIJ *B = (Mat_SeqBAIJ*)(baij->B->data); 156849ba73SBarry Smith PetscErrorCode ierr; 16b24ad042SBarry Smith PetscInt i,j,*aj = B->j,ec = 0,*garray; 17d0f46423SBarry Smith PetscInt bs = mat->rmap->bs,*stmp; 188016bdd1SSatish Balay IS from,to; 198016bdd1SSatish Balay Vec gvec; 20aa482453SBarry Smith #if defined (PETSC_USE_CTABLE) 210f5bd95cSBarry Smith PetscTable gid1_lid1; 220f5bd95cSBarry Smith PetscTablePosition tpos; 23b24ad042SBarry Smith PetscInt gid,lid; 246f531f54SSatish Balay #else 25b24ad042SBarry Smith PetscInt Nbs = baij->Nbs,*indices; 2673a2e727SSatish Balay #endif 278016bdd1SSatish Balay 283a40ed3dSBarry Smith PetscFunctionBegin; 2973a2e727SSatish Balay 30aa482453SBarry Smith #if defined (PETSC_USE_CTABLE) 3173a2e727SSatish Balay /* use a table - Mark Adams */ 32273d9f13SBarry Smith ierr = PetscTableCreate(B->mbs,&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 */ 390f5bd95cSBarry Smith ierr = PetscTableAdd(gid1_lid1,gid1,++ec);CHKERRQ(ierr); 4073a2e727SSatish Balay } 4173a2e727SSatish Balay } 4273a2e727SSatish Balay } 4373a2e727SSatish Balay /* form array of columns we need */ 44b24ad042SBarry Smith ierr = PetscMalloc((ec+1)*sizeof(PetscInt),&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++) { 540f5bd95cSBarry Smith ierr = PetscTableAdd(gid1_lid1,garray[i]+1,i+1);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; 6726283091SBarry Smith ierr = PetscLayoutSetUp((baij->B->cmap));CHKERRQ(ierr); 686bc0bbbfSBarry Smith ierr = PetscTableDestroy(&gid1_lid1);CHKERRQ(ierr); 6973a2e727SSatish Balay #else 70435da068SBarry Smith /* Make an array as long as the number of columns */ 71d9653453SSatish Balay /* mark those columns that are in baij->B */ 72b24ad042SBarry Smith ierr = PetscMalloc((Nbs+1)*sizeof(PetscInt),&indices);CHKERRQ(ierr); 73b24ad042SBarry Smith ierr = PetscMemzero(indices,Nbs*sizeof(PetscInt));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 */ 82b24ad042SBarry Smith ierr = PetscMalloc((ec+1)*sizeof(PetscInt),&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; 10326283091SBarry Smith ierr = PetscLayoutSetUp((baij->B->cmap));CHKERRQ(ierr); 104606d414cSSatish Balay ierr = PetscFree(indices);CHKERRQ(ierr); 10573a2e727SSatish Balay #endif 1068016bdd1SSatish Balay 1078016bdd1SSatish Balay /* create local vector that is used to scatter into */ 108029af93fSBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,ec*bs,&baij->lvec);CHKERRQ(ierr); 1098016bdd1SSatish Balay 110c16cb8f2SBarry Smith /* create two temporary index sets for building scatter-gather */ 111deff0451SBarry Smith ierr = ISCreateBlock(PETSC_COMM_SELF,bs,ec,garray,PETSC_COPY_VALUES,&from);CHKERRQ(ierr); 112c16cb8f2SBarry Smith 113b24ad042SBarry Smith ierr = PetscMalloc((ec+1)*sizeof(PetscInt),&stmp);CHKERRQ(ierr); 114e82e9f6bSBarry Smith for (i=0; i<ec; i++) { stmp[i] = i; } 115deff0451SBarry Smith ierr = ISCreateBlock(PETSC_COMM_SELF,bs,ec,stmp,PETSC_OWN_POINTER,&to);CHKERRQ(ierr); 1168016bdd1SSatish Balay 1178016bdd1SSatish Balay /* create temporary global vector to generate scatter context */ 118d0f46423SBarry Smith ierr = VecCreateMPIWithArray(((PetscObject)mat)->comm,mat->cmap->n,mat->cmap->N,PETSC_NULL,&gvec);CHKERRQ(ierr); 1198016bdd1SSatish Balay 120d9653453SSatish Balay ierr = VecScatterCreate(gvec,from,baij->lvec,to,&baij->Mvctx);CHKERRQ(ierr); 12190f02eecSBarry Smith 12252e6d16bSBarry Smith ierr = PetscLogObjectParent(mat,baij->Mvctx);CHKERRQ(ierr); 12352e6d16bSBarry Smith ierr = PetscLogObjectParent(mat,baij->lvec);CHKERRQ(ierr); 12452e6d16bSBarry Smith ierr = PetscLogObjectParent(mat,from);CHKERRQ(ierr); 12552e6d16bSBarry Smith ierr = PetscLogObjectParent(mat,to);CHKERRQ(ierr); 126d9653453SSatish Balay baij->garray = garray; 12752e6d16bSBarry Smith ierr = PetscLogObjectMemory(mat,(ec+1)*sizeof(PetscInt));CHKERRQ(ierr); 1286bf464f9SBarry Smith ierr = ISDestroy(&from);CHKERRQ(ierr); 1296bf464f9SBarry Smith ierr = ISDestroy(&to);CHKERRQ(ierr); 1306bf464f9SBarry Smith ierr = VecDestroy(&gvec);CHKERRQ(ierr); 1313a40ed3dSBarry Smith PetscFunctionReturn(0); 1328016bdd1SSatish Balay } 1338016bdd1SSatish Balay 1348016bdd1SSatish Balay /* 135d9653453SSatish Balay Takes the local part of an already assembled MPIBAIJ matrix 1368016bdd1SSatish Balay and disassembles it. This is to allow new nonzeros into the matrix 1378016bdd1SSatish Balay that require more communication in the matrix vector multiply. 1388016bdd1SSatish Balay Thus certain data-structures must be rebuilt. 1398016bdd1SSatish Balay 1408016bdd1SSatish Balay Kind of slow! But that's what application programmers get when 1418016bdd1SSatish Balay they are sloppy. 1428016bdd1SSatish Balay */ 1434a2ae208SSatish Balay #undef __FUNCT__ 1444a2ae208SSatish Balay #define __FUNCT__ "DisAssemble_MPIBAIJ" 145dfbe8321SBarry Smith PetscErrorCode DisAssemble_MPIBAIJ(Mat A) 1468016bdd1SSatish Balay { 147d9653453SSatish Balay Mat_MPIBAIJ *baij = (Mat_MPIBAIJ*)A->data; 148d9653453SSatish Balay Mat B = baij->B,Bnew; 149d9653453SSatish Balay Mat_SeqBAIJ *Bbaij = (Mat_SeqBAIJ*)B->data; 150dfbe8321SBarry Smith PetscErrorCode ierr; 151d0f46423SBarry Smith PetscInt i,j,mbs=Bbaij->mbs,n = A->cmap->N,col,*garray=baij->garray; 152d0f46423SBarry Smith PetscInt bs2 = baij->bs2,*nz,ec,m = A->rmap->n; 1533eda8832SBarry Smith MatScalar *a = Bbaij->a; 154dd6ea824SBarry Smith MatScalar *atmp; 15597e5c40aSBarry Smith 1568016bdd1SSatish Balay 1573a40ed3dSBarry Smith PetscFunctionBegin; 1588016bdd1SSatish Balay /* free stuff related to matrix-vec multiply */ 159b0a32e0cSBarry Smith ierr = VecGetSize(baij->lvec,&ec);CHKERRQ(ierr); /* needed for PetscLogObjectMemory below */ 1606bf464f9SBarry Smith ierr = VecDestroy(&baij->lvec);CHKERRQ(ierr); baij->lvec = 0; 1616bf464f9SBarry Smith ierr = VecScatterDestroy(&baij->Mvctx);CHKERRQ(ierr); baij->Mvctx = 0; 162d9653453SSatish Balay if (baij->colmap) { 163aa482453SBarry Smith #if defined (PETSC_USE_CTABLE) 1646bc0bbbfSBarry Smith ierr = PetscTableDestroy(&baij->colmap);CHKERRQ(ierr); 16548e59246SSatish Balay #else 166606d414cSSatish Balay ierr = PetscFree(baij->colmap);CHKERRQ(ierr); 16752e6d16bSBarry Smith ierr = PetscLogObjectMemory(A,-Bbaij->nbs*sizeof(PetscInt));CHKERRQ(ierr); 16848e59246SSatish Balay #endif 1698016bdd1SSatish Balay } 1708016bdd1SSatish Balay 1718016bdd1SSatish Balay /* make sure that B is assembled so we can access its values */ 1726d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1733eda8832SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1748016bdd1SSatish Balay 1758016bdd1SSatish Balay /* invent new B and copy stuff over */ 176b24ad042SBarry Smith ierr = PetscMalloc(mbs*sizeof(PetscInt),&nz);CHKERRQ(ierr); 177d9653453SSatish Balay for (i=0; i<mbs; i++) { 178d9653453SSatish Balay nz[i] = Bbaij->i[i+1]-Bbaij->i[i]; 1798016bdd1SSatish Balay } 1807adad957SLisandro Dalcin ierr = MatCreate(((PetscObject)B)->comm,&Bnew);CHKERRQ(ierr); 181f69a0ea3SMatthew Knepley ierr = MatSetSizes(Bnew,m,n,m,n);CHKERRQ(ierr); 1827adad957SLisandro Dalcin ierr = MatSetType(Bnew,((PetscObject)B)->type_name);CHKERRQ(ierr); 183d0f46423SBarry Smith ierr = MatSeqBAIJSetPreallocation(Bnew,B->rmap->bs,0,nz);CHKERRQ(ierr); 1844e0d8c25SBarry Smith ierr = MatSetOption(Bnew,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr); 185d9653453SSatish Balay 186bba1ac68SSatish Balay for (i=0; i<mbs; i++) { 187bba1ac68SSatish Balay for (j=Bbaij->i[i]; j<Bbaij->i[i+1]; j++) { 188bba1ac68SSatish Balay col = garray[Bbaij->j[j]]; 1893eda8832SBarry Smith atmp = a + j*bs2; 190bba1ac68SSatish Balay ierr = MatSetValuesBlocked_SeqBAIJ(Bnew,1,&i,1,&col,atmp,B->insertmode);CHKERRQ(ierr); 1918016bdd1SSatish Balay } 1928016bdd1SSatish Balay } 1934e0d8c25SBarry Smith ierr = MatSetOption(Bnew,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr); 194bba1ac68SSatish Balay 195bba1ac68SSatish Balay ierr = PetscFree(nz);CHKERRQ(ierr); 196606d414cSSatish Balay ierr = PetscFree(baij->garray);CHKERRQ(ierr); 19752e6d16bSBarry Smith ierr = PetscLogObjectMemory(A,-ec*sizeof(PetscInt));CHKERRQ(ierr); 1986bf464f9SBarry Smith ierr = MatDestroy(&B);CHKERRQ(ierr); 19952e6d16bSBarry Smith ierr = PetscLogObjectParent(A,Bnew);CHKERRQ(ierr); 200d9653453SSatish Balay baij->B = Bnew; 2018016bdd1SSatish Balay A->was_assembled = PETSC_FALSE; 2026a719282SBarry Smith A->assembled = PETSC_FALSE; 2033a40ed3dSBarry Smith PetscFunctionReturn(0); 2048016bdd1SSatish Balay } 2058016bdd1SSatish Balay 20604f1ad80SBarry Smith /* ugly stuff added for Glenn someday we should fix this up */ 20704f1ad80SBarry Smith 208b24ad042SBarry Smith static PetscInt *uglyrmapd = 0,*uglyrmapo = 0; /* mapping from the local ordering to the "diagonal" and "off-diagonal" 20904f1ad80SBarry Smith parts of the local matrix */ 2102cd6534aSBarry Smith static Vec uglydd = 0,uglyoo = 0; /* work vectors used to scale the two parts of the local matrix */ 21104f1ad80SBarry Smith 21204f1ad80SBarry Smith 213348d1a9dSSatish Balay #undef __FUNCT__ 214348d1a9dSSatish Balay #define __FUNCT__ "MatMPIBAIJDiagonalScaleLocalSetUp" 215dfbe8321SBarry Smith PetscErrorCode MatMPIBAIJDiagonalScaleLocalSetUp(Mat inA,Vec scale) 21604f1ad80SBarry Smith { 21704f1ad80SBarry Smith Mat_MPIBAIJ *ina = (Mat_MPIBAIJ*) inA->data; /*access private part of matrix */ 21804f1ad80SBarry Smith Mat_SeqBAIJ *B = (Mat_SeqBAIJ*)ina->B->data; 219dfbe8321SBarry Smith PetscErrorCode ierr; 220d0f46423SBarry Smith PetscInt bs = inA->rmap->bs,i,n,nt,j,cstart,cend,no,*garray = ina->garray,*lindices; 221b24ad042SBarry Smith PetscInt *r_rmapd,*r_rmapo; 22204f1ad80SBarry Smith 22304f1ad80SBarry Smith PetscFunctionBegin; 22404f1ad80SBarry Smith ierr = MatGetOwnershipRange(inA,&cstart,&cend);CHKERRQ(ierr); 22504f1ad80SBarry Smith ierr = MatGetSize(ina->A,PETSC_NULL,&n);CHKERRQ(ierr); 226*992144d0SBarry Smith ierr = PetscMalloc((inA->rmap->bmapping->n+1)*sizeof(PetscInt),&r_rmapd);CHKERRQ(ierr); 227*992144d0SBarry Smith ierr = PetscMemzero(r_rmapd,inA->rmap->bmapping->n*sizeof(PetscInt));CHKERRQ(ierr); 22804f1ad80SBarry Smith nt = 0; 229*992144d0SBarry Smith for (i=0; i<inA->rmap->bmapping->n; i++) { 230*992144d0SBarry Smith if (inA->rmap->bmapping->indices[i]*bs >= cstart && inA->rmap->bmapping->indices[i]*bs < cend) { 23104f1ad80SBarry Smith nt++; 232*992144d0SBarry Smith r_rmapd[i] = inA->rmap->bmapping->indices[i] + 1; 23304f1ad80SBarry Smith } 23404f1ad80SBarry Smith } 235e32f2f54SBarry Smith if (nt*bs != n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hmm nt*bs %D n %D",nt*bs,n); 236b24ad042SBarry Smith ierr = PetscMalloc((n+1)*sizeof(PetscInt),&uglyrmapd);CHKERRQ(ierr); 237*992144d0SBarry Smith for (i=0; i<inA->rmap->bmapping->n; i++) { 23804f1ad80SBarry Smith if (r_rmapd[i]){ 23904f1ad80SBarry Smith for (j=0; j<bs; j++) { 24004f1ad80SBarry Smith uglyrmapd[(r_rmapd[i]-1)*bs+j-cstart] = i*bs + j; 24104f1ad80SBarry Smith } 24204f1ad80SBarry Smith } 24304f1ad80SBarry Smith } 24404f1ad80SBarry Smith ierr = PetscFree(r_rmapd);CHKERRQ(ierr); 24504f1ad80SBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,n,&uglydd);CHKERRQ(ierr); 24604f1ad80SBarry Smith 247b24ad042SBarry Smith ierr = PetscMalloc((ina->Nbs+1)*sizeof(PetscInt),&lindices);CHKERRQ(ierr); 248b24ad042SBarry Smith ierr = PetscMemzero(lindices,ina->Nbs*sizeof(PetscInt));CHKERRQ(ierr); 24904f1ad80SBarry Smith for (i=0; i<B->nbs; i++) { 25004f1ad80SBarry Smith lindices[garray[i]] = i+1; 25104f1ad80SBarry Smith } 252*992144d0SBarry Smith no = inA->rmap->bmapping->n - nt; 253*992144d0SBarry Smith ierr = PetscMalloc((inA->rmap->bmapping->n+1)*sizeof(PetscInt),&r_rmapo);CHKERRQ(ierr); 254*992144d0SBarry Smith ierr = PetscMemzero(r_rmapo,inA->rmap->bmapping->n*sizeof(PetscInt));CHKERRQ(ierr); 25504f1ad80SBarry Smith nt = 0; 256*992144d0SBarry Smith for (i=0; i<inA->rmap->bmapping->n; i++) { 257*992144d0SBarry Smith if (lindices[inA->rmap->bmapping->indices[i]]) { 25804f1ad80SBarry Smith nt++; 259*992144d0SBarry Smith r_rmapo[i] = lindices[inA->rmap->bmapping->indices[i]]; 26004f1ad80SBarry Smith } 26104f1ad80SBarry Smith } 262e32f2f54SBarry Smith if (nt > no) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hmm nt %D no %D",nt,n); 26304f1ad80SBarry Smith ierr = PetscFree(lindices);CHKERRQ(ierr); 264b24ad042SBarry Smith ierr = PetscMalloc((nt*bs+1)*sizeof(PetscInt),&uglyrmapo);CHKERRQ(ierr); 265*992144d0SBarry Smith for (i=0; i<inA->rmap->bmapping->n; i++) { 26604f1ad80SBarry Smith if (r_rmapo[i]){ 26704f1ad80SBarry Smith for (j=0; j<bs; j++) { 26804f1ad80SBarry Smith uglyrmapo[(r_rmapo[i]-1)*bs+j] = i*bs + j; 26904f1ad80SBarry Smith } 27004f1ad80SBarry Smith } 27104f1ad80SBarry Smith } 27204f1ad80SBarry Smith ierr = PetscFree(r_rmapo);CHKERRQ(ierr); 27304f1ad80SBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,nt*bs,&uglyoo);CHKERRQ(ierr); 27404f1ad80SBarry Smith 27504f1ad80SBarry Smith PetscFunctionReturn(0); 27604f1ad80SBarry Smith } 27704f1ad80SBarry Smith 278348d1a9dSSatish Balay #undef __FUNCT__ 279348d1a9dSSatish Balay #define __FUNCT__ "MatMPIBAIJDiagonalScaleLocal" 2807087cfbeSBarry Smith PetscErrorCode MatMPIBAIJDiagonalScaleLocal(Mat A,Vec scale) 28104f1ad80SBarry Smith { 28292b32695SKris Buschelman /* This routine should really be abandoned as it duplicates MatDiagonalScaleLocal */ 2834ac538c5SBarry Smith PetscErrorCode ierr; 28492b32695SKris Buschelman 28592b32695SKris Buschelman PetscFunctionBegin; 2864ac538c5SBarry Smith ierr = PetscTryMethod(A,"MatDiagonalScaleLocal_C",(Mat,Vec),(A,scale));CHKERRQ(ierr); 28792b32695SKris Buschelman PetscFunctionReturn(0); 28892b32695SKris Buschelman } 28992b32695SKris Buschelman 29092b32695SKris Buschelman EXTERN_C_BEGIN 29192b32695SKris Buschelman #undef __FUNCT__ 29292b32695SKris Buschelman #define __FUNCT__ "MatDiagonalScaleLocal_MPIBAIJ" 2937087cfbeSBarry Smith PetscErrorCode MatDiagonalScaleLocal_MPIBAIJ(Mat A,Vec scale) 29492b32695SKris Buschelman { 29504f1ad80SBarry Smith Mat_MPIBAIJ *a = (Mat_MPIBAIJ*) A->data; /*access private part of matrix */ 296dfbe8321SBarry Smith PetscErrorCode ierr; 297b24ad042SBarry Smith PetscInt n,i; 29804f1ad80SBarry Smith PetscScalar *d,*o,*s; 29904f1ad80SBarry Smith 30004f1ad80SBarry Smith PetscFunctionBegin; 3012cd6534aSBarry Smith if (!uglyrmapd) { 3022cd6534aSBarry Smith ierr = MatMPIBAIJDiagonalScaleLocalSetUp(A,scale);CHKERRQ(ierr); 3032cd6534aSBarry Smith } 3042cd6534aSBarry Smith 3051ebc52fbSHong Zhang ierr = VecGetArray(scale,&s);CHKERRQ(ierr); 30604f1ad80SBarry Smith 30704f1ad80SBarry Smith ierr = VecGetLocalSize(uglydd,&n);CHKERRQ(ierr); 3081ebc52fbSHong Zhang ierr = VecGetArray(uglydd,&d);CHKERRQ(ierr); 30904f1ad80SBarry Smith for (i=0; i<n; i++) { 31004f1ad80SBarry Smith d[i] = s[uglyrmapd[i]]; /* copy "diagonal" (true local) portion of scale into dd vector */ 31104f1ad80SBarry Smith } 3121ebc52fbSHong Zhang ierr = VecRestoreArray(uglydd,&d);CHKERRQ(ierr); 31304f1ad80SBarry Smith /* column scale "diagonal" portion of local matrix */ 31404f1ad80SBarry Smith ierr = MatDiagonalScale(a->A,PETSC_NULL,uglydd);CHKERRQ(ierr); 31504f1ad80SBarry Smith 31604f1ad80SBarry Smith ierr = VecGetLocalSize(uglyoo,&n);CHKERRQ(ierr); 3171ebc52fbSHong Zhang ierr = VecGetArray(uglyoo,&o);CHKERRQ(ierr); 31804f1ad80SBarry Smith for (i=0; i<n; i++) { 31904f1ad80SBarry Smith o[i] = s[uglyrmapo[i]]; /* copy "off-diagonal" portion of scale into oo vector */ 32004f1ad80SBarry Smith } 3211ebc52fbSHong Zhang ierr = VecRestoreArray(scale,&s);CHKERRQ(ierr); 3221ebc52fbSHong Zhang ierr = VecRestoreArray(uglyoo,&o);CHKERRQ(ierr); 32304f1ad80SBarry Smith /* column scale "off-diagonal" portion of local matrix */ 32404f1ad80SBarry Smith ierr = MatDiagonalScale(a->B,PETSC_NULL,uglyoo);CHKERRQ(ierr); 32504f1ad80SBarry Smith 32604f1ad80SBarry Smith PetscFunctionReturn(0); 32704f1ad80SBarry Smith } 32892b32695SKris Buschelman EXTERN_C_END 3298016bdd1SSatish Balay 3303eda8832SBarry Smith 331