1*be1d678aSKris Buschelman #define PETSCMAT_DLL 2*be1d678aSKris Buschelman 38c79f6d3SBarry Smith /* 48c79f6d3SBarry Smith Support for the parallel AIJ matrix vector multiply 58c79f6d3SBarry Smith */ 670f55243SBarry Smith #include "src/mat/impls/aij/mpi/mpiaij.h" 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 23b1d57f15SBarry Smith PetscInt N = mat->N,*indices; 242066d6f7SSatish Balay #endif 252066d6f7SSatish Balay 263a40ed3dSBarry Smith PetscFunctionBegin; 272066d6f7SSatish Balay 28aa482453SBarry Smith #if defined (PETSC_USE_CTABLE) 29dc231df0SBarry Smith /* use a table - Mark Adams */ 30273d9f13SBarry Smith ierr = PetscTableCreate(aij->B->m,&gid1_lid1);CHKERRQ(ierr); 31273d9f13SBarry Smith for (i=0; i<aij->B->m; 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 */ 370f5bd95cSBarry Smith ierr = PetscTableAdd(gid1_lid1,gid1,++ec);CHKERRQ(ierr); 382066d6f7SSatish Balay } 392066d6f7SSatish Balay } 402066d6f7SSatish Balay } 412066d6f7SSatish Balay /* form array of columns we need */ 42b1d57f15SBarry Smith ierr = PetscMalloc((ec+1)*sizeof(PetscInt),&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++) { 530f5bd95cSBarry Smith ierr = PetscTableAdd(gid1_lid1,garray[i]+1,i+1);CHKERRQ(ierr); 542066d6f7SSatish Balay } 552066d6f7SSatish Balay /* compact out the extra columns in B */ 56273d9f13SBarry Smith for (i=0; i<aij->B->m; 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 } 64273d9f13SBarry Smith aij->B->n = aij->B->N = ec; 650f5bd95cSBarry Smith ierr = PetscTableDelete(gid1_lid1);CHKERRQ(ierr); 662066d6f7SSatish Balay /* Mark Adams */ 672066d6f7SSatish Balay #else 688c79f6d3SBarry Smith /* For the first stab we make an array as long as the number of columns */ 691eb62cbbSBarry Smith /* mark those columns that are in aij->B */ 70b1d57f15SBarry Smith ierr = PetscMalloc((N+1)*sizeof(PetscInt),&indices);CHKERRQ(ierr); 71b1d57f15SBarry Smith ierr = PetscMemzero(indices,N*sizeof(PetscInt));CHKERRQ(ierr); 72273d9f13SBarry Smith for (i=0; i<aij->B->m; 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 */ 80b1d57f15SBarry Smith ierr = PetscMalloc((ec+1)*sizeof(PetscInt),&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 */ 92273d9f13SBarry Smith for (i=0; i<aij->B->m; 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 } 97273d9f13SBarry Smith aij->B->n = aij->B->N = ec; 98606d414cSSatish Balay ierr = PetscFree(indices);CHKERRQ(ierr); 992066d6f7SSatish Balay #endif 1001eb62cbbSBarry Smith /* create local vector that is used to scatter into */ 101029af93fSBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,ec,&aij->lvec);CHKERRQ(ierr); 1021eb62cbbSBarry Smith 103d6dfbf8fSBarry Smith /* create two temporary Index sets for build scatter gather */ 104b9b97703SBarry Smith ierr = ISCreateGeneral(mat->comm,ec,garray,&from);CHKERRQ(ierr); 105029af93fSBarry Smith ierr = ISCreateStride(PETSC_COMM_SELF,ec,0,1,&to);CHKERRQ(ierr); 1061eb62cbbSBarry Smith 1071eb62cbbSBarry Smith /* create temporary global vector to generate scatter context */ 1081eb62cbbSBarry Smith /* this is inefficient, but otherwise we must do either 1091eb62cbbSBarry Smith 1) save garray until the first actual scatter when the vector is known or 1101eb62cbbSBarry Smith 2) have another way of generating a scatter context without a vector.*/ 111273d9f13SBarry Smith ierr = VecCreateMPI(mat->comm,mat->n,mat->N,&gvec);CHKERRQ(ierr); 1121eb62cbbSBarry Smith 1132d336d48SLois Curfman McInnes /* generate the scatter context */ 11408480c60SBarry Smith ierr = VecScatterCreate(gvec,from,aij->lvec,to,&aij->Mvctx);CHKERRQ(ierr); 11552e6d16bSBarry Smith ierr = PetscLogObjectParent(mat,aij->Mvctx);CHKERRQ(ierr); 11652e6d16bSBarry Smith ierr = PetscLogObjectParent(mat,aij->lvec);CHKERRQ(ierr); 11752e6d16bSBarry Smith ierr = PetscLogObjectParent(mat,from);CHKERRQ(ierr); 11852e6d16bSBarry Smith ierr = PetscLogObjectParent(mat,to);CHKERRQ(ierr); 1199e25ed09SBarry Smith aij->garray = garray; 12052e6d16bSBarry Smith ierr = PetscLogObjectMemory(mat,(ec+1)*sizeof(PetscInt));CHKERRQ(ierr); 12178b31e54SBarry Smith ierr = ISDestroy(from);CHKERRQ(ierr); 12278b31e54SBarry Smith ierr = ISDestroy(to);CHKERRQ(ierr); 123888f2ed8SSatish Balay ierr = VecDestroy(gvec);CHKERRQ(ierr); 1243a40ed3dSBarry Smith PetscFunctionReturn(0); 1258c79f6d3SBarry Smith } 1269e25ed09SBarry Smith 1279e25ed09SBarry Smith 1284a2ae208SSatish Balay #undef __FUNCT__ 1294a2ae208SSatish Balay #define __FUNCT__ "DisAssemble_MPIAIJ" 1302493cbb0SBarry Smith /* 1312493cbb0SBarry Smith Takes the local part of an already assembled MPIAIJ matrix 1322493cbb0SBarry Smith and disassembles it. This is to allow new nonzeros into the matrix 1332493cbb0SBarry Smith that require more communication in the matrix vector multiply. 1342493cbb0SBarry Smith Thus certain data-structures must be rebuilt. 1352493cbb0SBarry Smith 1362493cbb0SBarry Smith Kind of slow! But that's what application programmers get when 1372493cbb0SBarry Smith they are sloppy. 1382493cbb0SBarry Smith */ 139dfbe8321SBarry Smith PetscErrorCode DisAssemble_MPIAIJ(Mat A) 1402493cbb0SBarry Smith { 1412493cbb0SBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ*)A->data; 1422493cbb0SBarry Smith Mat B = aij->B,Bnew; 143ec8511deSBarry Smith Mat_SeqAIJ *Baij = (Mat_SeqAIJ*)B->data; 144dfbe8321SBarry Smith PetscErrorCode ierr; 14557809a77SBarry Smith PetscInt i,j,m = B->m,n = A->N,col,ct = 0,*garray = aij->garray,*nz,ec; 14687828ca2SBarry Smith PetscScalar v; 1472493cbb0SBarry Smith 1483a40ed3dSBarry Smith PetscFunctionBegin; 1492493cbb0SBarry Smith /* free stuff related to matrix-vec multiply */ 150b0a32e0cSBarry Smith ierr = VecGetSize(aij->lvec,&ec);CHKERRQ(ierr); /* needed for PetscLogObjectMemory below */ 1512493cbb0SBarry Smith ierr = VecDestroy(aij->lvec);CHKERRQ(ierr); aij->lvec = 0; 15208480c60SBarry Smith ierr = VecScatterDestroy(aij->Mvctx);CHKERRQ(ierr); aij->Mvctx = 0; 153464493b3SBarry Smith if (aij->colmap) { 154aa482453SBarry Smith #if defined (PETSC_USE_CTABLE) 1550f5bd95cSBarry Smith ierr = PetscTableDelete(aij->colmap);CHKERRQ(ierr); 1560f5bd95cSBarry Smith aij->colmap = 0; 1572066d6f7SSatish Balay #else 158606d414cSSatish Balay ierr = PetscFree(aij->colmap);CHKERRQ(ierr); 159606d414cSSatish Balay aij->colmap = 0; 16052e6d16bSBarry Smith ierr = PetscLogObjectMemory(A,-aij->B->n*sizeof(PetscInt));CHKERRQ(ierr); 1612066d6f7SSatish Balay #endif 162464493b3SBarry Smith } 1632493cbb0SBarry Smith 1642493cbb0SBarry Smith /* make sure that B is assembled so we can access its values */ 1656d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 166fe2f2677SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1672493cbb0SBarry Smith 1682493cbb0SBarry Smith /* invent new B and copy stuff over */ 169b1d57f15SBarry Smith ierr = PetscMalloc((m+1)*sizeof(PetscInt),&nz);CHKERRQ(ierr); 17048b35521SBarry Smith for (i=0; i<m; i++) { 17148b35521SBarry Smith nz[i] = Baij->i[i+1] - Baij->i[i]; 17248b35521SBarry Smith } 173f204ca49SKris Buschelman ierr = MatCreate(PETSC_COMM_SELF,m,n,m,n,&Bnew);CHKERRQ(ierr); 174f204ca49SKris Buschelman ierr = MatSetType(Bnew,B->type_name);CHKERRQ(ierr); 175f204ca49SKris Buschelman ierr = MatSeqAIJSetPreallocation(Bnew,0,nz);CHKERRQ(ierr); 176606d414cSSatish Balay ierr = PetscFree(nz);CHKERRQ(ierr); 1772493cbb0SBarry Smith for (i=0; i<m; i++) { 178bfec09a0SHong Zhang for (j=Baij->i[i]; j<Baij->i[i+1]; j++) { 179bfec09a0SHong Zhang col = garray[Baij->j[ct]]; 1802493cbb0SBarry Smith v = Baij->a[ct++]; 18183271157SBarry Smith ierr = MatSetValues(Bnew,1,&i,1,&col,&v,B->insertmode);CHKERRQ(ierr); 1822493cbb0SBarry Smith } 1832493cbb0SBarry Smith } 184606d414cSSatish Balay ierr = PetscFree(aij->garray);CHKERRQ(ierr); 185606d414cSSatish Balay aij->garray = 0; 18652e6d16bSBarry Smith ierr = PetscLogObjectMemory(A,-ec*sizeof(PetscInt));CHKERRQ(ierr); 1872493cbb0SBarry Smith ierr = MatDestroy(B);CHKERRQ(ierr); 18852e6d16bSBarry Smith ierr = PetscLogObjectParent(A,Bnew);CHKERRQ(ierr); 1892493cbb0SBarry Smith aij->B = Bnew; 190227d817aSBarry Smith A->was_assembled = PETSC_FALSE; 1913a40ed3dSBarry Smith PetscFunctionReturn(0); 1922493cbb0SBarry Smith } 1932493cbb0SBarry Smith 1942cd6534aSBarry Smith /* ugly stuff added for Glenn someday we should fix this up */ 1952cd6534aSBarry Smith 196b1d57f15SBarry Smith static PetscInt *auglyrmapd = 0,*auglyrmapo = 0; /* mapping from the local ordering to the "diagonal" and "off-diagonal" 1972cd6534aSBarry Smith parts of the local matrix */ 1982cd6534aSBarry Smith static Vec auglydd = 0,auglyoo = 0; /* work vectors used to scale the two parts of the local matrix */ 1992cd6534aSBarry Smith 2002cd6534aSBarry Smith 2012cd6534aSBarry Smith #undef __FUNCT__ 2022cd6534aSBarry Smith #define __FUNCT__ "MatMPIAIJDiagonalScaleLocalSetUp" 203dfbe8321SBarry Smith PetscErrorCode MatMPIAIJDiagonalScaleLocalSetUp(Mat inA,Vec scale) 2042cd6534aSBarry Smith { 2052cd6534aSBarry Smith Mat_MPIAIJ *ina = (Mat_MPIAIJ*) inA->data; /*access private part of matrix */ 206dfbe8321SBarry Smith PetscErrorCode ierr; 207b1d57f15SBarry Smith PetscInt i,n,nt,cstart,cend,no,*garray = ina->garray,*lindices; 208b1d57f15SBarry Smith PetscInt *r_rmapd,*r_rmapo; 2092cd6534aSBarry Smith 2102cd6534aSBarry Smith PetscFunctionBegin; 2112cd6534aSBarry Smith ierr = MatGetOwnershipRange(inA,&cstart,&cend);CHKERRQ(ierr); 2122cd6534aSBarry Smith ierr = MatGetSize(ina->A,PETSC_NULL,&n);CHKERRQ(ierr); 213b1d57f15SBarry Smith ierr = PetscMalloc((inA->mapping->n+1)*sizeof(PetscInt),&r_rmapd);CHKERRQ(ierr); 214b1d57f15SBarry Smith ierr = PetscMemzero(r_rmapd,inA->mapping->n*sizeof(PetscInt));CHKERRQ(ierr); 2152cd6534aSBarry Smith nt = 0; 2162cd6534aSBarry Smith for (i=0; i<inA->mapping->n; i++) { 2172cd6534aSBarry Smith if (inA->mapping->indices[i] >= cstart && inA->mapping->indices[i] < cend) { 2182cd6534aSBarry Smith nt++; 2192cd6534aSBarry Smith r_rmapd[i] = inA->mapping->indices[i] + 1; 2202cd6534aSBarry Smith } 2212cd6534aSBarry Smith } 22277431f27SBarry Smith if (nt != n) SETERRQ2(PETSC_ERR_PLIB,"Hmm nt %D n %D",nt,n); 223b1d57f15SBarry Smith ierr = PetscMalloc((n+1)*sizeof(PetscInt),&auglyrmapd);CHKERRQ(ierr); 2242cd6534aSBarry Smith for (i=0; i<inA->mapping->n; i++) { 2252cd6534aSBarry Smith if (r_rmapd[i]){ 2262cd6534aSBarry Smith auglyrmapd[(r_rmapd[i]-1)-cstart] = i; 2272cd6534aSBarry Smith } 2282cd6534aSBarry Smith } 2292cd6534aSBarry Smith ierr = PetscFree(r_rmapd);CHKERRQ(ierr); 2302cd6534aSBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,n,&auglydd);CHKERRQ(ierr); 2312cd6534aSBarry Smith 232b1d57f15SBarry Smith ierr = PetscMalloc((inA->N+1)*sizeof(PetscInt),&lindices);CHKERRQ(ierr); 233b1d57f15SBarry Smith ierr = PetscMemzero(lindices,inA->N*sizeof(PetscInt));CHKERRQ(ierr); 2342cd6534aSBarry Smith for (i=0; i<ina->B->n; i++) { 2352cd6534aSBarry Smith lindices[garray[i]] = i+1; 2362cd6534aSBarry Smith } 2372cd6534aSBarry Smith no = inA->mapping->n - nt; 238b1d57f15SBarry Smith ierr = PetscMalloc((inA->mapping->n+1)*sizeof(PetscInt),&r_rmapo);CHKERRQ(ierr); 239b1d57f15SBarry Smith ierr = PetscMemzero(r_rmapo,inA->mapping->n*sizeof(PetscInt));CHKERRQ(ierr); 2402cd6534aSBarry Smith nt = 0; 2412cd6534aSBarry Smith for (i=0; i<inA->mapping->n; i++) { 2422cd6534aSBarry Smith if (lindices[inA->mapping->indices[i]]) { 2432cd6534aSBarry Smith nt++; 2442cd6534aSBarry Smith r_rmapo[i] = lindices[inA->mapping->indices[i]]; 2452cd6534aSBarry Smith } 2462cd6534aSBarry Smith } 24777431f27SBarry Smith if (nt > no) SETERRQ2(PETSC_ERR_PLIB,"Hmm nt %D no %D",nt,n); 2482cd6534aSBarry Smith ierr = PetscFree(lindices);CHKERRQ(ierr); 249b1d57f15SBarry Smith ierr = PetscMalloc((nt+1)*sizeof(PetscInt),&auglyrmapo);CHKERRQ(ierr); 2502cd6534aSBarry Smith for (i=0; i<inA->mapping->n; i++) { 2512cd6534aSBarry Smith if (r_rmapo[i]){ 2522cd6534aSBarry Smith auglyrmapo[(r_rmapo[i]-1)] = i; 2532cd6534aSBarry Smith } 2542cd6534aSBarry Smith } 2552cd6534aSBarry Smith ierr = PetscFree(r_rmapo);CHKERRQ(ierr); 2562cd6534aSBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,nt,&auglyoo);CHKERRQ(ierr); 2572cd6534aSBarry Smith 2582cd6534aSBarry Smith PetscFunctionReturn(0); 2592cd6534aSBarry Smith } 2602cd6534aSBarry Smith 2612cd6534aSBarry Smith #undef __FUNCT__ 2622cd6534aSBarry Smith #define __FUNCT__ "MatMPIAIJDiagonalScaleLocal" 263dfbe8321SBarry Smith PetscErrorCode MatMPIAIJDiagonalScaleLocal(Mat A,Vec scale) 2642cd6534aSBarry Smith { 26592b32695SKris Buschelman /* This routine should really be abandoned as it duplicates MatDiagonalScaleLocal */ 266dfbe8321SBarry Smith PetscErrorCode ierr,(*f)(Mat,Vec); 26792b32695SKris Buschelman 26892b32695SKris Buschelman PetscFunctionBegin; 26992b32695SKris Buschelman ierr = PetscObjectQueryFunction((PetscObject)A,"MatDiagonalScaleLocal_C",(void (**)(void))&f);CHKERRQ(ierr); 27092b32695SKris Buschelman if (f) { 27192b32695SKris Buschelman ierr = (*f)(A,scale);CHKERRQ(ierr); 27292b32695SKris Buschelman } 27392b32695SKris Buschelman PetscFunctionReturn(0); 27492b32695SKris Buschelman } 27592b32695SKris Buschelman 27692b32695SKris Buschelman EXTERN_C_BEGIN 27792b32695SKris Buschelman #undef __FUNCT__ 27892b32695SKris Buschelman #define __FUNCT__ "MatDiagonalScaleLocal_MPIAIJ" 279*be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT 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; 2842cd6534aSBarry Smith PetscScalar *d,*o,*s; 2852cd6534aSBarry Smith 2862cd6534aSBarry Smith PetscFunctionBegin; 2872cd6534aSBarry Smith if (!auglyrmapd) { 2882cd6534aSBarry Smith ierr = MatMPIAIJDiagonalScaleLocalSetUp(A,scale);CHKERRQ(ierr); 2892cd6534aSBarry Smith } 2902cd6534aSBarry Smith 2911ebc52fbSHong Zhang ierr = VecGetArray(scale,&s);CHKERRQ(ierr); 2922cd6534aSBarry Smith 2932cd6534aSBarry Smith ierr = VecGetLocalSize(auglydd,&n);CHKERRQ(ierr); 2941ebc52fbSHong Zhang ierr = VecGetArray(auglydd,&d);CHKERRQ(ierr); 2952cd6534aSBarry Smith for (i=0; i<n; i++) { 2962cd6534aSBarry Smith d[i] = s[auglyrmapd[i]]; /* copy "diagonal" (true local) portion of scale into dd vector */ 2972cd6534aSBarry Smith } 2981ebc52fbSHong Zhang ierr = VecRestoreArray(auglydd,&d);CHKERRQ(ierr); 2992cd6534aSBarry Smith /* column scale "diagonal" portion of local matrix */ 3002cd6534aSBarry Smith ierr = MatDiagonalScale(a->A,PETSC_NULL,auglydd);CHKERRQ(ierr); 3012cd6534aSBarry Smith 3022cd6534aSBarry Smith ierr = VecGetLocalSize(auglyoo,&n);CHKERRQ(ierr); 3031ebc52fbSHong Zhang ierr = VecGetArray(auglyoo,&o);CHKERRQ(ierr); 3042cd6534aSBarry Smith for (i=0; i<n; i++) { 3052cd6534aSBarry Smith o[i] = s[auglyrmapo[i]]; /* copy "off-diagonal" portion of scale into oo vector */ 3062cd6534aSBarry Smith } 3071ebc52fbSHong Zhang ierr = VecRestoreArray(scale,&s);CHKERRQ(ierr); 3081ebc52fbSHong Zhang ierr = VecRestoreArray(auglyoo,&o);CHKERRQ(ierr); 3092cd6534aSBarry Smith /* column scale "off-diagonal" portion of local matrix */ 3102cd6534aSBarry Smith ierr = MatDiagonalScale(a->B,PETSC_NULL,auglyoo);CHKERRQ(ierr); 3112cd6534aSBarry Smith 3122cd6534aSBarry Smith PetscFunctionReturn(0); 3132cd6534aSBarry Smith } 31492b32695SKris Buschelman EXTERN_C_END 3152cd6534aSBarry Smith 31648b35521SBarry Smith 317