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; 14b3c64f9dSJunchao Zhang PetscInt i,j,*aj = B->j,*garray; 15b3c64f9dSJunchao Zhang PetscInt ec = 0; /* Number of nonzero external columns */ 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; 274b8d542aSHong Zhang if (!aij->garray) { 28*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(!aij->B,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing B mat"); 29aa482453SBarry Smith #if defined(PETSC_USE_CTABLE) 30c5bfad50SMark F. Adams /* use a table */ 31e23dfa41SBarry Smith ierr = PetscTableCreate(aij->B->rmap->n,mat->cmap->N+1,&gid1_lid1);CHKERRQ(ierr); 32d0f46423SBarry Smith for (i=0; i<aij->B->rmap->n; i++) { 332066d6f7SSatish Balay for (j=0; j<B->ilen[i]; j++) { 34b1d57f15SBarry Smith PetscInt data,gid1 = aj[B->i[i] + j] + 1; 350f5bd95cSBarry Smith ierr = PetscTableFind(gid1_lid1,gid1,&data);CHKERRQ(ierr); 36fa46199cSSatish Balay if (!data) { 372066d6f7SSatish Balay /* one based table */ 383861aac3SJed Brown ierr = PetscTableAdd(gid1_lid1,gid1,++ec,INSERT_VALUES);CHKERRQ(ierr); 392066d6f7SSatish Balay } 402066d6f7SSatish Balay } 412066d6f7SSatish Balay } 422066d6f7SSatish Balay /* form array of columns we need */ 43b3c64f9dSJunchao Zhang ierr = PetscMalloc1(ec,&garray);CHKERRQ(ierr); 440f5bd95cSBarry Smith ierr = PetscTableGetHeadPosition(gid1_lid1,&tpos);CHKERRQ(ierr); 452066d6f7SSatish Balay while (tpos) { 460f5bd95cSBarry Smith ierr = PetscTableGetNext(gid1_lid1,&tpos,&gid,&lid);CHKERRQ(ierr); 47b0a32e0cSBarry Smith gid--; 48b0a32e0cSBarry Smith lid--; 492066d6f7SSatish Balay garray[lid] = gid; 502066d6f7SSatish Balay } 510064e2bbSSatish Balay ierr = PetscSortInt(ec,garray);CHKERRQ(ierr); /* sort, and rebuild */ 520f5bd95cSBarry Smith ierr = PetscTableRemoveAll(gid1_lid1);CHKERRQ(ierr); 532066d6f7SSatish Balay for (i=0; i<ec; i++) { 543861aac3SJed Brown ierr = PetscTableAdd(gid1_lid1,garray[i]+1,i+1,INSERT_VALUES);CHKERRQ(ierr); 552066d6f7SSatish Balay } 562066d6f7SSatish Balay /* compact out the extra columns in B */ 57d0f46423SBarry Smith for (i=0; i<aij->B->rmap->n; i++) { 582066d6f7SSatish Balay for (j=0; j<B->ilen[i]; j++) { 59b1d57f15SBarry Smith PetscInt gid1 = aj[B->i[i] + j] + 1; 600f5bd95cSBarry Smith ierr = PetscTableFind(gid1_lid1,gid1,&lid);CHKERRQ(ierr); 61fa46199cSSatish Balay lid--; 62b3fb0a6cSMatthew Knepley aj[B->i[i] + j] = lid; 632066d6f7SSatish Balay } 642066d6f7SSatish Balay } 65ca5434daSLawrence Mitchell ierr = PetscLayoutDestroy(&aij->B->cmap);CHKERRQ(ierr); 66ca5434daSLawrence Mitchell ierr = PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)aij->B),ec,ec,1,&aij->B->cmap);CHKERRQ(ierr); 676bc0bbbfSBarry Smith ierr = PetscTableDestroy(&gid1_lid1);CHKERRQ(ierr); 682066d6f7SSatish Balay #else 6911285404SBarry Smith /* Make an array as long as the number of columns */ 701eb62cbbSBarry Smith /* mark those columns that are in aij->B */ 71b3c64f9dSJunchao Zhang ierr = PetscCalloc1(N,&indices);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 */ 80b3c64f9dSJunchao Zhang ierr = PetscMalloc1(ec,&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 } 97ca5434daSLawrence Mitchell ierr = PetscLayoutDestroy(&aij->B->cmap);CHKERRQ(ierr); 98ca5434daSLawrence Mitchell ierr = PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)aij->B),ec,ec,1,&aij->B->cmap);CHKERRQ(ierr); 99606d414cSSatish Balay ierr = PetscFree(indices);CHKERRQ(ierr); 1002066d6f7SSatish Balay #endif 1014b8d542aSHong Zhang } else { 1024b8d542aSHong Zhang garray = aij->garray; 1034b8d542aSHong Zhang } 1044b8d542aSHong Zhang 1054b8d542aSHong Zhang if (!aij->lvec) { 106*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(!aij->B,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing B mat"); 1077e8381f9SStefano Zampini ierr = MatCreateVecs(aij->B,&aij->lvec,NULL);CHKERRQ(ierr); 1084b8d542aSHong Zhang } 1097e8381f9SStefano Zampini ierr = VecGetSize(aij->lvec,&ec);CHKERRQ(ierr); 1101eb62cbbSBarry Smith 111d6dfbf8fSBarry Smith /* create two temporary Index sets for build scatter gather */ 112a3ebf921SJunchao Zhang ierr = ISCreateGeneral(PETSC_COMM_SELF,ec,garray,PETSC_COPY_VALUES,&from);CHKERRQ(ierr); 113029af93fSBarry Smith ierr = ISCreateStride(PETSC_COMM_SELF,ec,0,1,&to);CHKERRQ(ierr); 1141eb62cbbSBarry Smith 1151eb62cbbSBarry Smith /* create temporary global vector to generate scatter context */ 116b5eb4454SBarry Smith /* This does not allocate the array's memory so is efficient */ 117ce94432eSBarry Smith ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)mat),1,mat->cmap->n,mat->cmap->N,NULL,&gvec);CHKERRQ(ierr); 1181eb62cbbSBarry Smith 1192d336d48SLois Curfman McInnes /* generate the scatter context */ 1203f6a6bdaSHong Zhang ierr = VecScatterDestroy(&aij->Mvctx);CHKERRQ(ierr); 1219448b7f1SJunchao Zhang ierr = VecScatterCreate(gvec,from,aij->lvec,to,&aij->Mvctx);CHKERRQ(ierr); 1222f3b2168SJunchao Zhang ierr = VecScatterViewFromOptions(aij->Mvctx,(PetscObject)mat,"-matmult_vecscatter_view");CHKERRQ(ierr); 1233bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)aij->Mvctx);CHKERRQ(ierr); 1243bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)aij->lvec);CHKERRQ(ierr); 125ad1b2038SJunchao Zhang ierr = PetscLogObjectMemory((PetscObject)mat,ec*sizeof(PetscInt));CHKERRQ(ierr); 12667bb5161SHong Zhang aij->garray = garray; 1274b8d542aSHong Zhang 1283bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)from);CHKERRQ(ierr); 1293bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)to);CHKERRQ(ierr); 1302205254eSKarl Rupp 1316bf464f9SBarry Smith ierr = ISDestroy(&from);CHKERRQ(ierr); 1326bf464f9SBarry Smith ierr = ISDestroy(&to);CHKERRQ(ierr); 1336bf464f9SBarry Smith ierr = VecDestroy(&gvec);CHKERRQ(ierr); 1343a40ed3dSBarry Smith PetscFunctionReturn(0); 1358c79f6d3SBarry Smith } 1369e25ed09SBarry Smith 137fff043a9SJunchao Zhang /* Disassemble the off-diag portion of the MPIAIJXxx matrix. 138fff043a9SJunchao Zhang In other words, change the B from reduced format using local col ids 139fff043a9SJunchao Zhang to expanded format using global col ids, which will make it easier to 140fff043a9SJunchao Zhang insert new nonzeros (with global col ids) into the matrix. 141fff043a9SJunchao Zhang The off-diag B determines communication in the matrix vector multiply. 1422493cbb0SBarry Smith */ 143ab9863d7SBarry Smith PetscErrorCode MatDisAssemble_MPIAIJ(Mat A) 1442493cbb0SBarry Smith { 1452493cbb0SBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ*)A->data; 1462493cbb0SBarry Smith Mat B = aij->B,Bnew; 147ec8511deSBarry Smith Mat_SeqAIJ *Baij = (Mat_SeqAIJ*)B->data; 148dfbe8321SBarry Smith PetscErrorCode ierr; 149d0f46423SBarry Smith PetscInt i,j,m = B->rmap->n,n = A->cmap->N,col,ct = 0,*garray = aij->garray,*nz,ec; 15087828ca2SBarry Smith PetscScalar v; 151fff043a9SJunchao Zhang const PetscScalar *ba; 1522493cbb0SBarry Smith 1533a40ed3dSBarry Smith PetscFunctionBegin; 1542493cbb0SBarry Smith /* free stuff related to matrix-vec multiply */ 155b0a32e0cSBarry Smith ierr = VecGetSize(aij->lvec,&ec);CHKERRQ(ierr); /* needed for PetscLogObjectMemory below */ 1565e1f6667SBarry Smith ierr = VecDestroy(&aij->lvec);CHKERRQ(ierr); 157464493b3SBarry Smith if (aij->colmap) { 158aa482453SBarry Smith #if defined(PETSC_USE_CTABLE) 1596bc0bbbfSBarry Smith ierr = PetscTableDestroy(&aij->colmap);CHKERRQ(ierr); 1602066d6f7SSatish Balay #else 161606d414cSSatish Balay ierr = PetscFree(aij->colmap);CHKERRQ(ierr); 1623bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)A,-aij->B->cmap->n*sizeof(PetscInt));CHKERRQ(ierr); 1632066d6f7SSatish Balay #endif 164464493b3SBarry Smith } 1652493cbb0SBarry Smith 1662493cbb0SBarry Smith /* make sure that B is assembled so we can access its values */ 1676d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 168fe2f2677SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1692493cbb0SBarry Smith 1702493cbb0SBarry Smith /* invent new B and copy stuff over */ 171854ce69bSBarry Smith ierr = PetscMalloc1(m+1,&nz);CHKERRQ(ierr); 17248b35521SBarry Smith for (i=0; i<m; i++) { 17348b35521SBarry Smith nz[i] = Baij->i[i+1] - Baij->i[i]; 17448b35521SBarry Smith } 175f69a0ea3SMatthew Knepley ierr = MatCreate(PETSC_COMM_SELF,&Bnew);CHKERRQ(ierr); 176fff043a9SJunchao Zhang ierr = MatSetSizes(Bnew,m,n,m,n);CHKERRQ(ierr); /* Bnew now uses A->cmap->N as its col size */ 17733d57670SJed Brown ierr = MatSetBlockSizesFromMats(Bnew,A,A);CHKERRQ(ierr); 1787adad957SLisandro Dalcin ierr = MatSetType(Bnew,((PetscObject)B)->type_name);CHKERRQ(ierr); 179f204ca49SKris Buschelman ierr = MatSeqAIJSetPreallocation(Bnew,0,nz);CHKERRQ(ierr); 1802205254eSKarl Rupp 181b38c15b3SStefano Zampini if (Baij->nonew >= 0) { /* Inherit insertion error options (if positive). */ 182b38c15b3SStefano Zampini ((Mat_SeqAIJ*)Bnew->data)->nonew = Baij->nonew; 183b38c15b3SStefano Zampini } 184b38c15b3SStefano Zampini 18577341eacSDmitry Karpeev /* 18677341eacSDmitry Karpeev Ensure that B's nonzerostate is monotonically increasing. 18777341eacSDmitry Karpeev Or should this follow the MatSetValues() loop to preserve B's nonzerstate across a MatDisAssemble() call? 18877341eacSDmitry Karpeev */ 189f69fde56SShane Stafford Bnew->nonzerostate = B->nonzerostate; 1902205254eSKarl Rupp 191606d414cSSatish Balay ierr = PetscFree(nz);CHKERRQ(ierr); 192fff043a9SJunchao Zhang ierr = MatSeqAIJGetArrayRead(B,&ba);CHKERRQ(ierr); 1932493cbb0SBarry Smith for (i=0; i<m; i++) { 194bfec09a0SHong Zhang for (j=Baij->i[i]; j<Baij->i[i+1]; j++) { 195bfec09a0SHong Zhang col = garray[Baij->j[ct]]; 196fff043a9SJunchao Zhang v = ba[ct++]; 19783271157SBarry Smith ierr = MatSetValues(Bnew,1,&i,1,&col,&v,B->insertmode);CHKERRQ(ierr); 1982493cbb0SBarry Smith } 1992493cbb0SBarry Smith } 200fff043a9SJunchao Zhang ierr = MatSeqAIJRestoreArrayRead(B,&ba);CHKERRQ(ierr); 201fff043a9SJunchao Zhang 202606d414cSSatish Balay ierr = PetscFree(aij->garray);CHKERRQ(ierr); 2033bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)A,-ec*sizeof(PetscInt));CHKERRQ(ierr); 2046bf464f9SBarry Smith ierr = MatDestroy(&B);CHKERRQ(ierr); 2053bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)A,(PetscObject)Bnew);CHKERRQ(ierr); 2062205254eSKarl Rupp 2072493cbb0SBarry Smith aij->B = Bnew; 208227d817aSBarry Smith A->was_assembled = PETSC_FALSE; 2093a40ed3dSBarry Smith PetscFunctionReturn(0); 2102493cbb0SBarry Smith } 2112493cbb0SBarry Smith 2122cd6534aSBarry Smith /* ugly stuff added for Glenn someday we should fix this up */ 2132cd6534aSBarry Smith 214f4259b30SLisandro Dalcin static PetscInt *auglyrmapd = NULL,*auglyrmapo = NULL; /* mapping from the local ordering to the "diagonal" and "off-diagonal" parts of the local matrix */ 215f4259b30SLisandro Dalcin static Vec auglydd = NULL,auglyoo = NULL; /* work vectors used to scale the two parts of the local matrix */ 2162cd6534aSBarry Smith 217dfbe8321SBarry Smith PetscErrorCode MatMPIAIJDiagonalScaleLocalSetUp(Mat inA,Vec scale) 2182cd6534aSBarry Smith { 2192cd6534aSBarry Smith Mat_MPIAIJ *ina = (Mat_MPIAIJ*) inA->data; /*access private part of matrix */ 220dfbe8321SBarry Smith PetscErrorCode ierr; 221b1d57f15SBarry Smith PetscInt i,n,nt,cstart,cend,no,*garray = ina->garray,*lindices; 222b1d57f15SBarry Smith PetscInt *r_rmapd,*r_rmapo; 2232cd6534aSBarry Smith 2242cd6534aSBarry Smith PetscFunctionBegin; 2252cd6534aSBarry Smith ierr = MatGetOwnershipRange(inA,&cstart,&cend);CHKERRQ(ierr); 2260298fd71SBarry Smith ierr = MatGetSize(ina->A,NULL,&n);CHKERRQ(ierr); 2271795a4d1SJed Brown ierr = PetscCalloc1(inA->rmap->mapping->n+1,&r_rmapd);CHKERRQ(ierr); 2282cd6534aSBarry Smith nt = 0; 229992144d0SBarry Smith for (i=0; i<inA->rmap->mapping->n; i++) { 230992144d0SBarry Smith if (inA->rmap->mapping->indices[i] >= cstart && inA->rmap->mapping->indices[i] < cend) { 2312cd6534aSBarry Smith nt++; 232992144d0SBarry Smith r_rmapd[i] = inA->rmap->mapping->indices[i] + 1; 2332cd6534aSBarry Smith } 2342cd6534aSBarry Smith } 235*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(nt != n,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hmm nt %" PetscInt_FMT " n %" PetscInt_FMT,nt,n); 236854ce69bSBarry Smith ierr = PetscMalloc1(n+1,&auglyrmapd);CHKERRQ(ierr); 237992144d0SBarry Smith for (i=0; i<inA->rmap->mapping->n; i++) { 2382cd6534aSBarry Smith if (r_rmapd[i]) { 2392cd6534aSBarry Smith auglyrmapd[(r_rmapd[i]-1)-cstart] = i; 2402cd6534aSBarry Smith } 2412cd6534aSBarry Smith } 2422cd6534aSBarry Smith ierr = PetscFree(r_rmapd);CHKERRQ(ierr); 2432cd6534aSBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,n,&auglydd);CHKERRQ(ierr); 2442cd6534aSBarry Smith 2451795a4d1SJed Brown ierr = PetscCalloc1(inA->cmap->N+1,&lindices);CHKERRQ(ierr); 246d0f46423SBarry Smith for (i=0; i<ina->B->cmap->n; i++) { 2472cd6534aSBarry Smith lindices[garray[i]] = i+1; 2482cd6534aSBarry Smith } 249992144d0SBarry Smith no = inA->rmap->mapping->n - nt; 2501795a4d1SJed Brown ierr = PetscCalloc1(inA->rmap->mapping->n+1,&r_rmapo);CHKERRQ(ierr); 2512cd6534aSBarry Smith nt = 0; 252992144d0SBarry Smith for (i=0; i<inA->rmap->mapping->n; i++) { 253992144d0SBarry Smith if (lindices[inA->rmap->mapping->indices[i]]) { 2542cd6534aSBarry Smith nt++; 255992144d0SBarry Smith r_rmapo[i] = lindices[inA->rmap->mapping->indices[i]]; 2562cd6534aSBarry Smith } 2572cd6534aSBarry Smith } 258*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(nt > no,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hmm nt %" PetscInt_FMT " no %" PetscInt_FMT,nt,n); 2592cd6534aSBarry Smith ierr = PetscFree(lindices);CHKERRQ(ierr); 260854ce69bSBarry Smith ierr = PetscMalloc1(nt+1,&auglyrmapo);CHKERRQ(ierr); 261992144d0SBarry Smith for (i=0; i<inA->rmap->mapping->n; i++) { 2622cd6534aSBarry Smith if (r_rmapo[i]) { 2632cd6534aSBarry Smith auglyrmapo[(r_rmapo[i]-1)] = i; 2642cd6534aSBarry Smith } 2652cd6534aSBarry Smith } 2662cd6534aSBarry Smith ierr = PetscFree(r_rmapo);CHKERRQ(ierr); 2672cd6534aSBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,nt,&auglyoo);CHKERRQ(ierr); 2682cd6534aSBarry Smith PetscFunctionReturn(0); 2692cd6534aSBarry Smith } 2702cd6534aSBarry Smith 271dfbe8321SBarry Smith PetscErrorCode MatMPIAIJDiagonalScaleLocal(Mat A,Vec scale) 2722cd6534aSBarry Smith { 27392b32695SKris Buschelman /* This routine should really be abandoned as it duplicates MatDiagonalScaleLocal */ 2744ac538c5SBarry Smith PetscErrorCode ierr; 27592b32695SKris Buschelman 27692b32695SKris Buschelman PetscFunctionBegin; 2774ac538c5SBarry Smith ierr = PetscTryMethod(A,"MatDiagonalScaleLocal_C",(Mat,Vec),(A,scale));CHKERRQ(ierr); 27892b32695SKris Buschelman PetscFunctionReturn(0); 27992b32695SKris Buschelman } 28092b32695SKris Buschelman 2817087cfbeSBarry Smith PetscErrorCode MatDiagonalScaleLocal_MPIAIJ(Mat A,Vec scale) 28292b32695SKris Buschelman { 2832cd6534aSBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ*) A->data; /*access private part of matrix */ 284dfbe8321SBarry Smith PetscErrorCode ierr; 285b1d57f15SBarry Smith PetscInt n,i; 286bca11509SBarry Smith PetscScalar *d,*o; 287bca11509SBarry Smith const PetscScalar *s; 2882cd6534aSBarry Smith 2892cd6534aSBarry Smith PetscFunctionBegin; 2902cd6534aSBarry Smith if (!auglyrmapd) { 2912cd6534aSBarry Smith ierr = MatMPIAIJDiagonalScaleLocalSetUp(A,scale);CHKERRQ(ierr); 2922cd6534aSBarry Smith } 2932cd6534aSBarry Smith 294bca11509SBarry Smith ierr = VecGetArrayRead(scale,&s);CHKERRQ(ierr); 2952cd6534aSBarry Smith 2962cd6534aSBarry Smith ierr = VecGetLocalSize(auglydd,&n);CHKERRQ(ierr); 2971ebc52fbSHong Zhang ierr = VecGetArray(auglydd,&d);CHKERRQ(ierr); 2982cd6534aSBarry Smith for (i=0; i<n; i++) { 2992cd6534aSBarry Smith d[i] = s[auglyrmapd[i]]; /* copy "diagonal" (true local) portion of scale into dd vector */ 3002cd6534aSBarry Smith } 3011ebc52fbSHong Zhang ierr = VecRestoreArray(auglydd,&d);CHKERRQ(ierr); 3022cd6534aSBarry Smith /* column scale "diagonal" portion of local matrix */ 3030298fd71SBarry Smith ierr = MatDiagonalScale(a->A,NULL,auglydd);CHKERRQ(ierr); 3042cd6534aSBarry Smith 3052cd6534aSBarry Smith ierr = VecGetLocalSize(auglyoo,&n);CHKERRQ(ierr); 3061ebc52fbSHong Zhang ierr = VecGetArray(auglyoo,&o);CHKERRQ(ierr); 3072cd6534aSBarry Smith for (i=0; i<n; i++) { 3082cd6534aSBarry Smith o[i] = s[auglyrmapo[i]]; /* copy "off-diagonal" portion of scale into oo vector */ 3092cd6534aSBarry Smith } 310bca11509SBarry Smith ierr = VecRestoreArrayRead(scale,&s);CHKERRQ(ierr); 3111ebc52fbSHong Zhang ierr = VecRestoreArray(auglyoo,&o);CHKERRQ(ierr); 3122cd6534aSBarry Smith /* column scale "off-diagonal" portion of local matrix */ 3130298fd71SBarry Smith ierr = MatDiagonalScale(a->B,NULL,auglyoo);CHKERRQ(ierr); 3142cd6534aSBarry Smith PetscFunctionReturn(0); 3152cd6534aSBarry Smith } 316