xref: /petsc/src/mat/impls/baij/mpi/mmbaij.c (revision 3bb1ff401821b9e2ae019d3e61bc8ab4bd4e59d5)
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>
6186d4ecdSBarry Smith #include <petsc-private/isimpl.h>    /* needed because accesses data structure of ISLocalToGlobalMapping directly */
7bba1ac68SSatish Balay 
809573ac7SBarry 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;
30aa482453SBarry Smith #if defined(PETSC_USE_CTABLE)
3173a2e727SSatish Balay   /* use a table - Mark Adams */
32e23dfa41SBarry Smith   ierr = PetscTableCreate(B->mbs,baij->Nbs+1,&gid1_lid1);CHKERRQ(ierr);
3373a2e727SSatish Balay   for (i=0; i<B->mbs; i++) {
3473a2e727SSatish Balay     for (j=0; j<B->ilen[i]; j++) {
35b24ad042SBarry Smith       PetscInt data,gid1 = aj[B->i[i]+j] + 1;
360f5bd95cSBarry Smith       ierr = PetscTableFind(gid1_lid1,gid1,&data);CHKERRQ(ierr);
37fa46199cSSatish Balay       if (!data) {
3873a2e727SSatish Balay         /* one based table */
393861aac3SJed Brown         ierr = PetscTableAdd(gid1_lid1,gid1,++ec,INSERT_VALUES);CHKERRQ(ierr);
4073a2e727SSatish Balay       }
4173a2e727SSatish Balay     }
4273a2e727SSatish Balay   }
4373a2e727SSatish Balay   /* form array of columns we need */
44b24ad042SBarry Smith   ierr = PetscMalloc((ec+1)*sizeof(PetscInt),&garray);CHKERRQ(ierr);
450f5bd95cSBarry Smith   ierr = PetscTableGetHeadPosition(gid1_lid1,&tpos);CHKERRQ(ierr);
4673a2e727SSatish Balay   while (tpos) {
470f5bd95cSBarry Smith     ierr = PetscTableGetNext(gid1_lid1,&tpos,&gid,&lid);CHKERRQ(ierr);
4873a2e727SSatish Balay     gid--; lid--;
4973a2e727SSatish Balay     garray[lid] = gid;
5073a2e727SSatish Balay   }
510064e2bbSSatish Balay   ierr = PetscSortInt(ec,garray);CHKERRQ(ierr);
520f5bd95cSBarry Smith   ierr = PetscTableRemoveAll(gid1_lid1);CHKERRQ(ierr);
5373a2e727SSatish Balay   for (i=0; i<ec; i++) {
543861aac3SJed Brown     ierr = PetscTableAdd(gid1_lid1,garray[i]+1,i+1,INSERT_VALUES);CHKERRQ(ierr);
5573a2e727SSatish Balay   }
5673a2e727SSatish Balay   /* compact out the extra columns in B */
5773a2e727SSatish Balay   for (i=0; i<B->mbs; i++) {
5873a2e727SSatish Balay     for (j=0; j<B->ilen[i]; j++) {
59b24ad042SBarry Smith       PetscInt gid1 = aj[B->i[i] + j] + 1;
600f5bd95cSBarry Smith       ierr = PetscTableFind(gid1_lid1,gid1,&lid);CHKERRQ(ierr);
61fa46199cSSatish Balay       lid--;
6273a2e727SSatish Balay       aj[B->i[i]+j] = lid;
6373a2e727SSatish Balay     }
6473a2e727SSatish Balay   }
6573a2e727SSatish Balay   B->nbs           = ec;
66d0f46423SBarry Smith   baij->B->cmap->n = baij->B->cmap->N = ec*mat->rmap->bs;
6726fbe8dcSKarl Rupp 
6826283091SBarry Smith   ierr = PetscLayoutSetUp((baij->B->cmap));CHKERRQ(ierr);
696bc0bbbfSBarry Smith   ierr = PetscTableDestroy(&gid1_lid1);CHKERRQ(ierr);
7073a2e727SSatish Balay #else
71435da068SBarry Smith   /* Make an array as long as the number of columns */
72d9653453SSatish Balay   /* mark those columns that are in baij->B */
73b24ad042SBarry Smith   ierr = PetscMalloc((Nbs+1)*sizeof(PetscInt),&indices);CHKERRQ(ierr);
74b24ad042SBarry Smith   ierr = PetscMemzero(indices,Nbs*sizeof(PetscInt));CHKERRQ(ierr);
75d9653453SSatish Balay   for (i=0; i<B->mbs; i++) {
768016bdd1SSatish Balay     for (j=0; j<B->ilen[i]; j++) {
77d9653453SSatish Balay       if (!indices[aj[B->i[i] + j]]) ec++;
78d9653453SSatish Balay       indices[aj[B->i[i] + j]] = 1;
798016bdd1SSatish Balay     }
808016bdd1SSatish Balay   }
818016bdd1SSatish Balay 
828016bdd1SSatish Balay   /* form array of columns we need */
83b24ad042SBarry Smith   ierr = PetscMalloc((ec+1)*sizeof(PetscInt),&garray);CHKERRQ(ierr);
848016bdd1SSatish Balay   ec   = 0;
85d9653453SSatish Balay   for (i=0; i<Nbs; i++) {
860bdbc534SSatish Balay     if (indices[i]) {
870bdbc534SSatish Balay       garray[ec++] = i;
880bdbc534SSatish Balay     }
898016bdd1SSatish Balay   }
908016bdd1SSatish Balay 
918016bdd1SSatish Balay   /* make indices now point into garray */
928016bdd1SSatish Balay   for (i=0; i<ec; i++) {
93d9653453SSatish Balay     indices[garray[i]] = i;
948016bdd1SSatish Balay   }
958016bdd1SSatish Balay 
968016bdd1SSatish Balay   /* compact out the extra columns in B */
97d9653453SSatish Balay   for (i=0; i<B->mbs; i++) {
988016bdd1SSatish Balay     for (j=0; j<B->ilen[i]; j++) {
99d9653453SSatish Balay       aj[B->i[i] + j] = indices[aj[B->i[i] + j]];
1008016bdd1SSatish Balay     }
1018016bdd1SSatish Balay   }
102d9653453SSatish Balay   B->nbs           = ec;
103d0f46423SBarry Smith   baij->B->cmap->n = baij->B->cmap->N  = ec*mat->rmap->bs;
10426fbe8dcSKarl Rupp 
10526283091SBarry Smith   ierr = PetscLayoutSetUp((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 */
113deff0451SBarry Smith   ierr = ISCreateBlock(PETSC_COMM_SELF,bs,ec,garray,PETSC_COPY_VALUES,&from);CHKERRQ(ierr);
114c16cb8f2SBarry Smith 
115b24ad042SBarry Smith   ierr = PetscMalloc((ec+1)*sizeof(PetscInt),&stmp);CHKERRQ(ierr);
11626fbe8dcSKarl Rupp   for (i=0; i<ec; i++) stmp[i] = i;
117deff0451SBarry Smith   ierr = ISCreateBlock(PETSC_COMM_SELF,bs,ec,stmp,PETSC_OWN_POINTER,&to);CHKERRQ(ierr);
1188016bdd1SSatish Balay 
1198016bdd1SSatish Balay   /* create temporary global vector to generate scatter context */
120ce94432eSBarry Smith   ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)mat),1,mat->cmap->n,mat->cmap->N,NULL,&gvec);CHKERRQ(ierr);
1218016bdd1SSatish Balay 
122d9653453SSatish Balay   ierr = VecScatterCreate(gvec,from,baij->lvec,to,&baij->Mvctx);CHKERRQ(ierr);
12390f02eecSBarry Smith 
124*3bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)baij->Mvctx);CHKERRQ(ierr);
125*3bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)baij->lvec);CHKERRQ(ierr);
126*3bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)from);CHKERRQ(ierr);
127*3bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)to);CHKERRQ(ierr);
12826fbe8dcSKarl Rupp 
129d9653453SSatish Balay   baij->garray = garray;
13026fbe8dcSKarl Rupp 
131*3bb1ff40SBarry Smith   ierr = PetscLogObjectMemory((PetscObject)mat,(ec+1)*sizeof(PetscInt));CHKERRQ(ierr);
1326bf464f9SBarry Smith   ierr = ISDestroy(&from);CHKERRQ(ierr);
1336bf464f9SBarry Smith   ierr = ISDestroy(&to);CHKERRQ(ierr);
1346bf464f9SBarry Smith   ierr = VecDestroy(&gvec);CHKERRQ(ierr);
1353a40ed3dSBarry Smith   PetscFunctionReturn(0);
1368016bdd1SSatish Balay }
1378016bdd1SSatish Balay 
1388016bdd1SSatish Balay /*
139d9653453SSatish Balay      Takes the local part of an already assembled MPIBAIJ matrix
1408016bdd1SSatish Balay    and disassembles it. This is to allow new nonzeros into the matrix
1418016bdd1SSatish Balay    that require more communication in the matrix vector multiply.
1428016bdd1SSatish Balay    Thus certain data-structures must be rebuilt.
1438016bdd1SSatish Balay 
1448016bdd1SSatish Balay    Kind of slow! But that's what application programmers get when
1458016bdd1SSatish Balay    they are sloppy.
1468016bdd1SSatish Balay */
1474a2ae208SSatish Balay #undef __FUNCT__
148ab9863d7SBarry Smith #define __FUNCT__ "MatDisAssemble_MPIBAIJ"
149ab9863d7SBarry Smith PetscErrorCode MatDisAssemble_MPIBAIJ(Mat A)
1508016bdd1SSatish Balay {
151d9653453SSatish Balay   Mat_MPIBAIJ    *baij  = (Mat_MPIBAIJ*)A->data;
152d9653453SSatish Balay   Mat            B      = baij->B,Bnew;
153d9653453SSatish Balay   Mat_SeqBAIJ    *Bbaij = (Mat_SeqBAIJ*)B->data;
154dfbe8321SBarry Smith   PetscErrorCode ierr;
155d0f46423SBarry Smith   PetscInt       i,j,mbs=Bbaij->mbs,n = A->cmap->N,col,*garray=baij->garray;
156d0f46423SBarry Smith   PetscInt       bs2 = baij->bs2,*nz,ec,m = A->rmap->n;
1573eda8832SBarry Smith   MatScalar      *a  = Bbaij->a;
158dd6ea824SBarry Smith   MatScalar      *atmp;
15997e5c40aSBarry Smith 
1608016bdd1SSatish Balay 
1613a40ed3dSBarry Smith   PetscFunctionBegin;
1628016bdd1SSatish Balay   /* free stuff related to matrix-vec multiply */
163b0a32e0cSBarry Smith   ierr = VecGetSize(baij->lvec,&ec);CHKERRQ(ierr); /* needed for PetscLogObjectMemory below */
1646bf464f9SBarry Smith   ierr = VecDestroy(&baij->lvec);CHKERRQ(ierr); baij->lvec = 0;
1656bf464f9SBarry Smith   ierr = VecScatterDestroy(&baij->Mvctx);CHKERRQ(ierr); baij->Mvctx = 0;
166d9653453SSatish Balay   if (baij->colmap) {
167aa482453SBarry Smith #if defined(PETSC_USE_CTABLE)
1686bc0bbbfSBarry Smith     ierr = PetscTableDestroy(&baij->colmap);CHKERRQ(ierr);
16948e59246SSatish Balay #else
170606d414cSSatish Balay     ierr = PetscFree(baij->colmap);CHKERRQ(ierr);
171*3bb1ff40SBarry Smith     ierr = PetscLogObjectMemory((PetscObject)A,-Bbaij->nbs*sizeof(PetscInt));CHKERRQ(ierr);
17248e59246SSatish Balay #endif
1738016bdd1SSatish Balay   }
1748016bdd1SSatish Balay 
1758016bdd1SSatish Balay   /* make sure that B is assembled so we can access its values */
1766d4a8577SBarry Smith   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1773eda8832SBarry Smith   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1788016bdd1SSatish Balay 
1798016bdd1SSatish Balay   /* invent new B and copy stuff over */
180b24ad042SBarry Smith   ierr = PetscMalloc(mbs*sizeof(PetscInt),&nz);CHKERRQ(ierr);
181d9653453SSatish Balay   for (i=0; i<mbs; i++) {
182d9653453SSatish Balay     nz[i] = Bbaij->i[i+1]-Bbaij->i[i];
1838016bdd1SSatish Balay   }
184ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)B),&Bnew);CHKERRQ(ierr);
185f69a0ea3SMatthew Knepley   ierr = MatSetSizes(Bnew,m,n,m,n);CHKERRQ(ierr);
1867adad957SLisandro Dalcin   ierr = MatSetType(Bnew,((PetscObject)B)->type_name);CHKERRQ(ierr);
187d0f46423SBarry Smith   ierr = MatSeqBAIJSetPreallocation(Bnew,B->rmap->bs,0,nz);CHKERRQ(ierr);
18826fbe8dcSKarl Rupp 
1892576faa2SJed Brown   ((Mat_SeqBAIJ*)Bnew->data)->nonew = Bbaij->nonew; /* Inherit insertion error options. */
19026fbe8dcSKarl Rupp 
1914e0d8c25SBarry Smith   ierr = MatSetOption(Bnew,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr);
192d9653453SSatish Balay 
193bba1ac68SSatish Balay   for (i=0; i<mbs; i++) {
194bba1ac68SSatish Balay     for (j=Bbaij->i[i]; j<Bbaij->i[i+1]; j++) {
195bba1ac68SSatish Balay       col  = garray[Bbaij->j[j]];
1963eda8832SBarry Smith       atmp = a + j*bs2;
197bba1ac68SSatish Balay       ierr = MatSetValuesBlocked_SeqBAIJ(Bnew,1,&i,1,&col,atmp,B->insertmode);CHKERRQ(ierr);
1988016bdd1SSatish Balay     }
1998016bdd1SSatish Balay   }
2004e0d8c25SBarry Smith   ierr = MatSetOption(Bnew,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr);
201bba1ac68SSatish Balay 
202bba1ac68SSatish Balay   ierr = PetscFree(nz);CHKERRQ(ierr);
203606d414cSSatish Balay   ierr = PetscFree(baij->garray);CHKERRQ(ierr);
204*3bb1ff40SBarry Smith   ierr = PetscLogObjectMemory((PetscObject)A,-ec*sizeof(PetscInt));CHKERRQ(ierr);
2056bf464f9SBarry Smith   ierr = MatDestroy(&B);CHKERRQ(ierr);
206*3bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)A,(PetscObject)Bnew);CHKERRQ(ierr);
20726fbe8dcSKarl Rupp 
208d9653453SSatish Balay   baij->B          = Bnew;
2098016bdd1SSatish Balay   A->was_assembled = PETSC_FALSE;
2106a719282SBarry Smith   A->assembled     = PETSC_FALSE;
2113a40ed3dSBarry Smith   PetscFunctionReturn(0);
2128016bdd1SSatish Balay }
2138016bdd1SSatish Balay 
21404f1ad80SBarry Smith /*      ugly stuff added for Glenn someday we should fix this up */
21504f1ad80SBarry Smith 
21626fbe8dcSKarl Rupp static PetscInt *uglyrmapd = 0,*uglyrmapo = 0;  /* mapping from the local ordering to the "diagonal" and "off-diagonal" parts of the local matrix */
2172cd6534aSBarry Smith static Vec      uglydd     = 0,uglyoo     = 0;  /* work vectors used to scale the two parts of the local matrix */
21804f1ad80SBarry Smith 
21904f1ad80SBarry Smith 
220348d1a9dSSatish Balay #undef __FUNCT__
221348d1a9dSSatish Balay #define __FUNCT__ "MatMPIBAIJDiagonalScaleLocalSetUp"
222dfbe8321SBarry Smith PetscErrorCode MatMPIBAIJDiagonalScaleLocalSetUp(Mat inA,Vec scale)
22304f1ad80SBarry Smith {
22404f1ad80SBarry Smith   Mat_MPIBAIJ    *ina = (Mat_MPIBAIJ*) inA->data; /*access private part of matrix */
22504f1ad80SBarry Smith   Mat_SeqBAIJ    *B   = (Mat_SeqBAIJ*)ina->B->data;
226dfbe8321SBarry Smith   PetscErrorCode ierr;
227d0f46423SBarry Smith   PetscInt       bs = inA->rmap->bs,i,n,nt,j,cstart,cend,no,*garray = ina->garray,*lindices;
228b24ad042SBarry Smith   PetscInt       *r_rmapd,*r_rmapo;
22904f1ad80SBarry Smith 
23004f1ad80SBarry Smith   PetscFunctionBegin;
23104f1ad80SBarry Smith   ierr = MatGetOwnershipRange(inA,&cstart,&cend);CHKERRQ(ierr);
2320298fd71SBarry Smith   ierr = MatGetSize(ina->A,NULL,&n);CHKERRQ(ierr);
233992144d0SBarry Smith   ierr = PetscMalloc((inA->rmap->bmapping->n+1)*sizeof(PetscInt),&r_rmapd);CHKERRQ(ierr);
234992144d0SBarry Smith   ierr = PetscMemzero(r_rmapd,inA->rmap->bmapping->n*sizeof(PetscInt));CHKERRQ(ierr);
23504f1ad80SBarry Smith   nt   = 0;
236992144d0SBarry Smith   for (i=0; i<inA->rmap->bmapping->n; i++) {
237992144d0SBarry Smith     if (inA->rmap->bmapping->indices[i]*bs >= cstart && inA->rmap->bmapping->indices[i]*bs < cend) {
23804f1ad80SBarry Smith       nt++;
239992144d0SBarry Smith       r_rmapd[i] = inA->rmap->bmapping->indices[i] + 1;
24004f1ad80SBarry Smith     }
24104f1ad80SBarry Smith   }
242e32f2f54SBarry Smith   if (nt*bs != n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hmm nt*bs %D n %D",nt*bs,n);
243b24ad042SBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&uglyrmapd);CHKERRQ(ierr);
244992144d0SBarry Smith   for (i=0; i<inA->rmap->bmapping->n; i++) {
24504f1ad80SBarry Smith     if (r_rmapd[i]) {
24604f1ad80SBarry Smith       for (j=0; j<bs; j++) {
24704f1ad80SBarry Smith         uglyrmapd[(r_rmapd[i]-1)*bs+j-cstart] = i*bs + j;
24804f1ad80SBarry Smith       }
24904f1ad80SBarry Smith     }
25004f1ad80SBarry Smith   }
25104f1ad80SBarry Smith   ierr = PetscFree(r_rmapd);CHKERRQ(ierr);
25204f1ad80SBarry Smith   ierr = VecCreateSeq(PETSC_COMM_SELF,n,&uglydd);CHKERRQ(ierr);
25304f1ad80SBarry Smith 
254b24ad042SBarry Smith   ierr = PetscMalloc((ina->Nbs+1)*sizeof(PetscInt),&lindices);CHKERRQ(ierr);
255b24ad042SBarry Smith   ierr = PetscMemzero(lindices,ina->Nbs*sizeof(PetscInt));CHKERRQ(ierr);
25604f1ad80SBarry Smith   for (i=0; i<B->nbs; i++) {
25704f1ad80SBarry Smith     lindices[garray[i]] = i+1;
25804f1ad80SBarry Smith   }
259992144d0SBarry Smith   no   = inA->rmap->bmapping->n - nt;
260992144d0SBarry Smith   ierr = PetscMalloc((inA->rmap->bmapping->n+1)*sizeof(PetscInt),&r_rmapo);CHKERRQ(ierr);
261992144d0SBarry Smith   ierr = PetscMemzero(r_rmapo,inA->rmap->bmapping->n*sizeof(PetscInt));CHKERRQ(ierr);
26204f1ad80SBarry Smith   nt   = 0;
263992144d0SBarry Smith   for (i=0; i<inA->rmap->bmapping->n; i++) {
264992144d0SBarry Smith     if (lindices[inA->rmap->bmapping->indices[i]]) {
26504f1ad80SBarry Smith       nt++;
266992144d0SBarry Smith       r_rmapo[i] = lindices[inA->rmap->bmapping->indices[i]];
26704f1ad80SBarry Smith     }
26804f1ad80SBarry Smith   }
269e32f2f54SBarry Smith   if (nt > no) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hmm nt %D no %D",nt,n);
27004f1ad80SBarry Smith   ierr = PetscFree(lindices);CHKERRQ(ierr);
271b24ad042SBarry Smith   ierr = PetscMalloc((nt*bs+1)*sizeof(PetscInt),&uglyrmapo);CHKERRQ(ierr);
272992144d0SBarry Smith   for (i=0; i<inA->rmap->bmapping->n; i++) {
27304f1ad80SBarry Smith     if (r_rmapo[i]) {
27404f1ad80SBarry Smith       for (j=0; j<bs; j++) {
27504f1ad80SBarry Smith         uglyrmapo[(r_rmapo[i]-1)*bs+j] = i*bs + j;
27604f1ad80SBarry Smith       }
27704f1ad80SBarry Smith     }
27804f1ad80SBarry Smith   }
27904f1ad80SBarry Smith   ierr = PetscFree(r_rmapo);CHKERRQ(ierr);
28004f1ad80SBarry Smith   ierr = VecCreateSeq(PETSC_COMM_SELF,nt*bs,&uglyoo);CHKERRQ(ierr);
28104f1ad80SBarry Smith   PetscFunctionReturn(0);
28204f1ad80SBarry Smith }
28304f1ad80SBarry Smith 
284348d1a9dSSatish Balay #undef __FUNCT__
285348d1a9dSSatish Balay #define __FUNCT__ "MatMPIBAIJDiagonalScaleLocal"
2867087cfbeSBarry Smith PetscErrorCode  MatMPIBAIJDiagonalScaleLocal(Mat A,Vec scale)
28704f1ad80SBarry Smith {
28892b32695SKris Buschelman   /* This routine should really be abandoned as it duplicates MatDiagonalScaleLocal */
2894ac538c5SBarry Smith   PetscErrorCode ierr;
29092b32695SKris Buschelman 
29192b32695SKris Buschelman   PetscFunctionBegin;
2924ac538c5SBarry Smith   ierr = PetscTryMethod(A,"MatDiagonalScaleLocal_C",(Mat,Vec),(A,scale));CHKERRQ(ierr);
29392b32695SKris Buschelman   PetscFunctionReturn(0);
29492b32695SKris Buschelman }
29592b32695SKris Buschelman 
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 */
3190298fd71SBarry Smith   ierr = MatDiagonalScale(a->A,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 */
3290298fd71SBarry Smith   ierr = MatDiagonalScale(a->B,NULL,uglyoo);CHKERRQ(ierr);
33004f1ad80SBarry Smith   PetscFunctionReturn(0);
33104f1ad80SBarry Smith }
3328016bdd1SSatish Balay 
3333eda8832SBarry Smith 
334