xref: /petsc/src/mat/impls/sbaij/mpi/mmsbaij.c (revision dcca6d9d80ebd869fe6029bd05a3aa9faafef49e)
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;
22077aedafSJed Brown   PetscInt       *owners=sbaij->rangebs,*ec_owner,k;
23077aedafSJed Brown   const PetscInt *sowners;
2440781036SHong Zhang   PetscScalar    *ptr;
2540781036SHong Zhang 
2640781036SHong Zhang   PetscFunctionBegin;
276bf464f9SBarry Smith   ierr = VecScatterDestroy(&sbaij->sMvctx);CHKERRQ(ierr);
2840781036SHong Zhang 
2940781036SHong Zhang   /* For the first stab we make an array as long as the number of columns */
3040781036SHong Zhang   /* mark those columns that are in sbaij->B */
3174ed9c26SBarry Smith   ierr = PetscMalloc(Nbs*sizeof(PetscInt),&indices);CHKERRQ(ierr);
3213f74950SBarry Smith   ierr = PetscMemzero(indices,Nbs*sizeof(PetscInt));CHKERRQ(ierr);
3340781036SHong Zhang   for (i=0; i<mbs; i++) {
3440781036SHong Zhang     for (j=0; j<B->ilen[i]; j++) {
3540781036SHong Zhang       if (!indices[aj[B->i[i] + j]]) ec++;
3640781036SHong Zhang       indices[aj[B->i[i] + j]] = 1;
3740781036SHong Zhang     }
3840781036SHong Zhang   }
3940781036SHong Zhang 
4040781036SHong Zhang   /* form arrays of columns we need */
4174ed9c26SBarry Smith   ierr = PetscMalloc(ec*sizeof(PetscInt),&garray);CHKERRQ(ierr);
42*dcca6d9dSJed Brown   ierr = PetscMalloc2(2*ec,&sgarray,ec,&ec_owner);CHKERRQ(ierr);
4340781036SHong Zhang 
4440781036SHong Zhang   ec = 0;
4540781036SHong Zhang   for (j=0; j<size; j++) {
4640781036SHong Zhang     for (i=owners[j]; i<owners[j+1]; i++) {
4740781036SHong Zhang       if (indices[i]) {
4840781036SHong Zhang         garray[ec]   = i;
4940781036SHong Zhang         ec_owner[ec] = j;
5040781036SHong Zhang         ec++;
5140781036SHong Zhang       }
5240781036SHong Zhang     }
5340781036SHong Zhang   }
5440781036SHong Zhang 
5540781036SHong Zhang   /* make indices now point into garray */
56b424e231SHong Zhang   for (i=0; i<ec; i++) indices[garray[i]] = i;
5740781036SHong Zhang 
5840781036SHong Zhang   /* compact out the extra columns in B */
5940781036SHong Zhang   for (i=0; i<mbs; i++) {
6040781036SHong Zhang     for (j=0; j<B->ilen[i]; j++) aj[B->i[i] + j] = indices[aj[B->i[i] + j]];
6140781036SHong Zhang   }
6240781036SHong Zhang   B->nbs = ec;
6326fbe8dcSKarl Rupp 
64d0f46423SBarry Smith   sbaij->B->cmap->n = sbaij->B->cmap->N = ec*mat->rmap->bs;
6526fbe8dcSKarl Rupp 
6626283091SBarry Smith   ierr = PetscLayoutSetUp((sbaij->B->cmap));CHKERRQ(ierr);
6740781036SHong Zhang   ierr = PetscFree(indices);CHKERRQ(ierr);
6840781036SHong Zhang 
6940781036SHong Zhang   /* create local vector that is used to scatter into */
7040781036SHong Zhang   ierr = VecCreateSeq(PETSC_COMM_SELF,ec*bs,&sbaij->lvec);CHKERRQ(ierr);
7140781036SHong Zhang 
7240781036SHong Zhang   /* create two temporary index sets for building scatter-gather */
7374ed9c26SBarry Smith   ierr = PetscMalloc(2*ec*sizeof(PetscInt),&stmp);CHKERRQ(ierr);
74deff0451SBarry Smith   ierr = ISCreateBlock(PETSC_COMM_SELF,bs,ec,garray,PETSC_COPY_VALUES,&from);CHKERRQ(ierr);
7526fbe8dcSKarl Rupp   for (i=0; i<ec; i++) stmp[i] = i;
76deff0451SBarry Smith   ierr = ISCreateBlock(PETSC_COMM_SELF,bs,ec,stmp,PETSC_COPY_VALUES,&to);CHKERRQ(ierr);
7740781036SHong Zhang 
7893d1592dSHong Zhang   /* generate the scatter context
7993d1592dSHong Zhang      -- Mvctx and lvec are not used by MatMult_MPISBAIJ(), but usefule for some applications */
80ce94432eSBarry Smith   ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)mat),1,mat->cmap->n,mat->cmap->N,NULL,&gvec);CHKERRQ(ierr);
8140781036SHong Zhang   ierr = VecScatterCreate(gvec,from,sbaij->lvec,to,&sbaij->Mvctx);CHKERRQ(ierr);
826bf464f9SBarry Smith   ierr = VecDestroy(&gvec);CHKERRQ(ierr);
8340781036SHong Zhang 
8440781036SHong Zhang   sbaij->garray = garray;
8526fbe8dcSKarl Rupp 
863bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->Mvctx);CHKERRQ(ierr);
873bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->lvec);CHKERRQ(ierr);
883bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)from);CHKERRQ(ierr);
893bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)to);CHKERRQ(ierr);
9040781036SHong Zhang 
916bf464f9SBarry Smith   ierr = ISDestroy(&from);CHKERRQ(ierr);
926bf464f9SBarry Smith   ierr = ISDestroy(&to);CHKERRQ(ierr);
9340781036SHong Zhang 
9440781036SHong Zhang   /* create parallel vector that is used by SBAIJ matrix to scatter from/into */
9540781036SHong Zhang   lsize = (mbs + ec)*bs;
96ce94432eSBarry Smith   ierr  = VecCreateMPI(PetscObjectComm((PetscObject)mat),lsize,PETSC_DETERMINE,&sbaij->slvec0);CHKERRQ(ierr);
9740781036SHong Zhang   ierr  = VecDuplicate(sbaij->slvec0,&sbaij->slvec1);CHKERRQ(ierr);
9840781036SHong Zhang   ierr  = VecGetSize(sbaij->slvec0,&vec_size);CHKERRQ(ierr);
9940781036SHong Zhang 
100077aedafSJed Brown   ierr = VecGetOwnershipRanges(sbaij->slvec0,&sowners);CHKERRQ(ierr);
10140781036SHong Zhang 
10240781036SHong Zhang   /* x index in the IS sfrom */
10340781036SHong Zhang   for (i=0; i<ec; i++) {
10440781036SHong Zhang     j          = ec_owner[i];
10540781036SHong Zhang     sgarray[i] = garray[i] + (sowners[j]/bs - owners[j]);
10640781036SHong Zhang   }
10740781036SHong Zhang   /* b index in the IS sfrom */
10840781036SHong Zhang   k = sowners[rank]/bs + mbs;
10940781036SHong Zhang   for (i=ec,j=0; i< 2*ec; i++,j++) sgarray[i] = k + j;
110deff0451SBarry Smith   ierr = ISCreateBlock(PETSC_COMM_SELF,bs,2*ec,sgarray,PETSC_COPY_VALUES,&from);CHKERRQ(ierr);
11140781036SHong Zhang 
11240781036SHong Zhang   /* x index in the IS sto */
11340781036SHong Zhang   k = sowners[rank]/bs + mbs;
114e82e9f6bSBarry Smith   for (i=0; i<ec; i++) stmp[i] = (k + i);
11540781036SHong Zhang   /* b index in the IS sto */
116e82e9f6bSBarry Smith   for (i=ec; i<2*ec; i++) stmp[i] = sgarray[i-ec];
11740781036SHong Zhang 
118deff0451SBarry Smith   ierr = ISCreateBlock(PETSC_COMM_SELF,bs,2*ec,stmp,PETSC_COPY_VALUES,&to);CHKERRQ(ierr);
11940781036SHong Zhang 
12040781036SHong Zhang   ierr = VecScatterCreate(sbaij->slvec0,from,sbaij->slvec1,to,&sbaij->sMvctx);CHKERRQ(ierr);
12140781036SHong Zhang 
12240781036SHong Zhang   ierr = VecGetLocalSize(sbaij->slvec1,&nt);CHKERRQ(ierr);
1231ebc52fbSHong Zhang   ierr = VecGetArray(sbaij->slvec1,&ptr);CHKERRQ(ierr);
124778a2246SBarry Smith   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,bs*mbs,ptr,&sbaij->slvec1a);CHKERRQ(ierr);
125778a2246SBarry Smith   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,nt-bs*mbs,ptr+bs*mbs,&sbaij->slvec1b);CHKERRQ(ierr);
1261ebc52fbSHong Zhang   ierr = VecRestoreArray(sbaij->slvec1,&ptr);CHKERRQ(ierr);
12740781036SHong Zhang 
1281ebc52fbSHong Zhang   ierr = VecGetArray(sbaij->slvec0,&ptr);CHKERRQ(ierr);
129778a2246SBarry Smith   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,nt-bs*mbs,ptr+bs*mbs,&sbaij->slvec0b);CHKERRQ(ierr);
1301ebc52fbSHong Zhang   ierr = VecRestoreArray(sbaij->slvec0,&ptr);CHKERRQ(ierr);
13140781036SHong Zhang 
13240781036SHong Zhang   ierr = PetscFree(stmp);CHKERRQ(ierr);
133ce94432eSBarry Smith   ierr = MPI_Barrier(PetscObjectComm((PetscObject)mat));CHKERRQ(ierr);
13440781036SHong Zhang 
1353bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->sMvctx);CHKERRQ(ierr);
1363bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->slvec0);CHKERRQ(ierr);
1373bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->slvec1);CHKERRQ(ierr);
1383bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->slvec0b);CHKERRQ(ierr);
1393bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->slvec1a);CHKERRQ(ierr);
1403bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->slvec1b);CHKERRQ(ierr);
1413bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)from);CHKERRQ(ierr);
1423bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)to);CHKERRQ(ierr);
14340781036SHong Zhang 
1443bb1ff40SBarry Smith   ierr = PetscLogObjectMemory((PetscObject)mat,(ec+1)*sizeof(PetscInt));CHKERRQ(ierr);
1456bf464f9SBarry Smith   ierr = ISDestroy(&from);CHKERRQ(ierr);
1466bf464f9SBarry Smith   ierr = ISDestroy(&to);CHKERRQ(ierr);
14774ed9c26SBarry Smith   ierr = PetscFree2(sgarray,ec_owner);CHKERRQ(ierr);
14840781036SHong Zhang   PetscFunctionReturn(0);
14940781036SHong Zhang }
15040781036SHong Zhang 
15140781036SHong Zhang #undef __FUNCT__
15240781036SHong Zhang #define __FUNCT__ "MatSetUpMultiply_MPISBAIJ_2comm"
153dfbe8321SBarry Smith PetscErrorCode MatSetUpMultiply_MPISBAIJ_2comm(Mat mat)
15440781036SHong Zhang {
1552896befeSSatish Balay   Mat_MPISBAIJ   *baij = (Mat_MPISBAIJ*)mat->data;
156a30f8f8cSSatish Balay   Mat_SeqBAIJ    *B    = (Mat_SeqBAIJ*)(baij->B->data);
1576849ba73SBarry Smith   PetscErrorCode ierr;
15813f74950SBarry Smith   PetscInt       i,j,*aj = B->j,ec = 0,*garray;
159d0f46423SBarry Smith   PetscInt       bs = mat->rmap->bs,*stmp;
160a30f8f8cSSatish Balay   IS             from,to;
161a30f8f8cSSatish Balay   Vec            gvec;
162a30f8f8cSSatish Balay #if defined(PETSC_USE_CTABLE)
163a30f8f8cSSatish Balay   PetscTable         gid1_lid1;
164a30f8f8cSSatish Balay   PetscTablePosition tpos;
16513f74950SBarry Smith   PetscInt           gid,lid;
1666f531f54SSatish Balay #else
16713f74950SBarry Smith   PetscInt Nbs = baij->Nbs,*indices;
168a30f8f8cSSatish Balay #endif
169a30f8f8cSSatish Balay 
170a30f8f8cSSatish Balay   PetscFunctionBegin;
171a30f8f8cSSatish Balay #if defined(PETSC_USE_CTABLE)
172a30f8f8cSSatish Balay   /* use a table - Mark Adams */
173e23dfa41SBarry Smith   PetscTableCreate(B->mbs,baij->Nbs+1,&gid1_lid1);
174a30f8f8cSSatish Balay   for (i=0; i<B->mbs; i++) {
175a30f8f8cSSatish Balay     for (j=0; j<B->ilen[i]; j++) {
17613f74950SBarry Smith       PetscInt data,gid1 = aj[B->i[i]+j] + 1;
177a30f8f8cSSatish Balay       ierr = PetscTableFind(gid1_lid1,gid1,&data);CHKERRQ(ierr);
178a30f8f8cSSatish Balay       if (!data) {
179a30f8f8cSSatish Balay         /* one based table */
1803861aac3SJed Brown         ierr = PetscTableAdd(gid1_lid1,gid1,++ec,INSERT_VALUES);CHKERRQ(ierr);
181a30f8f8cSSatish Balay       }
182a30f8f8cSSatish Balay     }
183a30f8f8cSSatish Balay   }
184a30f8f8cSSatish Balay   /* form array of columns we need */
18574ed9c26SBarry Smith   ierr = PetscMalloc(ec*sizeof(PetscInt),&garray);CHKERRQ(ierr);
186a30f8f8cSSatish Balay   ierr = PetscTableGetHeadPosition(gid1_lid1,&tpos);CHKERRQ(ierr);
187a30f8f8cSSatish Balay   while (tpos) {
188a30f8f8cSSatish Balay     ierr = PetscTableGetNext(gid1_lid1,&tpos,&gid,&lid);CHKERRQ(ierr);
189a30f8f8cSSatish Balay     gid--; lid--;
190a30f8f8cSSatish Balay     garray[lid] = gid;
191a30f8f8cSSatish Balay   }
192a30f8f8cSSatish Balay   ierr = PetscSortInt(ec,garray);CHKERRQ(ierr);
193a30f8f8cSSatish Balay   ierr = PetscTableRemoveAll(gid1_lid1);CHKERRQ(ierr);
194a30f8f8cSSatish Balay   for (i=0; i<ec; i++) {
1953861aac3SJed Brown     ierr = PetscTableAdd(gid1_lid1,garray[i]+1,i+1,INSERT_VALUES);CHKERRQ(ierr);
196a30f8f8cSSatish Balay   }
197a30f8f8cSSatish Balay   /* compact out the extra columns in B */
198a30f8f8cSSatish Balay   for (i=0; i<B->mbs; i++) {
199a30f8f8cSSatish Balay     for (j=0; j<B->ilen[i]; j++) {
20013f74950SBarry Smith       PetscInt gid1 = aj[B->i[i] + j] + 1;
201a30f8f8cSSatish Balay       ierr = PetscTableFind(gid1_lid1,gid1,&lid);CHKERRQ(ierr);
202a30f8f8cSSatish Balay       lid--;
203a30f8f8cSSatish Balay       aj[B->i[i]+j] = lid;
204a30f8f8cSSatish Balay     }
205a30f8f8cSSatish Balay   }
206a30f8f8cSSatish Balay   B->nbs = ec;
20726fbe8dcSKarl Rupp 
208d0f46423SBarry Smith   baij->B->cmap->n = baij->B->cmap->N = ec*mat->rmap->bs;
20926fbe8dcSKarl Rupp 
21026283091SBarry Smith   ierr = PetscLayoutSetUp((baij->B->cmap));CHKERRQ(ierr);
2116bc0bbbfSBarry Smith   ierr = PetscTableDestroy(&gid1_lid1);CHKERRQ(ierr);
212a30f8f8cSSatish Balay #else
213a30f8f8cSSatish Balay   /* For the first stab we make an array as long as the number of columns */
214a30f8f8cSSatish Balay   /* mark those columns that are in baij->B */
21574ed9c26SBarry Smith   ierr = PetscMalloc(Nbs*sizeof(PetscInt),&indices);CHKERRQ(ierr);
21613f74950SBarry Smith   ierr = PetscMemzero(indices,Nbs*sizeof(PetscInt));CHKERRQ(ierr);
217a30f8f8cSSatish Balay   for (i=0; i<B->mbs; i++) {
218a30f8f8cSSatish Balay     for (j=0; j<B->ilen[i]; j++) {
219a30f8f8cSSatish Balay       if (!indices[aj[B->i[i] + j]]) ec++;
220a30f8f8cSSatish Balay       indices[aj[B->i[i] + j]] = 1;
221a30f8f8cSSatish Balay     }
222a30f8f8cSSatish Balay   }
223a30f8f8cSSatish Balay 
224a30f8f8cSSatish Balay   /* form array of columns we need */
22574ed9c26SBarry Smith   ierr = PetscMalloc(ec*sizeof(PetscInt),&garray);CHKERRQ(ierr);
226f3ef73ceSBarry Smith   ec = 0;
227f3ef73ceSBarry Smith   for (i=0; i<Nbs; i++) {
228f3ef73ceSBarry Smith     if (indices[i]) {
229f3ef73ceSBarry Smith       garray[ec++] = i;
230f3ef73ceSBarry Smith     }
231f3ef73ceSBarry Smith   }
232a30f8f8cSSatish Balay 
233a30f8f8cSSatish Balay   /* make indices now point into garray */
23426fbe8dcSKarl Rupp   for (i=0; i<ec; i++) indices[garray[i]] = i;
235a30f8f8cSSatish Balay 
236a30f8f8cSSatish Balay   /* compact out the extra columns in B */
237a30f8f8cSSatish Balay   for (i=0; i<B->mbs; i++) {
238a30f8f8cSSatish Balay     for (j=0; j<B->ilen[i]; j++) {
239a30f8f8cSSatish Balay       aj[B->i[i] + j] = indices[aj[B->i[i] + j]];
240a30f8f8cSSatish Balay     }
241a30f8f8cSSatish Balay   }
242a30f8f8cSSatish Balay   B->nbs = ec;
24326fbe8dcSKarl Rupp 
244d0f46423SBarry Smith   baij->B->cmap->n = ec*mat->rmap->bs;
24526fbe8dcSKarl Rupp 
246a30f8f8cSSatish Balay   ierr = PetscFree(indices);CHKERRQ(ierr);
247a30f8f8cSSatish Balay #endif
248633e10c7SHong Zhang 
249a30f8f8cSSatish Balay   /* create local vector that is used to scatter into */
250a30f8f8cSSatish Balay   ierr = VecCreateSeq(PETSC_COMM_SELF,ec*bs,&baij->lvec);CHKERRQ(ierr);
251a30f8f8cSSatish Balay 
252a30f8f8cSSatish Balay   /* create two temporary index sets for building scatter-gather */
253deff0451SBarry Smith   ierr = ISCreateBlock(PETSC_COMM_SELF,bs,ec,garray,PETSC_COPY_VALUES,&from);CHKERRQ(ierr);
254a30f8f8cSSatish Balay 
25574ed9c26SBarry Smith   ierr = PetscMalloc(ec*sizeof(PetscInt),&stmp);CHKERRQ(ierr);
25626fbe8dcSKarl Rupp   for (i=0; i<ec; i++) stmp[i] = i;
257deff0451SBarry Smith   ierr = ISCreateBlock(PETSC_COMM_SELF,bs,ec,stmp,PETSC_OWN_POINTER,&to);CHKERRQ(ierr);
258a30f8f8cSSatish Balay 
259a30f8f8cSSatish Balay   /* create temporary global vector to generate scatter context */
260ce94432eSBarry Smith   ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)mat),1,mat->cmap->n,mat->cmap->N,NULL,&gvec);CHKERRQ(ierr);
261a30f8f8cSSatish Balay   ierr = VecScatterCreate(gvec,from,baij->lvec,to,&baij->Mvctx);CHKERRQ(ierr);
2626bf464f9SBarry Smith   ierr = VecDestroy(&gvec);CHKERRQ(ierr);
263a30f8f8cSSatish Balay 
2643bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)baij->Mvctx);CHKERRQ(ierr);
2653bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)baij->lvec);CHKERRQ(ierr);
2663bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)from);CHKERRQ(ierr);
2673bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)to);CHKERRQ(ierr);
26826fbe8dcSKarl Rupp 
269a30f8f8cSSatish Balay   baij->garray = garray;
27026fbe8dcSKarl Rupp 
2713bb1ff40SBarry Smith   ierr = PetscLogObjectMemory((PetscObject)mat,(ec+1)*sizeof(PetscInt));CHKERRQ(ierr);
2726bf464f9SBarry Smith   ierr = ISDestroy(&from);CHKERRQ(ierr);
2736bf464f9SBarry Smith   ierr = ISDestroy(&to);CHKERRQ(ierr);
274a30f8f8cSSatish Balay   PetscFunctionReturn(0);
275a30f8f8cSSatish Balay }
276a30f8f8cSSatish Balay 
277a30f8f8cSSatish Balay /*
27801b2bd88SHong Zhang      Takes the local part of an already assembled MPISBAIJ matrix
279a30f8f8cSSatish Balay    and disassembles it. This is to allow new nonzeros into the matrix
280a30f8f8cSSatish Balay    that require more communication in the matrix vector multiply.
281a30f8f8cSSatish Balay    Thus certain data-structures must be rebuilt.
282a30f8f8cSSatish Balay 
283a30f8f8cSSatish Balay    Kind of slow! But that's what application programmers get when
284a30f8f8cSSatish Balay    they are sloppy.
285a30f8f8cSSatish Balay */
2864a2ae208SSatish Balay #undef __FUNCT__
287ab9863d7SBarry Smith #define __FUNCT__ "MatDisAssemble_MPISBAIJ"
288ab9863d7SBarry Smith PetscErrorCode MatDisAssemble_MPISBAIJ(Mat A)
289a30f8f8cSSatish Balay {
2902896befeSSatish Balay   Mat_MPISBAIJ   *baij  = (Mat_MPISBAIJ*)A->data;
291a30f8f8cSSatish Balay   Mat            B      = baij->B,Bnew;
292a30f8f8cSSatish Balay   Mat_SeqBAIJ    *Bbaij = (Mat_SeqBAIJ*)B->data;
293dfbe8321SBarry Smith   PetscErrorCode ierr;
294d0f46423SBarry Smith   PetscInt       i,j,mbs=Bbaij->mbs,n = A->cmap->N,col,*garray=baij->garray;
295d0f46423SBarry Smith   PetscInt       k,bs=A->rmap->bs,bs2=baij->bs2,*rvals,*nz,ec,m=A->rmap->n;
296a30f8f8cSSatish Balay   MatScalar      *a = Bbaij->a;
29787828ca2SBarry Smith   PetscScalar    *atmp;
298ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL_MAT_SINGLE)
29913f74950SBarry Smith   PetscInt l;
300a30f8f8cSSatish Balay #endif
301a30f8f8cSSatish Balay 
302a30f8f8cSSatish Balay   PetscFunctionBegin;
303ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL_MAT_SINGLE)
304d0f46423SBarry Smith   ierr = PetscMalloc(A->rmap->bs*sizeof(PetscScalar),&atmp);
30582502324SSatish Balay #endif
306a30f8f8cSSatish Balay   /* free stuff related to matrix-vec multiply */
307b0a32e0cSBarry Smith   ierr = VecGetSize(baij->lvec,&ec);CHKERRQ(ierr); /* needed for PetscLogObjectMemory below */
3086bf464f9SBarry Smith   ierr = VecDestroy(&baij->lvec);CHKERRQ(ierr);
3096bf464f9SBarry Smith   ierr = VecScatterDestroy(&baij->Mvctx);CHKERRQ(ierr);
31001b2bd88SHong Zhang 
3116bf464f9SBarry Smith   ierr = VecDestroy(&baij->slvec0);CHKERRQ(ierr);
3126bf464f9SBarry Smith   ierr = VecDestroy(&baij->slvec0b);CHKERRQ(ierr);
3136bf464f9SBarry Smith   ierr = VecDestroy(&baij->slvec1);CHKERRQ(ierr);
3146bf464f9SBarry Smith   ierr = VecDestroy(&baij->slvec1a);CHKERRQ(ierr);
3156bf464f9SBarry Smith   ierr = VecDestroy(&baij->slvec1b);CHKERRQ(ierr);
31601b2bd88SHong Zhang 
317a30f8f8cSSatish Balay   if (baij->colmap) {
318a30f8f8cSSatish Balay #if defined(PETSC_USE_CTABLE)
3196bc0bbbfSBarry Smith     ierr = PetscTableDestroy(&baij->colmap);CHKERRQ(ierr);
320a30f8f8cSSatish Balay #else
321a30f8f8cSSatish Balay     ierr = PetscFree(baij->colmap);CHKERRQ(ierr);
3223bb1ff40SBarry Smith     ierr = PetscLogObjectMemory((PetscObject)A,-Bbaij->nbs*sizeof(PetscInt));CHKERRQ(ierr);
323a30f8f8cSSatish Balay #endif
324a30f8f8cSSatish Balay   }
325a30f8f8cSSatish Balay 
326a30f8f8cSSatish Balay   /* make sure that B is assembled so we can access its values */
327a30f8f8cSSatish Balay   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
328a30f8f8cSSatish Balay   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
329a30f8f8cSSatish Balay 
330a30f8f8cSSatish Balay   /* invent new B and copy stuff over */
33113f74950SBarry Smith   ierr = PetscMalloc(mbs*sizeof(PetscInt),&nz);CHKERRQ(ierr);
332a30f8f8cSSatish Balay   for (i=0; i<mbs; i++) {
333a30f8f8cSSatish Balay     nz[i] = Bbaij->i[i+1]-Bbaij->i[i];
334a30f8f8cSSatish Balay   }
335f69a0ea3SMatthew Knepley   ierr = MatCreate(PETSC_COMM_SELF,&Bnew);CHKERRQ(ierr);
336f69a0ea3SMatthew Knepley   ierr = MatSetSizes(Bnew,m,n,m,n);CHKERRQ(ierr);
3377adad957SLisandro Dalcin   ierr = MatSetType(Bnew,((PetscObject)B)->type_name);CHKERRQ(ierr);
338d0f46423SBarry Smith   ierr = MatSeqBAIJSetPreallocation(Bnew,B->rmap->bs,0,nz);CHKERRQ(ierr);
33926fbe8dcSKarl Rupp 
3402576faa2SJed Brown   ((Mat_SeqSBAIJ*)Bnew->data)->nonew = Bbaij->nonew; /* Inherit insertion error options. */
34126fbe8dcSKarl Rupp 
342a30f8f8cSSatish Balay   ierr = PetscFree(nz);CHKERRQ(ierr);
343a30f8f8cSSatish Balay 
34413f74950SBarry Smith   ierr = PetscMalloc(bs*sizeof(PetscInt),&rvals);CHKERRQ(ierr);
345a30f8f8cSSatish Balay   for (i=0; i<mbs; i++) {
346a30f8f8cSSatish Balay     rvals[0] = bs*i;
34726fbe8dcSKarl Rupp     for (j=1; j<bs; j++) rvals[j] = rvals[j-1] + 1;
348a30f8f8cSSatish Balay     for (j=Bbaij->i[i]; j<Bbaij->i[i+1]; j++) {
349a30f8f8cSSatish Balay       col = garray[Bbaij->j[j]]*bs;
350a30f8f8cSSatish Balay       for (k=0; k<bs; k++) {
351ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL_MAT_SINGLE)
352a30f8f8cSSatish Balay         for (l=0; l<bs; l++) atmp[l] = a[j*bs2+l];
353a30f8f8cSSatish Balay #else
35471730473SSatish Balay         atmp = a+j*bs2 + k*bs;
355a30f8f8cSSatish Balay #endif
356c8407628SSatish Balay         ierr = MatSetValues_SeqSBAIJ(Bnew,bs,rvals,1,&col,atmp,B->insertmode);CHKERRQ(ierr);
357a30f8f8cSSatish Balay         col++;
358a30f8f8cSSatish Balay       }
359a30f8f8cSSatish Balay     }
360a30f8f8cSSatish Balay   }
361ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL_MAT_SINGLE)
362a30f8f8cSSatish Balay   ierr = PetscFree(atmp);CHKERRQ(ierr);
363a30f8f8cSSatish Balay #endif
364a30f8f8cSSatish Balay   ierr = PetscFree(baij->garray);CHKERRQ(ierr);
36526fbe8dcSKarl Rupp 
366a30f8f8cSSatish Balay   baij->garray = 0;
36726fbe8dcSKarl Rupp 
368a30f8f8cSSatish Balay   ierr = PetscFree(rvals);CHKERRQ(ierr);
3693bb1ff40SBarry Smith   ierr = PetscLogObjectMemory((PetscObject)A,-ec*sizeof(PetscInt));CHKERRQ(ierr);
3706bf464f9SBarry Smith   ierr = MatDestroy(&B);CHKERRQ(ierr);
3713bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)A,(PetscObject)Bnew);CHKERRQ(ierr);
37226fbe8dcSKarl Rupp 
373a30f8f8cSSatish Balay   baij->B = Bnew;
37426fbe8dcSKarl Rupp 
375a30f8f8cSSatish Balay   A->was_assembled = PETSC_FALSE;
376a30f8f8cSSatish Balay   PetscFunctionReturn(0);
377a30f8f8cSSatish Balay }
378a30f8f8cSSatish Balay 
379a30f8f8cSSatish Balay 
380a30f8f8cSSatish Balay 
381