xref: /petsc/src/mat/impls/aij/mpi/mmaij.c (revision 3bb1ff401821b9e2ae019d3e61bc8ab4bd4e59d5)
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 */
41b1d57f15SBarry Smith   ierr = PetscMalloc((ec+1)*sizeof(PetscInt),&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 */
70b1d57f15SBarry Smith   ierr = PetscMalloc((N+1)*sizeof(PetscInt),&indices);CHKERRQ(ierr);
71b1d57f15SBarry Smith   ierr = PetscMemzero(indices,N*sizeof(PetscInt));CHKERRQ(ierr);
72d0f46423SBarry Smith   for (i=0; i<aij->B->rmap->n; i++) {
73d6dfbf8fSBarry Smith     for (j=0; j<B->ilen[i]; j++) {
74bfec09a0SHong Zhang       if (!indices[aj[B->i[i] + j]]) ec++;
75bfec09a0SHong Zhang       indices[aj[B->i[i] + j]] = 1;
76416022c9SBarry Smith     }
771eb62cbbSBarry Smith   }
788c79f6d3SBarry Smith 
791eb62cbbSBarry Smith   /* form array of columns we need */
80b1d57f15SBarry Smith   ierr = PetscMalloc((ec+1)*sizeof(PetscInt),&garray);CHKERRQ(ierr);
811eb62cbbSBarry Smith   ec   = 0;
821eb62cbbSBarry Smith   for (i=0; i<N; i++) {
831eb62cbbSBarry Smith     if (indices[i]) garray[ec++] = i;
841eb62cbbSBarry Smith   }
851eb62cbbSBarry Smith 
861eb62cbbSBarry Smith   /* make indices now point into garray */
871eb62cbbSBarry Smith   for (i=0; i<ec; i++) {
88bfec09a0SHong Zhang     indices[garray[i]] = i;
891eb62cbbSBarry Smith   }
901eb62cbbSBarry Smith 
911eb62cbbSBarry Smith   /* compact out the extra columns in B */
92d0f46423SBarry Smith   for (i=0; i<aij->B->rmap->n; i++) {
93d6dfbf8fSBarry Smith     for (j=0; j<B->ilen[i]; j++) {
94bfec09a0SHong Zhang       aj[B->i[i] + j] = indices[aj[B->i[i] + j]];
951eb62cbbSBarry Smith     }
96d6dfbf8fSBarry Smith   }
97d0f46423SBarry Smith   aij->B->cmap->n = aij->B->cmap->N = ec;
982205254eSKarl Rupp 
9926283091SBarry Smith   ierr = PetscLayoutSetUp((aij->B->cmap));CHKERRQ(ierr);
100606d414cSSatish Balay   ierr = PetscFree(indices);CHKERRQ(ierr);
1012066d6f7SSatish Balay #endif
1021eb62cbbSBarry Smith   /* create local vector that is used to scatter into */
103029af93fSBarry Smith   ierr = VecCreateSeq(PETSC_COMM_SELF,ec,&aij->lvec);CHKERRQ(ierr);
1041eb62cbbSBarry Smith 
105d6dfbf8fSBarry Smith   /* create two temporary Index sets for build scatter gather */
10670b3c8c7SBarry Smith   ierr = ISCreateGeneral(((PetscObject)mat)->comm,ec,garray,PETSC_COPY_VALUES,&from);CHKERRQ(ierr);
107c5bfad50SMark F. Adams 
108029af93fSBarry Smith   ierr = ISCreateStride(PETSC_COMM_SELF,ec,0,1,&to);CHKERRQ(ierr);
1091eb62cbbSBarry Smith 
1101eb62cbbSBarry Smith   /* create temporary global vector to generate scatter context */
111b5eb4454SBarry Smith   /* This does not allocate the array's memory so is efficient */
112ce94432eSBarry Smith   ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)mat),1,mat->cmap->n,mat->cmap->N,NULL,&gvec);CHKERRQ(ierr);
1131eb62cbbSBarry Smith 
1142d336d48SLois Curfman McInnes   /* generate the scatter context */
11508480c60SBarry Smith   ierr = VecScatterCreate(gvec,from,aij->lvec,to,&aij->Mvctx);CHKERRQ(ierr);
116*3bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)aij->Mvctx);CHKERRQ(ierr);
117*3bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)aij->lvec);CHKERRQ(ierr);
118*3bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)from);CHKERRQ(ierr);
119*3bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)to);CHKERRQ(ierr);
1202205254eSKarl Rupp 
1219e25ed09SBarry Smith   aij->garray = garray;
1222205254eSKarl Rupp 
123*3bb1ff40SBarry Smith   ierr = PetscLogObjectMemory((PetscObject)mat,(ec+1)*sizeof(PetscInt));CHKERRQ(ierr);
1246bf464f9SBarry Smith   ierr = ISDestroy(&from);CHKERRQ(ierr);
1256bf464f9SBarry Smith   ierr = ISDestroy(&to);CHKERRQ(ierr);
1266bf464f9SBarry Smith   ierr = VecDestroy(&gvec);CHKERRQ(ierr);
1273a40ed3dSBarry Smith   PetscFunctionReturn(0);
1288c79f6d3SBarry Smith }
1299e25ed09SBarry Smith 
1309e25ed09SBarry Smith 
1314a2ae208SSatish Balay #undef __FUNCT__
132ab9863d7SBarry Smith #define __FUNCT__ "MatDisAssemble_MPIAIJ"
1332493cbb0SBarry Smith /*
1342493cbb0SBarry Smith      Takes the local part of an already assembled MPIAIJ matrix
1352493cbb0SBarry Smith    and disassembles it. This is to allow new nonzeros into the matrix
1362493cbb0SBarry Smith    that require more communication in the matrix vector multiply.
1372493cbb0SBarry Smith    Thus certain data-structures must be rebuilt.
1382493cbb0SBarry Smith 
1392493cbb0SBarry Smith    Kind of slow! But that's what application programmers get when
1402493cbb0SBarry Smith    they are sloppy.
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;
147dfbe8321SBarry Smith   PetscErrorCode ierr;
148d0f46423SBarry Smith   PetscInt       i,j,m = B->rmap->n,n = A->cmap->N,col,ct = 0,*garray = aij->garray,*nz,ec;
14987828ca2SBarry Smith   PetscScalar    v;
1502493cbb0SBarry Smith 
1513a40ed3dSBarry Smith   PetscFunctionBegin;
1522493cbb0SBarry Smith   /* free stuff related to matrix-vec multiply */
153b0a32e0cSBarry Smith   ierr = VecGetSize(aij->lvec,&ec);CHKERRQ(ierr); /* needed for PetscLogObjectMemory below */
1545e1f6667SBarry Smith   ierr = VecDestroy(&aij->lvec);CHKERRQ(ierr);
1555e1f6667SBarry Smith   ierr = VecScatterDestroy(&aij->Mvctx);CHKERRQ(ierr);
156464493b3SBarry Smith   if (aij->colmap) {
157aa482453SBarry Smith #if defined(PETSC_USE_CTABLE)
1586bc0bbbfSBarry Smith     ierr = PetscTableDestroy(&aij->colmap);CHKERRQ(ierr);
1592066d6f7SSatish Balay #else
160606d414cSSatish Balay     ierr = PetscFree(aij->colmap);CHKERRQ(ierr);
161*3bb1ff40SBarry Smith     ierr = PetscLogObjectMemory((PetscObject)A,-aij->B->cmap->n*sizeof(PetscInt));CHKERRQ(ierr);
1622066d6f7SSatish Balay #endif
163464493b3SBarry Smith   }
1642493cbb0SBarry Smith 
1652493cbb0SBarry Smith   /* make sure that B is assembled so we can access its values */
1666d4a8577SBarry Smith   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
167fe2f2677SBarry Smith   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1682493cbb0SBarry Smith 
1692493cbb0SBarry Smith   /* invent new B and copy stuff over */
170b1d57f15SBarry Smith   ierr = PetscMalloc((m+1)*sizeof(PetscInt),&nz);CHKERRQ(ierr);
17148b35521SBarry Smith   for (i=0; i<m; i++) {
17248b35521SBarry Smith     nz[i] = Baij->i[i+1] - Baij->i[i];
17348b35521SBarry Smith   }
174f69a0ea3SMatthew Knepley   ierr = MatCreate(PETSC_COMM_SELF,&Bnew);CHKERRQ(ierr);
175f69a0ea3SMatthew Knepley   ierr = MatSetSizes(Bnew,m,n,m,n);CHKERRQ(ierr);
176a2f3521dSMark F. Adams   ierr = MatSetBlockSizes(Bnew,A->rmap->bs,A->cmap->bs);CHKERRQ(ierr);
1777adad957SLisandro Dalcin   ierr = MatSetType(Bnew,((PetscObject)B)->type_name);CHKERRQ(ierr);
178f204ca49SKris Buschelman   ierr = MatSeqAIJSetPreallocation(Bnew,0,nz);CHKERRQ(ierr);
1792205254eSKarl Rupp 
1802576faa2SJed Brown   ((Mat_SeqAIJ*)Bnew->data)->nonew = Baij->nonew; /* Inherit insertion error options. */
1812205254eSKarl Rupp 
182606d414cSSatish Balay   ierr = PetscFree(nz);CHKERRQ(ierr);
1832493cbb0SBarry Smith   for (i=0; i<m; i++) {
184bfec09a0SHong Zhang     for (j=Baij->i[i]; j<Baij->i[i+1]; j++) {
185bfec09a0SHong Zhang       col  = garray[Baij->j[ct]];
1862493cbb0SBarry Smith       v    = Baij->a[ct++];
18783271157SBarry Smith       ierr = MatSetValues(Bnew,1,&i,1,&col,&v,B->insertmode);CHKERRQ(ierr);
1882493cbb0SBarry Smith     }
1892493cbb0SBarry Smith   }
190606d414cSSatish Balay   ierr = PetscFree(aij->garray);CHKERRQ(ierr);
191*3bb1ff40SBarry Smith   ierr = PetscLogObjectMemory((PetscObject)A,-ec*sizeof(PetscInt));CHKERRQ(ierr);
1926bf464f9SBarry Smith   ierr = MatDestroy(&B);CHKERRQ(ierr);
193*3bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)A,(PetscObject)Bnew);CHKERRQ(ierr);
1942205254eSKarl Rupp 
1952493cbb0SBarry Smith   aij->B           = Bnew;
196227d817aSBarry Smith   A->was_assembled = PETSC_FALSE;
1973a40ed3dSBarry Smith   PetscFunctionReturn(0);
1982493cbb0SBarry Smith }
1992493cbb0SBarry Smith 
2002cd6534aSBarry Smith /*      ugly stuff added for Glenn someday we should fix this up */
2012cd6534aSBarry Smith 
2022205254eSKarl Rupp static PetscInt *auglyrmapd = 0,*auglyrmapo = 0; /* mapping from the local ordering to the "diagonal" and "off-diagonal" parts of the local matrix */
2032cd6534aSBarry Smith static Vec auglydd          = 0,auglyoo     = 0; /* work vectors used to scale the two parts of the local matrix */
2042cd6534aSBarry Smith 
2052cd6534aSBarry Smith 
2062cd6534aSBarry Smith #undef __FUNCT__
2072cd6534aSBarry Smith #define __FUNCT__ "MatMPIAIJDiagonalScaleLocalSetUp"
208dfbe8321SBarry Smith PetscErrorCode MatMPIAIJDiagonalScaleLocalSetUp(Mat inA,Vec scale)
2092cd6534aSBarry Smith {
2102cd6534aSBarry Smith   Mat_MPIAIJ     *ina = (Mat_MPIAIJ*) inA->data; /*access private part of matrix */
211dfbe8321SBarry Smith   PetscErrorCode ierr;
212b1d57f15SBarry Smith   PetscInt       i,n,nt,cstart,cend,no,*garray = ina->garray,*lindices;
213b1d57f15SBarry Smith   PetscInt       *r_rmapd,*r_rmapo;
2142cd6534aSBarry Smith 
2152cd6534aSBarry Smith   PetscFunctionBegin;
2162cd6534aSBarry Smith   ierr = MatGetOwnershipRange(inA,&cstart,&cend);CHKERRQ(ierr);
2170298fd71SBarry Smith   ierr = MatGetSize(ina->A,NULL,&n);CHKERRQ(ierr);
218992144d0SBarry Smith   ierr = PetscMalloc((inA->rmap->mapping->n+1)*sizeof(PetscInt),&r_rmapd);CHKERRQ(ierr);
219992144d0SBarry Smith   ierr = PetscMemzero(r_rmapd,inA->rmap->mapping->n*sizeof(PetscInt));CHKERRQ(ierr);
2202cd6534aSBarry Smith   nt   = 0;
221992144d0SBarry Smith   for (i=0; i<inA->rmap->mapping->n; i++) {
222992144d0SBarry Smith     if (inA->rmap->mapping->indices[i] >= cstart && inA->rmap->mapping->indices[i] < cend) {
2232cd6534aSBarry Smith       nt++;
224992144d0SBarry Smith       r_rmapd[i] = inA->rmap->mapping->indices[i] + 1;
2252cd6534aSBarry Smith     }
2262cd6534aSBarry Smith   }
227e32f2f54SBarry Smith   if (nt != n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hmm nt %D n %D",nt,n);
228b1d57f15SBarry Smith   ierr = PetscMalloc((n+1)*sizeof(PetscInt),&auglyrmapd);CHKERRQ(ierr);
229992144d0SBarry Smith   for (i=0; i<inA->rmap->mapping->n; i++) {
2302cd6534aSBarry Smith     if (r_rmapd[i]) {
2312cd6534aSBarry Smith       auglyrmapd[(r_rmapd[i]-1)-cstart] = i;
2322cd6534aSBarry Smith     }
2332cd6534aSBarry Smith   }
2342cd6534aSBarry Smith   ierr = PetscFree(r_rmapd);CHKERRQ(ierr);
2352cd6534aSBarry Smith   ierr = VecCreateSeq(PETSC_COMM_SELF,n,&auglydd);CHKERRQ(ierr);
2362cd6534aSBarry Smith 
237d0f46423SBarry Smith   ierr = PetscMalloc((inA->cmap->N+1)*sizeof(PetscInt),&lindices);CHKERRQ(ierr);
238d0f46423SBarry Smith   ierr = PetscMemzero(lindices,inA->cmap->N*sizeof(PetscInt));CHKERRQ(ierr);
239d0f46423SBarry Smith   for (i=0; i<ina->B->cmap->n; i++) {
2402cd6534aSBarry Smith     lindices[garray[i]] = i+1;
2412cd6534aSBarry Smith   }
242992144d0SBarry Smith   no   = inA->rmap->mapping->n - nt;
243992144d0SBarry Smith   ierr = PetscMalloc((inA->rmap->mapping->n+1)*sizeof(PetscInt),&r_rmapo);CHKERRQ(ierr);
244992144d0SBarry Smith   ierr = PetscMemzero(r_rmapo,inA->rmap->mapping->n*sizeof(PetscInt));CHKERRQ(ierr);
2452cd6534aSBarry Smith   nt   = 0;
246992144d0SBarry Smith   for (i=0; i<inA->rmap->mapping->n; i++) {
247992144d0SBarry Smith     if (lindices[inA->rmap->mapping->indices[i]]) {
2482cd6534aSBarry Smith       nt++;
249992144d0SBarry Smith       r_rmapo[i] = lindices[inA->rmap->mapping->indices[i]];
2502cd6534aSBarry Smith     }
2512cd6534aSBarry Smith   }
252e32f2f54SBarry Smith   if (nt > no) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hmm nt %D no %D",nt,n);
2532cd6534aSBarry Smith   ierr = PetscFree(lindices);CHKERRQ(ierr);
254b1d57f15SBarry Smith   ierr = PetscMalloc((nt+1)*sizeof(PetscInt),&auglyrmapo);CHKERRQ(ierr);
255992144d0SBarry Smith   for (i=0; i<inA->rmap->mapping->n; i++) {
2562cd6534aSBarry Smith     if (r_rmapo[i]) {
2572cd6534aSBarry Smith       auglyrmapo[(r_rmapo[i]-1)] = i;
2582cd6534aSBarry Smith     }
2592cd6534aSBarry Smith   }
2602cd6534aSBarry Smith   ierr = PetscFree(r_rmapo);CHKERRQ(ierr);
2612cd6534aSBarry Smith   ierr = VecCreateSeq(PETSC_COMM_SELF,nt,&auglyoo);CHKERRQ(ierr);
2622cd6534aSBarry Smith   PetscFunctionReturn(0);
2632cd6534aSBarry Smith }
2642cd6534aSBarry Smith 
2652cd6534aSBarry Smith #undef __FUNCT__
2662cd6534aSBarry Smith #define __FUNCT__ "MatMPIAIJDiagonalScaleLocal"
267dfbe8321SBarry Smith PetscErrorCode MatMPIAIJDiagonalScaleLocal(Mat A,Vec scale)
2682cd6534aSBarry Smith {
26992b32695SKris Buschelman   /* This routine should really be abandoned as it duplicates MatDiagonalScaleLocal */
2704ac538c5SBarry Smith   PetscErrorCode ierr;
27192b32695SKris Buschelman 
27292b32695SKris Buschelman   PetscFunctionBegin;
2734ac538c5SBarry Smith   ierr = PetscTryMethod(A,"MatDiagonalScaleLocal_C",(Mat,Vec),(A,scale));CHKERRQ(ierr);
27492b32695SKris Buschelman   PetscFunctionReturn(0);
27592b32695SKris Buschelman }
27692b32695SKris Buschelman 
27792b32695SKris Buschelman #undef __FUNCT__
27892b32695SKris Buschelman #define __FUNCT__ "MatDiagonalScaleLocal_MPIAIJ"
2797087cfbeSBarry Smith PetscErrorCode  MatDiagonalScaleLocal_MPIAIJ(Mat A,Vec scale)
28092b32695SKris Buschelman {
2812cd6534aSBarry Smith   Mat_MPIAIJ     *a = (Mat_MPIAIJ*) A->data; /*access private part of matrix */
282dfbe8321SBarry Smith   PetscErrorCode ierr;
283b1d57f15SBarry Smith   PetscInt       n,i;
2842cd6534aSBarry Smith   PetscScalar    *d,*o,*s;
2852cd6534aSBarry Smith 
2862cd6534aSBarry Smith   PetscFunctionBegin;
2872cd6534aSBarry Smith   if (!auglyrmapd) {
2882cd6534aSBarry Smith     ierr = MatMPIAIJDiagonalScaleLocalSetUp(A,scale);CHKERRQ(ierr);
2892cd6534aSBarry Smith   }
2902cd6534aSBarry Smith 
2911ebc52fbSHong Zhang   ierr = VecGetArray(scale,&s);CHKERRQ(ierr);
2922cd6534aSBarry Smith 
2932cd6534aSBarry Smith   ierr = VecGetLocalSize(auglydd,&n);CHKERRQ(ierr);
2941ebc52fbSHong Zhang   ierr = VecGetArray(auglydd,&d);CHKERRQ(ierr);
2952cd6534aSBarry Smith   for (i=0; i<n; i++) {
2962cd6534aSBarry Smith     d[i] = s[auglyrmapd[i]]; /* copy "diagonal" (true local) portion of scale into dd vector */
2972cd6534aSBarry Smith   }
2981ebc52fbSHong Zhang   ierr = VecRestoreArray(auglydd,&d);CHKERRQ(ierr);
2992cd6534aSBarry Smith   /* column scale "diagonal" portion of local matrix */
3000298fd71SBarry Smith   ierr = MatDiagonalScale(a->A,NULL,auglydd);CHKERRQ(ierr);
3012cd6534aSBarry Smith 
3022cd6534aSBarry Smith   ierr = VecGetLocalSize(auglyoo,&n);CHKERRQ(ierr);
3031ebc52fbSHong Zhang   ierr = VecGetArray(auglyoo,&o);CHKERRQ(ierr);
3042cd6534aSBarry Smith   for (i=0; i<n; i++) {
3052cd6534aSBarry Smith     o[i] = s[auglyrmapo[i]]; /* copy "off-diagonal" portion of scale into oo vector */
3062cd6534aSBarry Smith   }
3071ebc52fbSHong Zhang   ierr = VecRestoreArray(scale,&s);CHKERRQ(ierr);
3081ebc52fbSHong Zhang   ierr = VecRestoreArray(auglyoo,&o);CHKERRQ(ierr);
3092cd6534aSBarry Smith   /* column scale "off-diagonal" portion of local matrix */
3100298fd71SBarry Smith   ierr = MatDiagonalScale(a->B,NULL,auglyoo);CHKERRQ(ierr);
3112cd6534aSBarry Smith   PetscFunctionReturn(0);
3122cd6534aSBarry Smith }
3132cd6534aSBarry Smith 
31448b35521SBarry Smith 
315