xref: /petsc/src/mat/impls/sbaij/mpi/mmsbaij.c (revision 6bc0bbbf38e6474be12d49c61632c43b1eeb2223)
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 
709573ac7SBarry Smith extern PetscErrorCode MatSetValues_SeqSBAIJ(Mat,PetscInt,const PetscInt [],PetscInt,const PetscInt [],const PetscScalar [],InsertMode);
87cff1425SSatish Balay 
9a30f8f8cSSatish Balay 
104a2ae208SSatish Balay #undef __FUNCT__
114a2ae208SSatish Balay #define __FUNCT__ "MatSetUpMultiply_MPISBAIJ"
12dfbe8321SBarry Smith PetscErrorCode MatSetUpMultiply_MPISBAIJ(Mat mat)
13a30f8f8cSSatish Balay {
1440781036SHong Zhang   Mat_MPISBAIJ   *sbaij = (Mat_MPISBAIJ*)mat->data;
1540781036SHong Zhang   Mat_SeqBAIJ    *B = (Mat_SeqBAIJ*)(sbaij->B->data);
166849ba73SBarry Smith   PetscErrorCode ierr;
1713f74950SBarry Smith   PetscInt       Nbs = sbaij->Nbs,i,j,*indices,*aj = B->j,ec = 0,*garray,*sgarray;
18d0f46423SBarry Smith   PetscInt       bs = mat->rmap->bs,*stmp,mbs=sbaij->mbs, vec_size,nt;
1940781036SHong Zhang   IS             from,to;
2040781036SHong Zhang   Vec            gvec;
2113f74950SBarry Smith   PetscMPIInt    rank=sbaij->rank,lsize,size=sbaij->size;
22899cda47SBarry Smith   PetscInt       *owners=sbaij->rangebs,*sowners,*ec_owner,k;
2340781036SHong Zhang   PetscScalar    *ptr;
2440781036SHong Zhang 
2540781036SHong Zhang   PetscFunctionBegin;
266bf464f9SBarry Smith   ierr = VecScatterDestroy(&sbaij->sMvctx);CHKERRQ(ierr);
2740781036SHong Zhang 
2840781036SHong Zhang   /* For the first stab we make an array as long as the number of columns */
2940781036SHong Zhang   /* mark those columns that are in sbaij->B */
3074ed9c26SBarry Smith   ierr = PetscMalloc(Nbs*sizeof(PetscInt),&indices);CHKERRQ(ierr);
3113f74950SBarry Smith   ierr = PetscMemzero(indices,Nbs*sizeof(PetscInt));CHKERRQ(ierr);
3240781036SHong Zhang   for (i=0; i<mbs; i++) {
3340781036SHong Zhang     for (j=0; j<B->ilen[i]; j++) {
3440781036SHong Zhang       if (!indices[aj[B->i[i] + j]]) ec++;
3540781036SHong Zhang       indices[aj[B->i[i] + j] ] = 1;
3640781036SHong Zhang     }
3740781036SHong Zhang   }
3840781036SHong Zhang 
3940781036SHong Zhang   /* form arrays of columns we need */
4074ed9c26SBarry Smith   ierr = PetscMalloc(ec*sizeof(PetscInt),&garray);CHKERRQ(ierr);
4174ed9c26SBarry Smith   ierr = PetscMalloc2(2*ec,PetscInt,&sgarray,ec,PetscInt,&ec_owner);CHKERRQ(ierr);
4240781036SHong Zhang 
4340781036SHong Zhang   ec = 0;
4440781036SHong Zhang   for (j=0; j<size; j++){
4540781036SHong Zhang     for (i=owners[j]; i<owners[j+1]; i++){
4640781036SHong Zhang       if (indices[i]) {
4740781036SHong Zhang         garray[ec]   = i;
4840781036SHong Zhang         ec_owner[ec] = j;
4940781036SHong Zhang         ec++;
5040781036SHong Zhang       }
5140781036SHong Zhang     }
5240781036SHong Zhang   }
5340781036SHong Zhang 
5440781036SHong Zhang   /* make indices now point into garray */
55b424e231SHong Zhang   for (i=0; i<ec; i++) indices[garray[i]] = i;
5640781036SHong Zhang 
5740781036SHong Zhang   /* compact out the extra columns in B */
5840781036SHong Zhang   for (i=0; i<mbs; i++) {
5940781036SHong Zhang     for (j=0; j<B->ilen[i]; j++) aj[B->i[i] + j] = indices[aj[B->i[i] + j]];
6040781036SHong Zhang   }
6140781036SHong Zhang   B->nbs      = ec;
62d0f46423SBarry Smith   sbaij->B->cmap->n = sbaij->B->cmap->N = ec*mat->rmap->bs;
6326283091SBarry Smith   ierr = PetscLayoutSetUp((sbaij->B->cmap));CHKERRQ(ierr);
6440781036SHong Zhang   ierr = PetscFree(indices);CHKERRQ(ierr);
6540781036SHong Zhang 
6640781036SHong Zhang   /* create local vector that is used to scatter into */
6740781036SHong Zhang   ierr = VecCreateSeq(PETSC_COMM_SELF,ec*bs,&sbaij->lvec);CHKERRQ(ierr);
6840781036SHong Zhang 
6940781036SHong Zhang   /* create two temporary index sets for building scatter-gather */
7074ed9c26SBarry Smith   ierr = PetscMalloc(2*ec*sizeof(PetscInt),&stmp);CHKERRQ(ierr);
71deff0451SBarry Smith   ierr = ISCreateBlock(PETSC_COMM_SELF,bs,ec,garray,PETSC_COPY_VALUES,&from);CHKERRQ(ierr);
72e82e9f6bSBarry Smith   for (i=0; i<ec; i++) { stmp[i] = i; }
73deff0451SBarry Smith   ierr = ISCreateBlock(PETSC_COMM_SELF,bs,ec,stmp,PETSC_COPY_VALUES,&to);CHKERRQ(ierr);
7440781036SHong Zhang 
7593d1592dSHong Zhang   /* generate the scatter context
7693d1592dSHong Zhang      -- Mvctx and lvec are not used by MatMult_MPISBAIJ(), but usefule for some applications */
77d0f46423SBarry Smith   ierr = VecCreateMPIWithArray(((PetscObject)mat)->comm,mat->cmap->n,mat->cmap->N,PETSC_NULL,&gvec);CHKERRQ(ierr);
7840781036SHong Zhang   ierr = VecScatterCreate(gvec,from,sbaij->lvec,to,&sbaij->Mvctx);CHKERRQ(ierr);
796bf464f9SBarry Smith   ierr = VecDestroy(&gvec);CHKERRQ(ierr);
8040781036SHong Zhang 
8140781036SHong Zhang   sbaij->garray = garray;
8252e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,sbaij->Mvctx);CHKERRQ(ierr);
8352e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,sbaij->lvec);CHKERRQ(ierr);
8452e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,from);CHKERRQ(ierr);
8552e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,to);CHKERRQ(ierr);
8640781036SHong Zhang 
876bf464f9SBarry Smith   ierr = ISDestroy(&from);CHKERRQ(ierr);
886bf464f9SBarry Smith   ierr = ISDestroy(&to);CHKERRQ(ierr);
8940781036SHong Zhang 
9040781036SHong Zhang   /* create parallel vector that is used by SBAIJ matrix to scatter from/into */
9140781036SHong Zhang   lsize = (mbs + ec)*bs;
927adad957SLisandro Dalcin   ierr = VecCreateMPI(((PetscObject)mat)->comm,lsize,PETSC_DETERMINE,&sbaij->slvec0);CHKERRQ(ierr);
9340781036SHong Zhang   ierr = VecDuplicate(sbaij->slvec0,&sbaij->slvec1);CHKERRQ(ierr);
9440781036SHong Zhang   ierr = VecGetSize(sbaij->slvec0,&vec_size);CHKERRQ(ierr);
9540781036SHong Zhang 
96d0f46423SBarry Smith   sowners = sbaij->slvec0->map->range;
9740781036SHong Zhang 
9840781036SHong Zhang   /* x index in the IS sfrom */
9940781036SHong Zhang   for (i=0; i<ec; i++) {
10040781036SHong Zhang     j = ec_owner[i];
10140781036SHong Zhang     sgarray[i]  = garray[i] + (sowners[j]/bs - owners[j]);
10240781036SHong Zhang   }
10340781036SHong Zhang   /* b index in the IS sfrom */
10440781036SHong Zhang   k = sowners[rank]/bs + mbs;
10540781036SHong Zhang   for (i=ec,j=0; i< 2*ec; i++,j++) sgarray[i] = k + j;
106deff0451SBarry Smith   ierr = ISCreateBlock(PETSC_COMM_SELF,bs,2*ec,sgarray,PETSC_COPY_VALUES,&from);CHKERRQ(ierr);
10740781036SHong Zhang 
10840781036SHong Zhang   /* x index in the IS sto */
10940781036SHong Zhang   k = sowners[rank]/bs + mbs;
110e82e9f6bSBarry Smith   for (i=0; i<ec; i++) stmp[i] = (k + i);
11140781036SHong Zhang   /* b index in the IS sto */
112e82e9f6bSBarry Smith   for (i=ec; i<2*ec; i++) stmp[i] = sgarray[i-ec];
11340781036SHong Zhang 
114deff0451SBarry Smith   ierr = ISCreateBlock(PETSC_COMM_SELF,bs,2*ec,stmp,PETSC_COPY_VALUES,&to);CHKERRQ(ierr);
11540781036SHong Zhang 
11640781036SHong Zhang   ierr = VecScatterCreate(sbaij->slvec0,from,sbaij->slvec1,to,&sbaij->sMvctx);CHKERRQ(ierr);
11740781036SHong Zhang 
11840781036SHong Zhang   ierr = VecGetLocalSize(sbaij->slvec1,&nt);CHKERRQ(ierr);
1191ebc52fbSHong Zhang   ierr = VecGetArray(sbaij->slvec1,&ptr);CHKERRQ(ierr);
12040781036SHong Zhang   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,bs*mbs,ptr,&sbaij->slvec1a);CHKERRQ(ierr);
12140781036SHong Zhang   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,nt-bs*mbs,ptr+bs*mbs,&sbaij->slvec1b);CHKERRQ(ierr);
1221ebc52fbSHong Zhang   ierr = VecRestoreArray(sbaij->slvec1,&ptr);CHKERRQ(ierr);
12340781036SHong Zhang 
1241ebc52fbSHong Zhang   ierr = VecGetArray(sbaij->slvec0,&ptr);CHKERRQ(ierr);
12540781036SHong Zhang   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,nt-bs*mbs,ptr+bs*mbs,&sbaij->slvec0b);CHKERRQ(ierr);
1261ebc52fbSHong Zhang   ierr = VecRestoreArray(sbaij->slvec0,&ptr);CHKERRQ(ierr);
12740781036SHong Zhang 
12840781036SHong Zhang   ierr = PetscFree(stmp);CHKERRQ(ierr);
1297adad957SLisandro Dalcin   ierr = MPI_Barrier(((PetscObject)mat)->comm);CHKERRQ(ierr);
13040781036SHong Zhang 
13152e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,sbaij->sMvctx);CHKERRQ(ierr);
13252e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,sbaij->slvec0);CHKERRQ(ierr);
13352e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,sbaij->slvec1);CHKERRQ(ierr);
13452e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,sbaij->slvec0b);CHKERRQ(ierr);
13552e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,sbaij->slvec1a);CHKERRQ(ierr);
13652e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,sbaij->slvec1b);CHKERRQ(ierr);
13752e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,from);CHKERRQ(ierr);
13852e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,to);CHKERRQ(ierr);
13940781036SHong Zhang 
14052e6d16bSBarry Smith   ierr = PetscLogObjectMemory(mat,(ec+1)*sizeof(PetscInt));CHKERRQ(ierr);
1416bf464f9SBarry Smith   ierr = ISDestroy(&from);CHKERRQ(ierr);
1426bf464f9SBarry Smith   ierr = ISDestroy(&to);CHKERRQ(ierr);
14374ed9c26SBarry Smith   ierr = PetscFree2(sgarray,ec_owner);CHKERRQ(ierr);
14440781036SHong Zhang   PetscFunctionReturn(0);
14540781036SHong Zhang }
14640781036SHong Zhang 
14740781036SHong Zhang #undef __FUNCT__
14840781036SHong Zhang #define __FUNCT__ "MatSetUpMultiply_MPISBAIJ_2comm"
149dfbe8321SBarry Smith PetscErrorCode MatSetUpMultiply_MPISBAIJ_2comm(Mat mat)
15040781036SHong Zhang {
1512896befeSSatish Balay   Mat_MPISBAIJ       *baij = (Mat_MPISBAIJ*)mat->data;
152a30f8f8cSSatish Balay   Mat_SeqBAIJ        *B = (Mat_SeqBAIJ*)(baij->B->data);
1536849ba73SBarry Smith   PetscErrorCode     ierr;
15413f74950SBarry Smith   PetscInt           i,j,*aj = B->j,ec = 0,*garray;
155d0f46423SBarry Smith   PetscInt           bs = mat->rmap->bs,*stmp;
156a30f8f8cSSatish Balay   IS                 from,to;
157a30f8f8cSSatish Balay   Vec                gvec;
158a30f8f8cSSatish Balay #if defined (PETSC_USE_CTABLE)
159a30f8f8cSSatish Balay   PetscTable         gid1_lid1;
160a30f8f8cSSatish Balay   PetscTablePosition tpos;
16113f74950SBarry Smith   PetscInt           gid,lid;
1626f531f54SSatish Balay #else
16313f74950SBarry Smith   PetscInt           Nbs = baij->Nbs,*indices;
164a30f8f8cSSatish Balay #endif
165a30f8f8cSSatish Balay 
166a30f8f8cSSatish Balay   PetscFunctionBegin;
167a30f8f8cSSatish Balay #if defined (PETSC_USE_CTABLE)
168a30f8f8cSSatish Balay   /* use a table - Mark Adams */
169a30f8f8cSSatish Balay   PetscTableCreate(B->mbs,&gid1_lid1);
170a30f8f8cSSatish Balay   for (i=0; i<B->mbs; i++) {
171a30f8f8cSSatish Balay     for (j=0; j<B->ilen[i]; j++) {
17213f74950SBarry Smith       PetscInt data,gid1 = aj[B->i[i]+j] + 1;
173a30f8f8cSSatish Balay       ierr = PetscTableFind(gid1_lid1,gid1,&data);CHKERRQ(ierr);
174a30f8f8cSSatish Balay       if (!data) {
175a30f8f8cSSatish Balay         /* one based table */
176a30f8f8cSSatish Balay         ierr = PetscTableAdd(gid1_lid1,gid1,++ec);CHKERRQ(ierr);
177a30f8f8cSSatish Balay       }
178a30f8f8cSSatish Balay     }
179a30f8f8cSSatish Balay   }
180a30f8f8cSSatish Balay   /* form array of columns we need */
18174ed9c26SBarry Smith   ierr = PetscMalloc(ec*sizeof(PetscInt),&garray);CHKERRQ(ierr);
182a30f8f8cSSatish Balay   ierr = PetscTableGetHeadPosition(gid1_lid1,&tpos);CHKERRQ(ierr);
183a30f8f8cSSatish Balay   while (tpos) {
184a30f8f8cSSatish Balay     ierr = PetscTableGetNext(gid1_lid1,&tpos,&gid,&lid);CHKERRQ(ierr);
185a30f8f8cSSatish Balay     gid--; lid--;
186a30f8f8cSSatish Balay     garray[lid] = gid;
187a30f8f8cSSatish Balay   }
188a30f8f8cSSatish Balay   ierr = PetscSortInt(ec,garray);CHKERRQ(ierr);
189a30f8f8cSSatish Balay   ierr = PetscTableRemoveAll(gid1_lid1);CHKERRQ(ierr);
190a30f8f8cSSatish Balay   for (i=0; i<ec; i++) {
191a30f8f8cSSatish Balay     ierr = PetscTableAdd(gid1_lid1,garray[i]+1,i+1);CHKERRQ(ierr);
192a30f8f8cSSatish Balay   }
193a30f8f8cSSatish Balay   /* compact out the extra columns in B */
194a30f8f8cSSatish Balay   for (i=0; i<B->mbs; i++) {
195a30f8f8cSSatish Balay     for (j=0; j<B->ilen[i]; j++) {
19613f74950SBarry Smith       PetscInt gid1 = aj[B->i[i] + j] + 1;
197a30f8f8cSSatish Balay       ierr = PetscTableFind(gid1_lid1,gid1,&lid);CHKERRQ(ierr);
198a30f8f8cSSatish Balay       lid --;
199a30f8f8cSSatish Balay       aj[B->i[i]+j] = lid;
200a30f8f8cSSatish Balay     }
201a30f8f8cSSatish Balay   }
202a30f8f8cSSatish Balay   B->nbs     = ec;
203d0f46423SBarry Smith   baij->B->cmap->n = baij->B->cmap->N = ec*mat->rmap->bs;
20426283091SBarry Smith   ierr = PetscLayoutSetUp((baij->B->cmap));CHKERRQ(ierr);
205*6bc0bbbfSBarry Smith   ierr = PetscTableDestroy(&gid1_lid1);CHKERRQ(ierr);
206a30f8f8cSSatish Balay #else
207a30f8f8cSSatish Balay   /* For the first stab we make an array as long as the number of columns */
208a30f8f8cSSatish Balay   /* mark those columns that are in baij->B */
20974ed9c26SBarry Smith   ierr = PetscMalloc(Nbs*sizeof(PetscInt),&indices);CHKERRQ(ierr);
21013f74950SBarry Smith   ierr = PetscMemzero(indices,Nbs*sizeof(PetscInt));CHKERRQ(ierr);
211a30f8f8cSSatish Balay   for (i=0; i<B->mbs; i++) {
212a30f8f8cSSatish Balay     for (j=0; j<B->ilen[i]; j++) {
213a30f8f8cSSatish Balay       if (!indices[aj[B->i[i] + j]]) ec++;
214a30f8f8cSSatish Balay       indices[aj[B->i[i] + j] ] = 1;
215a30f8f8cSSatish Balay     }
216a30f8f8cSSatish Balay   }
217a30f8f8cSSatish Balay 
218a30f8f8cSSatish Balay   /* form array of columns we need */
21974ed9c26SBarry Smith   ierr = PetscMalloc(ec*sizeof(PetscInt),&garray);CHKERRQ(ierr);
220f3ef73ceSBarry Smith   ec = 0;
221f3ef73ceSBarry Smith   for (i=0; i<Nbs; i++) {
222f3ef73ceSBarry Smith     if (indices[i]) {
223f3ef73ceSBarry Smith       garray[ec++] = i;
224f3ef73ceSBarry Smith     }
225f3ef73ceSBarry Smith   }
226a30f8f8cSSatish Balay 
227a30f8f8cSSatish Balay   /* make indices now point into garray */
228a30f8f8cSSatish Balay   for (i=0; i<ec; i++) {
229a30f8f8cSSatish Balay     indices[garray[i]] = i;
230a30f8f8cSSatish Balay   }
231a30f8f8cSSatish Balay 
232a30f8f8cSSatish Balay   /* compact out the extra columns in B */
233a30f8f8cSSatish Balay   for (i=0; i<B->mbs; i++) {
234a30f8f8cSSatish Balay     for (j=0; j<B->ilen[i]; j++) {
235a30f8f8cSSatish Balay       aj[B->i[i] + j] = indices[aj[B->i[i] + j]];
236a30f8f8cSSatish Balay     }
237a30f8f8cSSatish Balay   }
238a30f8f8cSSatish Balay   B->nbs       = ec;
239d0f46423SBarry Smith   baij->B->cmap->n   = ec*mat->rmap->bs;
240a30f8f8cSSatish Balay   ierr = PetscFree(indices);CHKERRQ(ierr);
241a30f8f8cSSatish Balay #endif
242633e10c7SHong Zhang 
243a30f8f8cSSatish Balay   /* create local vector that is used to scatter into */
244a30f8f8cSSatish Balay   ierr = VecCreateSeq(PETSC_COMM_SELF,ec*bs,&baij->lvec);CHKERRQ(ierr);
245a30f8f8cSSatish Balay 
246a30f8f8cSSatish Balay   /* create two temporary index sets for building scatter-gather */
247deff0451SBarry Smith   ierr = ISCreateBlock(PETSC_COMM_SELF,bs,ec,garray,PETSC_COPY_VALUES,&from);CHKERRQ(ierr);
248a30f8f8cSSatish Balay 
24974ed9c26SBarry Smith   ierr = PetscMalloc(ec*sizeof(PetscInt),&stmp);CHKERRQ(ierr);
250e82e9f6bSBarry Smith   for (i=0; i<ec; i++) { stmp[i] = i; }
251deff0451SBarry Smith   ierr = ISCreateBlock(PETSC_COMM_SELF,bs,ec,stmp,PETSC_OWN_POINTER,&to);CHKERRQ(ierr);
252a30f8f8cSSatish Balay 
253a30f8f8cSSatish Balay   /* create temporary global vector to generate scatter context */
254d0f46423SBarry Smith   ierr = VecCreateMPIWithArray(((PetscObject)mat)->comm,mat->cmap->n,mat->cmap->N,PETSC_NULL,&gvec);CHKERRQ(ierr);
255a30f8f8cSSatish Balay   ierr = VecScatterCreate(gvec,from,baij->lvec,to,&baij->Mvctx);CHKERRQ(ierr);
2566bf464f9SBarry Smith   ierr = VecDestroy(&gvec);CHKERRQ(ierr);
257a30f8f8cSSatish Balay 
25852e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,baij->Mvctx);CHKERRQ(ierr);
25952e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,baij->lvec);CHKERRQ(ierr);
26052e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,from);CHKERRQ(ierr);
26152e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,to);CHKERRQ(ierr);
262a30f8f8cSSatish Balay   baij->garray = garray;
26352e6d16bSBarry Smith   ierr = PetscLogObjectMemory(mat,(ec+1)*sizeof(PetscInt));CHKERRQ(ierr);
2646bf464f9SBarry Smith   ierr = ISDestroy(&from);CHKERRQ(ierr);
2656bf464f9SBarry Smith   ierr = ISDestroy(&to);CHKERRQ(ierr);
266a30f8f8cSSatish Balay   PetscFunctionReturn(0);
267a30f8f8cSSatish Balay }
268a30f8f8cSSatish Balay 
269a30f8f8cSSatish Balay /*
27001b2bd88SHong Zhang      Takes the local part of an already assembled MPISBAIJ matrix
271a30f8f8cSSatish Balay    and disassembles it. This is to allow new nonzeros into the matrix
272a30f8f8cSSatish Balay    that require more communication in the matrix vector multiply.
273a30f8f8cSSatish Balay    Thus certain data-structures must be rebuilt.
274a30f8f8cSSatish Balay 
275a30f8f8cSSatish Balay    Kind of slow! But that's what application programmers get when
276a30f8f8cSSatish Balay    they are sloppy.
277a30f8f8cSSatish Balay */
2784a2ae208SSatish Balay #undef __FUNCT__
2794a2ae208SSatish Balay #define __FUNCT__ "DisAssemble_MPISBAIJ"
280dfbe8321SBarry Smith PetscErrorCode DisAssemble_MPISBAIJ(Mat A)
281a30f8f8cSSatish Balay {
2822896befeSSatish Balay   Mat_MPISBAIJ   *baij = (Mat_MPISBAIJ*)A->data;
283a30f8f8cSSatish Balay   Mat            B = baij->B,Bnew;
284a30f8f8cSSatish Balay   Mat_SeqBAIJ    *Bbaij = (Mat_SeqBAIJ*)B->data;
285dfbe8321SBarry Smith   PetscErrorCode ierr;
286d0f46423SBarry Smith   PetscInt       i,j,mbs=Bbaij->mbs,n = A->cmap->N,col,*garray=baij->garray;
287d0f46423SBarry Smith   PetscInt       k,bs=A->rmap->bs,bs2=baij->bs2,*rvals,*nz,ec,m=A->rmap->n;
288a30f8f8cSSatish Balay   MatScalar      *a = Bbaij->a;
28987828ca2SBarry Smith   PetscScalar    *atmp;
290ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL_MAT_SINGLE)
29113f74950SBarry Smith   PetscInt       l;
292a30f8f8cSSatish Balay #endif
293a30f8f8cSSatish Balay 
294a30f8f8cSSatish Balay   PetscFunctionBegin;
295ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL_MAT_SINGLE)
296d0f46423SBarry Smith   ierr = PetscMalloc(A->rmap->bs*sizeof(PetscScalar),&atmp);
29782502324SSatish Balay #endif
298a30f8f8cSSatish Balay   /* free stuff related to matrix-vec multiply */
299b0a32e0cSBarry Smith   ierr = VecGetSize(baij->lvec,&ec);CHKERRQ(ierr); /* needed for PetscLogObjectMemory below */
3006bf464f9SBarry Smith   ierr = VecDestroy(&baij->lvec);CHKERRQ(ierr);
3016bf464f9SBarry Smith   ierr = VecScatterDestroy(&baij->Mvctx);CHKERRQ(ierr);
30201b2bd88SHong Zhang 
3036bf464f9SBarry Smith   ierr = VecDestroy(&baij->slvec0);CHKERRQ(ierr);
3046bf464f9SBarry Smith   ierr = VecDestroy(&baij->slvec0b);CHKERRQ(ierr);
3056bf464f9SBarry Smith   ierr = VecDestroy(&baij->slvec1);CHKERRQ(ierr);
3066bf464f9SBarry Smith   ierr = VecDestroy(&baij->slvec1a);CHKERRQ(ierr);
3076bf464f9SBarry Smith   ierr = VecDestroy(&baij->slvec1b);CHKERRQ(ierr);
30801b2bd88SHong Zhang 
309a30f8f8cSSatish Balay   if (baij->colmap) {
310a30f8f8cSSatish Balay #if defined (PETSC_USE_CTABLE)
311*6bc0bbbfSBarry Smith     ierr = PetscTableDestroy(&baij->colmap);CHKERRQ(ierr);
312a30f8f8cSSatish Balay #else
313a30f8f8cSSatish Balay     ierr = PetscFree(baij->colmap);CHKERRQ(ierr);
31452e6d16bSBarry Smith     ierr = PetscLogObjectMemory(A,-Bbaij->nbs*sizeof(PetscInt));CHKERRQ(ierr);
315a30f8f8cSSatish Balay #endif
316a30f8f8cSSatish Balay   }
317a30f8f8cSSatish Balay 
318a30f8f8cSSatish Balay   /* make sure that B is assembled so we can access its values */
319a30f8f8cSSatish Balay   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
320a30f8f8cSSatish Balay   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
321a30f8f8cSSatish Balay 
322a30f8f8cSSatish Balay   /* invent new B and copy stuff over */
32313f74950SBarry Smith   ierr = PetscMalloc(mbs*sizeof(PetscInt),&nz);CHKERRQ(ierr);
324a30f8f8cSSatish Balay   for (i=0; i<mbs; i++) {
325a30f8f8cSSatish Balay     nz[i] = Bbaij->i[i+1]-Bbaij->i[i];
326a30f8f8cSSatish Balay   }
327f69a0ea3SMatthew Knepley   ierr = MatCreate(PETSC_COMM_SELF,&Bnew);CHKERRQ(ierr);
328f69a0ea3SMatthew Knepley   ierr = MatSetSizes(Bnew,m,n,m,n);CHKERRQ(ierr);
3297adad957SLisandro Dalcin   ierr = MatSetType(Bnew,((PetscObject)B)->type_name);CHKERRQ(ierr);
330d0f46423SBarry Smith   ierr = MatSeqBAIJSetPreallocation(Bnew,B->rmap->bs,0,nz);CHKERRQ(ierr);
331a30f8f8cSSatish Balay   ierr = PetscFree(nz);CHKERRQ(ierr);
332a30f8f8cSSatish Balay 
33313f74950SBarry Smith   ierr = PetscMalloc(bs*sizeof(PetscInt),&rvals);CHKERRQ(ierr);
334a30f8f8cSSatish Balay   for (i=0; i<mbs; i++) {
335a30f8f8cSSatish Balay     rvals[0] = bs*i;
336a30f8f8cSSatish Balay     for (j=1; j<bs; j++) { rvals[j] = rvals[j-1] + 1; }
337a30f8f8cSSatish Balay     for (j=Bbaij->i[i]; j<Bbaij->i[i+1]; j++) {
338a30f8f8cSSatish Balay       col = garray[Bbaij->j[j]]*bs;
339a30f8f8cSSatish Balay       for (k=0; k<bs; k++) {
340ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL_MAT_SINGLE)
341a30f8f8cSSatish Balay         for (l=0; l<bs; l++) atmp[l] = a[j*bs2+l];
342a30f8f8cSSatish Balay #else
34371730473SSatish Balay         atmp = a+j*bs2 + k*bs;
344a30f8f8cSSatish Balay #endif
345c8407628SSatish Balay         ierr = MatSetValues_SeqSBAIJ(Bnew,bs,rvals,1,&col,atmp,B->insertmode);CHKERRQ(ierr);
346a30f8f8cSSatish Balay         col++;
347a30f8f8cSSatish Balay       }
348a30f8f8cSSatish Balay     }
349a30f8f8cSSatish Balay   }
350ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL_MAT_SINGLE)
351a30f8f8cSSatish Balay   ierr = PetscFree(atmp);CHKERRQ(ierr);
352a30f8f8cSSatish Balay #endif
353a30f8f8cSSatish Balay   ierr = PetscFree(baij->garray);CHKERRQ(ierr);
354a30f8f8cSSatish Balay   baij->garray = 0;
355a30f8f8cSSatish Balay   ierr = PetscFree(rvals);CHKERRQ(ierr);
35652e6d16bSBarry Smith   ierr = PetscLogObjectMemory(A,-ec*sizeof(PetscInt));CHKERRQ(ierr);
3576bf464f9SBarry Smith   ierr = MatDestroy(&B);CHKERRQ(ierr);
35852e6d16bSBarry Smith   ierr = PetscLogObjectParent(A,Bnew);CHKERRQ(ierr);
359a30f8f8cSSatish Balay   baij->B = Bnew;
360a30f8f8cSSatish Balay   A->was_assembled = PETSC_FALSE;
361a30f8f8cSSatish Balay   PetscFunctionReturn(0);
362a30f8f8cSSatish Balay }
363a30f8f8cSSatish Balay 
364a30f8f8cSSatish Balay 
365a30f8f8cSSatish Balay 
366