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 */ 41785e854fSJed Brown 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 } 63d0f46423SBarry Smith aij->B->cmap->n = aij->B->cmap->N = ec; 64cdce4254SBarry Smith aij->B->cmap->bs = 1; 652205254eSKarl Rupp 6626283091SBarry Smith ierr = PetscLayoutSetUp((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 */ 711795a4d1SJed Brown ierr = PetscCalloc1(N+1,&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 */ 80785e854fSJed Brown ierr = PetscMalloc1((ec+1),&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; 98cd0e7f71SBarry Smith aij->B->cmap->bs = 1; 992205254eSKarl Rupp 10026283091SBarry Smith ierr = PetscLayoutSetUp((aij->B->cmap));CHKERRQ(ierr); 101606d414cSSatish Balay ierr = PetscFree(indices);CHKERRQ(ierr); 1022066d6f7SSatish Balay #endif 1031eb62cbbSBarry Smith /* create local vector that is used to scatter into */ 104029af93fSBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,ec,&aij->lvec);CHKERRQ(ierr); 1051eb62cbbSBarry Smith 106d6dfbf8fSBarry Smith /* create two temporary Index sets for build scatter gather */ 10770b3c8c7SBarry Smith ierr = ISCreateGeneral(((PetscObject)mat)->comm,ec,garray,PETSC_COPY_VALUES,&from);CHKERRQ(ierr); 108c5bfad50SMark F. Adams 109029af93fSBarry Smith ierr = ISCreateStride(PETSC_COMM_SELF,ec,0,1,&to);CHKERRQ(ierr); 1101eb62cbbSBarry Smith 1111eb62cbbSBarry Smith /* create temporary global vector to generate scatter context */ 112b5eb4454SBarry Smith /* This does not allocate the array's memory so is efficient */ 113ce94432eSBarry Smith ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)mat),1,mat->cmap->n,mat->cmap->N,NULL,&gvec);CHKERRQ(ierr); 1141eb62cbbSBarry Smith 1152d336d48SLois Curfman McInnes /* generate the scatter context */ 11608480c60SBarry Smith ierr = VecScatterCreate(gvec,from,aij->lvec,to,&aij->Mvctx);CHKERRQ(ierr); 1173bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)aij->Mvctx);CHKERRQ(ierr); 1183bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)aij->lvec);CHKERRQ(ierr); 1193bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)from);CHKERRQ(ierr); 1203bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)to);CHKERRQ(ierr); 1212205254eSKarl Rupp 1229e25ed09SBarry Smith aij->garray = garray; 1232205254eSKarl Rupp 1243bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)mat,(ec+1)*sizeof(PetscInt));CHKERRQ(ierr); 1256bf464f9SBarry Smith ierr = ISDestroy(&from);CHKERRQ(ierr); 1266bf464f9SBarry Smith ierr = ISDestroy(&to);CHKERRQ(ierr); 1276bf464f9SBarry Smith ierr = VecDestroy(&gvec);CHKERRQ(ierr); 1283a40ed3dSBarry Smith PetscFunctionReturn(0); 1298c79f6d3SBarry Smith } 1309e25ed09SBarry Smith 1319e25ed09SBarry Smith 1324a2ae208SSatish Balay #undef __FUNCT__ 133ab9863d7SBarry Smith #define __FUNCT__ "MatDisAssemble_MPIAIJ" 1342493cbb0SBarry Smith /* 1352493cbb0SBarry Smith Takes the local part of an already assembled MPIAIJ matrix 1362493cbb0SBarry Smith and disassembles it. This is to allow new nonzeros into the matrix 1372493cbb0SBarry Smith that require more communication in the matrix vector multiply. 1382493cbb0SBarry Smith Thus certain data-structures must be rebuilt. 1392493cbb0SBarry Smith 1402493cbb0SBarry Smith Kind of slow! But that's what application programmers get when 1412493cbb0SBarry Smith they are sloppy. 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; 1512493cbb0SBarry Smith 1523a40ed3dSBarry Smith PetscFunctionBegin; 1532493cbb0SBarry Smith /* free stuff related to matrix-vec multiply */ 154b0a32e0cSBarry Smith ierr = VecGetSize(aij->lvec,&ec);CHKERRQ(ierr); /* needed for PetscLogObjectMemory below */ 1555e1f6667SBarry Smith ierr = VecDestroy(&aij->lvec);CHKERRQ(ierr); 1565e1f6667SBarry Smith ierr = VecScatterDestroy(&aij->Mvctx);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 */ 171785e854fSJed Brown 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 1812576faa2SJed Brown ((Mat_SeqAIJ*)Bnew->data)->nonew = Baij->nonew; /* Inherit insertion error options. */ 182*77341eacSDmitry Karpeev /* 183*77341eacSDmitry Karpeev Ensure that B's nonzerostate is monotonically increasing. 184*77341eacSDmitry Karpeev Or should this follow the MatSetValues() loop to preserve B's nonzerstate across a MatDisAssemble() call? 185*77341eacSDmitry Karpeev */ 186f69fde56SShane Stafford Bnew->nonzerostate = B->nonzerostate; 1872205254eSKarl Rupp 188606d414cSSatish Balay ierr = PetscFree(nz);CHKERRQ(ierr); 1892493cbb0SBarry Smith for (i=0; i<m; i++) { 190bfec09a0SHong Zhang for (j=Baij->i[i]; j<Baij->i[i+1]; j++) { 191bfec09a0SHong Zhang col = garray[Baij->j[ct]]; 1922493cbb0SBarry Smith v = Baij->a[ct++]; 19383271157SBarry Smith ierr = MatSetValues(Bnew,1,&i,1,&col,&v,B->insertmode);CHKERRQ(ierr); 1942493cbb0SBarry Smith } 1952493cbb0SBarry Smith } 196606d414cSSatish Balay ierr = PetscFree(aij->garray);CHKERRQ(ierr); 1973bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)A,-ec*sizeof(PetscInt));CHKERRQ(ierr); 1986bf464f9SBarry Smith ierr = MatDestroy(&B);CHKERRQ(ierr); 1993bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)A,(PetscObject)Bnew);CHKERRQ(ierr); 2002205254eSKarl Rupp 2012493cbb0SBarry Smith aij->B = Bnew; 202227d817aSBarry Smith A->was_assembled = PETSC_FALSE; 2033a40ed3dSBarry Smith PetscFunctionReturn(0); 2042493cbb0SBarry Smith } 2052493cbb0SBarry Smith 2062cd6534aSBarry Smith /* ugly stuff added for Glenn someday we should fix this up */ 2072cd6534aSBarry Smith 2082205254eSKarl Rupp static PetscInt *auglyrmapd = 0,*auglyrmapo = 0; /* mapping from the local ordering to the "diagonal" and "off-diagonal" parts of the local matrix */ 2092cd6534aSBarry Smith static Vec auglydd = 0,auglyoo = 0; /* work vectors used to scale the two parts of the local matrix */ 2102cd6534aSBarry Smith 2112cd6534aSBarry Smith 2122cd6534aSBarry Smith #undef __FUNCT__ 2132cd6534aSBarry Smith #define __FUNCT__ "MatMPIAIJDiagonalScaleLocalSetUp" 214dfbe8321SBarry Smith PetscErrorCode MatMPIAIJDiagonalScaleLocalSetUp(Mat inA,Vec scale) 2152cd6534aSBarry Smith { 2162cd6534aSBarry Smith Mat_MPIAIJ *ina = (Mat_MPIAIJ*) inA->data; /*access private part of matrix */ 217dfbe8321SBarry Smith PetscErrorCode ierr; 218b1d57f15SBarry Smith PetscInt i,n,nt,cstart,cend,no,*garray = ina->garray,*lindices; 219b1d57f15SBarry Smith PetscInt *r_rmapd,*r_rmapo; 2202cd6534aSBarry Smith 2212cd6534aSBarry Smith PetscFunctionBegin; 2222cd6534aSBarry Smith ierr = MatGetOwnershipRange(inA,&cstart,&cend);CHKERRQ(ierr); 2230298fd71SBarry Smith ierr = MatGetSize(ina->A,NULL,&n);CHKERRQ(ierr); 2241795a4d1SJed Brown ierr = PetscCalloc1(inA->rmap->mapping->n+1,&r_rmapd);CHKERRQ(ierr); 2252cd6534aSBarry Smith nt = 0; 226992144d0SBarry Smith for (i=0; i<inA->rmap->mapping->n; i++) { 227992144d0SBarry Smith if (inA->rmap->mapping->indices[i] >= cstart && inA->rmap->mapping->indices[i] < cend) { 2282cd6534aSBarry Smith nt++; 229992144d0SBarry Smith r_rmapd[i] = inA->rmap->mapping->indices[i] + 1; 2302cd6534aSBarry Smith } 2312cd6534aSBarry Smith } 232e32f2f54SBarry Smith if (nt != n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hmm nt %D n %D",nt,n); 233785e854fSJed Brown ierr = PetscMalloc1((n+1),&auglyrmapd);CHKERRQ(ierr); 234992144d0SBarry Smith for (i=0; i<inA->rmap->mapping->n; i++) { 2352cd6534aSBarry Smith if (r_rmapd[i]) { 2362cd6534aSBarry Smith auglyrmapd[(r_rmapd[i]-1)-cstart] = i; 2372cd6534aSBarry Smith } 2382cd6534aSBarry Smith } 2392cd6534aSBarry Smith ierr = PetscFree(r_rmapd);CHKERRQ(ierr); 2402cd6534aSBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,n,&auglydd);CHKERRQ(ierr); 2412cd6534aSBarry Smith 2421795a4d1SJed Brown ierr = PetscCalloc1(inA->cmap->N+1,&lindices);CHKERRQ(ierr); 243d0f46423SBarry Smith for (i=0; i<ina->B->cmap->n; i++) { 2442cd6534aSBarry Smith lindices[garray[i]] = i+1; 2452cd6534aSBarry Smith } 246992144d0SBarry Smith no = inA->rmap->mapping->n - nt; 2471795a4d1SJed Brown ierr = PetscCalloc1(inA->rmap->mapping->n+1,&r_rmapo);CHKERRQ(ierr); 2482cd6534aSBarry Smith nt = 0; 249992144d0SBarry Smith for (i=0; i<inA->rmap->mapping->n; i++) { 250992144d0SBarry Smith if (lindices[inA->rmap->mapping->indices[i]]) { 2512cd6534aSBarry Smith nt++; 252992144d0SBarry Smith r_rmapo[i] = lindices[inA->rmap->mapping->indices[i]]; 2532cd6534aSBarry Smith } 2542cd6534aSBarry Smith } 255e32f2f54SBarry Smith if (nt > no) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hmm nt %D no %D",nt,n); 2562cd6534aSBarry Smith ierr = PetscFree(lindices);CHKERRQ(ierr); 257785e854fSJed Brown ierr = PetscMalloc1((nt+1),&auglyrmapo);CHKERRQ(ierr); 258992144d0SBarry Smith for (i=0; i<inA->rmap->mapping->n; i++) { 2592cd6534aSBarry Smith if (r_rmapo[i]) { 2602cd6534aSBarry Smith auglyrmapo[(r_rmapo[i]-1)] = i; 2612cd6534aSBarry Smith } 2622cd6534aSBarry Smith } 2632cd6534aSBarry Smith ierr = PetscFree(r_rmapo);CHKERRQ(ierr); 2642cd6534aSBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,nt,&auglyoo);CHKERRQ(ierr); 2652cd6534aSBarry Smith PetscFunctionReturn(0); 2662cd6534aSBarry Smith } 2672cd6534aSBarry Smith 2682cd6534aSBarry Smith #undef __FUNCT__ 2692cd6534aSBarry Smith #define __FUNCT__ "MatMPIAIJDiagonalScaleLocal" 270dfbe8321SBarry Smith PetscErrorCode MatMPIAIJDiagonalScaleLocal(Mat A,Vec scale) 2712cd6534aSBarry Smith { 27292b32695SKris Buschelman /* This routine should really be abandoned as it duplicates MatDiagonalScaleLocal */ 2734ac538c5SBarry Smith PetscErrorCode ierr; 27492b32695SKris Buschelman 27592b32695SKris Buschelman PetscFunctionBegin; 2764ac538c5SBarry Smith ierr = PetscTryMethod(A,"MatDiagonalScaleLocal_C",(Mat,Vec),(A,scale));CHKERRQ(ierr); 27792b32695SKris Buschelman PetscFunctionReturn(0); 27892b32695SKris Buschelman } 27992b32695SKris Buschelman 28092b32695SKris Buschelman #undef __FUNCT__ 28192b32695SKris Buschelman #define __FUNCT__ "MatDiagonalScaleLocal_MPIAIJ" 2827087cfbeSBarry Smith PetscErrorCode MatDiagonalScaleLocal_MPIAIJ(Mat A,Vec scale) 28392b32695SKris Buschelman { 2842cd6534aSBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ*) A->data; /*access private part of matrix */ 285dfbe8321SBarry Smith PetscErrorCode ierr; 286b1d57f15SBarry Smith PetscInt n,i; 2872cd6534aSBarry Smith PetscScalar *d,*o,*s; 2882cd6534aSBarry Smith 2892cd6534aSBarry Smith PetscFunctionBegin; 2902cd6534aSBarry Smith if (!auglyrmapd) { 2912cd6534aSBarry Smith ierr = MatMPIAIJDiagonalScaleLocalSetUp(A,scale);CHKERRQ(ierr); 2922cd6534aSBarry Smith } 2932cd6534aSBarry Smith 2941ebc52fbSHong Zhang ierr = VecGetArray(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 } 3101ebc52fbSHong Zhang ierr = VecRestoreArray(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