1be1d678aSKris Buschelman #define PETSCMAT_DLL 2be1d678aSKris Buschelman 38c79f6d3SBarry Smith /* 48c79f6d3SBarry Smith Support for the parallel AIJ matrix vector multiply 58c79f6d3SBarry Smith */ 670f55243SBarry Smith #include "src/mat/impls/aij/mpi/mpiaij.h" 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; 1803919abeSBarry Smith PetscTruth 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 24b1d57f15SBarry Smith PetscInt N = mat->N,*indices; 252066d6f7SSatish Balay #endif 262066d6f7SSatish Balay 273a40ed3dSBarry Smith PetscFunctionBegin; 282066d6f7SSatish Balay 29aa482453SBarry Smith #if defined (PETSC_USE_CTABLE) 30dc231df0SBarry Smith /* use a table - Mark Adams */ 31273d9f13SBarry Smith ierr = PetscTableCreate(aij->B->m,&gid1_lid1);CHKERRQ(ierr); 32273d9f13SBarry Smith for (i=0; i<aij->B->m; i++) { 332066d6f7SSatish Balay for (j=0; j<B->ilen[i]; j++) { 34b1d57f15SBarry Smith PetscInt data,gid1 = aj[B->i[i] + j] + 1; 350f5bd95cSBarry Smith ierr = PetscTableFind(gid1_lid1,gid1,&data);CHKERRQ(ierr); 36fa46199cSSatish Balay if (!data) { 372066d6f7SSatish Balay /* one based table */ 380f5bd95cSBarry Smith ierr = PetscTableAdd(gid1_lid1,gid1,++ec);CHKERRQ(ierr); 392066d6f7SSatish Balay } 402066d6f7SSatish Balay } 412066d6f7SSatish Balay } 422066d6f7SSatish Balay /* form array of columns we need */ 43b1d57f15SBarry Smith ierr = PetscMalloc((ec+1)*sizeof(PetscInt),&garray);CHKERRQ(ierr); 440f5bd95cSBarry Smith ierr = PetscTableGetHeadPosition(gid1_lid1,&tpos);CHKERRQ(ierr); 452066d6f7SSatish Balay while (tpos) { 460f5bd95cSBarry Smith ierr = PetscTableGetNext(gid1_lid1,&tpos,&gid,&lid);CHKERRQ(ierr); 47b0a32e0cSBarry Smith gid--; 48b0a32e0cSBarry Smith lid--; 492066d6f7SSatish Balay garray[lid] = gid; 502066d6f7SSatish Balay } 510064e2bbSSatish Balay ierr = PetscSortInt(ec,garray);CHKERRQ(ierr); /* sort, and rebuild */ 520f5bd95cSBarry Smith ierr = PetscTableRemoveAll(gid1_lid1);CHKERRQ(ierr); 532066d6f7SSatish Balay for (i=0; i<ec; i++) { 540f5bd95cSBarry Smith ierr = PetscTableAdd(gid1_lid1,garray[i]+1,i+1);CHKERRQ(ierr); 552066d6f7SSatish Balay } 562066d6f7SSatish Balay /* compact out the extra columns in B */ 57273d9f13SBarry Smith for (i=0; i<aij->B->m; i++) { 582066d6f7SSatish Balay for (j=0; j<B->ilen[i]; j++) { 59b1d57f15SBarry Smith PetscInt gid1 = aj[B->i[i] + j] + 1; 600f5bd95cSBarry Smith ierr = PetscTableFind(gid1_lid1,gid1,&lid);CHKERRQ(ierr); 61fa46199cSSatish Balay lid --; 62b3fb0a6cSMatthew Knepley aj[B->i[i] + j] = lid; 632066d6f7SSatish Balay } 642066d6f7SSatish Balay } 65273d9f13SBarry Smith aij->B->n = aij->B->N = ec; 660f5bd95cSBarry Smith ierr = PetscTableDelete(gid1_lid1);CHKERRQ(ierr); 672066d6f7SSatish Balay /* Mark Adams */ 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); 73273d9f13SBarry Smith for (i=0; i<aij->B->m; 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 */ 93273d9f13SBarry Smith for (i=0; i<aij->B->m; 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 } 98273d9f13SBarry Smith aij->B->n = aij->B->N = ec; 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 */ 10503919abeSBarry Smith /* check for the special case where blocks are communicated for faster VecScatterXXX */ 10603919abeSBarry Smith useblockis = PETSC_TRUE; 10703919abeSBarry Smith if (mat->bs > 1) { 10803919abeSBarry Smith PetscInt bs = mat->bs,ibs,ga; 10903919abeSBarry Smith if (!(ec % bs)) { 11003919abeSBarry Smith for (i=0; i<ec/bs; i++) { 11103919abeSBarry Smith if ((ga = garray[ibs = i*bs]) % bs) { 11203919abeSBarry Smith useblockis = PETSC_FALSE; 11303919abeSBarry Smith break; 11403919abeSBarry Smith } 11503919abeSBarry Smith for (j=1; j<bs; j++) { 11603919abeSBarry Smith if (garray[ibs+j] != ga+j) { 11703919abeSBarry Smith useblockis = PETSC_FALSE; 11803919abeSBarry Smith break; 11903919abeSBarry Smith } 12003919abeSBarry Smith } 12103919abeSBarry Smith if (!useblockis) break; 12203919abeSBarry Smith } 12303919abeSBarry Smith } 12403919abeSBarry Smith } 12503919abeSBarry Smith if (useblockis) { 12603919abeSBarry Smith PetscInt *ga,bs = mat->bs,iec = ec/bs; 127*ae15b995SBarry Smith ierr = PetscInfo(mat,"Using block index set to define scatter\n"); 12803919abeSBarry Smith ierr = PetscMalloc((ec/mat->bs)*sizeof(PetscInt),&ga);CHKERRQ(ierr); 12903919abeSBarry Smith for (i=0; i<iec; i++) ga[i] = garray[i*bs]; 13003919abeSBarry Smith ierr = ISCreateBlock(mat->comm,bs,iec,ga,&from);CHKERRQ(ierr); 13103919abeSBarry Smith ierr = PetscFree(ga);CHKERRQ(ierr); 13203919abeSBarry Smith } else { 133b9b97703SBarry Smith ierr = ISCreateGeneral(mat->comm,ec,garray,&from);CHKERRQ(ierr); 13403919abeSBarry Smith } 135029af93fSBarry Smith ierr = ISCreateStride(PETSC_COMM_SELF,ec,0,1,&to);CHKERRQ(ierr); 1361eb62cbbSBarry Smith 1371eb62cbbSBarry Smith /* create temporary global vector to generate scatter context */ 1381eb62cbbSBarry Smith /* this is inefficient, but otherwise we must do either 1391eb62cbbSBarry Smith 1) save garray until the first actual scatter when the vector is known or 1401eb62cbbSBarry Smith 2) have another way of generating a scatter context without a vector.*/ 141273d9f13SBarry Smith ierr = VecCreateMPI(mat->comm,mat->n,mat->N,&gvec);CHKERRQ(ierr); 1421eb62cbbSBarry Smith 1432d336d48SLois Curfman McInnes /* generate the scatter context */ 14408480c60SBarry Smith ierr = VecScatterCreate(gvec,from,aij->lvec,to,&aij->Mvctx);CHKERRQ(ierr); 14552e6d16bSBarry Smith ierr = PetscLogObjectParent(mat,aij->Mvctx);CHKERRQ(ierr); 14652e6d16bSBarry Smith ierr = PetscLogObjectParent(mat,aij->lvec);CHKERRQ(ierr); 14752e6d16bSBarry Smith ierr = PetscLogObjectParent(mat,from);CHKERRQ(ierr); 14852e6d16bSBarry Smith ierr = PetscLogObjectParent(mat,to);CHKERRQ(ierr); 1499e25ed09SBarry Smith aij->garray = garray; 15052e6d16bSBarry Smith ierr = PetscLogObjectMemory(mat,(ec+1)*sizeof(PetscInt));CHKERRQ(ierr); 15178b31e54SBarry Smith ierr = ISDestroy(from);CHKERRQ(ierr); 15278b31e54SBarry Smith ierr = ISDestroy(to);CHKERRQ(ierr); 153888f2ed8SSatish Balay ierr = VecDestroy(gvec);CHKERRQ(ierr); 1543a40ed3dSBarry Smith PetscFunctionReturn(0); 1558c79f6d3SBarry Smith } 1569e25ed09SBarry Smith 1579e25ed09SBarry Smith 1584a2ae208SSatish Balay #undef __FUNCT__ 1594a2ae208SSatish Balay #define __FUNCT__ "DisAssemble_MPIAIJ" 1602493cbb0SBarry Smith /* 1612493cbb0SBarry Smith Takes the local part of an already assembled MPIAIJ matrix 1622493cbb0SBarry Smith and disassembles it. This is to allow new nonzeros into the matrix 1632493cbb0SBarry Smith that require more communication in the matrix vector multiply. 1642493cbb0SBarry Smith Thus certain data-structures must be rebuilt. 1652493cbb0SBarry Smith 1662493cbb0SBarry Smith Kind of slow! But that's what application programmers get when 1672493cbb0SBarry Smith they are sloppy. 1682493cbb0SBarry Smith */ 169dfbe8321SBarry Smith PetscErrorCode DisAssemble_MPIAIJ(Mat A) 1702493cbb0SBarry Smith { 1712493cbb0SBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ*)A->data; 1722493cbb0SBarry Smith Mat B = aij->B,Bnew; 173ec8511deSBarry Smith Mat_SeqAIJ *Baij = (Mat_SeqAIJ*)B->data; 174dfbe8321SBarry Smith PetscErrorCode ierr; 17557809a77SBarry Smith PetscInt i,j,m = B->m,n = A->N,col,ct = 0,*garray = aij->garray,*nz,ec; 17687828ca2SBarry Smith PetscScalar v; 1772493cbb0SBarry Smith 1783a40ed3dSBarry Smith PetscFunctionBegin; 1792493cbb0SBarry Smith /* free stuff related to matrix-vec multiply */ 180b0a32e0cSBarry Smith ierr = VecGetSize(aij->lvec,&ec);CHKERRQ(ierr); /* needed for PetscLogObjectMemory below */ 1812493cbb0SBarry Smith ierr = VecDestroy(aij->lvec);CHKERRQ(ierr); aij->lvec = 0; 18208480c60SBarry Smith ierr = VecScatterDestroy(aij->Mvctx);CHKERRQ(ierr); aij->Mvctx = 0; 183464493b3SBarry Smith if (aij->colmap) { 184aa482453SBarry Smith #if defined (PETSC_USE_CTABLE) 1850f5bd95cSBarry Smith ierr = PetscTableDelete(aij->colmap);CHKERRQ(ierr); 1860f5bd95cSBarry Smith aij->colmap = 0; 1872066d6f7SSatish Balay #else 188606d414cSSatish Balay ierr = PetscFree(aij->colmap);CHKERRQ(ierr); 189606d414cSSatish Balay aij->colmap = 0; 19052e6d16bSBarry Smith ierr = PetscLogObjectMemory(A,-aij->B->n*sizeof(PetscInt));CHKERRQ(ierr); 1912066d6f7SSatish Balay #endif 192464493b3SBarry Smith } 1932493cbb0SBarry Smith 1942493cbb0SBarry Smith /* make sure that B is assembled so we can access its values */ 1956d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 196fe2f2677SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1972493cbb0SBarry Smith 1982493cbb0SBarry Smith /* invent new B and copy stuff over */ 199b1d57f15SBarry Smith ierr = PetscMalloc((m+1)*sizeof(PetscInt),&nz);CHKERRQ(ierr); 20048b35521SBarry Smith for (i=0; i<m; i++) { 20148b35521SBarry Smith nz[i] = Baij->i[i+1] - Baij->i[i]; 20248b35521SBarry Smith } 203f69a0ea3SMatthew Knepley ierr = MatCreate(PETSC_COMM_SELF,&Bnew);CHKERRQ(ierr); 204f69a0ea3SMatthew Knepley ierr = MatSetSizes(Bnew,m,n,m,n);CHKERRQ(ierr); 205f204ca49SKris Buschelman ierr = MatSetType(Bnew,B->type_name);CHKERRQ(ierr); 206f204ca49SKris Buschelman ierr = MatSeqAIJSetPreallocation(Bnew,0,nz);CHKERRQ(ierr); 207606d414cSSatish Balay ierr = PetscFree(nz);CHKERRQ(ierr); 2082493cbb0SBarry Smith for (i=0; i<m; i++) { 209bfec09a0SHong Zhang for (j=Baij->i[i]; j<Baij->i[i+1]; j++) { 210bfec09a0SHong Zhang col = garray[Baij->j[ct]]; 2112493cbb0SBarry Smith v = Baij->a[ct++]; 21283271157SBarry Smith ierr = MatSetValues(Bnew,1,&i,1,&col,&v,B->insertmode);CHKERRQ(ierr); 2132493cbb0SBarry Smith } 2142493cbb0SBarry Smith } 215606d414cSSatish Balay ierr = PetscFree(aij->garray);CHKERRQ(ierr); 216606d414cSSatish Balay aij->garray = 0; 21752e6d16bSBarry Smith ierr = PetscLogObjectMemory(A,-ec*sizeof(PetscInt));CHKERRQ(ierr); 2182493cbb0SBarry Smith ierr = MatDestroy(B);CHKERRQ(ierr); 21952e6d16bSBarry Smith ierr = PetscLogObjectParent(A,Bnew);CHKERRQ(ierr); 2202493cbb0SBarry Smith aij->B = Bnew; 221227d817aSBarry Smith A->was_assembled = PETSC_FALSE; 2223a40ed3dSBarry Smith PetscFunctionReturn(0); 2232493cbb0SBarry Smith } 2242493cbb0SBarry Smith 2252cd6534aSBarry Smith /* ugly stuff added for Glenn someday we should fix this up */ 2262cd6534aSBarry Smith 227b1d57f15SBarry Smith static PetscInt *auglyrmapd = 0,*auglyrmapo = 0; /* mapping from the local ordering to the "diagonal" and "off-diagonal" 2282cd6534aSBarry Smith parts of the local matrix */ 2292cd6534aSBarry Smith static Vec auglydd = 0,auglyoo = 0; /* work vectors used to scale the two parts of the local matrix */ 2302cd6534aSBarry Smith 2312cd6534aSBarry Smith 2322cd6534aSBarry Smith #undef __FUNCT__ 2332cd6534aSBarry Smith #define __FUNCT__ "MatMPIAIJDiagonalScaleLocalSetUp" 234dfbe8321SBarry Smith PetscErrorCode MatMPIAIJDiagonalScaleLocalSetUp(Mat inA,Vec scale) 2352cd6534aSBarry Smith { 2362cd6534aSBarry Smith Mat_MPIAIJ *ina = (Mat_MPIAIJ*) inA->data; /*access private part of matrix */ 237dfbe8321SBarry Smith PetscErrorCode ierr; 238b1d57f15SBarry Smith PetscInt i,n,nt,cstart,cend,no,*garray = ina->garray,*lindices; 239b1d57f15SBarry Smith PetscInt *r_rmapd,*r_rmapo; 2402cd6534aSBarry Smith 2412cd6534aSBarry Smith PetscFunctionBegin; 2422cd6534aSBarry Smith ierr = MatGetOwnershipRange(inA,&cstart,&cend);CHKERRQ(ierr); 2432cd6534aSBarry Smith ierr = MatGetSize(ina->A,PETSC_NULL,&n);CHKERRQ(ierr); 244b1d57f15SBarry Smith ierr = PetscMalloc((inA->mapping->n+1)*sizeof(PetscInt),&r_rmapd);CHKERRQ(ierr); 245b1d57f15SBarry Smith ierr = PetscMemzero(r_rmapd,inA->mapping->n*sizeof(PetscInt));CHKERRQ(ierr); 2462cd6534aSBarry Smith nt = 0; 2472cd6534aSBarry Smith for (i=0; i<inA->mapping->n; i++) { 2482cd6534aSBarry Smith if (inA->mapping->indices[i] >= cstart && inA->mapping->indices[i] < cend) { 2492cd6534aSBarry Smith nt++; 2502cd6534aSBarry Smith r_rmapd[i] = inA->mapping->indices[i] + 1; 2512cd6534aSBarry Smith } 2522cd6534aSBarry Smith } 25377431f27SBarry Smith if (nt != n) SETERRQ2(PETSC_ERR_PLIB,"Hmm nt %D n %D",nt,n); 254b1d57f15SBarry Smith ierr = PetscMalloc((n+1)*sizeof(PetscInt),&auglyrmapd);CHKERRQ(ierr); 2552cd6534aSBarry Smith for (i=0; i<inA->mapping->n; i++) { 2562cd6534aSBarry Smith if (r_rmapd[i]){ 2572cd6534aSBarry Smith auglyrmapd[(r_rmapd[i]-1)-cstart] = i; 2582cd6534aSBarry Smith } 2592cd6534aSBarry Smith } 2602cd6534aSBarry Smith ierr = PetscFree(r_rmapd);CHKERRQ(ierr); 2612cd6534aSBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,n,&auglydd);CHKERRQ(ierr); 2622cd6534aSBarry Smith 263b1d57f15SBarry Smith ierr = PetscMalloc((inA->N+1)*sizeof(PetscInt),&lindices);CHKERRQ(ierr); 264b1d57f15SBarry Smith ierr = PetscMemzero(lindices,inA->N*sizeof(PetscInt));CHKERRQ(ierr); 2652cd6534aSBarry Smith for (i=0; i<ina->B->n; i++) { 2662cd6534aSBarry Smith lindices[garray[i]] = i+1; 2672cd6534aSBarry Smith } 2682cd6534aSBarry Smith no = inA->mapping->n - nt; 269b1d57f15SBarry Smith ierr = PetscMalloc((inA->mapping->n+1)*sizeof(PetscInt),&r_rmapo);CHKERRQ(ierr); 270b1d57f15SBarry Smith ierr = PetscMemzero(r_rmapo,inA->mapping->n*sizeof(PetscInt));CHKERRQ(ierr); 2712cd6534aSBarry Smith nt = 0; 2722cd6534aSBarry Smith for (i=0; i<inA->mapping->n; i++) { 2732cd6534aSBarry Smith if (lindices[inA->mapping->indices[i]]) { 2742cd6534aSBarry Smith nt++; 2752cd6534aSBarry Smith r_rmapo[i] = lindices[inA->mapping->indices[i]]; 2762cd6534aSBarry Smith } 2772cd6534aSBarry Smith } 27877431f27SBarry Smith if (nt > no) SETERRQ2(PETSC_ERR_PLIB,"Hmm nt %D no %D",nt,n); 2792cd6534aSBarry Smith ierr = PetscFree(lindices);CHKERRQ(ierr); 280b1d57f15SBarry Smith ierr = PetscMalloc((nt+1)*sizeof(PetscInt),&auglyrmapo);CHKERRQ(ierr); 2812cd6534aSBarry Smith for (i=0; i<inA->mapping->n; i++) { 2822cd6534aSBarry Smith if (r_rmapo[i]){ 2832cd6534aSBarry Smith auglyrmapo[(r_rmapo[i]-1)] = i; 2842cd6534aSBarry Smith } 2852cd6534aSBarry Smith } 2862cd6534aSBarry Smith ierr = PetscFree(r_rmapo);CHKERRQ(ierr); 2872cd6534aSBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,nt,&auglyoo);CHKERRQ(ierr); 2882cd6534aSBarry Smith 2892cd6534aSBarry Smith PetscFunctionReturn(0); 2902cd6534aSBarry Smith } 2912cd6534aSBarry Smith 2922cd6534aSBarry Smith #undef __FUNCT__ 2932cd6534aSBarry Smith #define __FUNCT__ "MatMPIAIJDiagonalScaleLocal" 294dfbe8321SBarry Smith PetscErrorCode MatMPIAIJDiagonalScaleLocal(Mat A,Vec scale) 2952cd6534aSBarry Smith { 29692b32695SKris Buschelman /* This routine should really be abandoned as it duplicates MatDiagonalScaleLocal */ 297dfbe8321SBarry Smith PetscErrorCode ierr,(*f)(Mat,Vec); 29892b32695SKris Buschelman 29992b32695SKris Buschelman PetscFunctionBegin; 30092b32695SKris Buschelman ierr = PetscObjectQueryFunction((PetscObject)A,"MatDiagonalScaleLocal_C",(void (**)(void))&f);CHKERRQ(ierr); 30192b32695SKris Buschelman if (f) { 30292b32695SKris Buschelman ierr = (*f)(A,scale);CHKERRQ(ierr); 30392b32695SKris Buschelman } 30492b32695SKris Buschelman PetscFunctionReturn(0); 30592b32695SKris Buschelman } 30692b32695SKris Buschelman 30792b32695SKris Buschelman EXTERN_C_BEGIN 30892b32695SKris Buschelman #undef __FUNCT__ 30992b32695SKris Buschelman #define __FUNCT__ "MatDiagonalScaleLocal_MPIAIJ" 310be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalScaleLocal_MPIAIJ(Mat A,Vec scale) 31192b32695SKris Buschelman { 3122cd6534aSBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ*) A->data; /*access private part of matrix */ 313dfbe8321SBarry Smith PetscErrorCode ierr; 314b1d57f15SBarry Smith PetscInt n,i; 3152cd6534aSBarry Smith PetscScalar *d,*o,*s; 3162cd6534aSBarry Smith 3172cd6534aSBarry Smith PetscFunctionBegin; 3182cd6534aSBarry Smith if (!auglyrmapd) { 3192cd6534aSBarry Smith ierr = MatMPIAIJDiagonalScaleLocalSetUp(A,scale);CHKERRQ(ierr); 3202cd6534aSBarry Smith } 3212cd6534aSBarry Smith 3221ebc52fbSHong Zhang ierr = VecGetArray(scale,&s);CHKERRQ(ierr); 3232cd6534aSBarry Smith 3242cd6534aSBarry Smith ierr = VecGetLocalSize(auglydd,&n);CHKERRQ(ierr); 3251ebc52fbSHong Zhang ierr = VecGetArray(auglydd,&d);CHKERRQ(ierr); 3262cd6534aSBarry Smith for (i=0; i<n; i++) { 3272cd6534aSBarry Smith d[i] = s[auglyrmapd[i]]; /* copy "diagonal" (true local) portion of scale into dd vector */ 3282cd6534aSBarry Smith } 3291ebc52fbSHong Zhang ierr = VecRestoreArray(auglydd,&d);CHKERRQ(ierr); 3302cd6534aSBarry Smith /* column scale "diagonal" portion of local matrix */ 3312cd6534aSBarry Smith ierr = MatDiagonalScale(a->A,PETSC_NULL,auglydd);CHKERRQ(ierr); 3322cd6534aSBarry Smith 3332cd6534aSBarry Smith ierr = VecGetLocalSize(auglyoo,&n);CHKERRQ(ierr); 3341ebc52fbSHong Zhang ierr = VecGetArray(auglyoo,&o);CHKERRQ(ierr); 3352cd6534aSBarry Smith for (i=0; i<n; i++) { 3362cd6534aSBarry Smith o[i] = s[auglyrmapo[i]]; /* copy "off-diagonal" portion of scale into oo vector */ 3372cd6534aSBarry Smith } 3381ebc52fbSHong Zhang ierr = VecRestoreArray(scale,&s);CHKERRQ(ierr); 3391ebc52fbSHong Zhang ierr = VecRestoreArray(auglyoo,&o);CHKERRQ(ierr); 3402cd6534aSBarry Smith /* column scale "off-diagonal" portion of local matrix */ 3412cd6534aSBarry Smith ierr = MatDiagonalScale(a->B,PETSC_NULL,auglyoo);CHKERRQ(ierr); 3422cd6534aSBarry Smith 3432cd6534aSBarry Smith PetscFunctionReturn(0); 3442cd6534aSBarry Smith } 34592b32695SKris Buschelman EXTERN_C_END 3462cd6534aSBarry Smith 34748b35521SBarry Smith 348