1be1d678aSKris Buschelman #define PETSCMAT_DLL 2be1d678aSKris Buschelman 38016bdd1SSatish Balay /* 4d9653453SSatish Balay Support for the parallel BAIJ matrix vector multiply 58016bdd1SSatish Balay */ 670f55243SBarry 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; 18899cda47SBarry 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; 676857c123SSatish Balay baij->B->cmap.n = baij->B->cmap.N = ec*mat->rmap.bs; 686148ca0dSBarry 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; 1046857c123SSatish Balay baij->B->cmap.n =baij->B->cmap.N = ec*mat->rmap.bs; 1056148ca0dSBarry 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 */ 1278016bdd1SSatish Balay /* this is inefficient, but otherwise we must do either 1288016bdd1SSatish Balay 1) save garray until the first actual scatter when the vector is known or 1298016bdd1SSatish Balay 2) have another way of generating a scatter context without a vector.*/ 130899cda47SBarry Smith ierr = VecCreateMPI(mat->comm,mat->cmap.n,mat->cmap.N,&gvec);CHKERRQ(ierr); 1318016bdd1SSatish Balay 132d9653453SSatish Balay ierr = VecScatterCreate(gvec,from,baij->lvec,to,&baij->Mvctx);CHKERRQ(ierr); 13390f02eecSBarry Smith 13452e6d16bSBarry Smith ierr = PetscLogObjectParent(mat,baij->Mvctx);CHKERRQ(ierr); 13552e6d16bSBarry Smith ierr = PetscLogObjectParent(mat,baij->lvec);CHKERRQ(ierr); 13652e6d16bSBarry Smith ierr = PetscLogObjectParent(mat,from);CHKERRQ(ierr); 13752e6d16bSBarry Smith ierr = PetscLogObjectParent(mat,to);CHKERRQ(ierr); 138d9653453SSatish Balay baij->garray = garray; 13952e6d16bSBarry Smith ierr = PetscLogObjectMemory(mat,(ec+1)*sizeof(PetscInt));CHKERRQ(ierr); 1408016bdd1SSatish Balay ierr = ISDestroy(from);CHKERRQ(ierr); 1418016bdd1SSatish Balay ierr = ISDestroy(to);CHKERRQ(ierr); 142888f2ed8SSatish Balay ierr = VecDestroy(gvec);CHKERRQ(ierr); 1433a40ed3dSBarry Smith PetscFunctionReturn(0); 1448016bdd1SSatish Balay } 1458016bdd1SSatish Balay 1468016bdd1SSatish Balay /* 147d9653453SSatish Balay Takes the local part of an already assembled MPIBAIJ matrix 1488016bdd1SSatish Balay and disassembles it. This is to allow new nonzeros into the matrix 1498016bdd1SSatish Balay that require more communication in the matrix vector multiply. 1508016bdd1SSatish Balay Thus certain data-structures must be rebuilt. 1518016bdd1SSatish Balay 1528016bdd1SSatish Balay Kind of slow! But that's what application programmers get when 1538016bdd1SSatish Balay they are sloppy. 1548016bdd1SSatish Balay */ 1554a2ae208SSatish Balay #undef __FUNCT__ 1564a2ae208SSatish Balay #define __FUNCT__ "DisAssemble_MPIBAIJ" 157dfbe8321SBarry Smith PetscErrorCode DisAssemble_MPIBAIJ(Mat A) 1588016bdd1SSatish Balay { 159d9653453SSatish Balay Mat_MPIBAIJ *baij = (Mat_MPIBAIJ*)A->data; 160d9653453SSatish Balay Mat B = baij->B,Bnew; 161d9653453SSatish Balay Mat_SeqBAIJ *Bbaij = (Mat_SeqBAIJ*)B->data; 162dfbe8321SBarry Smith PetscErrorCode ierr; 163899cda47SBarry Smith PetscInt i,j,mbs=Bbaij->mbs,n = A->cmap.n,col,*garray=baij->garray; 164899cda47SBarry Smith PetscInt bs2 = baij->bs2,*nz,ec,m = A->rmap.n; 1653eda8832SBarry Smith MatScalar *a = Bbaij->a; 16687828ca2SBarry Smith PetscScalar *atmp; 16776117effSSatish Balay #if defined(PETSC_USE_MAT_SINGLE) 168b24ad042SBarry Smith PetscInt k; 16976117effSSatish Balay #endif 1708016bdd1SSatish Balay 1713a40ed3dSBarry Smith PetscFunctionBegin; 1728016bdd1SSatish Balay /* free stuff related to matrix-vec multiply */ 173b0a32e0cSBarry Smith ierr = VecGetSize(baij->lvec,&ec);CHKERRQ(ierr); /* needed for PetscLogObjectMemory below */ 174d9653453SSatish Balay ierr = VecDestroy(baij->lvec);CHKERRQ(ierr); baij->lvec = 0; 175d9653453SSatish Balay ierr = VecScatterDestroy(baij->Mvctx);CHKERRQ(ierr); baij->Mvctx = 0; 176d9653453SSatish Balay if (baij->colmap) { 177aa482453SBarry Smith #if defined (PETSC_USE_CTABLE) 1789c666560SBarry Smith ierr = PetscTableDestroy(baij->colmap); baij->colmap = 0;CHKERRQ(ierr); 17948e59246SSatish Balay #else 180606d414cSSatish Balay ierr = PetscFree(baij->colmap);CHKERRQ(ierr); 181606d414cSSatish Balay baij->colmap = 0; 18252e6d16bSBarry Smith ierr = PetscLogObjectMemory(A,-Bbaij->nbs*sizeof(PetscInt));CHKERRQ(ierr); 18348e59246SSatish Balay #endif 1848016bdd1SSatish Balay } 1858016bdd1SSatish Balay 1868016bdd1SSatish Balay /* make sure that B is assembled so we can access its values */ 1876d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1883eda8832SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1898016bdd1SSatish Balay 1908016bdd1SSatish Balay /* invent new B and copy stuff over */ 191b24ad042SBarry Smith ierr = PetscMalloc(mbs*sizeof(PetscInt),&nz);CHKERRQ(ierr); 192d9653453SSatish Balay for (i=0; i<mbs; i++) { 193d9653453SSatish Balay nz[i] = Bbaij->i[i+1]-Bbaij->i[i]; 1948016bdd1SSatish Balay } 195f69a0ea3SMatthew Knepley ierr = MatCreate(B->comm,&Bnew);CHKERRQ(ierr); 196f69a0ea3SMatthew Knepley ierr = MatSetSizes(Bnew,m,n,m,n);CHKERRQ(ierr); 197f204ca49SKris Buschelman ierr = MatSetType(Bnew,B->type_name);CHKERRQ(ierr); 198899cda47SBarry Smith ierr = MatSeqBAIJSetPreallocation(Bnew,B->rmap.bs,0,nz);CHKERRQ(ierr); 199*4e0d8c25SBarry Smith ierr = MatSetOption(Bnew,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr); 200d9653453SSatish Balay 2015010ffefSSatish Balay #if defined(PETSC_USE_MAT_SINGLE) 202bba1ac68SSatish Balay ierr = PetscMalloc(bs2*sizeof(PetscScalar),&atmp);CHKERRQ(ierr); 203bba1ac68SSatish Balay #endif 204bba1ac68SSatish Balay for (i=0; i<mbs; i++) { 205bba1ac68SSatish Balay for (j=Bbaij->i[i]; j<Bbaij->i[i+1]; j++) { 206bba1ac68SSatish Balay col = garray[Bbaij->j[j]]; 207bba1ac68SSatish Balay #if defined(PETSC_USE_MAT_SINGLE) 208bba1ac68SSatish Balay for (k=0; k<bs2; k++) atmp[k] = a[j*bs2+k]; 2093eda8832SBarry Smith #else 2103eda8832SBarry Smith atmp = a + j*bs2; 2113eda8832SBarry Smith #endif 212bba1ac68SSatish Balay ierr = MatSetValuesBlocked_SeqBAIJ(Bnew,1,&i,1,&col,atmp,B->insertmode);CHKERRQ(ierr); 2138016bdd1SSatish Balay } 2148016bdd1SSatish Balay } 215*4e0d8c25SBarry Smith ierr = MatSetOption(Bnew,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr); 216bba1ac68SSatish Balay 2173eda8832SBarry Smith #if defined(PETSC_USE_MAT_SINGLE) 2183eda8832SBarry Smith ierr = PetscFree(atmp);CHKERRQ(ierr); 2193eda8832SBarry Smith #endif 220bba1ac68SSatish Balay 221bba1ac68SSatish Balay ierr = PetscFree(nz);CHKERRQ(ierr); 222606d414cSSatish Balay ierr = PetscFree(baij->garray);CHKERRQ(ierr); 223606d414cSSatish Balay baij->garray = 0; 22452e6d16bSBarry Smith ierr = PetscLogObjectMemory(A,-ec*sizeof(PetscInt));CHKERRQ(ierr); 2258016bdd1SSatish Balay ierr = MatDestroy(B);CHKERRQ(ierr); 22652e6d16bSBarry Smith ierr = PetscLogObjectParent(A,Bnew);CHKERRQ(ierr); 227d9653453SSatish Balay baij->B = Bnew; 2288016bdd1SSatish Balay A->was_assembled = PETSC_FALSE; 2293a40ed3dSBarry Smith PetscFunctionReturn(0); 2308016bdd1SSatish Balay } 2318016bdd1SSatish Balay 23204f1ad80SBarry Smith /* ugly stuff added for Glenn someday we should fix this up */ 23304f1ad80SBarry Smith 234b24ad042SBarry Smith static PetscInt *uglyrmapd = 0,*uglyrmapo = 0; /* mapping from the local ordering to the "diagonal" and "off-diagonal" 23504f1ad80SBarry Smith parts of the local matrix */ 2362cd6534aSBarry Smith static Vec uglydd = 0,uglyoo = 0; /* work vectors used to scale the two parts of the local matrix */ 23704f1ad80SBarry Smith 23804f1ad80SBarry Smith 239348d1a9dSSatish Balay #undef __FUNCT__ 240348d1a9dSSatish Balay #define __FUNCT__ "MatMPIBAIJDiagonalScaleLocalSetUp" 241dfbe8321SBarry Smith PetscErrorCode MatMPIBAIJDiagonalScaleLocalSetUp(Mat inA,Vec scale) 24204f1ad80SBarry Smith { 24304f1ad80SBarry Smith Mat_MPIBAIJ *ina = (Mat_MPIBAIJ*) inA->data; /*access private part of matrix */ 24404f1ad80SBarry Smith Mat_SeqBAIJ *B = (Mat_SeqBAIJ*)ina->B->data; 245dfbe8321SBarry Smith PetscErrorCode ierr; 246899cda47SBarry Smith PetscInt bs = inA->rmap.bs,i,n,nt,j,cstart,cend,no,*garray = ina->garray,*lindices; 247b24ad042SBarry Smith PetscInt *r_rmapd,*r_rmapo; 24804f1ad80SBarry Smith 24904f1ad80SBarry Smith PetscFunctionBegin; 25004f1ad80SBarry Smith ierr = MatGetOwnershipRange(inA,&cstart,&cend);CHKERRQ(ierr); 25104f1ad80SBarry Smith ierr = MatGetSize(ina->A,PETSC_NULL,&n);CHKERRQ(ierr); 252b24ad042SBarry Smith ierr = PetscMalloc((inA->bmapping->n+1)*sizeof(PetscInt),&r_rmapd);CHKERRQ(ierr); 253b24ad042SBarry Smith ierr = PetscMemzero(r_rmapd,inA->bmapping->n*sizeof(PetscInt));CHKERRQ(ierr); 25404f1ad80SBarry Smith nt = 0; 25504f1ad80SBarry Smith for (i=0; i<inA->bmapping->n; i++) { 25604f1ad80SBarry Smith if (inA->bmapping->indices[i]*bs >= cstart && inA->bmapping->indices[i]*bs < cend) { 25704f1ad80SBarry Smith nt++; 25804f1ad80SBarry Smith r_rmapd[i] = inA->bmapping->indices[i] + 1; 25904f1ad80SBarry Smith } 26004f1ad80SBarry Smith } 26177431f27SBarry Smith if (nt*bs != n) SETERRQ2(PETSC_ERR_PLIB,"Hmm nt*bs %D n %D",nt*bs,n); 262b24ad042SBarry Smith ierr = PetscMalloc((n+1)*sizeof(PetscInt),&uglyrmapd);CHKERRQ(ierr); 26304f1ad80SBarry Smith for (i=0; i<inA->bmapping->n; i++) { 26404f1ad80SBarry Smith if (r_rmapd[i]){ 26504f1ad80SBarry Smith for (j=0; j<bs; j++) { 26604f1ad80SBarry Smith uglyrmapd[(r_rmapd[i]-1)*bs+j-cstart] = i*bs + j; 26704f1ad80SBarry Smith } 26804f1ad80SBarry Smith } 26904f1ad80SBarry Smith } 27004f1ad80SBarry Smith ierr = PetscFree(r_rmapd);CHKERRQ(ierr); 27104f1ad80SBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,n,&uglydd);CHKERRQ(ierr); 27204f1ad80SBarry Smith 273b24ad042SBarry Smith ierr = PetscMalloc((ina->Nbs+1)*sizeof(PetscInt),&lindices);CHKERRQ(ierr); 274b24ad042SBarry Smith ierr = PetscMemzero(lindices,ina->Nbs*sizeof(PetscInt));CHKERRQ(ierr); 27504f1ad80SBarry Smith for (i=0; i<B->nbs; i++) { 27604f1ad80SBarry Smith lindices[garray[i]] = i+1; 27704f1ad80SBarry Smith } 27804f1ad80SBarry Smith no = inA->bmapping->n - nt; 279b24ad042SBarry Smith ierr = PetscMalloc((inA->bmapping->n+1)*sizeof(PetscInt),&r_rmapo);CHKERRQ(ierr); 280b24ad042SBarry Smith ierr = PetscMemzero(r_rmapo,inA->bmapping->n*sizeof(PetscInt));CHKERRQ(ierr); 28104f1ad80SBarry Smith nt = 0; 28204f1ad80SBarry Smith for (i=0; i<inA->bmapping->n; i++) { 28304f1ad80SBarry Smith if (lindices[inA->bmapping->indices[i]]) { 28404f1ad80SBarry Smith nt++; 28504f1ad80SBarry Smith r_rmapo[i] = lindices[inA->bmapping->indices[i]]; 28604f1ad80SBarry Smith } 28704f1ad80SBarry Smith } 28877431f27SBarry Smith if (nt > no) SETERRQ2(PETSC_ERR_PLIB,"Hmm nt %D no %D",nt,n); 28904f1ad80SBarry Smith ierr = PetscFree(lindices);CHKERRQ(ierr); 290b24ad042SBarry Smith ierr = PetscMalloc((nt*bs+1)*sizeof(PetscInt),&uglyrmapo);CHKERRQ(ierr); 29104f1ad80SBarry Smith for (i=0; i<inA->bmapping->n; i++) { 29204f1ad80SBarry Smith if (r_rmapo[i]){ 29304f1ad80SBarry Smith for (j=0; j<bs; j++) { 29404f1ad80SBarry Smith uglyrmapo[(r_rmapo[i]-1)*bs+j] = i*bs + j; 29504f1ad80SBarry Smith } 29604f1ad80SBarry Smith } 29704f1ad80SBarry Smith } 29804f1ad80SBarry Smith ierr = PetscFree(r_rmapo);CHKERRQ(ierr); 29904f1ad80SBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,nt*bs,&uglyoo);CHKERRQ(ierr); 30004f1ad80SBarry Smith 30104f1ad80SBarry Smith PetscFunctionReturn(0); 30204f1ad80SBarry Smith } 30304f1ad80SBarry Smith 304348d1a9dSSatish Balay #undef __FUNCT__ 305348d1a9dSSatish Balay #define __FUNCT__ "MatMPIBAIJDiagonalScaleLocal" 306be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJDiagonalScaleLocal(Mat A,Vec scale) 30704f1ad80SBarry Smith { 30892b32695SKris Buschelman /* This routine should really be abandoned as it duplicates MatDiagonalScaleLocal */ 309dfbe8321SBarry Smith PetscErrorCode ierr,(*f)(Mat,Vec); 31092b32695SKris Buschelman 31192b32695SKris Buschelman PetscFunctionBegin; 31292b32695SKris Buschelman ierr = PetscObjectQueryFunction((PetscObject)A,"MatDiagonalScaleLocal_C",(void (**)(void))&f);CHKERRQ(ierr); 31392b32695SKris Buschelman if (f) { 31492b32695SKris Buschelman ierr = (*f)(A,scale);CHKERRQ(ierr); 31592b32695SKris Buschelman } 31692b32695SKris Buschelman PetscFunctionReturn(0); 31792b32695SKris Buschelman } 31892b32695SKris Buschelman 31992b32695SKris Buschelman EXTERN_C_BEGIN 32092b32695SKris Buschelman #undef __FUNCT__ 32192b32695SKris Buschelman #define __FUNCT__ "MatDiagonalScaleLocal_MPIBAIJ" 322be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalScaleLocal_MPIBAIJ(Mat A,Vec scale) 32392b32695SKris Buschelman { 32404f1ad80SBarry Smith Mat_MPIBAIJ *a = (Mat_MPIBAIJ*) A->data; /*access private part of matrix */ 325dfbe8321SBarry Smith PetscErrorCode ierr; 326b24ad042SBarry Smith PetscInt n,i; 32704f1ad80SBarry Smith PetscScalar *d,*o,*s; 32804f1ad80SBarry Smith 32904f1ad80SBarry Smith PetscFunctionBegin; 3302cd6534aSBarry Smith if (!uglyrmapd) { 3312cd6534aSBarry Smith ierr = MatMPIBAIJDiagonalScaleLocalSetUp(A,scale);CHKERRQ(ierr); 3322cd6534aSBarry Smith } 3332cd6534aSBarry Smith 3341ebc52fbSHong Zhang ierr = VecGetArray(scale,&s);CHKERRQ(ierr); 33504f1ad80SBarry Smith 33604f1ad80SBarry Smith ierr = VecGetLocalSize(uglydd,&n);CHKERRQ(ierr); 3371ebc52fbSHong Zhang ierr = VecGetArray(uglydd,&d);CHKERRQ(ierr); 33804f1ad80SBarry Smith for (i=0; i<n; i++) { 33904f1ad80SBarry Smith d[i] = s[uglyrmapd[i]]; /* copy "diagonal" (true local) portion of scale into dd vector */ 34004f1ad80SBarry Smith } 3411ebc52fbSHong Zhang ierr = VecRestoreArray(uglydd,&d);CHKERRQ(ierr); 34204f1ad80SBarry Smith /* column scale "diagonal" portion of local matrix */ 34304f1ad80SBarry Smith ierr = MatDiagonalScale(a->A,PETSC_NULL,uglydd);CHKERRQ(ierr); 34404f1ad80SBarry Smith 34504f1ad80SBarry Smith ierr = VecGetLocalSize(uglyoo,&n);CHKERRQ(ierr); 3461ebc52fbSHong Zhang ierr = VecGetArray(uglyoo,&o);CHKERRQ(ierr); 34704f1ad80SBarry Smith for (i=0; i<n; i++) { 34804f1ad80SBarry Smith o[i] = s[uglyrmapo[i]]; /* copy "off-diagonal" portion of scale into oo vector */ 34904f1ad80SBarry Smith } 3501ebc52fbSHong Zhang ierr = VecRestoreArray(scale,&s);CHKERRQ(ierr); 3511ebc52fbSHong Zhang ierr = VecRestoreArray(uglyoo,&o);CHKERRQ(ierr); 35204f1ad80SBarry Smith /* column scale "off-diagonal" portion of local matrix */ 35304f1ad80SBarry Smith ierr = MatDiagonalScale(a->B,PETSC_NULL,uglyoo);CHKERRQ(ierr); 35404f1ad80SBarry Smith 35504f1ad80SBarry Smith PetscFunctionReturn(0); 35604f1ad80SBarry Smith } 35792b32695SKris Buschelman EXTERN_C_END 3588016bdd1SSatish Balay 3593eda8832SBarry Smith 360