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> 6*186d4ecdSBarry Smith #include <petsc-private/isimpl.h> /* needed because accesses data structure of ISLocalToGlobalMapping directly */ 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; 18ace3abfcSBarry Smith PetscBool useblockis; 19aa482453SBarry Smith #if defined(PETSC_USE_CTABLE) 200f5bd95cSBarry Smith PetscTable gid1_lid1; 210f5bd95cSBarry Smith PetscTablePosition tpos; 22b1d57f15SBarry Smith PetscInt gid,lid; 236f531f54SSatish Balay #else 24d0f46423SBarry Smith PetscInt N = mat->cmap->N,*indices; 252066d6f7SSatish Balay #endif 262066d6f7SSatish Balay 273a40ed3dSBarry Smith PetscFunctionBegin; 28aa482453SBarry Smith #if defined(PETSC_USE_CTABLE) 29c5bfad50SMark F. Adams /* use a table */ 30e23dfa41SBarry Smith ierr = PetscTableCreate(aij->B->rmap->n,mat->cmap->N+1,&gid1_lid1);CHKERRQ(ierr); 31d0f46423SBarry Smith for (i=0; i<aij->B->rmap->n; 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 */ 373861aac3SJed Brown ierr = PetscTableAdd(gid1_lid1,gid1,++ec,INSERT_VALUES);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++) { 533861aac3SJed Brown ierr = PetscTableAdd(gid1_lid1,garray[i]+1,i+1,INSERT_VALUES);CHKERRQ(ierr); 542066d6f7SSatish Balay } 552066d6f7SSatish Balay /* compact out the extra columns in B */ 56d0f46423SBarry Smith for (i=0; i<aij->B->rmap->n; 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 } 64d0f46423SBarry Smith aij->B->cmap->n = aij->B->cmap->N = ec; 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 */ 71b1d57f15SBarry Smith ierr = PetscMalloc((N+1)*sizeof(PetscInt),&indices);CHKERRQ(ierr); 72b1d57f15SBarry Smith ierr = PetscMemzero(indices,N*sizeof(PetscInt));CHKERRQ(ierr); 73d0f46423SBarry Smith for (i=0; i<aij->B->rmap->n; i++) { 74d6dfbf8fSBarry Smith for (j=0; j<B->ilen[i]; j++) { 75bfec09a0SHong Zhang if (!indices[aj[B->i[i] + j]]) ec++; 76bfec09a0SHong Zhang indices[aj[B->i[i] + j]] = 1; 77416022c9SBarry Smith } 781eb62cbbSBarry Smith } 798c79f6d3SBarry Smith 801eb62cbbSBarry Smith /* form array of columns we need */ 81b1d57f15SBarry Smith ierr = PetscMalloc((ec+1)*sizeof(PetscInt),&garray);CHKERRQ(ierr); 821eb62cbbSBarry Smith ec = 0; 831eb62cbbSBarry Smith for (i=0; i<N; i++) { 841eb62cbbSBarry Smith if (indices[i]) garray[ec++] = i; 851eb62cbbSBarry Smith } 861eb62cbbSBarry Smith 871eb62cbbSBarry Smith /* make indices now point into garray */ 881eb62cbbSBarry Smith for (i=0; i<ec; i++) { 89bfec09a0SHong Zhang indices[garray[i]] = i; 901eb62cbbSBarry Smith } 911eb62cbbSBarry Smith 921eb62cbbSBarry Smith /* compact out the extra columns in B */ 93d0f46423SBarry Smith for (i=0; i<aij->B->rmap->n; i++) { 94d6dfbf8fSBarry Smith for (j=0; j<B->ilen[i]; j++) { 95bfec09a0SHong Zhang aj[B->i[i] + j] = indices[aj[B->i[i] + j]]; 961eb62cbbSBarry Smith } 97d6dfbf8fSBarry Smith } 98d0f46423SBarry Smith aij->B->cmap->n = aij->B->cmap->N = ec; 992205254eSKarl Rupp 10026283091SBarry Smith ierr = PetscLayoutSetUp((aij->B->cmap));CHKERRQ(ierr); 101606d414cSSatish Balay ierr = PetscFree(indices);CHKERRQ(ierr); 1022066d6f7SSatish Balay #endif 1031eb62cbbSBarry Smith /* create local vector that is used to scatter into */ 104029af93fSBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,ec,&aij->lvec);CHKERRQ(ierr); 1051eb62cbbSBarry Smith 106d6dfbf8fSBarry Smith /* create two temporary Index sets for build scatter gather */ 10703919abeSBarry Smith /* check for the special case where blocks are communicated for faster VecScatterXXX */ 108a2f3521dSMark F. Adams useblockis = PETSC_FALSE; 109c5bfad50SMark F. Adams if (mat->cmap->bs > 1) { 110c5bfad50SMark F. Adams PetscInt bs = mat->cmap->bs,ibs,ga; 11103919abeSBarry Smith if (!(ec % bs)) { 1123ba34427SBarry Smith useblockis = PETSC_TRUE; 11303919abeSBarry Smith for (i=0; i<ec/bs; i++) { 11403919abeSBarry Smith if ((ga = garray[ibs = i*bs]) % bs) { 11503919abeSBarry Smith useblockis = PETSC_FALSE; 11603919abeSBarry Smith break; 11703919abeSBarry Smith } 11803919abeSBarry Smith for (j=1; j<bs; j++) { 11903919abeSBarry Smith if (garray[ibs+j] != ga+j) { 12003919abeSBarry Smith useblockis = PETSC_FALSE; 12103919abeSBarry Smith break; 12203919abeSBarry Smith } 12303919abeSBarry Smith } 12403919abeSBarry Smith if (!useblockis) break; 12503919abeSBarry Smith } 12603919abeSBarry Smith } 12703919abeSBarry Smith } 128c5bfad50SMark F. Adams #if defined(PETSC_USE_DEBUG) 129c5bfad50SMark F. Adams i = (PetscInt)useblockis; 130ce94432eSBarry Smith ierr = MPI_Allreduce(&i,&j,1,MPIU_INT,MPI_MIN,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr); 131c5bfad50SMark F. Adams if (j!=i) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Use of blocked not consistant (I am usning blocked)"); 132c5bfad50SMark F. Adams #endif 133c5bfad50SMark F. Adams 13403919abeSBarry Smith if (useblockis) { 135c5bfad50SMark F. Adams PetscInt *ga,bs = mat->cmap->bs,iec = ec/bs; 136c5bfad50SMark F. Adams if (ec%bs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"ec=%D bs=%D",ec,bs); 13722d28d08SBarry Smith ierr = PetscInfo(mat,"Using block index set to define scatter\n");CHKERRQ(ierr); 138c5bfad50SMark F. Adams ierr = PetscMalloc(iec*sizeof(PetscInt),&ga);CHKERRQ(ierr); 139e82e9f6bSBarry Smith for (i=0; i<iec; i++) ga[i] = garray[i*bs]/bs; 140ce94432eSBarry Smith ierr = ISCreateBlock(PetscObjectComm((PetscObject)mat),bs,iec,ga,PETSC_OWN_POINTER,&from);CHKERRQ(ierr); 14103919abeSBarry Smith } else { 142ce94432eSBarry Smith ierr = ISCreateGeneral(PetscObjectComm((PetscObject)mat),ec,garray,PETSC_COPY_VALUES,&from);CHKERRQ(ierr); 14303919abeSBarry Smith } 144c5bfad50SMark F. Adams 145029af93fSBarry Smith ierr = ISCreateStride(PETSC_COMM_SELF,ec,0,1,&to);CHKERRQ(ierr); 1461eb62cbbSBarry Smith 1471eb62cbbSBarry Smith /* create temporary global vector to generate scatter context */ 148b5eb4454SBarry Smith /* This does not allocate the array's memory so is efficient */ 149ce94432eSBarry Smith ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)mat),1,mat->cmap->n,mat->cmap->N,NULL,&gvec);CHKERRQ(ierr); 1501eb62cbbSBarry Smith 1512d336d48SLois Curfman McInnes /* generate the scatter context */ 15208480c60SBarry Smith ierr = VecScatterCreate(gvec,from,aij->lvec,to,&aij->Mvctx);CHKERRQ(ierr); 15352e6d16bSBarry Smith ierr = PetscLogObjectParent(mat,aij->Mvctx);CHKERRQ(ierr); 15452e6d16bSBarry Smith ierr = PetscLogObjectParent(mat,aij->lvec);CHKERRQ(ierr); 15552e6d16bSBarry Smith ierr = PetscLogObjectParent(mat,from);CHKERRQ(ierr); 15652e6d16bSBarry Smith ierr = PetscLogObjectParent(mat,to);CHKERRQ(ierr); 1572205254eSKarl Rupp 1589e25ed09SBarry Smith aij->garray = garray; 1592205254eSKarl Rupp 16052e6d16bSBarry Smith ierr = PetscLogObjectMemory(mat,(ec+1)*sizeof(PetscInt));CHKERRQ(ierr); 1616bf464f9SBarry Smith ierr = ISDestroy(&from);CHKERRQ(ierr); 1626bf464f9SBarry Smith ierr = ISDestroy(&to);CHKERRQ(ierr); 1636bf464f9SBarry Smith ierr = VecDestroy(&gvec);CHKERRQ(ierr); 1643a40ed3dSBarry Smith PetscFunctionReturn(0); 1658c79f6d3SBarry Smith } 1669e25ed09SBarry Smith 1679e25ed09SBarry Smith 1684a2ae208SSatish Balay #undef __FUNCT__ 169ab9863d7SBarry Smith #define __FUNCT__ "MatDisAssemble_MPIAIJ" 1702493cbb0SBarry Smith /* 1712493cbb0SBarry Smith Takes the local part of an already assembled MPIAIJ matrix 1722493cbb0SBarry Smith and disassembles it. This is to allow new nonzeros into the matrix 1732493cbb0SBarry Smith that require more communication in the matrix vector multiply. 1742493cbb0SBarry Smith Thus certain data-structures must be rebuilt. 1752493cbb0SBarry Smith 1762493cbb0SBarry Smith Kind of slow! But that's what application programmers get when 1772493cbb0SBarry Smith they are sloppy. 1782493cbb0SBarry Smith */ 179ab9863d7SBarry Smith PetscErrorCode MatDisAssemble_MPIAIJ(Mat A) 1802493cbb0SBarry Smith { 1812493cbb0SBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ*)A->data; 1822493cbb0SBarry Smith Mat B = aij->B,Bnew; 183ec8511deSBarry Smith Mat_SeqAIJ *Baij = (Mat_SeqAIJ*)B->data; 184dfbe8321SBarry Smith PetscErrorCode ierr; 185d0f46423SBarry Smith PetscInt i,j,m = B->rmap->n,n = A->cmap->N,col,ct = 0,*garray = aij->garray,*nz,ec; 18687828ca2SBarry Smith PetscScalar v; 1872493cbb0SBarry Smith 1883a40ed3dSBarry Smith PetscFunctionBegin; 1892493cbb0SBarry Smith /* free stuff related to matrix-vec multiply */ 190b0a32e0cSBarry Smith ierr = VecGetSize(aij->lvec,&ec);CHKERRQ(ierr); /* needed for PetscLogObjectMemory below */ 1915e1f6667SBarry Smith ierr = VecDestroy(&aij->lvec);CHKERRQ(ierr); 1925e1f6667SBarry Smith ierr = VecScatterDestroy(&aij->Mvctx);CHKERRQ(ierr); 193464493b3SBarry Smith if (aij->colmap) { 194aa482453SBarry Smith #if defined(PETSC_USE_CTABLE) 1956bc0bbbfSBarry Smith ierr = PetscTableDestroy(&aij->colmap);CHKERRQ(ierr); 1962066d6f7SSatish Balay #else 197606d414cSSatish Balay ierr = PetscFree(aij->colmap);CHKERRQ(ierr); 198d0f46423SBarry Smith ierr = PetscLogObjectMemory(A,-aij->B->cmap->n*sizeof(PetscInt));CHKERRQ(ierr); 1992066d6f7SSatish Balay #endif 200464493b3SBarry Smith } 2012493cbb0SBarry Smith 2022493cbb0SBarry Smith /* make sure that B is assembled so we can access its values */ 2036d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 204fe2f2677SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2052493cbb0SBarry Smith 2062493cbb0SBarry Smith /* invent new B and copy stuff over */ 207b1d57f15SBarry Smith ierr = PetscMalloc((m+1)*sizeof(PetscInt),&nz);CHKERRQ(ierr); 20848b35521SBarry Smith for (i=0; i<m; i++) { 20948b35521SBarry Smith nz[i] = Baij->i[i+1] - Baij->i[i]; 21048b35521SBarry Smith } 211f69a0ea3SMatthew Knepley ierr = MatCreate(PETSC_COMM_SELF,&Bnew);CHKERRQ(ierr); 212f69a0ea3SMatthew Knepley ierr = MatSetSizes(Bnew,m,n,m,n);CHKERRQ(ierr); 213a2f3521dSMark F. Adams ierr = MatSetBlockSizes(Bnew,A->rmap->bs,A->cmap->bs);CHKERRQ(ierr); 2147adad957SLisandro Dalcin ierr = MatSetType(Bnew,((PetscObject)B)->type_name);CHKERRQ(ierr); 215f204ca49SKris Buschelman ierr = MatSeqAIJSetPreallocation(Bnew,0,nz);CHKERRQ(ierr); 2162205254eSKarl Rupp 2172576faa2SJed Brown ((Mat_SeqAIJ*)Bnew->data)->nonew = Baij->nonew; /* Inherit insertion error options. */ 2182205254eSKarl Rupp 219606d414cSSatish Balay ierr = PetscFree(nz);CHKERRQ(ierr); 2202493cbb0SBarry Smith for (i=0; i<m; i++) { 221bfec09a0SHong Zhang for (j=Baij->i[i]; j<Baij->i[i+1]; j++) { 222bfec09a0SHong Zhang col = garray[Baij->j[ct]]; 2232493cbb0SBarry Smith v = Baij->a[ct++]; 22483271157SBarry Smith ierr = MatSetValues(Bnew,1,&i,1,&col,&v,B->insertmode);CHKERRQ(ierr); 2252493cbb0SBarry Smith } 2262493cbb0SBarry Smith } 227606d414cSSatish Balay ierr = PetscFree(aij->garray);CHKERRQ(ierr); 22852e6d16bSBarry Smith ierr = PetscLogObjectMemory(A,-ec*sizeof(PetscInt));CHKERRQ(ierr); 2296bf464f9SBarry Smith ierr = MatDestroy(&B);CHKERRQ(ierr); 23052e6d16bSBarry Smith ierr = PetscLogObjectParent(A,Bnew);CHKERRQ(ierr); 2312205254eSKarl Rupp 2322493cbb0SBarry Smith aij->B = Bnew; 233227d817aSBarry Smith A->was_assembled = PETSC_FALSE; 2343a40ed3dSBarry Smith PetscFunctionReturn(0); 2352493cbb0SBarry Smith } 2362493cbb0SBarry Smith 2372cd6534aSBarry Smith /* ugly stuff added for Glenn someday we should fix this up */ 2382cd6534aSBarry Smith 2392205254eSKarl Rupp static PetscInt *auglyrmapd = 0,*auglyrmapo = 0; /* mapping from the local ordering to the "diagonal" and "off-diagonal" parts of the local matrix */ 2402cd6534aSBarry Smith static Vec auglydd = 0,auglyoo = 0; /* work vectors used to scale the two parts of the local matrix */ 2412cd6534aSBarry Smith 2422cd6534aSBarry Smith 2432cd6534aSBarry Smith #undef __FUNCT__ 2442cd6534aSBarry Smith #define __FUNCT__ "MatMPIAIJDiagonalScaleLocalSetUp" 245dfbe8321SBarry Smith PetscErrorCode MatMPIAIJDiagonalScaleLocalSetUp(Mat inA,Vec scale) 2462cd6534aSBarry Smith { 2472cd6534aSBarry Smith Mat_MPIAIJ *ina = (Mat_MPIAIJ*) inA->data; /*access private part of matrix */ 248dfbe8321SBarry Smith PetscErrorCode ierr; 249b1d57f15SBarry Smith PetscInt i,n,nt,cstart,cend,no,*garray = ina->garray,*lindices; 250b1d57f15SBarry Smith PetscInt *r_rmapd,*r_rmapo; 2512cd6534aSBarry Smith 2522cd6534aSBarry Smith PetscFunctionBegin; 2532cd6534aSBarry Smith ierr = MatGetOwnershipRange(inA,&cstart,&cend);CHKERRQ(ierr); 2540298fd71SBarry Smith ierr = MatGetSize(ina->A,NULL,&n);CHKERRQ(ierr); 255992144d0SBarry Smith ierr = PetscMalloc((inA->rmap->mapping->n+1)*sizeof(PetscInt),&r_rmapd);CHKERRQ(ierr); 256992144d0SBarry Smith ierr = PetscMemzero(r_rmapd,inA->rmap->mapping->n*sizeof(PetscInt));CHKERRQ(ierr); 2572cd6534aSBarry Smith nt = 0; 258992144d0SBarry Smith for (i=0; i<inA->rmap->mapping->n; i++) { 259992144d0SBarry Smith if (inA->rmap->mapping->indices[i] >= cstart && inA->rmap->mapping->indices[i] < cend) { 2602cd6534aSBarry Smith nt++; 261992144d0SBarry Smith r_rmapd[i] = inA->rmap->mapping->indices[i] + 1; 2622cd6534aSBarry Smith } 2632cd6534aSBarry Smith } 264e32f2f54SBarry Smith if (nt != n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hmm nt %D n %D",nt,n); 265b1d57f15SBarry Smith ierr = PetscMalloc((n+1)*sizeof(PetscInt),&auglyrmapd);CHKERRQ(ierr); 266992144d0SBarry Smith for (i=0; i<inA->rmap->mapping->n; i++) { 2672cd6534aSBarry Smith if (r_rmapd[i]) { 2682cd6534aSBarry Smith auglyrmapd[(r_rmapd[i]-1)-cstart] = i; 2692cd6534aSBarry Smith } 2702cd6534aSBarry Smith } 2712cd6534aSBarry Smith ierr = PetscFree(r_rmapd);CHKERRQ(ierr); 2722cd6534aSBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,n,&auglydd);CHKERRQ(ierr); 2732cd6534aSBarry Smith 274d0f46423SBarry Smith ierr = PetscMalloc((inA->cmap->N+1)*sizeof(PetscInt),&lindices);CHKERRQ(ierr); 275d0f46423SBarry Smith ierr = PetscMemzero(lindices,inA->cmap->N*sizeof(PetscInt));CHKERRQ(ierr); 276d0f46423SBarry Smith for (i=0; i<ina->B->cmap->n; i++) { 2772cd6534aSBarry Smith lindices[garray[i]] = i+1; 2782cd6534aSBarry Smith } 279992144d0SBarry Smith no = inA->rmap->mapping->n - nt; 280992144d0SBarry Smith ierr = PetscMalloc((inA->rmap->mapping->n+1)*sizeof(PetscInt),&r_rmapo);CHKERRQ(ierr); 281992144d0SBarry Smith ierr = PetscMemzero(r_rmapo,inA->rmap->mapping->n*sizeof(PetscInt));CHKERRQ(ierr); 2822cd6534aSBarry Smith nt = 0; 283992144d0SBarry Smith for (i=0; i<inA->rmap->mapping->n; i++) { 284992144d0SBarry Smith if (lindices[inA->rmap->mapping->indices[i]]) { 2852cd6534aSBarry Smith nt++; 286992144d0SBarry Smith r_rmapo[i] = lindices[inA->rmap->mapping->indices[i]]; 2872cd6534aSBarry Smith } 2882cd6534aSBarry Smith } 289e32f2f54SBarry Smith if (nt > no) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hmm nt %D no %D",nt,n); 2902cd6534aSBarry Smith ierr = PetscFree(lindices);CHKERRQ(ierr); 291b1d57f15SBarry Smith ierr = PetscMalloc((nt+1)*sizeof(PetscInt),&auglyrmapo);CHKERRQ(ierr); 292992144d0SBarry Smith for (i=0; i<inA->rmap->mapping->n; i++) { 2932cd6534aSBarry Smith if (r_rmapo[i]) { 2942cd6534aSBarry Smith auglyrmapo[(r_rmapo[i]-1)] = i; 2952cd6534aSBarry Smith } 2962cd6534aSBarry Smith } 2972cd6534aSBarry Smith ierr = PetscFree(r_rmapo);CHKERRQ(ierr); 2982cd6534aSBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,nt,&auglyoo);CHKERRQ(ierr); 2992cd6534aSBarry Smith PetscFunctionReturn(0); 3002cd6534aSBarry Smith } 3012cd6534aSBarry Smith 3022cd6534aSBarry Smith #undef __FUNCT__ 3032cd6534aSBarry Smith #define __FUNCT__ "MatMPIAIJDiagonalScaleLocal" 304dfbe8321SBarry Smith PetscErrorCode MatMPIAIJDiagonalScaleLocal(Mat A,Vec scale) 3052cd6534aSBarry Smith { 30692b32695SKris Buschelman /* This routine should really be abandoned as it duplicates MatDiagonalScaleLocal */ 3074ac538c5SBarry Smith PetscErrorCode ierr; 30892b32695SKris Buschelman 30992b32695SKris Buschelman PetscFunctionBegin; 3104ac538c5SBarry Smith ierr = PetscTryMethod(A,"MatDiagonalScaleLocal_C",(Mat,Vec),(A,scale));CHKERRQ(ierr); 31192b32695SKris Buschelman PetscFunctionReturn(0); 31292b32695SKris Buschelman } 31392b32695SKris Buschelman 31492b32695SKris Buschelman EXTERN_C_BEGIN 31592b32695SKris Buschelman #undef __FUNCT__ 31692b32695SKris Buschelman #define __FUNCT__ "MatDiagonalScaleLocal_MPIAIJ" 3177087cfbeSBarry Smith PetscErrorCode MatDiagonalScaleLocal_MPIAIJ(Mat A,Vec scale) 31892b32695SKris Buschelman { 3192cd6534aSBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ*) A->data; /*access private part of matrix */ 320dfbe8321SBarry Smith PetscErrorCode ierr; 321b1d57f15SBarry Smith PetscInt n,i; 3222cd6534aSBarry Smith PetscScalar *d,*o,*s; 3232cd6534aSBarry Smith 3242cd6534aSBarry Smith PetscFunctionBegin; 3252cd6534aSBarry Smith if (!auglyrmapd) { 3262cd6534aSBarry Smith ierr = MatMPIAIJDiagonalScaleLocalSetUp(A,scale);CHKERRQ(ierr); 3272cd6534aSBarry Smith } 3282cd6534aSBarry Smith 3291ebc52fbSHong Zhang ierr = VecGetArray(scale,&s);CHKERRQ(ierr); 3302cd6534aSBarry Smith 3312cd6534aSBarry Smith ierr = VecGetLocalSize(auglydd,&n);CHKERRQ(ierr); 3321ebc52fbSHong Zhang ierr = VecGetArray(auglydd,&d);CHKERRQ(ierr); 3332cd6534aSBarry Smith for (i=0; i<n; i++) { 3342cd6534aSBarry Smith d[i] = s[auglyrmapd[i]]; /* copy "diagonal" (true local) portion of scale into dd vector */ 3352cd6534aSBarry Smith } 3361ebc52fbSHong Zhang ierr = VecRestoreArray(auglydd,&d);CHKERRQ(ierr); 3372cd6534aSBarry Smith /* column scale "diagonal" portion of local matrix */ 3380298fd71SBarry Smith ierr = MatDiagonalScale(a->A,NULL,auglydd);CHKERRQ(ierr); 3392cd6534aSBarry Smith 3402cd6534aSBarry Smith ierr = VecGetLocalSize(auglyoo,&n);CHKERRQ(ierr); 3411ebc52fbSHong Zhang ierr = VecGetArray(auglyoo,&o);CHKERRQ(ierr); 3422cd6534aSBarry Smith for (i=0; i<n; i++) { 3432cd6534aSBarry Smith o[i] = s[auglyrmapo[i]]; /* copy "off-diagonal" portion of scale into oo vector */ 3442cd6534aSBarry Smith } 3451ebc52fbSHong Zhang ierr = VecRestoreArray(scale,&s);CHKERRQ(ierr); 3461ebc52fbSHong Zhang ierr = VecRestoreArray(auglyoo,&o);CHKERRQ(ierr); 3472cd6534aSBarry Smith /* column scale "off-diagonal" portion of local matrix */ 3480298fd71SBarry Smith ierr = MatDiagonalScale(a->B,NULL,auglyoo);CHKERRQ(ierr); 3492cd6534aSBarry Smith PetscFunctionReturn(0); 3502cd6534aSBarry Smith } 35192b32695SKris Buschelman EXTERN_C_END 3522cd6534aSBarry Smith 35348b35521SBarry Smith 354