xref: /petsc/src/mat/impls/baij/mpi/mmbaij.c (revision 6857c123a4e6590fbaa83d0dbcd1ca4d08b61509)
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;
67*6857c123SSatish Balay   baij->B->cmap.n = baij->B->cmap.N = ec*mat->rmap.bs;
68*6857c123SSatish Balay   ierr = PetscMapInitialize(baij->B->comm,&(baij->B->cmap));CHKERRQ(ierr);
690f5bd95cSBarry Smith   ierr = PetscTableDelete(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;
104*6857c123SSatish Balay   baij->B->cmap.n =baij->B->cmap.N  = ec*mat->rmap.bs;
105*6857c123SSatish Balay   ierr = PetscMapInitialize(baij->B->comm,&(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 
1328016bdd1SSatish Balay   /* gnerate the scatter context */
133d9653453SSatish Balay   ierr = VecScatterCreate(gvec,from,baij->lvec,to,&baij->Mvctx);CHKERRQ(ierr);
13490f02eecSBarry Smith 
13590f02eecSBarry Smith   /*
13690f02eecSBarry Smith       Post the receives for the first matrix vector product. We sync-chronize after
13790f02eecSBarry Smith     this on the chance that the user immediately calls MatMult() after assemblying
13890f02eecSBarry Smith     the matrix.
13990f02eecSBarry Smith   */
14043a90d84SBarry Smith   ierr = VecScatterPostRecvs(gvec,baij->lvec,INSERT_VALUES,SCATTER_FORWARD,baij->Mvctx);CHKERRQ(ierr);
141ca161407SBarry Smith   ierr = MPI_Barrier(mat->comm);CHKERRQ(ierr);
14290f02eecSBarry Smith 
14352e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,baij->Mvctx);CHKERRQ(ierr);
14452e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,baij->lvec);CHKERRQ(ierr);
14552e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,from);CHKERRQ(ierr);
14652e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,to);CHKERRQ(ierr);
147d9653453SSatish Balay   baij->garray = garray;
14852e6d16bSBarry Smith   ierr = PetscLogObjectMemory(mat,(ec+1)*sizeof(PetscInt));CHKERRQ(ierr);
1498016bdd1SSatish Balay   ierr = ISDestroy(from);CHKERRQ(ierr);
1508016bdd1SSatish Balay   ierr = ISDestroy(to);CHKERRQ(ierr);
151888f2ed8SSatish Balay   ierr = VecDestroy(gvec);CHKERRQ(ierr);
1523a40ed3dSBarry Smith   PetscFunctionReturn(0);
1538016bdd1SSatish Balay }
1548016bdd1SSatish Balay 
1558016bdd1SSatish Balay /*
156d9653453SSatish Balay      Takes the local part of an already assembled MPIBAIJ matrix
1578016bdd1SSatish Balay    and disassembles it. This is to allow new nonzeros into the matrix
1588016bdd1SSatish Balay    that require more communication in the matrix vector multiply.
1598016bdd1SSatish Balay    Thus certain data-structures must be rebuilt.
1608016bdd1SSatish Balay 
1618016bdd1SSatish Balay    Kind of slow! But that's what application programmers get when
1628016bdd1SSatish Balay    they are sloppy.
1638016bdd1SSatish Balay */
1644a2ae208SSatish Balay #undef __FUNCT__
1654a2ae208SSatish Balay #define __FUNCT__ "DisAssemble_MPIBAIJ"
166dfbe8321SBarry Smith PetscErrorCode DisAssemble_MPIBAIJ(Mat A)
1678016bdd1SSatish Balay {
168d9653453SSatish Balay   Mat_MPIBAIJ    *baij = (Mat_MPIBAIJ*)A->data;
169d9653453SSatish Balay   Mat            B = baij->B,Bnew;
170d9653453SSatish Balay   Mat_SeqBAIJ    *Bbaij = (Mat_SeqBAIJ*)B->data;
171dfbe8321SBarry Smith   PetscErrorCode ierr;
172899cda47SBarry Smith   PetscInt       i,j,mbs=Bbaij->mbs,n = A->cmap.n,col,*garray=baij->garray;
173899cda47SBarry Smith   PetscInt       bs2 = baij->bs2,*nz,ec,m = A->rmap.n;
1743eda8832SBarry Smith   MatScalar      *a = Bbaij->a;
17587828ca2SBarry Smith   PetscScalar    *atmp;
17676117effSSatish Balay #if defined(PETSC_USE_MAT_SINGLE)
177b24ad042SBarry Smith   PetscInt       k;
17876117effSSatish Balay #endif
1798016bdd1SSatish Balay 
1803a40ed3dSBarry Smith   PetscFunctionBegin;
1818016bdd1SSatish Balay   /* free stuff related to matrix-vec multiply */
182b0a32e0cSBarry Smith   ierr = VecGetSize(baij->lvec,&ec);CHKERRQ(ierr); /* needed for PetscLogObjectMemory below */
183d9653453SSatish Balay   ierr = VecDestroy(baij->lvec);CHKERRQ(ierr); baij->lvec = 0;
184d9653453SSatish Balay   ierr = VecScatterDestroy(baij->Mvctx);CHKERRQ(ierr); baij->Mvctx = 0;
185d9653453SSatish Balay   if (baij->colmap) {
186aa482453SBarry Smith #if defined (PETSC_USE_CTABLE)
1870f5bd95cSBarry Smith     ierr = PetscTableDelete(baij->colmap); baij->colmap = 0;CHKERRQ(ierr);
18848e59246SSatish Balay #else
189606d414cSSatish Balay     ierr = PetscFree(baij->colmap);CHKERRQ(ierr);
190606d414cSSatish Balay     baij->colmap = 0;
19152e6d16bSBarry Smith     ierr = PetscLogObjectMemory(A,-Bbaij->nbs*sizeof(PetscInt));CHKERRQ(ierr);
19248e59246SSatish Balay #endif
1938016bdd1SSatish Balay   }
1948016bdd1SSatish Balay 
1958016bdd1SSatish Balay   /* make sure that B is assembled so we can access its values */
1966d4a8577SBarry Smith   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1973eda8832SBarry Smith   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1988016bdd1SSatish Balay 
1998016bdd1SSatish Balay   /* invent new B and copy stuff over */
200b24ad042SBarry Smith   ierr = PetscMalloc(mbs*sizeof(PetscInt),&nz);CHKERRQ(ierr);
201d9653453SSatish Balay   for (i=0; i<mbs; i++) {
202d9653453SSatish Balay     nz[i] = Bbaij->i[i+1]-Bbaij->i[i];
2038016bdd1SSatish Balay   }
204f69a0ea3SMatthew Knepley   ierr = MatCreate(B->comm,&Bnew);CHKERRQ(ierr);
205f69a0ea3SMatthew Knepley   ierr = MatSetSizes(Bnew,m,n,m,n);CHKERRQ(ierr);
206f204ca49SKris Buschelman   ierr = MatSetType(Bnew,B->type_name);CHKERRQ(ierr);
207899cda47SBarry Smith   ierr = MatSeqBAIJSetPreallocation(Bnew,B->rmap.bs,0,nz);CHKERRQ(ierr);
208bba1ac68SSatish Balay   ierr = MatSetOption(Bnew,MAT_COLUMN_ORIENTED);CHKERRQ(ierr);
209d9653453SSatish Balay 
2105010ffefSSatish Balay #if defined(PETSC_USE_MAT_SINGLE)
211bba1ac68SSatish Balay   ierr = PetscMalloc(bs2*sizeof(PetscScalar),&atmp);CHKERRQ(ierr);
212bba1ac68SSatish Balay #endif
213bba1ac68SSatish Balay     for (i=0; i<mbs; i++) {
214bba1ac68SSatish Balay       for (j=Bbaij->i[i]; j<Bbaij->i[i+1]; j++) {
215bba1ac68SSatish Balay         col  = garray[Bbaij->j[j]];
216bba1ac68SSatish Balay #if defined(PETSC_USE_MAT_SINGLE)
217bba1ac68SSatish Balay         for (k=0; k<bs2; k++) atmp[k] = a[j*bs2+k];
2183eda8832SBarry Smith #else
2193eda8832SBarry Smith         atmp = a + j*bs2;
2203eda8832SBarry Smith #endif
221bba1ac68SSatish Balay         ierr = MatSetValuesBlocked_SeqBAIJ(Bnew,1,&i,1,&col,atmp,B->insertmode);CHKERRQ(ierr);
2228016bdd1SSatish Balay       }
2238016bdd1SSatish Balay     }
224bba1ac68SSatish Balay   ierr = MatSetOption(Bnew,MAT_ROW_ORIENTED);CHKERRQ(ierr);
225bba1ac68SSatish Balay 
2263eda8832SBarry Smith #if defined(PETSC_USE_MAT_SINGLE)
2273eda8832SBarry Smith   ierr = PetscFree(atmp);CHKERRQ(ierr);
2283eda8832SBarry Smith #endif
229bba1ac68SSatish Balay 
230bba1ac68SSatish Balay   ierr = PetscFree(nz);CHKERRQ(ierr);
231606d414cSSatish Balay   ierr = PetscFree(baij->garray);CHKERRQ(ierr);
232606d414cSSatish Balay   baij->garray = 0;
23352e6d16bSBarry Smith   ierr = PetscLogObjectMemory(A,-ec*sizeof(PetscInt));CHKERRQ(ierr);
2348016bdd1SSatish Balay   ierr = MatDestroy(B);CHKERRQ(ierr);
23552e6d16bSBarry Smith   ierr = PetscLogObjectParent(A,Bnew);CHKERRQ(ierr);
236d9653453SSatish Balay   baij->B = Bnew;
2378016bdd1SSatish Balay   A->was_assembled = PETSC_FALSE;
2383a40ed3dSBarry Smith   PetscFunctionReturn(0);
2398016bdd1SSatish Balay }
2408016bdd1SSatish Balay 
24104f1ad80SBarry Smith /*      ugly stuff added for Glenn someday we should fix this up */
24204f1ad80SBarry Smith 
243b24ad042SBarry Smith static PetscInt *uglyrmapd = 0,*uglyrmapo = 0;  /* mapping from the local ordering to the "diagonal" and "off-diagonal"
24404f1ad80SBarry Smith                                       parts of the local matrix */
2452cd6534aSBarry Smith static Vec uglydd = 0,uglyoo = 0;   /* work vectors used to scale the two parts of the local matrix */
24604f1ad80SBarry Smith 
24704f1ad80SBarry Smith 
248348d1a9dSSatish Balay #undef __FUNCT__
249348d1a9dSSatish Balay #define __FUNCT__ "MatMPIBAIJDiagonalScaleLocalSetUp"
250dfbe8321SBarry Smith PetscErrorCode MatMPIBAIJDiagonalScaleLocalSetUp(Mat inA,Vec scale)
25104f1ad80SBarry Smith {
25204f1ad80SBarry Smith   Mat_MPIBAIJ    *ina = (Mat_MPIBAIJ*) inA->data; /*access private part of matrix */
25304f1ad80SBarry Smith   Mat_SeqBAIJ    *B = (Mat_SeqBAIJ*)ina->B->data;
254dfbe8321SBarry Smith   PetscErrorCode ierr;
255899cda47SBarry Smith   PetscInt       bs = inA->rmap.bs,i,n,nt,j,cstart,cend,no,*garray = ina->garray,*lindices;
256b24ad042SBarry Smith   PetscInt       *r_rmapd,*r_rmapo;
25704f1ad80SBarry Smith 
25804f1ad80SBarry Smith   PetscFunctionBegin;
25904f1ad80SBarry Smith   ierr = MatGetOwnershipRange(inA,&cstart,&cend);CHKERRQ(ierr);
26004f1ad80SBarry Smith   ierr = MatGetSize(ina->A,PETSC_NULL,&n);CHKERRQ(ierr);
261b24ad042SBarry Smith   ierr = PetscMalloc((inA->bmapping->n+1)*sizeof(PetscInt),&r_rmapd);CHKERRQ(ierr);
262b24ad042SBarry Smith   ierr = PetscMemzero(r_rmapd,inA->bmapping->n*sizeof(PetscInt));CHKERRQ(ierr);
26304f1ad80SBarry Smith   nt   = 0;
26404f1ad80SBarry Smith   for (i=0; i<inA->bmapping->n; i++) {
26504f1ad80SBarry Smith     if (inA->bmapping->indices[i]*bs >= cstart && inA->bmapping->indices[i]*bs < cend) {
26604f1ad80SBarry Smith       nt++;
26704f1ad80SBarry Smith       r_rmapd[i] = inA->bmapping->indices[i] + 1;
26804f1ad80SBarry Smith     }
26904f1ad80SBarry Smith   }
27077431f27SBarry Smith   if (nt*bs != n) SETERRQ2(PETSC_ERR_PLIB,"Hmm nt*bs %D n %D",nt*bs,n);
271b24ad042SBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&uglyrmapd);CHKERRQ(ierr);
27204f1ad80SBarry Smith   for (i=0; i<inA->bmapping->n; i++) {
27304f1ad80SBarry Smith     if (r_rmapd[i]){
27404f1ad80SBarry Smith       for (j=0; j<bs; j++) {
27504f1ad80SBarry Smith         uglyrmapd[(r_rmapd[i]-1)*bs+j-cstart] = i*bs + j;
27604f1ad80SBarry Smith       }
27704f1ad80SBarry Smith     }
27804f1ad80SBarry Smith   }
27904f1ad80SBarry Smith   ierr = PetscFree(r_rmapd);CHKERRQ(ierr);
28004f1ad80SBarry Smith   ierr = VecCreateSeq(PETSC_COMM_SELF,n,&uglydd);CHKERRQ(ierr);
28104f1ad80SBarry Smith 
282b24ad042SBarry Smith   ierr = PetscMalloc((ina->Nbs+1)*sizeof(PetscInt),&lindices);CHKERRQ(ierr);
283b24ad042SBarry Smith   ierr = PetscMemzero(lindices,ina->Nbs*sizeof(PetscInt));CHKERRQ(ierr);
28404f1ad80SBarry Smith   for (i=0; i<B->nbs; i++) {
28504f1ad80SBarry Smith     lindices[garray[i]] = i+1;
28604f1ad80SBarry Smith   }
28704f1ad80SBarry Smith   no   = inA->bmapping->n - nt;
288b24ad042SBarry Smith   ierr = PetscMalloc((inA->bmapping->n+1)*sizeof(PetscInt),&r_rmapo);CHKERRQ(ierr);
289b24ad042SBarry Smith   ierr = PetscMemzero(r_rmapo,inA->bmapping->n*sizeof(PetscInt));CHKERRQ(ierr);
29004f1ad80SBarry Smith   nt   = 0;
29104f1ad80SBarry Smith   for (i=0; i<inA->bmapping->n; i++) {
29204f1ad80SBarry Smith     if (lindices[inA->bmapping->indices[i]]) {
29304f1ad80SBarry Smith       nt++;
29404f1ad80SBarry Smith       r_rmapo[i] = lindices[inA->bmapping->indices[i]];
29504f1ad80SBarry Smith     }
29604f1ad80SBarry Smith   }
29777431f27SBarry Smith   if (nt > no) SETERRQ2(PETSC_ERR_PLIB,"Hmm nt %D no %D",nt,n);
29804f1ad80SBarry Smith   ierr = PetscFree(lindices);CHKERRQ(ierr);
299b24ad042SBarry Smith   ierr = PetscMalloc((nt*bs+1)*sizeof(PetscInt),&uglyrmapo);CHKERRQ(ierr);
30004f1ad80SBarry Smith   for (i=0; i<inA->bmapping->n; i++) {
30104f1ad80SBarry Smith     if (r_rmapo[i]){
30204f1ad80SBarry Smith       for (j=0; j<bs; j++) {
30304f1ad80SBarry Smith         uglyrmapo[(r_rmapo[i]-1)*bs+j] = i*bs + j;
30404f1ad80SBarry Smith       }
30504f1ad80SBarry Smith     }
30604f1ad80SBarry Smith   }
30704f1ad80SBarry Smith   ierr = PetscFree(r_rmapo);CHKERRQ(ierr);
30804f1ad80SBarry Smith   ierr = VecCreateSeq(PETSC_COMM_SELF,nt*bs,&uglyoo);CHKERRQ(ierr);
30904f1ad80SBarry Smith 
31004f1ad80SBarry Smith   PetscFunctionReturn(0);
31104f1ad80SBarry Smith }
31204f1ad80SBarry Smith 
313348d1a9dSSatish Balay #undef __FUNCT__
314348d1a9dSSatish Balay #define __FUNCT__ "MatMPIBAIJDiagonalScaleLocal"
315be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatMPIBAIJDiagonalScaleLocal(Mat A,Vec scale)
31604f1ad80SBarry Smith {
31792b32695SKris Buschelman   /* This routine should really be abandoned as it duplicates MatDiagonalScaleLocal */
318dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(Mat,Vec);
31992b32695SKris Buschelman 
32092b32695SKris Buschelman   PetscFunctionBegin;
32192b32695SKris Buschelman   ierr = PetscObjectQueryFunction((PetscObject)A,"MatDiagonalScaleLocal_C",(void (**)(void))&f);CHKERRQ(ierr);
32292b32695SKris Buschelman   if (f) {
32392b32695SKris Buschelman     ierr = (*f)(A,scale);CHKERRQ(ierr);
32492b32695SKris Buschelman   }
32592b32695SKris Buschelman   PetscFunctionReturn(0);
32692b32695SKris Buschelman }
32792b32695SKris Buschelman 
32892b32695SKris Buschelman EXTERN_C_BEGIN
32992b32695SKris Buschelman #undef __FUNCT__
33092b32695SKris Buschelman #define __FUNCT__ "MatDiagonalScaleLocal_MPIBAIJ"
331be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalScaleLocal_MPIBAIJ(Mat A,Vec scale)
33292b32695SKris Buschelman {
33304f1ad80SBarry Smith   Mat_MPIBAIJ    *a = (Mat_MPIBAIJ*) A->data; /*access private part of matrix */
334dfbe8321SBarry Smith   PetscErrorCode ierr;
335b24ad042SBarry Smith   PetscInt       n,i;
33604f1ad80SBarry Smith   PetscScalar    *d,*o,*s;
33704f1ad80SBarry Smith 
33804f1ad80SBarry Smith   PetscFunctionBegin;
3392cd6534aSBarry Smith   if (!uglyrmapd) {
3402cd6534aSBarry Smith     ierr = MatMPIBAIJDiagonalScaleLocalSetUp(A,scale);CHKERRQ(ierr);
3412cd6534aSBarry Smith   }
3422cd6534aSBarry Smith 
3431ebc52fbSHong Zhang   ierr = VecGetArray(scale,&s);CHKERRQ(ierr);
34404f1ad80SBarry Smith 
34504f1ad80SBarry Smith   ierr = VecGetLocalSize(uglydd,&n);CHKERRQ(ierr);
3461ebc52fbSHong Zhang   ierr = VecGetArray(uglydd,&d);CHKERRQ(ierr);
34704f1ad80SBarry Smith   for (i=0; i<n; i++) {
34804f1ad80SBarry Smith     d[i] = s[uglyrmapd[i]]; /* copy "diagonal" (true local) portion of scale into dd vector */
34904f1ad80SBarry Smith   }
3501ebc52fbSHong Zhang   ierr = VecRestoreArray(uglydd,&d);CHKERRQ(ierr);
35104f1ad80SBarry Smith   /* column scale "diagonal" portion of local matrix */
35204f1ad80SBarry Smith   ierr = MatDiagonalScale(a->A,PETSC_NULL,uglydd);CHKERRQ(ierr);
35304f1ad80SBarry Smith 
35404f1ad80SBarry Smith   ierr = VecGetLocalSize(uglyoo,&n);CHKERRQ(ierr);
3551ebc52fbSHong Zhang   ierr = VecGetArray(uglyoo,&o);CHKERRQ(ierr);
35604f1ad80SBarry Smith   for (i=0; i<n; i++) {
35704f1ad80SBarry Smith     o[i] = s[uglyrmapo[i]]; /* copy "off-diagonal" portion of scale into oo vector */
35804f1ad80SBarry Smith   }
3591ebc52fbSHong Zhang   ierr = VecRestoreArray(scale,&s);CHKERRQ(ierr);
3601ebc52fbSHong Zhang   ierr = VecRestoreArray(uglyoo,&o);CHKERRQ(ierr);
36104f1ad80SBarry Smith   /* column scale "off-diagonal" portion of local matrix */
36204f1ad80SBarry Smith   ierr = MatDiagonalScale(a->B,PETSC_NULL,uglyoo);CHKERRQ(ierr);
36304f1ad80SBarry Smith 
36404f1ad80SBarry Smith   PetscFunctionReturn(0);
36504f1ad80SBarry Smith }
36692b32695SKris Buschelman EXTERN_C_END
3678016bdd1SSatish Balay 
3683eda8832SBarry Smith 
369