xref: /petsc/src/mat/impls/baij/mpi/mmbaij.c (revision cac4c232dc4f93991e342196e27ef7b0655dac7b)
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);
12b24ad042SBarry Smith   PetscInt       i,j,*aj = B->j,ec = 0,*garray;
13d0f46423SBarry Smith   PetscInt       bs = mat->rmap->bs,*stmp;
148016bdd1SSatish Balay   IS             from,to;
158016bdd1SSatish Balay   Vec            gvec;
16aa482453SBarry Smith #if defined(PETSC_USE_CTABLE)
170f5bd95cSBarry Smith   PetscTable         gid1_lid1;
180f5bd95cSBarry Smith   PetscTablePosition tpos;
19b24ad042SBarry Smith   PetscInt           gid,lid;
206f531f54SSatish Balay #else
21b24ad042SBarry Smith   PetscInt Nbs = baij->Nbs,*indices;
2273a2e727SSatish Balay #endif
238016bdd1SSatish Balay 
243a40ed3dSBarry Smith   PetscFunctionBegin;
25aa482453SBarry Smith #if defined(PETSC_USE_CTABLE)
2673a2e727SSatish Balay   /* use a table - Mark Adams */
279566063dSJacob Faibussowitsch   PetscCall(PetscTableCreate(B->mbs,baij->Nbs+1,&gid1_lid1));
2873a2e727SSatish Balay   for (i=0; i<B->mbs; i++) {
2973a2e727SSatish Balay     for (j=0; j<B->ilen[i]; j++) {
30b24ad042SBarry Smith       PetscInt data,gid1 = aj[B->i[i]+j] + 1;
319566063dSJacob Faibussowitsch       PetscCall(PetscTableFind(gid1_lid1,gid1,&data));
32fa46199cSSatish Balay       if (!data) {
3373a2e727SSatish Balay         /* one based table */
349566063dSJacob Faibussowitsch         PetscCall(PetscTableAdd(gid1_lid1,gid1,++ec,INSERT_VALUES));
3573a2e727SSatish Balay       }
3673a2e727SSatish Balay     }
3773a2e727SSatish Balay   }
3873a2e727SSatish Balay   /* form array of columns we need */
399566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(ec,&garray));
409566063dSJacob Faibussowitsch   PetscCall(PetscTableGetHeadPosition(gid1_lid1,&tpos));
4173a2e727SSatish Balay   while (tpos) {
429566063dSJacob Faibussowitsch     PetscCall(PetscTableGetNext(gid1_lid1,&tpos,&gid,&lid));
4373a2e727SSatish Balay     gid--; lid--;
4473a2e727SSatish Balay     garray[lid] = gid;
4573a2e727SSatish Balay   }
469566063dSJacob Faibussowitsch   PetscCall(PetscSortInt(ec,garray));
479566063dSJacob Faibussowitsch   PetscCall(PetscTableRemoveAll(gid1_lid1));
4873a2e727SSatish Balay   for (i=0; i<ec; i++) {
499566063dSJacob Faibussowitsch     PetscCall(PetscTableAdd(gid1_lid1,garray[i]+1,i+1,INSERT_VALUES));
5073a2e727SSatish Balay   }
5173a2e727SSatish Balay   /* compact out the extra columns in B */
5273a2e727SSatish Balay   for (i=0; i<B->mbs; i++) {
5373a2e727SSatish Balay     for (j=0; j<B->ilen[i]; j++) {
54b24ad042SBarry Smith       PetscInt gid1 = aj[B->i[i] + j] + 1;
559566063dSJacob Faibussowitsch       PetscCall(PetscTableFind(gid1_lid1,gid1,&lid));
56fa46199cSSatish Balay       lid--;
5773a2e727SSatish Balay       aj[B->i[i]+j] = lid;
5873a2e727SSatish Balay     }
5973a2e727SSatish Balay   }
6073a2e727SSatish Balay   B->nbs           = ec;
619566063dSJacob Faibussowitsch   PetscCall(PetscLayoutDestroy(&baij->B->cmap));
629566063dSJacob Faibussowitsch   PetscCall(PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)baij->B),ec*mat->rmap->bs,ec*mat->rmap->bs,mat->rmap->bs,&baij->B->cmap));
639566063dSJacob Faibussowitsch   PetscCall(PetscTableDestroy(&gid1_lid1));
6473a2e727SSatish Balay #else
65435da068SBarry Smith   /* Make an array as long as the number of columns */
66d9653453SSatish Balay   /* mark those columns that are in baij->B */
679566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(Nbs,&indices));
68d9653453SSatish Balay   for (i=0; i<B->mbs; i++) {
698016bdd1SSatish Balay     for (j=0; j<B->ilen[i]; j++) {
70d9653453SSatish Balay       if (!indices[aj[B->i[i] + j]]) ec++;
71d9653453SSatish Balay       indices[aj[B->i[i] + j]] = 1;
728016bdd1SSatish Balay     }
738016bdd1SSatish Balay   }
748016bdd1SSatish Balay 
758016bdd1SSatish Balay   /* form array of columns we need */
769566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(ec,&garray));
778016bdd1SSatish Balay   ec   = 0;
78d9653453SSatish Balay   for (i=0; i<Nbs; i++) {
790bdbc534SSatish Balay     if (indices[i]) {
800bdbc534SSatish Balay       garray[ec++] = i;
810bdbc534SSatish Balay     }
828016bdd1SSatish Balay   }
838016bdd1SSatish Balay 
848016bdd1SSatish Balay   /* make indices now point into garray */
858016bdd1SSatish Balay   for (i=0; i<ec; i++) {
86d9653453SSatish Balay     indices[garray[i]] = i;
878016bdd1SSatish Balay   }
888016bdd1SSatish Balay 
898016bdd1SSatish Balay   /* compact out the extra columns in B */
90d9653453SSatish Balay   for (i=0; i<B->mbs; i++) {
918016bdd1SSatish Balay     for (j=0; j<B->ilen[i]; j++) {
92d9653453SSatish Balay       aj[B->i[i] + j] = indices[aj[B->i[i] + j]];
938016bdd1SSatish Balay     }
948016bdd1SSatish Balay   }
95d9653453SSatish Balay   B->nbs           = ec;
969566063dSJacob Faibussowitsch   PetscCall(PetscLayoutDestroy(&baij->B->cmap));
979566063dSJacob Faibussowitsch   PetscCall(PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)baij->B),ec*mat->rmap->bs,ec*mat->rmap->bs,mat->rmap->bs,&baij->B->cmap));
989566063dSJacob Faibussowitsch   PetscCall(PetscFree(indices));
9973a2e727SSatish Balay #endif
1008016bdd1SSatish Balay 
1018016bdd1SSatish Balay   /* create local vector that is used to scatter into */
1029566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF,ec*bs,&baij->lvec));
1038016bdd1SSatish Balay 
104c16cb8f2SBarry Smith   /* create two temporary index sets for building scatter-gather */
1059566063dSJacob Faibussowitsch   PetscCall(ISCreateBlock(PETSC_COMM_SELF,bs,ec,garray,PETSC_COPY_VALUES,&from));
106c16cb8f2SBarry Smith 
1079566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(ec,&stmp));
10826fbe8dcSKarl Rupp   for (i=0; i<ec; i++) stmp[i] = i;
1099566063dSJacob Faibussowitsch   PetscCall(ISCreateBlock(PETSC_COMM_SELF,bs,ec,stmp,PETSC_OWN_POINTER,&to));
1108016bdd1SSatish Balay 
1118016bdd1SSatish Balay   /* create temporary global vector to generate scatter context */
1129566063dSJacob Faibussowitsch   PetscCall(VecCreateMPIWithArray(PetscObjectComm((PetscObject)mat),1,mat->cmap->n,mat->cmap->N,NULL,&gvec));
1138016bdd1SSatish Balay 
1149566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(gvec,from,baij->lvec,to,&baij->Mvctx));
1159566063dSJacob Faibussowitsch   PetscCall(VecScatterViewFromOptions(baij->Mvctx,(PetscObject)mat,"-matmult_vecscatter_view"));
11690f02eecSBarry Smith 
1179566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectParent((PetscObject)mat,(PetscObject)baij->Mvctx));
1189566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectParent((PetscObject)mat,(PetscObject)baij->lvec));
1199566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectParent((PetscObject)mat,(PetscObject)from));
1209566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectParent((PetscObject)mat,(PetscObject)to));
12126fbe8dcSKarl Rupp 
122d9653453SSatish Balay   baij->garray = garray;
12326fbe8dcSKarl Rupp 
1249566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectMemory((PetscObject)mat,ec*sizeof(PetscInt)));
1259566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&from));
1269566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&to));
1279566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&gvec));
1283a40ed3dSBarry Smith   PetscFunctionReturn(0);
1298016bdd1SSatish Balay }
1308016bdd1SSatish Balay 
1318016bdd1SSatish Balay /*
132d9653453SSatish Balay      Takes the local part of an already assembled MPIBAIJ matrix
1338016bdd1SSatish Balay    and disassembles it. This is to allow new nonzeros into the matrix
1348016bdd1SSatish Balay    that require more communication in the matrix vector multiply.
1358016bdd1SSatish Balay    Thus certain data-structures must be rebuilt.
1368016bdd1SSatish Balay 
1378016bdd1SSatish Balay    Kind of slow! But that's what application programmers get when
1388016bdd1SSatish Balay    they are sloppy.
1398016bdd1SSatish Balay */
140ab9863d7SBarry Smith PetscErrorCode MatDisAssemble_MPIBAIJ(Mat A)
1418016bdd1SSatish Balay {
142d9653453SSatish Balay   Mat_MPIBAIJ    *baij  = (Mat_MPIBAIJ*)A->data;
143d9653453SSatish Balay   Mat            B      = baij->B,Bnew;
144d9653453SSatish Balay   Mat_SeqBAIJ    *Bbaij = (Mat_SeqBAIJ*)B->data;
145d0f46423SBarry Smith   PetscInt       i,j,mbs=Bbaij->mbs,n = A->cmap->N,col,*garray=baij->garray;
146d0f46423SBarry Smith   PetscInt       bs2 = baij->bs2,*nz,ec,m = A->rmap->n;
1473eda8832SBarry Smith   MatScalar      *a  = Bbaij->a;
148dd6ea824SBarry Smith   MatScalar      *atmp;
14997e5c40aSBarry Smith 
1503a40ed3dSBarry Smith   PetscFunctionBegin;
1518016bdd1SSatish Balay   /* free stuff related to matrix-vec multiply */
1529566063dSJacob Faibussowitsch   PetscCall(VecGetSize(baij->lvec,&ec)); /* needed for PetscLogObjectMemory below */
1539566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&baij->lvec)); baij->lvec = NULL;
1549566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&baij->Mvctx)); baij->Mvctx = NULL;
155d9653453SSatish Balay   if (baij->colmap) {
156aa482453SBarry Smith #if defined(PETSC_USE_CTABLE)
1579566063dSJacob Faibussowitsch     PetscCall(PetscTableDestroy(&baij->colmap));
15848e59246SSatish Balay #else
1599566063dSJacob Faibussowitsch     PetscCall(PetscFree(baij->colmap));
1609566063dSJacob Faibussowitsch     PetscCall(PetscLogObjectMemory((PetscObject)A,-Bbaij->nbs*sizeof(PetscInt)));
16148e59246SSatish Balay #endif
1628016bdd1SSatish Balay   }
1638016bdd1SSatish Balay 
1648016bdd1SSatish Balay   /* make sure that B is assembled so we can access its values */
1659566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
1669566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
1678016bdd1SSatish Balay 
1688016bdd1SSatish Balay   /* invent new B and copy stuff over */
1699566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(mbs,&nz));
170d9653453SSatish Balay   for (i=0; i<mbs; i++) {
171d9653453SSatish Balay     nz[i] = Bbaij->i[i+1]-Bbaij->i[i];
1728016bdd1SSatish Balay   }
1739566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)B),&Bnew));
1749566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(Bnew,m,n,m,n));
1759566063dSJacob Faibussowitsch   PetscCall(MatSetType(Bnew,((PetscObject)B)->type_name));
1769566063dSJacob Faibussowitsch   PetscCall(MatSeqBAIJSetPreallocation(Bnew,B->rmap->bs,0,nz));
177b38c15b3SStefano Zampini   if (Bbaij->nonew >= 0) { /* Inherit insertion error options (if positive). */
178b38c15b3SStefano Zampini     ((Mat_SeqBAIJ*)Bnew->data)->nonew = Bbaij->nonew;
179b38c15b3SStefano Zampini   }
18026fbe8dcSKarl Rupp 
1819566063dSJacob Faibussowitsch   PetscCall(MatSetOption(Bnew,MAT_ROW_ORIENTED,PETSC_FALSE));
18277341eacSDmitry Karpeev   /*
18377341eacSDmitry Karpeev    Ensure that B's nonzerostate is monotonically increasing.
18477341eacSDmitry Karpeev    Or should this follow the MatSetValuesBlocked() loop to preserve B's nonzerstate across a MatDisAssemble() call?
18577341eacSDmitry Karpeev    */
18677341eacSDmitry Karpeev   Bnew->nonzerostate = B->nonzerostate;
187d9653453SSatish Balay 
188bba1ac68SSatish Balay   for (i=0; i<mbs; i++) {
189bba1ac68SSatish Balay     for (j=Bbaij->i[i]; j<Bbaij->i[i+1]; j++) {
190bba1ac68SSatish Balay       col  = garray[Bbaij->j[j]];
1913eda8832SBarry Smith       atmp = a + j*bs2;
1929566063dSJacob Faibussowitsch       PetscCall(MatSetValuesBlocked_SeqBAIJ(Bnew,1,&i,1,&col,atmp,B->insertmode));
1938016bdd1SSatish Balay     }
1948016bdd1SSatish Balay   }
1959566063dSJacob Faibussowitsch   PetscCall(MatSetOption(Bnew,MAT_ROW_ORIENTED,PETSC_TRUE));
196bba1ac68SSatish Balay 
1979566063dSJacob Faibussowitsch   PetscCall(PetscFree(nz));
1989566063dSJacob Faibussowitsch   PetscCall(PetscFree(baij->garray));
1999566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectMemory((PetscObject)A,-ec*sizeof(PetscInt)));
2009566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
2019566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectParent((PetscObject)A,(PetscObject)Bnew));
20226fbe8dcSKarl Rupp 
203d9653453SSatish Balay   baij->B          = Bnew;
2048016bdd1SSatish Balay   A->was_assembled = PETSC_FALSE;
2056a719282SBarry Smith   A->assembled     = PETSC_FALSE;
2063a40ed3dSBarry Smith   PetscFunctionReturn(0);
2078016bdd1SSatish Balay }
2088016bdd1SSatish Balay 
20904f1ad80SBarry Smith /*      ugly stuff added for Glenn someday we should fix this up */
21004f1ad80SBarry Smith 
211f4259b30SLisandro Dalcin static PetscInt *uglyrmapd = NULL,*uglyrmapo = NULL;  /* mapping from the local ordering to the "diagonal" and "off-diagonal" parts of the local matrix */
212f4259b30SLisandro Dalcin static Vec      uglydd     = NULL,uglyoo     = NULL;  /* work vectors used to scale the two parts of the local matrix */
21304f1ad80SBarry Smith 
214dfbe8321SBarry Smith PetscErrorCode MatMPIBAIJDiagonalScaleLocalSetUp(Mat inA,Vec scale)
21504f1ad80SBarry Smith {
21604f1ad80SBarry Smith   Mat_MPIBAIJ    *ina = (Mat_MPIBAIJ*) inA->data; /*access private part of matrix */
21704f1ad80SBarry Smith   Mat_SeqBAIJ    *B   = (Mat_SeqBAIJ*)ina->B->data;
218d0f46423SBarry Smith   PetscInt       bs = inA->rmap->bs,i,n,nt,j,cstart,cend,no,*garray = ina->garray,*lindices;
219b24ad042SBarry Smith   PetscInt       *r_rmapd,*r_rmapo;
22004f1ad80SBarry Smith 
22104f1ad80SBarry Smith   PetscFunctionBegin;
2229566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(inA,&cstart,&cend));
2239566063dSJacob Faibussowitsch   PetscCall(MatGetSize(ina->A,NULL,&n));
2249566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(inA->rmap->mapping->n+1,&r_rmapd));
22504f1ad80SBarry Smith   nt   = 0;
22645b6f7e9SBarry Smith   for (i=0; i<inA->rmap->mapping->n; i++) {
22745b6f7e9SBarry Smith     if (inA->rmap->mapping->indices[i]*bs >= cstart && inA->rmap->mapping->indices[i]*bs < cend) {
22804f1ad80SBarry Smith       nt++;
22945b6f7e9SBarry Smith       r_rmapd[i] = inA->rmap->mapping->indices[i] + 1;
23004f1ad80SBarry Smith     }
23104f1ad80SBarry Smith   }
2322c71b3e2SJacob Faibussowitsch   PetscCheckFalse(nt*bs != n,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hmm nt*bs %" PetscInt_FMT " n %" PetscInt_FMT,nt*bs,n);
2339566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n+1,&uglyrmapd));
23445b6f7e9SBarry Smith   for (i=0; i<inA->rmap->mapping->n; i++) {
23504f1ad80SBarry Smith     if (r_rmapd[i]) {
23604f1ad80SBarry Smith       for (j=0; j<bs; j++) {
23704f1ad80SBarry Smith         uglyrmapd[(r_rmapd[i]-1)*bs+j-cstart] = i*bs + j;
23804f1ad80SBarry Smith       }
23904f1ad80SBarry Smith     }
24004f1ad80SBarry Smith   }
2419566063dSJacob Faibussowitsch   PetscCall(PetscFree(r_rmapd));
2429566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF,n,&uglydd));
24304f1ad80SBarry Smith 
2449566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(ina->Nbs+1,&lindices));
24504f1ad80SBarry Smith   for (i=0; i<B->nbs; i++) {
24604f1ad80SBarry Smith     lindices[garray[i]] = i+1;
24704f1ad80SBarry Smith   }
24845b6f7e9SBarry Smith   no   = inA->rmap->mapping->n - nt;
2499566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(inA->rmap->mapping->n+1,&r_rmapo));
25004f1ad80SBarry Smith   nt   = 0;
25145b6f7e9SBarry Smith   for (i=0; i<inA->rmap->mapping->n; i++) {
25245b6f7e9SBarry Smith     if (lindices[inA->rmap->mapping->indices[i]]) {
25304f1ad80SBarry Smith       nt++;
25445b6f7e9SBarry Smith       r_rmapo[i] = lindices[inA->rmap->mapping->indices[i]];
25504f1ad80SBarry Smith     }
25604f1ad80SBarry Smith   }
2572c71b3e2SJacob Faibussowitsch   PetscCheckFalse(nt > no,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hmm nt %" PetscInt_FMT " no %" PetscInt_FMT,nt,n);
2589566063dSJacob Faibussowitsch   PetscCall(PetscFree(lindices));
2599566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nt*bs+1,&uglyrmapo));
26045b6f7e9SBarry Smith   for (i=0; i<inA->rmap->mapping->n; i++) {
26104f1ad80SBarry Smith     if (r_rmapo[i]) {
26204f1ad80SBarry Smith       for (j=0; j<bs; j++) {
26304f1ad80SBarry Smith         uglyrmapo[(r_rmapo[i]-1)*bs+j] = i*bs + j;
26404f1ad80SBarry Smith       }
26504f1ad80SBarry Smith     }
26604f1ad80SBarry Smith   }
2679566063dSJacob Faibussowitsch   PetscCall(PetscFree(r_rmapo));
2689566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF,nt*bs,&uglyoo));
26904f1ad80SBarry Smith   PetscFunctionReturn(0);
27004f1ad80SBarry Smith }
27104f1ad80SBarry Smith 
2727087cfbeSBarry Smith PetscErrorCode  MatMPIBAIJDiagonalScaleLocal(Mat A,Vec scale)
27304f1ad80SBarry Smith {
27492b32695SKris Buschelman   /* This routine should really be abandoned as it duplicates MatDiagonalScaleLocal */
27592b32695SKris Buschelman 
27692b32695SKris Buschelman   PetscFunctionBegin;
277*cac4c232SBarry Smith   PetscTryMethod(A,"MatDiagonalScaleLocal_C",(Mat,Vec),(A,scale));
27892b32695SKris Buschelman   PetscFunctionReturn(0);
27992b32695SKris Buschelman }
28092b32695SKris Buschelman 
2817087cfbeSBarry Smith PetscErrorCode  MatDiagonalScaleLocal_MPIBAIJ(Mat A,Vec scale)
28292b32695SKris Buschelman {
28304f1ad80SBarry Smith   Mat_MPIBAIJ       *a = (Mat_MPIBAIJ*) A->data; /*access private part of matrix */
284b24ad042SBarry Smith   PetscInt          n,i;
285bca11509SBarry Smith   PetscScalar       *d,*o;
286bca11509SBarry Smith   const PetscScalar *s;
28704f1ad80SBarry Smith 
28804f1ad80SBarry Smith   PetscFunctionBegin;
2892cd6534aSBarry Smith   if (!uglyrmapd) {
2909566063dSJacob Faibussowitsch     PetscCall(MatMPIBAIJDiagonalScaleLocalSetUp(A,scale));
2912cd6534aSBarry Smith   }
2922cd6534aSBarry Smith 
2939566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(scale,&s));
29404f1ad80SBarry Smith 
2959566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(uglydd,&n));
2969566063dSJacob Faibussowitsch   PetscCall(VecGetArray(uglydd,&d));
29704f1ad80SBarry Smith   for (i=0; i<n; i++) {
29804f1ad80SBarry Smith     d[i] = s[uglyrmapd[i]]; /* copy "diagonal" (true local) portion of scale into dd vector */
29904f1ad80SBarry Smith   }
3009566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(uglydd,&d));
30104f1ad80SBarry Smith   /* column scale "diagonal" portion of local matrix */
3029566063dSJacob Faibussowitsch   PetscCall(MatDiagonalScale(a->A,NULL,uglydd));
30304f1ad80SBarry Smith 
3049566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(uglyoo,&n));
3059566063dSJacob Faibussowitsch   PetscCall(VecGetArray(uglyoo,&o));
30604f1ad80SBarry Smith   for (i=0; i<n; i++) {
30704f1ad80SBarry Smith     o[i] = s[uglyrmapo[i]]; /* copy "off-diagonal" portion of scale into oo vector */
30804f1ad80SBarry Smith   }
3099566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(scale,&s));
3109566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(uglyoo,&o));
31104f1ad80SBarry Smith   /* column scale "off-diagonal" portion of local matrix */
3129566063dSJacob Faibussowitsch   PetscCall(MatDiagonalScale(a->B,NULL,uglyoo));
31304f1ad80SBarry Smith   PetscFunctionReturn(0);
31404f1ad80SBarry Smith }
315