18c79f6d3SBarry Smith /* 28c79f6d3SBarry Smith Support for the parallel AIJ matrix vector multiply 38c79f6d3SBarry Smith */ 470f55243SBarry Smith #include "src/mat/impls/aij/mpi/mpiaij.h" 58c79f6d3SBarry Smith 64a2ae208SSatish Balay #undef __FUNCT__ 74a2ae208SSatish Balay #define __FUNCT__ "MatSetUpMultiply_MPIAIJ" 8dfbe8321SBarry Smith PetscErrorCode MatSetUpMultiply_MPIAIJ(Mat mat) 98c79f6d3SBarry Smith { 1044a69424SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ*)mat->data; 11ec8511deSBarry Smith Mat_SeqAIJ *B = (Mat_SeqAIJ*)(aij->B->data); 12*6849ba73SBarry Smith PetscErrorCode ierr; 13*6849ba73SBarry Smith int i,j,*aj = B->j,ec = 0,*garray; 141eb62cbbSBarry Smith IS from,to; 151eb62cbbSBarry Smith Vec gvec; 16aa482453SBarry Smith #if defined (PETSC_USE_CTABLE) 170f5bd95cSBarry Smith PetscTable gid1_lid1; 180f5bd95cSBarry Smith PetscTablePosition tpos; 192066d6f7SSatish Balay int gid,lid; 206f531f54SSatish Balay #else 216f531f54SSatish Balay int N = mat->N,*indices; 226f531f54SSatish Balay 232066d6f7SSatish Balay #endif 242066d6f7SSatish Balay 253a40ed3dSBarry Smith PetscFunctionBegin; 262066d6f7SSatish Balay 27aa482453SBarry Smith #if defined (PETSC_USE_CTABLE) 282066d6f7SSatish Balay /* use a table - Mark Adams (this has not been tested with "shift") */ 29273d9f13SBarry Smith ierr = PetscTableCreate(aij->B->m,&gid1_lid1);CHKERRQ(ierr); 30273d9f13SBarry Smith for (i=0; i<aij->B->m; i++) { 312066d6f7SSatish Balay for (j=0; j<B->ilen[i]; j++) { 32b3fb0a6cSMatthew Knepley int 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 */ 360f5bd95cSBarry Smith ierr = PetscTableAdd(gid1_lid1,gid1,++ec);CHKERRQ(ierr); 372066d6f7SSatish Balay } 382066d6f7SSatish Balay } 392066d6f7SSatish Balay } 402066d6f7SSatish Balay /* form array of columns we need */ 41b0a32e0cSBarry Smith ierr = PetscMalloc((ec+1)*sizeof(int),&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++) { 520f5bd95cSBarry Smith ierr = PetscTableAdd(gid1_lid1,garray[i]+1,i+1);CHKERRQ(ierr); 532066d6f7SSatish Balay } 542066d6f7SSatish Balay /* compact out the extra columns in B */ 55273d9f13SBarry Smith for (i=0; i<aij->B->m; i++) { 562066d6f7SSatish Balay for (j=0; j<B->ilen[i]; j++) { 57b3fb0a6cSMatthew Knepley int 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 } 63273d9f13SBarry Smith aij->B->n = aij->B->N = ec; 640f5bd95cSBarry Smith ierr = PetscTableDelete(gid1_lid1);CHKERRQ(ierr); 652066d6f7SSatish Balay /* Mark Adams */ 662066d6f7SSatish Balay #else 678c79f6d3SBarry Smith /* For the first stab we make an array as long as the number of columns */ 681eb62cbbSBarry Smith /* mark those columns that are in aij->B */ 69b0a32e0cSBarry Smith ierr = PetscMalloc((N+1)*sizeof(int),&indices);CHKERRQ(ierr); 70549d3d68SSatish Balay ierr = PetscMemzero(indices,N*sizeof(int));CHKERRQ(ierr); 71273d9f13SBarry Smith for (i=0; i<aij->B->m; 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 */ 79b0a32e0cSBarry Smith ierr = PetscMalloc((ec+1)*sizeof(int),&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 */ 91273d9f13SBarry Smith for (i=0; i<aij->B->m; 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 } 96273d9f13SBarry Smith aij->B->n = aij->B->N = ec; 97606d414cSSatish Balay ierr = PetscFree(indices);CHKERRQ(ierr); 982066d6f7SSatish Balay #endif 991eb62cbbSBarry Smith /* create local vector that is used to scatter into */ 100029af93fSBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,ec,&aij->lvec);CHKERRQ(ierr); 1011eb62cbbSBarry Smith 102d6dfbf8fSBarry Smith /* create two temporary Index sets for build scatter gather */ 103b9b97703SBarry Smith ierr = ISCreateGeneral(mat->comm,ec,garray,&from);CHKERRQ(ierr); 104029af93fSBarry Smith ierr = ISCreateStride(PETSC_COMM_SELF,ec,0,1,&to);CHKERRQ(ierr); 1051eb62cbbSBarry Smith 1061eb62cbbSBarry Smith /* create temporary global vector to generate scatter context */ 1071eb62cbbSBarry Smith /* this is inefficient, but otherwise we must do either 1081eb62cbbSBarry Smith 1) save garray until the first actual scatter when the vector is known or 1091eb62cbbSBarry Smith 2) have another way of generating a scatter context without a vector.*/ 110273d9f13SBarry Smith ierr = VecCreateMPI(mat->comm,mat->n,mat->N,&gvec);CHKERRQ(ierr); 1111eb62cbbSBarry Smith 1122d336d48SLois Curfman McInnes /* generate the scatter context */ 11308480c60SBarry Smith ierr = VecScatterCreate(gvec,from,aij->lvec,to,&aij->Mvctx);CHKERRQ(ierr); 114b0a32e0cSBarry Smith PetscLogObjectParent(mat,aij->Mvctx); 115b0a32e0cSBarry Smith PetscLogObjectParent(mat,aij->lvec); 116b0a32e0cSBarry Smith PetscLogObjectParent(mat,from); 117b0a32e0cSBarry Smith PetscLogObjectParent(mat,to); 1189e25ed09SBarry Smith aij->garray = garray; 119b0a32e0cSBarry Smith PetscLogObjectMemory(mat,(ec+1)*sizeof(int)); 12078b31e54SBarry Smith ierr = ISDestroy(from);CHKERRQ(ierr); 12178b31e54SBarry Smith ierr = ISDestroy(to);CHKERRQ(ierr); 122888f2ed8SSatish Balay ierr = VecDestroy(gvec);CHKERRQ(ierr); 1233a40ed3dSBarry Smith PetscFunctionReturn(0); 1248c79f6d3SBarry Smith } 1259e25ed09SBarry Smith 1269e25ed09SBarry Smith 1274a2ae208SSatish Balay #undef __FUNCT__ 1284a2ae208SSatish Balay #define __FUNCT__ "DisAssemble_MPIAIJ" 1292493cbb0SBarry Smith /* 1302493cbb0SBarry Smith Takes the local part of an already assembled MPIAIJ matrix 1312493cbb0SBarry Smith and disassembles it. This is to allow new nonzeros into the matrix 1322493cbb0SBarry Smith that require more communication in the matrix vector multiply. 1332493cbb0SBarry Smith Thus certain data-structures must be rebuilt. 1342493cbb0SBarry Smith 1352493cbb0SBarry Smith Kind of slow! But that's what application programmers get when 1362493cbb0SBarry Smith they are sloppy. 1372493cbb0SBarry Smith */ 138dfbe8321SBarry Smith PetscErrorCode DisAssemble_MPIAIJ(Mat A) 1392493cbb0SBarry Smith { 1402493cbb0SBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ*)A->data; 1412493cbb0SBarry Smith Mat B = aij->B,Bnew; 142ec8511deSBarry Smith Mat_SeqAIJ *Baij = (Mat_SeqAIJ*)B->data; 143dfbe8321SBarry Smith PetscErrorCode ierr; 144dfbe8321SBarry Smith int i,j,m = B->m,n = A->N,col,ct = 0,*garray = aij->garray; 145bfec09a0SHong Zhang int *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; 160b0a32e0cSBarry Smith PetscLogObjectMemory(A,-aij->B->n*sizeof(int)); 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 */ 169b0a32e0cSBarry Smith ierr = PetscMalloc((m+1)*sizeof(int),&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; 186b0a32e0cSBarry Smith PetscLogObjectMemory(A,-ec*sizeof(int)); 1872493cbb0SBarry Smith ierr = MatDestroy(B);CHKERRQ(ierr); 188b0a32e0cSBarry Smith PetscLogObjectParent(A,Bnew); 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 1962cd6534aSBarry Smith static int *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; 207dfbe8321SBarry Smith int i,n,nt,cstart,cend,no,*garray = ina->garray,*lindices; 2082cd6534aSBarry Smith int *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); 2132cd6534aSBarry Smith ierr = PetscMalloc((inA->mapping->n+1)*sizeof(int),&r_rmapd);CHKERRQ(ierr); 2142cd6534aSBarry Smith ierr = PetscMemzero(r_rmapd,inA->mapping->n*sizeof(int));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 } 2222cd6534aSBarry Smith if (nt != n) SETERRQ2(1,"Hmm nt %d n %d",nt,n); 2232cd6534aSBarry Smith ierr = PetscMalloc((n+1)*sizeof(int),&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 2322cd6534aSBarry Smith ierr = PetscMalloc((inA->N+1)*sizeof(int),&lindices);CHKERRQ(ierr); 2332cd6534aSBarry Smith ierr = PetscMemzero(lindices,inA->N*sizeof(int));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; 2382cd6534aSBarry Smith ierr = PetscMalloc((inA->mapping->n+1)*sizeof(int),&r_rmapo);CHKERRQ(ierr); 2392cd6534aSBarry Smith ierr = PetscMemzero(r_rmapo,inA->mapping->n*sizeof(int));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 } 2472cd6534aSBarry Smith if (nt > no) SETERRQ2(1,"Hmm nt %d no %d",nt,n); 2482cd6534aSBarry Smith ierr = PetscFree(lindices);CHKERRQ(ierr); 2492cd6534aSBarry Smith ierr = PetscMalloc((nt+1)*sizeof(int),&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" 279dfbe8321SBarry 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; 283dfbe8321SBarry Smith int 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