xref: /petsc/src/mat/impls/sbaij/mpi/mmsbaij.c (revision 4c500f23f41a3f3a0e9b2cac83fcce6b0e919e07)
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;
1213f74950SBarry Smith   PetscInt       Nbs = sbaij->Nbs,i,j,*indices,*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;
1613f74950SBarry Smith   PetscMPIInt    rank   =sbaij->rank,lsize,size=sbaij->size;
17077aedafSJed Brown   PetscInt       *owners=sbaij->rangebs,*ec_owner,k;
18077aedafSJed Brown   const PetscInt *sowners;
1940781036SHong Zhang   PetscScalar    *ptr;
2040781036SHong Zhang 
2140781036SHong Zhang   PetscFunctionBegin;
226bf464f9SBarry Smith   ierr = VecScatterDestroy(&sbaij->sMvctx);CHKERRQ(ierr);
2340781036SHong Zhang 
2440781036SHong Zhang   /* For the first stab we make an array as long as the number of columns */
2540781036SHong Zhang   /* mark those columns that are in sbaij->B */
261795a4d1SJed Brown   ierr = PetscCalloc1(Nbs,&indices);CHKERRQ(ierr);
2740781036SHong Zhang   for (i=0; i<mbs; i++) {
2840781036SHong Zhang     for (j=0; j<B->ilen[i]; j++) {
2940781036SHong Zhang       if (!indices[aj[B->i[i] + j]]) ec++;
3040781036SHong Zhang       indices[aj[B->i[i] + j]] = 1;
3140781036SHong Zhang     }
3240781036SHong Zhang   }
3340781036SHong Zhang 
3440781036SHong Zhang   /* form arrays of columns we need */
35785e854fSJed Brown   ierr = PetscMalloc1(ec,&garray);CHKERRQ(ierr);
36dcca6d9dSJed Brown   ierr = PetscMalloc2(2*ec,&sgarray,ec,&ec_owner);CHKERRQ(ierr);
3740781036SHong Zhang 
3840781036SHong Zhang   ec = 0;
3940781036SHong Zhang   for (j=0; j<size; j++) {
4040781036SHong Zhang     for (i=owners[j]; i<owners[j+1]; i++) {
4140781036SHong Zhang       if (indices[i]) {
4240781036SHong Zhang         garray[ec]   = i;
4340781036SHong Zhang         ec_owner[ec] = j;
4440781036SHong Zhang         ec++;
4540781036SHong Zhang       }
4640781036SHong Zhang     }
4740781036SHong Zhang   }
4840781036SHong Zhang 
4940781036SHong Zhang   /* make indices now point into garray */
50b424e231SHong Zhang   for (i=0; i<ec; i++) indices[garray[i]] = i;
5140781036SHong Zhang 
5240781036SHong Zhang   /* compact out the extra columns in B */
5340781036SHong Zhang   for (i=0; i<mbs; i++) {
5440781036SHong Zhang     for (j=0; j<B->ilen[i]; j++) aj[B->i[i] + j] = indices[aj[B->i[i] + j]];
5540781036SHong Zhang   }
5640781036SHong Zhang   B->nbs = ec;
5726fbe8dcSKarl Rupp 
58d0f46423SBarry Smith   sbaij->B->cmap->n = sbaij->B->cmap->N = ec*mat->rmap->bs;
5926fbe8dcSKarl Rupp 
6026283091SBarry Smith   ierr = PetscLayoutSetUp((sbaij->B->cmap));CHKERRQ(ierr);
6140781036SHong Zhang   ierr = PetscFree(indices);CHKERRQ(ierr);
6240781036SHong Zhang 
6340781036SHong Zhang   /* create local vector that is used to scatter into */
6440781036SHong Zhang   ierr = VecCreateSeq(PETSC_COMM_SELF,ec*bs,&sbaij->lvec);CHKERRQ(ierr);
6540781036SHong Zhang 
6640781036SHong Zhang   /* create two temporary index sets for building scatter-gather */
67785e854fSJed Brown   ierr = PetscMalloc1(2*ec,&stmp);CHKERRQ(ierr);
68deff0451SBarry Smith   ierr = ISCreateBlock(PETSC_COMM_SELF,bs,ec,garray,PETSC_COPY_VALUES,&from);CHKERRQ(ierr);
6926fbe8dcSKarl Rupp   for (i=0; i<ec; i++) stmp[i] = i;
70deff0451SBarry Smith   ierr = ISCreateBlock(PETSC_COMM_SELF,bs,ec,stmp,PETSC_COPY_VALUES,&to);CHKERRQ(ierr);
7140781036SHong Zhang 
7293d1592dSHong Zhang   /* generate the scatter context
73*4c500f23SPierre Jolivet      -- Mvctx and lvec are not used by MatMult_MPISBAIJ(), but usefull for some applications */
74ce94432eSBarry Smith   ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)mat),1,mat->cmap->n,mat->cmap->N,NULL,&gvec);CHKERRQ(ierr);
759448b7f1SJunchao Zhang   ierr = VecScatterCreate(gvec,from,sbaij->lvec,to,&sbaij->Mvctx);CHKERRQ(ierr);
766bf464f9SBarry Smith   ierr = VecDestroy(&gvec);CHKERRQ(ierr);
7740781036SHong Zhang 
7840781036SHong Zhang   sbaij->garray = garray;
7926fbe8dcSKarl Rupp 
803bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->Mvctx);CHKERRQ(ierr);
813bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->lvec);CHKERRQ(ierr);
823bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)from);CHKERRQ(ierr);
833bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)to);CHKERRQ(ierr);
8440781036SHong Zhang 
856bf464f9SBarry Smith   ierr = ISDestroy(&from);CHKERRQ(ierr);
866bf464f9SBarry Smith   ierr = ISDestroy(&to);CHKERRQ(ierr);
8740781036SHong Zhang 
8840781036SHong Zhang   /* create parallel vector that is used by SBAIJ matrix to scatter from/into */
8940781036SHong Zhang   lsize = (mbs + ec)*bs;
90ce94432eSBarry Smith   ierr  = VecCreateMPI(PetscObjectComm((PetscObject)mat),lsize,PETSC_DETERMINE,&sbaij->slvec0);CHKERRQ(ierr);
9140781036SHong Zhang   ierr  = VecDuplicate(sbaij->slvec0,&sbaij->slvec1);CHKERRQ(ierr);
9240781036SHong Zhang   ierr  = VecGetSize(sbaij->slvec0,&vec_size);CHKERRQ(ierr);
9340781036SHong Zhang 
94077aedafSJed Brown   ierr = VecGetOwnershipRanges(sbaij->slvec0,&sowners);CHKERRQ(ierr);
9540781036SHong Zhang 
9640781036SHong Zhang   /* x index in the IS sfrom */
9740781036SHong Zhang   for (i=0; i<ec; i++) {
9840781036SHong Zhang     j          = ec_owner[i];
9940781036SHong Zhang     sgarray[i] = garray[i] + (sowners[j]/bs - owners[j]);
10040781036SHong Zhang   }
10140781036SHong Zhang   /* b index in the IS sfrom */
10240781036SHong Zhang   k = sowners[rank]/bs + mbs;
10340781036SHong Zhang   for (i=ec,j=0; i< 2*ec; i++,j++) sgarray[i] = k + j;
104deff0451SBarry Smith   ierr = ISCreateBlock(PETSC_COMM_SELF,bs,2*ec,sgarray,PETSC_COPY_VALUES,&from);CHKERRQ(ierr);
10540781036SHong Zhang 
10640781036SHong Zhang   /* x index in the IS sto */
10740781036SHong Zhang   k = sowners[rank]/bs + mbs;
108e82e9f6bSBarry Smith   for (i=0; i<ec; i++) stmp[i] = (k + i);
10940781036SHong Zhang   /* b index in the IS sto */
110e82e9f6bSBarry Smith   for (i=ec; i<2*ec; i++) stmp[i] = sgarray[i-ec];
11140781036SHong Zhang 
112deff0451SBarry Smith   ierr = ISCreateBlock(PETSC_COMM_SELF,bs,2*ec,stmp,PETSC_COPY_VALUES,&to);CHKERRQ(ierr);
11340781036SHong Zhang 
1149448b7f1SJunchao Zhang   ierr = VecScatterCreate(sbaij->slvec0,from,sbaij->slvec1,to,&sbaij->sMvctx);CHKERRQ(ierr);
11540781036SHong Zhang 
11640781036SHong Zhang   ierr = VecGetLocalSize(sbaij->slvec1,&nt);CHKERRQ(ierr);
1171ebc52fbSHong Zhang   ierr = VecGetArray(sbaij->slvec1,&ptr);CHKERRQ(ierr);
118778a2246SBarry Smith   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,bs*mbs,ptr,&sbaij->slvec1a);CHKERRQ(ierr);
119778a2246SBarry Smith   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,nt-bs*mbs,ptr+bs*mbs,&sbaij->slvec1b);CHKERRQ(ierr);
1201ebc52fbSHong Zhang   ierr = VecRestoreArray(sbaij->slvec1,&ptr);CHKERRQ(ierr);
12140781036SHong Zhang 
1221ebc52fbSHong Zhang   ierr = VecGetArray(sbaij->slvec0,&ptr);CHKERRQ(ierr);
123778a2246SBarry Smith   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,nt-bs*mbs,ptr+bs*mbs,&sbaij->slvec0b);CHKERRQ(ierr);
1241ebc52fbSHong Zhang   ierr = VecRestoreArray(sbaij->slvec0,&ptr);CHKERRQ(ierr);
12540781036SHong Zhang 
12640781036SHong Zhang   ierr = PetscFree(stmp);CHKERRQ(ierr);
127ce94432eSBarry Smith   ierr = MPI_Barrier(PetscObjectComm((PetscObject)mat));CHKERRQ(ierr);
12840781036SHong Zhang 
1293bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->sMvctx);CHKERRQ(ierr);
1303bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->slvec0);CHKERRQ(ierr);
1313bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->slvec1);CHKERRQ(ierr);
1323bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->slvec0b);CHKERRQ(ierr);
1333bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->slvec1a);CHKERRQ(ierr);
1343bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->slvec1b);CHKERRQ(ierr);
1353bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)from);CHKERRQ(ierr);
1363bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)to);CHKERRQ(ierr);
13740781036SHong Zhang 
1383bb1ff40SBarry Smith   ierr = PetscLogObjectMemory((PetscObject)mat,(ec+1)*sizeof(PetscInt));CHKERRQ(ierr);
1396bf464f9SBarry Smith   ierr = ISDestroy(&from);CHKERRQ(ierr);
1406bf464f9SBarry Smith   ierr = ISDestroy(&to);CHKERRQ(ierr);
14174ed9c26SBarry Smith   ierr = PetscFree2(sgarray,ec_owner);CHKERRQ(ierr);
14240781036SHong Zhang   PetscFunctionReturn(0);
14340781036SHong Zhang }
14440781036SHong Zhang 
145dfbe8321SBarry Smith PetscErrorCode MatSetUpMultiply_MPISBAIJ_2comm(Mat mat)
14640781036SHong Zhang {
1472896befeSSatish Balay   Mat_MPISBAIJ   *baij = (Mat_MPISBAIJ*)mat->data;
148a30f8f8cSSatish Balay   Mat_SeqBAIJ    *B    = (Mat_SeqBAIJ*)(baij->B->data);
1496849ba73SBarry Smith   PetscErrorCode ierr;
15013f74950SBarry Smith   PetscInt       i,j,*aj = B->j,ec = 0,*garray;
151d0f46423SBarry Smith   PetscInt       bs = mat->rmap->bs,*stmp;
152a30f8f8cSSatish Balay   IS             from,to;
153a30f8f8cSSatish Balay   Vec            gvec;
154a30f8f8cSSatish Balay #if defined(PETSC_USE_CTABLE)
155a30f8f8cSSatish Balay   PetscTable         gid1_lid1;
156a30f8f8cSSatish Balay   PetscTablePosition tpos;
15713f74950SBarry Smith   PetscInt           gid,lid;
1586f531f54SSatish Balay #else
15913f74950SBarry Smith   PetscInt Nbs = baij->Nbs,*indices;
160a30f8f8cSSatish Balay #endif
161a30f8f8cSSatish Balay 
162a30f8f8cSSatish Balay   PetscFunctionBegin;
163a30f8f8cSSatish Balay #if defined(PETSC_USE_CTABLE)
164a30f8f8cSSatish Balay   /* use a table - Mark Adams */
165e23dfa41SBarry Smith   PetscTableCreate(B->mbs,baij->Nbs+1,&gid1_lid1);
166a30f8f8cSSatish Balay   for (i=0; i<B->mbs; i++) {
167a30f8f8cSSatish Balay     for (j=0; j<B->ilen[i]; j++) {
16813f74950SBarry Smith       PetscInt data,gid1 = aj[B->i[i]+j] + 1;
169a30f8f8cSSatish Balay       ierr = PetscTableFind(gid1_lid1,gid1,&data);CHKERRQ(ierr);
170a30f8f8cSSatish Balay       if (!data) {
171a30f8f8cSSatish Balay         /* one based table */
1723861aac3SJed Brown         ierr = PetscTableAdd(gid1_lid1,gid1,++ec,INSERT_VALUES);CHKERRQ(ierr);
173a30f8f8cSSatish Balay       }
174a30f8f8cSSatish Balay     }
175a30f8f8cSSatish Balay   }
176a30f8f8cSSatish Balay   /* form array of columns we need */
177785e854fSJed Brown   ierr = PetscMalloc1(ec,&garray);CHKERRQ(ierr);
178a30f8f8cSSatish Balay   ierr = PetscTableGetHeadPosition(gid1_lid1,&tpos);CHKERRQ(ierr);
179a30f8f8cSSatish Balay   while (tpos) {
180a30f8f8cSSatish Balay     ierr = PetscTableGetNext(gid1_lid1,&tpos,&gid,&lid);CHKERRQ(ierr);
181a30f8f8cSSatish Balay     gid--; lid--;
182a30f8f8cSSatish Balay     garray[lid] = gid;
183a30f8f8cSSatish Balay   }
184a30f8f8cSSatish Balay   ierr = PetscSortInt(ec,garray);CHKERRQ(ierr);
185a30f8f8cSSatish Balay   ierr = PetscTableRemoveAll(gid1_lid1);CHKERRQ(ierr);
186a30f8f8cSSatish Balay   for (i=0; i<ec; i++) {
1873861aac3SJed Brown     ierr = PetscTableAdd(gid1_lid1,garray[i]+1,i+1,INSERT_VALUES);CHKERRQ(ierr);
188a30f8f8cSSatish Balay   }
189a30f8f8cSSatish Balay   /* compact out the extra columns in B */
190a30f8f8cSSatish Balay   for (i=0; i<B->mbs; i++) {
191a30f8f8cSSatish Balay     for (j=0; j<B->ilen[i]; j++) {
19213f74950SBarry Smith       PetscInt gid1 = aj[B->i[i] + j] + 1;
193a30f8f8cSSatish Balay       ierr = PetscTableFind(gid1_lid1,gid1,&lid);CHKERRQ(ierr);
194a30f8f8cSSatish Balay       lid--;
195a30f8f8cSSatish Balay       aj[B->i[i]+j] = lid;
196a30f8f8cSSatish Balay     }
197a30f8f8cSSatish Balay   }
198a30f8f8cSSatish Balay   B->nbs = ec;
19926fbe8dcSKarl Rupp 
200d0f46423SBarry Smith   baij->B->cmap->n = baij->B->cmap->N = ec*mat->rmap->bs;
20126fbe8dcSKarl Rupp 
20226283091SBarry Smith   ierr = PetscLayoutSetUp((baij->B->cmap));CHKERRQ(ierr);
2036bc0bbbfSBarry Smith   ierr = PetscTableDestroy(&gid1_lid1);CHKERRQ(ierr);
204a30f8f8cSSatish Balay #else
205a30f8f8cSSatish Balay   /* For the first stab we make an array as long as the number of columns */
206a30f8f8cSSatish Balay   /* mark those columns that are in baij->B */
2071795a4d1SJed Brown   ierr = PetscCalloc1(Nbs,&indices);CHKERRQ(ierr);
208a30f8f8cSSatish Balay   for (i=0; i<B->mbs; i++) {
209a30f8f8cSSatish Balay     for (j=0; j<B->ilen[i]; j++) {
210a30f8f8cSSatish Balay       if (!indices[aj[B->i[i] + j]]) ec++;
211a30f8f8cSSatish Balay       indices[aj[B->i[i] + j]] = 1;
212a30f8f8cSSatish Balay     }
213a30f8f8cSSatish Balay   }
214a30f8f8cSSatish Balay 
215a30f8f8cSSatish Balay   /* form array of columns we need */
216785e854fSJed Brown   ierr = PetscMalloc1(ec,&garray);CHKERRQ(ierr);
217f3ef73ceSBarry Smith   ec = 0;
218f3ef73ceSBarry Smith   for (i=0; i<Nbs; i++) {
219f3ef73ceSBarry Smith     if (indices[i]) {
220f3ef73ceSBarry Smith       garray[ec++] = i;
221f3ef73ceSBarry Smith     }
222f3ef73ceSBarry Smith   }
223a30f8f8cSSatish Balay 
224a30f8f8cSSatish Balay   /* make indices now point into garray */
22526fbe8dcSKarl Rupp   for (i=0; i<ec; i++) indices[garray[i]] = i;
226a30f8f8cSSatish Balay 
227a30f8f8cSSatish Balay   /* compact out the extra columns in B */
228a30f8f8cSSatish Balay   for (i=0; i<B->mbs; i++) {
229a30f8f8cSSatish Balay     for (j=0; j<B->ilen[i]; j++) {
230a30f8f8cSSatish Balay       aj[B->i[i] + j] = indices[aj[B->i[i] + j]];
231a30f8f8cSSatish Balay     }
232a30f8f8cSSatish Balay   }
233a30f8f8cSSatish Balay   B->nbs = ec;
23426fbe8dcSKarl Rupp 
235d0f46423SBarry Smith   baij->B->cmap->n = ec*mat->rmap->bs;
23626fbe8dcSKarl Rupp 
237a30f8f8cSSatish Balay   ierr = PetscFree(indices);CHKERRQ(ierr);
238a30f8f8cSSatish Balay #endif
239633e10c7SHong Zhang 
240a30f8f8cSSatish Balay   /* create local vector that is used to scatter into */
241a30f8f8cSSatish Balay   ierr = VecCreateSeq(PETSC_COMM_SELF,ec*bs,&baij->lvec);CHKERRQ(ierr);
242a30f8f8cSSatish Balay 
243a30f8f8cSSatish Balay   /* create two temporary index sets for building scatter-gather */
244deff0451SBarry Smith   ierr = ISCreateBlock(PETSC_COMM_SELF,bs,ec,garray,PETSC_COPY_VALUES,&from);CHKERRQ(ierr);
245a30f8f8cSSatish Balay 
246785e854fSJed Brown   ierr = PetscMalloc1(ec,&stmp);CHKERRQ(ierr);
24726fbe8dcSKarl Rupp   for (i=0; i<ec; i++) stmp[i] = i;
248deff0451SBarry Smith   ierr = ISCreateBlock(PETSC_COMM_SELF,bs,ec,stmp,PETSC_OWN_POINTER,&to);CHKERRQ(ierr);
249a30f8f8cSSatish Balay 
250a30f8f8cSSatish Balay   /* create temporary global vector to generate scatter context */
251ce94432eSBarry Smith   ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)mat),1,mat->cmap->n,mat->cmap->N,NULL,&gvec);CHKERRQ(ierr);
2529448b7f1SJunchao Zhang   ierr = VecScatterCreate(gvec,from,baij->lvec,to,&baij->Mvctx);CHKERRQ(ierr);
2536bf464f9SBarry Smith   ierr = VecDestroy(&gvec);CHKERRQ(ierr);
254a30f8f8cSSatish Balay 
2553bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)baij->Mvctx);CHKERRQ(ierr);
2563bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)baij->lvec);CHKERRQ(ierr);
2573bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)from);CHKERRQ(ierr);
2583bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)to);CHKERRQ(ierr);
25926fbe8dcSKarl Rupp 
260a30f8f8cSSatish Balay   baij->garray = garray;
26126fbe8dcSKarl Rupp 
2623bb1ff40SBarry Smith   ierr = PetscLogObjectMemory((PetscObject)mat,(ec+1)*sizeof(PetscInt));CHKERRQ(ierr);
2636bf464f9SBarry Smith   ierr = ISDestroy(&from);CHKERRQ(ierr);
2646bf464f9SBarry Smith   ierr = ISDestroy(&to);CHKERRQ(ierr);
265a30f8f8cSSatish Balay   PetscFunctionReturn(0);
266a30f8f8cSSatish Balay }
267a30f8f8cSSatish Balay 
268a30f8f8cSSatish Balay /*
26901b2bd88SHong Zhang      Takes the local part of an already assembled MPISBAIJ matrix
270a30f8f8cSSatish Balay    and disassembles it. This is to allow new nonzeros into the matrix
271a30f8f8cSSatish Balay    that require more communication in the matrix vector multiply.
272a30f8f8cSSatish Balay    Thus certain data-structures must be rebuilt.
273a30f8f8cSSatish Balay 
274a30f8f8cSSatish Balay    Kind of slow! But that's what application programmers get when
275a30f8f8cSSatish Balay    they are sloppy.
276a30f8f8cSSatish Balay */
277ab9863d7SBarry Smith PetscErrorCode MatDisAssemble_MPISBAIJ(Mat A)
278a30f8f8cSSatish Balay {
2792896befeSSatish Balay   Mat_MPISBAIJ   *baij  = (Mat_MPISBAIJ*)A->data;
280a30f8f8cSSatish Balay   Mat            B      = baij->B,Bnew;
281a30f8f8cSSatish Balay   Mat_SeqBAIJ    *Bbaij = (Mat_SeqBAIJ*)B->data;
282dfbe8321SBarry Smith   PetscErrorCode ierr;
283d0f46423SBarry Smith   PetscInt       i,j,mbs=Bbaij->mbs,n = A->cmap->N,col,*garray=baij->garray;
284d0f46423SBarry Smith   PetscInt       k,bs=A->rmap->bs,bs2=baij->bs2,*rvals,*nz,ec,m=A->rmap->n;
285a30f8f8cSSatish Balay   MatScalar      *a = Bbaij->a;
28687828ca2SBarry Smith   PetscScalar    *atmp;
287ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL_MAT_SINGLE)
28813f74950SBarry Smith   PetscInt l;
289a30f8f8cSSatish Balay #endif
290a30f8f8cSSatish Balay 
291a30f8f8cSSatish Balay   PetscFunctionBegin;
292ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL_MAT_SINGLE)
293785e854fSJed Brown   ierr = PetscMalloc1(A->rmap->bs,&atmp);
29482502324SSatish Balay #endif
295a30f8f8cSSatish Balay   /* free stuff related to matrix-vec multiply */
296b0a32e0cSBarry Smith   ierr = VecGetSize(baij->lvec,&ec);CHKERRQ(ierr); /* needed for PetscLogObjectMemory below */
2976bf464f9SBarry Smith   ierr = VecDestroy(&baij->lvec);CHKERRQ(ierr);
2986bf464f9SBarry Smith   ierr = VecScatterDestroy(&baij->Mvctx);CHKERRQ(ierr);
29901b2bd88SHong Zhang 
3006bf464f9SBarry Smith   ierr = VecDestroy(&baij->slvec0);CHKERRQ(ierr);
3016bf464f9SBarry Smith   ierr = VecDestroy(&baij->slvec0b);CHKERRQ(ierr);
3026bf464f9SBarry Smith   ierr = VecDestroy(&baij->slvec1);CHKERRQ(ierr);
3036bf464f9SBarry Smith   ierr = VecDestroy(&baij->slvec1a);CHKERRQ(ierr);
3046bf464f9SBarry Smith   ierr = VecDestroy(&baij->slvec1b);CHKERRQ(ierr);
30501b2bd88SHong Zhang 
306a30f8f8cSSatish Balay   if (baij->colmap) {
307a30f8f8cSSatish Balay #if defined(PETSC_USE_CTABLE)
3086bc0bbbfSBarry Smith     ierr = PetscTableDestroy(&baij->colmap);CHKERRQ(ierr);
309a30f8f8cSSatish Balay #else
310a30f8f8cSSatish Balay     ierr = PetscFree(baij->colmap);CHKERRQ(ierr);
3113bb1ff40SBarry Smith     ierr = PetscLogObjectMemory((PetscObject)A,-Bbaij->nbs*sizeof(PetscInt));CHKERRQ(ierr);
312a30f8f8cSSatish Balay #endif
313a30f8f8cSSatish Balay   }
314a30f8f8cSSatish Balay 
315a30f8f8cSSatish Balay   /* make sure that B is assembled so we can access its values */
316a30f8f8cSSatish Balay   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
317a30f8f8cSSatish Balay   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
318a30f8f8cSSatish Balay 
319a30f8f8cSSatish Balay   /* invent new B and copy stuff over */
320785e854fSJed Brown   ierr = PetscMalloc1(mbs,&nz);CHKERRQ(ierr);
321a30f8f8cSSatish Balay   for (i=0; i<mbs; i++) {
322a30f8f8cSSatish Balay     nz[i] = Bbaij->i[i+1]-Bbaij->i[i];
323a30f8f8cSSatish Balay   }
324f69a0ea3SMatthew Knepley   ierr = MatCreate(PETSC_COMM_SELF,&Bnew);CHKERRQ(ierr);
325f69a0ea3SMatthew Knepley   ierr = MatSetSizes(Bnew,m,n,m,n);CHKERRQ(ierr);
3267adad957SLisandro Dalcin   ierr = MatSetType(Bnew,((PetscObject)B)->type_name);CHKERRQ(ierr);
327d0f46423SBarry Smith   ierr = MatSeqBAIJSetPreallocation(Bnew,B->rmap->bs,0,nz);CHKERRQ(ierr);
328a30f8f8cSSatish Balay   ierr = PetscFree(nz);CHKERRQ(ierr);
329a30f8f8cSSatish Balay 
330b38c15b3SStefano Zampini   if (Bbaij->nonew >= 0) { /* Inherit insertion error options (if positive). */
331b38c15b3SStefano Zampini     ((Mat_SeqSBAIJ*)Bnew->data)->nonew = Bbaij->nonew;
332b38c15b3SStefano Zampini   }
333b38c15b3SStefano Zampini 
33477341eacSDmitry Karpeev   /*
33577341eacSDmitry Karpeev    Ensure that B's nonzerostate is monotonically increasing.
33677341eacSDmitry Karpeev    Or should this follow the MatSetValues() loop to preserve B's nonzerstate across a MatDisAssemble() call?
33777341eacSDmitry Karpeev    */
33877341eacSDmitry Karpeev   Bnew->nonzerostate = B->nonzerostate;
339785e854fSJed Brown   ierr = PetscMalloc1(bs,&rvals);CHKERRQ(ierr);
340a30f8f8cSSatish Balay   for (i=0; i<mbs; i++) {
341a30f8f8cSSatish Balay     rvals[0] = bs*i;
34226fbe8dcSKarl Rupp     for (j=1; j<bs; j++) rvals[j] = rvals[j-1] + 1;
343a30f8f8cSSatish Balay     for (j=Bbaij->i[i]; j<Bbaij->i[i+1]; j++) {
344a30f8f8cSSatish Balay       col = garray[Bbaij->j[j]]*bs;
345a30f8f8cSSatish Balay       for (k=0; k<bs; k++) {
346ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL_MAT_SINGLE)
347a30f8f8cSSatish Balay         for (l=0; l<bs; l++) atmp[l] = a[j*bs2+l];
348a30f8f8cSSatish Balay #else
34971730473SSatish Balay         atmp = a+j*bs2 + k*bs;
350a30f8f8cSSatish Balay #endif
351c8407628SSatish Balay         ierr = MatSetValues_SeqSBAIJ(Bnew,bs,rvals,1,&col,atmp,B->insertmode);CHKERRQ(ierr);
352a30f8f8cSSatish Balay         col++;
353a30f8f8cSSatish Balay       }
354a30f8f8cSSatish Balay     }
355a30f8f8cSSatish Balay   }
356ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL_MAT_SINGLE)
357a30f8f8cSSatish Balay   ierr = PetscFree(atmp);CHKERRQ(ierr);
358a30f8f8cSSatish Balay #endif
359a30f8f8cSSatish Balay   ierr = PetscFree(baij->garray);CHKERRQ(ierr);
36026fbe8dcSKarl Rupp 
361a30f8f8cSSatish Balay   baij->garray = 0;
36226fbe8dcSKarl Rupp 
363a30f8f8cSSatish Balay   ierr = PetscFree(rvals);CHKERRQ(ierr);
3643bb1ff40SBarry Smith   ierr = PetscLogObjectMemory((PetscObject)A,-ec*sizeof(PetscInt));CHKERRQ(ierr);
3656bf464f9SBarry Smith   ierr = MatDestroy(&B);CHKERRQ(ierr);
3663bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)A,(PetscObject)Bnew);CHKERRQ(ierr);
36726fbe8dcSKarl Rupp 
368a30f8f8cSSatish Balay   baij->B = Bnew;
36926fbe8dcSKarl Rupp 
370a30f8f8cSSatish Balay   A->was_assembled = PETSC_FALSE;
371a30f8f8cSSatish Balay   PetscFunctionReturn(0);
372a30f8f8cSSatish Balay }
373