xref: /petsc/src/mat/impls/sbaij/mpi/mmsbaij.c (revision 7adad9577d7f34aa90d02da3975546c0571d04a7)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
2a30f8f8cSSatish Balay 
3a30f8f8cSSatish Balay /*
4a30f8f8cSSatish Balay    Support for the parallel SBAIJ matrix vector multiply
5a30f8f8cSSatish Balay */
62896befeSSatish Balay #include "src/mat/impls/sbaij/mpi/mpisbaij.h"
77cff1425SSatish Balay 
813f74950SBarry Smith EXTERN PetscErrorCode MatSetValues_SeqSBAIJ(Mat,PetscInt,const PetscInt [],PetscInt,const PetscInt [],const PetscScalar [],InsertMode);
97cff1425SSatish Balay 
10a30f8f8cSSatish Balay 
114a2ae208SSatish Balay #undef __FUNCT__
124a2ae208SSatish Balay #define __FUNCT__ "MatSetUpMultiply_MPISBAIJ"
13dfbe8321SBarry Smith PetscErrorCode MatSetUpMultiply_MPISBAIJ(Mat mat)
14a30f8f8cSSatish Balay {
1540781036SHong Zhang   Mat_MPISBAIJ   *sbaij = (Mat_MPISBAIJ*)mat->data;
1640781036SHong Zhang   Mat_SeqBAIJ    *B = (Mat_SeqBAIJ*)(sbaij->B->data);
176849ba73SBarry Smith   PetscErrorCode ierr;
1813f74950SBarry Smith   PetscInt       Nbs = sbaij->Nbs,i,j,*indices,*aj = B->j,ec = 0,*garray,*sgarray;
19899cda47SBarry Smith   PetscInt       bs = mat->rmap.bs,*stmp,mbs=sbaij->mbs, vec_size,nt;
2040781036SHong Zhang   IS             from,to;
2140781036SHong Zhang   Vec            gvec;
2213f74950SBarry Smith   PetscMPIInt    rank=sbaij->rank,lsize,size=sbaij->size;
23899cda47SBarry Smith   PetscInt       *owners=sbaij->rangebs,*sowners,*ec_owner,k;
2440781036SHong Zhang   PetscScalar    *ptr;
2540781036SHong Zhang 
2640781036SHong Zhang   PetscFunctionBegin;
2701b2bd88SHong Zhang   if (sbaij->sMvctx) {
2801b2bd88SHong Zhang     /* This two lines should be in DisAssemble_MPISBAIJ(). Don't know why it causes crash there? */
2901b2bd88SHong Zhang     ierr = VecScatterDestroy(sbaij->sMvctx);CHKERRQ(ierr);
3001b2bd88SHong Zhang     sbaij->sMvctx = 0;
3140781036SHong Zhang   }
3240781036SHong Zhang 
3340781036SHong Zhang   /* For the first stab we make an array as long as the number of columns */
3440781036SHong Zhang   /* mark those columns that are in sbaij->B */
3513f74950SBarry Smith   ierr = PetscMalloc((Nbs+1)*sizeof(PetscInt),&indices);CHKERRQ(ierr);
3613f74950SBarry Smith   ierr = PetscMemzero(indices,Nbs*sizeof(PetscInt));CHKERRQ(ierr);
3740781036SHong Zhang   for (i=0; i<mbs; i++) {
3840781036SHong Zhang     for (j=0; j<B->ilen[i]; j++) {
3940781036SHong Zhang       if (!indices[aj[B->i[i] + j]]) ec++;
4040781036SHong Zhang       indices[aj[B->i[i] + j] ] = 1;
4140781036SHong Zhang     }
4240781036SHong Zhang   }
4340781036SHong Zhang 
4440781036SHong Zhang   /* form arrays of columns we need */
4513f74950SBarry Smith   ierr = PetscMalloc((ec+1)*sizeof(PetscInt),&garray);CHKERRQ(ierr);
4613f74950SBarry Smith   ierr = PetscMalloc((3*ec+1)*sizeof(PetscInt),&sgarray);CHKERRQ(ierr);
4740781036SHong Zhang   ec_owner = sgarray + 2*ec;
4840781036SHong Zhang 
4940781036SHong Zhang   ec = 0;
5040781036SHong Zhang   for (j=0; j<size; j++){
5140781036SHong Zhang     for (i=owners[j]; i<owners[j+1]; i++){
5240781036SHong Zhang       if (indices[i]) {
5340781036SHong Zhang         garray[ec]   = i;
5440781036SHong Zhang         ec_owner[ec] = j;
5540781036SHong Zhang         ec++;
5640781036SHong Zhang       }
5740781036SHong Zhang     }
5840781036SHong Zhang   }
5940781036SHong Zhang 
6040781036SHong Zhang   /* make indices now point into garray */
61b424e231SHong Zhang   for (i=0; i<ec; i++) indices[garray[i]] = i;
6240781036SHong Zhang 
6340781036SHong Zhang   /* compact out the extra columns in B */
6440781036SHong Zhang   for (i=0; i<mbs; i++) {
6540781036SHong Zhang     for (j=0; j<B->ilen[i]; j++) aj[B->i[i] + j] = indices[aj[B->i[i] + j]];
6640781036SHong Zhang   }
6740781036SHong Zhang   B->nbs      = ec;
686857c123SSatish Balay   sbaij->B->cmap.n = sbaij->B->cmap.N = ec*mat->rmap.bs;
696148ca0dSBarry Smith   ierr = PetscMapSetUp(&(sbaij->B->cmap));CHKERRQ(ierr);
7040781036SHong Zhang   ierr = PetscFree(indices);CHKERRQ(ierr);
7140781036SHong Zhang 
7240781036SHong Zhang   /* create local vector that is used to scatter into */
7340781036SHong Zhang   ierr = VecCreateSeq(PETSC_COMM_SELF,ec*bs,&sbaij->lvec);CHKERRQ(ierr);
7440781036SHong Zhang 
7540781036SHong Zhang   /* create two temporary index sets for building scatter-gather */
7613f74950SBarry Smith   ierr = PetscMalloc((2*ec+1)*sizeof(PetscInt),&stmp);CHKERRQ(ierr);
7740781036SHong Zhang   for (i=0; i<ec; i++) stmp[i] = bs*garray[i];
7840781036SHong Zhang   ierr = ISCreateBlock(PETSC_COMM_SELF,bs,ec,stmp,&from);CHKERRQ(ierr);
7940781036SHong Zhang 
8040781036SHong Zhang   for (i=0; i<ec; i++) { stmp[i] = bs*i; }
8140781036SHong Zhang   ierr = ISCreateBlock(PETSC_COMM_SELF,bs,ec,stmp,&to);CHKERRQ(ierr);
8240781036SHong Zhang 
8393d1592dSHong Zhang   /* generate the scatter context
8493d1592dSHong Zhang      -- Mvctx and lvec are not used by MatMult_MPISBAIJ(), but usefule for some applications */
85*7adad957SLisandro Dalcin   ierr = VecCreateMPI(((PetscObject)mat)->comm,mat->cmap.n,mat->cmap.N,&gvec);CHKERRQ(ierr);
8640781036SHong Zhang   ierr = VecScatterCreate(gvec,from,sbaij->lvec,to,&sbaij->Mvctx);CHKERRQ(ierr);
8740781036SHong Zhang 
8840781036SHong Zhang   sbaij->garray = garray;
8952e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,sbaij->Mvctx);CHKERRQ(ierr);
9052e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,sbaij->lvec);CHKERRQ(ierr);
9152e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,from);CHKERRQ(ierr);
9252e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,to);CHKERRQ(ierr);
9340781036SHong Zhang 
9440781036SHong Zhang   ierr = ISDestroy(from);CHKERRQ(ierr);
9540781036SHong Zhang   ierr = ISDestroy(to);CHKERRQ(ierr);
9640781036SHong Zhang 
9740781036SHong Zhang   /* create parallel vector that is used by SBAIJ matrix to scatter from/into */
9840781036SHong Zhang   lsize = (mbs + ec)*bs;
99*7adad957SLisandro Dalcin   ierr = VecCreateMPI(((PetscObject)mat)->comm,lsize,PETSC_DETERMINE,&sbaij->slvec0);CHKERRQ(ierr);
10040781036SHong Zhang   ierr = VecDuplicate(sbaij->slvec0,&sbaij->slvec1);CHKERRQ(ierr);
10140781036SHong Zhang   ierr = VecGetSize(sbaij->slvec0,&vec_size);CHKERRQ(ierr);
10240781036SHong Zhang 
103bba805e6SBarry Smith   sowners = sbaij->slvec0->map.range;
10440781036SHong Zhang 
10540781036SHong Zhang   /* x index in the IS sfrom */
10640781036SHong Zhang   for (i=0; i<ec; i++) {
10740781036SHong Zhang     j = ec_owner[i];
10840781036SHong Zhang     sgarray[i]  = garray[i] + (sowners[j]/bs - owners[j]);
10940781036SHong Zhang   }
11040781036SHong Zhang   /* b index in the IS sfrom */
11140781036SHong Zhang   k = sowners[rank]/bs + mbs;
11240781036SHong Zhang   for (i=ec,j=0; i< 2*ec; i++,j++) sgarray[i] = k + j;
11340781036SHong Zhang 
11440781036SHong Zhang   for (i=0; i<2*ec; i++) stmp[i] = bs*sgarray[i];
11540781036SHong Zhang   ierr = ISCreateBlock(PETSC_COMM_SELF,bs,2*ec,stmp,&from);CHKERRQ(ierr);
11640781036SHong Zhang 
11740781036SHong Zhang   /* x index in the IS sto */
11840781036SHong Zhang   k = sowners[rank]/bs + mbs;
11940781036SHong Zhang   for (i=0; i<ec; i++) stmp[i] = bs*(k + i);
12040781036SHong Zhang   /* b index in the IS sto */
12140781036SHong Zhang   for (i=ec; i<2*ec; i++) stmp[i] = bs*sgarray[i-ec];
12240781036SHong Zhang 
12340781036SHong Zhang   ierr = ISCreateBlock(PETSC_COMM_SELF,bs,2*ec,stmp,&to);CHKERRQ(ierr);
12440781036SHong Zhang 
12540781036SHong Zhang   ierr = VecScatterCreate(sbaij->slvec0,from,sbaij->slvec1,to,&sbaij->sMvctx);CHKERRQ(ierr);
12640781036SHong Zhang 
12740781036SHong Zhang   ierr = VecGetLocalSize(sbaij->slvec1,&nt);CHKERRQ(ierr);
1281ebc52fbSHong Zhang   ierr = VecGetArray(sbaij->slvec1,&ptr);CHKERRQ(ierr);
12940781036SHong Zhang   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,bs*mbs,ptr,&sbaij->slvec1a);CHKERRQ(ierr);
13040781036SHong Zhang   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,nt-bs*mbs,ptr+bs*mbs,&sbaij->slvec1b);CHKERRQ(ierr);
1311ebc52fbSHong Zhang   ierr = VecRestoreArray(sbaij->slvec1,&ptr);CHKERRQ(ierr);
13240781036SHong Zhang 
1331ebc52fbSHong Zhang   ierr = VecGetArray(sbaij->slvec0,&ptr);CHKERRQ(ierr);
13440781036SHong Zhang   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,nt-bs*mbs,ptr+bs*mbs,&sbaij->slvec0b);CHKERRQ(ierr);
1351ebc52fbSHong Zhang   ierr = VecRestoreArray(sbaij->slvec0,&ptr);CHKERRQ(ierr);
13640781036SHong Zhang 
13740781036SHong Zhang   ierr = PetscFree(stmp);CHKERRQ(ierr);
138*7adad957SLisandro Dalcin   ierr = MPI_Barrier(((PetscObject)mat)->comm);CHKERRQ(ierr);
13940781036SHong Zhang 
14052e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,sbaij->sMvctx);CHKERRQ(ierr);
14152e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,sbaij->slvec0);CHKERRQ(ierr);
14252e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,sbaij->slvec1);CHKERRQ(ierr);
14352e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,sbaij->slvec0b);CHKERRQ(ierr);
14452e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,sbaij->slvec1a);CHKERRQ(ierr);
14552e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,sbaij->slvec1b);CHKERRQ(ierr);
14652e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,from);CHKERRQ(ierr);
14752e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,to);CHKERRQ(ierr);
14840781036SHong Zhang 
14952e6d16bSBarry Smith   ierr = PetscLogObjectMemory(mat,(ec+1)*sizeof(PetscInt));CHKERRQ(ierr);
15040781036SHong Zhang   ierr = ISDestroy(from);CHKERRQ(ierr);
15140781036SHong Zhang   ierr = ISDestroy(to);CHKERRQ(ierr);
15240781036SHong Zhang   ierr = VecDestroy(gvec);CHKERRQ(ierr);
15340781036SHong Zhang   ierr = PetscFree(sgarray);CHKERRQ(ierr);
15440781036SHong Zhang   PetscFunctionReturn(0);
15540781036SHong Zhang }
15640781036SHong Zhang 
15740781036SHong Zhang #undef __FUNCT__
15840781036SHong Zhang #define __FUNCT__ "MatSetUpMultiply_MPISBAIJ_2comm"
159dfbe8321SBarry Smith PetscErrorCode MatSetUpMultiply_MPISBAIJ_2comm(Mat mat)
16040781036SHong Zhang {
1612896befeSSatish Balay   Mat_MPISBAIJ       *baij = (Mat_MPISBAIJ*)mat->data;
162a30f8f8cSSatish Balay   Mat_SeqBAIJ        *B = (Mat_SeqBAIJ*)(baij->B->data);
1636849ba73SBarry Smith   PetscErrorCode     ierr;
16413f74950SBarry Smith   PetscInt           i,j,*aj = B->j,ec = 0,*garray;
165899cda47SBarry Smith   PetscInt           bs = mat->rmap.bs,*stmp;
166a30f8f8cSSatish Balay   IS                 from,to;
167a30f8f8cSSatish Balay   Vec                gvec;
168a30f8f8cSSatish Balay #if defined (PETSC_USE_CTABLE)
169a30f8f8cSSatish Balay   PetscTable         gid1_lid1;
170a30f8f8cSSatish Balay   PetscTablePosition tpos;
17113f74950SBarry Smith   PetscInt           gid,lid;
1726f531f54SSatish Balay #else
17313f74950SBarry Smith   PetscInt           Nbs = baij->Nbs,*indices;
174a30f8f8cSSatish Balay #endif
175a30f8f8cSSatish Balay 
176a30f8f8cSSatish Balay   PetscFunctionBegin;
177a30f8f8cSSatish Balay #if defined (PETSC_USE_CTABLE)
178a30f8f8cSSatish Balay   /* use a table - Mark Adams */
179a30f8f8cSSatish Balay   PetscTableCreate(B->mbs,&gid1_lid1);
180a30f8f8cSSatish Balay   for (i=0; i<B->mbs; i++) {
181a30f8f8cSSatish Balay     for (j=0; j<B->ilen[i]; j++) {
18213f74950SBarry Smith       PetscInt data,gid1 = aj[B->i[i]+j] + 1;
183a30f8f8cSSatish Balay       ierr = PetscTableFind(gid1_lid1,gid1,&data);CHKERRQ(ierr);
184a30f8f8cSSatish Balay       if (!data) {
185a30f8f8cSSatish Balay         /* one based table */
186a30f8f8cSSatish Balay         ierr = PetscTableAdd(gid1_lid1,gid1,++ec);CHKERRQ(ierr);
187a30f8f8cSSatish Balay       }
188a30f8f8cSSatish Balay     }
189a30f8f8cSSatish Balay   }
190a30f8f8cSSatish Balay   /* form array of columns we need */
19113f74950SBarry Smith   ierr = PetscMalloc((ec+1)*sizeof(PetscInt),&garray);CHKERRQ(ierr);
192a30f8f8cSSatish Balay   ierr = PetscTableGetHeadPosition(gid1_lid1,&tpos);CHKERRQ(ierr);
193a30f8f8cSSatish Balay   while (tpos) {
194a30f8f8cSSatish Balay     ierr = PetscTableGetNext(gid1_lid1,&tpos,&gid,&lid);CHKERRQ(ierr);
195a30f8f8cSSatish Balay     gid--; lid--;
196a30f8f8cSSatish Balay     garray[lid] = gid;
197a30f8f8cSSatish Balay   }
198a30f8f8cSSatish Balay   ierr = PetscSortInt(ec,garray);CHKERRQ(ierr);
199a30f8f8cSSatish Balay   ierr = PetscTableRemoveAll(gid1_lid1);CHKERRQ(ierr);
200a30f8f8cSSatish Balay   for (i=0; i<ec; i++) {
201a30f8f8cSSatish Balay     ierr = PetscTableAdd(gid1_lid1,garray[i]+1,i+1);CHKERRQ(ierr);
202a30f8f8cSSatish Balay   }
203a30f8f8cSSatish Balay   /* compact out the extra columns in B */
204a30f8f8cSSatish Balay   for (i=0; i<B->mbs; i++) {
205a30f8f8cSSatish Balay     for (j=0; j<B->ilen[i]; j++) {
20613f74950SBarry Smith       PetscInt gid1 = aj[B->i[i] + j] + 1;
207a30f8f8cSSatish Balay       ierr = PetscTableFind(gid1_lid1,gid1,&lid);CHKERRQ(ierr);
208a30f8f8cSSatish Balay       lid --;
209a30f8f8cSSatish Balay       aj[B->i[i]+j] = lid;
210a30f8f8cSSatish Balay     }
211a30f8f8cSSatish Balay   }
212a30f8f8cSSatish Balay   B->nbs     = ec;
2136857c123SSatish Balay   baij->B->cmap.n = baij->B->cmap.N = ec*mat->rmap.bs;
2146148ca0dSBarry Smith   ierr = PetscMapSetUp(&(baij->B->cmap));CHKERRQ(ierr);
2159c666560SBarry Smith   ierr = PetscTableDestroy(gid1_lid1);CHKERRQ(ierr);
216a30f8f8cSSatish Balay   /* Mark Adams */
217a30f8f8cSSatish Balay #else
218a30f8f8cSSatish Balay   /* For the first stab we make an array as long as the number of columns */
219a30f8f8cSSatish Balay   /* mark those columns that are in baij->B */
22013f74950SBarry Smith   ierr = PetscMalloc((Nbs+1)*sizeof(PetscInt),&indices);CHKERRQ(ierr);
22113f74950SBarry Smith   ierr = PetscMemzero(indices,Nbs*sizeof(PetscInt));CHKERRQ(ierr);
222a30f8f8cSSatish Balay   for (i=0; i<B->mbs; i++) {
223a30f8f8cSSatish Balay     for (j=0; j<B->ilen[i]; j++) {
224a30f8f8cSSatish Balay       if (!indices[aj[B->i[i] + j]]) ec++;
225a30f8f8cSSatish Balay       indices[aj[B->i[i] + j] ] = 1;
226a30f8f8cSSatish Balay     }
227a30f8f8cSSatish Balay   }
228a30f8f8cSSatish Balay 
229a30f8f8cSSatish Balay   /* form array of columns we need */
23013f74950SBarry Smith   ierr = PetscMalloc((ec+1)*sizeof(PetscInt),&garray);CHKERRQ(ierr);
231f3ef73ceSBarry Smith   ec = 0;
232f3ef73ceSBarry Smith   for (i=0; i<Nbs; i++) {
233f3ef73ceSBarry Smith     if (indices[i]) {
234f3ef73ceSBarry Smith       garray[ec++] = i;
235f3ef73ceSBarry Smith     }
236f3ef73ceSBarry Smith   }
237a30f8f8cSSatish Balay 
238a30f8f8cSSatish Balay   /* make indices now point into garray */
239a30f8f8cSSatish Balay   for (i=0; i<ec; i++) {
240a30f8f8cSSatish Balay     indices[garray[i]] = i;
241a30f8f8cSSatish Balay   }
242a30f8f8cSSatish Balay 
243a30f8f8cSSatish Balay   /* compact out the extra columns in B */
244a30f8f8cSSatish Balay   for (i=0; i<B->mbs; i++) {
245a30f8f8cSSatish Balay     for (j=0; j<B->ilen[i]; j++) {
246a30f8f8cSSatish Balay       aj[B->i[i] + j] = indices[aj[B->i[i] + j]];
247a30f8f8cSSatish Balay     }
248a30f8f8cSSatish Balay   }
249a30f8f8cSSatish Balay   B->nbs       = ec;
250899cda47SBarry Smith   baij->B->cmap.n   = ec*mat->rmap.bs;
251a30f8f8cSSatish Balay   ierr = PetscFree(indices);CHKERRQ(ierr);
252a30f8f8cSSatish Balay #endif
253633e10c7SHong Zhang 
254a30f8f8cSSatish Balay   /* create local vector that is used to scatter into */
255a30f8f8cSSatish Balay   ierr = VecCreateSeq(PETSC_COMM_SELF,ec*bs,&baij->lvec);CHKERRQ(ierr);
256a30f8f8cSSatish Balay 
257a30f8f8cSSatish Balay   /* create two temporary index sets for building scatter-gather */
258eb7adc28SSatish Balay   for (i=0; i<ec; i++) {
259a30f8f8cSSatish Balay     garray[i] = bs*garray[i];
260a30f8f8cSSatish Balay   }
261a30f8f8cSSatish Balay   ierr = ISCreateBlock(PETSC_COMM_SELF,bs,ec,garray,&from);CHKERRQ(ierr);
26286e15631SSatish Balay   for (i=0; i<ec; i++) {
263a30f8f8cSSatish Balay     garray[i] = garray[i]/bs;
264a30f8f8cSSatish Balay   }
265a30f8f8cSSatish Balay 
26613f74950SBarry Smith   ierr = PetscMalloc((ec+1)*sizeof(PetscInt),&stmp);CHKERRQ(ierr);
267a30f8f8cSSatish Balay   for (i=0; i<ec; i++) { stmp[i] = bs*i; }
268a30f8f8cSSatish Balay   ierr = ISCreateBlock(PETSC_COMM_SELF,bs,ec,stmp,&to);CHKERRQ(ierr);
269a30f8f8cSSatish Balay   ierr = PetscFree(stmp);CHKERRQ(ierr);
270a30f8f8cSSatish Balay 
271a30f8f8cSSatish Balay   /* create temporary global vector to generate scatter context */
272a30f8f8cSSatish Balay   /* this is inefficient, but otherwise we must do either
273a30f8f8cSSatish Balay      1) save garray until the first actual scatter when the vector is known or
274a30f8f8cSSatish Balay      2) have another way of generating a scatter context without a vector.*/
275*7adad957SLisandro Dalcin   ierr = VecCreateMPI(((PetscObject)mat)->comm,mat->cmap.n,mat->cmap.N,&gvec);CHKERRQ(ierr);
276a30f8f8cSSatish Balay 
277a30f8f8cSSatish Balay   ierr = VecScatterCreate(gvec,from,baij->lvec,to,&baij->Mvctx);CHKERRQ(ierr);
278a30f8f8cSSatish Balay 
27952e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,baij->Mvctx);CHKERRQ(ierr);
28052e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,baij->lvec);CHKERRQ(ierr);
28152e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,from);CHKERRQ(ierr);
28252e6d16bSBarry Smith   ierr = PetscLogObjectParent(mat,to);CHKERRQ(ierr);
283a30f8f8cSSatish Balay   baij->garray = garray;
28452e6d16bSBarry Smith   ierr = PetscLogObjectMemory(mat,(ec+1)*sizeof(PetscInt));CHKERRQ(ierr);
285a30f8f8cSSatish Balay   ierr = ISDestroy(from);CHKERRQ(ierr);
286a30f8f8cSSatish Balay   ierr = ISDestroy(to);CHKERRQ(ierr);
287a30f8f8cSSatish Balay   ierr = VecDestroy(gvec);CHKERRQ(ierr);
288a30f8f8cSSatish Balay   PetscFunctionReturn(0);
289a30f8f8cSSatish Balay }
290a30f8f8cSSatish Balay 
291a30f8f8cSSatish Balay /*
29201b2bd88SHong Zhang      Takes the local part of an already assembled MPISBAIJ matrix
293a30f8f8cSSatish Balay    and disassembles it. This is to allow new nonzeros into the matrix
294a30f8f8cSSatish Balay    that require more communication in the matrix vector multiply.
295a30f8f8cSSatish Balay    Thus certain data-structures must be rebuilt.
296a30f8f8cSSatish Balay 
297a30f8f8cSSatish Balay    Kind of slow! But that's what application programmers get when
298a30f8f8cSSatish Balay    they are sloppy.
299a30f8f8cSSatish Balay */
3004a2ae208SSatish Balay #undef __FUNCT__
3014a2ae208SSatish Balay #define __FUNCT__ "DisAssemble_MPISBAIJ"
302dfbe8321SBarry Smith PetscErrorCode DisAssemble_MPISBAIJ(Mat A)
303a30f8f8cSSatish Balay {
3042896befeSSatish Balay   Mat_MPISBAIJ   *baij = (Mat_MPISBAIJ*)A->data;
305a30f8f8cSSatish Balay   Mat            B = baij->B,Bnew;
306a30f8f8cSSatish Balay   Mat_SeqBAIJ    *Bbaij = (Mat_SeqBAIJ*)B->data;
307dfbe8321SBarry Smith   PetscErrorCode ierr;
308899cda47SBarry Smith   PetscInt       i,j,mbs=Bbaij->mbs,n = A->cmap.n,col,*garray=baij->garray;
309899cda47SBarry Smith   PetscInt       k,bs=A->rmap.bs,bs2=baij->bs2,*rvals,*nz,ec,m=A->rmap.n;
310a30f8f8cSSatish Balay   MatScalar      *a = Bbaij->a;
31187828ca2SBarry Smith   PetscScalar    *atmp;
31282502324SSatish Balay #if defined(PETSC_USE_MAT_SINGLE)
31313f74950SBarry Smith   PetscInt       l;
314a30f8f8cSSatish Balay #endif
315a30f8f8cSSatish Balay 
316a30f8f8cSSatish Balay   PetscFunctionBegin;
31782502324SSatish Balay #if defined(PETSC_USE_MAT_SINGLE)
318899cda47SBarry Smith   ierr = PetscMalloc(A->rmap.bs*sizeof(PetscScalar),&atmp);
31982502324SSatish Balay #endif
320a30f8f8cSSatish Balay   /* free stuff related to matrix-vec multiply */
321b0a32e0cSBarry Smith   ierr = VecGetSize(baij->lvec,&ec);CHKERRQ(ierr); /* needed for PetscLogObjectMemory below */
32201b2bd88SHong Zhang   ierr = VecDestroy(baij->lvec);CHKERRQ(ierr);
32301b2bd88SHong Zhang   baij->lvec = 0;
32401b2bd88SHong Zhang   ierr = VecScatterDestroy(baij->Mvctx);CHKERRQ(ierr);
32501b2bd88SHong Zhang   baij->Mvctx = 0;
32601b2bd88SHong Zhang 
32701b2bd88SHong Zhang   ierr = VecDestroy(baij->slvec0);CHKERRQ(ierr);
32801b2bd88SHong Zhang   ierr = VecDestroy(baij->slvec0b);CHKERRQ(ierr);
32901b2bd88SHong Zhang   baij->slvec0 = 0;
33001b2bd88SHong Zhang   ierr = VecDestroy(baij->slvec1);CHKERRQ(ierr);
33101b2bd88SHong Zhang   ierr = VecDestroy(baij->slvec1a);CHKERRQ(ierr);
33201b2bd88SHong Zhang   ierr = VecDestroy(baij->slvec1b);CHKERRQ(ierr);
33301b2bd88SHong Zhang   baij->slvec1 = 0;
33401b2bd88SHong Zhang 
335a30f8f8cSSatish Balay   if (baij->colmap) {
336a30f8f8cSSatish Balay #if defined (PETSC_USE_CTABLE)
3379c666560SBarry Smith     ierr = PetscTableDestroy(baij->colmap); baij->colmap = 0;CHKERRQ(ierr);
338a30f8f8cSSatish Balay #else
339a30f8f8cSSatish Balay     ierr = PetscFree(baij->colmap);CHKERRQ(ierr);
340a30f8f8cSSatish Balay     baij->colmap = 0;
34152e6d16bSBarry Smith     ierr = PetscLogObjectMemory(A,-Bbaij->nbs*sizeof(PetscInt));CHKERRQ(ierr);
342a30f8f8cSSatish Balay #endif
343a30f8f8cSSatish Balay   }
344a30f8f8cSSatish Balay 
345a30f8f8cSSatish Balay   /* make sure that B is assembled so we can access its values */
346a30f8f8cSSatish Balay   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
347a30f8f8cSSatish Balay   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
348a30f8f8cSSatish Balay 
349a30f8f8cSSatish Balay   /* invent new B and copy stuff over */
35013f74950SBarry Smith   ierr = PetscMalloc(mbs*sizeof(PetscInt),&nz);CHKERRQ(ierr);
351a30f8f8cSSatish Balay   for (i=0; i<mbs; i++) {
352a30f8f8cSSatish Balay     nz[i] = Bbaij->i[i+1]-Bbaij->i[i];
353a30f8f8cSSatish Balay   }
354f69a0ea3SMatthew Knepley   ierr = MatCreate(PETSC_COMM_SELF,&Bnew);CHKERRQ(ierr);
355f69a0ea3SMatthew Knepley   ierr = MatSetSizes(Bnew,m,n,m,n);CHKERRQ(ierr);
356*7adad957SLisandro Dalcin   ierr = MatSetType(Bnew,((PetscObject)B)->type_name);CHKERRQ(ierr);
357899cda47SBarry Smith   ierr = MatSeqBAIJSetPreallocation(Bnew,B->rmap.bs,0,nz);CHKERRQ(ierr);
358a30f8f8cSSatish Balay   ierr = PetscFree(nz);CHKERRQ(ierr);
359a30f8f8cSSatish Balay 
36013f74950SBarry Smith   ierr = PetscMalloc(bs*sizeof(PetscInt),&rvals);CHKERRQ(ierr);
361a30f8f8cSSatish Balay   for (i=0; i<mbs; i++) {
362a30f8f8cSSatish Balay     rvals[0] = bs*i;
363a30f8f8cSSatish Balay     for (j=1; j<bs; j++) { rvals[j] = rvals[j-1] + 1; }
364a30f8f8cSSatish Balay     for (j=Bbaij->i[i]; j<Bbaij->i[i+1]; j++) {
365a30f8f8cSSatish Balay       col = garray[Bbaij->j[j]]*bs;
366a30f8f8cSSatish Balay       for (k=0; k<bs; k++) {
367a30f8f8cSSatish Balay #if defined(PETSC_USE_MAT_SINGLE)
368a30f8f8cSSatish Balay         for (l=0; l<bs; l++) atmp[l] = a[j*bs2+l];
369a30f8f8cSSatish Balay #else
37071730473SSatish Balay         atmp = a+j*bs2 + k*bs;
371a30f8f8cSSatish Balay #endif
372c8407628SSatish Balay         ierr = MatSetValues_SeqSBAIJ(Bnew,bs,rvals,1,&col,atmp,B->insertmode);CHKERRQ(ierr);
373a30f8f8cSSatish Balay         col++;
374a30f8f8cSSatish Balay       }
375a30f8f8cSSatish Balay     }
376a30f8f8cSSatish Balay   }
377a30f8f8cSSatish Balay #if defined(PETSC_USE_MAT_SINGLE)
378a30f8f8cSSatish Balay   ierr = PetscFree(atmp);CHKERRQ(ierr);
379a30f8f8cSSatish Balay #endif
380a30f8f8cSSatish Balay   ierr = PetscFree(baij->garray);CHKERRQ(ierr);
381a30f8f8cSSatish Balay   baij->garray = 0;
382a30f8f8cSSatish Balay   ierr = PetscFree(rvals);CHKERRQ(ierr);
38352e6d16bSBarry Smith   ierr = PetscLogObjectMemory(A,-ec*sizeof(PetscInt));CHKERRQ(ierr);
384a30f8f8cSSatish Balay   ierr = MatDestroy(B);CHKERRQ(ierr);
38552e6d16bSBarry Smith   ierr = PetscLogObjectParent(A,Bnew);CHKERRQ(ierr);
386a30f8f8cSSatish Balay   baij->B = Bnew;
387a30f8f8cSSatish Balay   A->was_assembled = PETSC_FALSE;
388a30f8f8cSSatish Balay   PetscFunctionReturn(0);
389a30f8f8cSSatish Balay }
390a30f8f8cSSatish Balay 
391a30f8f8cSSatish Balay 
392a30f8f8cSSatish Balay 
393