1bba1ac68SSatish Balay /*$Id: mmbaij.c,v 1.46 2001/09/25 00:31:36 balay Exp $*/ 28016bdd1SSatish Balay 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 87cff1425SSatish Balay EXTERN int MatSetValuesBlocked_SeqBAIJ(Mat,int,const int[],int,const int[],const PetscScalar[],InsertMode); 98016bdd1SSatish Balay 104a2ae208SSatish Balay #undef __FUNCT__ 114a2ae208SSatish Balay #define __FUNCT__ "MatSetUpMultiply_MPIBAIJ" 12d9653453SSatish Balay int 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); 16d9653453SSatish Balay int Nbs = baij->Nbs,i,j,*indices,*aj = B->j,ierr,ec = 0,*garray; 173d3cf644SBarry Smith int bs = baij->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; 2373a2e727SSatish Balay int gid,lid; 2473a2e727SSatish Balay #endif 258016bdd1SSatish Balay 263a40ed3dSBarry Smith PetscFunctionBegin; 2773a2e727SSatish Balay 28aa482453SBarry Smith #if defined (PETSC_USE_CTABLE) 2973a2e727SSatish Balay /* use a table - Mark Adams */ 30273d9f13SBarry Smith ierr = PetscTableCreate(B->mbs,&gid1_lid1);CHKERRQ(ierr); 3173a2e727SSatish Balay for (i=0; i<B->mbs; i++) { 3273a2e727SSatish Balay for (j=0; j<B->ilen[i]; j++) { 33fa46199cSSatish Balay int data,gid1 = aj[B->i[i]+j] + 1; 340f5bd95cSBarry Smith ierr = PetscTableFind(gid1_lid1,gid1,&data) ;CHKERRQ(ierr); 35fa46199cSSatish Balay if (!data) { 3673a2e727SSatish Balay /* one based table */ 370f5bd95cSBarry Smith ierr = PetscTableAdd(gid1_lid1,gid1,++ec);CHKERRQ(ierr); 3873a2e727SSatish Balay } 3973a2e727SSatish Balay } 4073a2e727SSatish Balay } 4173a2e727SSatish Balay /* form array of columns we need */ 42b0a32e0cSBarry Smith ierr = PetscMalloc((ec+1)*sizeof(int),&garray);CHKERRQ(ierr); 430f5bd95cSBarry Smith ierr = PetscTableGetHeadPosition(gid1_lid1,&tpos);CHKERRQ(ierr); 4473a2e727SSatish Balay while (tpos) { 450f5bd95cSBarry Smith ierr = PetscTableGetNext(gid1_lid1,&tpos,&gid,&lid);CHKERRQ(ierr); 4673a2e727SSatish Balay gid--; lid--; 4773a2e727SSatish Balay garray[lid] = gid; 4873a2e727SSatish Balay } 490064e2bbSSatish Balay ierr = PetscSortInt(ec,garray);CHKERRQ(ierr); 500f5bd95cSBarry Smith ierr = PetscTableRemoveAll(gid1_lid1);CHKERRQ(ierr); 5173a2e727SSatish Balay for (i=0; i<ec; i++) { 520f5bd95cSBarry Smith ierr = PetscTableAdd(gid1_lid1,garray[i]+1,i+1);CHKERRQ(ierr); 5373a2e727SSatish Balay } 5473a2e727SSatish Balay /* compact out the extra columns in B */ 5573a2e727SSatish Balay for (i=0; i<B->mbs; i++) { 5673a2e727SSatish Balay for (j=0; j<B->ilen[i]; j++) { 5773a2e727SSatish Balay int gid1 = aj[B->i[i] + j] + 1; 580f5bd95cSBarry Smith ierr = PetscTableFind(gid1_lid1,gid1,&lid);CHKERRQ(ierr); 59fa46199cSSatish Balay lid --; 6073a2e727SSatish Balay aj[B->i[i]+j] = lid; 6173a2e727SSatish Balay } 6273a2e727SSatish Balay } 6373a2e727SSatish Balay B->nbs = ec; 64273d9f13SBarry Smith baij->B->n = ec*B->bs; 650f5bd95cSBarry Smith ierr = PetscTableDelete(gid1_lid1);CHKERRQ(ierr); 6673a2e727SSatish Balay /* Mark Adams */ 6773a2e727SSatish Balay #else 68435da068SBarry Smith /* Make an array as long as the number of columns */ 69d9653453SSatish Balay /* mark those columns that are in baij->B */ 70b0a32e0cSBarry Smith ierr = PetscMalloc((Nbs+1)*sizeof(int),&indices);CHKERRQ(ierr); 71549d3d68SSatish Balay ierr = PetscMemzero(indices,Nbs*sizeof(int));CHKERRQ(ierr); 72d9653453SSatish Balay for (i=0; i<B->mbs; i++) { 738016bdd1SSatish Balay for (j=0; j<B->ilen[i]; j++) { 74d9653453SSatish Balay if (!indices[aj[B->i[i] + j]]) ec++; 75d9653453SSatish Balay indices[aj[B->i[i] + j]] = 1; 768016bdd1SSatish Balay } 778016bdd1SSatish Balay } 788016bdd1SSatish Balay 798016bdd1SSatish Balay /* form array of columns we need */ 80b0a32e0cSBarry Smith ierr = PetscMalloc((ec+1)*sizeof(int),&garray);CHKERRQ(ierr); 818016bdd1SSatish Balay ec = 0; 82d9653453SSatish Balay for (i=0; i<Nbs; i++) { 830bdbc534SSatish Balay if (indices[i]) { 840bdbc534SSatish Balay garray[ec++] = i; 850bdbc534SSatish Balay } 868016bdd1SSatish Balay } 878016bdd1SSatish Balay 888016bdd1SSatish Balay /* make indices now point into garray */ 898016bdd1SSatish Balay for (i=0; i<ec; i++) { 90d9653453SSatish Balay indices[garray[i]] = i; 918016bdd1SSatish Balay } 928016bdd1SSatish Balay 938016bdd1SSatish Balay /* compact out the extra columns in B */ 94d9653453SSatish Balay for (i=0; i<B->mbs; i++) { 958016bdd1SSatish Balay for (j=0; j<B->ilen[i]; j++) { 96d9653453SSatish Balay aj[B->i[i] + j] = indices[aj[B->i[i] + j]]; 978016bdd1SSatish Balay } 988016bdd1SSatish Balay } 99d9653453SSatish Balay B->nbs = ec; 100273d9f13SBarry Smith baij->B->n = ec*B->bs; 101606d414cSSatish Balay ierr = PetscFree(indices);CHKERRQ(ierr); 10273a2e727SSatish Balay #endif 1038016bdd1SSatish Balay 1048016bdd1SSatish Balay /* create local vector that is used to scatter into */ 105029af93fSBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,ec*bs,&baij->lvec);CHKERRQ(ierr); 1068016bdd1SSatish Balay 107c16cb8f2SBarry Smith /* create two temporary index sets for building scatter-gather */ 1083d3cf644SBarry Smith for (i=0; i<ec; i++) { 109c16cb8f2SBarry Smith garray[i] = bs*garray[i]; 110c16cb8f2SBarry Smith } 111029af93fSBarry Smith ierr = ISCreateBlock(PETSC_COMM_SELF,bs,ec,garray,&from);CHKERRQ(ierr); 1123d3cf644SBarry Smith for (i=0; i<ec; i++) { 113c16cb8f2SBarry Smith garray[i] = garray[i]/bs; 114c16cb8f2SBarry Smith } 115c16cb8f2SBarry Smith 116b0a32e0cSBarry Smith ierr = PetscMalloc((ec+1)*sizeof(int),&stmp);CHKERRQ(ierr); 117537820f0SBarry Smith for (i=0; i<ec; i++) { stmp[i] = bs*i; } 118029af93fSBarry Smith ierr = ISCreateBlock(PETSC_COMM_SELF,bs,ec,stmp,&to);CHKERRQ(ierr); 119606d414cSSatish Balay ierr = PetscFree(stmp);CHKERRQ(ierr); 1208016bdd1SSatish Balay 1218016bdd1SSatish Balay /* create temporary global vector to generate scatter context */ 1228016bdd1SSatish Balay /* this is inefficient, but otherwise we must do either 1238016bdd1SSatish Balay 1) save garray until the first actual scatter when the vector is known or 1248016bdd1SSatish Balay 2) have another way of generating a scatter context without a vector.*/ 125273d9f13SBarry Smith ierr = VecCreateMPI(mat->comm,mat->n,mat->N,&gvec);CHKERRQ(ierr); 1268016bdd1SSatish Balay 1278016bdd1SSatish Balay /* gnerate the scatter context */ 128d9653453SSatish Balay ierr = VecScatterCreate(gvec,from,baij->lvec,to,&baij->Mvctx);CHKERRQ(ierr); 12990f02eecSBarry Smith 13090f02eecSBarry Smith /* 13190f02eecSBarry Smith Post the receives for the first matrix vector product. We sync-chronize after 13290f02eecSBarry Smith this on the chance that the user immediately calls MatMult() after assemblying 13390f02eecSBarry Smith the matrix. 13490f02eecSBarry Smith */ 13543a90d84SBarry Smith ierr = VecScatterPostRecvs(gvec,baij->lvec,INSERT_VALUES,SCATTER_FORWARD,baij->Mvctx);CHKERRQ(ierr); 136ca161407SBarry Smith ierr = MPI_Barrier(mat->comm);CHKERRQ(ierr); 13790f02eecSBarry Smith 138b0a32e0cSBarry Smith PetscLogObjectParent(mat,baij->Mvctx); 139b0a32e0cSBarry Smith PetscLogObjectParent(mat,baij->lvec); 140b0a32e0cSBarry Smith PetscLogObjectParent(mat,from); 141b0a32e0cSBarry Smith PetscLogObjectParent(mat,to); 142d9653453SSatish Balay baij->garray = garray; 143b0a32e0cSBarry Smith PetscLogObjectMemory(mat,(ec+1)*sizeof(int)); 1448016bdd1SSatish Balay ierr = ISDestroy(from);CHKERRQ(ierr); 1458016bdd1SSatish Balay ierr = ISDestroy(to);CHKERRQ(ierr); 146888f2ed8SSatish Balay ierr = VecDestroy(gvec);CHKERRQ(ierr); 1473a40ed3dSBarry Smith PetscFunctionReturn(0); 1488016bdd1SSatish Balay } 1498016bdd1SSatish Balay 1508016bdd1SSatish Balay /* 151d9653453SSatish Balay Takes the local part of an already assembled MPIBAIJ matrix 1528016bdd1SSatish Balay and disassembles it. This is to allow new nonzeros into the matrix 1538016bdd1SSatish Balay that require more communication in the matrix vector multiply. 1548016bdd1SSatish Balay Thus certain data-structures must be rebuilt. 1558016bdd1SSatish Balay 1568016bdd1SSatish Balay Kind of slow! But that's what application programmers get when 1578016bdd1SSatish Balay they are sloppy. 1588016bdd1SSatish Balay */ 1594a2ae208SSatish Balay #undef __FUNCT__ 1604a2ae208SSatish Balay #define __FUNCT__ "DisAssemble_MPIBAIJ" 161d9653453SSatish Balay int DisAssemble_MPIBAIJ(Mat A) 1628016bdd1SSatish Balay { 163d9653453SSatish Balay Mat_MPIBAIJ *baij = (Mat_MPIBAIJ*)A->data; 164d9653453SSatish Balay Mat B = baij->B,Bnew; 165d9653453SSatish Balay Mat_SeqBAIJ *Bbaij = (Mat_SeqBAIJ*)B->data; 166273d9f13SBarry Smith int ierr,i,j,mbs=Bbaij->mbs,n = A->N,col,*garray=baij->garray; 167bba1ac68SSatish Balay int bs2 = baij->bs2,*nz,ec,m = A->m; 1683eda8832SBarry Smith MatScalar *a = Bbaij->a; 16987828ca2SBarry Smith PetscScalar *atmp; 17076117effSSatish Balay #if defined(PETSC_USE_MAT_SINGLE) 17176117effSSatish Balay int k; 17276117effSSatish Balay #endif 1738016bdd1SSatish Balay 1743a40ed3dSBarry Smith PetscFunctionBegin; 1758016bdd1SSatish Balay /* free stuff related to matrix-vec multiply */ 176b0a32e0cSBarry Smith ierr = VecGetSize(baij->lvec,&ec);CHKERRQ(ierr); /* needed for PetscLogObjectMemory below */ 177d9653453SSatish Balay ierr = VecDestroy(baij->lvec);CHKERRQ(ierr); baij->lvec = 0; 178d9653453SSatish Balay ierr = VecScatterDestroy(baij->Mvctx);CHKERRQ(ierr); baij->Mvctx = 0; 179d9653453SSatish Balay if (baij->colmap) { 180aa482453SBarry Smith #if defined (PETSC_USE_CTABLE) 1810f5bd95cSBarry Smith ierr = PetscTableDelete(baij->colmap); baij->colmap = 0;CHKERRQ(ierr); 18248e59246SSatish Balay #else 183606d414cSSatish Balay ierr = PetscFree(baij->colmap);CHKERRQ(ierr); 184606d414cSSatish Balay baij->colmap = 0; 185b0a32e0cSBarry Smith PetscLogObjectMemory(A,-Bbaij->nbs*sizeof(int)); 18648e59246SSatish Balay #endif 1878016bdd1SSatish Balay } 1888016bdd1SSatish Balay 1898016bdd1SSatish Balay /* make sure that B is assembled so we can access its values */ 1906d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1913eda8832SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1928016bdd1SSatish Balay 1938016bdd1SSatish Balay /* invent new B and copy stuff over */ 19482502324SSatish Balay ierr = PetscMalloc(mbs*sizeof(int),&nz);CHKERRQ(ierr); 195d9653453SSatish Balay for (i=0; i<mbs; i++) { 196d9653453SSatish Balay nz[i] = Bbaij->i[i+1]-Bbaij->i[i]; 1978016bdd1SSatish Balay } 198*f204ca49SKris Buschelman ierr = MatCreate(B->comm,m,n,m,n,&Bnew);CHKERRQ(ierr); 199*f204ca49SKris Buschelman ierr = MatSetType(Bnew,B->type_name);CHKERRQ(ierr); 200*f204ca49SKris Buschelman ierr = MatSeqBAIJSetPreallocation(Bnew,baij->bs,0,nz);CHKERRQ(ierr); 201bba1ac68SSatish Balay ierr = MatSetOption(Bnew,MAT_COLUMN_ORIENTED);CHKERRQ(ierr); 202d9653453SSatish Balay 2035010ffefSSatish Balay #if defined(PETSC_USE_MAT_SINGLE) 204bba1ac68SSatish Balay ierr = PetscMalloc(bs2*sizeof(PetscScalar),&atmp);CHKERRQ(ierr); 205bba1ac68SSatish Balay #endif 206bba1ac68SSatish Balay for (i=0; i<mbs; i++) { 207bba1ac68SSatish Balay for (j=Bbaij->i[i]; j<Bbaij->i[i+1]; j++) { 208bba1ac68SSatish Balay col = garray[Bbaij->j[j]]; 209bba1ac68SSatish Balay #if defined(PETSC_USE_MAT_SINGLE) 210bba1ac68SSatish Balay for (k=0; k<bs2; k++) atmp[k] = a[j*bs2+k]; 2113eda8832SBarry Smith #else 2123eda8832SBarry Smith atmp = a + j*bs2; 2133eda8832SBarry Smith #endif 214bba1ac68SSatish Balay ierr = MatSetValuesBlocked_SeqBAIJ(Bnew,1,&i,1,&col,atmp,B->insertmode);CHKERRQ(ierr); 2158016bdd1SSatish Balay } 2168016bdd1SSatish Balay } 217bba1ac68SSatish Balay ierr = MatSetOption(Bnew,MAT_ROW_ORIENTED);CHKERRQ(ierr); 218bba1ac68SSatish Balay 2193eda8832SBarry Smith #if defined(PETSC_USE_MAT_SINGLE) 2203eda8832SBarry Smith ierr = PetscFree(atmp);CHKERRQ(ierr); 2213eda8832SBarry Smith #endif 222bba1ac68SSatish Balay 223bba1ac68SSatish Balay ierr = PetscFree(nz);CHKERRQ(ierr); 224606d414cSSatish Balay ierr = PetscFree(baij->garray);CHKERRQ(ierr); 225606d414cSSatish Balay baij->garray = 0; 226b0a32e0cSBarry Smith PetscLogObjectMemory(A,-ec*sizeof(int)); 2278016bdd1SSatish Balay ierr = MatDestroy(B);CHKERRQ(ierr); 228b0a32e0cSBarry Smith PetscLogObjectParent(A,Bnew); 229d9653453SSatish Balay baij->B = Bnew; 2308016bdd1SSatish Balay A->was_assembled = PETSC_FALSE; 2313a40ed3dSBarry Smith PetscFunctionReturn(0); 2328016bdd1SSatish Balay } 2338016bdd1SSatish Balay 23404f1ad80SBarry Smith /* ugly stuff added for Glenn someday we should fix this up */ 23504f1ad80SBarry Smith 2362cd6534aSBarry Smith static int *uglyrmapd = 0,*uglyrmapo = 0; /* mapping from the local ordering to the "diagonal" and "off-diagonal" 23704f1ad80SBarry Smith parts of the local matrix */ 2382cd6534aSBarry Smith static Vec uglydd = 0,uglyoo = 0; /* work vectors used to scale the two parts of the local matrix */ 23904f1ad80SBarry Smith 24004f1ad80SBarry Smith 241348d1a9dSSatish Balay #undef __FUNCT__ 242348d1a9dSSatish Balay #define __FUNCT__ "MatMPIBAIJDiagonalScaleLocalSetUp" 24304f1ad80SBarry Smith int MatMPIBAIJDiagonalScaleLocalSetUp(Mat inA,Vec scale) 24404f1ad80SBarry Smith { 24504f1ad80SBarry Smith Mat_MPIBAIJ *ina = (Mat_MPIBAIJ*) inA->data; /*access private part of matrix */ 24604f1ad80SBarry Smith Mat_SeqBAIJ *A = (Mat_SeqBAIJ*)ina->A->data; 24704f1ad80SBarry Smith Mat_SeqBAIJ *B = (Mat_SeqBAIJ*)ina->B->data; 24804f1ad80SBarry Smith int ierr,bs = A->bs,i,n,nt,j,cstart,cend,no,*garray = ina->garray,*lindices; 24904f1ad80SBarry Smith int *r_rmapd,*r_rmapo; 25004f1ad80SBarry Smith 25104f1ad80SBarry Smith PetscFunctionBegin; 25204f1ad80SBarry Smith ierr = MatGetOwnershipRange(inA,&cstart,&cend);CHKERRQ(ierr); 25304f1ad80SBarry Smith ierr = MatGetSize(ina->A,PETSC_NULL,&n);CHKERRQ(ierr); 25404f1ad80SBarry Smith ierr = PetscMalloc((inA->bmapping->n+1)*sizeof(int),&r_rmapd);CHKERRQ(ierr); 25504f1ad80SBarry Smith ierr = PetscMemzero(r_rmapd,inA->bmapping->n*sizeof(int));CHKERRQ(ierr); 25604f1ad80SBarry Smith nt = 0; 25704f1ad80SBarry Smith for (i=0; i<inA->bmapping->n; i++) { 25804f1ad80SBarry Smith if (inA->bmapping->indices[i]*bs >= cstart && inA->bmapping->indices[i]*bs < cend) { 25904f1ad80SBarry Smith nt++; 26004f1ad80SBarry Smith r_rmapd[i] = inA->bmapping->indices[i] + 1; 26104f1ad80SBarry Smith } 26204f1ad80SBarry Smith } 26304f1ad80SBarry Smith if (nt*bs != n) SETERRQ2(1,"Hmm nt*bs %d n %d",nt*bs,n); 26404f1ad80SBarry Smith ierr = PetscMalloc((n+1)*sizeof(int),&uglyrmapd);CHKERRQ(ierr); 26504f1ad80SBarry Smith for (i=0; i<inA->bmapping->n; i++) { 26604f1ad80SBarry Smith if (r_rmapd[i]){ 26704f1ad80SBarry Smith for (j=0; j<bs; j++) { 26804f1ad80SBarry Smith uglyrmapd[(r_rmapd[i]-1)*bs+j-cstart] = i*bs + j; 26904f1ad80SBarry Smith } 27004f1ad80SBarry Smith } 27104f1ad80SBarry Smith } 27204f1ad80SBarry Smith ierr = PetscFree(r_rmapd);CHKERRQ(ierr); 27304f1ad80SBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,n,&uglydd);CHKERRQ(ierr); 27404f1ad80SBarry Smith 27504f1ad80SBarry Smith ierr = PetscMalloc((ina->Nbs+1)*sizeof(int),&lindices);CHKERRQ(ierr); 27604f1ad80SBarry Smith ierr = PetscMemzero(lindices,ina->Nbs*sizeof(int));CHKERRQ(ierr); 27704f1ad80SBarry Smith for (i=0; i<B->nbs; i++) { 27804f1ad80SBarry Smith lindices[garray[i]] = i+1; 27904f1ad80SBarry Smith } 28004f1ad80SBarry Smith no = inA->bmapping->n - nt; 28104f1ad80SBarry Smith ierr = PetscMalloc((inA->bmapping->n+1)*sizeof(int),&r_rmapo);CHKERRQ(ierr); 28204f1ad80SBarry Smith ierr = PetscMemzero(r_rmapo,inA->bmapping->n*sizeof(int));CHKERRQ(ierr); 28304f1ad80SBarry Smith nt = 0; 28404f1ad80SBarry Smith for (i=0; i<inA->bmapping->n; i++) { 28504f1ad80SBarry Smith if (lindices[inA->bmapping->indices[i]]) { 28604f1ad80SBarry Smith nt++; 28704f1ad80SBarry Smith r_rmapo[i] = lindices[inA->bmapping->indices[i]]; 28804f1ad80SBarry Smith } 28904f1ad80SBarry Smith } 29004f1ad80SBarry Smith if (nt > no) SETERRQ2(1,"Hmm nt %d no %d",nt,n); 29104f1ad80SBarry Smith ierr = PetscFree(lindices);CHKERRQ(ierr); 29204f1ad80SBarry Smith ierr = PetscMalloc((nt*bs+1)*sizeof(int),&uglyrmapo);CHKERRQ(ierr); 29304f1ad80SBarry Smith for (i=0; i<inA->bmapping->n; i++) { 29404f1ad80SBarry Smith if (r_rmapo[i]){ 29504f1ad80SBarry Smith for (j=0; j<bs; j++) { 29604f1ad80SBarry Smith uglyrmapo[(r_rmapo[i]-1)*bs+j] = i*bs + j; 29704f1ad80SBarry Smith } 29804f1ad80SBarry Smith } 29904f1ad80SBarry Smith } 30004f1ad80SBarry Smith ierr = PetscFree(r_rmapo);CHKERRQ(ierr); 30104f1ad80SBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,nt*bs,&uglyoo);CHKERRQ(ierr); 30204f1ad80SBarry Smith 30304f1ad80SBarry Smith PetscFunctionReturn(0); 30404f1ad80SBarry Smith } 30504f1ad80SBarry Smith 306348d1a9dSSatish Balay #undef __FUNCT__ 307348d1a9dSSatish Balay #define __FUNCT__ "MatMPIBAIJDiagonalScaleLocal" 30804f1ad80SBarry Smith int MatMPIBAIJDiagonalScaleLocal(Mat A,Vec scale) 30904f1ad80SBarry Smith { 31092b32695SKris Buschelman /* This routine should really be abandoned as it duplicates MatDiagonalScaleLocal */ 31192b32695SKris Buschelman int ierr,(*f)(Mat,Vec); 31292b32695SKris Buschelman 31392b32695SKris Buschelman PetscFunctionBegin; 31492b32695SKris Buschelman ierr = PetscObjectQueryFunction((PetscObject)A,"MatDiagonalScaleLocal_C",(void (**)(void))&f);CHKERRQ(ierr); 31592b32695SKris Buschelman if (f) { 31692b32695SKris Buschelman ierr = (*f)(A,scale);CHKERRQ(ierr); 31792b32695SKris Buschelman } 31892b32695SKris Buschelman PetscFunctionReturn(0); 31992b32695SKris Buschelman } 32092b32695SKris Buschelman 32192b32695SKris Buschelman EXTERN_C_BEGIN 32292b32695SKris Buschelman #undef __FUNCT__ 32392b32695SKris Buschelman #define __FUNCT__ "MatDiagonalScaleLocal_MPIBAIJ" 32492b32695SKris Buschelman int MatDiagonalScaleLocal_MPIBAIJ(Mat A,Vec scale) 32592b32695SKris Buschelman { 32604f1ad80SBarry Smith Mat_MPIBAIJ *a = (Mat_MPIBAIJ*) A->data; /*access private part of matrix */ 32704f1ad80SBarry Smith int ierr,n,i; 32804f1ad80SBarry Smith PetscScalar *d,*o,*s; 32904f1ad80SBarry Smith 33004f1ad80SBarry Smith PetscFunctionBegin; 3312cd6534aSBarry Smith if (!uglyrmapd) { 3322cd6534aSBarry Smith ierr = MatMPIBAIJDiagonalScaleLocalSetUp(A,scale);CHKERRQ(ierr); 3332cd6534aSBarry Smith } 3342cd6534aSBarry Smith 335b1d4fb26SBarry Smith ierr = VecGetArrayFast(scale,&s);CHKERRQ(ierr); 33604f1ad80SBarry Smith 33704f1ad80SBarry Smith ierr = VecGetLocalSize(uglydd,&n);CHKERRQ(ierr); 338b1d4fb26SBarry Smith ierr = VecGetArrayFast(uglydd,&d);CHKERRQ(ierr); 33904f1ad80SBarry Smith for (i=0; i<n; i++) { 34004f1ad80SBarry Smith d[i] = s[uglyrmapd[i]]; /* copy "diagonal" (true local) portion of scale into dd vector */ 34104f1ad80SBarry Smith } 342b1d4fb26SBarry Smith ierr = VecRestoreArrayFast(uglydd,&d);CHKERRQ(ierr); 34304f1ad80SBarry Smith /* column scale "diagonal" portion of local matrix */ 34404f1ad80SBarry Smith ierr = MatDiagonalScale(a->A,PETSC_NULL,uglydd);CHKERRQ(ierr); 34504f1ad80SBarry Smith 34604f1ad80SBarry Smith ierr = VecGetLocalSize(uglyoo,&n);CHKERRQ(ierr); 347b1d4fb26SBarry Smith ierr = VecGetArrayFast(uglyoo,&o);CHKERRQ(ierr); 34804f1ad80SBarry Smith for (i=0; i<n; i++) { 34904f1ad80SBarry Smith o[i] = s[uglyrmapo[i]]; /* copy "off-diagonal" portion of scale into oo vector */ 35004f1ad80SBarry Smith } 351b1d4fb26SBarry Smith ierr = VecRestoreArrayFast(scale,&s);CHKERRQ(ierr); 352b1d4fb26SBarry Smith ierr = VecRestoreArrayFast(uglyoo,&o);CHKERRQ(ierr); 35304f1ad80SBarry Smith /* column scale "off-diagonal" portion of local matrix */ 35404f1ad80SBarry Smith ierr = MatDiagonalScale(a->B,PETSC_NULL,uglyoo);CHKERRQ(ierr); 35504f1ad80SBarry Smith 35604f1ad80SBarry Smith PetscFunctionReturn(0); 35704f1ad80SBarry Smith } 35892b32695SKris Buschelman EXTERN_C_END 3598016bdd1SSatish Balay 3603eda8832SBarry Smith 361