1be1d678aSKris Buschelman #define PETSCMAT_DLL 2be1d678aSKris Buschelman 38016bdd1SSatish Balay /* 4d9653453SSatish Balay Support for the parallel BAIJ matrix vector multiply 58016bdd1SSatish Balay */ 6*7c4f633dSBarry Smith #include "../src/mat/impls/baij/mpi/mpibaij.h" 7bba1ac68SSatish Balay 8b24ad042SBarry 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; 3073a2e727SSatish Balay 31aa482453SBarry Smith #if defined (PETSC_USE_CTABLE) 3273a2e727SSatish Balay /* use a table - Mark Adams */ 33273d9f13SBarry Smith ierr = PetscTableCreate(B->mbs,&gid1_lid1);CHKERRQ(ierr); 3473a2e727SSatish Balay for (i=0; i<B->mbs; i++) { 3573a2e727SSatish Balay for (j=0; j<B->ilen[i]; j++) { 36b24ad042SBarry Smith PetscInt data,gid1 = aj[B->i[i]+j] + 1; 370f5bd95cSBarry Smith ierr = PetscTableFind(gid1_lid1,gid1,&data);CHKERRQ(ierr); 38fa46199cSSatish Balay if (!data) { 3973a2e727SSatish Balay /* one based table */ 400f5bd95cSBarry Smith ierr = PetscTableAdd(gid1_lid1,gid1,++ec);CHKERRQ(ierr); 4173a2e727SSatish Balay } 4273a2e727SSatish Balay } 4373a2e727SSatish Balay } 4473a2e727SSatish Balay /* form array of columns we need */ 45b24ad042SBarry Smith ierr = PetscMalloc((ec+1)*sizeof(PetscInt),&garray);CHKERRQ(ierr); 460f5bd95cSBarry Smith ierr = PetscTableGetHeadPosition(gid1_lid1,&tpos);CHKERRQ(ierr); 4773a2e727SSatish Balay while (tpos) { 480f5bd95cSBarry Smith ierr = PetscTableGetNext(gid1_lid1,&tpos,&gid,&lid);CHKERRQ(ierr); 4973a2e727SSatish Balay gid--; lid--; 5073a2e727SSatish Balay garray[lid] = gid; 5173a2e727SSatish Balay } 520064e2bbSSatish Balay ierr = PetscSortInt(ec,garray);CHKERRQ(ierr); 530f5bd95cSBarry Smith ierr = PetscTableRemoveAll(gid1_lid1);CHKERRQ(ierr); 5473a2e727SSatish Balay for (i=0; i<ec; i++) { 550f5bd95cSBarry Smith ierr = PetscTableAdd(gid1_lid1,garray[i]+1,i+1);CHKERRQ(ierr); 5673a2e727SSatish Balay } 5773a2e727SSatish Balay /* compact out the extra columns in B */ 5873a2e727SSatish Balay for (i=0; i<B->mbs; i++) { 5973a2e727SSatish Balay for (j=0; j<B->ilen[i]; j++) { 60b24ad042SBarry Smith PetscInt gid1 = aj[B->i[i] + j] + 1; 610f5bd95cSBarry Smith ierr = PetscTableFind(gid1_lid1,gid1,&lid);CHKERRQ(ierr); 62fa46199cSSatish Balay lid --; 6373a2e727SSatish Balay aj[B->i[i]+j] = lid; 6473a2e727SSatish Balay } 6573a2e727SSatish Balay } 6673a2e727SSatish Balay B->nbs = ec; 67d0f46423SBarry Smith baij->B->cmap->n = baij->B->cmap->N = ec*mat->rmap->bs; 68d0f46423SBarry Smith ierr = PetscMapSetUp((baij->B->cmap));CHKERRQ(ierr); 699c666560SBarry Smith ierr = PetscTableDestroy(gid1_lid1);CHKERRQ(ierr); 7073a2e727SSatish Balay /* Mark Adams */ 7173a2e727SSatish Balay #else 72435da068SBarry Smith /* Make an array as long as the number of columns */ 73d9653453SSatish Balay /* mark those columns that are in baij->B */ 74b24ad042SBarry Smith ierr = PetscMalloc((Nbs+1)*sizeof(PetscInt),&indices);CHKERRQ(ierr); 75b24ad042SBarry Smith ierr = PetscMemzero(indices,Nbs*sizeof(PetscInt));CHKERRQ(ierr); 76d9653453SSatish Balay for (i=0; i<B->mbs; i++) { 778016bdd1SSatish Balay for (j=0; j<B->ilen[i]; j++) { 78d9653453SSatish Balay if (!indices[aj[B->i[i] + j]]) ec++; 79d9653453SSatish Balay indices[aj[B->i[i] + j]] = 1; 808016bdd1SSatish Balay } 818016bdd1SSatish Balay } 828016bdd1SSatish Balay 838016bdd1SSatish Balay /* form array of columns we need */ 84b24ad042SBarry Smith ierr = PetscMalloc((ec+1)*sizeof(PetscInt),&garray);CHKERRQ(ierr); 858016bdd1SSatish Balay ec = 0; 86d9653453SSatish Balay for (i=0; i<Nbs; i++) { 870bdbc534SSatish Balay if (indices[i]) { 880bdbc534SSatish Balay garray[ec++] = i; 890bdbc534SSatish Balay } 908016bdd1SSatish Balay } 918016bdd1SSatish Balay 928016bdd1SSatish Balay /* make indices now point into garray */ 938016bdd1SSatish Balay for (i=0; i<ec; i++) { 94d9653453SSatish Balay indices[garray[i]] = i; 958016bdd1SSatish Balay } 968016bdd1SSatish Balay 978016bdd1SSatish Balay /* compact out the extra columns in B */ 98d9653453SSatish Balay for (i=0; i<B->mbs; i++) { 998016bdd1SSatish Balay for (j=0; j<B->ilen[i]; j++) { 100d9653453SSatish Balay aj[B->i[i] + j] = indices[aj[B->i[i] + j]]; 1018016bdd1SSatish Balay } 1028016bdd1SSatish Balay } 103d9653453SSatish Balay B->nbs = ec; 104d0f46423SBarry Smith baij->B->cmap->n =baij->B->cmap->N = ec*mat->rmap->bs; 105d0f46423SBarry Smith ierr = PetscMapSetUp((baij->B->cmap));CHKERRQ(ierr); 106606d414cSSatish Balay ierr = PetscFree(indices);CHKERRQ(ierr); 10773a2e727SSatish Balay #endif 1088016bdd1SSatish Balay 1098016bdd1SSatish Balay /* create local vector that is used to scatter into */ 110029af93fSBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,ec*bs,&baij->lvec);CHKERRQ(ierr); 1118016bdd1SSatish Balay 112c16cb8f2SBarry Smith /* create two temporary index sets for building scatter-gather */ 1133d3cf644SBarry Smith for (i=0; i<ec; i++) { 114c16cb8f2SBarry Smith garray[i] = bs*garray[i]; 115c16cb8f2SBarry Smith } 116029af93fSBarry Smith ierr = ISCreateBlock(PETSC_COMM_SELF,bs,ec,garray,&from);CHKERRQ(ierr); 1173d3cf644SBarry Smith for (i=0; i<ec; i++) { 118c16cb8f2SBarry Smith garray[i] = garray[i]/bs; 119c16cb8f2SBarry Smith } 120c16cb8f2SBarry Smith 121b24ad042SBarry Smith ierr = PetscMalloc((ec+1)*sizeof(PetscInt),&stmp);CHKERRQ(ierr); 122537820f0SBarry Smith for (i=0; i<ec; i++) { stmp[i] = bs*i; } 123029af93fSBarry Smith ierr = ISCreateBlock(PETSC_COMM_SELF,bs,ec,stmp,&to);CHKERRQ(ierr); 124606d414cSSatish Balay ierr = PetscFree(stmp);CHKERRQ(ierr); 1258016bdd1SSatish Balay 1268016bdd1SSatish Balay /* create temporary global vector to generate scatter context */ 127d0f46423SBarry Smith ierr = VecCreateMPIWithArray(((PetscObject)mat)->comm,mat->cmap->n,mat->cmap->N,PETSC_NULL,&gvec);CHKERRQ(ierr); 1288016bdd1SSatish Balay 129d9653453SSatish Balay ierr = VecScatterCreate(gvec,from,baij->lvec,to,&baij->Mvctx);CHKERRQ(ierr); 13090f02eecSBarry Smith 13152e6d16bSBarry Smith ierr = PetscLogObjectParent(mat,baij->Mvctx);CHKERRQ(ierr); 13252e6d16bSBarry Smith ierr = PetscLogObjectParent(mat,baij->lvec);CHKERRQ(ierr); 13352e6d16bSBarry Smith ierr = PetscLogObjectParent(mat,from);CHKERRQ(ierr); 13452e6d16bSBarry Smith ierr = PetscLogObjectParent(mat,to);CHKERRQ(ierr); 135d9653453SSatish Balay baij->garray = garray; 13652e6d16bSBarry Smith ierr = PetscLogObjectMemory(mat,(ec+1)*sizeof(PetscInt));CHKERRQ(ierr); 1378016bdd1SSatish Balay ierr = ISDestroy(from);CHKERRQ(ierr); 1388016bdd1SSatish Balay ierr = ISDestroy(to);CHKERRQ(ierr); 139888f2ed8SSatish Balay ierr = VecDestroy(gvec);CHKERRQ(ierr); 1403a40ed3dSBarry Smith PetscFunctionReturn(0); 1418016bdd1SSatish Balay } 1428016bdd1SSatish Balay 1438016bdd1SSatish Balay /* 144d9653453SSatish Balay Takes the local part of an already assembled MPIBAIJ matrix 1458016bdd1SSatish Balay and disassembles it. This is to allow new nonzeros into the matrix 1468016bdd1SSatish Balay that require more communication in the matrix vector multiply. 1478016bdd1SSatish Balay Thus certain data-structures must be rebuilt. 1488016bdd1SSatish Balay 1498016bdd1SSatish Balay Kind of slow! But that's what application programmers get when 1508016bdd1SSatish Balay they are sloppy. 1518016bdd1SSatish Balay */ 1524a2ae208SSatish Balay #undef __FUNCT__ 1534a2ae208SSatish Balay #define __FUNCT__ "DisAssemble_MPIBAIJ" 154dfbe8321SBarry Smith PetscErrorCode DisAssemble_MPIBAIJ(Mat A) 1558016bdd1SSatish Balay { 156d9653453SSatish Balay Mat_MPIBAIJ *baij = (Mat_MPIBAIJ*)A->data; 157d9653453SSatish Balay Mat B = baij->B,Bnew; 158d9653453SSatish Balay Mat_SeqBAIJ *Bbaij = (Mat_SeqBAIJ*)B->data; 159dfbe8321SBarry Smith PetscErrorCode ierr; 160d0f46423SBarry Smith PetscInt i,j,mbs=Bbaij->mbs,n = A->cmap->N,col,*garray=baij->garray; 161d0f46423SBarry Smith PetscInt bs2 = baij->bs2,*nz,ec,m = A->rmap->n; 1623eda8832SBarry Smith MatScalar *a = Bbaij->a; 163dd6ea824SBarry Smith MatScalar *atmp; 16497e5c40aSBarry Smith 1658016bdd1SSatish Balay 1663a40ed3dSBarry Smith PetscFunctionBegin; 1678016bdd1SSatish Balay /* free stuff related to matrix-vec multiply */ 168b0a32e0cSBarry Smith ierr = VecGetSize(baij->lvec,&ec);CHKERRQ(ierr); /* needed for PetscLogObjectMemory below */ 169d9653453SSatish Balay ierr = VecDestroy(baij->lvec);CHKERRQ(ierr); baij->lvec = 0; 170d9653453SSatish Balay ierr = VecScatterDestroy(baij->Mvctx);CHKERRQ(ierr); baij->Mvctx = 0; 171d9653453SSatish Balay if (baij->colmap) { 172aa482453SBarry Smith #if defined (PETSC_USE_CTABLE) 1739c666560SBarry Smith ierr = PetscTableDestroy(baij->colmap); baij->colmap = 0;CHKERRQ(ierr); 17448e59246SSatish Balay #else 175606d414cSSatish Balay ierr = PetscFree(baij->colmap);CHKERRQ(ierr); 176606d414cSSatish Balay baij->colmap = 0; 17752e6d16bSBarry Smith ierr = PetscLogObjectMemory(A,-Bbaij->nbs*sizeof(PetscInt));CHKERRQ(ierr); 17848e59246SSatish Balay #endif 1798016bdd1SSatish Balay } 1808016bdd1SSatish Balay 1818016bdd1SSatish Balay /* make sure that B is assembled so we can access its values */ 1826d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1833eda8832SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1848016bdd1SSatish Balay 1858016bdd1SSatish Balay /* invent new B and copy stuff over */ 186b24ad042SBarry Smith ierr = PetscMalloc(mbs*sizeof(PetscInt),&nz);CHKERRQ(ierr); 187d9653453SSatish Balay for (i=0; i<mbs; i++) { 188d9653453SSatish Balay nz[i] = Bbaij->i[i+1]-Bbaij->i[i]; 1898016bdd1SSatish Balay } 1907adad957SLisandro Dalcin ierr = MatCreate(((PetscObject)B)->comm,&Bnew);CHKERRQ(ierr); 191f69a0ea3SMatthew Knepley ierr = MatSetSizes(Bnew,m,n,m,n);CHKERRQ(ierr); 1927adad957SLisandro Dalcin ierr = MatSetType(Bnew,((PetscObject)B)->type_name);CHKERRQ(ierr); 193d0f46423SBarry Smith ierr = MatSeqBAIJSetPreallocation(Bnew,B->rmap->bs,0,nz);CHKERRQ(ierr); 1944e0d8c25SBarry Smith ierr = MatSetOption(Bnew,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr); 195d9653453SSatish Balay 196bba1ac68SSatish Balay for (i=0; i<mbs; i++) { 197bba1ac68SSatish Balay for (j=Bbaij->i[i]; j<Bbaij->i[i+1]; j++) { 198bba1ac68SSatish Balay col = garray[Bbaij->j[j]]; 1993eda8832SBarry Smith atmp = a + j*bs2; 200bba1ac68SSatish Balay ierr = MatSetValuesBlocked_SeqBAIJ(Bnew,1,&i,1,&col,atmp,B->insertmode);CHKERRQ(ierr); 2018016bdd1SSatish Balay } 2028016bdd1SSatish Balay } 2034e0d8c25SBarry Smith ierr = MatSetOption(Bnew,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr); 204bba1ac68SSatish Balay 205bba1ac68SSatish Balay ierr = PetscFree(nz);CHKERRQ(ierr); 206606d414cSSatish Balay ierr = PetscFree(baij->garray);CHKERRQ(ierr); 207606d414cSSatish Balay baij->garray = 0; 20852e6d16bSBarry Smith ierr = PetscLogObjectMemory(A,-ec*sizeof(PetscInt));CHKERRQ(ierr); 2098016bdd1SSatish Balay ierr = MatDestroy(B);CHKERRQ(ierr); 21052e6d16bSBarry Smith ierr = PetscLogObjectParent(A,Bnew);CHKERRQ(ierr); 211d9653453SSatish Balay baij->B = Bnew; 2128016bdd1SSatish Balay A->was_assembled = PETSC_FALSE; 2133a40ed3dSBarry Smith PetscFunctionReturn(0); 2148016bdd1SSatish Balay } 2158016bdd1SSatish Balay 21604f1ad80SBarry Smith /* ugly stuff added for Glenn someday we should fix this up */ 21704f1ad80SBarry Smith 218b24ad042SBarry Smith static PetscInt *uglyrmapd = 0,*uglyrmapo = 0; /* mapping from the local ordering to the "diagonal" and "off-diagonal" 21904f1ad80SBarry Smith parts of the local matrix */ 2202cd6534aSBarry Smith static Vec uglydd = 0,uglyoo = 0; /* work vectors used to scale the two parts of the local matrix */ 22104f1ad80SBarry Smith 22204f1ad80SBarry Smith 223348d1a9dSSatish Balay #undef __FUNCT__ 224348d1a9dSSatish Balay #define __FUNCT__ "MatMPIBAIJDiagonalScaleLocalSetUp" 225dfbe8321SBarry Smith PetscErrorCode MatMPIBAIJDiagonalScaleLocalSetUp(Mat inA,Vec scale) 22604f1ad80SBarry Smith { 22704f1ad80SBarry Smith Mat_MPIBAIJ *ina = (Mat_MPIBAIJ*) inA->data; /*access private part of matrix */ 22804f1ad80SBarry Smith Mat_SeqBAIJ *B = (Mat_SeqBAIJ*)ina->B->data; 229dfbe8321SBarry Smith PetscErrorCode ierr; 230d0f46423SBarry Smith PetscInt bs = inA->rmap->bs,i,n,nt,j,cstart,cend,no,*garray = ina->garray,*lindices; 231b24ad042SBarry Smith PetscInt *r_rmapd,*r_rmapo; 23204f1ad80SBarry Smith 23304f1ad80SBarry Smith PetscFunctionBegin; 23404f1ad80SBarry Smith ierr = MatGetOwnershipRange(inA,&cstart,&cend);CHKERRQ(ierr); 23504f1ad80SBarry Smith ierr = MatGetSize(ina->A,PETSC_NULL,&n);CHKERRQ(ierr); 236b24ad042SBarry Smith ierr = PetscMalloc((inA->bmapping->n+1)*sizeof(PetscInt),&r_rmapd);CHKERRQ(ierr); 237b24ad042SBarry Smith ierr = PetscMemzero(r_rmapd,inA->bmapping->n*sizeof(PetscInt));CHKERRQ(ierr); 23804f1ad80SBarry Smith nt = 0; 23904f1ad80SBarry Smith for (i=0; i<inA->bmapping->n; i++) { 24004f1ad80SBarry Smith if (inA->bmapping->indices[i]*bs >= cstart && inA->bmapping->indices[i]*bs < cend) { 24104f1ad80SBarry Smith nt++; 24204f1ad80SBarry Smith r_rmapd[i] = inA->bmapping->indices[i] + 1; 24304f1ad80SBarry Smith } 24404f1ad80SBarry Smith } 24577431f27SBarry Smith if (nt*bs != n) SETERRQ2(PETSC_ERR_PLIB,"Hmm nt*bs %D n %D",nt*bs,n); 246b24ad042SBarry Smith ierr = PetscMalloc((n+1)*sizeof(PetscInt),&uglyrmapd);CHKERRQ(ierr); 24704f1ad80SBarry Smith for (i=0; i<inA->bmapping->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 257b24ad042SBarry Smith ierr = PetscMalloc((ina->Nbs+1)*sizeof(PetscInt),&lindices);CHKERRQ(ierr); 258b24ad042SBarry Smith ierr = PetscMemzero(lindices,ina->Nbs*sizeof(PetscInt));CHKERRQ(ierr); 25904f1ad80SBarry Smith for (i=0; i<B->nbs; i++) { 26004f1ad80SBarry Smith lindices[garray[i]] = i+1; 26104f1ad80SBarry Smith } 26204f1ad80SBarry Smith no = inA->bmapping->n - nt; 263b24ad042SBarry Smith ierr = PetscMalloc((inA->bmapping->n+1)*sizeof(PetscInt),&r_rmapo);CHKERRQ(ierr); 264b24ad042SBarry Smith ierr = PetscMemzero(r_rmapo,inA->bmapping->n*sizeof(PetscInt));CHKERRQ(ierr); 26504f1ad80SBarry Smith nt = 0; 26604f1ad80SBarry Smith for (i=0; i<inA->bmapping->n; i++) { 26704f1ad80SBarry Smith if (lindices[inA->bmapping->indices[i]]) { 26804f1ad80SBarry Smith nt++; 26904f1ad80SBarry Smith r_rmapo[i] = lindices[inA->bmapping->indices[i]]; 27004f1ad80SBarry Smith } 27104f1ad80SBarry Smith } 27277431f27SBarry Smith if (nt > no) SETERRQ2(PETSC_ERR_PLIB,"Hmm nt %D no %D",nt,n); 27304f1ad80SBarry Smith ierr = PetscFree(lindices);CHKERRQ(ierr); 274b24ad042SBarry Smith ierr = PetscMalloc((nt*bs+1)*sizeof(PetscInt),&uglyrmapo);CHKERRQ(ierr); 27504f1ad80SBarry Smith for (i=0; i<inA->bmapping->n; i++) { 27604f1ad80SBarry Smith if (r_rmapo[i]){ 27704f1ad80SBarry Smith for (j=0; j<bs; j++) { 27804f1ad80SBarry Smith uglyrmapo[(r_rmapo[i]-1)*bs+j] = i*bs + j; 27904f1ad80SBarry Smith } 28004f1ad80SBarry Smith } 28104f1ad80SBarry Smith } 28204f1ad80SBarry Smith ierr = PetscFree(r_rmapo);CHKERRQ(ierr); 28304f1ad80SBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,nt*bs,&uglyoo);CHKERRQ(ierr); 28404f1ad80SBarry Smith 28504f1ad80SBarry Smith PetscFunctionReturn(0); 28604f1ad80SBarry Smith } 28704f1ad80SBarry Smith 288348d1a9dSSatish Balay #undef __FUNCT__ 289348d1a9dSSatish Balay #define __FUNCT__ "MatMPIBAIJDiagonalScaleLocal" 290be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJDiagonalScaleLocal(Mat A,Vec scale) 29104f1ad80SBarry Smith { 29292b32695SKris Buschelman /* This routine should really be abandoned as it duplicates MatDiagonalScaleLocal */ 293dfbe8321SBarry Smith PetscErrorCode ierr,(*f)(Mat,Vec); 29492b32695SKris Buschelman 29592b32695SKris Buschelman PetscFunctionBegin; 29692b32695SKris Buschelman ierr = PetscObjectQueryFunction((PetscObject)A,"MatDiagonalScaleLocal_C",(void (**)(void))&f);CHKERRQ(ierr); 29792b32695SKris Buschelman if (f) { 29892b32695SKris Buschelman ierr = (*f)(A,scale);CHKERRQ(ierr); 29992b32695SKris Buschelman } 30092b32695SKris Buschelman PetscFunctionReturn(0); 30192b32695SKris Buschelman } 30292b32695SKris Buschelman 30392b32695SKris Buschelman EXTERN_C_BEGIN 30492b32695SKris Buschelman #undef __FUNCT__ 30592b32695SKris Buschelman #define __FUNCT__ "MatDiagonalScaleLocal_MPIBAIJ" 306be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalScaleLocal_MPIBAIJ(Mat A,Vec scale) 30792b32695SKris Buschelman { 30804f1ad80SBarry Smith Mat_MPIBAIJ *a = (Mat_MPIBAIJ*) A->data; /*access private part of matrix */ 309dfbe8321SBarry Smith PetscErrorCode ierr; 310b24ad042SBarry Smith PetscInt n,i; 31104f1ad80SBarry Smith PetscScalar *d,*o,*s; 31204f1ad80SBarry Smith 31304f1ad80SBarry Smith PetscFunctionBegin; 3142cd6534aSBarry Smith if (!uglyrmapd) { 3152cd6534aSBarry Smith ierr = MatMPIBAIJDiagonalScaleLocalSetUp(A,scale);CHKERRQ(ierr); 3162cd6534aSBarry Smith } 3172cd6534aSBarry Smith 3181ebc52fbSHong Zhang ierr = VecGetArray(scale,&s);CHKERRQ(ierr); 31904f1ad80SBarry Smith 32004f1ad80SBarry Smith ierr = VecGetLocalSize(uglydd,&n);CHKERRQ(ierr); 3211ebc52fbSHong Zhang ierr = VecGetArray(uglydd,&d);CHKERRQ(ierr); 32204f1ad80SBarry Smith for (i=0; i<n; i++) { 32304f1ad80SBarry Smith d[i] = s[uglyrmapd[i]]; /* copy "diagonal" (true local) portion of scale into dd vector */ 32404f1ad80SBarry Smith } 3251ebc52fbSHong Zhang ierr = VecRestoreArray(uglydd,&d);CHKERRQ(ierr); 32604f1ad80SBarry Smith /* column scale "diagonal" portion of local matrix */ 32704f1ad80SBarry Smith ierr = MatDiagonalScale(a->A,PETSC_NULL,uglydd);CHKERRQ(ierr); 32804f1ad80SBarry Smith 32904f1ad80SBarry Smith ierr = VecGetLocalSize(uglyoo,&n);CHKERRQ(ierr); 3301ebc52fbSHong Zhang ierr = VecGetArray(uglyoo,&o);CHKERRQ(ierr); 33104f1ad80SBarry Smith for (i=0; i<n; i++) { 33204f1ad80SBarry Smith o[i] = s[uglyrmapo[i]]; /* copy "off-diagonal" portion of scale into oo vector */ 33304f1ad80SBarry Smith } 3341ebc52fbSHong Zhang ierr = VecRestoreArray(scale,&s);CHKERRQ(ierr); 3351ebc52fbSHong Zhang ierr = VecRestoreArray(uglyoo,&o);CHKERRQ(ierr); 33604f1ad80SBarry Smith /* column scale "off-diagonal" portion of local matrix */ 33704f1ad80SBarry Smith ierr = MatDiagonalScale(a->B,PETSC_NULL,uglyoo);CHKERRQ(ierr); 33804f1ad80SBarry Smith 33904f1ad80SBarry Smith PetscFunctionReturn(0); 34004f1ad80SBarry Smith } 34192b32695SKris Buschelman EXTERN_C_END 3428016bdd1SSatish Balay 3433eda8832SBarry Smith 344