xref: /petsc/src/mat/impls/sbaij/mpi/mmsbaij.c (revision a11b0bdaf4192ea9b3ddaf9d3c9adfe6dc29abe9)
1a30f8f8cSSatish Balay 
2a30f8f8cSSatish Balay /*
3a30f8f8cSSatish Balay    Support for the parallel SBAIJ matrix vector multiply
4a30f8f8cSSatish Balay */
5c6db04a5SJed Brown #include <../src/mat/impls/sbaij/mpi/mpisbaij.h>
67cff1425SSatish Balay 
7dfbe8321SBarry Smith PetscErrorCode MatSetUpMultiply_MPISBAIJ(Mat mat)
8a30f8f8cSSatish Balay {
940781036SHong Zhang   Mat_MPISBAIJ   *sbaij = (Mat_MPISBAIJ*)mat->data;
1040781036SHong Zhang   Mat_SeqBAIJ    *B     = (Mat_SeqBAIJ*)(sbaij->B->data);
116849ba73SBarry Smith   PetscErrorCode ierr;
12*a11b0bdaSJunchao Zhang   PetscInt       Nbs = sbaij->Nbs,i,j,*aj = B->j,ec = 0,*garray,*sgarray;
13d0f46423SBarry Smith   PetscInt       bs  = mat->rmap->bs,*stmp,mbs=sbaij->mbs, vec_size,nt;
1440781036SHong Zhang   IS             from,to;
1540781036SHong Zhang   Vec            gvec;
16*a11b0bdaSJunchao Zhang   PetscMPIInt    rank = sbaij->rank,lsize;
17077aedafSJed Brown   PetscInt       *owners = sbaij->rangebs,*ec_owner,k;
18077aedafSJed Brown   const PetscInt *sowners;
1940781036SHong Zhang   PetscScalar    *ptr;
20*a11b0bdaSJunchao Zhang #if defined(PETSC_USE_CTABLE)
21*a11b0bdaSJunchao Zhang   PetscTable         gid1_lid1; /* one-based gid to lid table */
22*a11b0bdaSJunchao Zhang   PetscTablePosition tpos;
23*a11b0bdaSJunchao Zhang   PetscInt           gid,lid;
24*a11b0bdaSJunchao Zhang #else
25*a11b0bdaSJunchao Zhang   PetscInt           *indices;
26*a11b0bdaSJunchao Zhang #endif
2740781036SHong Zhang 
2840781036SHong Zhang   PetscFunctionBegin;
29*a11b0bdaSJunchao Zhang #if defined(PETSC_USE_CTABLE)
30*a11b0bdaSJunchao Zhang   ierr = PetscTableCreate(mbs,Nbs+1,&gid1_lid1);CHKERRQ(ierr);
31*a11b0bdaSJunchao Zhang   for (i=0; i<B->mbs; i++) {
32*a11b0bdaSJunchao Zhang     for (j=0; j<B->ilen[i]; j++) {
33*a11b0bdaSJunchao Zhang       PetscInt data,gid1 = aj[B->i[i]+j] + 1;
34*a11b0bdaSJunchao Zhang       ierr = PetscTableFind(gid1_lid1,gid1,&data);CHKERRQ(ierr);
35*a11b0bdaSJunchao Zhang       if (!data) {ierr = PetscTableAdd(gid1_lid1,gid1,++ec,INSERT_VALUES);CHKERRQ(ierr);}
36*a11b0bdaSJunchao Zhang     }
37*a11b0bdaSJunchao Zhang   }
38*a11b0bdaSJunchao Zhang   /* form array of columns we need */
39*a11b0bdaSJunchao Zhang   ierr = PetscMalloc1(ec,&garray);CHKERRQ(ierr);
40*a11b0bdaSJunchao Zhang   ierr = PetscTableGetHeadPosition(gid1_lid1,&tpos);CHKERRQ(ierr);
41*a11b0bdaSJunchao Zhang   while (tpos) {
42*a11b0bdaSJunchao Zhang     ierr = PetscTableGetNext(gid1_lid1,&tpos,&gid,&lid);CHKERRQ(ierr);
43*a11b0bdaSJunchao Zhang     gid--; lid--;
44*a11b0bdaSJunchao Zhang     garray[lid] = gid;
45*a11b0bdaSJunchao Zhang   }
46*a11b0bdaSJunchao Zhang   ierr = PetscSortInt(ec,garray);CHKERRQ(ierr);
47*a11b0bdaSJunchao Zhang   ierr = PetscTableRemoveAll(gid1_lid1);CHKERRQ(ierr);
48*a11b0bdaSJunchao Zhang   for (i=0; i<ec; i++) {ierr = PetscTableAdd(gid1_lid1,garray[i]+1,i+1,INSERT_VALUES);CHKERRQ(ierr);}
49*a11b0bdaSJunchao Zhang   /* compact out the extra columns in B */
50*a11b0bdaSJunchao Zhang   for (i=0; i<B->mbs; i++) {
51*a11b0bdaSJunchao Zhang     for (j=0; j<B->ilen[i]; j++) {
52*a11b0bdaSJunchao Zhang       PetscInt gid1 = aj[B->i[i] + j] + 1;
53*a11b0bdaSJunchao Zhang       ierr = PetscTableFind(gid1_lid1,gid1,&lid);CHKERRQ(ierr);
54*a11b0bdaSJunchao Zhang       lid--;
55*a11b0bdaSJunchao Zhang       aj[B->i[i]+j] = lid;
56*a11b0bdaSJunchao Zhang     }
57*a11b0bdaSJunchao Zhang   }
58*a11b0bdaSJunchao Zhang   ierr = PetscTableDestroy(&gid1_lid1);CHKERRQ(ierr);
59*a11b0bdaSJunchao Zhang   ierr = PetscMalloc2(2*ec,&sgarray,ec,&ec_owner);CHKERRQ(ierr);
60*a11b0bdaSJunchao Zhang   for (i=j=0; i<ec; i++) {
61*a11b0bdaSJunchao Zhang     while (garray[i]>=owners[j+1]) j++;
62*a11b0bdaSJunchao Zhang     ec_owner[i] = j;
63*a11b0bdaSJunchao Zhang   }
64*a11b0bdaSJunchao Zhang #else
6540781036SHong Zhang   /* For the first stab we make an array as long as the number of columns */
6640781036SHong Zhang   /* mark those columns that are in sbaij->B */
671795a4d1SJed Brown   ierr = PetscCalloc1(Nbs,&indices);CHKERRQ(ierr);
6840781036SHong Zhang   for (i=0; i<mbs; i++) {
6940781036SHong Zhang     for (j=0; j<B->ilen[i]; j++) {
7040781036SHong Zhang       if (!indices[aj[B->i[i] + j]]) ec++;
7140781036SHong Zhang       indices[aj[B->i[i] + j]] = 1;
7240781036SHong Zhang     }
7340781036SHong Zhang   }
7440781036SHong Zhang 
7540781036SHong Zhang   /* form arrays of columns we need */
76785e854fSJed Brown   ierr = PetscMalloc1(ec,&garray);CHKERRQ(ierr);
77dcca6d9dSJed Brown   ierr = PetscMalloc2(2*ec,&sgarray,ec,&ec_owner);CHKERRQ(ierr);
7840781036SHong Zhang 
7940781036SHong Zhang   ec = 0;
80*a11b0bdaSJunchao Zhang   for (j=0; j<sbaij->size; j++) {
8140781036SHong Zhang     for (i=owners[j]; i<owners[j+1]; i++) {
8240781036SHong Zhang       if (indices[i]) {
8340781036SHong Zhang         garray[ec]   = i;
8440781036SHong Zhang         ec_owner[ec] = j;
8540781036SHong Zhang         ec++;
8640781036SHong Zhang       }
8740781036SHong Zhang     }
8840781036SHong Zhang   }
8940781036SHong Zhang 
9040781036SHong Zhang   /* make indices now point into garray */
91b424e231SHong Zhang   for (i=0; i<ec; i++) indices[garray[i]] = i;
9240781036SHong Zhang 
9340781036SHong Zhang   /* compact out the extra columns in B */
9440781036SHong Zhang   for (i=0; i<mbs; i++) {
9540781036SHong Zhang     for (j=0; j<B->ilen[i]; j++) aj[B->i[i] + j] = indices[aj[B->i[i] + j]];
9640781036SHong Zhang   }
97*a11b0bdaSJunchao Zhang   ierr = PetscFree(indices);CHKERRQ(ierr);
98*a11b0bdaSJunchao Zhang #endif
9940781036SHong Zhang   B->nbs = ec;
100ca5434daSLawrence Mitchell   ierr = PetscLayoutDestroy(&sbaij->B->cmap);CHKERRQ(ierr);
101ca5434daSLawrence Mitchell   ierr = PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)sbaij->B),ec*mat->rmap->bs,ec*mat->rmap->bs,mat->rmap->bs,&sbaij->B->cmap);CHKERRQ(ierr);
10240781036SHong Zhang 
103*a11b0bdaSJunchao Zhang   ierr = VecScatterDestroy(&sbaij->sMvctx);CHKERRQ(ierr);
10440781036SHong Zhang   /* create local vector that is used to scatter into */
10540781036SHong Zhang   ierr = VecCreateSeq(PETSC_COMM_SELF,ec*bs,&sbaij->lvec);CHKERRQ(ierr);
10640781036SHong Zhang 
10740781036SHong Zhang   /* create two temporary index sets for building scatter-gather */
108785e854fSJed Brown   ierr = PetscMalloc1(2*ec,&stmp);CHKERRQ(ierr);
109deff0451SBarry Smith   ierr = ISCreateBlock(PETSC_COMM_SELF,bs,ec,garray,PETSC_COPY_VALUES,&from);CHKERRQ(ierr);
11026fbe8dcSKarl Rupp   for (i=0; i<ec; i++) stmp[i] = i;
111deff0451SBarry Smith   ierr = ISCreateBlock(PETSC_COMM_SELF,bs,ec,stmp,PETSC_COPY_VALUES,&to);CHKERRQ(ierr);
11240781036SHong Zhang 
11393d1592dSHong Zhang   /* generate the scatter context
1144c500f23SPierre Jolivet      -- Mvctx and lvec are not used by MatMult_MPISBAIJ(), but usefull for some applications */
115ce94432eSBarry Smith   ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)mat),1,mat->cmap->n,mat->cmap->N,NULL,&gvec);CHKERRQ(ierr);
1169448b7f1SJunchao Zhang   ierr = VecScatterCreate(gvec,from,sbaij->lvec,to,&sbaij->Mvctx);CHKERRQ(ierr);
1176bf464f9SBarry Smith   ierr = VecDestroy(&gvec);CHKERRQ(ierr);
11840781036SHong Zhang 
11940781036SHong Zhang   sbaij->garray = garray;
12026fbe8dcSKarl Rupp 
1213bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->Mvctx);CHKERRQ(ierr);
1223bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->lvec);CHKERRQ(ierr);
1233bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)from);CHKERRQ(ierr);
1243bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)to);CHKERRQ(ierr);
12540781036SHong Zhang 
1266bf464f9SBarry Smith   ierr = ISDestroy(&from);CHKERRQ(ierr);
1276bf464f9SBarry Smith   ierr = ISDestroy(&to);CHKERRQ(ierr);
12840781036SHong Zhang 
12940781036SHong Zhang   /* create parallel vector that is used by SBAIJ matrix to scatter from/into */
13040781036SHong Zhang   lsize = (mbs + ec)*bs;
131ce94432eSBarry Smith   ierr  = VecCreateMPI(PetscObjectComm((PetscObject)mat),lsize,PETSC_DETERMINE,&sbaij->slvec0);CHKERRQ(ierr);
13240781036SHong Zhang   ierr  = VecDuplicate(sbaij->slvec0,&sbaij->slvec1);CHKERRQ(ierr);
13340781036SHong Zhang   ierr  = VecGetSize(sbaij->slvec0,&vec_size);CHKERRQ(ierr);
13440781036SHong Zhang 
135077aedafSJed Brown   ierr = VecGetOwnershipRanges(sbaij->slvec0,&sowners);CHKERRQ(ierr);
13640781036SHong Zhang 
13740781036SHong Zhang   /* x index in the IS sfrom */
13840781036SHong Zhang   for (i=0; i<ec; i++) {
13940781036SHong Zhang     j          = ec_owner[i];
14040781036SHong Zhang     sgarray[i] = garray[i] + (sowners[j]/bs - owners[j]);
14140781036SHong Zhang   }
14240781036SHong Zhang   /* b index in the IS sfrom */
14340781036SHong Zhang   k = sowners[rank]/bs + mbs;
14440781036SHong Zhang   for (i=ec,j=0; i< 2*ec; i++,j++) sgarray[i] = k + j;
145deff0451SBarry Smith   ierr = ISCreateBlock(PETSC_COMM_SELF,bs,2*ec,sgarray,PETSC_COPY_VALUES,&from);CHKERRQ(ierr);
14640781036SHong Zhang 
14740781036SHong Zhang   /* x index in the IS sto */
14840781036SHong Zhang   k = sowners[rank]/bs + mbs;
149e82e9f6bSBarry Smith   for (i=0; i<ec; i++) stmp[i] = (k + i);
15040781036SHong Zhang   /* b index in the IS sto */
151e82e9f6bSBarry Smith   for (i=ec; i<2*ec; i++) stmp[i] = sgarray[i-ec];
15240781036SHong Zhang 
153deff0451SBarry Smith   ierr = ISCreateBlock(PETSC_COMM_SELF,bs,2*ec,stmp,PETSC_COPY_VALUES,&to);CHKERRQ(ierr);
15440781036SHong Zhang 
1559448b7f1SJunchao Zhang   ierr = VecScatterCreate(sbaij->slvec0,from,sbaij->slvec1,to,&sbaij->sMvctx);CHKERRQ(ierr);
15640781036SHong Zhang 
15740781036SHong Zhang   ierr = VecGetLocalSize(sbaij->slvec1,&nt);CHKERRQ(ierr);
1581ebc52fbSHong Zhang   ierr = VecGetArray(sbaij->slvec1,&ptr);CHKERRQ(ierr);
159778a2246SBarry Smith   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,bs*mbs,ptr,&sbaij->slvec1a);CHKERRQ(ierr);
160778a2246SBarry Smith   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,nt-bs*mbs,ptr+bs*mbs,&sbaij->slvec1b);CHKERRQ(ierr);
1611ebc52fbSHong Zhang   ierr = VecRestoreArray(sbaij->slvec1,&ptr);CHKERRQ(ierr);
16240781036SHong Zhang 
1631ebc52fbSHong Zhang   ierr = VecGetArray(sbaij->slvec0,&ptr);CHKERRQ(ierr);
164778a2246SBarry Smith   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,nt-bs*mbs,ptr+bs*mbs,&sbaij->slvec0b);CHKERRQ(ierr);
1651ebc52fbSHong Zhang   ierr = VecRestoreArray(sbaij->slvec0,&ptr);CHKERRQ(ierr);
16640781036SHong Zhang 
16740781036SHong Zhang   ierr = PetscFree(stmp);CHKERRQ(ierr);
168ce94432eSBarry Smith   ierr = MPI_Barrier(PetscObjectComm((PetscObject)mat));CHKERRQ(ierr);
16940781036SHong Zhang 
1703bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->sMvctx);CHKERRQ(ierr);
1713bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->slvec0);CHKERRQ(ierr);
1723bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->slvec1);CHKERRQ(ierr);
1733bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->slvec0b);CHKERRQ(ierr);
1743bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->slvec1a);CHKERRQ(ierr);
1753bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->slvec1b);CHKERRQ(ierr);
1763bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)from);CHKERRQ(ierr);
1773bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)to);CHKERRQ(ierr);
17840781036SHong Zhang 
1793bb1ff40SBarry Smith   ierr = PetscLogObjectMemory((PetscObject)mat,(ec+1)*sizeof(PetscInt));CHKERRQ(ierr);
1806bf464f9SBarry Smith   ierr = ISDestroy(&from);CHKERRQ(ierr);
1816bf464f9SBarry Smith   ierr = ISDestroy(&to);CHKERRQ(ierr);
18274ed9c26SBarry Smith   ierr = PetscFree2(sgarray,ec_owner);CHKERRQ(ierr);
18340781036SHong Zhang   PetscFunctionReturn(0);
18440781036SHong Zhang }
18540781036SHong Zhang 
186a30f8f8cSSatish Balay /*
18701b2bd88SHong Zhang      Takes the local part of an already assembled MPISBAIJ matrix
188a30f8f8cSSatish Balay    and disassembles it. This is to allow new nonzeros into the matrix
189a30f8f8cSSatish Balay    that require more communication in the matrix vector multiply.
190a30f8f8cSSatish Balay    Thus certain data-structures must be rebuilt.
191a30f8f8cSSatish Balay 
192a30f8f8cSSatish Balay    Kind of slow! But that's what application programmers get when
193a30f8f8cSSatish Balay    they are sloppy.
194a30f8f8cSSatish Balay */
195ab9863d7SBarry Smith PetscErrorCode MatDisAssemble_MPISBAIJ(Mat A)
196a30f8f8cSSatish Balay {
1972896befeSSatish Balay   Mat_MPISBAIJ   *baij  = (Mat_MPISBAIJ*)A->data;
198a30f8f8cSSatish Balay   Mat            B      = baij->B,Bnew;
199a30f8f8cSSatish Balay   Mat_SeqBAIJ    *Bbaij = (Mat_SeqBAIJ*)B->data;
200dfbe8321SBarry Smith   PetscErrorCode ierr;
201d0f46423SBarry Smith   PetscInt       i,j,mbs=Bbaij->mbs,n = A->cmap->N,col,*garray=baij->garray;
202d0f46423SBarry Smith   PetscInt       k,bs=A->rmap->bs,bs2=baij->bs2,*rvals,*nz,ec,m=A->rmap->n;
203a30f8f8cSSatish Balay   MatScalar      *a = Bbaij->a;
20487828ca2SBarry Smith   PetscScalar    *atmp;
205ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL_MAT_SINGLE)
20613f74950SBarry Smith   PetscInt l;
207a30f8f8cSSatish Balay #endif
208a30f8f8cSSatish Balay 
209a30f8f8cSSatish Balay   PetscFunctionBegin;
210ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL_MAT_SINGLE)
211785e854fSJed Brown   ierr = PetscMalloc1(A->rmap->bs,&atmp);
21282502324SSatish Balay #endif
213a30f8f8cSSatish Balay   /* free stuff related to matrix-vec multiply */
214b0a32e0cSBarry Smith   ierr = VecGetSize(baij->lvec,&ec);CHKERRQ(ierr); /* needed for PetscLogObjectMemory below */
2156bf464f9SBarry Smith   ierr = VecDestroy(&baij->lvec);CHKERRQ(ierr);
2166bf464f9SBarry Smith   ierr = VecScatterDestroy(&baij->Mvctx);CHKERRQ(ierr);
21701b2bd88SHong Zhang 
2186bf464f9SBarry Smith   ierr = VecDestroy(&baij->slvec0);CHKERRQ(ierr);
2196bf464f9SBarry Smith   ierr = VecDestroy(&baij->slvec0b);CHKERRQ(ierr);
2206bf464f9SBarry Smith   ierr = VecDestroy(&baij->slvec1);CHKERRQ(ierr);
2216bf464f9SBarry Smith   ierr = VecDestroy(&baij->slvec1a);CHKERRQ(ierr);
2226bf464f9SBarry Smith   ierr = VecDestroy(&baij->slvec1b);CHKERRQ(ierr);
22301b2bd88SHong Zhang 
224a30f8f8cSSatish Balay   if (baij->colmap) {
225a30f8f8cSSatish Balay #if defined(PETSC_USE_CTABLE)
2266bc0bbbfSBarry Smith     ierr = PetscTableDestroy(&baij->colmap);CHKERRQ(ierr);
227a30f8f8cSSatish Balay #else
228a30f8f8cSSatish Balay     ierr = PetscFree(baij->colmap);CHKERRQ(ierr);
2293bb1ff40SBarry Smith     ierr = PetscLogObjectMemory((PetscObject)A,-Bbaij->nbs*sizeof(PetscInt));CHKERRQ(ierr);
230a30f8f8cSSatish Balay #endif
231a30f8f8cSSatish Balay   }
232a30f8f8cSSatish Balay 
233a30f8f8cSSatish Balay   /* make sure that B is assembled so we can access its values */
234a30f8f8cSSatish Balay   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
235a30f8f8cSSatish Balay   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
236a30f8f8cSSatish Balay 
237a30f8f8cSSatish Balay   /* invent new B and copy stuff over */
238785e854fSJed Brown   ierr = PetscMalloc1(mbs,&nz);CHKERRQ(ierr);
239a30f8f8cSSatish Balay   for (i=0; i<mbs; i++) {
240a30f8f8cSSatish Balay     nz[i] = Bbaij->i[i+1]-Bbaij->i[i];
241a30f8f8cSSatish Balay   }
242f69a0ea3SMatthew Knepley   ierr = MatCreate(PETSC_COMM_SELF,&Bnew);CHKERRQ(ierr);
243f69a0ea3SMatthew Knepley   ierr = MatSetSizes(Bnew,m,n,m,n);CHKERRQ(ierr);
2447adad957SLisandro Dalcin   ierr = MatSetType(Bnew,((PetscObject)B)->type_name);CHKERRQ(ierr);
245d0f46423SBarry Smith   ierr = MatSeqBAIJSetPreallocation(Bnew,B->rmap->bs,0,nz);CHKERRQ(ierr);
246a30f8f8cSSatish Balay   ierr = PetscFree(nz);CHKERRQ(ierr);
247a30f8f8cSSatish Balay 
248b38c15b3SStefano Zampini   if (Bbaij->nonew >= 0) { /* Inherit insertion error options (if positive). */
249b38c15b3SStefano Zampini     ((Mat_SeqSBAIJ*)Bnew->data)->nonew = Bbaij->nonew;
250b38c15b3SStefano Zampini   }
251b38c15b3SStefano Zampini 
25277341eacSDmitry Karpeev   /*
25377341eacSDmitry Karpeev    Ensure that B's nonzerostate is monotonically increasing.
25477341eacSDmitry Karpeev    Or should this follow the MatSetValues() loop to preserve B's nonzerstate across a MatDisAssemble() call?
25577341eacSDmitry Karpeev    */
25677341eacSDmitry Karpeev   Bnew->nonzerostate = B->nonzerostate;
257785e854fSJed Brown   ierr = PetscMalloc1(bs,&rvals);CHKERRQ(ierr);
258a30f8f8cSSatish Balay   for (i=0; i<mbs; i++) {
259a30f8f8cSSatish Balay     rvals[0] = bs*i;
26026fbe8dcSKarl Rupp     for (j=1; j<bs; j++) rvals[j] = rvals[j-1] + 1;
261a30f8f8cSSatish Balay     for (j=Bbaij->i[i]; j<Bbaij->i[i+1]; j++) {
262a30f8f8cSSatish Balay       col = garray[Bbaij->j[j]]*bs;
263a30f8f8cSSatish Balay       for (k=0; k<bs; k++) {
264ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL_MAT_SINGLE)
265a30f8f8cSSatish Balay         for (l=0; l<bs; l++) atmp[l] = a[j*bs2+l];
266a30f8f8cSSatish Balay #else
26771730473SSatish Balay         atmp = a+j*bs2 + k*bs;
268a30f8f8cSSatish Balay #endif
269c8407628SSatish Balay         ierr = MatSetValues_SeqSBAIJ(Bnew,bs,rvals,1,&col,atmp,B->insertmode);CHKERRQ(ierr);
270a30f8f8cSSatish Balay         col++;
271a30f8f8cSSatish Balay       }
272a30f8f8cSSatish Balay     }
273a30f8f8cSSatish Balay   }
274ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL_MAT_SINGLE)
275a30f8f8cSSatish Balay   ierr = PetscFree(atmp);CHKERRQ(ierr);
276a30f8f8cSSatish Balay #endif
277a30f8f8cSSatish Balay   ierr = PetscFree(baij->garray);CHKERRQ(ierr);
27826fbe8dcSKarl Rupp 
279a30f8f8cSSatish Balay   baij->garray = 0;
28026fbe8dcSKarl Rupp 
281a30f8f8cSSatish Balay   ierr = PetscFree(rvals);CHKERRQ(ierr);
2823bb1ff40SBarry Smith   ierr = PetscLogObjectMemory((PetscObject)A,-ec*sizeof(PetscInt));CHKERRQ(ierr);
2836bf464f9SBarry Smith   ierr = MatDestroy(&B);CHKERRQ(ierr);
2843bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)A,(PetscObject)Bnew);CHKERRQ(ierr);
28526fbe8dcSKarl Rupp 
286a30f8f8cSSatish Balay   baij->B = Bnew;
28726fbe8dcSKarl Rupp 
288a30f8f8cSSatish Balay   A->was_assembled = PETSC_FALSE;
289a30f8f8cSSatish Balay   PetscFunctionReturn(0);
290a30f8f8cSSatish Balay }
291