173f4d377SMatthew Knepley /*$Id: mmaij.c,v 1.59 2001/08/07 03:02:49 balay Exp $*/ 28c79f6d3SBarry Smith 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" 1044a69424SLois Curfman McInnes int 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); 14273d9f13SBarry Smith int N = mat->N,i,j,*indices,*aj = B->j,ierr,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; 202066d6f7SSatish Balay int gid,lid; 212066d6f7SSatish Balay #endif 222066d6f7SSatish Balay 233a40ed3dSBarry Smith PetscFunctionBegin; 242066d6f7SSatish Balay 25aa482453SBarry Smith #if defined (PETSC_USE_CTABLE) 262066d6f7SSatish Balay /* use a table - Mark Adams (this has not been tested with "shift") */ 27273d9f13SBarry Smith ierr = PetscTableCreate(aij->B->m,&gid1_lid1);CHKERRQ(ierr); 28273d9f13SBarry Smith for (i=0; i<aij->B->m; i++) { 292066d6f7SSatish Balay for (j=0; j<B->ilen[i]; j++) { 30b3fb0a6cSMatthew Knepley int 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 */ 340f5bd95cSBarry Smith ierr = PetscTableAdd(gid1_lid1,gid1,++ec);CHKERRQ(ierr); 352066d6f7SSatish Balay } 362066d6f7SSatish Balay } 372066d6f7SSatish Balay } 382066d6f7SSatish Balay /* form array of columns we need */ 39b0a32e0cSBarry Smith ierr = PetscMalloc((ec+1)*sizeof(int),&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++) { 500f5bd95cSBarry Smith ierr = PetscTableAdd(gid1_lid1,garray[i]+1,i+1);CHKERRQ(ierr); 512066d6f7SSatish Balay } 522066d6f7SSatish Balay /* compact out the extra columns in B */ 53273d9f13SBarry Smith for (i=0; i<aij->B->m; i++) { 542066d6f7SSatish Balay for (j=0; j<B->ilen[i]; j++) { 55b3fb0a6cSMatthew Knepley int 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 } 61273d9f13SBarry Smith aij->B->n = aij->B->N = ec; 620f5bd95cSBarry Smith ierr = PetscTableDelete(gid1_lid1);CHKERRQ(ierr); 632066d6f7SSatish Balay /* Mark Adams */ 642066d6f7SSatish Balay #else 658c79f6d3SBarry Smith /* For the first stab we make an array as long as the number of columns */ 661eb62cbbSBarry Smith /* mark those columns that are in aij->B */ 67b0a32e0cSBarry Smith ierr = PetscMalloc((N+1)*sizeof(int),&indices);CHKERRQ(ierr); 68549d3d68SSatish Balay ierr = PetscMemzero(indices,N*sizeof(int));CHKERRQ(ierr); 69273d9f13SBarry Smith for (i=0; i<aij->B->m; i++) { 70d6dfbf8fSBarry Smith for (j=0; j<B->ilen[i]; j++) { 71bfec09a0SHong Zhang if (!indices[aj[B->i[i] + j] ]) ec++; 72bfec09a0SHong Zhang indices[aj[B->i[i] + j] ] = 1; 73416022c9SBarry Smith } 741eb62cbbSBarry Smith } 758c79f6d3SBarry Smith 761eb62cbbSBarry Smith /* form array of columns we need */ 77b0a32e0cSBarry Smith ierr = PetscMalloc((ec+1)*sizeof(int),&garray);CHKERRQ(ierr); 781eb62cbbSBarry Smith ec = 0; 791eb62cbbSBarry Smith for (i=0; i<N; i++) { 801eb62cbbSBarry Smith if (indices[i]) garray[ec++] = i; 811eb62cbbSBarry Smith } 821eb62cbbSBarry Smith 831eb62cbbSBarry Smith /* make indices now point into garray */ 841eb62cbbSBarry Smith for (i=0; i<ec; i++) { 85bfec09a0SHong Zhang indices[garray[i]] = i; 861eb62cbbSBarry Smith } 871eb62cbbSBarry Smith 881eb62cbbSBarry Smith /* compact out the extra columns in B */ 89273d9f13SBarry Smith for (i=0; i<aij->B->m; i++) { 90d6dfbf8fSBarry Smith for (j=0; j<B->ilen[i]; j++) { 91bfec09a0SHong Zhang aj[B->i[i] + j] = indices[aj[B->i[i] + j]]; 921eb62cbbSBarry Smith } 93d6dfbf8fSBarry Smith } 94273d9f13SBarry Smith aij->B->n = aij->B->N = ec; 95606d414cSSatish Balay ierr = PetscFree(indices);CHKERRQ(ierr); 962066d6f7SSatish Balay #endif 971eb62cbbSBarry Smith /* create local vector that is used to scatter into */ 98029af93fSBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,ec,&aij->lvec);CHKERRQ(ierr); 991eb62cbbSBarry Smith 100d6dfbf8fSBarry Smith /* create two temporary Index sets for build scatter gather */ 101b9b97703SBarry Smith ierr = ISCreateGeneral(mat->comm,ec,garray,&from);CHKERRQ(ierr); 102029af93fSBarry Smith ierr = ISCreateStride(PETSC_COMM_SELF,ec,0,1,&to);CHKERRQ(ierr); 1031eb62cbbSBarry Smith 1041eb62cbbSBarry Smith /* create temporary global vector to generate scatter context */ 1051eb62cbbSBarry Smith /* this is inefficient, but otherwise we must do either 1061eb62cbbSBarry Smith 1) save garray until the first actual scatter when the vector is known or 1071eb62cbbSBarry Smith 2) have another way of generating a scatter context without a vector.*/ 108273d9f13SBarry Smith ierr = VecCreateMPI(mat->comm,mat->n,mat->N,&gvec);CHKERRQ(ierr); 1091eb62cbbSBarry Smith 1102d336d48SLois Curfman McInnes /* generate the scatter context */ 11108480c60SBarry Smith ierr = VecScatterCreate(gvec,from,aij->lvec,to,&aij->Mvctx);CHKERRQ(ierr); 112b0a32e0cSBarry Smith PetscLogObjectParent(mat,aij->Mvctx); 113b0a32e0cSBarry Smith PetscLogObjectParent(mat,aij->lvec); 114b0a32e0cSBarry Smith PetscLogObjectParent(mat,from); 115b0a32e0cSBarry Smith PetscLogObjectParent(mat,to); 1169e25ed09SBarry Smith aij->garray = garray; 117b0a32e0cSBarry Smith PetscLogObjectMemory(mat,(ec+1)*sizeof(int)); 11878b31e54SBarry Smith ierr = ISDestroy(from);CHKERRQ(ierr); 11978b31e54SBarry Smith ierr = ISDestroy(to);CHKERRQ(ierr); 120888f2ed8SSatish Balay ierr = VecDestroy(gvec);CHKERRQ(ierr); 1213a40ed3dSBarry Smith PetscFunctionReturn(0); 1228c79f6d3SBarry Smith } 1239e25ed09SBarry Smith 1249e25ed09SBarry Smith 1254a2ae208SSatish Balay #undef __FUNCT__ 1264a2ae208SSatish Balay #define __FUNCT__ "DisAssemble_MPIAIJ" 1272493cbb0SBarry Smith /* 1282493cbb0SBarry Smith Takes the local part of an already assembled MPIAIJ matrix 1292493cbb0SBarry Smith and disassembles it. This is to allow new nonzeros into the matrix 1302493cbb0SBarry Smith that require more communication in the matrix vector multiply. 1312493cbb0SBarry Smith Thus certain data-structures must be rebuilt. 1322493cbb0SBarry Smith 1332493cbb0SBarry Smith Kind of slow! But that's what application programmers get when 1342493cbb0SBarry Smith they are sloppy. 1352493cbb0SBarry Smith */ 1362493cbb0SBarry Smith int DisAssemble_MPIAIJ(Mat A) 1372493cbb0SBarry Smith { 1382493cbb0SBarry Smith Mat_MPIAIJ *aij = (Mat_MPIAIJ*)A->data; 1392493cbb0SBarry Smith Mat B = aij->B,Bnew; 140ec8511deSBarry Smith Mat_SeqAIJ *Baij = (Mat_SeqAIJ*)B->data; 141273d9f13SBarry Smith int ierr,i,j,m = B->m,n = A->N,col,ct = 0,*garray = aij->garray; 142bfec09a0SHong Zhang int *nz,ec; 14387828ca2SBarry Smith PetscScalar v; 1442493cbb0SBarry Smith 1453a40ed3dSBarry Smith PetscFunctionBegin; 1462493cbb0SBarry Smith /* free stuff related to matrix-vec multiply */ 147b0a32e0cSBarry Smith ierr = VecGetSize(aij->lvec,&ec);CHKERRQ(ierr); /* needed for PetscLogObjectMemory below */ 1482493cbb0SBarry Smith ierr = VecDestroy(aij->lvec);CHKERRQ(ierr); aij->lvec = 0; 14908480c60SBarry Smith ierr = VecScatterDestroy(aij->Mvctx);CHKERRQ(ierr); aij->Mvctx = 0; 150464493b3SBarry Smith if (aij->colmap) { 151aa482453SBarry Smith #if defined (PETSC_USE_CTABLE) 1520f5bd95cSBarry Smith ierr = PetscTableDelete(aij->colmap);CHKERRQ(ierr); 1530f5bd95cSBarry Smith aij->colmap = 0; 1542066d6f7SSatish Balay #else 155606d414cSSatish Balay ierr = PetscFree(aij->colmap);CHKERRQ(ierr); 156606d414cSSatish Balay aij->colmap = 0; 157b0a32e0cSBarry Smith PetscLogObjectMemory(A,-aij->B->n*sizeof(int)); 1582066d6f7SSatish Balay #endif 159464493b3SBarry Smith } 1602493cbb0SBarry Smith 1612493cbb0SBarry Smith /* make sure that B is assembled so we can access its values */ 1626d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 163fe2f2677SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1642493cbb0SBarry Smith 1652493cbb0SBarry Smith /* invent new B and copy stuff over */ 166b0a32e0cSBarry Smith ierr = PetscMalloc((m+1)*sizeof(int),&nz);CHKERRQ(ierr); 16748b35521SBarry Smith for (i=0; i<m; i++) { 16848b35521SBarry Smith nz[i] = Baij->i[i+1] - Baij->i[i]; 16948b35521SBarry Smith } 170*f204ca49SKris Buschelman ierr = MatCreate(PETSC_COMM_SELF,m,n,m,n,&Bnew);CHKERRQ(ierr); 171*f204ca49SKris Buschelman ierr = MatSetType(Bnew,B->type_name);CHKERRQ(ierr); 172*f204ca49SKris Buschelman ierr = MatSeqAIJSetPreallocation(Bnew,0,nz);CHKERRQ(ierr); 173606d414cSSatish Balay ierr = PetscFree(nz);CHKERRQ(ierr); 1742493cbb0SBarry Smith for (i=0; i<m; i++) { 175bfec09a0SHong Zhang for (j=Baij->i[i]; j<Baij->i[i+1]; j++) { 176bfec09a0SHong Zhang col = garray[Baij->j[ct]]; 1772493cbb0SBarry Smith v = Baij->a[ct++]; 17883271157SBarry Smith ierr = MatSetValues(Bnew,1,&i,1,&col,&v,B->insertmode);CHKERRQ(ierr); 1792493cbb0SBarry Smith } 1802493cbb0SBarry Smith } 181606d414cSSatish Balay ierr = PetscFree(aij->garray);CHKERRQ(ierr); 182606d414cSSatish Balay aij->garray = 0; 183b0a32e0cSBarry Smith PetscLogObjectMemory(A,-ec*sizeof(int)); 1842493cbb0SBarry Smith ierr = MatDestroy(B);CHKERRQ(ierr); 185b0a32e0cSBarry Smith PetscLogObjectParent(A,Bnew); 1862493cbb0SBarry Smith aij->B = Bnew; 187227d817aSBarry Smith A->was_assembled = PETSC_FALSE; 1883a40ed3dSBarry Smith PetscFunctionReturn(0); 1892493cbb0SBarry Smith } 1902493cbb0SBarry Smith 1912cd6534aSBarry Smith /* ugly stuff added for Glenn someday we should fix this up */ 1922cd6534aSBarry Smith 1932cd6534aSBarry Smith static int *auglyrmapd = 0,*auglyrmapo = 0; /* mapping from the local ordering to the "diagonal" and "off-diagonal" 1942cd6534aSBarry Smith parts of the local matrix */ 1952cd6534aSBarry Smith static Vec auglydd = 0,auglyoo = 0; /* work vectors used to scale the two parts of the local matrix */ 1962cd6534aSBarry Smith 1972cd6534aSBarry Smith 1982cd6534aSBarry Smith #undef __FUNCT__ 1992cd6534aSBarry Smith #define __FUNCT__ "MatMPIAIJDiagonalScaleLocalSetUp" 2002cd6534aSBarry Smith int MatMPIAIJDiagonalScaleLocalSetUp(Mat inA,Vec scale) 2012cd6534aSBarry Smith { 2022cd6534aSBarry Smith Mat_MPIAIJ *ina = (Mat_MPIAIJ*) inA->data; /*access private part of matrix */ 2032cd6534aSBarry Smith int ierr,i,n,nt,cstart,cend,no,*garray = ina->garray,*lindices; 2042cd6534aSBarry Smith int *r_rmapd,*r_rmapo; 2052cd6534aSBarry Smith 2062cd6534aSBarry Smith PetscFunctionBegin; 2072cd6534aSBarry Smith ierr = MatGetOwnershipRange(inA,&cstart,&cend);CHKERRQ(ierr); 2082cd6534aSBarry Smith ierr = MatGetSize(ina->A,PETSC_NULL,&n);CHKERRQ(ierr); 2092cd6534aSBarry Smith ierr = PetscMalloc((inA->mapping->n+1)*sizeof(int),&r_rmapd);CHKERRQ(ierr); 2102cd6534aSBarry Smith ierr = PetscMemzero(r_rmapd,inA->mapping->n*sizeof(int));CHKERRQ(ierr); 2112cd6534aSBarry Smith nt = 0; 2122cd6534aSBarry Smith for (i=0; i<inA->mapping->n; i++) { 2132cd6534aSBarry Smith if (inA->mapping->indices[i] >= cstart && inA->mapping->indices[i] < cend) { 2142cd6534aSBarry Smith nt++; 2152cd6534aSBarry Smith r_rmapd[i] = inA->mapping->indices[i] + 1; 2162cd6534aSBarry Smith } 2172cd6534aSBarry Smith } 2182cd6534aSBarry Smith if (nt != n) SETERRQ2(1,"Hmm nt %d n %d",nt,n); 2192cd6534aSBarry Smith ierr = PetscMalloc((n+1)*sizeof(int),&auglyrmapd);CHKERRQ(ierr); 2202cd6534aSBarry Smith for (i=0; i<inA->mapping->n; i++) { 2212cd6534aSBarry Smith if (r_rmapd[i]){ 2222cd6534aSBarry Smith auglyrmapd[(r_rmapd[i]-1)-cstart] = i; 2232cd6534aSBarry Smith } 2242cd6534aSBarry Smith } 2252cd6534aSBarry Smith ierr = PetscFree(r_rmapd);CHKERRQ(ierr); 2262cd6534aSBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,n,&auglydd);CHKERRQ(ierr); 2272cd6534aSBarry Smith 2282cd6534aSBarry Smith ierr = PetscMalloc((inA->N+1)*sizeof(int),&lindices);CHKERRQ(ierr); 2292cd6534aSBarry Smith ierr = PetscMemzero(lindices,inA->N*sizeof(int));CHKERRQ(ierr); 2302cd6534aSBarry Smith for (i=0; i<ina->B->n; i++) { 2312cd6534aSBarry Smith lindices[garray[i]] = i+1; 2322cd6534aSBarry Smith } 2332cd6534aSBarry Smith no = inA->mapping->n - nt; 2342cd6534aSBarry Smith ierr = PetscMalloc((inA->mapping->n+1)*sizeof(int),&r_rmapo);CHKERRQ(ierr); 2352cd6534aSBarry Smith ierr = PetscMemzero(r_rmapo,inA->mapping->n*sizeof(int));CHKERRQ(ierr); 2362cd6534aSBarry Smith nt = 0; 2372cd6534aSBarry Smith for (i=0; i<inA->mapping->n; i++) { 2382cd6534aSBarry Smith if (lindices[inA->mapping->indices[i]]) { 2392cd6534aSBarry Smith nt++; 2402cd6534aSBarry Smith r_rmapo[i] = lindices[inA->mapping->indices[i]]; 2412cd6534aSBarry Smith } 2422cd6534aSBarry Smith } 2432cd6534aSBarry Smith if (nt > no) SETERRQ2(1,"Hmm nt %d no %d",nt,n); 2442cd6534aSBarry Smith ierr = PetscFree(lindices);CHKERRQ(ierr); 2452cd6534aSBarry Smith ierr = PetscMalloc((nt+1)*sizeof(int),&auglyrmapo);CHKERRQ(ierr); 2462cd6534aSBarry Smith for (i=0; i<inA->mapping->n; i++) { 2472cd6534aSBarry Smith if (r_rmapo[i]){ 2482cd6534aSBarry Smith auglyrmapo[(r_rmapo[i]-1)] = i; 2492cd6534aSBarry Smith } 2502cd6534aSBarry Smith } 2512cd6534aSBarry Smith ierr = PetscFree(r_rmapo);CHKERRQ(ierr); 2522cd6534aSBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,nt,&auglyoo);CHKERRQ(ierr); 2532cd6534aSBarry Smith 2542cd6534aSBarry Smith PetscFunctionReturn(0); 2552cd6534aSBarry Smith } 2562cd6534aSBarry Smith 2572cd6534aSBarry Smith #undef __FUNCT__ 2582cd6534aSBarry Smith #define __FUNCT__ "MatMPIAIJDiagonalScaleLocal" 2592cd6534aSBarry Smith int MatMPIAIJDiagonalScaleLocal(Mat A,Vec scale) 2602cd6534aSBarry Smith { 26192b32695SKris Buschelman /* This routine should really be abandoned as it duplicates MatDiagonalScaleLocal */ 26292b32695SKris Buschelman int ierr,(*f)(Mat,Vec); 26392b32695SKris Buschelman 26492b32695SKris Buschelman PetscFunctionBegin; 26592b32695SKris Buschelman ierr = PetscObjectQueryFunction((PetscObject)A,"MatDiagonalScaleLocal_C",(void (**)(void))&f);CHKERRQ(ierr); 26692b32695SKris Buschelman if (f) { 26792b32695SKris Buschelman ierr = (*f)(A,scale);CHKERRQ(ierr); 26892b32695SKris Buschelman } 26992b32695SKris Buschelman PetscFunctionReturn(0); 27092b32695SKris Buschelman } 27192b32695SKris Buschelman 27292b32695SKris Buschelman EXTERN_C_BEGIN 27392b32695SKris Buschelman #undef __FUNCT__ 27492b32695SKris Buschelman #define __FUNCT__ "MatDiagonalScaleLocal_MPIAIJ" 27592b32695SKris Buschelman int MatDiagonalScaleLocal_MPIAIJ(Mat A,Vec scale) 27692b32695SKris Buschelman { 2772cd6534aSBarry Smith Mat_MPIAIJ *a = (Mat_MPIAIJ*) A->data; /*access private part of matrix */ 2782cd6534aSBarry Smith int ierr,n,i; 2792cd6534aSBarry Smith PetscScalar *d,*o,*s; 2802cd6534aSBarry Smith 2812cd6534aSBarry Smith PetscFunctionBegin; 2822cd6534aSBarry Smith if (!auglyrmapd) { 2832cd6534aSBarry Smith ierr = MatMPIAIJDiagonalScaleLocalSetUp(A,scale);CHKERRQ(ierr); 2842cd6534aSBarry Smith } 2852cd6534aSBarry Smith 286b1d4fb26SBarry Smith ierr = VecGetArrayFast(scale,&s);CHKERRQ(ierr); 2872cd6534aSBarry Smith 2882cd6534aSBarry Smith ierr = VecGetLocalSize(auglydd,&n);CHKERRQ(ierr); 289b1d4fb26SBarry Smith ierr = VecGetArrayFast(auglydd,&d);CHKERRQ(ierr); 2902cd6534aSBarry Smith for (i=0; i<n; i++) { 2912cd6534aSBarry Smith d[i] = s[auglyrmapd[i]]; /* copy "diagonal" (true local) portion of scale into dd vector */ 2922cd6534aSBarry Smith } 293b1d4fb26SBarry Smith ierr = VecRestoreArrayFast(auglydd,&d);CHKERRQ(ierr); 2942cd6534aSBarry Smith /* column scale "diagonal" portion of local matrix */ 2952cd6534aSBarry Smith ierr = MatDiagonalScale(a->A,PETSC_NULL,auglydd);CHKERRQ(ierr); 2962cd6534aSBarry Smith 2972cd6534aSBarry Smith ierr = VecGetLocalSize(auglyoo,&n);CHKERRQ(ierr); 298b1d4fb26SBarry Smith ierr = VecGetArrayFast(auglyoo,&o);CHKERRQ(ierr); 2992cd6534aSBarry Smith for (i=0; i<n; i++) { 3002cd6534aSBarry Smith o[i] = s[auglyrmapo[i]]; /* copy "off-diagonal" portion of scale into oo vector */ 3012cd6534aSBarry Smith } 302b1d4fb26SBarry Smith ierr = VecRestoreArrayFast(scale,&s);CHKERRQ(ierr); 303b1d4fb26SBarry Smith ierr = VecRestoreArrayFast(auglyoo,&o);CHKERRQ(ierr); 3042cd6534aSBarry Smith /* column scale "off-diagonal" portion of local matrix */ 3052cd6534aSBarry Smith ierr = MatDiagonalScale(a->B,PETSC_NULL,auglyoo);CHKERRQ(ierr); 3062cd6534aSBarry Smith 3072cd6534aSBarry Smith PetscFunctionReturn(0); 3082cd6534aSBarry Smith } 30992b32695SKris Buschelman EXTERN_C_END 3102cd6534aSBarry Smith 31148b35521SBarry Smith 312