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> 6aeda0f58SHong Zhang #include <petsc/private/vecimpl.h> 7af0996ceSBarry Smith #include <petsc/private/isimpl.h> /* needed because accesses data structure of ISLocalToGlobalMapping directly */ 88c79f6d3SBarry Smith 9dfbe8321SBarry Smith PetscErrorCode MatSetUpMultiply_MPIAIJ(Mat mat) 108c79f6d3SBarry Smith { 1144a69424SLois Curfman McInnes Mat_MPIAIJ *aij = (Mat_MPIAIJ*)mat->data; 12ec8511deSBarry Smith Mat_SeqAIJ *B = (Mat_SeqAIJ*)(aij->B->data); 136849ba73SBarry Smith PetscErrorCode ierr; 14b1d57f15SBarry Smith PetscInt i,j,*aj = B->j,ec = 0,*garray; 151eb62cbbSBarry Smith IS from,to; 161eb62cbbSBarry Smith Vec gvec; 17aa482453SBarry Smith #if defined(PETSC_USE_CTABLE) 180f5bd95cSBarry Smith PetscTable gid1_lid1; 190f5bd95cSBarry Smith PetscTablePosition tpos; 20b1d57f15SBarry Smith PetscInt gid,lid; 216f531f54SSatish Balay #else 22d0f46423SBarry Smith PetscInt N = mat->cmap->N,*indices; 232066d6f7SSatish Balay #endif 242066d6f7SSatish Balay 253a40ed3dSBarry Smith PetscFunctionBegin; 264b8d542aSHong Zhang if (!aij->garray) { 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 */ 41854ce69bSBarry Smith 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 */ 80854ce69bSBarry Smith 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 1034b8d542aSHong Zhang } else { 1044b8d542aSHong Zhang garray = aij->garray; 1054b8d542aSHong Zhang } 1064b8d542aSHong Zhang 1074b8d542aSHong Zhang if (!aij->lvec) { 1081eb62cbbSBarry Smith /* create local vector that is used to scatter into */ 109029af93fSBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,ec,&aij->lvec);CHKERRQ(ierr); 1104b8d542aSHong Zhang } else { 1114b8d542aSHong Zhang ierr = VecGetSize(aij->lvec,&ec);CHKERRQ(ierr); 1124b8d542aSHong Zhang } 1131eb62cbbSBarry Smith 114d6dfbf8fSBarry Smith /* create two temporary Index sets for build scatter gather */ 11570b3c8c7SBarry Smith ierr = ISCreateGeneral(((PetscObject)mat)->comm,ec,garray,PETSC_COPY_VALUES,&from);CHKERRQ(ierr); 116c5bfad50SMark F. Adams 117029af93fSBarry Smith ierr = ISCreateStride(PETSC_COMM_SELF,ec,0,1,&to);CHKERRQ(ierr); 1181eb62cbbSBarry Smith 1191eb62cbbSBarry Smith /* create temporary global vector to generate scatter context */ 120b5eb4454SBarry Smith /* This does not allocate the array's memory so is efficient */ 121ce94432eSBarry Smith ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)mat),1,mat->cmap->n,mat->cmap->N,NULL,&gvec);CHKERRQ(ierr); 1221eb62cbbSBarry Smith 1232d336d48SLois Curfman McInnes /* generate the scatter context */ 124a78d8160SHong Zhang if (aij->Mvctx_mpi1_flg) { 12501ad2aeeSHong Zhang ierr = VecScatterDestroy(&aij->Mvctx_mpi1);CHKERRQ(ierr); 126*803a1b88SHong Zhang ierr = VecScatterCreate(gvec,from,aij->lvec,to,&aij->Mvctx_mpi1);CHKERRQ(ierr); 127*803a1b88SHong Zhang ierr = VecScatterSetType(aij->Mvctx_mpi1,VECSCATTERMPI1);CHKERRQ(ierr); 1284b8d542aSHong Zhang ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)aij->Mvctx_mpi1);CHKERRQ(ierr); 1294b8d542aSHong Zhang } else { 1303f6a6bdaSHong Zhang ierr = VecScatterDestroy(&aij->Mvctx);CHKERRQ(ierr); 13108480c60SBarry Smith ierr = VecScatterCreate(gvec,from,aij->lvec,to,&aij->Mvctx);CHKERRQ(ierr); 1323bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)aij->Mvctx);CHKERRQ(ierr); 1333bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)aij->lvec);CHKERRQ(ierr); 1344b8d542aSHong Zhang ierr = PetscLogObjectMemory((PetscObject)mat,(ec+1)*sizeof(PetscInt));CHKERRQ(ierr); 1354b8d542aSHong Zhang } 13667bb5161SHong Zhang aij->garray = garray; 1374b8d542aSHong Zhang 1383bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)from);CHKERRQ(ierr); 1393bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)to);CHKERRQ(ierr); 1402205254eSKarl Rupp 1416bf464f9SBarry Smith ierr = ISDestroy(&from);CHKERRQ(ierr); 1426bf464f9SBarry Smith ierr = ISDestroy(&to);CHKERRQ(ierr); 1436bf464f9SBarry Smith ierr = VecDestroy(&gvec);CHKERRQ(ierr); 1443a40ed3dSBarry Smith PetscFunctionReturn(0); 1458c79f6d3SBarry Smith } 1469e25ed09SBarry Smith 1472493cbb0SBarry Smith /* 1482493cbb0SBarry Smith Takes the local part of an already assembled MPIAIJ matrix 1492493cbb0SBarry Smith and disassembles it. This is to allow new nonzeros into the matrix 1502493cbb0SBarry Smith that require more communication in the matrix vector multiply. 1512493cbb0SBarry Smith Thus certain data-structures must be rebuilt. 1522493cbb0SBarry Smith 1532493cbb0SBarry Smith Kind of slow! But that's what application programmers get when 1542493cbb0SBarry Smith they are sloppy. 1552493cbb0SBarry Smith */ 156ab9863d7SBarry Smith PetscErrorCode MatDisAssemble_MPIAIJ(Mat A) 1572493cbb0SBarry Smith { 1582493cbb0SBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ*)A->data; 1592493cbb0SBarry Smith Mat B = aij->B,Bnew; 160ec8511deSBarry Smith Mat_SeqAIJ *Baij = (Mat_SeqAIJ*)B->data; 161dfbe8321SBarry Smith PetscErrorCode ierr; 162d0f46423SBarry Smith PetscInt i,j,m = B->rmap->n,n = A->cmap->N,col,ct = 0,*garray = aij->garray,*nz,ec; 16387828ca2SBarry Smith PetscScalar v; 1642493cbb0SBarry Smith 1653a40ed3dSBarry Smith PetscFunctionBegin; 1662493cbb0SBarry Smith /* free stuff related to matrix-vec multiply */ 167b0a32e0cSBarry Smith ierr = VecGetSize(aij->lvec,&ec);CHKERRQ(ierr); /* needed for PetscLogObjectMemory below */ 1685e1f6667SBarry Smith ierr = VecDestroy(&aij->lvec);CHKERRQ(ierr); 169464493b3SBarry Smith if (aij->colmap) { 170aa482453SBarry Smith #if defined(PETSC_USE_CTABLE) 1716bc0bbbfSBarry Smith ierr = PetscTableDestroy(&aij->colmap);CHKERRQ(ierr); 1722066d6f7SSatish Balay #else 173606d414cSSatish Balay ierr = PetscFree(aij->colmap);CHKERRQ(ierr); 1743bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)A,-aij->B->cmap->n*sizeof(PetscInt));CHKERRQ(ierr); 1752066d6f7SSatish Balay #endif 176464493b3SBarry Smith } 1772493cbb0SBarry Smith 1782493cbb0SBarry Smith /* make sure that B is assembled so we can access its values */ 1796d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 180fe2f2677SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1812493cbb0SBarry Smith 1822493cbb0SBarry Smith /* invent new B and copy stuff over */ 183854ce69bSBarry Smith ierr = PetscMalloc1(m+1,&nz);CHKERRQ(ierr); 18448b35521SBarry Smith for (i=0; i<m; i++) { 18548b35521SBarry Smith nz[i] = Baij->i[i+1] - Baij->i[i]; 18648b35521SBarry Smith } 187f69a0ea3SMatthew Knepley ierr = MatCreate(PETSC_COMM_SELF,&Bnew);CHKERRQ(ierr); 188f69a0ea3SMatthew Knepley ierr = MatSetSizes(Bnew,m,n,m,n);CHKERRQ(ierr); 18933d57670SJed Brown ierr = MatSetBlockSizesFromMats(Bnew,A,A);CHKERRQ(ierr); 1907adad957SLisandro Dalcin ierr = MatSetType(Bnew,((PetscObject)B)->type_name);CHKERRQ(ierr); 191f204ca49SKris Buschelman ierr = MatSeqAIJSetPreallocation(Bnew,0,nz);CHKERRQ(ierr); 1922205254eSKarl Rupp 1932576faa2SJed Brown ((Mat_SeqAIJ*)Bnew->data)->nonew = Baij->nonew; /* Inherit insertion error options. */ 19477341eacSDmitry Karpeev /* 19577341eacSDmitry Karpeev Ensure that B's nonzerostate is monotonically increasing. 19677341eacSDmitry Karpeev Or should this follow the MatSetValues() loop to preserve B's nonzerstate across a MatDisAssemble() call? 19777341eacSDmitry Karpeev */ 198f69fde56SShane Stafford Bnew->nonzerostate = B->nonzerostate; 1992205254eSKarl Rupp 200606d414cSSatish Balay ierr = PetscFree(nz);CHKERRQ(ierr); 2012493cbb0SBarry Smith for (i=0; i<m; i++) { 202bfec09a0SHong Zhang for (j=Baij->i[i]; j<Baij->i[i+1]; j++) { 203bfec09a0SHong Zhang col = garray[Baij->j[ct]]; 2042493cbb0SBarry Smith v = Baij->a[ct++]; 20583271157SBarry Smith ierr = MatSetValues(Bnew,1,&i,1,&col,&v,B->insertmode);CHKERRQ(ierr); 2062493cbb0SBarry Smith } 2072493cbb0SBarry Smith } 208606d414cSSatish Balay ierr = PetscFree(aij->garray);CHKERRQ(ierr); 2093bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)A,-ec*sizeof(PetscInt));CHKERRQ(ierr); 2106bf464f9SBarry Smith ierr = MatDestroy(&B);CHKERRQ(ierr); 2113bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)A,(PetscObject)Bnew);CHKERRQ(ierr); 2122205254eSKarl Rupp 2132493cbb0SBarry Smith aij->B = Bnew; 214227d817aSBarry Smith A->was_assembled = PETSC_FALSE; 2153a40ed3dSBarry Smith PetscFunctionReturn(0); 2162493cbb0SBarry Smith } 2172493cbb0SBarry Smith 2182cd6534aSBarry Smith /* ugly stuff added for Glenn someday we should fix this up */ 2192cd6534aSBarry Smith 2202205254eSKarl Rupp static PetscInt *auglyrmapd = 0,*auglyrmapo = 0; /* mapping from the local ordering to the "diagonal" and "off-diagonal" parts of the local matrix */ 2212cd6534aSBarry Smith static Vec auglydd = 0,auglyoo = 0; /* work vectors used to scale the two parts of the local matrix */ 2222cd6534aSBarry Smith 2232cd6534aSBarry Smith 224dfbe8321SBarry Smith PetscErrorCode MatMPIAIJDiagonalScaleLocalSetUp(Mat inA,Vec scale) 2252cd6534aSBarry Smith { 2262cd6534aSBarry Smith Mat_MPIAIJ *ina = (Mat_MPIAIJ*) inA->data; /*access private part of matrix */ 227dfbe8321SBarry Smith PetscErrorCode ierr; 228b1d57f15SBarry Smith PetscInt i,n,nt,cstart,cend,no,*garray = ina->garray,*lindices; 229b1d57f15SBarry Smith PetscInt *r_rmapd,*r_rmapo; 2302cd6534aSBarry Smith 2312cd6534aSBarry Smith PetscFunctionBegin; 2322cd6534aSBarry Smith ierr = MatGetOwnershipRange(inA,&cstart,&cend);CHKERRQ(ierr); 2330298fd71SBarry Smith ierr = MatGetSize(ina->A,NULL,&n);CHKERRQ(ierr); 2341795a4d1SJed Brown ierr = PetscCalloc1(inA->rmap->mapping->n+1,&r_rmapd);CHKERRQ(ierr); 2352cd6534aSBarry Smith nt = 0; 236992144d0SBarry Smith for (i=0; i<inA->rmap->mapping->n; i++) { 237992144d0SBarry Smith if (inA->rmap->mapping->indices[i] >= cstart && inA->rmap->mapping->indices[i] < cend) { 2382cd6534aSBarry Smith nt++; 239992144d0SBarry Smith r_rmapd[i] = inA->rmap->mapping->indices[i] + 1; 2402cd6534aSBarry Smith } 2412cd6534aSBarry Smith } 242e32f2f54SBarry Smith if (nt != n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hmm nt %D n %D",nt,n); 243854ce69bSBarry Smith ierr = PetscMalloc1(n+1,&auglyrmapd);CHKERRQ(ierr); 244992144d0SBarry Smith for (i=0; i<inA->rmap->mapping->n; i++) { 2452cd6534aSBarry Smith if (r_rmapd[i]) { 2462cd6534aSBarry Smith auglyrmapd[(r_rmapd[i]-1)-cstart] = i; 2472cd6534aSBarry Smith } 2482cd6534aSBarry Smith } 2492cd6534aSBarry Smith ierr = PetscFree(r_rmapd);CHKERRQ(ierr); 2502cd6534aSBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,n,&auglydd);CHKERRQ(ierr); 2512cd6534aSBarry Smith 2521795a4d1SJed Brown ierr = PetscCalloc1(inA->cmap->N+1,&lindices);CHKERRQ(ierr); 253d0f46423SBarry Smith for (i=0; i<ina->B->cmap->n; i++) { 2542cd6534aSBarry Smith lindices[garray[i]] = i+1; 2552cd6534aSBarry Smith } 256992144d0SBarry Smith no = inA->rmap->mapping->n - nt; 2571795a4d1SJed Brown ierr = PetscCalloc1(inA->rmap->mapping->n+1,&r_rmapo);CHKERRQ(ierr); 2582cd6534aSBarry Smith nt = 0; 259992144d0SBarry Smith for (i=0; i<inA->rmap->mapping->n; i++) { 260992144d0SBarry Smith if (lindices[inA->rmap->mapping->indices[i]]) { 2612cd6534aSBarry Smith nt++; 262992144d0SBarry Smith r_rmapo[i] = lindices[inA->rmap->mapping->indices[i]]; 2632cd6534aSBarry Smith } 2642cd6534aSBarry Smith } 265e32f2f54SBarry Smith if (nt > no) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hmm nt %D no %D",nt,n); 2662cd6534aSBarry Smith ierr = PetscFree(lindices);CHKERRQ(ierr); 267854ce69bSBarry Smith ierr = PetscMalloc1(nt+1,&auglyrmapo);CHKERRQ(ierr); 268992144d0SBarry Smith for (i=0; i<inA->rmap->mapping->n; i++) { 2692cd6534aSBarry Smith if (r_rmapo[i]) { 2702cd6534aSBarry Smith auglyrmapo[(r_rmapo[i]-1)] = i; 2712cd6534aSBarry Smith } 2722cd6534aSBarry Smith } 2732cd6534aSBarry Smith ierr = PetscFree(r_rmapo);CHKERRQ(ierr); 2742cd6534aSBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,nt,&auglyoo);CHKERRQ(ierr); 2752cd6534aSBarry Smith PetscFunctionReturn(0); 2762cd6534aSBarry Smith } 2772cd6534aSBarry Smith 278dfbe8321SBarry Smith PetscErrorCode MatMPIAIJDiagonalScaleLocal(Mat A,Vec scale) 2792cd6534aSBarry Smith { 28092b32695SKris Buschelman /* This routine should really be abandoned as it duplicates MatDiagonalScaleLocal */ 2814ac538c5SBarry Smith PetscErrorCode ierr; 28292b32695SKris Buschelman 28392b32695SKris Buschelman PetscFunctionBegin; 2844ac538c5SBarry Smith ierr = PetscTryMethod(A,"MatDiagonalScaleLocal_C",(Mat,Vec),(A,scale));CHKERRQ(ierr); 28592b32695SKris Buschelman PetscFunctionReturn(0); 28692b32695SKris Buschelman } 28792b32695SKris Buschelman 2887087cfbeSBarry Smith PetscErrorCode MatDiagonalScaleLocal_MPIAIJ(Mat A,Vec scale) 28992b32695SKris Buschelman { 2902cd6534aSBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ*) A->data; /*access private part of matrix */ 291dfbe8321SBarry Smith PetscErrorCode ierr; 292b1d57f15SBarry Smith PetscInt n,i; 293bca11509SBarry Smith PetscScalar *d,*o; 294bca11509SBarry Smith const PetscScalar *s; 2952cd6534aSBarry Smith 2962cd6534aSBarry Smith PetscFunctionBegin; 2972cd6534aSBarry Smith if (!auglyrmapd) { 2982cd6534aSBarry Smith ierr = MatMPIAIJDiagonalScaleLocalSetUp(A,scale);CHKERRQ(ierr); 2992cd6534aSBarry Smith } 3002cd6534aSBarry Smith 301bca11509SBarry Smith ierr = VecGetArrayRead(scale,&s);CHKERRQ(ierr); 3022cd6534aSBarry Smith 3032cd6534aSBarry Smith ierr = VecGetLocalSize(auglydd,&n);CHKERRQ(ierr); 3041ebc52fbSHong Zhang ierr = VecGetArray(auglydd,&d);CHKERRQ(ierr); 3052cd6534aSBarry Smith for (i=0; i<n; i++) { 3062cd6534aSBarry Smith d[i] = s[auglyrmapd[i]]; /* copy "diagonal" (true local) portion of scale into dd vector */ 3072cd6534aSBarry Smith } 3081ebc52fbSHong Zhang ierr = VecRestoreArray(auglydd,&d);CHKERRQ(ierr); 3092cd6534aSBarry Smith /* column scale "diagonal" portion of local matrix */ 3100298fd71SBarry Smith ierr = MatDiagonalScale(a->A,NULL,auglydd);CHKERRQ(ierr); 3112cd6534aSBarry Smith 3122cd6534aSBarry Smith ierr = VecGetLocalSize(auglyoo,&n);CHKERRQ(ierr); 3131ebc52fbSHong Zhang ierr = VecGetArray(auglyoo,&o);CHKERRQ(ierr); 3142cd6534aSBarry Smith for (i=0; i<n; i++) { 3152cd6534aSBarry Smith o[i] = s[auglyrmapo[i]]; /* copy "off-diagonal" portion of scale into oo vector */ 3162cd6534aSBarry Smith } 317bca11509SBarry Smith ierr = VecRestoreArrayRead(scale,&s);CHKERRQ(ierr); 3181ebc52fbSHong Zhang ierr = VecRestoreArray(auglyoo,&o);CHKERRQ(ierr); 3192cd6534aSBarry Smith /* column scale "off-diagonal" portion of local matrix */ 3200298fd71SBarry Smith ierr = MatDiagonalScale(a->B,NULL,auglyoo);CHKERRQ(ierr); 3212cd6534aSBarry Smith PetscFunctionReturn(0); 3222cd6534aSBarry Smith } 323