xref: /petsc/src/mat/impls/aij/mpi/mmaij.c (revision 3f6a6bda4d6eca75b1f8884d49cdb7fdd0dd3ec2)
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>
6af0996ceSBarry Smith #include <petsc/private/isimpl.h>    /* needed because accesses data structure of ISLocalToGlobalMapping directly */
78c79f6d3SBarry Smith 
8dfbe8321SBarry Smith PetscErrorCode MatSetUpMultiply_MPIAIJ(Mat mat)
98c79f6d3SBarry Smith {
1044a69424SLois Curfman McInnes   Mat_MPIAIJ     *aij = (Mat_MPIAIJ*)mat->data;
11ec8511deSBarry Smith   Mat_SeqAIJ     *B   = (Mat_SeqAIJ*)(aij->B->data);
126849ba73SBarry Smith   PetscErrorCode ierr;
13b1d57f15SBarry Smith   PetscInt       i,j,*aj = B->j,ec = 0,*garray;
141eb62cbbSBarry Smith   IS             from,to;
151eb62cbbSBarry Smith   Vec            gvec;
16aa482453SBarry Smith #if defined(PETSC_USE_CTABLE)
170f5bd95cSBarry Smith   PetscTable         gid1_lid1;
180f5bd95cSBarry Smith   PetscTablePosition tpos;
19b1d57f15SBarry Smith   PetscInt           gid,lid;
206f531f54SSatish Balay #else
21d0f46423SBarry Smith   PetscInt N = mat->cmap->N,*indices;
222066d6f7SSatish Balay #endif
232066d6f7SSatish Balay 
243a40ed3dSBarry Smith   PetscFunctionBegin;
254b8d542aSHong Zhang   if (!aij->garray) {
26aa482453SBarry Smith #if defined(PETSC_USE_CTABLE)
27c5bfad50SMark F. Adams     /* use a table */
28e23dfa41SBarry Smith     ierr = PetscTableCreate(aij->B->rmap->n,mat->cmap->N+1,&gid1_lid1);CHKERRQ(ierr);
29d0f46423SBarry Smith     for (i=0; i<aij->B->rmap->n; i++) {
302066d6f7SSatish Balay       for (j=0; j<B->ilen[i]; j++) {
31b1d57f15SBarry Smith         PetscInt data,gid1 = aj[B->i[i] + j] + 1;
320f5bd95cSBarry Smith         ierr = PetscTableFind(gid1_lid1,gid1,&data);CHKERRQ(ierr);
33fa46199cSSatish Balay         if (!data) {
342066d6f7SSatish Balay           /* one based table */
353861aac3SJed Brown           ierr = PetscTableAdd(gid1_lid1,gid1,++ec,INSERT_VALUES);CHKERRQ(ierr);
362066d6f7SSatish Balay         }
372066d6f7SSatish Balay       }
382066d6f7SSatish Balay     }
392066d6f7SSatish 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);
422066d6f7SSatish Balay     while (tpos) {
430f5bd95cSBarry Smith       ierr = PetscTableGetNext(gid1_lid1,&tpos,&gid,&lid);CHKERRQ(ierr);
44b0a32e0cSBarry Smith       gid--;
45b0a32e0cSBarry Smith       lid--;
462066d6f7SSatish Balay       garray[lid] = gid;
472066d6f7SSatish Balay     }
480064e2bbSSatish Balay     ierr = PetscSortInt(ec,garray);CHKERRQ(ierr); /* sort, and rebuild */
490f5bd95cSBarry Smith     ierr = PetscTableRemoveAll(gid1_lid1);CHKERRQ(ierr);
502066d6f7SSatish Balay     for (i=0; i<ec; i++) {
513861aac3SJed Brown       ierr = PetscTableAdd(gid1_lid1,garray[i]+1,i+1,INSERT_VALUES);CHKERRQ(ierr);
522066d6f7SSatish Balay     }
532066d6f7SSatish Balay     /* compact out the extra columns in B */
54d0f46423SBarry Smith     for (i=0; i<aij->B->rmap->n; i++) {
552066d6f7SSatish Balay       for (j=0; j<B->ilen[i]; j++) {
56b1d57f15SBarry Smith         PetscInt gid1 = aj[B->i[i] + j] + 1;
570f5bd95cSBarry Smith         ierr = PetscTableFind(gid1_lid1,gid1,&lid);CHKERRQ(ierr);
58fa46199cSSatish Balay         lid--;
59b3fb0a6cSMatthew Knepley         aj[B->i[i] + j] = lid;
602066d6f7SSatish Balay       }
612066d6f7SSatish Balay     }
62d0f46423SBarry Smith     aij->B->cmap->n = aij->B->cmap->N = ec;
63cdce4254SBarry Smith     aij->B->cmap->bs = 1;
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 */
701795a4d1SJed 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 */
79854ce69bSBarry Smith     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;
97cd0e7f71SBarry Smith     aij->B->cmap->bs = 1;
982205254eSKarl Rupp 
9926283091SBarry Smith     ierr = PetscLayoutSetUp((aij->B->cmap));CHKERRQ(ierr);
100606d414cSSatish Balay     ierr = PetscFree(indices);CHKERRQ(ierr);
1012066d6f7SSatish Balay #endif
1024b8d542aSHong Zhang   } else {
1034b8d542aSHong Zhang     garray = aij->garray;
1044b8d542aSHong Zhang   }
1054b8d542aSHong Zhang 
1064b8d542aSHong Zhang   if (!aij->lvec) {
1071eb62cbbSBarry Smith     /* create local vector that is used to scatter into */
108029af93fSBarry Smith     ierr = VecCreateSeq(PETSC_COMM_SELF,ec,&aij->lvec);CHKERRQ(ierr);
1094b8d542aSHong Zhang   } else {
1104b8d542aSHong Zhang     ierr = VecGetSize(aij->lvec,&ec);CHKERRQ(ierr);
1114b8d542aSHong Zhang   }
1121eb62cbbSBarry Smith 
113d6dfbf8fSBarry Smith   /* create two temporary Index sets for build scatter gather */
11470b3c8c7SBarry Smith   ierr = ISCreateGeneral(((PetscObject)mat)->comm,ec,garray,PETSC_COPY_VALUES,&from);CHKERRQ(ierr);
115c5bfad50SMark F. Adams 
116029af93fSBarry Smith   ierr = ISCreateStride(PETSC_COMM_SELF,ec,0,1,&to);CHKERRQ(ierr);
1171eb62cbbSBarry Smith 
1181eb62cbbSBarry Smith   /* create temporary global vector to generate scatter context */
119b5eb4454SBarry Smith   /* This does not allocate the array's memory so is efficient */
120ce94432eSBarry Smith   ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)mat),1,mat->cmap->n,mat->cmap->N,NULL,&gvec);CHKERRQ(ierr);
1211eb62cbbSBarry Smith 
1222d336d48SLois Curfman McInnes   /* generate the scatter context */
12301ad2aeeSHong Zhang   if (aij->Mvctx_mpi1_flg || mat->mpi1) {
12401ad2aeeSHong Zhang     if (aij->Mvctx_mpi1) {
12501ad2aeeSHong Zhang       /* Must destroy existing Mvctx_mpi1, then recreate it. See src/ksp/ksp/examples/tutorials/network/runex1_nest_2 */
12601ad2aeeSHong Zhang       ierr = VecScatterDestroy(&aij->Mvctx_mpi1);CHKERRQ(ierr);
12701ad2aeeSHong Zhang     }
1284b8d542aSHong Zhang     ierr = VecScatterCreateMPI1(gvec,from,aij->lvec,to,&aij->Mvctx_mpi1);CHKERRQ(ierr);
1294b8d542aSHong Zhang     ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)aij->Mvctx_mpi1);CHKERRQ(ierr);
1304b8d542aSHong Zhang   } else {
131*3f6a6bdaSHong Zhang     if (aij->Mvctx) {
132*3f6a6bdaSHong Zhang       ierr = VecScatterDestroy(&aij->Mvctx);CHKERRQ(ierr);
133*3f6a6bdaSHong Zhang     }
134*3f6a6bdaSHong Zhang 
13508480c60SBarry Smith     ierr = VecScatterCreate(gvec,from,aij->lvec,to,&aij->Mvctx);CHKERRQ(ierr);
1363bb1ff40SBarry Smith     ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)aij->Mvctx);CHKERRQ(ierr);
1373bb1ff40SBarry Smith     ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)aij->lvec);CHKERRQ(ierr);
1384b8d542aSHong Zhang     ierr = PetscLogObjectMemory((PetscObject)mat,(ec+1)*sizeof(PetscInt));CHKERRQ(ierr);
1394b8d542aSHong Zhang   }
14067bb5161SHong Zhang   aij->garray = garray;
1414b8d542aSHong Zhang 
1423bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)from);CHKERRQ(ierr);
1433bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)to);CHKERRQ(ierr);
1442205254eSKarl Rupp 
1456bf464f9SBarry Smith   ierr = ISDestroy(&from);CHKERRQ(ierr);
1466bf464f9SBarry Smith   ierr = ISDestroy(&to);CHKERRQ(ierr);
1476bf464f9SBarry Smith   ierr = VecDestroy(&gvec);CHKERRQ(ierr);
1483a40ed3dSBarry Smith   PetscFunctionReturn(0);
1498c79f6d3SBarry Smith }
1509e25ed09SBarry Smith 
1512493cbb0SBarry Smith /*
1522493cbb0SBarry Smith      Takes the local part of an already assembled MPIAIJ matrix
1532493cbb0SBarry Smith    and disassembles it. This is to allow new nonzeros into the matrix
1542493cbb0SBarry Smith    that require more communication in the matrix vector multiply.
1552493cbb0SBarry Smith    Thus certain data-structures must be rebuilt.
1562493cbb0SBarry Smith 
1572493cbb0SBarry Smith    Kind of slow! But that's what application programmers get when
1582493cbb0SBarry Smith    they are sloppy.
1592493cbb0SBarry Smith */
160ab9863d7SBarry Smith PetscErrorCode MatDisAssemble_MPIAIJ(Mat A)
1612493cbb0SBarry Smith {
1622493cbb0SBarry Smith   Mat_MPIAIJ     *aij  = (Mat_MPIAIJ*)A->data;
1632493cbb0SBarry Smith   Mat            B     = aij->B,Bnew;
164ec8511deSBarry Smith   Mat_SeqAIJ     *Baij = (Mat_SeqAIJ*)B->data;
165dfbe8321SBarry Smith   PetscErrorCode ierr;
166d0f46423SBarry Smith   PetscInt       i,j,m = B->rmap->n,n = A->cmap->N,col,ct = 0,*garray = aij->garray,*nz,ec;
16787828ca2SBarry Smith   PetscScalar    v;
1682493cbb0SBarry Smith 
1693a40ed3dSBarry Smith   PetscFunctionBegin;
1702493cbb0SBarry Smith   /* free stuff related to matrix-vec multiply */
171b0a32e0cSBarry Smith   ierr = VecGetSize(aij->lvec,&ec);CHKERRQ(ierr); /* needed for PetscLogObjectMemory below */
1725e1f6667SBarry Smith   ierr = VecDestroy(&aij->lvec);CHKERRQ(ierr);
173464493b3SBarry Smith   if (aij->colmap) {
174aa482453SBarry Smith #if defined(PETSC_USE_CTABLE)
1756bc0bbbfSBarry Smith     ierr = PetscTableDestroy(&aij->colmap);CHKERRQ(ierr);
1762066d6f7SSatish Balay #else
177606d414cSSatish Balay     ierr = PetscFree(aij->colmap);CHKERRQ(ierr);
1783bb1ff40SBarry Smith     ierr = PetscLogObjectMemory((PetscObject)A,-aij->B->cmap->n*sizeof(PetscInt));CHKERRQ(ierr);
1792066d6f7SSatish Balay #endif
180464493b3SBarry Smith   }
1812493cbb0SBarry Smith 
1822493cbb0SBarry Smith   /* make sure that B is assembled so we can access its values */
1836d4a8577SBarry Smith   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
184fe2f2677SBarry Smith   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1852493cbb0SBarry Smith 
1862493cbb0SBarry Smith   /* invent new B and copy stuff over */
187854ce69bSBarry Smith   ierr = PetscMalloc1(m+1,&nz);CHKERRQ(ierr);
18848b35521SBarry Smith   for (i=0; i<m; i++) {
18948b35521SBarry Smith     nz[i] = Baij->i[i+1] - Baij->i[i];
19048b35521SBarry Smith   }
191f69a0ea3SMatthew Knepley   ierr = MatCreate(PETSC_COMM_SELF,&Bnew);CHKERRQ(ierr);
192f69a0ea3SMatthew Knepley   ierr = MatSetSizes(Bnew,m,n,m,n);CHKERRQ(ierr);
19333d57670SJed Brown   ierr = MatSetBlockSizesFromMats(Bnew,A,A);CHKERRQ(ierr);
1947adad957SLisandro Dalcin   ierr = MatSetType(Bnew,((PetscObject)B)->type_name);CHKERRQ(ierr);
195f204ca49SKris Buschelman   ierr = MatSeqAIJSetPreallocation(Bnew,0,nz);CHKERRQ(ierr);
1962205254eSKarl Rupp 
1972576faa2SJed Brown   ((Mat_SeqAIJ*)Bnew->data)->nonew = Baij->nonew; /* Inherit insertion error options. */
19877341eacSDmitry Karpeev   /*
19977341eacSDmitry Karpeev    Ensure that B's nonzerostate is monotonically increasing.
20077341eacSDmitry Karpeev    Or should this follow the MatSetValues() loop to preserve B's nonzerstate across a MatDisAssemble() call?
20177341eacSDmitry Karpeev    */
202f69fde56SShane Stafford   Bnew->nonzerostate = B->nonzerostate;
2032205254eSKarl Rupp 
204606d414cSSatish Balay   ierr = PetscFree(nz);CHKERRQ(ierr);
2052493cbb0SBarry Smith   for (i=0; i<m; i++) {
206bfec09a0SHong Zhang     for (j=Baij->i[i]; j<Baij->i[i+1]; j++) {
207bfec09a0SHong Zhang       col  = garray[Baij->j[ct]];
2082493cbb0SBarry Smith       v    = Baij->a[ct++];
20983271157SBarry Smith       ierr = MatSetValues(Bnew,1,&i,1,&col,&v,B->insertmode);CHKERRQ(ierr);
2102493cbb0SBarry Smith     }
2112493cbb0SBarry Smith   }
212606d414cSSatish Balay   ierr = PetscFree(aij->garray);CHKERRQ(ierr);
2133bb1ff40SBarry Smith   ierr = PetscLogObjectMemory((PetscObject)A,-ec*sizeof(PetscInt));CHKERRQ(ierr);
2146bf464f9SBarry Smith   ierr = MatDestroy(&B);CHKERRQ(ierr);
2153bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)A,(PetscObject)Bnew);CHKERRQ(ierr);
2162205254eSKarl Rupp 
2172493cbb0SBarry Smith   aij->B           = Bnew;
218227d817aSBarry Smith   A->was_assembled = PETSC_FALSE;
2193a40ed3dSBarry Smith   PetscFunctionReturn(0);
2202493cbb0SBarry Smith }
2212493cbb0SBarry Smith 
2222cd6534aSBarry Smith /*      ugly stuff added for Glenn someday we should fix this up */
2232cd6534aSBarry Smith 
2242205254eSKarl Rupp static PetscInt *auglyrmapd = 0,*auglyrmapo = 0; /* mapping from the local ordering to the "diagonal" and "off-diagonal" parts of the local matrix */
2252cd6534aSBarry Smith static Vec auglydd          = 0,auglyoo     = 0; /* work vectors used to scale the two parts of the local matrix */
2262cd6534aSBarry Smith 
2272cd6534aSBarry Smith 
228dfbe8321SBarry Smith PetscErrorCode MatMPIAIJDiagonalScaleLocalSetUp(Mat inA,Vec scale)
2292cd6534aSBarry Smith {
2302cd6534aSBarry Smith   Mat_MPIAIJ     *ina = (Mat_MPIAIJ*) inA->data; /*access private part of matrix */
231dfbe8321SBarry Smith   PetscErrorCode ierr;
232b1d57f15SBarry Smith   PetscInt       i,n,nt,cstart,cend,no,*garray = ina->garray,*lindices;
233b1d57f15SBarry Smith   PetscInt       *r_rmapd,*r_rmapo;
2342cd6534aSBarry Smith 
2352cd6534aSBarry Smith   PetscFunctionBegin;
2362cd6534aSBarry Smith   ierr = MatGetOwnershipRange(inA,&cstart,&cend);CHKERRQ(ierr);
2370298fd71SBarry Smith   ierr = MatGetSize(ina->A,NULL,&n);CHKERRQ(ierr);
2381795a4d1SJed Brown   ierr = PetscCalloc1(inA->rmap->mapping->n+1,&r_rmapd);CHKERRQ(ierr);
2392cd6534aSBarry Smith   nt   = 0;
240992144d0SBarry Smith   for (i=0; i<inA->rmap->mapping->n; i++) {
241992144d0SBarry Smith     if (inA->rmap->mapping->indices[i] >= cstart && inA->rmap->mapping->indices[i] < cend) {
2422cd6534aSBarry Smith       nt++;
243992144d0SBarry Smith       r_rmapd[i] = inA->rmap->mapping->indices[i] + 1;
2442cd6534aSBarry Smith     }
2452cd6534aSBarry Smith   }
246e32f2f54SBarry Smith   if (nt != n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hmm nt %D n %D",nt,n);
247854ce69bSBarry Smith   ierr = PetscMalloc1(n+1,&auglyrmapd);CHKERRQ(ierr);
248992144d0SBarry Smith   for (i=0; i<inA->rmap->mapping->n; i++) {
2492cd6534aSBarry Smith     if (r_rmapd[i]) {
2502cd6534aSBarry Smith       auglyrmapd[(r_rmapd[i]-1)-cstart] = i;
2512cd6534aSBarry Smith     }
2522cd6534aSBarry Smith   }
2532cd6534aSBarry Smith   ierr = PetscFree(r_rmapd);CHKERRQ(ierr);
2542cd6534aSBarry Smith   ierr = VecCreateSeq(PETSC_COMM_SELF,n,&auglydd);CHKERRQ(ierr);
2552cd6534aSBarry Smith 
2561795a4d1SJed Brown   ierr = PetscCalloc1(inA->cmap->N+1,&lindices);CHKERRQ(ierr);
257d0f46423SBarry Smith   for (i=0; i<ina->B->cmap->n; i++) {
2582cd6534aSBarry Smith     lindices[garray[i]] = i+1;
2592cd6534aSBarry Smith   }
260992144d0SBarry Smith   no   = inA->rmap->mapping->n - nt;
2611795a4d1SJed Brown   ierr = PetscCalloc1(inA->rmap->mapping->n+1,&r_rmapo);CHKERRQ(ierr);
2622cd6534aSBarry Smith   nt   = 0;
263992144d0SBarry Smith   for (i=0; i<inA->rmap->mapping->n; i++) {
264992144d0SBarry Smith     if (lindices[inA->rmap->mapping->indices[i]]) {
2652cd6534aSBarry Smith       nt++;
266992144d0SBarry Smith       r_rmapo[i] = lindices[inA->rmap->mapping->indices[i]];
2672cd6534aSBarry Smith     }
2682cd6534aSBarry Smith   }
269e32f2f54SBarry Smith   if (nt > no) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hmm nt %D no %D",nt,n);
2702cd6534aSBarry Smith   ierr = PetscFree(lindices);CHKERRQ(ierr);
271854ce69bSBarry Smith   ierr = PetscMalloc1(nt+1,&auglyrmapo);CHKERRQ(ierr);
272992144d0SBarry Smith   for (i=0; i<inA->rmap->mapping->n; i++) {
2732cd6534aSBarry Smith     if (r_rmapo[i]) {
2742cd6534aSBarry Smith       auglyrmapo[(r_rmapo[i]-1)] = i;
2752cd6534aSBarry Smith     }
2762cd6534aSBarry Smith   }
2772cd6534aSBarry Smith   ierr = PetscFree(r_rmapo);CHKERRQ(ierr);
2782cd6534aSBarry Smith   ierr = VecCreateSeq(PETSC_COMM_SELF,nt,&auglyoo);CHKERRQ(ierr);
2792cd6534aSBarry Smith   PetscFunctionReturn(0);
2802cd6534aSBarry Smith }
2812cd6534aSBarry Smith 
282dfbe8321SBarry Smith PetscErrorCode MatMPIAIJDiagonalScaleLocal(Mat A,Vec scale)
2832cd6534aSBarry Smith {
28492b32695SKris Buschelman   /* This routine should really be abandoned as it duplicates MatDiagonalScaleLocal */
2854ac538c5SBarry Smith   PetscErrorCode ierr;
28692b32695SKris Buschelman 
28792b32695SKris Buschelman   PetscFunctionBegin;
2884ac538c5SBarry Smith   ierr = PetscTryMethod(A,"MatDiagonalScaleLocal_C",(Mat,Vec),(A,scale));CHKERRQ(ierr);
28992b32695SKris Buschelman   PetscFunctionReturn(0);
29092b32695SKris Buschelman }
29192b32695SKris Buschelman 
2927087cfbeSBarry Smith PetscErrorCode  MatDiagonalScaleLocal_MPIAIJ(Mat A,Vec scale)
29392b32695SKris Buschelman {
2942cd6534aSBarry Smith   Mat_MPIAIJ        *a = (Mat_MPIAIJ*) A->data; /*access private part of matrix */
295dfbe8321SBarry Smith   PetscErrorCode    ierr;
296b1d57f15SBarry Smith   PetscInt          n,i;
297bca11509SBarry Smith   PetscScalar       *d,*o;
298bca11509SBarry Smith   const PetscScalar *s;
2992cd6534aSBarry Smith 
3002cd6534aSBarry Smith   PetscFunctionBegin;
3012cd6534aSBarry Smith   if (!auglyrmapd) {
3022cd6534aSBarry Smith     ierr = MatMPIAIJDiagonalScaleLocalSetUp(A,scale);CHKERRQ(ierr);
3032cd6534aSBarry Smith   }
3042cd6534aSBarry Smith 
305bca11509SBarry Smith   ierr = VecGetArrayRead(scale,&s);CHKERRQ(ierr);
3062cd6534aSBarry Smith 
3072cd6534aSBarry Smith   ierr = VecGetLocalSize(auglydd,&n);CHKERRQ(ierr);
3081ebc52fbSHong Zhang   ierr = VecGetArray(auglydd,&d);CHKERRQ(ierr);
3092cd6534aSBarry Smith   for (i=0; i<n; i++) {
3102cd6534aSBarry Smith     d[i] = s[auglyrmapd[i]]; /* copy "diagonal" (true local) portion of scale into dd vector */
3112cd6534aSBarry Smith   }
3121ebc52fbSHong Zhang   ierr = VecRestoreArray(auglydd,&d);CHKERRQ(ierr);
3132cd6534aSBarry Smith   /* column scale "diagonal" portion of local matrix */
3140298fd71SBarry Smith   ierr = MatDiagonalScale(a->A,NULL,auglydd);CHKERRQ(ierr);
3152cd6534aSBarry Smith 
3162cd6534aSBarry Smith   ierr = VecGetLocalSize(auglyoo,&n);CHKERRQ(ierr);
3171ebc52fbSHong Zhang   ierr = VecGetArray(auglyoo,&o);CHKERRQ(ierr);
3182cd6534aSBarry Smith   for (i=0; i<n; i++) {
3192cd6534aSBarry Smith     o[i] = s[auglyrmapo[i]]; /* copy "off-diagonal" portion of scale into oo vector */
3202cd6534aSBarry Smith   }
321bca11509SBarry Smith   ierr = VecRestoreArrayRead(scale,&s);CHKERRQ(ierr);
3221ebc52fbSHong Zhang   ierr = VecRestoreArray(auglyoo,&o);CHKERRQ(ierr);
3232cd6534aSBarry Smith   /* column scale "off-diagonal" portion of local matrix */
3240298fd71SBarry Smith   ierr = MatDiagonalScale(a->B,NULL,auglyoo);CHKERRQ(ierr);
3252cd6534aSBarry Smith   PetscFunctionReturn(0);
3262cd6534aSBarry Smith }
327