xref: /petsc/src/mat/impls/aij/mpi/mmaij.c (revision 1795a4d16c893ec2fc06bbbc6c5ce592a2de75d4)
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>
6186d4ecdSBarry Smith #include <petsc-private/isimpl.h>    /* needed because accesses data structure of ISLocalToGlobalMapping directly */
78c79f6d3SBarry Smith 
84a2ae208SSatish Balay #undef __FUNCT__
94a2ae208SSatish Balay #define __FUNCT__ "MatSetUpMultiply_MPIAIJ"
10dfbe8321SBarry Smith PetscErrorCode MatSetUpMultiply_MPIAIJ(Mat mat)
118c79f6d3SBarry Smith {
1244a69424SLois Curfman McInnes   Mat_MPIAIJ     *aij = (Mat_MPIAIJ*)mat->data;
13ec8511deSBarry Smith   Mat_SeqAIJ     *B   = (Mat_SeqAIJ*)(aij->B->data);
146849ba73SBarry Smith   PetscErrorCode ierr;
15b1d57f15SBarry Smith   PetscInt       i,j,*aj = B->j,ec = 0,*garray;
161eb62cbbSBarry Smith   IS             from,to;
171eb62cbbSBarry Smith   Vec            gvec;
18aa482453SBarry Smith #if defined(PETSC_USE_CTABLE)
190f5bd95cSBarry Smith   PetscTable         gid1_lid1;
200f5bd95cSBarry Smith   PetscTablePosition tpos;
21b1d57f15SBarry Smith   PetscInt           gid,lid;
226f531f54SSatish Balay #else
23d0f46423SBarry Smith   PetscInt N = mat->cmap->N,*indices;
242066d6f7SSatish Balay #endif
252066d6f7SSatish Balay 
263a40ed3dSBarry Smith   PetscFunctionBegin;
27aa482453SBarry Smith #if defined(PETSC_USE_CTABLE)
28c5bfad50SMark F. Adams   /* use a table */
29e23dfa41SBarry Smith   ierr = PetscTableCreate(aij->B->rmap->n,mat->cmap->N+1,&gid1_lid1);CHKERRQ(ierr);
30d0f46423SBarry Smith   for (i=0; i<aij->B->rmap->n; i++) {
312066d6f7SSatish Balay     for (j=0; j<B->ilen[i]; j++) {
32b1d57f15SBarry Smith       PetscInt data,gid1 = aj[B->i[i] + j] + 1;
330f5bd95cSBarry Smith       ierr = PetscTableFind(gid1_lid1,gid1,&data);CHKERRQ(ierr);
34fa46199cSSatish Balay       if (!data) {
352066d6f7SSatish Balay         /* one based table */
363861aac3SJed Brown         ierr = PetscTableAdd(gid1_lid1,gid1,++ec,INSERT_VALUES);CHKERRQ(ierr);
372066d6f7SSatish Balay       }
382066d6f7SSatish Balay     }
392066d6f7SSatish Balay   }
402066d6f7SSatish Balay   /* form array of columns we need */
41785e854fSJed Brown   ierr = PetscMalloc1((ec+1),&garray);CHKERRQ(ierr);
420f5bd95cSBarry Smith   ierr = PetscTableGetHeadPosition(gid1_lid1,&tpos);CHKERRQ(ierr);
432066d6f7SSatish Balay   while (tpos) {
440f5bd95cSBarry Smith     ierr = PetscTableGetNext(gid1_lid1,&tpos,&gid,&lid);CHKERRQ(ierr);
45b0a32e0cSBarry Smith     gid--;
46b0a32e0cSBarry Smith     lid--;
472066d6f7SSatish Balay     garray[lid] = gid;
482066d6f7SSatish Balay   }
490064e2bbSSatish Balay   ierr = PetscSortInt(ec,garray);CHKERRQ(ierr); /* sort, and rebuild */
500f5bd95cSBarry Smith   ierr = PetscTableRemoveAll(gid1_lid1);CHKERRQ(ierr);
512066d6f7SSatish Balay   for (i=0; i<ec; i++) {
523861aac3SJed Brown     ierr = PetscTableAdd(gid1_lid1,garray[i]+1,i+1,INSERT_VALUES);CHKERRQ(ierr);
532066d6f7SSatish Balay   }
542066d6f7SSatish Balay   /* compact out the extra columns in B */
55d0f46423SBarry Smith   for (i=0; i<aij->B->rmap->n; i++) {
562066d6f7SSatish Balay     for (j=0; j<B->ilen[i]; j++) {
57b1d57f15SBarry Smith       PetscInt gid1 = aj[B->i[i] + j] + 1;
580f5bd95cSBarry Smith       ierr = PetscTableFind(gid1_lid1,gid1,&lid);CHKERRQ(ierr);
59fa46199cSSatish Balay       lid--;
60b3fb0a6cSMatthew Knepley       aj[B->i[i] + j] = lid;
612066d6f7SSatish Balay     }
622066d6f7SSatish Balay   }
63d0f46423SBarry Smith   aij->B->cmap->n = aij->B->cmap->N = ec;
642205254eSKarl Rupp 
6526283091SBarry Smith   ierr = PetscLayoutSetUp((aij->B->cmap));CHKERRQ(ierr);
666bc0bbbfSBarry Smith   ierr = PetscTableDestroy(&gid1_lid1);CHKERRQ(ierr);
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 */
70*1795a4d1SJed Brown   ierr = PetscCalloc1(N+1,&indices);CHKERRQ(ierr);
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 */
79785e854fSJed Brown   ierr = PetscMalloc1((ec+1),&garray);CHKERRQ(ierr);
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   }
96d0f46423SBarry Smith   aij->B->cmap->n = aij->B->cmap->N = ec;
972205254eSKarl Rupp 
9826283091SBarry Smith   ierr = PetscLayoutSetUp((aij->B->cmap));CHKERRQ(ierr);
99606d414cSSatish Balay   ierr = PetscFree(indices);CHKERRQ(ierr);
1002066d6f7SSatish Balay #endif
1011eb62cbbSBarry Smith   /* create local vector that is used to scatter into */
102029af93fSBarry Smith   ierr = VecCreateSeq(PETSC_COMM_SELF,ec,&aij->lvec);CHKERRQ(ierr);
1031eb62cbbSBarry Smith 
104d6dfbf8fSBarry Smith   /* create two temporary Index sets for build scatter gather */
10570b3c8c7SBarry Smith   ierr = ISCreateGeneral(((PetscObject)mat)->comm,ec,garray,PETSC_COPY_VALUES,&from);CHKERRQ(ierr);
106c5bfad50SMark F. Adams 
107029af93fSBarry Smith   ierr = ISCreateStride(PETSC_COMM_SELF,ec,0,1,&to);CHKERRQ(ierr);
1081eb62cbbSBarry Smith 
1091eb62cbbSBarry Smith   /* create temporary global vector to generate scatter context */
110b5eb4454SBarry Smith   /* This does not allocate the array's memory so is efficient */
111ce94432eSBarry Smith   ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)mat),1,mat->cmap->n,mat->cmap->N,NULL,&gvec);CHKERRQ(ierr);
1121eb62cbbSBarry Smith 
1132d336d48SLois Curfman McInnes   /* generate the scatter context */
11408480c60SBarry Smith   ierr = VecScatterCreate(gvec,from,aij->lvec,to,&aij->Mvctx);CHKERRQ(ierr);
1153bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)aij->Mvctx);CHKERRQ(ierr);
1163bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)aij->lvec);CHKERRQ(ierr);
1173bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)from);CHKERRQ(ierr);
1183bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)to);CHKERRQ(ierr);
1192205254eSKarl Rupp 
1209e25ed09SBarry Smith   aij->garray = garray;
1212205254eSKarl Rupp 
1223bb1ff40SBarry Smith   ierr = PetscLogObjectMemory((PetscObject)mat,(ec+1)*sizeof(PetscInt));CHKERRQ(ierr);
1236bf464f9SBarry Smith   ierr = ISDestroy(&from);CHKERRQ(ierr);
1246bf464f9SBarry Smith   ierr = ISDestroy(&to);CHKERRQ(ierr);
1256bf464f9SBarry Smith   ierr = VecDestroy(&gvec);CHKERRQ(ierr);
1263a40ed3dSBarry Smith   PetscFunctionReturn(0);
1278c79f6d3SBarry Smith }
1289e25ed09SBarry Smith 
1299e25ed09SBarry Smith 
1304a2ae208SSatish Balay #undef __FUNCT__
131ab9863d7SBarry Smith #define __FUNCT__ "MatDisAssemble_MPIAIJ"
1322493cbb0SBarry Smith /*
1332493cbb0SBarry Smith      Takes the local part of an already assembled MPIAIJ matrix
1342493cbb0SBarry Smith    and disassembles it. This is to allow new nonzeros into the matrix
1352493cbb0SBarry Smith    that require more communication in the matrix vector multiply.
1362493cbb0SBarry Smith    Thus certain data-structures must be rebuilt.
1372493cbb0SBarry Smith 
1382493cbb0SBarry Smith    Kind of slow! But that's what application programmers get when
1392493cbb0SBarry Smith    they are sloppy.
1402493cbb0SBarry Smith */
141ab9863d7SBarry Smith PetscErrorCode MatDisAssemble_MPIAIJ(Mat A)
1422493cbb0SBarry Smith {
1432493cbb0SBarry Smith   Mat_MPIAIJ     *aij  = (Mat_MPIAIJ*)A->data;
1442493cbb0SBarry Smith   Mat            B     = aij->B,Bnew;
145ec8511deSBarry Smith   Mat_SeqAIJ     *Baij = (Mat_SeqAIJ*)B->data;
146dfbe8321SBarry Smith   PetscErrorCode ierr;
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;
1492493cbb0SBarry Smith 
1503a40ed3dSBarry Smith   PetscFunctionBegin;
1512493cbb0SBarry Smith   /* free stuff related to matrix-vec multiply */
152b0a32e0cSBarry Smith   ierr = VecGetSize(aij->lvec,&ec);CHKERRQ(ierr); /* needed for PetscLogObjectMemory below */
1535e1f6667SBarry Smith   ierr = VecDestroy(&aij->lvec);CHKERRQ(ierr);
1545e1f6667SBarry Smith   ierr = VecScatterDestroy(&aij->Mvctx);CHKERRQ(ierr);
155464493b3SBarry Smith   if (aij->colmap) {
156aa482453SBarry Smith #if defined(PETSC_USE_CTABLE)
1576bc0bbbfSBarry Smith     ierr = PetscTableDestroy(&aij->colmap);CHKERRQ(ierr);
1582066d6f7SSatish Balay #else
159606d414cSSatish Balay     ierr = PetscFree(aij->colmap);CHKERRQ(ierr);
1603bb1ff40SBarry Smith     ierr = PetscLogObjectMemory((PetscObject)A,-aij->B->cmap->n*sizeof(PetscInt));CHKERRQ(ierr);
1612066d6f7SSatish Balay #endif
162464493b3SBarry Smith   }
1632493cbb0SBarry Smith 
1642493cbb0SBarry Smith   /* make sure that B is assembled so we can access its values */
1656d4a8577SBarry Smith   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
166fe2f2677SBarry Smith   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1672493cbb0SBarry Smith 
1682493cbb0SBarry Smith   /* invent new B and copy stuff over */
169785e854fSJed Brown   ierr = PetscMalloc1((m+1),&nz);CHKERRQ(ierr);
17048b35521SBarry Smith   for (i=0; i<m; i++) {
17148b35521SBarry Smith     nz[i] = Baij->i[i+1] - Baij->i[i];
17248b35521SBarry Smith   }
173f69a0ea3SMatthew Knepley   ierr = MatCreate(PETSC_COMM_SELF,&Bnew);CHKERRQ(ierr);
174f69a0ea3SMatthew Knepley   ierr = MatSetSizes(Bnew,m,n,m,n);CHKERRQ(ierr);
175a2f3521dSMark F. Adams   ierr = MatSetBlockSizes(Bnew,A->rmap->bs,A->cmap->bs);CHKERRQ(ierr);
1767adad957SLisandro Dalcin   ierr = MatSetType(Bnew,((PetscObject)B)->type_name);CHKERRQ(ierr);
177f204ca49SKris Buschelman   ierr = MatSeqAIJSetPreallocation(Bnew,0,nz);CHKERRQ(ierr);
1782205254eSKarl Rupp 
1792576faa2SJed Brown   ((Mat_SeqAIJ*)Bnew->data)->nonew = Baij->nonew; /* Inherit insertion error options. */
1802205254eSKarl Rupp 
181606d414cSSatish Balay   ierr = PetscFree(nz);CHKERRQ(ierr);
1822493cbb0SBarry Smith   for (i=0; i<m; i++) {
183bfec09a0SHong Zhang     for (j=Baij->i[i]; j<Baij->i[i+1]; j++) {
184bfec09a0SHong Zhang       col  = garray[Baij->j[ct]];
1852493cbb0SBarry Smith       v    = Baij->a[ct++];
18683271157SBarry Smith       ierr = MatSetValues(Bnew,1,&i,1,&col,&v,B->insertmode);CHKERRQ(ierr);
1872493cbb0SBarry Smith     }
1882493cbb0SBarry Smith   }
189606d414cSSatish Balay   ierr = PetscFree(aij->garray);CHKERRQ(ierr);
1903bb1ff40SBarry Smith   ierr = PetscLogObjectMemory((PetscObject)A,-ec*sizeof(PetscInt));CHKERRQ(ierr);
1916bf464f9SBarry Smith   ierr = MatDestroy(&B);CHKERRQ(ierr);
1923bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)A,(PetscObject)Bnew);CHKERRQ(ierr);
1932205254eSKarl Rupp 
1942493cbb0SBarry Smith   aij->B           = Bnew;
195227d817aSBarry Smith   A->was_assembled = PETSC_FALSE;
1963a40ed3dSBarry Smith   PetscFunctionReturn(0);
1972493cbb0SBarry Smith }
1982493cbb0SBarry Smith 
1992cd6534aSBarry Smith /*      ugly stuff added for Glenn someday we should fix this up */
2002cd6534aSBarry Smith 
2012205254eSKarl Rupp static PetscInt *auglyrmapd = 0,*auglyrmapo = 0; /* mapping from the local ordering to the "diagonal" and "off-diagonal" parts of the local matrix */
2022cd6534aSBarry Smith static Vec auglydd          = 0,auglyoo     = 0; /* work vectors used to scale the two parts of the local matrix */
2032cd6534aSBarry Smith 
2042cd6534aSBarry Smith 
2052cd6534aSBarry Smith #undef __FUNCT__
2062cd6534aSBarry Smith #define __FUNCT__ "MatMPIAIJDiagonalScaleLocalSetUp"
207dfbe8321SBarry Smith PetscErrorCode MatMPIAIJDiagonalScaleLocalSetUp(Mat inA,Vec scale)
2082cd6534aSBarry Smith {
2092cd6534aSBarry Smith   Mat_MPIAIJ     *ina = (Mat_MPIAIJ*) inA->data; /*access private part of matrix */
210dfbe8321SBarry Smith   PetscErrorCode ierr;
211b1d57f15SBarry Smith   PetscInt       i,n,nt,cstart,cend,no,*garray = ina->garray,*lindices;
212b1d57f15SBarry Smith   PetscInt       *r_rmapd,*r_rmapo;
2132cd6534aSBarry Smith 
2142cd6534aSBarry Smith   PetscFunctionBegin;
2152cd6534aSBarry Smith   ierr = MatGetOwnershipRange(inA,&cstart,&cend);CHKERRQ(ierr);
2160298fd71SBarry Smith   ierr = MatGetSize(ina->A,NULL,&n);CHKERRQ(ierr);
217*1795a4d1SJed Brown   ierr = PetscCalloc1(inA->rmap->mapping->n+1,&r_rmapd);CHKERRQ(ierr);
2182cd6534aSBarry Smith   nt   = 0;
219992144d0SBarry Smith   for (i=0; i<inA->rmap->mapping->n; i++) {
220992144d0SBarry Smith     if (inA->rmap->mapping->indices[i] >= cstart && inA->rmap->mapping->indices[i] < cend) {
2212cd6534aSBarry Smith       nt++;
222992144d0SBarry Smith       r_rmapd[i] = inA->rmap->mapping->indices[i] + 1;
2232cd6534aSBarry Smith     }
2242cd6534aSBarry Smith   }
225e32f2f54SBarry Smith   if (nt != n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hmm nt %D n %D",nt,n);
226785e854fSJed Brown   ierr = PetscMalloc1((n+1),&auglyrmapd);CHKERRQ(ierr);
227992144d0SBarry Smith   for (i=0; i<inA->rmap->mapping->n; i++) {
2282cd6534aSBarry Smith     if (r_rmapd[i]) {
2292cd6534aSBarry Smith       auglyrmapd[(r_rmapd[i]-1)-cstart] = i;
2302cd6534aSBarry Smith     }
2312cd6534aSBarry Smith   }
2322cd6534aSBarry Smith   ierr = PetscFree(r_rmapd);CHKERRQ(ierr);
2332cd6534aSBarry Smith   ierr = VecCreateSeq(PETSC_COMM_SELF,n,&auglydd);CHKERRQ(ierr);
2342cd6534aSBarry Smith 
235*1795a4d1SJed Brown   ierr = PetscCalloc1(inA->cmap->N+1,&lindices);CHKERRQ(ierr);
236d0f46423SBarry Smith   for (i=0; i<ina->B->cmap->n; i++) {
2372cd6534aSBarry Smith     lindices[garray[i]] = i+1;
2382cd6534aSBarry Smith   }
239992144d0SBarry Smith   no   = inA->rmap->mapping->n - nt;
240*1795a4d1SJed Brown   ierr = PetscCalloc1(inA->rmap->mapping->n+1,&r_rmapo);CHKERRQ(ierr);
2412cd6534aSBarry Smith   nt   = 0;
242992144d0SBarry Smith   for (i=0; i<inA->rmap->mapping->n; i++) {
243992144d0SBarry Smith     if (lindices[inA->rmap->mapping->indices[i]]) {
2442cd6534aSBarry Smith       nt++;
245992144d0SBarry Smith       r_rmapo[i] = lindices[inA->rmap->mapping->indices[i]];
2462cd6534aSBarry Smith     }
2472cd6534aSBarry Smith   }
248e32f2f54SBarry Smith   if (nt > no) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hmm nt %D no %D",nt,n);
2492cd6534aSBarry Smith   ierr = PetscFree(lindices);CHKERRQ(ierr);
250785e854fSJed Brown   ierr = PetscMalloc1((nt+1),&auglyrmapo);CHKERRQ(ierr);
251992144d0SBarry Smith   for (i=0; i<inA->rmap->mapping->n; i++) {
2522cd6534aSBarry Smith     if (r_rmapo[i]) {
2532cd6534aSBarry Smith       auglyrmapo[(r_rmapo[i]-1)] = i;
2542cd6534aSBarry Smith     }
2552cd6534aSBarry Smith   }
2562cd6534aSBarry Smith   ierr = PetscFree(r_rmapo);CHKERRQ(ierr);
2572cd6534aSBarry Smith   ierr = VecCreateSeq(PETSC_COMM_SELF,nt,&auglyoo);CHKERRQ(ierr);
2582cd6534aSBarry Smith   PetscFunctionReturn(0);
2592cd6534aSBarry Smith }
2602cd6534aSBarry Smith 
2612cd6534aSBarry Smith #undef __FUNCT__
2622cd6534aSBarry Smith #define __FUNCT__ "MatMPIAIJDiagonalScaleLocal"
263dfbe8321SBarry Smith PetscErrorCode MatMPIAIJDiagonalScaleLocal(Mat A,Vec scale)
2642cd6534aSBarry Smith {
26592b32695SKris Buschelman   /* This routine should really be abandoned as it duplicates MatDiagonalScaleLocal */
2664ac538c5SBarry Smith   PetscErrorCode ierr;
26792b32695SKris Buschelman 
26892b32695SKris Buschelman   PetscFunctionBegin;
2694ac538c5SBarry Smith   ierr = PetscTryMethod(A,"MatDiagonalScaleLocal_C",(Mat,Vec),(A,scale));CHKERRQ(ierr);
27092b32695SKris Buschelman   PetscFunctionReturn(0);
27192b32695SKris Buschelman }
27292b32695SKris Buschelman 
27392b32695SKris Buschelman #undef __FUNCT__
27492b32695SKris Buschelman #define __FUNCT__ "MatDiagonalScaleLocal_MPIAIJ"
2757087cfbeSBarry Smith PetscErrorCode  MatDiagonalScaleLocal_MPIAIJ(Mat A,Vec scale)
27692b32695SKris Buschelman {
2772cd6534aSBarry Smith   Mat_MPIAIJ     *a = (Mat_MPIAIJ*) A->data; /*access private part of matrix */
278dfbe8321SBarry Smith   PetscErrorCode ierr;
279b1d57f15SBarry Smith   PetscInt       n,i;
2802cd6534aSBarry Smith   PetscScalar    *d,*o,*s;
2812cd6534aSBarry Smith 
2822cd6534aSBarry Smith   PetscFunctionBegin;
2832cd6534aSBarry Smith   if (!auglyrmapd) {
2842cd6534aSBarry Smith     ierr = MatMPIAIJDiagonalScaleLocalSetUp(A,scale);CHKERRQ(ierr);
2852cd6534aSBarry Smith   }
2862cd6534aSBarry Smith 
2871ebc52fbSHong Zhang   ierr = VecGetArray(scale,&s);CHKERRQ(ierr);
2882cd6534aSBarry Smith 
2892cd6534aSBarry Smith   ierr = VecGetLocalSize(auglydd,&n);CHKERRQ(ierr);
2901ebc52fbSHong Zhang   ierr = VecGetArray(auglydd,&d);CHKERRQ(ierr);
2912cd6534aSBarry Smith   for (i=0; i<n; i++) {
2922cd6534aSBarry Smith     d[i] = s[auglyrmapd[i]]; /* copy "diagonal" (true local) portion of scale into dd vector */
2932cd6534aSBarry Smith   }
2941ebc52fbSHong Zhang   ierr = VecRestoreArray(auglydd,&d);CHKERRQ(ierr);
2952cd6534aSBarry Smith   /* column scale "diagonal" portion of local matrix */
2960298fd71SBarry Smith   ierr = MatDiagonalScale(a->A,NULL,auglydd);CHKERRQ(ierr);
2972cd6534aSBarry Smith 
2982cd6534aSBarry Smith   ierr = VecGetLocalSize(auglyoo,&n);CHKERRQ(ierr);
2991ebc52fbSHong Zhang   ierr = VecGetArray(auglyoo,&o);CHKERRQ(ierr);
3002cd6534aSBarry Smith   for (i=0; i<n; i++) {
3012cd6534aSBarry Smith     o[i] = s[auglyrmapo[i]]; /* copy "off-diagonal" portion of scale into oo vector */
3022cd6534aSBarry Smith   }
3031ebc52fbSHong Zhang   ierr = VecRestoreArray(scale,&s);CHKERRQ(ierr);
3041ebc52fbSHong Zhang   ierr = VecRestoreArray(auglyoo,&o);CHKERRQ(ierr);
3052cd6534aSBarry Smith   /* column scale "off-diagonal" portion of local matrix */
3060298fd71SBarry Smith   ierr = MatDiagonalScale(a->B,NULL,auglyoo);CHKERRQ(ierr);
3072cd6534aSBarry Smith   PetscFunctionReturn(0);
3082cd6534aSBarry Smith }
3092cd6534aSBarry Smith 
31048b35521SBarry Smith 
311