xref: /petsc/src/mat/impls/aij/mpi/mmaij.c (revision cac4c232dc4f93991e342196e27ef7b0655dac7b)
1be1d678aSKris Buschelman 
28c79f6d3SBarry Smith /*
38c79f6d3SBarry Smith    Support for the parallel AIJ matrix vector multiply
48c79f6d3SBarry Smith */
5c6db04a5SJed Brown #include <../src/mat/impls/aij/mpi/mpiaij.h>
6aeda0f58SHong Zhang #include <petsc/private/vecimpl.h>
7af0996ceSBarry Smith #include <petsc/private/isimpl.h>    /* needed because accesses data structure of ISLocalToGlobalMapping directly */
88c79f6d3SBarry Smith 
9dfbe8321SBarry Smith PetscErrorCode MatSetUpMultiply_MPIAIJ(Mat mat)
108c79f6d3SBarry Smith {
1144a69424SLois Curfman McInnes   Mat_MPIAIJ         *aij = (Mat_MPIAIJ*)mat->data;
12ec8511deSBarry Smith   Mat_SeqAIJ         *B   = (Mat_SeqAIJ*)(aij->B->data);
13b3c64f9dSJunchao Zhang   PetscInt           i,j,*aj = B->j,*garray;
14b3c64f9dSJunchao Zhang   PetscInt           ec = 0; /* Number of nonzero external columns */
151eb62cbbSBarry Smith   IS                 from,to;
161eb62cbbSBarry Smith   Vec                gvec;
17aa482453SBarry Smith #if defined(PETSC_USE_CTABLE)
180f5bd95cSBarry Smith   PetscTable         gid1_lid1;
190f5bd95cSBarry Smith   PetscTablePosition tpos;
20b1d57f15SBarry Smith   PetscInt           gid,lid;
216f531f54SSatish Balay #else
22d0f46423SBarry Smith   PetscInt           N = mat->cmap->N,*indices;
232066d6f7SSatish Balay #endif
242066d6f7SSatish Balay 
253a40ed3dSBarry Smith   PetscFunctionBegin;
264b8d542aSHong Zhang   if (!aij->garray) {
2728b400f6SJacob Faibussowitsch     PetscCheck(aij->B,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing B mat");
28aa482453SBarry Smith #if defined(PETSC_USE_CTABLE)
29c5bfad50SMark F. Adams     /* use a table */
309566063dSJacob Faibussowitsch     PetscCall(PetscTableCreate(aij->B->rmap->n,mat->cmap->N+1,&gid1_lid1));
31d0f46423SBarry Smith     for (i=0; i<aij->B->rmap->n; i++) {
322066d6f7SSatish Balay       for (j=0; j<B->ilen[i]; j++) {
33b1d57f15SBarry Smith         PetscInt data,gid1 = aj[B->i[i] + j] + 1;
349566063dSJacob Faibussowitsch         PetscCall(PetscTableFind(gid1_lid1,gid1,&data));
35fa46199cSSatish Balay         if (!data) {
362066d6f7SSatish Balay           /* one based table */
379566063dSJacob Faibussowitsch           PetscCall(PetscTableAdd(gid1_lid1,gid1,++ec,INSERT_VALUES));
382066d6f7SSatish Balay         }
392066d6f7SSatish Balay       }
402066d6f7SSatish Balay     }
412066d6f7SSatish Balay     /* form array of columns we need */
429566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(ec,&garray));
439566063dSJacob Faibussowitsch     PetscCall(PetscTableGetHeadPosition(gid1_lid1,&tpos));
442066d6f7SSatish Balay     while (tpos) {
459566063dSJacob Faibussowitsch       PetscCall(PetscTableGetNext(gid1_lid1,&tpos,&gid,&lid));
46b0a32e0cSBarry Smith       gid--;
47b0a32e0cSBarry Smith       lid--;
482066d6f7SSatish Balay       garray[lid] = gid;
492066d6f7SSatish Balay     }
509566063dSJacob Faibussowitsch     PetscCall(PetscSortInt(ec,garray)); /* sort, and rebuild */
519566063dSJacob Faibussowitsch     PetscCall(PetscTableRemoveAll(gid1_lid1));
522066d6f7SSatish Balay     for (i=0; i<ec; i++) {
539566063dSJacob Faibussowitsch       PetscCall(PetscTableAdd(gid1_lid1,garray[i]+1,i+1,INSERT_VALUES));
542066d6f7SSatish Balay     }
552066d6f7SSatish Balay     /* compact out the extra columns in B */
56d0f46423SBarry Smith     for (i=0; i<aij->B->rmap->n; i++) {
572066d6f7SSatish Balay       for (j=0; j<B->ilen[i]; j++) {
58b1d57f15SBarry Smith         PetscInt gid1 = aj[B->i[i] + j] + 1;
599566063dSJacob Faibussowitsch         PetscCall(PetscTableFind(gid1_lid1,gid1,&lid));
60fa46199cSSatish Balay         lid--;
61b3fb0a6cSMatthew Knepley         aj[B->i[i] + j] = lid;
622066d6f7SSatish Balay       }
632066d6f7SSatish Balay     }
649566063dSJacob Faibussowitsch     PetscCall(PetscLayoutDestroy(&aij->B->cmap));
659566063dSJacob Faibussowitsch     PetscCall(PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)aij->B),ec,ec,1,&aij->B->cmap));
669566063dSJacob Faibussowitsch     PetscCall(PetscTableDestroy(&gid1_lid1));
672066d6f7SSatish Balay #else
6811285404SBarry Smith     /* Make an array as long as the number of columns */
691eb62cbbSBarry Smith     /* mark those columns that are in aij->B */
709566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(N,&indices));
71d0f46423SBarry Smith     for (i=0; i<aij->B->rmap->n; i++) {
72d6dfbf8fSBarry Smith       for (j=0; j<B->ilen[i]; j++) {
73bfec09a0SHong Zhang         if (!indices[aj[B->i[i] + j]]) ec++;
74bfec09a0SHong Zhang         indices[aj[B->i[i] + j]] = 1;
75416022c9SBarry Smith       }
761eb62cbbSBarry Smith     }
778c79f6d3SBarry Smith 
781eb62cbbSBarry Smith     /* form array of columns we need */
799566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(ec,&garray));
801eb62cbbSBarry Smith     ec   = 0;
811eb62cbbSBarry Smith     for (i=0; i<N; i++) {
821eb62cbbSBarry Smith       if (indices[i]) garray[ec++] = i;
831eb62cbbSBarry Smith     }
841eb62cbbSBarry Smith 
851eb62cbbSBarry Smith     /* make indices now point into garray */
861eb62cbbSBarry Smith     for (i=0; i<ec; i++) {
87bfec09a0SHong Zhang       indices[garray[i]] = i;
881eb62cbbSBarry Smith     }
891eb62cbbSBarry Smith 
901eb62cbbSBarry Smith     /* compact out the extra columns in B */
91d0f46423SBarry Smith     for (i=0; i<aij->B->rmap->n; i++) {
92d6dfbf8fSBarry Smith       for (j=0; j<B->ilen[i]; j++) {
93bfec09a0SHong Zhang         aj[B->i[i] + j] = indices[aj[B->i[i] + j]];
941eb62cbbSBarry Smith       }
95d6dfbf8fSBarry Smith     }
969566063dSJacob Faibussowitsch     PetscCall(PetscLayoutDestroy(&aij->B->cmap));
979566063dSJacob Faibussowitsch     PetscCall(PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)aij->B),ec,ec,1,&aij->B->cmap));
989566063dSJacob Faibussowitsch     PetscCall(PetscFree(indices));
992066d6f7SSatish Balay #endif
1004b8d542aSHong Zhang   } else {
1014b8d542aSHong Zhang     garray = aij->garray;
1024b8d542aSHong Zhang   }
1034b8d542aSHong Zhang 
1044b8d542aSHong Zhang   if (!aij->lvec) {
10528b400f6SJacob Faibussowitsch     PetscCheck(aij->B,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing B mat");
1069566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(aij->B,&aij->lvec,NULL));
1074b8d542aSHong Zhang   }
1089566063dSJacob Faibussowitsch   PetscCall(VecGetSize(aij->lvec,&ec));
1091eb62cbbSBarry Smith 
110d6dfbf8fSBarry Smith   /* create two temporary Index sets for build scatter gather */
1119566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF,ec,garray,PETSC_COPY_VALUES,&from));
1129566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(PETSC_COMM_SELF,ec,0,1,&to));
1131eb62cbbSBarry Smith 
1141eb62cbbSBarry Smith   /* create temporary global vector to generate scatter context */
115b5eb4454SBarry Smith   /* This does not allocate the array's memory so is efficient */
1169566063dSJacob Faibussowitsch   PetscCall(VecCreateMPIWithArray(PetscObjectComm((PetscObject)mat),1,mat->cmap->n,mat->cmap->N,NULL,&gvec));
1171eb62cbbSBarry Smith 
1182d336d48SLois Curfman McInnes   /* generate the scatter context */
1199566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&aij->Mvctx));
1209566063dSJacob Faibussowitsch   PetscCall(VecScatterCreate(gvec,from,aij->lvec,to,&aij->Mvctx));
1219566063dSJacob Faibussowitsch   PetscCall(VecScatterViewFromOptions(aij->Mvctx,(PetscObject)mat,"-matmult_vecscatter_view"));
1229566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectParent((PetscObject)mat,(PetscObject)aij->Mvctx));
1239566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectParent((PetscObject)mat,(PetscObject)aij->lvec));
1249566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectMemory((PetscObject)mat,ec*sizeof(PetscInt)));
12567bb5161SHong Zhang   aij->garray = garray;
1264b8d542aSHong Zhang 
1279566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectParent((PetscObject)mat,(PetscObject)from));
1289566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectParent((PetscObject)mat,(PetscObject)to));
1292205254eSKarl Rupp 
1309566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&from));
1319566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&to));
1329566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&gvec));
1333a40ed3dSBarry Smith   PetscFunctionReturn(0);
1348c79f6d3SBarry Smith }
1359e25ed09SBarry Smith 
136fff043a9SJunchao Zhang /* Disassemble the off-diag portion of the MPIAIJXxx matrix.
137fff043a9SJunchao Zhang    In other words, change the B from reduced format using local col ids
138fff043a9SJunchao Zhang    to expanded format using global col ids, which will make it easier to
139fff043a9SJunchao Zhang    insert new nonzeros (with global col ids) into the matrix.
140fff043a9SJunchao Zhang    The off-diag B determines communication in the matrix vector multiply.
1412493cbb0SBarry Smith */
142ab9863d7SBarry Smith PetscErrorCode MatDisAssemble_MPIAIJ(Mat A)
1432493cbb0SBarry Smith {
1442493cbb0SBarry Smith   Mat_MPIAIJ        *aij  = (Mat_MPIAIJ*)A->data;
1452493cbb0SBarry Smith   Mat               B     = aij->B,Bnew;
146ec8511deSBarry Smith   Mat_SeqAIJ        *Baij = (Mat_SeqAIJ*)B->data;
147d0f46423SBarry Smith   PetscInt          i,j,m = B->rmap->n,n = A->cmap->N,col,ct = 0,*garray = aij->garray,*nz,ec;
14887828ca2SBarry Smith   PetscScalar       v;
149fff043a9SJunchao Zhang   const PetscScalar *ba;
1502493cbb0SBarry Smith 
1513a40ed3dSBarry Smith   PetscFunctionBegin;
1522493cbb0SBarry Smith   /* free stuff related to matrix-vec multiply */
1539566063dSJacob Faibussowitsch   PetscCall(VecGetSize(aij->lvec,&ec)); /* needed for PetscLogObjectMemory below */
1549566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&aij->lvec));
155464493b3SBarry Smith   if (aij->colmap) {
156aa482453SBarry Smith #if defined(PETSC_USE_CTABLE)
1579566063dSJacob Faibussowitsch     PetscCall(PetscTableDestroy(&aij->colmap));
1582066d6f7SSatish Balay #else
1599566063dSJacob Faibussowitsch     PetscCall(PetscFree(aij->colmap));
1609566063dSJacob Faibussowitsch     PetscCall(PetscLogObjectMemory((PetscObject)A,-aij->B->cmap->n*sizeof(PetscInt)));
1612066d6f7SSatish Balay #endif
162464493b3SBarry Smith   }
1632493cbb0SBarry Smith 
1642493cbb0SBarry Smith   /* 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));
1672493cbb0SBarry Smith 
1682493cbb0SBarry Smith   /* invent new B and copy stuff over */
1699566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(m+1,&nz));
17048b35521SBarry Smith   for (i=0; i<m; i++) {
17148b35521SBarry Smith     nz[i] = Baij->i[i+1] - Baij->i[i];
17248b35521SBarry Smith   }
1739566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF,&Bnew));
1749566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(Bnew,m,n,m,n)); /* Bnew now uses A->cmap->N as its col size */
1759566063dSJacob Faibussowitsch   PetscCall(MatSetBlockSizesFromMats(Bnew,A,A));
1769566063dSJacob Faibussowitsch   PetscCall(MatSetType(Bnew,((PetscObject)B)->type_name));
1779566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(Bnew,0,nz));
1782205254eSKarl Rupp 
179b38c15b3SStefano Zampini   if (Baij->nonew >= 0) { /* Inherit insertion error options (if positive). */
180b38c15b3SStefano Zampini     ((Mat_SeqAIJ*)Bnew->data)->nonew = Baij->nonew;
181b38c15b3SStefano Zampini   }
182b38c15b3SStefano Zampini 
18377341eacSDmitry Karpeev   /*
18477341eacSDmitry Karpeev    Ensure that B's nonzerostate is monotonically increasing.
18577341eacSDmitry Karpeev    Or should this follow the MatSetValues() loop to preserve B's nonzerstate across a MatDisAssemble() call?
18677341eacSDmitry Karpeev    */
187f69fde56SShane Stafford   Bnew->nonzerostate = B->nonzerostate;
1882205254eSKarl Rupp 
1899566063dSJacob Faibussowitsch   PetscCall(PetscFree(nz));
1909566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJGetArrayRead(B,&ba));
1912493cbb0SBarry Smith   for (i=0; i<m; i++) {
192bfec09a0SHong Zhang     for (j=Baij->i[i]; j<Baij->i[i+1]; j++) {
193bfec09a0SHong Zhang       col  = garray[Baij->j[ct]];
194fff043a9SJunchao Zhang       v    = ba[ct++];
1959566063dSJacob Faibussowitsch       PetscCall(MatSetValues(Bnew,1,&i,1,&col,&v,B->insertmode));
1962493cbb0SBarry Smith     }
1972493cbb0SBarry Smith   }
1989566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJRestoreArrayRead(B,&ba));
199fff043a9SJunchao Zhang 
2009566063dSJacob Faibussowitsch   PetscCall(PetscFree(aij->garray));
2019566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectMemory((PetscObject)A,-ec*sizeof(PetscInt)));
2029566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
2039566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectParent((PetscObject)A,(PetscObject)Bnew));
2042205254eSKarl Rupp 
2052493cbb0SBarry Smith   aij->B           = Bnew;
206227d817aSBarry Smith   A->was_assembled = PETSC_FALSE;
2073a40ed3dSBarry Smith   PetscFunctionReturn(0);
2082493cbb0SBarry Smith }
2092493cbb0SBarry Smith 
2102cd6534aSBarry Smith /*      ugly stuff added for Glenn someday we should fix this up */
2112cd6534aSBarry Smith 
212f4259b30SLisandro Dalcin static PetscInt *auglyrmapd = NULL,*auglyrmapo = NULL; /* mapping from the local ordering to the "diagonal" and "off-diagonal" parts of the local matrix */
213f4259b30SLisandro Dalcin static Vec auglydd          = NULL,auglyoo     = NULL; /* work vectors used to scale the two parts of the local matrix */
2142cd6534aSBarry Smith 
215dfbe8321SBarry Smith PetscErrorCode MatMPIAIJDiagonalScaleLocalSetUp(Mat inA,Vec scale)
2162cd6534aSBarry Smith {
2172cd6534aSBarry Smith   Mat_MPIAIJ     *ina = (Mat_MPIAIJ*) inA->data; /*access private part of matrix */
218b1d57f15SBarry Smith   PetscInt       i,n,nt,cstart,cend,no,*garray = ina->garray,*lindices;
219b1d57f15SBarry Smith   PetscInt       *r_rmapd,*r_rmapo;
2202cd6534aSBarry Smith 
2212cd6534aSBarry 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));
2252cd6534aSBarry Smith   nt   = 0;
226992144d0SBarry Smith   for (i=0; i<inA->rmap->mapping->n; i++) {
227992144d0SBarry Smith     if (inA->rmap->mapping->indices[i] >= cstart && inA->rmap->mapping->indices[i] < cend) {
2282cd6534aSBarry Smith       nt++;
229992144d0SBarry Smith       r_rmapd[i] = inA->rmap->mapping->indices[i] + 1;
2302cd6534aSBarry Smith     }
2312cd6534aSBarry Smith   }
2322c71b3e2SJacob Faibussowitsch   PetscCheckFalse(nt != n,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hmm nt %" PetscInt_FMT " n %" PetscInt_FMT,nt,n);
2339566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n+1,&auglyrmapd));
234992144d0SBarry Smith   for (i=0; i<inA->rmap->mapping->n; i++) {
2352cd6534aSBarry Smith     if (r_rmapd[i]) {
2362cd6534aSBarry Smith       auglyrmapd[(r_rmapd[i]-1)-cstart] = i;
2372cd6534aSBarry Smith     }
2382cd6534aSBarry Smith   }
2399566063dSJacob Faibussowitsch   PetscCall(PetscFree(r_rmapd));
2409566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF,n,&auglydd));
2412cd6534aSBarry Smith 
2429566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(inA->cmap->N+1,&lindices));
243d0f46423SBarry Smith   for (i=0; i<ina->B->cmap->n; i++) {
2442cd6534aSBarry Smith     lindices[garray[i]] = i+1;
2452cd6534aSBarry Smith   }
246992144d0SBarry Smith   no   = inA->rmap->mapping->n - nt;
2479566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(inA->rmap->mapping->n+1,&r_rmapo));
2482cd6534aSBarry Smith   nt   = 0;
249992144d0SBarry Smith   for (i=0; i<inA->rmap->mapping->n; i++) {
250992144d0SBarry Smith     if (lindices[inA->rmap->mapping->indices[i]]) {
2512cd6534aSBarry Smith       nt++;
252992144d0SBarry Smith       r_rmapo[i] = lindices[inA->rmap->mapping->indices[i]];
2532cd6534aSBarry Smith     }
2542cd6534aSBarry Smith   }
2552c71b3e2SJacob Faibussowitsch   PetscCheckFalse(nt > no,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hmm nt %" PetscInt_FMT " no %" PetscInt_FMT,nt,n);
2569566063dSJacob Faibussowitsch   PetscCall(PetscFree(lindices));
2579566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nt+1,&auglyrmapo));
258992144d0SBarry Smith   for (i=0; i<inA->rmap->mapping->n; i++) {
2592cd6534aSBarry Smith     if (r_rmapo[i]) {
2602cd6534aSBarry Smith       auglyrmapo[(r_rmapo[i]-1)] = i;
2612cd6534aSBarry Smith     }
2622cd6534aSBarry Smith   }
2639566063dSJacob Faibussowitsch   PetscCall(PetscFree(r_rmapo));
2649566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF,nt,&auglyoo));
2652cd6534aSBarry Smith   PetscFunctionReturn(0);
2662cd6534aSBarry Smith }
2672cd6534aSBarry Smith 
268dfbe8321SBarry Smith PetscErrorCode MatMPIAIJDiagonalScaleLocal(Mat A,Vec scale)
2692cd6534aSBarry Smith {
27092b32695SKris Buschelman   /* This routine should really be abandoned as it duplicates MatDiagonalScaleLocal */
27192b32695SKris Buschelman 
27292b32695SKris Buschelman   PetscFunctionBegin;
273*cac4c232SBarry Smith   PetscTryMethod(A,"MatDiagonalScaleLocal_C",(Mat,Vec),(A,scale));
27492b32695SKris Buschelman   PetscFunctionReturn(0);
27592b32695SKris Buschelman }
27692b32695SKris Buschelman 
2777087cfbeSBarry Smith PetscErrorCode  MatDiagonalScaleLocal_MPIAIJ(Mat A,Vec scale)
27892b32695SKris Buschelman {
2792cd6534aSBarry Smith   Mat_MPIAIJ        *a = (Mat_MPIAIJ*) A->data; /*access private part of matrix */
280b1d57f15SBarry Smith   PetscInt          n,i;
281bca11509SBarry Smith   PetscScalar       *d,*o;
282bca11509SBarry Smith   const PetscScalar *s;
2832cd6534aSBarry Smith 
2842cd6534aSBarry Smith   PetscFunctionBegin;
2852cd6534aSBarry Smith   if (!auglyrmapd) {
2869566063dSJacob Faibussowitsch     PetscCall(MatMPIAIJDiagonalScaleLocalSetUp(A,scale));
2872cd6534aSBarry Smith   }
2882cd6534aSBarry Smith 
2899566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(scale,&s));
2902cd6534aSBarry Smith 
2919566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(auglydd,&n));
2929566063dSJacob Faibussowitsch   PetscCall(VecGetArray(auglydd,&d));
2932cd6534aSBarry Smith   for (i=0; i<n; i++) {
2942cd6534aSBarry Smith     d[i] = s[auglyrmapd[i]]; /* copy "diagonal" (true local) portion of scale into dd vector */
2952cd6534aSBarry Smith   }
2969566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(auglydd,&d));
2972cd6534aSBarry Smith   /* column scale "diagonal" portion of local matrix */
2989566063dSJacob Faibussowitsch   PetscCall(MatDiagonalScale(a->A,NULL,auglydd));
2992cd6534aSBarry Smith 
3009566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(auglyoo,&n));
3019566063dSJacob Faibussowitsch   PetscCall(VecGetArray(auglyoo,&o));
3022cd6534aSBarry Smith   for (i=0; i<n; i++) {
3032cd6534aSBarry Smith     o[i] = s[auglyrmapo[i]]; /* copy "off-diagonal" portion of scale into oo vector */
3042cd6534aSBarry Smith   }
3059566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(scale,&s));
3069566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(auglyoo,&o));
3072cd6534aSBarry Smith   /* column scale "off-diagonal" portion of local matrix */
3089566063dSJacob Faibussowitsch   PetscCall(MatDiagonalScale(a->B,NULL,auglyoo));
3092cd6534aSBarry Smith   PetscFunctionReturn(0);
3102cd6534aSBarry Smith }
311