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); 136849ba73SBarry Smith PetscErrorCode ierr; 14b1d57f15SBarry Smith PetscInt i,j,*aj = B->j,ec = 0,*garray; 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) { 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 */ 41854ce69bSBarry Smith 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 } 63*ca5434daSLawrence Mitchell ierr = PetscLayoutDestroy(&aij->B->cmap);CHKERRQ(ierr); 64*ca5434daSLawrence Mitchell ierr = PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)aij->B),ec,ec,1,&aij->B->cmap);CHKERRQ(ierr); 656bc0bbbfSBarry Smith ierr = PetscTableDestroy(&gid1_lid1);CHKERRQ(ierr); 662066d6f7SSatish Balay #else 6711285404SBarry Smith /* Make an array as long as the number of columns */ 681eb62cbbSBarry Smith /* mark those columns that are in aij->B */ 691795a4d1SJed Brown ierr = PetscCalloc1(N+1,&indices);CHKERRQ(ierr); 70d0f46423SBarry Smith for (i=0; i<aij->B->rmap->n; i++) { 71d6dfbf8fSBarry Smith for (j=0; j<B->ilen[i]; j++) { 72bfec09a0SHong Zhang if (!indices[aj[B->i[i] + j]]) ec++; 73bfec09a0SHong Zhang indices[aj[B->i[i] + j]] = 1; 74416022c9SBarry Smith } 751eb62cbbSBarry Smith } 768c79f6d3SBarry Smith 771eb62cbbSBarry Smith /* form array of columns we need */ 78854ce69bSBarry Smith ierr = PetscMalloc1(ec+1,&garray);CHKERRQ(ierr); 791eb62cbbSBarry Smith ec = 0; 801eb62cbbSBarry Smith for (i=0; i<N; i++) { 811eb62cbbSBarry Smith if (indices[i]) garray[ec++] = i; 821eb62cbbSBarry Smith } 831eb62cbbSBarry Smith 841eb62cbbSBarry Smith /* make indices now point into garray */ 851eb62cbbSBarry Smith for (i=0; i<ec; i++) { 86bfec09a0SHong Zhang indices[garray[i]] = i; 871eb62cbbSBarry Smith } 881eb62cbbSBarry Smith 891eb62cbbSBarry Smith /* compact out the extra columns in B */ 90d0f46423SBarry Smith for (i=0; i<aij->B->rmap->n; i++) { 91d6dfbf8fSBarry Smith for (j=0; j<B->ilen[i]; j++) { 92bfec09a0SHong Zhang aj[B->i[i] + j] = indices[aj[B->i[i] + j]]; 931eb62cbbSBarry Smith } 94d6dfbf8fSBarry Smith } 95*ca5434daSLawrence Mitchell ierr = PetscLayoutDestroy(&aij->B->cmap);CHKERRQ(ierr); 96*ca5434daSLawrence Mitchell ierr = PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)aij->B),ec,ec,1,&aij->B->cmap);CHKERRQ(ierr); 97606d414cSSatish Balay ierr = PetscFree(indices);CHKERRQ(ierr); 982066d6f7SSatish Balay #endif 994b8d542aSHong Zhang } else { 1004b8d542aSHong Zhang garray = aij->garray; 1014b8d542aSHong Zhang } 1024b8d542aSHong Zhang 1034b8d542aSHong Zhang if (!aij->lvec) { 1041eb62cbbSBarry Smith /* create local vector that is used to scatter into */ 105029af93fSBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,ec,&aij->lvec);CHKERRQ(ierr); 1064b8d542aSHong Zhang } else { 1074b8d542aSHong Zhang ierr = VecGetSize(aij->lvec,&ec);CHKERRQ(ierr); 1084b8d542aSHong Zhang } 1091eb62cbbSBarry Smith 110d6dfbf8fSBarry Smith /* create two temporary Index sets for build scatter gather */ 111a3ebf921SJunchao Zhang ierr = ISCreateGeneral(PETSC_COMM_SELF,ec,garray,PETSC_COPY_VALUES,&from);CHKERRQ(ierr); 112029af93fSBarry Smith ierr = ISCreateStride(PETSC_COMM_SELF,ec,0,1,&to);CHKERRQ(ierr); 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 */ 116ce94432eSBarry Smith ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)mat),1,mat->cmap->n,mat->cmap->N,NULL,&gvec);CHKERRQ(ierr); 1171eb62cbbSBarry Smith 1182d336d48SLois Curfman McInnes /* generate the scatter context */ 119a78d8160SHong Zhang if (aij->Mvctx_mpi1_flg) { 12001ad2aeeSHong Zhang ierr = VecScatterDestroy(&aij->Mvctx_mpi1);CHKERRQ(ierr); 1219448b7f1SJunchao Zhang ierr = VecScatterCreate(gvec,from,aij->lvec,to,&aij->Mvctx_mpi1);CHKERRQ(ierr); 122803a1b88SHong Zhang ierr = VecScatterSetType(aij->Mvctx_mpi1,VECSCATTERMPI1);CHKERRQ(ierr); 1234b8d542aSHong Zhang ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)aij->Mvctx_mpi1);CHKERRQ(ierr); 1244b8d542aSHong Zhang } else { 1253f6a6bdaSHong Zhang ierr = VecScatterDestroy(&aij->Mvctx);CHKERRQ(ierr); 1269448b7f1SJunchao Zhang ierr = VecScatterCreate(gvec,from,aij->lvec,to,&aij->Mvctx);CHKERRQ(ierr); 1273bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)aij->Mvctx);CHKERRQ(ierr); 1283bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)aij->lvec);CHKERRQ(ierr); 1294b8d542aSHong Zhang ierr = PetscLogObjectMemory((PetscObject)mat,(ec+1)*sizeof(PetscInt));CHKERRQ(ierr); 1304b8d542aSHong Zhang } 13167bb5161SHong Zhang aij->garray = garray; 1324b8d542aSHong Zhang 1333bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)from);CHKERRQ(ierr); 1343bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)to);CHKERRQ(ierr); 1352205254eSKarl Rupp 1366bf464f9SBarry Smith ierr = ISDestroy(&from);CHKERRQ(ierr); 1376bf464f9SBarry Smith ierr = ISDestroy(&to);CHKERRQ(ierr); 1386bf464f9SBarry Smith ierr = VecDestroy(&gvec);CHKERRQ(ierr); 1393a40ed3dSBarry Smith PetscFunctionReturn(0); 1408c79f6d3SBarry Smith } 1419e25ed09SBarry Smith 1422493cbb0SBarry Smith /* 1432493cbb0SBarry Smith Takes the local part of an already assembled MPIAIJ matrix 1442493cbb0SBarry Smith and disassembles it. This is to allow new nonzeros into the matrix 1452493cbb0SBarry Smith that require more communication in the matrix vector multiply. 1462493cbb0SBarry Smith Thus certain data-structures must be rebuilt. 1472493cbb0SBarry Smith 1482493cbb0SBarry Smith Kind of slow! But that's what application programmers get when 1492493cbb0SBarry Smith they are sloppy. 1502493cbb0SBarry Smith */ 151ab9863d7SBarry Smith PetscErrorCode MatDisAssemble_MPIAIJ(Mat A) 1522493cbb0SBarry Smith { 1532493cbb0SBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ*)A->data; 1542493cbb0SBarry Smith Mat B = aij->B,Bnew; 155ec8511deSBarry Smith Mat_SeqAIJ *Baij = (Mat_SeqAIJ*)B->data; 156dfbe8321SBarry Smith PetscErrorCode ierr; 157d0f46423SBarry Smith PetscInt i,j,m = B->rmap->n,n = A->cmap->N,col,ct = 0,*garray = aij->garray,*nz,ec; 15887828ca2SBarry Smith PetscScalar v; 1592493cbb0SBarry Smith 1603a40ed3dSBarry Smith PetscFunctionBegin; 1612493cbb0SBarry Smith /* free stuff related to matrix-vec multiply */ 162b0a32e0cSBarry Smith ierr = VecGetSize(aij->lvec,&ec);CHKERRQ(ierr); /* needed for PetscLogObjectMemory below */ 1635e1f6667SBarry Smith ierr = VecDestroy(&aij->lvec);CHKERRQ(ierr); 164464493b3SBarry Smith if (aij->colmap) { 165aa482453SBarry Smith #if defined(PETSC_USE_CTABLE) 1666bc0bbbfSBarry Smith ierr = PetscTableDestroy(&aij->colmap);CHKERRQ(ierr); 1672066d6f7SSatish Balay #else 168606d414cSSatish Balay ierr = PetscFree(aij->colmap);CHKERRQ(ierr); 1693bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)A,-aij->B->cmap->n*sizeof(PetscInt));CHKERRQ(ierr); 1702066d6f7SSatish Balay #endif 171464493b3SBarry Smith } 1722493cbb0SBarry Smith 1732493cbb0SBarry Smith /* make sure that B is assembled so we can access its values */ 1746d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 175fe2f2677SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1762493cbb0SBarry Smith 1772493cbb0SBarry Smith /* invent new B and copy stuff over */ 178854ce69bSBarry Smith ierr = PetscMalloc1(m+1,&nz);CHKERRQ(ierr); 17948b35521SBarry Smith for (i=0; i<m; i++) { 18048b35521SBarry Smith nz[i] = Baij->i[i+1] - Baij->i[i]; 18148b35521SBarry Smith } 182f69a0ea3SMatthew Knepley ierr = MatCreate(PETSC_COMM_SELF,&Bnew);CHKERRQ(ierr); 183f69a0ea3SMatthew Knepley ierr = MatSetSizes(Bnew,m,n,m,n);CHKERRQ(ierr); 18433d57670SJed Brown ierr = MatSetBlockSizesFromMats(Bnew,A,A);CHKERRQ(ierr); 1857adad957SLisandro Dalcin ierr = MatSetType(Bnew,((PetscObject)B)->type_name);CHKERRQ(ierr); 186f204ca49SKris Buschelman ierr = MatSeqAIJSetPreallocation(Bnew,0,nz);CHKERRQ(ierr); 1872205254eSKarl Rupp 188b38c15b3SStefano Zampini if (Baij->nonew >= 0) { /* Inherit insertion error options (if positive). */ 189b38c15b3SStefano Zampini ((Mat_SeqAIJ*)Bnew->data)->nonew = Baij->nonew; 190b38c15b3SStefano Zampini } 191b38c15b3SStefano Zampini 19277341eacSDmitry Karpeev /* 19377341eacSDmitry Karpeev Ensure that B's nonzerostate is monotonically increasing. 19477341eacSDmitry Karpeev Or should this follow the MatSetValues() loop to preserve B's nonzerstate across a MatDisAssemble() call? 19577341eacSDmitry Karpeev */ 196f69fde56SShane Stafford Bnew->nonzerostate = B->nonzerostate; 1972205254eSKarl Rupp 198606d414cSSatish Balay ierr = PetscFree(nz);CHKERRQ(ierr); 1992493cbb0SBarry Smith for (i=0; i<m; i++) { 200bfec09a0SHong Zhang for (j=Baij->i[i]; j<Baij->i[i+1]; j++) { 201bfec09a0SHong Zhang col = garray[Baij->j[ct]]; 2022493cbb0SBarry Smith v = Baij->a[ct++]; 20383271157SBarry Smith ierr = MatSetValues(Bnew,1,&i,1,&col,&v,B->insertmode);CHKERRQ(ierr); 2042493cbb0SBarry Smith } 2052493cbb0SBarry Smith } 206606d414cSSatish Balay ierr = PetscFree(aij->garray);CHKERRQ(ierr); 2073bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)A,-ec*sizeof(PetscInt));CHKERRQ(ierr); 2086bf464f9SBarry Smith ierr = MatDestroy(&B);CHKERRQ(ierr); 2093bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)A,(PetscObject)Bnew);CHKERRQ(ierr); 2102205254eSKarl Rupp 2112493cbb0SBarry Smith aij->B = Bnew; 212227d817aSBarry Smith A->was_assembled = PETSC_FALSE; 2133a40ed3dSBarry Smith PetscFunctionReturn(0); 2142493cbb0SBarry Smith } 2152493cbb0SBarry Smith 2162cd6534aSBarry Smith /* ugly stuff added for Glenn someday we should fix this up */ 2172cd6534aSBarry Smith 2182205254eSKarl Rupp static PetscInt *auglyrmapd = 0,*auglyrmapo = 0; /* mapping from the local ordering to the "diagonal" and "off-diagonal" parts of the local matrix */ 2192cd6534aSBarry Smith static Vec auglydd = 0,auglyoo = 0; /* work vectors used to scale the two parts of the local matrix */ 2202cd6534aSBarry Smith 2212cd6534aSBarry Smith 222dfbe8321SBarry Smith PetscErrorCode MatMPIAIJDiagonalScaleLocalSetUp(Mat inA,Vec scale) 2232cd6534aSBarry Smith { 2242cd6534aSBarry Smith Mat_MPIAIJ *ina = (Mat_MPIAIJ*) inA->data; /*access private part of matrix */ 225dfbe8321SBarry Smith PetscErrorCode ierr; 226b1d57f15SBarry Smith PetscInt i,n,nt,cstart,cend,no,*garray = ina->garray,*lindices; 227b1d57f15SBarry Smith PetscInt *r_rmapd,*r_rmapo; 2282cd6534aSBarry Smith 2292cd6534aSBarry Smith PetscFunctionBegin; 2302cd6534aSBarry Smith ierr = MatGetOwnershipRange(inA,&cstart,&cend);CHKERRQ(ierr); 2310298fd71SBarry Smith ierr = MatGetSize(ina->A,NULL,&n);CHKERRQ(ierr); 2321795a4d1SJed Brown ierr = PetscCalloc1(inA->rmap->mapping->n+1,&r_rmapd);CHKERRQ(ierr); 2332cd6534aSBarry Smith nt = 0; 234992144d0SBarry Smith for (i=0; i<inA->rmap->mapping->n; i++) { 235992144d0SBarry Smith if (inA->rmap->mapping->indices[i] >= cstart && inA->rmap->mapping->indices[i] < cend) { 2362cd6534aSBarry Smith nt++; 237992144d0SBarry Smith r_rmapd[i] = inA->rmap->mapping->indices[i] + 1; 2382cd6534aSBarry Smith } 2392cd6534aSBarry Smith } 240e32f2f54SBarry Smith if (nt != n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hmm nt %D n %D",nt,n); 241854ce69bSBarry Smith ierr = PetscMalloc1(n+1,&auglyrmapd);CHKERRQ(ierr); 242992144d0SBarry Smith for (i=0; i<inA->rmap->mapping->n; i++) { 2432cd6534aSBarry Smith if (r_rmapd[i]) { 2442cd6534aSBarry Smith auglyrmapd[(r_rmapd[i]-1)-cstart] = i; 2452cd6534aSBarry Smith } 2462cd6534aSBarry Smith } 2472cd6534aSBarry Smith ierr = PetscFree(r_rmapd);CHKERRQ(ierr); 2482cd6534aSBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,n,&auglydd);CHKERRQ(ierr); 2492cd6534aSBarry Smith 2501795a4d1SJed Brown ierr = PetscCalloc1(inA->cmap->N+1,&lindices);CHKERRQ(ierr); 251d0f46423SBarry Smith for (i=0; i<ina->B->cmap->n; i++) { 2522cd6534aSBarry Smith lindices[garray[i]] = i+1; 2532cd6534aSBarry Smith } 254992144d0SBarry Smith no = inA->rmap->mapping->n - nt; 2551795a4d1SJed Brown ierr = PetscCalloc1(inA->rmap->mapping->n+1,&r_rmapo);CHKERRQ(ierr); 2562cd6534aSBarry Smith nt = 0; 257992144d0SBarry Smith for (i=0; i<inA->rmap->mapping->n; i++) { 258992144d0SBarry Smith if (lindices[inA->rmap->mapping->indices[i]]) { 2592cd6534aSBarry Smith nt++; 260992144d0SBarry Smith r_rmapo[i] = lindices[inA->rmap->mapping->indices[i]]; 2612cd6534aSBarry Smith } 2622cd6534aSBarry Smith } 263e32f2f54SBarry Smith if (nt > no) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hmm nt %D no %D",nt,n); 2642cd6534aSBarry Smith ierr = PetscFree(lindices);CHKERRQ(ierr); 265854ce69bSBarry Smith ierr = PetscMalloc1(nt+1,&auglyrmapo);CHKERRQ(ierr); 266992144d0SBarry Smith for (i=0; i<inA->rmap->mapping->n; i++) { 2672cd6534aSBarry Smith if (r_rmapo[i]) { 2682cd6534aSBarry Smith auglyrmapo[(r_rmapo[i]-1)] = i; 2692cd6534aSBarry Smith } 2702cd6534aSBarry Smith } 2712cd6534aSBarry Smith ierr = PetscFree(r_rmapo);CHKERRQ(ierr); 2722cd6534aSBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,nt,&auglyoo);CHKERRQ(ierr); 2732cd6534aSBarry Smith PetscFunctionReturn(0); 2742cd6534aSBarry Smith } 2752cd6534aSBarry Smith 276dfbe8321SBarry Smith PetscErrorCode MatMPIAIJDiagonalScaleLocal(Mat A,Vec scale) 2772cd6534aSBarry Smith { 27892b32695SKris Buschelman /* This routine should really be abandoned as it duplicates MatDiagonalScaleLocal */ 2794ac538c5SBarry Smith PetscErrorCode ierr; 28092b32695SKris Buschelman 28192b32695SKris Buschelman PetscFunctionBegin; 2824ac538c5SBarry Smith ierr = PetscTryMethod(A,"MatDiagonalScaleLocal_C",(Mat,Vec),(A,scale));CHKERRQ(ierr); 28392b32695SKris Buschelman PetscFunctionReturn(0); 28492b32695SKris Buschelman } 28592b32695SKris Buschelman 2867087cfbeSBarry Smith PetscErrorCode MatDiagonalScaleLocal_MPIAIJ(Mat A,Vec scale) 28792b32695SKris Buschelman { 2882cd6534aSBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ*) A->data; /*access private part of matrix */ 289dfbe8321SBarry Smith PetscErrorCode ierr; 290b1d57f15SBarry Smith PetscInt n,i; 291bca11509SBarry Smith PetscScalar *d,*o; 292bca11509SBarry Smith const PetscScalar *s; 2932cd6534aSBarry Smith 2942cd6534aSBarry Smith PetscFunctionBegin; 2952cd6534aSBarry Smith if (!auglyrmapd) { 2962cd6534aSBarry Smith ierr = MatMPIAIJDiagonalScaleLocalSetUp(A,scale);CHKERRQ(ierr); 2972cd6534aSBarry Smith } 2982cd6534aSBarry Smith 299bca11509SBarry Smith ierr = VecGetArrayRead(scale,&s);CHKERRQ(ierr); 3002cd6534aSBarry Smith 3012cd6534aSBarry Smith ierr = VecGetLocalSize(auglydd,&n);CHKERRQ(ierr); 3021ebc52fbSHong Zhang ierr = VecGetArray(auglydd,&d);CHKERRQ(ierr); 3032cd6534aSBarry Smith for (i=0; i<n; i++) { 3042cd6534aSBarry Smith d[i] = s[auglyrmapd[i]]; /* copy "diagonal" (true local) portion of scale into dd vector */ 3052cd6534aSBarry Smith } 3061ebc52fbSHong Zhang ierr = VecRestoreArray(auglydd,&d);CHKERRQ(ierr); 3072cd6534aSBarry Smith /* column scale "diagonal" portion of local matrix */ 3080298fd71SBarry Smith ierr = MatDiagonalScale(a->A,NULL,auglydd);CHKERRQ(ierr); 3092cd6534aSBarry Smith 3102cd6534aSBarry Smith ierr = VecGetLocalSize(auglyoo,&n);CHKERRQ(ierr); 3111ebc52fbSHong Zhang ierr = VecGetArray(auglyoo,&o);CHKERRQ(ierr); 3122cd6534aSBarry Smith for (i=0; i<n; i++) { 3132cd6534aSBarry Smith o[i] = s[auglyrmapo[i]]; /* copy "off-diagonal" portion of scale into oo vector */ 3142cd6534aSBarry Smith } 315bca11509SBarry Smith ierr = VecRestoreArrayRead(scale,&s);CHKERRQ(ierr); 3161ebc52fbSHong Zhang ierr = VecRestoreArray(auglyoo,&o);CHKERRQ(ierr); 3172cd6534aSBarry Smith /* column scale "off-diagonal" portion of local matrix */ 3180298fd71SBarry Smith ierr = MatDiagonalScale(a->B,NULL,auglyoo);CHKERRQ(ierr); 3192cd6534aSBarry Smith PetscFunctionReturn(0); 3202cd6534aSBarry Smith } 321