xref: /petsc/src/mat/impls/baij/mpi/mmbaij.c (revision bca11509f5effeca090a0d43c32015d8daa7ad4e)
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>
6af0996ceSBarry Smith #include <petsc/private/isimpl.h>    /* needed because accesses data structure of ISLocalToGlobalMapping directly */
7bba1ac68SSatish Balay 
8dfbe8321SBarry Smith PetscErrorCode MatSetUpMultiply_MPIBAIJ(Mat mat)
98016bdd1SSatish Balay {
10d9653453SSatish Balay   Mat_MPIBAIJ    *baij = (Mat_MPIBAIJ*)mat->data;
11d9653453SSatish Balay   Mat_SeqBAIJ    *B    = (Mat_SeqBAIJ*)(baij->B->data);
126849ba73SBarry Smith   PetscErrorCode ierr;
13b24ad042SBarry Smith   PetscInt       i,j,*aj = B->j,ec = 0,*garray;
14d0f46423SBarry Smith   PetscInt       bs = mat->rmap->bs,*stmp;
158016bdd1SSatish Balay   IS             from,to;
168016bdd1SSatish Balay   Vec            gvec;
17aa482453SBarry Smith #if defined(PETSC_USE_CTABLE)
180f5bd95cSBarry Smith   PetscTable         gid1_lid1;
190f5bd95cSBarry Smith   PetscTablePosition tpos;
20b24ad042SBarry Smith   PetscInt           gid,lid;
216f531f54SSatish Balay #else
22b24ad042SBarry Smith   PetscInt Nbs = baij->Nbs,*indices;
2373a2e727SSatish Balay #endif
248016bdd1SSatish Balay 
253a40ed3dSBarry Smith   PetscFunctionBegin;
26aa482453SBarry Smith #if defined(PETSC_USE_CTABLE)
2773a2e727SSatish Balay   /* use a table - Mark Adams */
28e23dfa41SBarry Smith   ierr = PetscTableCreate(B->mbs,baij->Nbs+1,&gid1_lid1);CHKERRQ(ierr);
2973a2e727SSatish Balay   for (i=0; i<B->mbs; i++) {
3073a2e727SSatish Balay     for (j=0; j<B->ilen[i]; j++) {
31b24ad042SBarry Smith       PetscInt data,gid1 = aj[B->i[i]+j] + 1;
320f5bd95cSBarry Smith       ierr = PetscTableFind(gid1_lid1,gid1,&data);CHKERRQ(ierr);
33fa46199cSSatish Balay       if (!data) {
3473a2e727SSatish Balay         /* one based table */
353861aac3SJed Brown         ierr = PetscTableAdd(gid1_lid1,gid1,++ec,INSERT_VALUES);CHKERRQ(ierr);
3673a2e727SSatish Balay       }
3773a2e727SSatish Balay     }
3873a2e727SSatish Balay   }
3973a2e727SSatish Balay   /* form array of columns we need */
40854ce69bSBarry Smith   ierr = PetscMalloc1(ec+1,&garray);CHKERRQ(ierr);
410f5bd95cSBarry Smith   ierr = PetscTableGetHeadPosition(gid1_lid1,&tpos);CHKERRQ(ierr);
4273a2e727SSatish Balay   while (tpos) {
430f5bd95cSBarry Smith     ierr = PetscTableGetNext(gid1_lid1,&tpos,&gid,&lid);CHKERRQ(ierr);
4473a2e727SSatish Balay     gid--; lid--;
4573a2e727SSatish Balay     garray[lid] = gid;
4673a2e727SSatish Balay   }
470064e2bbSSatish Balay   ierr = PetscSortInt(ec,garray);CHKERRQ(ierr);
480f5bd95cSBarry Smith   ierr = PetscTableRemoveAll(gid1_lid1);CHKERRQ(ierr);
4973a2e727SSatish Balay   for (i=0; i<ec; i++) {
503861aac3SJed Brown     ierr = PetscTableAdd(gid1_lid1,garray[i]+1,i+1,INSERT_VALUES);CHKERRQ(ierr);
5173a2e727SSatish Balay   }
5273a2e727SSatish Balay   /* compact out the extra columns in B */
5373a2e727SSatish Balay   for (i=0; i<B->mbs; i++) {
5473a2e727SSatish Balay     for (j=0; j<B->ilen[i]; j++) {
55b24ad042SBarry Smith       PetscInt gid1 = aj[B->i[i] + j] + 1;
560f5bd95cSBarry Smith       ierr = PetscTableFind(gid1_lid1,gid1,&lid);CHKERRQ(ierr);
57fa46199cSSatish Balay       lid--;
5873a2e727SSatish Balay       aj[B->i[i]+j] = lid;
5973a2e727SSatish Balay     }
6073a2e727SSatish Balay   }
6173a2e727SSatish Balay   B->nbs           = ec;
62d0f46423SBarry Smith   baij->B->cmap->n = baij->B->cmap->N = ec*mat->rmap->bs;
6326fbe8dcSKarl Rupp 
6426283091SBarry Smith   ierr = PetscLayoutSetUp((baij->B->cmap));CHKERRQ(ierr);
656bc0bbbfSBarry Smith   ierr = PetscTableDestroy(&gid1_lid1);CHKERRQ(ierr);
6673a2e727SSatish Balay #else
67435da068SBarry Smith   /* Make an array as long as the number of columns */
68d9653453SSatish Balay   /* mark those columns that are in baij->B */
691795a4d1SJed Brown   ierr = PetscCalloc1(Nbs+1,&indices);CHKERRQ(ierr);
70d9653453SSatish Balay   for (i=0; i<B->mbs; i++) {
718016bdd1SSatish Balay     for (j=0; j<B->ilen[i]; j++) {
72d9653453SSatish Balay       if (!indices[aj[B->i[i] + j]]) ec++;
73d9653453SSatish Balay       indices[aj[B->i[i] + j]] = 1;
748016bdd1SSatish Balay     }
758016bdd1SSatish Balay   }
768016bdd1SSatish Balay 
778016bdd1SSatish Balay   /* form array of columns we need */
78854ce69bSBarry Smith   ierr = PetscMalloc1(ec+1,&garray);CHKERRQ(ierr);
798016bdd1SSatish Balay   ec   = 0;
80d9653453SSatish Balay   for (i=0; i<Nbs; i++) {
810bdbc534SSatish Balay     if (indices[i]) {
820bdbc534SSatish Balay       garray[ec++] = i;
830bdbc534SSatish Balay     }
848016bdd1SSatish Balay   }
858016bdd1SSatish Balay 
868016bdd1SSatish Balay   /* make indices now point into garray */
878016bdd1SSatish Balay   for (i=0; i<ec; i++) {
88d9653453SSatish Balay     indices[garray[i]] = i;
898016bdd1SSatish Balay   }
908016bdd1SSatish Balay 
918016bdd1SSatish Balay   /* compact out the extra columns in B */
92d9653453SSatish Balay   for (i=0; i<B->mbs; i++) {
938016bdd1SSatish Balay     for (j=0; j<B->ilen[i]; j++) {
94d9653453SSatish Balay       aj[B->i[i] + j] = indices[aj[B->i[i] + j]];
958016bdd1SSatish Balay     }
968016bdd1SSatish Balay   }
97d9653453SSatish Balay   B->nbs           = ec;
98d0f46423SBarry Smith   baij->B->cmap->n = baij->B->cmap->N  = ec*mat->rmap->bs;
9926fbe8dcSKarl Rupp 
10026283091SBarry Smith   ierr = PetscLayoutSetUp((baij->B->cmap));CHKERRQ(ierr);
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 */
108deff0451SBarry Smith   ierr = ISCreateBlock(PETSC_COMM_SELF,bs,ec,garray,PETSC_COPY_VALUES,&from);CHKERRQ(ierr);
109c16cb8f2SBarry Smith 
110854ce69bSBarry Smith   ierr = PetscMalloc1(ec+1,&stmp);CHKERRQ(ierr);
11126fbe8dcSKarl Rupp   for (i=0; i<ec; i++) stmp[i] = i;
112deff0451SBarry Smith   ierr = ISCreateBlock(PETSC_COMM_SELF,bs,ec,stmp,PETSC_OWN_POINTER,&to);CHKERRQ(ierr);
1138016bdd1SSatish Balay 
1148016bdd1SSatish Balay   /* create temporary global vector to generate scatter context */
115ce94432eSBarry Smith   ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)mat),1,mat->cmap->n,mat->cmap->N,NULL,&gvec);CHKERRQ(ierr);
1168016bdd1SSatish Balay 
117d9653453SSatish Balay   ierr = VecScatterCreate(gvec,from,baij->lvec,to,&baij->Mvctx);CHKERRQ(ierr);
11890f02eecSBarry Smith 
1193bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)baij->Mvctx);CHKERRQ(ierr);
1203bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)baij->lvec);CHKERRQ(ierr);
1213bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)from);CHKERRQ(ierr);
1223bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)to);CHKERRQ(ierr);
12326fbe8dcSKarl Rupp 
124d9653453SSatish Balay   baij->garray = garray;
12526fbe8dcSKarl Rupp 
1263bb1ff40SBarry Smith   ierr = PetscLogObjectMemory((PetscObject)mat,(ec+1)*sizeof(PetscInt));CHKERRQ(ierr);
1276bf464f9SBarry Smith   ierr = ISDestroy(&from);CHKERRQ(ierr);
1286bf464f9SBarry Smith   ierr = ISDestroy(&to);CHKERRQ(ierr);
1296bf464f9SBarry Smith   ierr = VecDestroy(&gvec);CHKERRQ(ierr);
1303a40ed3dSBarry Smith   PetscFunctionReturn(0);
1318016bdd1SSatish Balay }
1328016bdd1SSatish Balay 
1338016bdd1SSatish Balay /*
134d9653453SSatish Balay      Takes the local part of an already assembled MPIBAIJ matrix
1358016bdd1SSatish Balay    and disassembles it. This is to allow new nonzeros into the matrix
1368016bdd1SSatish Balay    that require more communication in the matrix vector multiply.
1378016bdd1SSatish Balay    Thus certain data-structures must be rebuilt.
1388016bdd1SSatish Balay 
1398016bdd1SSatish Balay    Kind of slow! But that's what application programmers get when
1408016bdd1SSatish Balay    they are sloppy.
1418016bdd1SSatish Balay */
142ab9863d7SBarry Smith PetscErrorCode MatDisAssemble_MPIBAIJ(Mat A)
1438016bdd1SSatish Balay {
144d9653453SSatish Balay   Mat_MPIBAIJ    *baij  = (Mat_MPIBAIJ*)A->data;
145d9653453SSatish Balay   Mat            B      = baij->B,Bnew;
146d9653453SSatish Balay   Mat_SeqBAIJ    *Bbaij = (Mat_SeqBAIJ*)B->data;
147dfbe8321SBarry Smith   PetscErrorCode ierr;
148d0f46423SBarry Smith   PetscInt       i,j,mbs=Bbaij->mbs,n = A->cmap->N,col,*garray=baij->garray;
149d0f46423SBarry Smith   PetscInt       bs2 = baij->bs2,*nz,ec,m = A->rmap->n;
1503eda8832SBarry Smith   MatScalar      *a  = Bbaij->a;
151dd6ea824SBarry Smith   MatScalar      *atmp;
15297e5c40aSBarry Smith 
1538016bdd1SSatish Balay 
1543a40ed3dSBarry Smith   PetscFunctionBegin;
1558016bdd1SSatish Balay   /* free stuff related to matrix-vec multiply */
156b0a32e0cSBarry Smith   ierr = VecGetSize(baij->lvec,&ec);CHKERRQ(ierr); /* needed for PetscLogObjectMemory below */
1576bf464f9SBarry Smith   ierr = VecDestroy(&baij->lvec);CHKERRQ(ierr); baij->lvec = 0;
1586bf464f9SBarry Smith   ierr = VecScatterDestroy(&baij->Mvctx);CHKERRQ(ierr); baij->Mvctx = 0;
159d9653453SSatish Balay   if (baij->colmap) {
160aa482453SBarry Smith #if defined(PETSC_USE_CTABLE)
1616bc0bbbfSBarry Smith     ierr = PetscTableDestroy(&baij->colmap);CHKERRQ(ierr);
16248e59246SSatish Balay #else
163606d414cSSatish Balay     ierr = PetscFree(baij->colmap);CHKERRQ(ierr);
1643bb1ff40SBarry Smith     ierr = PetscLogObjectMemory((PetscObject)A,-Bbaij->nbs*sizeof(PetscInt));CHKERRQ(ierr);
16548e59246SSatish Balay #endif
1668016bdd1SSatish Balay   }
1678016bdd1SSatish Balay 
1688016bdd1SSatish Balay   /* make sure that B is assembled so we can access its values */
1696d4a8577SBarry Smith   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1703eda8832SBarry Smith   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1718016bdd1SSatish Balay 
1728016bdd1SSatish Balay   /* invent new B and copy stuff over */
173785e854fSJed Brown   ierr = PetscMalloc1(mbs,&nz);CHKERRQ(ierr);
174d9653453SSatish Balay   for (i=0; i<mbs; i++) {
175d9653453SSatish Balay     nz[i] = Bbaij->i[i+1]-Bbaij->i[i];
1768016bdd1SSatish Balay   }
177ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)B),&Bnew);CHKERRQ(ierr);
178f69a0ea3SMatthew Knepley   ierr = MatSetSizes(Bnew,m,n,m,n);CHKERRQ(ierr);
1797adad957SLisandro Dalcin   ierr = MatSetType(Bnew,((PetscObject)B)->type_name);CHKERRQ(ierr);
180d0f46423SBarry Smith   ierr = MatSeqBAIJSetPreallocation(Bnew,B->rmap->bs,0,nz);CHKERRQ(ierr);
18126fbe8dcSKarl Rupp 
1822576faa2SJed Brown   ((Mat_SeqBAIJ*)Bnew->data)->nonew = Bbaij->nonew; /* Inherit insertion error options. */
18326fbe8dcSKarl Rupp 
1844e0d8c25SBarry Smith   ierr = MatSetOption(Bnew,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr);
18577341eacSDmitry Karpeev   /*
18677341eacSDmitry Karpeev    Ensure that B's nonzerostate is monotonically increasing.
18777341eacSDmitry Karpeev    Or should this follow the MatSetValuesBlocked() loop to preserve B's nonzerstate across a MatDisAssemble() call?
18877341eacSDmitry Karpeev    */
18977341eacSDmitry Karpeev   Bnew->nonzerostate = B->nonzerostate;
190d9653453SSatish Balay 
191bba1ac68SSatish Balay   for (i=0; i<mbs; i++) {
192bba1ac68SSatish Balay     for (j=Bbaij->i[i]; j<Bbaij->i[i+1]; j++) {
193bba1ac68SSatish Balay       col  = garray[Bbaij->j[j]];
1943eda8832SBarry Smith       atmp = a + j*bs2;
195bba1ac68SSatish Balay       ierr = MatSetValuesBlocked_SeqBAIJ(Bnew,1,&i,1,&col,atmp,B->insertmode);CHKERRQ(ierr);
1968016bdd1SSatish Balay     }
1978016bdd1SSatish Balay   }
1984e0d8c25SBarry Smith   ierr = MatSetOption(Bnew,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr);
199bba1ac68SSatish Balay 
200bba1ac68SSatish Balay   ierr = PetscFree(nz);CHKERRQ(ierr);
201606d414cSSatish Balay   ierr = PetscFree(baij->garray);CHKERRQ(ierr);
2023bb1ff40SBarry Smith   ierr = PetscLogObjectMemory((PetscObject)A,-ec*sizeof(PetscInt));CHKERRQ(ierr);
2036bf464f9SBarry Smith   ierr = MatDestroy(&B);CHKERRQ(ierr);
2043bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)A,(PetscObject)Bnew);CHKERRQ(ierr);
20526fbe8dcSKarl Rupp 
206d9653453SSatish Balay   baij->B          = Bnew;
2078016bdd1SSatish Balay   A->was_assembled = PETSC_FALSE;
2086a719282SBarry Smith   A->assembled     = PETSC_FALSE;
2093a40ed3dSBarry Smith   PetscFunctionReturn(0);
2108016bdd1SSatish Balay }
2118016bdd1SSatish Balay 
21204f1ad80SBarry Smith /*      ugly stuff added for Glenn someday we should fix this up */
21304f1ad80SBarry Smith 
21426fbe8dcSKarl Rupp static PetscInt *uglyrmapd = 0,*uglyrmapo = 0;  /* mapping from the local ordering to the "diagonal" and "off-diagonal" parts of the local matrix */
2152cd6534aSBarry Smith static Vec      uglydd     = 0,uglyoo     = 0;  /* work vectors used to scale the two parts of the local matrix */
21604f1ad80SBarry Smith 
21704f1ad80SBarry Smith 
218dfbe8321SBarry Smith PetscErrorCode MatMPIBAIJDiagonalScaleLocalSetUp(Mat inA,Vec scale)
21904f1ad80SBarry Smith {
22004f1ad80SBarry Smith   Mat_MPIBAIJ    *ina = (Mat_MPIBAIJ*) inA->data; /*access private part of matrix */
22104f1ad80SBarry Smith   Mat_SeqBAIJ    *B   = (Mat_SeqBAIJ*)ina->B->data;
222dfbe8321SBarry Smith   PetscErrorCode ierr;
223d0f46423SBarry Smith   PetscInt       bs = inA->rmap->bs,i,n,nt,j,cstart,cend,no,*garray = ina->garray,*lindices;
224b24ad042SBarry Smith   PetscInt       *r_rmapd,*r_rmapo;
22504f1ad80SBarry Smith 
22604f1ad80SBarry Smith   PetscFunctionBegin;
22704f1ad80SBarry Smith   ierr = MatGetOwnershipRange(inA,&cstart,&cend);CHKERRQ(ierr);
2280298fd71SBarry Smith   ierr = MatGetSize(ina->A,NULL,&n);CHKERRQ(ierr);
229854ce69bSBarry Smith   ierr = PetscCalloc1(inA->rmap->mapping->n+1,&r_rmapd);CHKERRQ(ierr);
23004f1ad80SBarry Smith   nt   = 0;
23145b6f7e9SBarry Smith   for (i=0; i<inA->rmap->mapping->n; i++) {
23245b6f7e9SBarry Smith     if (inA->rmap->mapping->indices[i]*bs >= cstart && inA->rmap->mapping->indices[i]*bs < cend) {
23304f1ad80SBarry Smith       nt++;
23445b6f7e9SBarry Smith       r_rmapd[i] = inA->rmap->mapping->indices[i] + 1;
23504f1ad80SBarry Smith     }
23604f1ad80SBarry Smith   }
237e32f2f54SBarry Smith   if (nt*bs != n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hmm nt*bs %D n %D",nt*bs,n);
238854ce69bSBarry Smith   ierr = PetscMalloc1(n+1,&uglyrmapd);CHKERRQ(ierr);
23945b6f7e9SBarry Smith   for (i=0; i<inA->rmap->mapping->n; i++) {
24004f1ad80SBarry Smith     if (r_rmapd[i]) {
24104f1ad80SBarry Smith       for (j=0; j<bs; j++) {
24204f1ad80SBarry Smith         uglyrmapd[(r_rmapd[i]-1)*bs+j-cstart] = i*bs + j;
24304f1ad80SBarry Smith       }
24404f1ad80SBarry Smith     }
24504f1ad80SBarry Smith   }
24604f1ad80SBarry Smith   ierr = PetscFree(r_rmapd);CHKERRQ(ierr);
24704f1ad80SBarry Smith   ierr = VecCreateSeq(PETSC_COMM_SELF,n,&uglydd);CHKERRQ(ierr);
24804f1ad80SBarry Smith 
249854ce69bSBarry Smith   ierr = PetscCalloc1(ina->Nbs+1,&lindices);CHKERRQ(ierr);
25004f1ad80SBarry Smith   for (i=0; i<B->nbs; i++) {
25104f1ad80SBarry Smith     lindices[garray[i]] = i+1;
25204f1ad80SBarry Smith   }
25345b6f7e9SBarry Smith   no   = inA->rmap->mapping->n - nt;
254854ce69bSBarry Smith   ierr = PetscCalloc1(inA->rmap->mapping->n+1,&r_rmapo);CHKERRQ(ierr);
25504f1ad80SBarry Smith   nt   = 0;
25645b6f7e9SBarry Smith   for (i=0; i<inA->rmap->mapping->n; i++) {
25745b6f7e9SBarry Smith     if (lindices[inA->rmap->mapping->indices[i]]) {
25804f1ad80SBarry Smith       nt++;
25945b6f7e9SBarry Smith       r_rmapo[i] = lindices[inA->rmap->mapping->indices[i]];
26004f1ad80SBarry Smith     }
26104f1ad80SBarry Smith   }
262e32f2f54SBarry Smith   if (nt > no) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hmm nt %D no %D",nt,n);
26304f1ad80SBarry Smith   ierr = PetscFree(lindices);CHKERRQ(ierr);
264854ce69bSBarry Smith   ierr = PetscMalloc1(nt*bs+1,&uglyrmapo);CHKERRQ(ierr);
26545b6f7e9SBarry Smith   for (i=0; i<inA->rmap->mapping->n; i++) {
26604f1ad80SBarry Smith     if (r_rmapo[i]) {
26704f1ad80SBarry Smith       for (j=0; j<bs; j++) {
26804f1ad80SBarry Smith         uglyrmapo[(r_rmapo[i]-1)*bs+j] = i*bs + j;
26904f1ad80SBarry Smith       }
27004f1ad80SBarry Smith     }
27104f1ad80SBarry Smith   }
27204f1ad80SBarry Smith   ierr = PetscFree(r_rmapo);CHKERRQ(ierr);
27304f1ad80SBarry Smith   ierr = VecCreateSeq(PETSC_COMM_SELF,nt*bs,&uglyoo);CHKERRQ(ierr);
27404f1ad80SBarry Smith   PetscFunctionReturn(0);
27504f1ad80SBarry Smith }
27604f1ad80SBarry Smith 
2777087cfbeSBarry Smith PetscErrorCode  MatMPIBAIJDiagonalScaleLocal(Mat A,Vec scale)
27804f1ad80SBarry Smith {
27992b32695SKris Buschelman   /* This routine should really be abandoned as it duplicates MatDiagonalScaleLocal */
2804ac538c5SBarry Smith   PetscErrorCode ierr;
28192b32695SKris Buschelman 
28292b32695SKris Buschelman   PetscFunctionBegin;
2834ac538c5SBarry Smith   ierr = PetscTryMethod(A,"MatDiagonalScaleLocal_C",(Mat,Vec),(A,scale));CHKERRQ(ierr);
28492b32695SKris Buschelman   PetscFunctionReturn(0);
28592b32695SKris Buschelman }
28692b32695SKris Buschelman 
2877087cfbeSBarry Smith PetscErrorCode  MatDiagonalScaleLocal_MPIBAIJ(Mat A,Vec scale)
28892b32695SKris Buschelman {
28904f1ad80SBarry Smith   Mat_MPIBAIJ       *a = (Mat_MPIBAIJ*) A->data; /*access private part of matrix */
290dfbe8321SBarry Smith   PetscErrorCode    ierr;
291b24ad042SBarry Smith   PetscInt          n,i;
292*bca11509SBarry Smith   PetscScalar       *d,*o;
293*bca11509SBarry Smith   const PetscScalar *s;
29404f1ad80SBarry Smith 
29504f1ad80SBarry Smith   PetscFunctionBegin;
2962cd6534aSBarry Smith   if (!uglyrmapd) {
2972cd6534aSBarry Smith     ierr = MatMPIBAIJDiagonalScaleLocalSetUp(A,scale);CHKERRQ(ierr);
2982cd6534aSBarry Smith   }
2992cd6534aSBarry Smith 
300*bca11509SBarry Smith   ierr = VecGetArrayRead(scale,&s);CHKERRQ(ierr);
30104f1ad80SBarry Smith 
30204f1ad80SBarry Smith   ierr = VecGetLocalSize(uglydd,&n);CHKERRQ(ierr);
3031ebc52fbSHong Zhang   ierr = VecGetArray(uglydd,&d);CHKERRQ(ierr);
30404f1ad80SBarry Smith   for (i=0; i<n; i++) {
30504f1ad80SBarry Smith     d[i] = s[uglyrmapd[i]]; /* copy "diagonal" (true local) portion of scale into dd vector */
30604f1ad80SBarry Smith   }
3071ebc52fbSHong Zhang   ierr = VecRestoreArray(uglydd,&d);CHKERRQ(ierr);
30804f1ad80SBarry Smith   /* column scale "diagonal" portion of local matrix */
3090298fd71SBarry Smith   ierr = MatDiagonalScale(a->A,NULL,uglydd);CHKERRQ(ierr);
31004f1ad80SBarry Smith 
31104f1ad80SBarry Smith   ierr = VecGetLocalSize(uglyoo,&n);CHKERRQ(ierr);
3121ebc52fbSHong Zhang   ierr = VecGetArray(uglyoo,&o);CHKERRQ(ierr);
31304f1ad80SBarry Smith   for (i=0; i<n; i++) {
31404f1ad80SBarry Smith     o[i] = s[uglyrmapo[i]]; /* copy "off-diagonal" portion of scale into oo vector */
31504f1ad80SBarry Smith   }
316*bca11509SBarry Smith   ierr = VecRestoreArrayRead(scale,&s);CHKERRQ(ierr);
3171ebc52fbSHong Zhang   ierr = VecRestoreArray(uglyoo,&o);CHKERRQ(ierr);
31804f1ad80SBarry Smith   /* column scale "off-diagonal" portion of local matrix */
3190298fd71SBarry Smith   ierr = MatDiagonalScale(a->B,NULL,uglyoo);CHKERRQ(ierr);
32004f1ad80SBarry Smith   PetscFunctionReturn(0);
32104f1ad80SBarry Smith }
322