xref: /petsc/src/mat/impls/sbaij/mpi/mmsbaij.c (revision 6849ba73f22fecb8f92ef896a42e4e8bd4cd6965)
1a30f8f8cSSatish Balay 
2a30f8f8cSSatish Balay /*
3a30f8f8cSSatish Balay    Support for the parallel SBAIJ matrix vector multiply
4a30f8f8cSSatish Balay */
52896befeSSatish Balay #include "src/mat/impls/sbaij/mpi/mpisbaij.h"
67cff1425SSatish Balay 
7dfbe8321SBarry Smith EXTERN PetscErrorCode MatSetValues_SeqSBAIJ(Mat,int,const int [],int,const int [],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);
16*6849ba73SBarry Smith   PetscErrorCode ierr;
17*6849ba73SBarry Smith   int                Nbs = sbaij->Nbs,i,j,*indices,*aj = B->j,ec = 0,*garray,*sgarray;
1840781036SHong Zhang   int                bs = sbaij->bs,*stmp,mbs=sbaij->mbs, vec_size,nt;
1940781036SHong Zhang   IS                 from,to;
2040781036SHong Zhang   Vec                gvec;
2140781036SHong Zhang   int                rank=sbaij->rank,lsize,size=sbaij->size;
22b424e231SHong Zhang   int                *owners=sbaij->rowners,*sowners,*ec_owner,k;
2340781036SHong Zhang   PetscMap           vecmap;
2440781036SHong Zhang   PetscScalar        *ptr;
2540781036SHong Zhang 
2640781036SHong Zhang   PetscFunctionBegin;
2740781036SHong Zhang   if (sbaij->lvec) {
2840781036SHong Zhang     ierr = VecDestroy(sbaij->lvec);CHKERRQ(ierr);
2940781036SHong Zhang     sbaij->lvec = 0;
3040781036SHong Zhang   }
3140781036SHong Zhang   if (sbaij->Mvctx) {
3240781036SHong Zhang     ierr = VecScatterDestroy(sbaij->Mvctx);CHKERRQ(ierr);
3340781036SHong Zhang     sbaij->Mvctx = 0;
3440781036SHong Zhang   }
3540781036SHong Zhang 
3640781036SHong Zhang   /* For the first stab we make an array as long as the number of columns */
3740781036SHong Zhang   /* mark those columns that are in sbaij->B */
3840781036SHong Zhang   ierr = PetscMalloc((Nbs+1)*sizeof(int),&indices);CHKERRQ(ierr);
3940781036SHong Zhang   ierr = PetscMemzero(indices,Nbs*sizeof(int));CHKERRQ(ierr);
4040781036SHong Zhang   for (i=0; i<mbs; i++) {
4140781036SHong Zhang     for (j=0; j<B->ilen[i]; j++) {
4240781036SHong Zhang       if (!indices[aj[B->i[i] + j]]) ec++;
4340781036SHong Zhang       indices[aj[B->i[i] + j] ] = 1;
4440781036SHong Zhang     }
4540781036SHong Zhang   }
4640781036SHong Zhang 
4740781036SHong Zhang   /* form arrays of columns we need */
4840781036SHong Zhang   ierr = PetscMalloc((ec+1)*sizeof(int),&garray);CHKERRQ(ierr);
4940781036SHong Zhang   ierr = PetscMalloc((3*ec+1)*sizeof(int),&sgarray);CHKERRQ(ierr);
5040781036SHong Zhang   ec_owner = sgarray + 2*ec;
5140781036SHong Zhang 
5240781036SHong Zhang   ec = 0;
5340781036SHong Zhang   for (j=0; j<size; j++){
5440781036SHong Zhang     for (i=owners[j]; i<owners[j+1]; i++){
5540781036SHong Zhang       if (indices[i]) {
5640781036SHong Zhang         garray[ec]   = i;
5740781036SHong Zhang         ec_owner[ec] = j;
5840781036SHong Zhang         ec++;
5940781036SHong Zhang       }
6040781036SHong Zhang     }
6140781036SHong Zhang   }
6240781036SHong Zhang 
6340781036SHong Zhang   /* make indices now point into garray */
64b424e231SHong Zhang   for (i=0; i<ec; i++) indices[garray[i]] = i;
6540781036SHong Zhang 
6640781036SHong Zhang   /* compact out the extra columns in B */
6740781036SHong Zhang   for (i=0; i<mbs; i++) {
6840781036SHong Zhang     for (j=0; j<B->ilen[i]; j++) aj[B->i[i] + j] = indices[aj[B->i[i] + j]];
6940781036SHong Zhang   }
7040781036SHong Zhang   B->nbs      = ec;
7140781036SHong Zhang   sbaij->B->n = ec*B->bs;
7240781036SHong Zhang   ierr = PetscFree(indices);CHKERRQ(ierr);
7340781036SHong Zhang 
7440781036SHong Zhang   /* create local vector that is used to scatter into */
7540781036SHong Zhang   ierr = VecCreateSeq(PETSC_COMM_SELF,ec*bs,&sbaij->lvec);CHKERRQ(ierr);
7640781036SHong Zhang 
7740781036SHong Zhang   /* create two temporary index sets for building scatter-gather */
7840781036SHong Zhang   ierr = PetscMalloc((2*ec+1)*sizeof(int),&stmp);CHKERRQ(ierr);
7940781036SHong Zhang   for (i=0; i<ec; i++) stmp[i] = bs*garray[i];
8040781036SHong Zhang   ierr = ISCreateBlock(PETSC_COMM_SELF,bs,ec,stmp,&from);CHKERRQ(ierr);
8140781036SHong Zhang 
8240781036SHong Zhang   for (i=0; i<ec; i++) { stmp[i] = bs*i; }
8340781036SHong Zhang   ierr = ISCreateBlock(PETSC_COMM_SELF,bs,ec,stmp,&to);CHKERRQ(ierr);
8440781036SHong Zhang 
8593d1592dSHong Zhang   /* generate the scatter context
8693d1592dSHong Zhang      -- Mvctx and lvec are not used by MatMult_MPISBAIJ(), but usefule for some applications */
8740781036SHong Zhang   ierr = VecCreateMPI(mat->comm,mat->n,mat->N,&gvec);CHKERRQ(ierr);
8840781036SHong Zhang   ierr = VecScatterCreate(gvec,from,sbaij->lvec,to,&sbaij->Mvctx);CHKERRQ(ierr);
8940781036SHong Zhang   ierr = VecScatterPostRecvs(gvec,sbaij->lvec,INSERT_VALUES,SCATTER_FORWARD,sbaij->Mvctx);CHKERRQ(ierr);
9040781036SHong Zhang 
9140781036SHong Zhang   sbaij->garray = garray;
9240781036SHong Zhang   PetscLogObjectParent(mat,sbaij->Mvctx);
9340781036SHong Zhang   PetscLogObjectParent(mat,sbaij->lvec);
9440781036SHong Zhang   PetscLogObjectParent(mat,from);
9540781036SHong Zhang   PetscLogObjectParent(mat,to);
9640781036SHong Zhang 
9740781036SHong Zhang   ierr = ISDestroy(from);CHKERRQ(ierr);
9840781036SHong Zhang   ierr = ISDestroy(to);CHKERRQ(ierr);
9940781036SHong Zhang 
10040781036SHong Zhang   /* create parallel vector that is used by SBAIJ matrix to scatter from/into */
10140781036SHong Zhang   lsize = (mbs + ec)*bs;
10240781036SHong Zhang   ierr = VecCreateMPI(mat->comm,lsize,PETSC_DETERMINE,&sbaij->slvec0);CHKERRQ(ierr);
10340781036SHong Zhang   ierr = VecDuplicate(sbaij->slvec0,&sbaij->slvec1);CHKERRQ(ierr);
10440781036SHong Zhang   ierr = VecGetSize(sbaij->slvec0,&vec_size);CHKERRQ(ierr);
10540781036SHong Zhang 
10640781036SHong Zhang   ierr = VecGetPetscMap(sbaij->slvec0,&vecmap);CHKERRQ(ierr);
10740781036SHong Zhang   ierr = PetscMapGetGlobalRange(vecmap,&sowners);CHKERRQ(ierr);
10840781036SHong Zhang 
10940781036SHong Zhang   /* x index in the IS sfrom */
11040781036SHong Zhang   for (i=0; i<ec; i++) {
11140781036SHong Zhang     j = ec_owner[i];
11240781036SHong Zhang     sgarray[i]  = garray[i] + (sowners[j]/bs - owners[j]);
11340781036SHong Zhang   }
11440781036SHong Zhang   /* b index in the IS sfrom */
11540781036SHong Zhang   k = sowners[rank]/bs + mbs;
11640781036SHong Zhang   for (i=ec,j=0; i< 2*ec; i++,j++) sgarray[i] = k + j;
11740781036SHong Zhang 
11840781036SHong Zhang   for (i=0; i<2*ec; i++) stmp[i] = bs*sgarray[i];
11940781036SHong Zhang   ierr = ISCreateBlock(PETSC_COMM_SELF,bs,2*ec,stmp,&from);CHKERRQ(ierr);
12040781036SHong Zhang 
12140781036SHong Zhang   /* x index in the IS sto */
12240781036SHong Zhang   k = sowners[rank]/bs + mbs;
12340781036SHong Zhang   for (i=0; i<ec; i++) stmp[i] = bs*(k + i);
12440781036SHong Zhang   /* b index in the IS sto */
12540781036SHong Zhang   for (i=ec; i<2*ec; i++) stmp[i] = bs*sgarray[i-ec];
12640781036SHong Zhang 
12740781036SHong Zhang   ierr = ISCreateBlock(PETSC_COMM_SELF,bs,2*ec,stmp,&to);CHKERRQ(ierr);
12840781036SHong Zhang 
12940781036SHong Zhang   /* gnerate the SBAIJ scatter context */
13040781036SHong Zhang   ierr = VecScatterCreate(sbaij->slvec0,from,sbaij->slvec1,to,&sbaij->sMvctx);CHKERRQ(ierr);
13140781036SHong Zhang 
13240781036SHong Zhang    /*
13340781036SHong Zhang       Post the receives for the first matrix vector product. We sync-chronize after
13440781036SHong Zhang     this on the chance that the user immediately calls MatMult() after assemblying
13540781036SHong Zhang     the matrix.
13640781036SHong Zhang   */
13740781036SHong Zhang   ierr = VecScatterPostRecvs(sbaij->slvec0,sbaij->slvec1,INSERT_VALUES,SCATTER_FORWARD,sbaij->sMvctx);CHKERRQ(ierr);
13840781036SHong Zhang 
13940781036SHong Zhang   ierr = VecGetLocalSize(sbaij->slvec1,&nt);CHKERRQ(ierr);
1401ebc52fbSHong Zhang   ierr = VecGetArray(sbaij->slvec1,&ptr);CHKERRQ(ierr);
14140781036SHong Zhang   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,bs*mbs,ptr,&sbaij->slvec1a);CHKERRQ(ierr);
14240781036SHong Zhang   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,nt-bs*mbs,ptr+bs*mbs,&sbaij->slvec1b);CHKERRQ(ierr);
1431ebc52fbSHong Zhang   ierr = VecRestoreArray(sbaij->slvec1,&ptr);CHKERRQ(ierr);
14440781036SHong Zhang 
1451ebc52fbSHong Zhang   ierr = VecGetArray(sbaij->slvec0,&ptr);CHKERRQ(ierr);
14640781036SHong Zhang   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,nt-bs*mbs,ptr+bs*mbs,&sbaij->slvec0b);CHKERRQ(ierr);
1471ebc52fbSHong Zhang   ierr = VecRestoreArray(sbaij->slvec0,&ptr);CHKERRQ(ierr);
14840781036SHong Zhang 
14940781036SHong Zhang   ierr = PetscFree(stmp);CHKERRQ(ierr);
15040781036SHong Zhang   ierr = MPI_Barrier(mat->comm);CHKERRQ(ierr);
15140781036SHong Zhang 
15240781036SHong Zhang   PetscLogObjectParent(mat,sbaij->sMvctx);
15340781036SHong Zhang   PetscLogObjectParent(mat,sbaij->slvec0);
15440781036SHong Zhang   PetscLogObjectParent(mat,sbaij->slvec1);
15540781036SHong Zhang   PetscLogObjectParent(mat,sbaij->slvec0b);
15640781036SHong Zhang   PetscLogObjectParent(mat,sbaij->slvec1a);
15740781036SHong Zhang   PetscLogObjectParent(mat,sbaij->slvec1b);
15840781036SHong Zhang   PetscLogObjectParent(mat,from);
15940781036SHong Zhang   PetscLogObjectParent(mat,to);
16040781036SHong Zhang 
16140781036SHong Zhang   PetscLogObjectMemory(mat,(ec+1)*sizeof(int));
16240781036SHong Zhang   ierr = ISDestroy(from);CHKERRQ(ierr);
16340781036SHong Zhang   ierr = ISDestroy(to);CHKERRQ(ierr);
16440781036SHong Zhang   ierr = VecDestroy(gvec);CHKERRQ(ierr);
16540781036SHong Zhang   ierr = PetscFree(sgarray);CHKERRQ(ierr);
16640781036SHong Zhang 
16740781036SHong Zhang   PetscFunctionReturn(0);
16840781036SHong Zhang }
16940781036SHong Zhang 
17040781036SHong Zhang #undef __FUNCT__
17140781036SHong Zhang #define __FUNCT__ "MatSetUpMultiply_MPISBAIJ_2comm"
172dfbe8321SBarry Smith PetscErrorCode MatSetUpMultiply_MPISBAIJ_2comm(Mat mat)
17340781036SHong Zhang {
1742896befeSSatish Balay   Mat_MPISBAIJ       *baij = (Mat_MPISBAIJ*)mat->data;
175a30f8f8cSSatish Balay   Mat_SeqBAIJ        *B = (Mat_SeqBAIJ*)(baij->B->data);
176*6849ba73SBarry Smith   PetscErrorCode ierr;
177*6849ba73SBarry Smith   int                i,j,*aj = B->j,ec = 0,*garray;
17886e15631SSatish Balay   int                bs = baij->bs,*stmp;
179a30f8f8cSSatish Balay   IS                 from,to;
180a30f8f8cSSatish Balay   Vec                gvec;
181a30f8f8cSSatish Balay #if defined (PETSC_USE_CTABLE)
182a30f8f8cSSatish Balay   PetscTable         gid1_lid1;
183a30f8f8cSSatish Balay   PetscTablePosition tpos;
184a30f8f8cSSatish Balay   int                gid,lid;
1856f531f54SSatish Balay #else
1866f531f54SSatish Balay   int                Nbs = baij->Nbs,*indices;
187a30f8f8cSSatish Balay #endif
188a30f8f8cSSatish Balay 
189a30f8f8cSSatish Balay   PetscFunctionBegin;
190a30f8f8cSSatish Balay 
191a30f8f8cSSatish Balay #if defined (PETSC_USE_CTABLE)
192a30f8f8cSSatish Balay   /* use a table - Mark Adams */
193a30f8f8cSSatish Balay   PetscTableCreate(B->mbs,&gid1_lid1);
194a30f8f8cSSatish Balay   for (i=0; i<B->mbs; i++) {
195a30f8f8cSSatish Balay     for (j=0; j<B->ilen[i]; j++) {
196a30f8f8cSSatish Balay       int data,gid1 = aj[B->i[i]+j] + 1;
197a30f8f8cSSatish Balay       ierr = PetscTableFind(gid1_lid1,gid1,&data);CHKERRQ(ierr);
198a30f8f8cSSatish Balay       if (!data) {
199a30f8f8cSSatish Balay         /* one based table */
200a30f8f8cSSatish Balay         ierr = PetscTableAdd(gid1_lid1,gid1,++ec);CHKERRQ(ierr);
201a30f8f8cSSatish Balay       }
202a30f8f8cSSatish Balay     }
203a30f8f8cSSatish Balay   }
204a30f8f8cSSatish Balay   /* form array of columns we need */
205b0a32e0cSBarry Smith   ierr = PetscMalloc((ec+1)*sizeof(int),&garray);CHKERRQ(ierr);
206a30f8f8cSSatish Balay   ierr = PetscTableGetHeadPosition(gid1_lid1,&tpos);CHKERRQ(ierr);
207a30f8f8cSSatish Balay   while (tpos) {
208a30f8f8cSSatish Balay     ierr = PetscTableGetNext(gid1_lid1,&tpos,&gid,&lid);CHKERRQ(ierr);
209a30f8f8cSSatish Balay     gid--; lid--;
210a30f8f8cSSatish Balay     garray[lid] = gid;
211a30f8f8cSSatish Balay   }
212a30f8f8cSSatish Balay   ierr = PetscSortInt(ec,garray);CHKERRQ(ierr);
213a30f8f8cSSatish Balay   ierr = PetscTableRemoveAll(gid1_lid1);CHKERRQ(ierr);
214a30f8f8cSSatish Balay   for (i=0; i<ec; i++) {
215a30f8f8cSSatish Balay     ierr = PetscTableAdd(gid1_lid1,garray[i]+1,i+1);CHKERRQ(ierr);
216a30f8f8cSSatish Balay   }
217a30f8f8cSSatish Balay   /* compact out the extra columns in B */
218a30f8f8cSSatish Balay   for (i=0; i<B->mbs; i++) {
219a30f8f8cSSatish Balay     for (j=0; j<B->ilen[i]; j++) {
220a30f8f8cSSatish Balay       int gid1 = aj[B->i[i] + j] + 1;
221a30f8f8cSSatish Balay       ierr = PetscTableFind(gid1_lid1,gid1,&lid);CHKERRQ(ierr);
222a30f8f8cSSatish Balay       lid --;
223a30f8f8cSSatish Balay       aj[B->i[i]+j] = lid;
224a30f8f8cSSatish Balay     }
225a30f8f8cSSatish Balay   }
226a30f8f8cSSatish Balay   B->nbs     = ec;
227273d9f13SBarry Smith   baij->B->n = ec*B->bs;
228a30f8f8cSSatish Balay   ierr = PetscTableDelete(gid1_lid1);CHKERRQ(ierr);
229a30f8f8cSSatish Balay   /* Mark Adams */
230a30f8f8cSSatish Balay #else
231a30f8f8cSSatish Balay   /* For the first stab we make an array as long as the number of columns */
232a30f8f8cSSatish Balay   /* mark those columns that are in baij->B */
233b0a32e0cSBarry Smith   ierr = PetscMalloc((Nbs+1)*sizeof(int),&indices);CHKERRQ(ierr);
234a30f8f8cSSatish Balay   ierr = PetscMemzero(indices,Nbs*sizeof(int));CHKERRQ(ierr);
235a30f8f8cSSatish Balay   for (i=0; i<B->mbs; i++) {
236a30f8f8cSSatish Balay     for (j=0; j<B->ilen[i]; j++) {
237a30f8f8cSSatish Balay       if (!indices[aj[B->i[i] + j]]) ec++;
238a30f8f8cSSatish Balay       indices[aj[B->i[i] + j] ] = 1;
239a30f8f8cSSatish Balay     }
240a30f8f8cSSatish Balay   }
241a30f8f8cSSatish Balay 
242a30f8f8cSSatish Balay   /* form array of columns we need */
243b0a32e0cSBarry Smith   ierr = PetscMalloc((ec+1)*sizeof(int),&garray);CHKERRQ(ierr);
244f3ef73ceSBarry Smith   ec = 0;
245f3ef73ceSBarry Smith   for (i=0; i<Nbs; i++) {
246f3ef73ceSBarry Smith     if (indices[i]) {
247f3ef73ceSBarry Smith       garray[ec++] = i;
248f3ef73ceSBarry Smith     }
249f3ef73ceSBarry Smith   }
250a30f8f8cSSatish Balay 
251a30f8f8cSSatish Balay   /* make indices now point into garray */
252a30f8f8cSSatish Balay   for (i=0; i<ec; i++) {
253a30f8f8cSSatish Balay     indices[garray[i]] = i;
254a30f8f8cSSatish Balay   }
255a30f8f8cSSatish Balay 
256a30f8f8cSSatish Balay   /* compact out the extra columns in B */
257a30f8f8cSSatish Balay   for (i=0; i<B->mbs; i++) {
258a30f8f8cSSatish Balay     for (j=0; j<B->ilen[i]; j++) {
259a30f8f8cSSatish Balay       aj[B->i[i] + j] = indices[aj[B->i[i] + j]];
260a30f8f8cSSatish Balay     }
261a30f8f8cSSatish Balay   }
262a30f8f8cSSatish Balay   B->nbs       = ec;
263273d9f13SBarry Smith   baij->B->n   = ec*B->bs;
264a30f8f8cSSatish Balay   ierr = PetscFree(indices);CHKERRQ(ierr);
265a30f8f8cSSatish Balay #endif
266633e10c7SHong Zhang 
267a30f8f8cSSatish Balay   /* create local vector that is used to scatter into */
268a30f8f8cSSatish Balay   ierr = VecCreateSeq(PETSC_COMM_SELF,ec*bs,&baij->lvec);CHKERRQ(ierr);
269a30f8f8cSSatish Balay 
270a30f8f8cSSatish Balay   /* create two temporary index sets for building scatter-gather */
271eb7adc28SSatish Balay   for (i=0; i<ec; i++) {
272a30f8f8cSSatish Balay     garray[i] = bs*garray[i];
273a30f8f8cSSatish Balay   }
274a30f8f8cSSatish Balay   ierr = ISCreateBlock(PETSC_COMM_SELF,bs,ec,garray,&from);CHKERRQ(ierr);
27586e15631SSatish Balay   for (i=0; i<ec; i++) {
276a30f8f8cSSatish Balay     garray[i] = garray[i]/bs;
277a30f8f8cSSatish Balay   }
278a30f8f8cSSatish Balay 
279b0a32e0cSBarry Smith   ierr = PetscMalloc((ec+1)*sizeof(int),&stmp);CHKERRQ(ierr);
280a30f8f8cSSatish Balay   for (i=0; i<ec; i++) { stmp[i] = bs*i; }
281a30f8f8cSSatish Balay   ierr = ISCreateBlock(PETSC_COMM_SELF,bs,ec,stmp,&to);CHKERRQ(ierr);
282a30f8f8cSSatish Balay   ierr = PetscFree(stmp);CHKERRQ(ierr);
283a30f8f8cSSatish Balay 
284a30f8f8cSSatish Balay   /* create temporary global vector to generate scatter context */
285a30f8f8cSSatish Balay   /* this is inefficient, but otherwise we must do either
286a30f8f8cSSatish Balay      1) save garray until the first actual scatter when the vector is known or
287a30f8f8cSSatish Balay      2) have another way of generating a scatter context without a vector.*/
288273d9f13SBarry Smith   ierr = VecCreateMPI(mat->comm,mat->n,mat->N,&gvec);CHKERRQ(ierr);
289a30f8f8cSSatish Balay 
290a30f8f8cSSatish Balay   /* gnerate the scatter context */
291a30f8f8cSSatish Balay   ierr = VecScatterCreate(gvec,from,baij->lvec,to,&baij->Mvctx);CHKERRQ(ierr);
292a30f8f8cSSatish Balay 
293a30f8f8cSSatish Balay   /*
294a30f8f8cSSatish Balay       Post the receives for the first matrix vector product. We sync-chronize after
295a30f8f8cSSatish Balay     this on the chance that the user immediately calls MatMult() after assemblying
296a30f8f8cSSatish Balay     the matrix.
297a30f8f8cSSatish Balay   */
298a30f8f8cSSatish Balay   ierr = VecScatterPostRecvs(gvec,baij->lvec,INSERT_VALUES,SCATTER_FORWARD,baij->Mvctx);CHKERRQ(ierr);
299a30f8f8cSSatish Balay   ierr = MPI_Barrier(mat->comm);CHKERRQ(ierr);
300a30f8f8cSSatish Balay 
301b0a32e0cSBarry Smith   PetscLogObjectParent(mat,baij->Mvctx);
302b0a32e0cSBarry Smith   PetscLogObjectParent(mat,baij->lvec);
303b0a32e0cSBarry Smith   PetscLogObjectParent(mat,from);
304b0a32e0cSBarry Smith   PetscLogObjectParent(mat,to);
305a30f8f8cSSatish Balay   baij->garray = garray;
306b0a32e0cSBarry Smith   PetscLogObjectMemory(mat,(ec+1)*sizeof(int));
307a30f8f8cSSatish Balay   ierr = ISDestroy(from);CHKERRQ(ierr);
308a30f8f8cSSatish Balay   ierr = ISDestroy(to);CHKERRQ(ierr);
309a30f8f8cSSatish Balay   ierr = VecDestroy(gvec);CHKERRQ(ierr);
310a30f8f8cSSatish Balay   PetscFunctionReturn(0);
311a30f8f8cSSatish Balay }
312a30f8f8cSSatish Balay 
313a30f8f8cSSatish Balay 
314a30f8f8cSSatish Balay /*
315a30f8f8cSSatish Balay      Takes the local part of an already assembled MPIBAIJ matrix
316a30f8f8cSSatish Balay    and disassembles it. This is to allow new nonzeros into the matrix
317a30f8f8cSSatish Balay    that require more communication in the matrix vector multiply.
318a30f8f8cSSatish Balay    Thus certain data-structures must be rebuilt.
319a30f8f8cSSatish Balay 
320a30f8f8cSSatish Balay    Kind of slow! But that's what application programmers get when
321a30f8f8cSSatish Balay    they are sloppy.
322a30f8f8cSSatish Balay */
3234a2ae208SSatish Balay #undef __FUNCT__
3244a2ae208SSatish Balay #define __FUNCT__ "DisAssemble_MPISBAIJ"
325dfbe8321SBarry Smith PetscErrorCode DisAssemble_MPISBAIJ(Mat A)
326a30f8f8cSSatish Balay {
3272896befeSSatish Balay   Mat_MPISBAIJ   *baij = (Mat_MPISBAIJ*)A->data;
328a30f8f8cSSatish Balay   Mat           B = baij->B,Bnew;
329a30f8f8cSSatish Balay   Mat_SeqBAIJ   *Bbaij = (Mat_SeqBAIJ*)B->data;
330dfbe8321SBarry Smith   PetscErrorCode ierr;
331dfbe8321SBarry Smith   int i,j,mbs=Bbaij->mbs,n = A->N,col,*garray=baij->garray;
332273d9f13SBarry Smith   int           k,bs=baij->bs,bs2=baij->bs2,*rvals,*nz,ec,m=A->m;
333a30f8f8cSSatish Balay   MatScalar     *a = Bbaij->a;
33487828ca2SBarry Smith   PetscScalar   *atmp;
33582502324SSatish Balay #if defined(PETSC_USE_MAT_SINGLE)
33682502324SSatish Balay   int           l;
337a30f8f8cSSatish Balay #endif
338a30f8f8cSSatish Balay 
339a30f8f8cSSatish Balay   PetscFunctionBegin;
34082502324SSatish Balay #if defined(PETSC_USE_MAT_SINGLE)
34187828ca2SBarry Smith   ierr = PetscMalloc(baij->bs*sizeof(PetscScalar),&atmp);
34282502324SSatish Balay #endif
343a30f8f8cSSatish Balay   /* free stuff related to matrix-vec multiply */
344b0a32e0cSBarry Smith   ierr = VecGetSize(baij->lvec,&ec);CHKERRQ(ierr); /* needed for PetscLogObjectMemory below */
345a30f8f8cSSatish Balay   ierr = VecDestroy(baij->lvec);CHKERRQ(ierr); baij->lvec = 0;
346a30f8f8cSSatish Balay   ierr = VecScatterDestroy(baij->Mvctx);CHKERRQ(ierr); baij->Mvctx = 0;
347a30f8f8cSSatish Balay   if (baij->colmap) {
348a30f8f8cSSatish Balay #if defined (PETSC_USE_CTABLE)
349a30f8f8cSSatish Balay     ierr = PetscTableDelete(baij->colmap); baij->colmap = 0;CHKERRQ(ierr);
350a30f8f8cSSatish Balay #else
351a30f8f8cSSatish Balay     ierr = PetscFree(baij->colmap);CHKERRQ(ierr);
352a30f8f8cSSatish Balay     baij->colmap = 0;
353b0a32e0cSBarry Smith     PetscLogObjectMemory(A,-Bbaij->nbs*sizeof(int));
354a30f8f8cSSatish Balay #endif
355a30f8f8cSSatish Balay   }
356a30f8f8cSSatish Balay 
357a30f8f8cSSatish Balay   /* make sure that B is assembled so we can access its values */
358a30f8f8cSSatish Balay   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
359a30f8f8cSSatish Balay   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
360a30f8f8cSSatish Balay 
361a30f8f8cSSatish Balay   /* invent new B and copy stuff over */
36282502324SSatish Balay   ierr = PetscMalloc(mbs*sizeof(int),&nz);CHKERRQ(ierr);
363a30f8f8cSSatish Balay   for (i=0; i<mbs; i++) {
364a30f8f8cSSatish Balay     nz[i] = Bbaij->i[i+1]-Bbaij->i[i];
365a30f8f8cSSatish Balay   }
366be5d1d56SKris Buschelman   ierr = MatCreate(PETSC_COMM_SELF,m,n,m,n,&Bnew);CHKERRQ(ierr);
367be5d1d56SKris Buschelman   ierr = MatSetType(Bnew,B->type_name);CHKERRQ(ierr);
368be5d1d56SKris Buschelman   ierr = MatSeqBAIJSetPreallocation(Bnew,baij->bs,0,nz);CHKERRQ(ierr);
369a30f8f8cSSatish Balay   ierr = PetscFree(nz);CHKERRQ(ierr);
370a30f8f8cSSatish Balay 
37182502324SSatish Balay   ierr = PetscMalloc(bs*sizeof(int),&rvals);CHKERRQ(ierr);
372a30f8f8cSSatish Balay   for (i=0; i<mbs; i++) {
373a30f8f8cSSatish Balay     rvals[0] = bs*i;
374a30f8f8cSSatish Balay     for (j=1; j<bs; j++) { rvals[j] = rvals[j-1] + 1; }
375a30f8f8cSSatish Balay     for (j=Bbaij->i[i]; j<Bbaij->i[i+1]; j++) {
376a30f8f8cSSatish Balay       col = garray[Bbaij->j[j]]*bs;
377a30f8f8cSSatish Balay       for (k=0; k<bs; k++) {
378a30f8f8cSSatish Balay #if defined(PETSC_USE_MAT_SINGLE)
379a30f8f8cSSatish Balay         for (l=0; l<bs; l++) atmp[l] = a[j*bs2+l];
380a30f8f8cSSatish Balay #else
38171730473SSatish Balay         atmp = a+j*bs2 + k*bs;
382a30f8f8cSSatish Balay #endif
383c8407628SSatish Balay         ierr = MatSetValues_SeqSBAIJ(Bnew,bs,rvals,1,&col,atmp,B->insertmode);CHKERRQ(ierr);
384a30f8f8cSSatish Balay         col++;
385a30f8f8cSSatish Balay       }
386a30f8f8cSSatish Balay     }
387a30f8f8cSSatish Balay   }
388a30f8f8cSSatish Balay #if defined(PETSC_USE_MAT_SINGLE)
389a30f8f8cSSatish Balay   ierr = PetscFree(atmp);CHKERRQ(ierr);
390a30f8f8cSSatish Balay #endif
391a30f8f8cSSatish Balay   ierr = PetscFree(baij->garray);CHKERRQ(ierr);
392a30f8f8cSSatish Balay   baij->garray = 0;
393a30f8f8cSSatish Balay   ierr = PetscFree(rvals);CHKERRQ(ierr);
394b0a32e0cSBarry Smith   PetscLogObjectMemory(A,-ec*sizeof(int));
395a30f8f8cSSatish Balay   ierr = MatDestroy(B);CHKERRQ(ierr);
396b0a32e0cSBarry Smith   PetscLogObjectParent(A,Bnew);
397a30f8f8cSSatish Balay   baij->B = Bnew;
398a30f8f8cSSatish Balay   A->was_assembled = PETSC_FALSE;
399a30f8f8cSSatish Balay   PetscFunctionReturn(0);
400a30f8f8cSSatish Balay }
401a30f8f8cSSatish Balay 
402a30f8f8cSSatish Balay 
403a30f8f8cSSatish Balay 
404