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> 6af0996ceSBarry Smith #include <petsc/private/isimpl.h> /* needed because accesses data structure of ISLocalToGlobalMapping directly */ 78c79f6d3SBarry Smith 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; 13b1d57f15SBarry 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; 19b1d57f15SBarry Smith PetscInt gid,lid; 206f531f54SSatish Balay #else 21d0f46423SBarry Smith PetscInt N = mat->cmap->N,*indices; 222066d6f7SSatish Balay #endif 232066d6f7SSatish Balay 243a40ed3dSBarry Smith PetscFunctionBegin; 25aa482453SBarry Smith #if defined(PETSC_USE_CTABLE) 26c5bfad50SMark F. Adams /* use a table */ 27e23dfa41SBarry Smith ierr = PetscTableCreate(aij->B->rmap->n,mat->cmap->N+1,&gid1_lid1);CHKERRQ(ierr); 28d0f46423SBarry Smith for (i=0; i<aij->B->rmap->n; i++) { 292066d6f7SSatish Balay for (j=0; j<B->ilen[i]; j++) { 30b1d57f15SBarry Smith PetscInt data,gid1 = aj[B->i[i] + j] + 1; 310f5bd95cSBarry Smith ierr = PetscTableFind(gid1_lid1,gid1,&data);CHKERRQ(ierr); 32fa46199cSSatish Balay if (!data) { 332066d6f7SSatish Balay /* one based table */ 343861aac3SJed Brown ierr = PetscTableAdd(gid1_lid1,gid1,++ec,INSERT_VALUES);CHKERRQ(ierr); 352066d6f7SSatish Balay } 362066d6f7SSatish Balay } 372066d6f7SSatish Balay } 382066d6f7SSatish Balay /* form array of columns we need */ 39854ce69bSBarry Smith ierr = PetscMalloc1(ec+1,&garray);CHKERRQ(ierr); 400f5bd95cSBarry Smith ierr = PetscTableGetHeadPosition(gid1_lid1,&tpos);CHKERRQ(ierr); 412066d6f7SSatish Balay while (tpos) { 420f5bd95cSBarry Smith ierr = PetscTableGetNext(gid1_lid1,&tpos,&gid,&lid);CHKERRQ(ierr); 43b0a32e0cSBarry Smith gid--; 44b0a32e0cSBarry Smith lid--; 452066d6f7SSatish Balay garray[lid] = gid; 462066d6f7SSatish Balay } 470064e2bbSSatish Balay ierr = PetscSortInt(ec,garray);CHKERRQ(ierr); /* sort, and rebuild */ 480f5bd95cSBarry Smith ierr = PetscTableRemoveAll(gid1_lid1);CHKERRQ(ierr); 492066d6f7SSatish Balay for (i=0; i<ec; i++) { 503861aac3SJed Brown ierr = PetscTableAdd(gid1_lid1,garray[i]+1,i+1,INSERT_VALUES);CHKERRQ(ierr); 512066d6f7SSatish Balay } 522066d6f7SSatish Balay /* compact out the extra columns in B */ 53d0f46423SBarry Smith for (i=0; i<aij->B->rmap->n; i++) { 542066d6f7SSatish Balay for (j=0; j<B->ilen[i]; j++) { 55b1d57f15SBarry Smith PetscInt gid1 = aj[B->i[i] + j] + 1; 560f5bd95cSBarry Smith ierr = PetscTableFind(gid1_lid1,gid1,&lid);CHKERRQ(ierr); 57fa46199cSSatish Balay lid--; 58b3fb0a6cSMatthew Knepley aj[B->i[i] + j] = lid; 592066d6f7SSatish Balay } 602066d6f7SSatish Balay } 61d0f46423SBarry Smith aij->B->cmap->n = aij->B->cmap->N = ec; 62cdce4254SBarry Smith aij->B->cmap->bs = 1; 632205254eSKarl Rupp 6426283091SBarry Smith ierr = PetscLayoutSetUp((aij->B->cmap));CHKERRQ(ierr); 656bc0bbbfSBarry Smith ierr = PetscTableDestroy(&gid1_lid1);CHKERRQ(ierr); 662066d6f7SSatish Balay #else 6711285404SBarry Smith /* Make an array as long as the number of columns */ 681eb62cbbSBarry Smith /* mark those columns that are in aij->B */ 691795a4d1SJed Brown ierr = PetscCalloc1(N+1,&indices);CHKERRQ(ierr); 70d0f46423SBarry Smith for (i=0; i<aij->B->rmap->n; 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 */ 78854ce69bSBarry Smith ierr = PetscMalloc1(ec+1,&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 */ 90d0f46423SBarry Smith for (i=0; i<aij->B->rmap->n; 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 } 95d0f46423SBarry Smith aij->B->cmap->n = aij->B->cmap->N = ec; 96cd0e7f71SBarry Smith aij->B->cmap->bs = 1; 972205254eSKarl Rupp 9826283091SBarry Smith ierr = PetscLayoutSetUp((aij->B->cmap));CHKERRQ(ierr); 99606d414cSSatish Balay ierr = PetscFree(indices);CHKERRQ(ierr); 1002066d6f7SSatish Balay #endif 1011eb62cbbSBarry Smith /* create local vector that is used to scatter into */ 102029af93fSBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,ec,&aij->lvec);CHKERRQ(ierr); 1031eb62cbbSBarry Smith 104d6dfbf8fSBarry Smith /* create two temporary Index sets for build scatter gather */ 10570b3c8c7SBarry Smith ierr = ISCreateGeneral(((PetscObject)mat)->comm,ec,garray,PETSC_COPY_VALUES,&from);CHKERRQ(ierr); 106c5bfad50SMark F. Adams 107029af93fSBarry Smith ierr = ISCreateStride(PETSC_COMM_SELF,ec,0,1,&to);CHKERRQ(ierr); 1081eb62cbbSBarry Smith 1091eb62cbbSBarry Smith /* create temporary global vector to generate scatter context */ 110b5eb4454SBarry Smith /* This does not allocate the array's memory so is efficient */ 111ce94432eSBarry Smith ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)mat),1,mat->cmap->n,mat->cmap->N,NULL,&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); 1153bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)aij->Mvctx);CHKERRQ(ierr); 1163bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)aij->lvec);CHKERRQ(ierr); 1173bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)from);CHKERRQ(ierr); 1183bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)to);CHKERRQ(ierr); 1192205254eSKarl Rupp 1209e25ed09SBarry Smith aij->garray = garray; 1212205254eSKarl Rupp 1223bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)mat,(ec+1)*sizeof(PetscInt));CHKERRQ(ierr); 1236bf464f9SBarry Smith ierr = ISDestroy(&from);CHKERRQ(ierr); 1246bf464f9SBarry Smith ierr = ISDestroy(&to);CHKERRQ(ierr); 1256bf464f9SBarry Smith ierr = VecDestroy(&gvec);CHKERRQ(ierr); 1263a40ed3dSBarry Smith PetscFunctionReturn(0); 1278c79f6d3SBarry Smith } 1289e25ed09SBarry Smith 1299e25ed09SBarry Smith 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 */ 139ab9863d7SBarry Smith PetscErrorCode MatDisAssemble_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; 145d0f46423SBarry Smith PetscInt i,j,m = B->rmap->n,n = A->cmap->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 */ 1515e1f6667SBarry Smith ierr = VecDestroy(&aij->lvec);CHKERRQ(ierr); 1525e1f6667SBarry Smith ierr = VecScatterDestroy(&aij->Mvctx);CHKERRQ(ierr); 153464493b3SBarry Smith if (aij->colmap) { 154aa482453SBarry Smith #if defined(PETSC_USE_CTABLE) 1556bc0bbbfSBarry Smith ierr = PetscTableDestroy(&aij->colmap);CHKERRQ(ierr); 1562066d6f7SSatish Balay #else 157606d414cSSatish Balay ierr = PetscFree(aij->colmap);CHKERRQ(ierr); 1583bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)A,-aij->B->cmap->n*sizeof(PetscInt));CHKERRQ(ierr); 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 */ 167854ce69bSBarry Smith ierr = PetscMalloc1(m+1,&nz);CHKERRQ(ierr); 16848b35521SBarry Smith for (i=0; i<m; i++) { 16948b35521SBarry Smith nz[i] = Baij->i[i+1] - Baij->i[i]; 17048b35521SBarry Smith } 171f69a0ea3SMatthew Knepley ierr = MatCreate(PETSC_COMM_SELF,&Bnew);CHKERRQ(ierr); 172f69a0ea3SMatthew Knepley ierr = MatSetSizes(Bnew,m,n,m,n);CHKERRQ(ierr); 17333d57670SJed Brown ierr = MatSetBlockSizesFromMats(Bnew,A,A);CHKERRQ(ierr); 1747adad957SLisandro Dalcin ierr = MatSetType(Bnew,((PetscObject)B)->type_name);CHKERRQ(ierr); 175f204ca49SKris Buschelman ierr = MatSeqAIJSetPreallocation(Bnew,0,nz);CHKERRQ(ierr); 1762205254eSKarl Rupp 1772576faa2SJed Brown ((Mat_SeqAIJ*)Bnew->data)->nonew = Baij->nonew; /* Inherit insertion error options. */ 17877341eacSDmitry Karpeev /* 17977341eacSDmitry Karpeev Ensure that B's nonzerostate is monotonically increasing. 18077341eacSDmitry Karpeev Or should this follow the MatSetValues() loop to preserve B's nonzerstate across a MatDisAssemble() call? 18177341eacSDmitry Karpeev */ 182f69fde56SShane Stafford Bnew->nonzerostate = B->nonzerostate; 1832205254eSKarl Rupp 184606d414cSSatish Balay ierr = PetscFree(nz);CHKERRQ(ierr); 1852493cbb0SBarry Smith for (i=0; i<m; i++) { 186bfec09a0SHong Zhang for (j=Baij->i[i]; j<Baij->i[i+1]; j++) { 187bfec09a0SHong Zhang col = garray[Baij->j[ct]]; 1882493cbb0SBarry Smith v = Baij->a[ct++]; 18983271157SBarry Smith ierr = MatSetValues(Bnew,1,&i,1,&col,&v,B->insertmode);CHKERRQ(ierr); 1902493cbb0SBarry Smith } 1912493cbb0SBarry Smith } 192606d414cSSatish Balay ierr = PetscFree(aij->garray);CHKERRQ(ierr); 1933bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)A,-ec*sizeof(PetscInt));CHKERRQ(ierr); 1946bf464f9SBarry Smith ierr = MatDestroy(&B);CHKERRQ(ierr); 1953bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)A,(PetscObject)Bnew);CHKERRQ(ierr); 1962205254eSKarl Rupp 1972493cbb0SBarry Smith aij->B = Bnew; 198227d817aSBarry Smith A->was_assembled = PETSC_FALSE; 1993a40ed3dSBarry Smith PetscFunctionReturn(0); 2002493cbb0SBarry Smith } 2012493cbb0SBarry Smith 2022cd6534aSBarry Smith /* ugly stuff added for Glenn someday we should fix this up */ 2032cd6534aSBarry Smith 2042205254eSKarl Rupp static PetscInt *auglyrmapd = 0,*auglyrmapo = 0; /* mapping from the local ordering to the "diagonal" and "off-diagonal" parts of the local matrix */ 2052cd6534aSBarry Smith static Vec auglydd = 0,auglyoo = 0; /* work vectors used to scale the two parts of the local matrix */ 2062cd6534aSBarry Smith 2072cd6534aSBarry Smith 208dfbe8321SBarry Smith PetscErrorCode MatMPIAIJDiagonalScaleLocalSetUp(Mat inA,Vec scale) 2092cd6534aSBarry Smith { 2102cd6534aSBarry Smith Mat_MPIAIJ *ina = (Mat_MPIAIJ*) inA->data; /*access private part of matrix */ 211dfbe8321SBarry Smith PetscErrorCode ierr; 212b1d57f15SBarry Smith PetscInt i,n,nt,cstart,cend,no,*garray = ina->garray,*lindices; 213b1d57f15SBarry Smith PetscInt *r_rmapd,*r_rmapo; 2142cd6534aSBarry Smith 2152cd6534aSBarry Smith PetscFunctionBegin; 2162cd6534aSBarry Smith ierr = MatGetOwnershipRange(inA,&cstart,&cend);CHKERRQ(ierr); 2170298fd71SBarry Smith ierr = MatGetSize(ina->A,NULL,&n);CHKERRQ(ierr); 2181795a4d1SJed Brown ierr = PetscCalloc1(inA->rmap->mapping->n+1,&r_rmapd);CHKERRQ(ierr); 2192cd6534aSBarry Smith nt = 0; 220992144d0SBarry Smith for (i=0; i<inA->rmap->mapping->n; i++) { 221992144d0SBarry Smith if (inA->rmap->mapping->indices[i] >= cstart && inA->rmap->mapping->indices[i] < cend) { 2222cd6534aSBarry Smith nt++; 223992144d0SBarry Smith r_rmapd[i] = inA->rmap->mapping->indices[i] + 1; 2242cd6534aSBarry Smith } 2252cd6534aSBarry Smith } 226e32f2f54SBarry Smith if (nt != n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hmm nt %D n %D",nt,n); 227854ce69bSBarry Smith ierr = PetscMalloc1(n+1,&auglyrmapd);CHKERRQ(ierr); 228992144d0SBarry Smith for (i=0; i<inA->rmap->mapping->n; i++) { 2292cd6534aSBarry Smith if (r_rmapd[i]) { 2302cd6534aSBarry Smith auglyrmapd[(r_rmapd[i]-1)-cstart] = i; 2312cd6534aSBarry Smith } 2322cd6534aSBarry Smith } 2332cd6534aSBarry Smith ierr = PetscFree(r_rmapd);CHKERRQ(ierr); 2342cd6534aSBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,n,&auglydd);CHKERRQ(ierr); 2352cd6534aSBarry Smith 2361795a4d1SJed Brown ierr = PetscCalloc1(inA->cmap->N+1,&lindices);CHKERRQ(ierr); 237d0f46423SBarry Smith for (i=0; i<ina->B->cmap->n; i++) { 2382cd6534aSBarry Smith lindices[garray[i]] = i+1; 2392cd6534aSBarry Smith } 240992144d0SBarry Smith no = inA->rmap->mapping->n - nt; 2411795a4d1SJed Brown ierr = PetscCalloc1(inA->rmap->mapping->n+1,&r_rmapo);CHKERRQ(ierr); 2422cd6534aSBarry Smith nt = 0; 243992144d0SBarry Smith for (i=0; i<inA->rmap->mapping->n; i++) { 244992144d0SBarry Smith if (lindices[inA->rmap->mapping->indices[i]]) { 2452cd6534aSBarry Smith nt++; 246992144d0SBarry Smith r_rmapo[i] = lindices[inA->rmap->mapping->indices[i]]; 2472cd6534aSBarry Smith } 2482cd6534aSBarry Smith } 249e32f2f54SBarry Smith if (nt > no) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hmm nt %D no %D",nt,n); 2502cd6534aSBarry Smith ierr = PetscFree(lindices);CHKERRQ(ierr); 251854ce69bSBarry Smith ierr = PetscMalloc1(nt+1,&auglyrmapo);CHKERRQ(ierr); 252992144d0SBarry Smith for (i=0; i<inA->rmap->mapping->n; i++) { 2532cd6534aSBarry Smith if (r_rmapo[i]) { 2542cd6534aSBarry Smith auglyrmapo[(r_rmapo[i]-1)] = i; 2552cd6534aSBarry Smith } 2562cd6534aSBarry Smith } 2572cd6534aSBarry Smith ierr = PetscFree(r_rmapo);CHKERRQ(ierr); 2582cd6534aSBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,nt,&auglyoo);CHKERRQ(ierr); 2592cd6534aSBarry Smith PetscFunctionReturn(0); 2602cd6534aSBarry Smith } 2612cd6534aSBarry Smith 262dfbe8321SBarry Smith PetscErrorCode MatMPIAIJDiagonalScaleLocal(Mat A,Vec scale) 2632cd6534aSBarry Smith { 26492b32695SKris Buschelman /* This routine should really be abandoned as it duplicates MatDiagonalScaleLocal */ 2654ac538c5SBarry Smith PetscErrorCode ierr; 26692b32695SKris Buschelman 26792b32695SKris Buschelman PetscFunctionBegin; 2684ac538c5SBarry Smith ierr = PetscTryMethod(A,"MatDiagonalScaleLocal_C",(Mat,Vec),(A,scale));CHKERRQ(ierr); 26992b32695SKris Buschelman PetscFunctionReturn(0); 27092b32695SKris Buschelman } 27192b32695SKris Buschelman 2727087cfbeSBarry Smith PetscErrorCode MatDiagonalScaleLocal_MPIAIJ(Mat A,Vec scale) 27392b32695SKris Buschelman { 2742cd6534aSBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ*) A->data; /*access private part of matrix */ 275dfbe8321SBarry Smith PetscErrorCode ierr; 276b1d57f15SBarry Smith PetscInt n,i; 277*bca11509SBarry Smith PetscScalar *d,*o; 278*bca11509SBarry Smith const PetscScalar *s; 2792cd6534aSBarry Smith 2802cd6534aSBarry Smith PetscFunctionBegin; 2812cd6534aSBarry Smith if (!auglyrmapd) { 2822cd6534aSBarry Smith ierr = MatMPIAIJDiagonalScaleLocalSetUp(A,scale);CHKERRQ(ierr); 2832cd6534aSBarry Smith } 2842cd6534aSBarry Smith 285*bca11509SBarry Smith ierr = VecGetArrayRead(scale,&s);CHKERRQ(ierr); 2862cd6534aSBarry Smith 2872cd6534aSBarry Smith ierr = VecGetLocalSize(auglydd,&n);CHKERRQ(ierr); 2881ebc52fbSHong Zhang ierr = VecGetArray(auglydd,&d);CHKERRQ(ierr); 2892cd6534aSBarry Smith for (i=0; i<n; i++) { 2902cd6534aSBarry Smith d[i] = s[auglyrmapd[i]]; /* copy "diagonal" (true local) portion of scale into dd vector */ 2912cd6534aSBarry Smith } 2921ebc52fbSHong Zhang ierr = VecRestoreArray(auglydd,&d);CHKERRQ(ierr); 2932cd6534aSBarry Smith /* column scale "diagonal" portion of local matrix */ 2940298fd71SBarry Smith ierr = MatDiagonalScale(a->A,NULL,auglydd);CHKERRQ(ierr); 2952cd6534aSBarry Smith 2962cd6534aSBarry Smith ierr = VecGetLocalSize(auglyoo,&n);CHKERRQ(ierr); 2971ebc52fbSHong Zhang ierr = VecGetArray(auglyoo,&o);CHKERRQ(ierr); 2982cd6534aSBarry Smith for (i=0; i<n; i++) { 2992cd6534aSBarry Smith o[i] = s[auglyrmapo[i]]; /* copy "off-diagonal" portion of scale into oo vector */ 3002cd6534aSBarry Smith } 301*bca11509SBarry Smith ierr = VecRestoreArrayRead(scale,&s);CHKERRQ(ierr); 3021ebc52fbSHong Zhang ierr = VecRestoreArray(auglyoo,&o);CHKERRQ(ierr); 3032cd6534aSBarry Smith /* column scale "off-diagonal" portion of local matrix */ 3040298fd71SBarry Smith ierr = MatDiagonalScale(a->B,NULL,auglyoo);CHKERRQ(ierr); 3052cd6534aSBarry Smith PetscFunctionReturn(0); 3062cd6534aSBarry Smith } 307