xref: /petsc/src/mat/impls/baij/mpi/mmbaij.c (revision 26fbe8dc1a3c99fb8dddfa572c8c6b3b4ce3ca53)
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;
29aa482453SBarry Smith #if defined(PETSC_USE_CTABLE)
3073a2e727SSatish Balay   /* use a table - Mark Adams */
31e23dfa41SBarry Smith   ierr = PetscTableCreate(B->mbs,baij->Nbs+1,&gid1_lid1);CHKERRQ(ierr);
3273a2e727SSatish Balay   for (i=0; i<B->mbs; i++) {
3373a2e727SSatish Balay     for (j=0; j<B->ilen[i]; j++) {
34b24ad042SBarry Smith       PetscInt data,gid1 = aj[B->i[i]+j] + 1;
350f5bd95cSBarry Smith       ierr = PetscTableFind(gid1_lid1,gid1,&data);CHKERRQ(ierr);
36fa46199cSSatish Balay       if (!data) {
3773a2e727SSatish Balay         /* one based table */
383861aac3SJed Brown         ierr = PetscTableAdd(gid1_lid1,gid1,++ec,INSERT_VALUES);CHKERRQ(ierr);
3973a2e727SSatish Balay       }
4073a2e727SSatish Balay     }
4173a2e727SSatish Balay   }
4273a2e727SSatish Balay   /* form array of columns we need */
43b24ad042SBarry Smith   ierr = PetscMalloc((ec+1)*sizeof(PetscInt),&garray);CHKERRQ(ierr);
440f5bd95cSBarry Smith   ierr = PetscTableGetHeadPosition(gid1_lid1,&tpos);CHKERRQ(ierr);
4573a2e727SSatish Balay   while (tpos) {
460f5bd95cSBarry Smith     ierr = PetscTableGetNext(gid1_lid1,&tpos,&gid,&lid);CHKERRQ(ierr);
4773a2e727SSatish Balay     gid--; lid--;
4873a2e727SSatish Balay     garray[lid] = gid;
4973a2e727SSatish Balay   }
500064e2bbSSatish Balay   ierr = PetscSortInt(ec,garray);CHKERRQ(ierr);
510f5bd95cSBarry Smith   ierr = PetscTableRemoveAll(gid1_lid1);CHKERRQ(ierr);
5273a2e727SSatish Balay   for (i=0; i<ec; i++) {
533861aac3SJed Brown     ierr = PetscTableAdd(gid1_lid1,garray[i]+1,i+1,INSERT_VALUES);CHKERRQ(ierr);
5473a2e727SSatish Balay   }
5573a2e727SSatish Balay   /* compact out the extra columns in B */
5673a2e727SSatish Balay   for (i=0; i<B->mbs; i++) {
5773a2e727SSatish Balay     for (j=0; j<B->ilen[i]; j++) {
58b24ad042SBarry Smith       PetscInt gid1 = aj[B->i[i] + j] + 1;
590f5bd95cSBarry Smith       ierr = PetscTableFind(gid1_lid1,gid1,&lid);CHKERRQ(ierr);
60fa46199cSSatish Balay       lid--;
6173a2e727SSatish Balay       aj[B->i[i]+j] = lid;
6273a2e727SSatish Balay     }
6373a2e727SSatish Balay   }
6473a2e727SSatish Balay   B->nbs           = ec;
65d0f46423SBarry Smith   baij->B->cmap->n = baij->B->cmap->N = ec*mat->rmap->bs;
66*26fbe8dcSKarl Rupp 
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;
103*26fbe8dcSKarl Rupp 
10426283091SBarry Smith   ierr = PetscLayoutSetUp((baij->B->cmap));CHKERRQ(ierr);
105606d414cSSatish Balay   ierr = PetscFree(indices);CHKERRQ(ierr);
10673a2e727SSatish Balay #endif
1078016bdd1SSatish Balay 
1088016bdd1SSatish Balay   /* create local vector that is used to scatter into */
109029af93fSBarry Smith   ierr = VecCreateSeq(PETSC_COMM_SELF,ec*bs,&baij->lvec);CHKERRQ(ierr);
1108016bdd1SSatish Balay 
111c16cb8f2SBarry Smith   /* create two temporary index sets for building scatter-gather */
112deff0451SBarry Smith   ierr = ISCreateBlock(PETSC_COMM_SELF,bs,ec,garray,PETSC_COPY_VALUES,&from);CHKERRQ(ierr);
113c16cb8f2SBarry Smith 
114b24ad042SBarry Smith   ierr = PetscMalloc((ec+1)*sizeof(PetscInt),&stmp);CHKERRQ(ierr);
115*26fbe8dcSKarl Rupp   for (i=0; i<ec; i++) stmp[i] = i;
116deff0451SBarry Smith   ierr = ISCreateBlock(PETSC_COMM_SELF,bs,ec,stmp,PETSC_OWN_POINTER,&to);CHKERRQ(ierr);
1178016bdd1SSatish Balay 
1188016bdd1SSatish Balay   /* create temporary global vector to generate scatter context */
119778a2246SBarry Smith   ierr = VecCreateMPIWithArray(((PetscObject)mat)->comm,1,mat->cmap->n,mat->cmap->N,PETSC_NULL,&gvec);CHKERRQ(ierr);
1208016bdd1SSatish Balay 
121d9653453SSatish Balay   ierr = VecScatterCreate(gvec,from,baij->lvec,to,&baij->Mvctx);CHKERRQ(ierr);
12290f02eecSBarry Smith 
12352e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,baij->Mvctx);CHKERRQ(ierr);
12452e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,baij->lvec);CHKERRQ(ierr);
12552e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,from);CHKERRQ(ierr);
12652e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,to);CHKERRQ(ierr);
127*26fbe8dcSKarl Rupp 
128d9653453SSatish Balay   baij->garray = garray;
129*26fbe8dcSKarl Rupp 
13052e6d16bSBarry Smith   ierr = PetscLogObjectMemory(mat,(ec+1)*sizeof(PetscInt));CHKERRQ(ierr);
1316bf464f9SBarry Smith   ierr = ISDestroy(&from);CHKERRQ(ierr);
1326bf464f9SBarry Smith   ierr = ISDestroy(&to);CHKERRQ(ierr);
1336bf464f9SBarry Smith   ierr = VecDestroy(&gvec);CHKERRQ(ierr);
1343a40ed3dSBarry Smith   PetscFunctionReturn(0);
1358016bdd1SSatish Balay }
1368016bdd1SSatish Balay 
1378016bdd1SSatish Balay /*
138d9653453SSatish Balay      Takes the local part of an already assembled MPIBAIJ matrix
1398016bdd1SSatish Balay    and disassembles it. This is to allow new nonzeros into the matrix
1408016bdd1SSatish Balay    that require more communication in the matrix vector multiply.
1418016bdd1SSatish Balay    Thus certain data-structures must be rebuilt.
1428016bdd1SSatish Balay 
1438016bdd1SSatish Balay    Kind of slow! But that's what application programmers get when
1448016bdd1SSatish Balay    they are sloppy.
1458016bdd1SSatish Balay */
1464a2ae208SSatish Balay #undef __FUNCT__
147ab9863d7SBarry Smith #define __FUNCT__ "MatDisAssemble_MPIBAIJ"
148ab9863d7SBarry Smith PetscErrorCode MatDisAssemble_MPIBAIJ(Mat A)
1498016bdd1SSatish Balay {
150d9653453SSatish Balay   Mat_MPIBAIJ    *baij  = (Mat_MPIBAIJ*)A->data;
151d9653453SSatish Balay   Mat            B      = baij->B,Bnew;
152d9653453SSatish Balay   Mat_SeqBAIJ    *Bbaij = (Mat_SeqBAIJ*)B->data;
153dfbe8321SBarry Smith   PetscErrorCode ierr;
154d0f46423SBarry Smith   PetscInt       i,j,mbs=Bbaij->mbs,n = A->cmap->N,col,*garray=baij->garray;
155d0f46423SBarry Smith   PetscInt       bs2 = baij->bs2,*nz,ec,m = A->rmap->n;
1563eda8832SBarry Smith   MatScalar      *a  = Bbaij->a;
157dd6ea824SBarry Smith   MatScalar      *atmp;
15897e5c40aSBarry Smith 
1598016bdd1SSatish Balay 
1603a40ed3dSBarry Smith   PetscFunctionBegin;
1618016bdd1SSatish Balay   /* free stuff related to matrix-vec multiply */
162b0a32e0cSBarry Smith   ierr = VecGetSize(baij->lvec,&ec);CHKERRQ(ierr); /* needed for PetscLogObjectMemory below */
1636bf464f9SBarry Smith   ierr = VecDestroy(&baij->lvec);CHKERRQ(ierr); baij->lvec = 0;
1646bf464f9SBarry Smith   ierr = VecScatterDestroy(&baij->Mvctx);CHKERRQ(ierr); baij->Mvctx = 0;
165d9653453SSatish Balay   if (baij->colmap) {
166aa482453SBarry Smith #if defined(PETSC_USE_CTABLE)
1676bc0bbbfSBarry Smith     ierr = PetscTableDestroy(&baij->colmap);CHKERRQ(ierr);
16848e59246SSatish Balay #else
169606d414cSSatish Balay     ierr = PetscFree(baij->colmap);CHKERRQ(ierr);
17052e6d16bSBarry Smith     ierr = PetscLogObjectMemory(A,-Bbaij->nbs*sizeof(PetscInt));CHKERRQ(ierr);
17148e59246SSatish Balay #endif
1728016bdd1SSatish Balay   }
1738016bdd1SSatish Balay 
1748016bdd1SSatish Balay   /* make sure that B is assembled so we can access its values */
1756d4a8577SBarry Smith   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1763eda8832SBarry Smith   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1778016bdd1SSatish Balay 
1788016bdd1SSatish Balay   /* invent new B and copy stuff over */
179b24ad042SBarry Smith   ierr = PetscMalloc(mbs*sizeof(PetscInt),&nz);CHKERRQ(ierr);
180d9653453SSatish Balay   for (i=0; i<mbs; i++) {
181d9653453SSatish Balay     nz[i] = Bbaij->i[i+1]-Bbaij->i[i];
1828016bdd1SSatish Balay   }
1837adad957SLisandro Dalcin   ierr = MatCreate(((PetscObject)B)->comm,&Bnew);CHKERRQ(ierr);
184f69a0ea3SMatthew Knepley   ierr = MatSetSizes(Bnew,m,n,m,n);CHKERRQ(ierr);
1857adad957SLisandro Dalcin   ierr = MatSetType(Bnew,((PetscObject)B)->type_name);CHKERRQ(ierr);
186d0f46423SBarry Smith   ierr = MatSeqBAIJSetPreallocation(Bnew,B->rmap->bs,0,nz);CHKERRQ(ierr);
187*26fbe8dcSKarl Rupp 
1882576faa2SJed Brown   ((Mat_SeqBAIJ*)Bnew->data)->nonew = Bbaij->nonew; /* Inherit insertion error options. */
189*26fbe8dcSKarl Rupp 
1904e0d8c25SBarry Smith   ierr = MatSetOption(Bnew,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr);
191d9653453SSatish Balay 
192bba1ac68SSatish Balay   for (i=0; i<mbs; i++) {
193bba1ac68SSatish Balay     for (j=Bbaij->i[i]; j<Bbaij->i[i+1]; j++) {
194bba1ac68SSatish Balay       col  = garray[Bbaij->j[j]];
1953eda8832SBarry Smith       atmp = a + j*bs2;
196bba1ac68SSatish Balay       ierr = MatSetValuesBlocked_SeqBAIJ(Bnew,1,&i,1,&col,atmp,B->insertmode);CHKERRQ(ierr);
1978016bdd1SSatish Balay     }
1988016bdd1SSatish Balay   }
1994e0d8c25SBarry Smith   ierr = MatSetOption(Bnew,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr);
200bba1ac68SSatish Balay 
201bba1ac68SSatish Balay   ierr = PetscFree(nz);CHKERRQ(ierr);
202606d414cSSatish Balay   ierr = PetscFree(baij->garray);CHKERRQ(ierr);
20352e6d16bSBarry Smith   ierr = PetscLogObjectMemory(A,-ec*sizeof(PetscInt));CHKERRQ(ierr);
2046bf464f9SBarry Smith   ierr = MatDestroy(&B);CHKERRQ(ierr);
20552e6d16bSBarry Smith   ierr = PetscLogObjectParent(A,Bnew);CHKERRQ(ierr);
206*26fbe8dcSKarl Rupp 
207d9653453SSatish Balay   baij->B          = Bnew;
2088016bdd1SSatish Balay   A->was_assembled = PETSC_FALSE;
2096a719282SBarry Smith   A->assembled     = PETSC_FALSE;
2103a40ed3dSBarry Smith   PetscFunctionReturn(0);
2118016bdd1SSatish Balay }
2128016bdd1SSatish Balay 
21304f1ad80SBarry Smith /*      ugly stuff added for Glenn someday we should fix this up */
21404f1ad80SBarry Smith 
215*26fbe8dcSKarl Rupp static PetscInt *uglyrmapd = 0,*uglyrmapo = 0;  /* mapping from the local ordering to the "diagonal" and "off-diagonal" parts of the local matrix */
2162cd6534aSBarry Smith static Vec      uglydd     = 0,uglyoo     = 0;  /* work vectors used to scale the two parts of the local matrix */
21704f1ad80SBarry Smith 
21804f1ad80SBarry Smith 
219348d1a9dSSatish Balay #undef __FUNCT__
220348d1a9dSSatish Balay #define __FUNCT__ "MatMPIBAIJDiagonalScaleLocalSetUp"
221dfbe8321SBarry Smith PetscErrorCode MatMPIBAIJDiagonalScaleLocalSetUp(Mat inA,Vec scale)
22204f1ad80SBarry Smith {
22304f1ad80SBarry Smith   Mat_MPIBAIJ    *ina = (Mat_MPIBAIJ*) inA->data; /*access private part of matrix */
22404f1ad80SBarry Smith   Mat_SeqBAIJ    *B   = (Mat_SeqBAIJ*)ina->B->data;
225dfbe8321SBarry Smith   PetscErrorCode ierr;
226d0f46423SBarry Smith   PetscInt       bs = inA->rmap->bs,i,n,nt,j,cstart,cend,no,*garray = ina->garray,*lindices;
227b24ad042SBarry Smith   PetscInt       *r_rmapd,*r_rmapo;
22804f1ad80SBarry Smith 
22904f1ad80SBarry Smith   PetscFunctionBegin;
23004f1ad80SBarry Smith   ierr = MatGetOwnershipRange(inA,&cstart,&cend);CHKERRQ(ierr);
23104f1ad80SBarry Smith   ierr = MatGetSize(ina->A,PETSC_NULL,&n);CHKERRQ(ierr);
232992144d0SBarry Smith   ierr = PetscMalloc((inA->rmap->bmapping->n+1)*sizeof(PetscInt),&r_rmapd);CHKERRQ(ierr);
233992144d0SBarry Smith   ierr = PetscMemzero(r_rmapd,inA->rmap->bmapping->n*sizeof(PetscInt));CHKERRQ(ierr);
23404f1ad80SBarry Smith   nt   = 0;
235992144d0SBarry Smith   for (i=0; i<inA->rmap->bmapping->n; i++) {
236992144d0SBarry Smith     if (inA->rmap->bmapping->indices[i]*bs >= cstart && inA->rmap->bmapping->indices[i]*bs < cend) {
23704f1ad80SBarry Smith       nt++;
238992144d0SBarry Smith       r_rmapd[i] = inA->rmap->bmapping->indices[i] + 1;
23904f1ad80SBarry Smith     }
24004f1ad80SBarry Smith   }
241e32f2f54SBarry Smith   if (nt*bs != n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hmm nt*bs %D n %D",nt*bs,n);
242b24ad042SBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&uglyrmapd);CHKERRQ(ierr);
243992144d0SBarry Smith   for (i=0; i<inA->rmap->bmapping->n; i++) {
24404f1ad80SBarry Smith     if (r_rmapd[i]) {
24504f1ad80SBarry Smith       for (j=0; j<bs; j++) {
24604f1ad80SBarry Smith         uglyrmapd[(r_rmapd[i]-1)*bs+j-cstart] = i*bs + j;
24704f1ad80SBarry Smith       }
24804f1ad80SBarry Smith     }
24904f1ad80SBarry Smith   }
25004f1ad80SBarry Smith   ierr = PetscFree(r_rmapd);CHKERRQ(ierr);
25104f1ad80SBarry Smith   ierr = VecCreateSeq(PETSC_COMM_SELF,n,&uglydd);CHKERRQ(ierr);
25204f1ad80SBarry Smith 
253b24ad042SBarry Smith   ierr = PetscMalloc((ina->Nbs+1)*sizeof(PetscInt),&lindices);CHKERRQ(ierr);
254b24ad042SBarry Smith   ierr = PetscMemzero(lindices,ina->Nbs*sizeof(PetscInt));CHKERRQ(ierr);
25504f1ad80SBarry Smith   for (i=0; i<B->nbs; i++) {
25604f1ad80SBarry Smith     lindices[garray[i]] = i+1;
25704f1ad80SBarry Smith   }
258992144d0SBarry Smith   no   = inA->rmap->bmapping->n - nt;
259992144d0SBarry Smith   ierr = PetscMalloc((inA->rmap->bmapping->n+1)*sizeof(PetscInt),&r_rmapo);CHKERRQ(ierr);
260992144d0SBarry Smith   ierr = PetscMemzero(r_rmapo,inA->rmap->bmapping->n*sizeof(PetscInt));CHKERRQ(ierr);
26104f1ad80SBarry Smith   nt   = 0;
262992144d0SBarry Smith   for (i=0; i<inA->rmap->bmapping->n; i++) {
263992144d0SBarry Smith     if (lindices[inA->rmap->bmapping->indices[i]]) {
26404f1ad80SBarry Smith       nt++;
265992144d0SBarry Smith       r_rmapo[i] = lindices[inA->rmap->bmapping->indices[i]];
26604f1ad80SBarry Smith     }
26704f1ad80SBarry Smith   }
268e32f2f54SBarry Smith   if (nt > no) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hmm nt %D no %D",nt,n);
26904f1ad80SBarry Smith   ierr = PetscFree(lindices);CHKERRQ(ierr);
270b24ad042SBarry Smith   ierr = PetscMalloc((nt*bs+1)*sizeof(PetscInt),&uglyrmapo);CHKERRQ(ierr);
271992144d0SBarry Smith   for (i=0; i<inA->rmap->bmapping->n; i++) {
27204f1ad80SBarry Smith     if (r_rmapo[i]) {
27304f1ad80SBarry Smith       for (j=0; j<bs; j++) {
27404f1ad80SBarry Smith         uglyrmapo[(r_rmapo[i]-1)*bs+j] = i*bs + j;
27504f1ad80SBarry Smith       }
27604f1ad80SBarry Smith     }
27704f1ad80SBarry Smith   }
27804f1ad80SBarry Smith   ierr = PetscFree(r_rmapo);CHKERRQ(ierr);
27904f1ad80SBarry Smith   ierr = VecCreateSeq(PETSC_COMM_SELF,nt*bs,&uglyoo);CHKERRQ(ierr);
28004f1ad80SBarry Smith   PetscFunctionReturn(0);
28104f1ad80SBarry Smith }
28204f1ad80SBarry Smith 
283348d1a9dSSatish Balay #undef __FUNCT__
284348d1a9dSSatish Balay #define __FUNCT__ "MatMPIBAIJDiagonalScaleLocal"
2857087cfbeSBarry Smith PetscErrorCode  MatMPIBAIJDiagonalScaleLocal(Mat A,Vec scale)
28604f1ad80SBarry Smith {
28792b32695SKris Buschelman   /* This routine should really be abandoned as it duplicates MatDiagonalScaleLocal */
2884ac538c5SBarry Smith   PetscErrorCode ierr;
28992b32695SKris Buschelman 
29092b32695SKris Buschelman   PetscFunctionBegin;
2914ac538c5SBarry Smith   ierr = PetscTryMethod(A,"MatDiagonalScaleLocal_C",(Mat,Vec),(A,scale));CHKERRQ(ierr);
29292b32695SKris Buschelman   PetscFunctionReturn(0);
29392b32695SKris Buschelman }
29492b32695SKris Buschelman 
29592b32695SKris Buschelman EXTERN_C_BEGIN
29692b32695SKris Buschelman #undef __FUNCT__
29792b32695SKris Buschelman #define __FUNCT__ "MatDiagonalScaleLocal_MPIBAIJ"
2987087cfbeSBarry Smith PetscErrorCode  MatDiagonalScaleLocal_MPIBAIJ(Mat A,Vec scale)
29992b32695SKris Buschelman {
30004f1ad80SBarry Smith   Mat_MPIBAIJ    *a = (Mat_MPIBAIJ*) A->data; /*access private part of matrix */
301dfbe8321SBarry Smith   PetscErrorCode ierr;
302b24ad042SBarry Smith   PetscInt       n,i;
30304f1ad80SBarry Smith   PetscScalar    *d,*o,*s;
30404f1ad80SBarry Smith 
30504f1ad80SBarry Smith   PetscFunctionBegin;
3062cd6534aSBarry Smith   if (!uglyrmapd) {
3072cd6534aSBarry Smith     ierr = MatMPIBAIJDiagonalScaleLocalSetUp(A,scale);CHKERRQ(ierr);
3082cd6534aSBarry Smith   }
3092cd6534aSBarry Smith 
3101ebc52fbSHong Zhang   ierr = VecGetArray(scale,&s);CHKERRQ(ierr);
31104f1ad80SBarry Smith 
31204f1ad80SBarry Smith   ierr = VecGetLocalSize(uglydd,&n);CHKERRQ(ierr);
3131ebc52fbSHong Zhang   ierr = VecGetArray(uglydd,&d);CHKERRQ(ierr);
31404f1ad80SBarry Smith   for (i=0; i<n; i++) {
31504f1ad80SBarry Smith     d[i] = s[uglyrmapd[i]]; /* copy "diagonal" (true local) portion of scale into dd vector */
31604f1ad80SBarry Smith   }
3171ebc52fbSHong Zhang   ierr = VecRestoreArray(uglydd,&d);CHKERRQ(ierr);
31804f1ad80SBarry Smith   /* column scale "diagonal" portion of local matrix */
31904f1ad80SBarry Smith   ierr = MatDiagonalScale(a->A,PETSC_NULL,uglydd);CHKERRQ(ierr);
32004f1ad80SBarry Smith 
32104f1ad80SBarry Smith   ierr = VecGetLocalSize(uglyoo,&n);CHKERRQ(ierr);
3221ebc52fbSHong Zhang   ierr = VecGetArray(uglyoo,&o);CHKERRQ(ierr);
32304f1ad80SBarry Smith   for (i=0; i<n; i++) {
32404f1ad80SBarry Smith     o[i] = s[uglyrmapo[i]]; /* copy "off-diagonal" portion of scale into oo vector */
32504f1ad80SBarry Smith   }
3261ebc52fbSHong Zhang   ierr = VecRestoreArray(scale,&s);CHKERRQ(ierr);
3271ebc52fbSHong Zhang   ierr = VecRestoreArray(uglyoo,&o);CHKERRQ(ierr);
32804f1ad80SBarry Smith   /* column scale "off-diagonal" portion of local matrix */
32904f1ad80SBarry Smith   ierr = MatDiagonalScale(a->B,PETSC_NULL,uglyoo);CHKERRQ(ierr);
33004f1ad80SBarry Smith   PetscFunctionReturn(0);
33104f1ad80SBarry Smith }
33292b32695SKris Buschelman EXTERN_C_END
3338016bdd1SSatish Balay 
3343eda8832SBarry Smith 
335