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> 68c79f6d3SBarry Smith 74a2ae208SSatish Balay #undef __FUNCT__ 84a2ae208SSatish Balay #define __FUNCT__ "MatSetUpMultiply_MPIAIJ" 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; 17ace3abfcSBarry Smith PetscBool useblockis; 18aa482453SBarry Smith #if defined(PETSC_USE_CTABLE) 190f5bd95cSBarry Smith PetscTable gid1_lid1; 200f5bd95cSBarry Smith PetscTablePosition tpos; 21b1d57f15SBarry Smith PetscInt gid,lid; 226f531f54SSatish Balay #else 23d0f46423SBarry Smith PetscInt N = mat->cmap->N,*indices; 242066d6f7SSatish Balay #endif 252066d6f7SSatish Balay 263a40ed3dSBarry Smith PetscFunctionBegin; 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 */ 41b1d57f15SBarry Smith ierr = PetscMalloc((ec+1)*sizeof(PetscInt),&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; 64*2205254eSKarl Rupp 6526283091SBarry Smith ierr = PetscLayoutSetUp((aij->B->cmap));CHKERRQ(ierr); 666bc0bbbfSBarry Smith ierr = PetscTableDestroy(&gid1_lid1);CHKERRQ(ierr); 672066d6f7SSatish Balay #else 6811285404SBarry Smith /* Make an array as long as the number of columns */ 691eb62cbbSBarry Smith /* mark those columns that are in aij->B */ 70b1d57f15SBarry Smith ierr = PetscMalloc((N+1)*sizeof(PetscInt),&indices);CHKERRQ(ierr); 71b1d57f15SBarry Smith ierr = PetscMemzero(indices,N*sizeof(PetscInt));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 */ 80b1d57f15SBarry Smith ierr = PetscMalloc((ec+1)*sizeof(PetscInt),&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; 98*2205254eSKarl Rupp 9926283091SBarry Smith ierr = PetscLayoutSetUp((aij->B->cmap));CHKERRQ(ierr); 100606d414cSSatish Balay ierr = PetscFree(indices);CHKERRQ(ierr); 1012066d6f7SSatish Balay #endif 1021eb62cbbSBarry Smith /* create local vector that is used to scatter into */ 103029af93fSBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,ec,&aij->lvec);CHKERRQ(ierr); 1041eb62cbbSBarry Smith 105d6dfbf8fSBarry Smith /* create two temporary Index sets for build scatter gather */ 10603919abeSBarry Smith /* check for the special case where blocks are communicated for faster VecScatterXXX */ 107a2f3521dSMark F. Adams useblockis = PETSC_FALSE; 108c5bfad50SMark F. Adams if (mat->cmap->bs > 1) { 109c5bfad50SMark F. Adams PetscInt bs = mat->cmap->bs,ibs,ga; 11003919abeSBarry Smith if (!(ec % bs)) { 1113ba34427SBarry Smith useblockis = PETSC_TRUE; 11203919abeSBarry Smith for (i=0; i<ec/bs; i++) { 11303919abeSBarry Smith if ((ga = garray[ibs = i*bs]) % bs) { 11403919abeSBarry Smith useblockis = PETSC_FALSE; 11503919abeSBarry Smith break; 11603919abeSBarry Smith } 11703919abeSBarry Smith for (j=1; j<bs; j++) { 11803919abeSBarry Smith if (garray[ibs+j] != ga+j) { 11903919abeSBarry Smith useblockis = PETSC_FALSE; 12003919abeSBarry Smith break; 12103919abeSBarry Smith } 12203919abeSBarry Smith } 12303919abeSBarry Smith if (!useblockis) break; 12403919abeSBarry Smith } 12503919abeSBarry Smith } 12603919abeSBarry Smith } 127c5bfad50SMark F. Adams #if defined(PETSC_USE_DEBUG) 128c5bfad50SMark F. Adams i = (PetscInt)useblockis; 129c5bfad50SMark F. Adams ierr = MPI_Allreduce(&i,&j,1,MPIU_INT,MPI_MIN,((PetscObject)mat)->comm);CHKERRQ(ierr); 130c5bfad50SMark F. Adams if (j!=i) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Use of blocked not consistant (I am usning blocked)"); 131c5bfad50SMark F. Adams #endif 132c5bfad50SMark F. Adams 13303919abeSBarry Smith if (useblockis) { 134c5bfad50SMark F. Adams PetscInt *ga,bs = mat->cmap->bs,iec = ec/bs; 135c5bfad50SMark F. Adams if (ec%bs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"ec=%D bs=%D",ec,bs); 13622d28d08SBarry Smith ierr = PetscInfo(mat,"Using block index set to define scatter\n");CHKERRQ(ierr); 137c5bfad50SMark F. Adams ierr = PetscMalloc(iec*sizeof(PetscInt),&ga);CHKERRQ(ierr); 138e82e9f6bSBarry Smith for (i=0; i<iec; i++) ga[i] = garray[i*bs]/bs; 139deff0451SBarry Smith ierr = ISCreateBlock(((PetscObject)mat)->comm,bs,iec,ga,PETSC_OWN_POINTER,&from);CHKERRQ(ierr); 14003919abeSBarry Smith } else { 14170b3c8c7SBarry Smith ierr = ISCreateGeneral(((PetscObject)mat)->comm,ec,garray,PETSC_COPY_VALUES,&from);CHKERRQ(ierr); 14203919abeSBarry Smith } 143c5bfad50SMark F. Adams 144029af93fSBarry Smith ierr = ISCreateStride(PETSC_COMM_SELF,ec,0,1,&to);CHKERRQ(ierr); 1451eb62cbbSBarry Smith 1461eb62cbbSBarry Smith /* create temporary global vector to generate scatter context */ 147b5eb4454SBarry Smith /* This does not allocate the array's memory so is efficient */ 148778a2246SBarry Smith ierr = VecCreateMPIWithArray(((PetscObject)mat)->comm,1,mat->cmap->n,mat->cmap->N,PETSC_NULL,&gvec);CHKERRQ(ierr); 1491eb62cbbSBarry Smith 1502d336d48SLois Curfman McInnes /* generate the scatter context */ 15108480c60SBarry Smith ierr = VecScatterCreate(gvec,from,aij->lvec,to,&aij->Mvctx);CHKERRQ(ierr); 15252e6d16bSBarry Smith ierr = PetscLogObjectParent(mat,aij->Mvctx);CHKERRQ(ierr); 15352e6d16bSBarry Smith ierr = PetscLogObjectParent(mat,aij->lvec);CHKERRQ(ierr); 15452e6d16bSBarry Smith ierr = PetscLogObjectParent(mat,from);CHKERRQ(ierr); 15552e6d16bSBarry Smith ierr = PetscLogObjectParent(mat,to);CHKERRQ(ierr); 156*2205254eSKarl Rupp 1579e25ed09SBarry Smith aij->garray = garray; 158*2205254eSKarl Rupp 15952e6d16bSBarry Smith ierr = PetscLogObjectMemory(mat,(ec+1)*sizeof(PetscInt));CHKERRQ(ierr); 1606bf464f9SBarry Smith ierr = ISDestroy(&from);CHKERRQ(ierr); 1616bf464f9SBarry Smith ierr = ISDestroy(&to);CHKERRQ(ierr); 1626bf464f9SBarry Smith ierr = VecDestroy(&gvec);CHKERRQ(ierr); 1633a40ed3dSBarry Smith PetscFunctionReturn(0); 1648c79f6d3SBarry Smith } 1659e25ed09SBarry Smith 1669e25ed09SBarry Smith 1674a2ae208SSatish Balay #undef __FUNCT__ 168ab9863d7SBarry Smith #define __FUNCT__ "MatDisAssemble_MPIAIJ" 1692493cbb0SBarry Smith /* 1702493cbb0SBarry Smith Takes the local part of an already assembled MPIAIJ matrix 1712493cbb0SBarry Smith and disassembles it. This is to allow new nonzeros into the matrix 1722493cbb0SBarry Smith that require more communication in the matrix vector multiply. 1732493cbb0SBarry Smith Thus certain data-structures must be rebuilt. 1742493cbb0SBarry Smith 1752493cbb0SBarry Smith Kind of slow! But that's what application programmers get when 1762493cbb0SBarry Smith they are sloppy. 1772493cbb0SBarry Smith */ 178ab9863d7SBarry Smith PetscErrorCode MatDisAssemble_MPIAIJ(Mat A) 1792493cbb0SBarry Smith { 1802493cbb0SBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ*)A->data; 1812493cbb0SBarry Smith Mat B = aij->B,Bnew; 182ec8511deSBarry Smith Mat_SeqAIJ *Baij = (Mat_SeqAIJ*)B->data; 183dfbe8321SBarry Smith PetscErrorCode ierr; 184d0f46423SBarry Smith PetscInt i,j,m = B->rmap->n,n = A->cmap->N,col,ct = 0,*garray = aij->garray,*nz,ec; 18587828ca2SBarry Smith PetscScalar v; 1862493cbb0SBarry Smith 1873a40ed3dSBarry Smith PetscFunctionBegin; 1882493cbb0SBarry Smith /* free stuff related to matrix-vec multiply */ 189b0a32e0cSBarry Smith ierr = VecGetSize(aij->lvec,&ec);CHKERRQ(ierr); /* needed for PetscLogObjectMemory below */ 1905e1f6667SBarry Smith ierr = VecDestroy(&aij->lvec);CHKERRQ(ierr); 1915e1f6667SBarry Smith ierr = VecScatterDestroy(&aij->Mvctx);CHKERRQ(ierr); 192464493b3SBarry Smith if (aij->colmap) { 193aa482453SBarry Smith #if defined(PETSC_USE_CTABLE) 1946bc0bbbfSBarry Smith ierr = PetscTableDestroy(&aij->colmap);CHKERRQ(ierr); 1952066d6f7SSatish Balay #else 196606d414cSSatish Balay ierr = PetscFree(aij->colmap);CHKERRQ(ierr); 197d0f46423SBarry Smith ierr = PetscLogObjectMemory(A,-aij->B->cmap->n*sizeof(PetscInt));CHKERRQ(ierr); 1982066d6f7SSatish Balay #endif 199464493b3SBarry Smith } 2002493cbb0SBarry Smith 2012493cbb0SBarry Smith /* make sure that B is assembled so we can access its values */ 2026d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 203fe2f2677SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2042493cbb0SBarry Smith 2052493cbb0SBarry Smith /* invent new B and copy stuff over */ 206b1d57f15SBarry Smith ierr = PetscMalloc((m+1)*sizeof(PetscInt),&nz);CHKERRQ(ierr); 20748b35521SBarry Smith for (i=0; i<m; i++) { 20848b35521SBarry Smith nz[i] = Baij->i[i+1] - Baij->i[i]; 20948b35521SBarry Smith } 210f69a0ea3SMatthew Knepley ierr = MatCreate(PETSC_COMM_SELF,&Bnew);CHKERRQ(ierr); 211f69a0ea3SMatthew Knepley ierr = MatSetSizes(Bnew,m,n,m,n);CHKERRQ(ierr); 212a2f3521dSMark F. Adams ierr = MatSetBlockSizes(Bnew,A->rmap->bs,A->cmap->bs);CHKERRQ(ierr); 2137adad957SLisandro Dalcin ierr = MatSetType(Bnew,((PetscObject)B)->type_name);CHKERRQ(ierr); 214f204ca49SKris Buschelman ierr = MatSeqAIJSetPreallocation(Bnew,0,nz);CHKERRQ(ierr); 215*2205254eSKarl Rupp 2162576faa2SJed Brown ((Mat_SeqAIJ*)Bnew->data)->nonew = Baij->nonew; /* Inherit insertion error options. */ 217*2205254eSKarl Rupp 218606d414cSSatish Balay ierr = PetscFree(nz);CHKERRQ(ierr); 2192493cbb0SBarry Smith for (i=0; i<m; i++) { 220bfec09a0SHong Zhang for (j=Baij->i[i]; j<Baij->i[i+1]; j++) { 221bfec09a0SHong Zhang col = garray[Baij->j[ct]]; 2222493cbb0SBarry Smith v = Baij->a[ct++]; 22383271157SBarry Smith ierr = MatSetValues(Bnew,1,&i,1,&col,&v,B->insertmode);CHKERRQ(ierr); 2242493cbb0SBarry Smith } 2252493cbb0SBarry Smith } 226606d414cSSatish Balay ierr = PetscFree(aij->garray);CHKERRQ(ierr); 22752e6d16bSBarry Smith ierr = PetscLogObjectMemory(A,-ec*sizeof(PetscInt));CHKERRQ(ierr); 2286bf464f9SBarry Smith ierr = MatDestroy(&B);CHKERRQ(ierr); 22952e6d16bSBarry Smith ierr = PetscLogObjectParent(A,Bnew);CHKERRQ(ierr); 230*2205254eSKarl Rupp 2312493cbb0SBarry Smith aij->B = Bnew; 232227d817aSBarry Smith A->was_assembled = PETSC_FALSE; 2333a40ed3dSBarry Smith PetscFunctionReturn(0); 2342493cbb0SBarry Smith } 2352493cbb0SBarry Smith 2362cd6534aSBarry Smith /* ugly stuff added for Glenn someday we should fix this up */ 2372cd6534aSBarry Smith 238*2205254eSKarl Rupp static PetscInt *auglyrmapd = 0,*auglyrmapo = 0; /* mapping from the local ordering to the "diagonal" and "off-diagonal" parts of the local matrix */ 2392cd6534aSBarry Smith static Vec auglydd = 0,auglyoo = 0; /* work vectors used to scale the two parts of the local matrix */ 2402cd6534aSBarry Smith 2412cd6534aSBarry Smith 2422cd6534aSBarry Smith #undef __FUNCT__ 2432cd6534aSBarry Smith #define __FUNCT__ "MatMPIAIJDiagonalScaleLocalSetUp" 244dfbe8321SBarry Smith PetscErrorCode MatMPIAIJDiagonalScaleLocalSetUp(Mat inA,Vec scale) 2452cd6534aSBarry Smith { 2462cd6534aSBarry Smith Mat_MPIAIJ *ina = (Mat_MPIAIJ*) inA->data; /*access private part of matrix */ 247dfbe8321SBarry Smith PetscErrorCode ierr; 248b1d57f15SBarry Smith PetscInt i,n,nt,cstart,cend,no,*garray = ina->garray,*lindices; 249b1d57f15SBarry Smith PetscInt *r_rmapd,*r_rmapo; 2502cd6534aSBarry Smith 2512cd6534aSBarry Smith PetscFunctionBegin; 2522cd6534aSBarry Smith ierr = MatGetOwnershipRange(inA,&cstart,&cend);CHKERRQ(ierr); 2532cd6534aSBarry Smith ierr = MatGetSize(ina->A,PETSC_NULL,&n);CHKERRQ(ierr); 254992144d0SBarry Smith ierr = PetscMalloc((inA->rmap->mapping->n+1)*sizeof(PetscInt),&r_rmapd);CHKERRQ(ierr); 255992144d0SBarry Smith ierr = PetscMemzero(r_rmapd,inA->rmap->mapping->n*sizeof(PetscInt));CHKERRQ(ierr); 2562cd6534aSBarry Smith nt = 0; 257992144d0SBarry Smith for (i=0; i<inA->rmap->mapping->n; i++) { 258992144d0SBarry Smith if (inA->rmap->mapping->indices[i] >= cstart && inA->rmap->mapping->indices[i] < cend) { 2592cd6534aSBarry Smith nt++; 260992144d0SBarry Smith r_rmapd[i] = inA->rmap->mapping->indices[i] + 1; 2612cd6534aSBarry Smith } 2622cd6534aSBarry Smith } 263e32f2f54SBarry Smith if (nt != n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hmm nt %D n %D",nt,n); 264b1d57f15SBarry Smith ierr = PetscMalloc((n+1)*sizeof(PetscInt),&auglyrmapd);CHKERRQ(ierr); 265992144d0SBarry Smith for (i=0; i<inA->rmap->mapping->n; i++) { 2662cd6534aSBarry Smith if (r_rmapd[i]) { 2672cd6534aSBarry Smith auglyrmapd[(r_rmapd[i]-1)-cstart] = i; 2682cd6534aSBarry Smith } 2692cd6534aSBarry Smith } 2702cd6534aSBarry Smith ierr = PetscFree(r_rmapd);CHKERRQ(ierr); 2712cd6534aSBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,n,&auglydd);CHKERRQ(ierr); 2722cd6534aSBarry Smith 273d0f46423SBarry Smith ierr = PetscMalloc((inA->cmap->N+1)*sizeof(PetscInt),&lindices);CHKERRQ(ierr); 274d0f46423SBarry Smith ierr = PetscMemzero(lindices,inA->cmap->N*sizeof(PetscInt));CHKERRQ(ierr); 275d0f46423SBarry Smith for (i=0; i<ina->B->cmap->n; i++) { 2762cd6534aSBarry Smith lindices[garray[i]] = i+1; 2772cd6534aSBarry Smith } 278992144d0SBarry Smith no = inA->rmap->mapping->n - nt; 279992144d0SBarry Smith ierr = PetscMalloc((inA->rmap->mapping->n+1)*sizeof(PetscInt),&r_rmapo);CHKERRQ(ierr); 280992144d0SBarry Smith ierr = PetscMemzero(r_rmapo,inA->rmap->mapping->n*sizeof(PetscInt));CHKERRQ(ierr); 2812cd6534aSBarry Smith nt = 0; 282992144d0SBarry Smith for (i=0; i<inA->rmap->mapping->n; i++) { 283992144d0SBarry Smith if (lindices[inA->rmap->mapping->indices[i]]) { 2842cd6534aSBarry Smith nt++; 285992144d0SBarry Smith r_rmapo[i] = lindices[inA->rmap->mapping->indices[i]]; 2862cd6534aSBarry Smith } 2872cd6534aSBarry Smith } 288e32f2f54SBarry Smith if (nt > no) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hmm nt %D no %D",nt,n); 2892cd6534aSBarry Smith ierr = PetscFree(lindices);CHKERRQ(ierr); 290b1d57f15SBarry Smith ierr = PetscMalloc((nt+1)*sizeof(PetscInt),&auglyrmapo);CHKERRQ(ierr); 291992144d0SBarry Smith for (i=0; i<inA->rmap->mapping->n; i++) { 2922cd6534aSBarry Smith if (r_rmapo[i]) { 2932cd6534aSBarry Smith auglyrmapo[(r_rmapo[i]-1)] = i; 2942cd6534aSBarry Smith } 2952cd6534aSBarry Smith } 2962cd6534aSBarry Smith ierr = PetscFree(r_rmapo);CHKERRQ(ierr); 2972cd6534aSBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,nt,&auglyoo);CHKERRQ(ierr); 2982cd6534aSBarry Smith PetscFunctionReturn(0); 2992cd6534aSBarry Smith } 3002cd6534aSBarry Smith 3012cd6534aSBarry Smith #undef __FUNCT__ 3022cd6534aSBarry Smith #define __FUNCT__ "MatMPIAIJDiagonalScaleLocal" 303dfbe8321SBarry Smith PetscErrorCode MatMPIAIJDiagonalScaleLocal(Mat A,Vec scale) 3042cd6534aSBarry Smith { 30592b32695SKris Buschelman /* This routine should really be abandoned as it duplicates MatDiagonalScaleLocal */ 3064ac538c5SBarry Smith PetscErrorCode ierr; 30792b32695SKris Buschelman 30892b32695SKris Buschelman PetscFunctionBegin; 3094ac538c5SBarry Smith ierr = PetscTryMethod(A,"MatDiagonalScaleLocal_C",(Mat,Vec),(A,scale));CHKERRQ(ierr); 31092b32695SKris Buschelman PetscFunctionReturn(0); 31192b32695SKris Buschelman } 31292b32695SKris Buschelman 31392b32695SKris Buschelman EXTERN_C_BEGIN 31492b32695SKris Buschelman #undef __FUNCT__ 31592b32695SKris Buschelman #define __FUNCT__ "MatDiagonalScaleLocal_MPIAIJ" 3167087cfbeSBarry Smith PetscErrorCode MatDiagonalScaleLocal_MPIAIJ(Mat A,Vec scale) 31792b32695SKris Buschelman { 3182cd6534aSBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ*) A->data; /*access private part of matrix */ 319dfbe8321SBarry Smith PetscErrorCode ierr; 320b1d57f15SBarry Smith PetscInt n,i; 3212cd6534aSBarry Smith PetscScalar *d,*o,*s; 3222cd6534aSBarry Smith 3232cd6534aSBarry Smith PetscFunctionBegin; 3242cd6534aSBarry Smith if (!auglyrmapd) { 3252cd6534aSBarry Smith ierr = MatMPIAIJDiagonalScaleLocalSetUp(A,scale);CHKERRQ(ierr); 3262cd6534aSBarry Smith } 3272cd6534aSBarry Smith 3281ebc52fbSHong Zhang ierr = VecGetArray(scale,&s);CHKERRQ(ierr); 3292cd6534aSBarry Smith 3302cd6534aSBarry Smith ierr = VecGetLocalSize(auglydd,&n);CHKERRQ(ierr); 3311ebc52fbSHong Zhang ierr = VecGetArray(auglydd,&d);CHKERRQ(ierr); 3322cd6534aSBarry Smith for (i=0; i<n; i++) { 3332cd6534aSBarry Smith d[i] = s[auglyrmapd[i]]; /* copy "diagonal" (true local) portion of scale into dd vector */ 3342cd6534aSBarry Smith } 3351ebc52fbSHong Zhang ierr = VecRestoreArray(auglydd,&d);CHKERRQ(ierr); 3362cd6534aSBarry Smith /* column scale "diagonal" portion of local matrix */ 3372cd6534aSBarry Smith ierr = MatDiagonalScale(a->A,PETSC_NULL,auglydd);CHKERRQ(ierr); 3382cd6534aSBarry Smith 3392cd6534aSBarry Smith ierr = VecGetLocalSize(auglyoo,&n);CHKERRQ(ierr); 3401ebc52fbSHong Zhang ierr = VecGetArray(auglyoo,&o);CHKERRQ(ierr); 3412cd6534aSBarry Smith for (i=0; i<n; i++) { 3422cd6534aSBarry Smith o[i] = s[auglyrmapo[i]]; /* copy "off-diagonal" portion of scale into oo vector */ 3432cd6534aSBarry Smith } 3441ebc52fbSHong Zhang ierr = VecRestoreArray(scale,&s);CHKERRQ(ierr); 3451ebc52fbSHong Zhang ierr = VecRestoreArray(auglyoo,&o);CHKERRQ(ierr); 3462cd6534aSBarry Smith /* column scale "off-diagonal" portion of local matrix */ 3472cd6534aSBarry Smith ierr = MatDiagonalScale(a->B,PETSC_NULL,auglyoo);CHKERRQ(ierr); 3482cd6534aSBarry Smith PetscFunctionReturn(0); 3492cd6534aSBarry Smith } 35092b32695SKris Buschelman EXTERN_C_END 3512cd6534aSBarry Smith 35248b35521SBarry Smith 353