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) { 27*7e8381f9SStefano Zampini if (!aij->B) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing B mat"); 28aa482453SBarry Smith #if defined(PETSC_USE_CTABLE) 29c5bfad50SMark F. Adams /* use a table */ 30e23dfa41SBarry Smith ierr = PetscTableCreate(aij->B->rmap->n,mat->cmap->N+1,&gid1_lid1);CHKERRQ(ierr); 31d0f46423SBarry Smith for (i=0; i<aij->B->rmap->n; i++) { 322066d6f7SSatish Balay for (j=0; j<B->ilen[i]; j++) { 33b1d57f15SBarry Smith PetscInt data,gid1 = aj[B->i[i] + j] + 1; 340f5bd95cSBarry Smith ierr = PetscTableFind(gid1_lid1,gid1,&data);CHKERRQ(ierr); 35fa46199cSSatish Balay if (!data) { 362066d6f7SSatish Balay /* one based table */ 373861aac3SJed Brown ierr = PetscTableAdd(gid1_lid1,gid1,++ec,INSERT_VALUES);CHKERRQ(ierr); 382066d6f7SSatish Balay } 392066d6f7SSatish Balay } 402066d6f7SSatish Balay } 412066d6f7SSatish Balay /* form array of columns we need */ 42854ce69bSBarry Smith ierr = PetscMalloc1(ec+1,&garray);CHKERRQ(ierr); 430f5bd95cSBarry Smith ierr = PetscTableGetHeadPosition(gid1_lid1,&tpos);CHKERRQ(ierr); 442066d6f7SSatish Balay while (tpos) { 450f5bd95cSBarry Smith ierr = PetscTableGetNext(gid1_lid1,&tpos,&gid,&lid);CHKERRQ(ierr); 46b0a32e0cSBarry Smith gid--; 47b0a32e0cSBarry Smith lid--; 482066d6f7SSatish Balay garray[lid] = gid; 492066d6f7SSatish Balay } 500064e2bbSSatish Balay ierr = PetscSortInt(ec,garray);CHKERRQ(ierr); /* sort, and rebuild */ 510f5bd95cSBarry Smith ierr = PetscTableRemoveAll(gid1_lid1);CHKERRQ(ierr); 522066d6f7SSatish Balay for (i=0; i<ec; i++) { 533861aac3SJed Brown ierr = PetscTableAdd(gid1_lid1,garray[i]+1,i+1,INSERT_VALUES);CHKERRQ(ierr); 542066d6f7SSatish Balay } 552066d6f7SSatish Balay /* compact out the extra columns in B */ 56d0f46423SBarry Smith for (i=0; i<aij->B->rmap->n; i++) { 572066d6f7SSatish Balay for (j=0; j<B->ilen[i]; j++) { 58b1d57f15SBarry Smith PetscInt gid1 = aj[B->i[i] + j] + 1; 590f5bd95cSBarry Smith ierr = PetscTableFind(gid1_lid1,gid1,&lid);CHKERRQ(ierr); 60fa46199cSSatish Balay lid--; 61b3fb0a6cSMatthew Knepley aj[B->i[i] + j] = lid; 622066d6f7SSatish Balay } 632066d6f7SSatish Balay } 64ca5434daSLawrence Mitchell ierr = PetscLayoutDestroy(&aij->B->cmap);CHKERRQ(ierr); 65ca5434daSLawrence Mitchell ierr = PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)aij->B),ec,ec,1,&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 } 96ca5434daSLawrence Mitchell ierr = PetscLayoutDestroy(&aij->B->cmap);CHKERRQ(ierr); 97ca5434daSLawrence Mitchell ierr = PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)aij->B),ec,ec,1,&aij->B->cmap);CHKERRQ(ierr); 98606d414cSSatish Balay ierr = PetscFree(indices);CHKERRQ(ierr); 992066d6f7SSatish Balay #endif 1004b8d542aSHong Zhang } else { 1014b8d542aSHong Zhang garray = aij->garray; 1024b8d542aSHong Zhang } 1034b8d542aSHong Zhang 1044b8d542aSHong Zhang if (!aij->lvec) { 105*7e8381f9SStefano Zampini if (!aij->B) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing B mat"); 106*7e8381f9SStefano Zampini ierr = MatCreateVecs(aij->B,&aij->lvec,NULL);CHKERRQ(ierr); 1074b8d542aSHong Zhang } 108*7e8381f9SStefano Zampini ierr = VecGetSize(aij->lvec,&ec);CHKERRQ(ierr); 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 */ 1193f6a6bdaSHong Zhang ierr = VecScatterDestroy(&aij->Mvctx);CHKERRQ(ierr); 1209448b7f1SJunchao Zhang ierr = VecScatterCreate(gvec,from,aij->lvec,to,&aij->Mvctx);CHKERRQ(ierr); 1213bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)aij->Mvctx);CHKERRQ(ierr); 1223bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)aij->lvec);CHKERRQ(ierr); 1234b8d542aSHong Zhang ierr = PetscLogObjectMemory((PetscObject)mat,(ec+1)*sizeof(PetscInt));CHKERRQ(ierr); 12467bb5161SHong Zhang aij->garray = garray; 1254b8d542aSHong Zhang 1263bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)from);CHKERRQ(ierr); 1273bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)to);CHKERRQ(ierr); 1282205254eSKarl Rupp 1296bf464f9SBarry Smith ierr = ISDestroy(&from);CHKERRQ(ierr); 1306bf464f9SBarry Smith ierr = ISDestroy(&to);CHKERRQ(ierr); 1316bf464f9SBarry Smith ierr = VecDestroy(&gvec);CHKERRQ(ierr); 1323a40ed3dSBarry Smith PetscFunctionReturn(0); 1338c79f6d3SBarry Smith } 1349e25ed09SBarry Smith 1352493cbb0SBarry Smith /* 1362493cbb0SBarry Smith Takes the local part of an already assembled MPIAIJ matrix 1372493cbb0SBarry Smith and disassembles it. This is to allow new nonzeros into the matrix 1382493cbb0SBarry Smith that require more communication in the matrix vector multiply. 1392493cbb0SBarry Smith Thus certain data-structures must be rebuilt. 1402493cbb0SBarry Smith 1412493cbb0SBarry Smith Kind of slow! But that's what application programmers get when 1422493cbb0SBarry Smith they are sloppy. 1432493cbb0SBarry Smith */ 144ab9863d7SBarry Smith PetscErrorCode MatDisAssemble_MPIAIJ(Mat A) 1452493cbb0SBarry Smith { 1462493cbb0SBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ*)A->data; 1472493cbb0SBarry Smith Mat B = aij->B,Bnew; 148ec8511deSBarry Smith Mat_SeqAIJ *Baij = (Mat_SeqAIJ*)B->data; 149dfbe8321SBarry Smith PetscErrorCode ierr; 150d0f46423SBarry Smith PetscInt i,j,m = B->rmap->n,n = A->cmap->N,col,ct = 0,*garray = aij->garray,*nz,ec; 15187828ca2SBarry Smith PetscScalar v; 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); 176f69a0ea3SMatthew Knepley ierr = MatSetSizes(Bnew,m,n,m,n);CHKERRQ(ierr); 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); 1922493cbb0SBarry Smith for (i=0; i<m; i++) { 193bfec09a0SHong Zhang for (j=Baij->i[i]; j<Baij->i[i+1]; j++) { 194bfec09a0SHong Zhang col = garray[Baij->j[ct]]; 1952493cbb0SBarry Smith v = Baij->a[ct++]; 19683271157SBarry Smith ierr = MatSetValues(Bnew,1,&i,1,&col,&v,B->insertmode);CHKERRQ(ierr); 1972493cbb0SBarry Smith } 1982493cbb0SBarry Smith } 199606d414cSSatish Balay ierr = PetscFree(aij->garray);CHKERRQ(ierr); 2003bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)A,-ec*sizeof(PetscInt));CHKERRQ(ierr); 2016bf464f9SBarry Smith ierr = MatDestroy(&B);CHKERRQ(ierr); 2023bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)A,(PetscObject)Bnew);CHKERRQ(ierr); 2032205254eSKarl Rupp 2042493cbb0SBarry Smith aij->B = Bnew; 205227d817aSBarry Smith A->was_assembled = PETSC_FALSE; 2063a40ed3dSBarry Smith PetscFunctionReturn(0); 2072493cbb0SBarry Smith } 2082493cbb0SBarry Smith 2092cd6534aSBarry Smith /* ugly stuff added for Glenn someday we should fix this up */ 2102cd6534aSBarry Smith 211f4259b30SLisandro Dalcin static PetscInt *auglyrmapd = NULL,*auglyrmapo = NULL; /* mapping from the local ordering to the "diagonal" and "off-diagonal" parts of the local matrix */ 212f4259b30SLisandro Dalcin static Vec auglydd = NULL,auglyoo = NULL; /* work vectors used to scale the two parts of the local matrix */ 2132cd6534aSBarry Smith 2142cd6534aSBarry Smith 215dfbe8321SBarry Smith PetscErrorCode MatMPIAIJDiagonalScaleLocalSetUp(Mat inA,Vec scale) 2162cd6534aSBarry Smith { 2172cd6534aSBarry Smith Mat_MPIAIJ *ina = (Mat_MPIAIJ*) inA->data; /*access private part of matrix */ 218dfbe8321SBarry Smith PetscErrorCode ierr; 219b1d57f15SBarry Smith PetscInt i,n,nt,cstart,cend,no,*garray = ina->garray,*lindices; 220b1d57f15SBarry Smith PetscInt *r_rmapd,*r_rmapo; 2212cd6534aSBarry Smith 2222cd6534aSBarry Smith PetscFunctionBegin; 2232cd6534aSBarry Smith ierr = MatGetOwnershipRange(inA,&cstart,&cend);CHKERRQ(ierr); 2240298fd71SBarry Smith ierr = MatGetSize(ina->A,NULL,&n);CHKERRQ(ierr); 2251795a4d1SJed Brown ierr = PetscCalloc1(inA->rmap->mapping->n+1,&r_rmapd);CHKERRQ(ierr); 2262cd6534aSBarry Smith nt = 0; 227992144d0SBarry Smith for (i=0; i<inA->rmap->mapping->n; i++) { 228992144d0SBarry Smith if (inA->rmap->mapping->indices[i] >= cstart && inA->rmap->mapping->indices[i] < cend) { 2292cd6534aSBarry Smith nt++; 230992144d0SBarry Smith r_rmapd[i] = inA->rmap->mapping->indices[i] + 1; 2312cd6534aSBarry Smith } 2322cd6534aSBarry Smith } 233e32f2f54SBarry Smith if (nt != n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hmm nt %D n %D",nt,n); 234854ce69bSBarry Smith ierr = PetscMalloc1(n+1,&auglyrmapd);CHKERRQ(ierr); 235992144d0SBarry Smith for (i=0; i<inA->rmap->mapping->n; i++) { 2362cd6534aSBarry Smith if (r_rmapd[i]) { 2372cd6534aSBarry Smith auglyrmapd[(r_rmapd[i]-1)-cstart] = i; 2382cd6534aSBarry Smith } 2392cd6534aSBarry Smith } 2402cd6534aSBarry Smith ierr = PetscFree(r_rmapd);CHKERRQ(ierr); 2412cd6534aSBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,n,&auglydd);CHKERRQ(ierr); 2422cd6534aSBarry Smith 2431795a4d1SJed Brown ierr = PetscCalloc1(inA->cmap->N+1,&lindices);CHKERRQ(ierr); 244d0f46423SBarry Smith for (i=0; i<ina->B->cmap->n; i++) { 2452cd6534aSBarry Smith lindices[garray[i]] = i+1; 2462cd6534aSBarry Smith } 247992144d0SBarry Smith no = inA->rmap->mapping->n - nt; 2481795a4d1SJed Brown ierr = PetscCalloc1(inA->rmap->mapping->n+1,&r_rmapo);CHKERRQ(ierr); 2492cd6534aSBarry Smith nt = 0; 250992144d0SBarry Smith for (i=0; i<inA->rmap->mapping->n; i++) { 251992144d0SBarry Smith if (lindices[inA->rmap->mapping->indices[i]]) { 2522cd6534aSBarry Smith nt++; 253992144d0SBarry Smith r_rmapo[i] = lindices[inA->rmap->mapping->indices[i]]; 2542cd6534aSBarry Smith } 2552cd6534aSBarry Smith } 256e32f2f54SBarry Smith if (nt > no) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hmm nt %D no %D",nt,n); 2572cd6534aSBarry Smith ierr = PetscFree(lindices);CHKERRQ(ierr); 258854ce69bSBarry Smith ierr = PetscMalloc1(nt+1,&auglyrmapo);CHKERRQ(ierr); 259992144d0SBarry Smith for (i=0; i<inA->rmap->mapping->n; i++) { 2602cd6534aSBarry Smith if (r_rmapo[i]) { 2612cd6534aSBarry Smith auglyrmapo[(r_rmapo[i]-1)] = i; 2622cd6534aSBarry Smith } 2632cd6534aSBarry Smith } 2642cd6534aSBarry Smith ierr = PetscFree(r_rmapo);CHKERRQ(ierr); 2652cd6534aSBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,nt,&auglyoo);CHKERRQ(ierr); 2662cd6534aSBarry Smith PetscFunctionReturn(0); 2672cd6534aSBarry Smith } 2682cd6534aSBarry Smith 269dfbe8321SBarry Smith PetscErrorCode MatMPIAIJDiagonalScaleLocal(Mat A,Vec scale) 2702cd6534aSBarry Smith { 27192b32695SKris Buschelman /* This routine should really be abandoned as it duplicates MatDiagonalScaleLocal */ 2724ac538c5SBarry Smith PetscErrorCode ierr; 27392b32695SKris Buschelman 27492b32695SKris Buschelman PetscFunctionBegin; 2754ac538c5SBarry Smith ierr = PetscTryMethod(A,"MatDiagonalScaleLocal_C",(Mat,Vec),(A,scale));CHKERRQ(ierr); 27692b32695SKris Buschelman PetscFunctionReturn(0); 27792b32695SKris Buschelman } 27892b32695SKris Buschelman 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; 284bca11509SBarry Smith PetscScalar *d,*o; 285bca11509SBarry Smith const PetscScalar *s; 2862cd6534aSBarry Smith 2872cd6534aSBarry Smith PetscFunctionBegin; 2882cd6534aSBarry Smith if (!auglyrmapd) { 2892cd6534aSBarry Smith ierr = MatMPIAIJDiagonalScaleLocalSetUp(A,scale);CHKERRQ(ierr); 2902cd6534aSBarry Smith } 2912cd6534aSBarry Smith 292bca11509SBarry Smith ierr = VecGetArrayRead(scale,&s);CHKERRQ(ierr); 2932cd6534aSBarry Smith 2942cd6534aSBarry Smith ierr = VecGetLocalSize(auglydd,&n);CHKERRQ(ierr); 2951ebc52fbSHong Zhang ierr = VecGetArray(auglydd,&d);CHKERRQ(ierr); 2962cd6534aSBarry Smith for (i=0; i<n; i++) { 2972cd6534aSBarry Smith d[i] = s[auglyrmapd[i]]; /* copy "diagonal" (true local) portion of scale into dd vector */ 2982cd6534aSBarry Smith } 2991ebc52fbSHong Zhang ierr = VecRestoreArray(auglydd,&d);CHKERRQ(ierr); 3002cd6534aSBarry Smith /* column scale "diagonal" portion of local matrix */ 3010298fd71SBarry Smith ierr = MatDiagonalScale(a->A,NULL,auglydd);CHKERRQ(ierr); 3022cd6534aSBarry Smith 3032cd6534aSBarry Smith ierr = VecGetLocalSize(auglyoo,&n);CHKERRQ(ierr); 3041ebc52fbSHong Zhang ierr = VecGetArray(auglyoo,&o);CHKERRQ(ierr); 3052cd6534aSBarry Smith for (i=0; i<n; i++) { 3062cd6534aSBarry Smith o[i] = s[auglyrmapo[i]]; /* copy "off-diagonal" portion of scale into oo vector */ 3072cd6534aSBarry Smith } 308bca11509SBarry Smith ierr = VecRestoreArrayRead(scale,&s);CHKERRQ(ierr); 3091ebc52fbSHong Zhang ierr = VecRestoreArray(auglyoo,&o);CHKERRQ(ierr); 3102cd6534aSBarry Smith /* column scale "off-diagonal" portion of local matrix */ 3110298fd71SBarry Smith ierr = MatDiagonalScale(a->B,NULL,auglyoo);CHKERRQ(ierr); 3122cd6534aSBarry Smith PetscFunctionReturn(0); 3132cd6534aSBarry Smith } 314