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); 126849ba73SBarry Smith PetscErrorCode ierr; 13*b1d57f15SBarry Smith PetscInt 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; 19*b1d57f15SBarry Smith PetscInt gid,lid; 206f531f54SSatish Balay #else 21*b1d57f15SBarry Smith PetscInt N = mat->N,*indices; 222066d6f7SSatish Balay #endif 232066d6f7SSatish Balay 243a40ed3dSBarry Smith PetscFunctionBegin; 252066d6f7SSatish Balay 26aa482453SBarry Smith #if defined (PETSC_USE_CTABLE) 272066d6f7SSatish Balay /* use a table - Mark Adams (this has not been tested with "shift") */ 28273d9f13SBarry Smith ierr = PetscTableCreate(aij->B->m,&gid1_lid1);CHKERRQ(ierr); 29273d9f13SBarry Smith for (i=0; i<aij->B->m; i++) { 302066d6f7SSatish Balay for (j=0; j<B->ilen[i]; j++) { 31*b1d57f15SBarry Smith PetscInt data,gid1 = aj[B->i[i] + j] + 1; 320f5bd95cSBarry Smith ierr = PetscTableFind(gid1_lid1,gid1,&data);CHKERRQ(ierr); 33fa46199cSSatish Balay if (!data) { 342066d6f7SSatish Balay /* one based table */ 350f5bd95cSBarry Smith ierr = PetscTableAdd(gid1_lid1,gid1,++ec);CHKERRQ(ierr); 362066d6f7SSatish Balay } 372066d6f7SSatish Balay } 382066d6f7SSatish Balay } 392066d6f7SSatish Balay /* form array of columns we need */ 40*b1d57f15SBarry Smith ierr = PetscMalloc((ec+1)*sizeof(PetscInt),&garray);CHKERRQ(ierr); 410f5bd95cSBarry Smith ierr = PetscTableGetHeadPosition(gid1_lid1,&tpos);CHKERRQ(ierr); 422066d6f7SSatish Balay while (tpos) { 430f5bd95cSBarry Smith ierr = PetscTableGetNext(gid1_lid1,&tpos,&gid,&lid);CHKERRQ(ierr); 44b0a32e0cSBarry Smith gid--; 45b0a32e0cSBarry Smith lid--; 462066d6f7SSatish Balay garray[lid] = gid; 472066d6f7SSatish Balay } 480064e2bbSSatish Balay ierr = PetscSortInt(ec,garray);CHKERRQ(ierr); /* sort, and rebuild */ 490f5bd95cSBarry Smith ierr = PetscTableRemoveAll(gid1_lid1);CHKERRQ(ierr); 502066d6f7SSatish Balay for (i=0; i<ec; i++) { 510f5bd95cSBarry Smith ierr = PetscTableAdd(gid1_lid1,garray[i]+1,i+1);CHKERRQ(ierr); 522066d6f7SSatish Balay } 532066d6f7SSatish Balay /* compact out the extra columns in B */ 54273d9f13SBarry Smith for (i=0; i<aij->B->m; i++) { 552066d6f7SSatish Balay for (j=0; j<B->ilen[i]; j++) { 56*b1d57f15SBarry Smith PetscInt gid1 = aj[B->i[i] + j] + 1; 570f5bd95cSBarry Smith ierr = PetscTableFind(gid1_lid1,gid1,&lid);CHKERRQ(ierr); 58fa46199cSSatish Balay lid --; 59b3fb0a6cSMatthew Knepley aj[B->i[i] + j] = lid; 602066d6f7SSatish Balay } 612066d6f7SSatish Balay } 62273d9f13SBarry Smith aij->B->n = aij->B->N = ec; 630f5bd95cSBarry Smith ierr = PetscTableDelete(gid1_lid1);CHKERRQ(ierr); 642066d6f7SSatish Balay /* Mark Adams */ 652066d6f7SSatish Balay #else 668c79f6d3SBarry Smith /* For the first stab we make an array as long as the number of columns */ 671eb62cbbSBarry Smith /* mark those columns that are in aij->B */ 68*b1d57f15SBarry Smith ierr = PetscMalloc((N+1)*sizeof(PetscInt),&indices);CHKERRQ(ierr); 69*b1d57f15SBarry Smith ierr = PetscMemzero(indices,N*sizeof(PetscInt));CHKERRQ(ierr); 70273d9f13SBarry Smith for (i=0; i<aij->B->m; i++) { 71d6dfbf8fSBarry Smith for (j=0; j<B->ilen[i]; j++) { 72bfec09a0SHong Zhang if (!indices[aj[B->i[i] + j] ]) ec++; 73bfec09a0SHong Zhang indices[aj[B->i[i] + j] ] = 1; 74416022c9SBarry Smith } 751eb62cbbSBarry Smith } 768c79f6d3SBarry Smith 771eb62cbbSBarry Smith /* form array of columns we need */ 78*b1d57f15SBarry Smith ierr = PetscMalloc((ec+1)*sizeof(PetscInt),&garray);CHKERRQ(ierr); 791eb62cbbSBarry Smith ec = 0; 801eb62cbbSBarry Smith for (i=0; i<N; i++) { 811eb62cbbSBarry Smith if (indices[i]) garray[ec++] = i; 821eb62cbbSBarry Smith } 831eb62cbbSBarry Smith 841eb62cbbSBarry Smith /* make indices now point into garray */ 851eb62cbbSBarry Smith for (i=0; i<ec; i++) { 86bfec09a0SHong Zhang indices[garray[i]] = i; 871eb62cbbSBarry Smith } 881eb62cbbSBarry Smith 891eb62cbbSBarry Smith /* compact out the extra columns in B */ 90273d9f13SBarry Smith for (i=0; i<aij->B->m; i++) { 91d6dfbf8fSBarry Smith for (j=0; j<B->ilen[i]; j++) { 92bfec09a0SHong Zhang aj[B->i[i] + j] = indices[aj[B->i[i] + j]]; 931eb62cbbSBarry Smith } 94d6dfbf8fSBarry Smith } 95273d9f13SBarry Smith aij->B->n = aij->B->N = ec; 96606d414cSSatish Balay ierr = PetscFree(indices);CHKERRQ(ierr); 972066d6f7SSatish Balay #endif 981eb62cbbSBarry Smith /* create local vector that is used to scatter into */ 99029af93fSBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,ec,&aij->lvec);CHKERRQ(ierr); 1001eb62cbbSBarry Smith 101d6dfbf8fSBarry Smith /* create two temporary Index sets for build scatter gather */ 102b9b97703SBarry Smith ierr = ISCreateGeneral(mat->comm,ec,garray,&from);CHKERRQ(ierr); 103029af93fSBarry Smith ierr = ISCreateStride(PETSC_COMM_SELF,ec,0,1,&to);CHKERRQ(ierr); 1041eb62cbbSBarry Smith 1051eb62cbbSBarry Smith /* create temporary global vector to generate scatter context */ 1061eb62cbbSBarry Smith /* this is inefficient, but otherwise we must do either 1071eb62cbbSBarry Smith 1) save garray until the first actual scatter when the vector is known or 1081eb62cbbSBarry Smith 2) have another way of generating a scatter context without a vector.*/ 109273d9f13SBarry Smith ierr = VecCreateMPI(mat->comm,mat->n,mat->N,&gvec);CHKERRQ(ierr); 1101eb62cbbSBarry Smith 1112d336d48SLois Curfman McInnes /* generate the scatter context */ 11208480c60SBarry Smith ierr = VecScatterCreate(gvec,from,aij->lvec,to,&aij->Mvctx);CHKERRQ(ierr); 113b0a32e0cSBarry Smith PetscLogObjectParent(mat,aij->Mvctx); 114b0a32e0cSBarry Smith PetscLogObjectParent(mat,aij->lvec); 115b0a32e0cSBarry Smith PetscLogObjectParent(mat,from); 116b0a32e0cSBarry Smith PetscLogObjectParent(mat,to); 1179e25ed09SBarry Smith aij->garray = garray; 118*b1d57f15SBarry Smith PetscLogObjectMemory(mat,(ec+1)*sizeof(PetscInt)); 11978b31e54SBarry Smith ierr = ISDestroy(from);CHKERRQ(ierr); 12078b31e54SBarry Smith ierr = ISDestroy(to);CHKERRQ(ierr); 121888f2ed8SSatish Balay ierr = VecDestroy(gvec);CHKERRQ(ierr); 1223a40ed3dSBarry Smith PetscFunctionReturn(0); 1238c79f6d3SBarry Smith } 1249e25ed09SBarry Smith 1259e25ed09SBarry Smith 1264a2ae208SSatish Balay #undef __FUNCT__ 1274a2ae208SSatish Balay #define __FUNCT__ "DisAssemble_MPIAIJ" 1282493cbb0SBarry Smith /* 1292493cbb0SBarry Smith Takes the local part of an already assembled MPIAIJ matrix 1302493cbb0SBarry Smith and disassembles it. This is to allow new nonzeros into the matrix 1312493cbb0SBarry Smith that require more communication in the matrix vector multiply. 1322493cbb0SBarry Smith Thus certain data-structures must be rebuilt. 1332493cbb0SBarry Smith 1342493cbb0SBarry Smith Kind of slow! But that's what application programmers get when 1352493cbb0SBarry Smith they are sloppy. 1362493cbb0SBarry Smith */ 137dfbe8321SBarry Smith PetscErrorCode DisAssemble_MPIAIJ(Mat A) 1382493cbb0SBarry Smith { 1392493cbb0SBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ*)A->data; 1402493cbb0SBarry Smith Mat B = aij->B,Bnew; 141ec8511deSBarry Smith Mat_SeqAIJ *Baij = (Mat_SeqAIJ*)B->data; 142dfbe8321SBarry Smith PetscErrorCode ierr; 14357809a77SBarry Smith PetscInt i,j,m = B->m,n = A->N,col,ct = 0,*garray = aij->garray,*nz,ec; 14487828ca2SBarry Smith PetscScalar v; 1452493cbb0SBarry Smith 1463a40ed3dSBarry Smith PetscFunctionBegin; 1472493cbb0SBarry Smith /* free stuff related to matrix-vec multiply */ 148b0a32e0cSBarry Smith ierr = VecGetSize(aij->lvec,&ec);CHKERRQ(ierr); /* needed for PetscLogObjectMemory below */ 1492493cbb0SBarry Smith ierr = VecDestroy(aij->lvec);CHKERRQ(ierr); aij->lvec = 0; 15008480c60SBarry Smith ierr = VecScatterDestroy(aij->Mvctx);CHKERRQ(ierr); aij->Mvctx = 0; 151464493b3SBarry Smith if (aij->colmap) { 152aa482453SBarry Smith #if defined (PETSC_USE_CTABLE) 1530f5bd95cSBarry Smith ierr = PetscTableDelete(aij->colmap);CHKERRQ(ierr); 1540f5bd95cSBarry Smith aij->colmap = 0; 1552066d6f7SSatish Balay #else 156606d414cSSatish Balay ierr = PetscFree(aij->colmap);CHKERRQ(ierr); 157606d414cSSatish Balay aij->colmap = 0; 158*b1d57f15SBarry Smith PetscLogObjectMemory(A,-aij->B->n*sizeof(PetscInt)); 1592066d6f7SSatish Balay #endif 160464493b3SBarry Smith } 1612493cbb0SBarry Smith 1622493cbb0SBarry Smith /* make sure that B is assembled so we can access its values */ 1636d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 164fe2f2677SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1652493cbb0SBarry Smith 1662493cbb0SBarry Smith /* invent new B and copy stuff over */ 167*b1d57f15SBarry Smith ierr = PetscMalloc((m+1)*sizeof(PetscInt),&nz);CHKERRQ(ierr); 16848b35521SBarry Smith for (i=0; i<m; i++) { 16948b35521SBarry Smith nz[i] = Baij->i[i+1] - Baij->i[i]; 17048b35521SBarry Smith } 171f204ca49SKris Buschelman ierr = MatCreate(PETSC_COMM_SELF,m,n,m,n,&Bnew);CHKERRQ(ierr); 172f204ca49SKris Buschelman ierr = MatSetType(Bnew,B->type_name);CHKERRQ(ierr); 173f204ca49SKris Buschelman ierr = MatSeqAIJSetPreallocation(Bnew,0,nz);CHKERRQ(ierr); 174606d414cSSatish Balay ierr = PetscFree(nz);CHKERRQ(ierr); 1752493cbb0SBarry Smith for (i=0; i<m; i++) { 176bfec09a0SHong Zhang for (j=Baij->i[i]; j<Baij->i[i+1]; j++) { 177bfec09a0SHong Zhang col = garray[Baij->j[ct]]; 1782493cbb0SBarry Smith v = Baij->a[ct++]; 17983271157SBarry Smith ierr = MatSetValues(Bnew,1,&i,1,&col,&v,B->insertmode);CHKERRQ(ierr); 1802493cbb0SBarry Smith } 1812493cbb0SBarry Smith } 182606d414cSSatish Balay ierr = PetscFree(aij->garray);CHKERRQ(ierr); 183606d414cSSatish Balay aij->garray = 0; 184*b1d57f15SBarry Smith PetscLogObjectMemory(A,-ec*sizeof(PetscInt)); 1852493cbb0SBarry Smith ierr = MatDestroy(B);CHKERRQ(ierr); 186b0a32e0cSBarry Smith PetscLogObjectParent(A,Bnew); 1872493cbb0SBarry Smith aij->B = Bnew; 188227d817aSBarry Smith A->was_assembled = PETSC_FALSE; 1893a40ed3dSBarry Smith PetscFunctionReturn(0); 1902493cbb0SBarry Smith } 1912493cbb0SBarry Smith 1922cd6534aSBarry Smith /* ugly stuff added for Glenn someday we should fix this up */ 1932cd6534aSBarry Smith 194*b1d57f15SBarry Smith static PetscInt *auglyrmapd = 0,*auglyrmapo = 0; /* mapping from the local ordering to the "diagonal" and "off-diagonal" 1952cd6534aSBarry Smith parts of the local matrix */ 1962cd6534aSBarry Smith static Vec auglydd = 0,auglyoo = 0; /* work vectors used to scale the two parts of the local matrix */ 1972cd6534aSBarry Smith 1982cd6534aSBarry Smith 1992cd6534aSBarry Smith #undef __FUNCT__ 2002cd6534aSBarry Smith #define __FUNCT__ "MatMPIAIJDiagonalScaleLocalSetUp" 201dfbe8321SBarry Smith PetscErrorCode MatMPIAIJDiagonalScaleLocalSetUp(Mat inA,Vec scale) 2022cd6534aSBarry Smith { 2032cd6534aSBarry Smith Mat_MPIAIJ *ina = (Mat_MPIAIJ*) inA->data; /*access private part of matrix */ 204dfbe8321SBarry Smith PetscErrorCode ierr; 205*b1d57f15SBarry Smith PetscInt i,n,nt,cstart,cend,no,*garray = ina->garray,*lindices; 206*b1d57f15SBarry Smith PetscInt *r_rmapd,*r_rmapo; 2072cd6534aSBarry Smith 2082cd6534aSBarry Smith PetscFunctionBegin; 2092cd6534aSBarry Smith ierr = MatGetOwnershipRange(inA,&cstart,&cend);CHKERRQ(ierr); 2102cd6534aSBarry Smith ierr = MatGetSize(ina->A,PETSC_NULL,&n);CHKERRQ(ierr); 211*b1d57f15SBarry Smith ierr = PetscMalloc((inA->mapping->n+1)*sizeof(PetscInt),&r_rmapd);CHKERRQ(ierr); 212*b1d57f15SBarry Smith ierr = PetscMemzero(r_rmapd,inA->mapping->n*sizeof(PetscInt));CHKERRQ(ierr); 2132cd6534aSBarry Smith nt = 0; 2142cd6534aSBarry Smith for (i=0; i<inA->mapping->n; i++) { 2152cd6534aSBarry Smith if (inA->mapping->indices[i] >= cstart && inA->mapping->indices[i] < cend) { 2162cd6534aSBarry Smith nt++; 2172cd6534aSBarry Smith r_rmapd[i] = inA->mapping->indices[i] + 1; 2182cd6534aSBarry Smith } 2192cd6534aSBarry Smith } 2202cd6534aSBarry Smith if (nt != n) SETERRQ2(1,"Hmm nt %d n %d",nt,n); 221*b1d57f15SBarry Smith ierr = PetscMalloc((n+1)*sizeof(PetscInt),&auglyrmapd);CHKERRQ(ierr); 2222cd6534aSBarry Smith for (i=0; i<inA->mapping->n; i++) { 2232cd6534aSBarry Smith if (r_rmapd[i]){ 2242cd6534aSBarry Smith auglyrmapd[(r_rmapd[i]-1)-cstart] = i; 2252cd6534aSBarry Smith } 2262cd6534aSBarry Smith } 2272cd6534aSBarry Smith ierr = PetscFree(r_rmapd);CHKERRQ(ierr); 2282cd6534aSBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,n,&auglydd);CHKERRQ(ierr); 2292cd6534aSBarry Smith 230*b1d57f15SBarry Smith ierr = PetscMalloc((inA->N+1)*sizeof(PetscInt),&lindices);CHKERRQ(ierr); 231*b1d57f15SBarry Smith ierr = PetscMemzero(lindices,inA->N*sizeof(PetscInt));CHKERRQ(ierr); 2322cd6534aSBarry Smith for (i=0; i<ina->B->n; i++) { 2332cd6534aSBarry Smith lindices[garray[i]] = i+1; 2342cd6534aSBarry Smith } 2352cd6534aSBarry Smith no = inA->mapping->n - nt; 236*b1d57f15SBarry Smith ierr = PetscMalloc((inA->mapping->n+1)*sizeof(PetscInt),&r_rmapo);CHKERRQ(ierr); 237*b1d57f15SBarry Smith ierr = PetscMemzero(r_rmapo,inA->mapping->n*sizeof(PetscInt));CHKERRQ(ierr); 2382cd6534aSBarry Smith nt = 0; 2392cd6534aSBarry Smith for (i=0; i<inA->mapping->n; i++) { 2402cd6534aSBarry Smith if (lindices[inA->mapping->indices[i]]) { 2412cd6534aSBarry Smith nt++; 2422cd6534aSBarry Smith r_rmapo[i] = lindices[inA->mapping->indices[i]]; 2432cd6534aSBarry Smith } 2442cd6534aSBarry Smith } 2452cd6534aSBarry Smith if (nt > no) SETERRQ2(1,"Hmm nt %d no %d",nt,n); 2462cd6534aSBarry Smith ierr = PetscFree(lindices);CHKERRQ(ierr); 247*b1d57f15SBarry Smith ierr = PetscMalloc((nt+1)*sizeof(PetscInt),&auglyrmapo);CHKERRQ(ierr); 2482cd6534aSBarry Smith for (i=0; i<inA->mapping->n; i++) { 2492cd6534aSBarry Smith if (r_rmapo[i]){ 2502cd6534aSBarry Smith auglyrmapo[(r_rmapo[i]-1)] = i; 2512cd6534aSBarry Smith } 2522cd6534aSBarry Smith } 2532cd6534aSBarry Smith ierr = PetscFree(r_rmapo);CHKERRQ(ierr); 2542cd6534aSBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,nt,&auglyoo);CHKERRQ(ierr); 2552cd6534aSBarry Smith 2562cd6534aSBarry Smith PetscFunctionReturn(0); 2572cd6534aSBarry Smith } 2582cd6534aSBarry Smith 2592cd6534aSBarry Smith #undef __FUNCT__ 2602cd6534aSBarry Smith #define __FUNCT__ "MatMPIAIJDiagonalScaleLocal" 261dfbe8321SBarry Smith PetscErrorCode MatMPIAIJDiagonalScaleLocal(Mat A,Vec scale) 2622cd6534aSBarry Smith { 26392b32695SKris Buschelman /* This routine should really be abandoned as it duplicates MatDiagonalScaleLocal */ 264dfbe8321SBarry Smith PetscErrorCode ierr,(*f)(Mat,Vec); 26592b32695SKris Buschelman 26692b32695SKris Buschelman PetscFunctionBegin; 26792b32695SKris Buschelman ierr = PetscObjectQueryFunction((PetscObject)A,"MatDiagonalScaleLocal_C",(void (**)(void))&f);CHKERRQ(ierr); 26892b32695SKris Buschelman if (f) { 26992b32695SKris Buschelman ierr = (*f)(A,scale);CHKERRQ(ierr); 27092b32695SKris Buschelman } 27192b32695SKris Buschelman PetscFunctionReturn(0); 27292b32695SKris Buschelman } 27392b32695SKris Buschelman 27492b32695SKris Buschelman EXTERN_C_BEGIN 27592b32695SKris Buschelman #undef __FUNCT__ 27692b32695SKris Buschelman #define __FUNCT__ "MatDiagonalScaleLocal_MPIAIJ" 277dfbe8321SBarry Smith PetscErrorCode MatDiagonalScaleLocal_MPIAIJ(Mat A,Vec scale) 27892b32695SKris Buschelman { 2792cd6534aSBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ*) A->data; /*access private part of matrix */ 280dfbe8321SBarry Smith PetscErrorCode ierr; 281*b1d57f15SBarry Smith PetscInt n,i; 2822cd6534aSBarry Smith PetscScalar *d,*o,*s; 2832cd6534aSBarry Smith 2842cd6534aSBarry Smith PetscFunctionBegin; 2852cd6534aSBarry Smith if (!auglyrmapd) { 2862cd6534aSBarry Smith ierr = MatMPIAIJDiagonalScaleLocalSetUp(A,scale);CHKERRQ(ierr); 2872cd6534aSBarry Smith } 2882cd6534aSBarry Smith 2891ebc52fbSHong Zhang ierr = VecGetArray(scale,&s);CHKERRQ(ierr); 2902cd6534aSBarry Smith 2912cd6534aSBarry Smith ierr = VecGetLocalSize(auglydd,&n);CHKERRQ(ierr); 2921ebc52fbSHong Zhang ierr = VecGetArray(auglydd,&d);CHKERRQ(ierr); 2932cd6534aSBarry Smith for (i=0; i<n; i++) { 2942cd6534aSBarry Smith d[i] = s[auglyrmapd[i]]; /* copy "diagonal" (true local) portion of scale into dd vector */ 2952cd6534aSBarry Smith } 2961ebc52fbSHong Zhang ierr = VecRestoreArray(auglydd,&d);CHKERRQ(ierr); 2972cd6534aSBarry Smith /* column scale "diagonal" portion of local matrix */ 2982cd6534aSBarry Smith ierr = MatDiagonalScale(a->A,PETSC_NULL,auglydd);CHKERRQ(ierr); 2992cd6534aSBarry Smith 3002cd6534aSBarry Smith ierr = VecGetLocalSize(auglyoo,&n);CHKERRQ(ierr); 3011ebc52fbSHong Zhang ierr = VecGetArray(auglyoo,&o);CHKERRQ(ierr); 3022cd6534aSBarry Smith for (i=0; i<n; i++) { 3032cd6534aSBarry Smith o[i] = s[auglyrmapo[i]]; /* copy "off-diagonal" portion of scale into oo vector */ 3042cd6534aSBarry Smith } 3051ebc52fbSHong Zhang ierr = VecRestoreArray(scale,&s);CHKERRQ(ierr); 3061ebc52fbSHong Zhang ierr = VecRestoreArray(auglyoo,&o);CHKERRQ(ierr); 3072cd6534aSBarry Smith /* column scale "off-diagonal" portion of local matrix */ 3082cd6534aSBarry Smith ierr = MatDiagonalScale(a->B,PETSC_NULL,auglyoo);CHKERRQ(ierr); 3092cd6534aSBarry Smith 3102cd6534aSBarry Smith PetscFunctionReturn(0); 3112cd6534aSBarry Smith } 31292b32695SKris Buschelman EXTERN_C_END 3132cd6534aSBarry Smith 31448b35521SBarry Smith 315