xref: /petsc/src/mat/impls/sbaij/mpi/mmsbaij.c (revision 1795a4d16c893ec2fc06bbbc6c5ce592a2de75d4)
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 */
31*1795a4d1SJed Brown   ierr = PetscCalloc1(Nbs,&indices);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 */
40785e854fSJed Brown   ierr = PetscMalloc1(ec,&garray);CHKERRQ(ierr);
41dcca6d9dSJed Brown   ierr = PetscMalloc2(2*ec,&sgarray,ec,&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;
6226fbe8dcSKarl Rupp 
63d0f46423SBarry Smith   sbaij->B->cmap->n = sbaij->B->cmap->N = ec*mat->rmap->bs;
6426fbe8dcSKarl Rupp 
6526283091SBarry Smith   ierr = PetscLayoutSetUp((sbaij->B->cmap));CHKERRQ(ierr);
6640781036SHong Zhang   ierr = PetscFree(indices);CHKERRQ(ierr);
6740781036SHong Zhang 
6840781036SHong Zhang   /* create local vector that is used to scatter into */
6940781036SHong Zhang   ierr = VecCreateSeq(PETSC_COMM_SELF,ec*bs,&sbaij->lvec);CHKERRQ(ierr);
7040781036SHong Zhang 
7140781036SHong Zhang   /* create two temporary index sets for building scatter-gather */
72785e854fSJed Brown   ierr = PetscMalloc1(2*ec,&stmp);CHKERRQ(ierr);
73deff0451SBarry Smith   ierr = ISCreateBlock(PETSC_COMM_SELF,bs,ec,garray,PETSC_COPY_VALUES,&from);CHKERRQ(ierr);
7426fbe8dcSKarl Rupp   for (i=0; i<ec; i++) stmp[i] = i;
75deff0451SBarry Smith   ierr = ISCreateBlock(PETSC_COMM_SELF,bs,ec,stmp,PETSC_COPY_VALUES,&to);CHKERRQ(ierr);
7640781036SHong Zhang 
7793d1592dSHong Zhang   /* generate the scatter context
7893d1592dSHong Zhang      -- Mvctx and lvec are not used by MatMult_MPISBAIJ(), but usefule for some applications */
79ce94432eSBarry Smith   ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)mat),1,mat->cmap->n,mat->cmap->N,NULL,&gvec);CHKERRQ(ierr);
8040781036SHong Zhang   ierr = VecScatterCreate(gvec,from,sbaij->lvec,to,&sbaij->Mvctx);CHKERRQ(ierr);
816bf464f9SBarry Smith   ierr = VecDestroy(&gvec);CHKERRQ(ierr);
8240781036SHong Zhang 
8340781036SHong Zhang   sbaij->garray = garray;
8426fbe8dcSKarl Rupp 
853bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->Mvctx);CHKERRQ(ierr);
863bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->lvec);CHKERRQ(ierr);
873bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)from);CHKERRQ(ierr);
883bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)to);CHKERRQ(ierr);
8940781036SHong Zhang 
906bf464f9SBarry Smith   ierr = ISDestroy(&from);CHKERRQ(ierr);
916bf464f9SBarry Smith   ierr = ISDestroy(&to);CHKERRQ(ierr);
9240781036SHong Zhang 
9340781036SHong Zhang   /* create parallel vector that is used by SBAIJ matrix to scatter from/into */
9440781036SHong Zhang   lsize = (mbs + ec)*bs;
95ce94432eSBarry Smith   ierr  = VecCreateMPI(PetscObjectComm((PetscObject)mat),lsize,PETSC_DETERMINE,&sbaij->slvec0);CHKERRQ(ierr);
9640781036SHong Zhang   ierr  = VecDuplicate(sbaij->slvec0,&sbaij->slvec1);CHKERRQ(ierr);
9740781036SHong Zhang   ierr  = VecGetSize(sbaij->slvec0,&vec_size);CHKERRQ(ierr);
9840781036SHong Zhang 
99077aedafSJed Brown   ierr = VecGetOwnershipRanges(sbaij->slvec0,&sowners);CHKERRQ(ierr);
10040781036SHong Zhang 
10140781036SHong Zhang   /* x index in the IS sfrom */
10240781036SHong Zhang   for (i=0; i<ec; i++) {
10340781036SHong Zhang     j          = ec_owner[i];
10440781036SHong Zhang     sgarray[i] = garray[i] + (sowners[j]/bs - owners[j]);
10540781036SHong Zhang   }
10640781036SHong Zhang   /* b index in the IS sfrom */
10740781036SHong Zhang   k = sowners[rank]/bs + mbs;
10840781036SHong Zhang   for (i=ec,j=0; i< 2*ec; i++,j++) sgarray[i] = k + j;
109deff0451SBarry Smith   ierr = ISCreateBlock(PETSC_COMM_SELF,bs,2*ec,sgarray,PETSC_COPY_VALUES,&from);CHKERRQ(ierr);
11040781036SHong Zhang 
11140781036SHong Zhang   /* x index in the IS sto */
11240781036SHong Zhang   k = sowners[rank]/bs + mbs;
113e82e9f6bSBarry Smith   for (i=0; i<ec; i++) stmp[i] = (k + i);
11440781036SHong Zhang   /* b index in the IS sto */
115e82e9f6bSBarry Smith   for (i=ec; i<2*ec; i++) stmp[i] = sgarray[i-ec];
11640781036SHong Zhang 
117deff0451SBarry Smith   ierr = ISCreateBlock(PETSC_COMM_SELF,bs,2*ec,stmp,PETSC_COPY_VALUES,&to);CHKERRQ(ierr);
11840781036SHong Zhang 
11940781036SHong Zhang   ierr = VecScatterCreate(sbaij->slvec0,from,sbaij->slvec1,to,&sbaij->sMvctx);CHKERRQ(ierr);
12040781036SHong Zhang 
12140781036SHong Zhang   ierr = VecGetLocalSize(sbaij->slvec1,&nt);CHKERRQ(ierr);
1221ebc52fbSHong Zhang   ierr = VecGetArray(sbaij->slvec1,&ptr);CHKERRQ(ierr);
123778a2246SBarry Smith   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,bs*mbs,ptr,&sbaij->slvec1a);CHKERRQ(ierr);
124778a2246SBarry Smith   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,nt-bs*mbs,ptr+bs*mbs,&sbaij->slvec1b);CHKERRQ(ierr);
1251ebc52fbSHong Zhang   ierr = VecRestoreArray(sbaij->slvec1,&ptr);CHKERRQ(ierr);
12640781036SHong Zhang 
1271ebc52fbSHong Zhang   ierr = VecGetArray(sbaij->slvec0,&ptr);CHKERRQ(ierr);
128778a2246SBarry Smith   ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,nt-bs*mbs,ptr+bs*mbs,&sbaij->slvec0b);CHKERRQ(ierr);
1291ebc52fbSHong Zhang   ierr = VecRestoreArray(sbaij->slvec0,&ptr);CHKERRQ(ierr);
13040781036SHong Zhang 
13140781036SHong Zhang   ierr = PetscFree(stmp);CHKERRQ(ierr);
132ce94432eSBarry Smith   ierr = MPI_Barrier(PetscObjectComm((PetscObject)mat));CHKERRQ(ierr);
13340781036SHong Zhang 
1343bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->sMvctx);CHKERRQ(ierr);
1353bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->slvec0);CHKERRQ(ierr);
1363bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->slvec1);CHKERRQ(ierr);
1373bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->slvec0b);CHKERRQ(ierr);
1383bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->slvec1a);CHKERRQ(ierr);
1393bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->slvec1b);CHKERRQ(ierr);
1403bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)from);CHKERRQ(ierr);
1413bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)to);CHKERRQ(ierr);
14240781036SHong Zhang 
1433bb1ff40SBarry Smith   ierr = PetscLogObjectMemory((PetscObject)mat,(ec+1)*sizeof(PetscInt));CHKERRQ(ierr);
1446bf464f9SBarry Smith   ierr = ISDestroy(&from);CHKERRQ(ierr);
1456bf464f9SBarry Smith   ierr = ISDestroy(&to);CHKERRQ(ierr);
14674ed9c26SBarry Smith   ierr = PetscFree2(sgarray,ec_owner);CHKERRQ(ierr);
14740781036SHong Zhang   PetscFunctionReturn(0);
14840781036SHong Zhang }
14940781036SHong Zhang 
15040781036SHong Zhang #undef __FUNCT__
15140781036SHong Zhang #define __FUNCT__ "MatSetUpMultiply_MPISBAIJ_2comm"
152dfbe8321SBarry Smith PetscErrorCode MatSetUpMultiply_MPISBAIJ_2comm(Mat mat)
15340781036SHong Zhang {
1542896befeSSatish Balay   Mat_MPISBAIJ   *baij = (Mat_MPISBAIJ*)mat->data;
155a30f8f8cSSatish Balay   Mat_SeqBAIJ    *B    = (Mat_SeqBAIJ*)(baij->B->data);
1566849ba73SBarry Smith   PetscErrorCode ierr;
15713f74950SBarry Smith   PetscInt       i,j,*aj = B->j,ec = 0,*garray;
158d0f46423SBarry Smith   PetscInt       bs = mat->rmap->bs,*stmp;
159a30f8f8cSSatish Balay   IS             from,to;
160a30f8f8cSSatish Balay   Vec            gvec;
161a30f8f8cSSatish Balay #if defined(PETSC_USE_CTABLE)
162a30f8f8cSSatish Balay   PetscTable         gid1_lid1;
163a30f8f8cSSatish Balay   PetscTablePosition tpos;
16413f74950SBarry Smith   PetscInt           gid,lid;
1656f531f54SSatish Balay #else
16613f74950SBarry Smith   PetscInt Nbs = baij->Nbs,*indices;
167a30f8f8cSSatish Balay #endif
168a30f8f8cSSatish Balay 
169a30f8f8cSSatish Balay   PetscFunctionBegin;
170a30f8f8cSSatish Balay #if defined(PETSC_USE_CTABLE)
171a30f8f8cSSatish Balay   /* use a table - Mark Adams */
172e23dfa41SBarry Smith   PetscTableCreate(B->mbs,baij->Nbs+1,&gid1_lid1);
173a30f8f8cSSatish Balay   for (i=0; i<B->mbs; i++) {
174a30f8f8cSSatish Balay     for (j=0; j<B->ilen[i]; j++) {
17513f74950SBarry Smith       PetscInt data,gid1 = aj[B->i[i]+j] + 1;
176a30f8f8cSSatish Balay       ierr = PetscTableFind(gid1_lid1,gid1,&data);CHKERRQ(ierr);
177a30f8f8cSSatish Balay       if (!data) {
178a30f8f8cSSatish Balay         /* one based table */
1793861aac3SJed Brown         ierr = PetscTableAdd(gid1_lid1,gid1,++ec,INSERT_VALUES);CHKERRQ(ierr);
180a30f8f8cSSatish Balay       }
181a30f8f8cSSatish Balay     }
182a30f8f8cSSatish Balay   }
183a30f8f8cSSatish Balay   /* form array of columns we need */
184785e854fSJed Brown   ierr = PetscMalloc1(ec,&garray);CHKERRQ(ierr);
185a30f8f8cSSatish Balay   ierr = PetscTableGetHeadPosition(gid1_lid1,&tpos);CHKERRQ(ierr);
186a30f8f8cSSatish Balay   while (tpos) {
187a30f8f8cSSatish Balay     ierr = PetscTableGetNext(gid1_lid1,&tpos,&gid,&lid);CHKERRQ(ierr);
188a30f8f8cSSatish Balay     gid--; lid--;
189a30f8f8cSSatish Balay     garray[lid] = gid;
190a30f8f8cSSatish Balay   }
191a30f8f8cSSatish Balay   ierr = PetscSortInt(ec,garray);CHKERRQ(ierr);
192a30f8f8cSSatish Balay   ierr = PetscTableRemoveAll(gid1_lid1);CHKERRQ(ierr);
193a30f8f8cSSatish Balay   for (i=0; i<ec; i++) {
1943861aac3SJed Brown     ierr = PetscTableAdd(gid1_lid1,garray[i]+1,i+1,INSERT_VALUES);CHKERRQ(ierr);
195a30f8f8cSSatish Balay   }
196a30f8f8cSSatish Balay   /* compact out the extra columns in B */
197a30f8f8cSSatish Balay   for (i=0; i<B->mbs; i++) {
198a30f8f8cSSatish Balay     for (j=0; j<B->ilen[i]; j++) {
19913f74950SBarry Smith       PetscInt gid1 = aj[B->i[i] + j] + 1;
200a30f8f8cSSatish Balay       ierr = PetscTableFind(gid1_lid1,gid1,&lid);CHKERRQ(ierr);
201a30f8f8cSSatish Balay       lid--;
202a30f8f8cSSatish Balay       aj[B->i[i]+j] = lid;
203a30f8f8cSSatish Balay     }
204a30f8f8cSSatish Balay   }
205a30f8f8cSSatish Balay   B->nbs = ec;
20626fbe8dcSKarl Rupp 
207d0f46423SBarry Smith   baij->B->cmap->n = baij->B->cmap->N = ec*mat->rmap->bs;
20826fbe8dcSKarl Rupp 
20926283091SBarry Smith   ierr = PetscLayoutSetUp((baij->B->cmap));CHKERRQ(ierr);
2106bc0bbbfSBarry Smith   ierr = PetscTableDestroy(&gid1_lid1);CHKERRQ(ierr);
211a30f8f8cSSatish Balay #else
212a30f8f8cSSatish Balay   /* For the first stab we make an array as long as the number of columns */
213a30f8f8cSSatish Balay   /* mark those columns that are in baij->B */
214*1795a4d1SJed Brown   ierr = PetscCalloc1(Nbs,&indices);CHKERRQ(ierr);
215a30f8f8cSSatish Balay   for (i=0; i<B->mbs; i++) {
216a30f8f8cSSatish Balay     for (j=0; j<B->ilen[i]; j++) {
217a30f8f8cSSatish Balay       if (!indices[aj[B->i[i] + j]]) ec++;
218a30f8f8cSSatish Balay       indices[aj[B->i[i] + j]] = 1;
219a30f8f8cSSatish Balay     }
220a30f8f8cSSatish Balay   }
221a30f8f8cSSatish Balay 
222a30f8f8cSSatish Balay   /* form array of columns we need */
223785e854fSJed Brown   ierr = PetscMalloc1(ec,&garray);CHKERRQ(ierr);
224f3ef73ceSBarry Smith   ec = 0;
225f3ef73ceSBarry Smith   for (i=0; i<Nbs; i++) {
226f3ef73ceSBarry Smith     if (indices[i]) {
227f3ef73ceSBarry Smith       garray[ec++] = i;
228f3ef73ceSBarry Smith     }
229f3ef73ceSBarry Smith   }
230a30f8f8cSSatish Balay 
231a30f8f8cSSatish Balay   /* make indices now point into garray */
23226fbe8dcSKarl Rupp   for (i=0; i<ec; i++) indices[garray[i]] = i;
233a30f8f8cSSatish Balay 
234a30f8f8cSSatish Balay   /* compact out the extra columns in B */
235a30f8f8cSSatish Balay   for (i=0; i<B->mbs; i++) {
236a30f8f8cSSatish Balay     for (j=0; j<B->ilen[i]; j++) {
237a30f8f8cSSatish Balay       aj[B->i[i] + j] = indices[aj[B->i[i] + j]];
238a30f8f8cSSatish Balay     }
239a30f8f8cSSatish Balay   }
240a30f8f8cSSatish Balay   B->nbs = ec;
24126fbe8dcSKarl Rupp 
242d0f46423SBarry Smith   baij->B->cmap->n = ec*mat->rmap->bs;
24326fbe8dcSKarl Rupp 
244a30f8f8cSSatish Balay   ierr = PetscFree(indices);CHKERRQ(ierr);
245a30f8f8cSSatish Balay #endif
246633e10c7SHong Zhang 
247a30f8f8cSSatish Balay   /* create local vector that is used to scatter into */
248a30f8f8cSSatish Balay   ierr = VecCreateSeq(PETSC_COMM_SELF,ec*bs,&baij->lvec);CHKERRQ(ierr);
249a30f8f8cSSatish Balay 
250a30f8f8cSSatish Balay   /* create two temporary index sets for building scatter-gather */
251deff0451SBarry Smith   ierr = ISCreateBlock(PETSC_COMM_SELF,bs,ec,garray,PETSC_COPY_VALUES,&from);CHKERRQ(ierr);
252a30f8f8cSSatish Balay 
253785e854fSJed Brown   ierr = PetscMalloc1(ec,&stmp);CHKERRQ(ierr);
25426fbe8dcSKarl Rupp   for (i=0; i<ec; i++) stmp[i] = i;
255deff0451SBarry Smith   ierr = ISCreateBlock(PETSC_COMM_SELF,bs,ec,stmp,PETSC_OWN_POINTER,&to);CHKERRQ(ierr);
256a30f8f8cSSatish Balay 
257a30f8f8cSSatish Balay   /* create temporary global vector to generate scatter context */
258ce94432eSBarry Smith   ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)mat),1,mat->cmap->n,mat->cmap->N,NULL,&gvec);CHKERRQ(ierr);
259a30f8f8cSSatish Balay   ierr = VecScatterCreate(gvec,from,baij->lvec,to,&baij->Mvctx);CHKERRQ(ierr);
2606bf464f9SBarry Smith   ierr = VecDestroy(&gvec);CHKERRQ(ierr);
261a30f8f8cSSatish Balay 
2623bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)baij->Mvctx);CHKERRQ(ierr);
2633bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)baij->lvec);CHKERRQ(ierr);
2643bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)from);CHKERRQ(ierr);
2653bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)mat,(PetscObject)to);CHKERRQ(ierr);
26626fbe8dcSKarl Rupp 
267a30f8f8cSSatish Balay   baij->garray = garray;
26826fbe8dcSKarl Rupp 
2693bb1ff40SBarry Smith   ierr = PetscLogObjectMemory((PetscObject)mat,(ec+1)*sizeof(PetscInt));CHKERRQ(ierr);
2706bf464f9SBarry Smith   ierr = ISDestroy(&from);CHKERRQ(ierr);
2716bf464f9SBarry Smith   ierr = ISDestroy(&to);CHKERRQ(ierr);
272a30f8f8cSSatish Balay   PetscFunctionReturn(0);
273a30f8f8cSSatish Balay }
274a30f8f8cSSatish Balay 
275a30f8f8cSSatish Balay /*
27601b2bd88SHong Zhang      Takes the local part of an already assembled MPISBAIJ matrix
277a30f8f8cSSatish Balay    and disassembles it. This is to allow new nonzeros into the matrix
278a30f8f8cSSatish Balay    that require more communication in the matrix vector multiply.
279a30f8f8cSSatish Balay    Thus certain data-structures must be rebuilt.
280a30f8f8cSSatish Balay 
281a30f8f8cSSatish Balay    Kind of slow! But that's what application programmers get when
282a30f8f8cSSatish Balay    they are sloppy.
283a30f8f8cSSatish Balay */
2844a2ae208SSatish Balay #undef __FUNCT__
285ab9863d7SBarry Smith #define __FUNCT__ "MatDisAssemble_MPISBAIJ"
286ab9863d7SBarry Smith PetscErrorCode MatDisAssemble_MPISBAIJ(Mat A)
287a30f8f8cSSatish Balay {
2882896befeSSatish Balay   Mat_MPISBAIJ   *baij  = (Mat_MPISBAIJ*)A->data;
289a30f8f8cSSatish Balay   Mat            B      = baij->B,Bnew;
290a30f8f8cSSatish Balay   Mat_SeqBAIJ    *Bbaij = (Mat_SeqBAIJ*)B->data;
291dfbe8321SBarry Smith   PetscErrorCode ierr;
292d0f46423SBarry Smith   PetscInt       i,j,mbs=Bbaij->mbs,n = A->cmap->N,col,*garray=baij->garray;
293d0f46423SBarry Smith   PetscInt       k,bs=A->rmap->bs,bs2=baij->bs2,*rvals,*nz,ec,m=A->rmap->n;
294a30f8f8cSSatish Balay   MatScalar      *a = Bbaij->a;
29587828ca2SBarry Smith   PetscScalar    *atmp;
296ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL_MAT_SINGLE)
29713f74950SBarry Smith   PetscInt l;
298a30f8f8cSSatish Balay #endif
299a30f8f8cSSatish Balay 
300a30f8f8cSSatish Balay   PetscFunctionBegin;
301ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL_MAT_SINGLE)
302785e854fSJed Brown   ierr = PetscMalloc1(A->rmap->bs,&atmp);
30382502324SSatish Balay #endif
304a30f8f8cSSatish Balay   /* free stuff related to matrix-vec multiply */
305b0a32e0cSBarry Smith   ierr = VecGetSize(baij->lvec,&ec);CHKERRQ(ierr); /* needed for PetscLogObjectMemory below */
3066bf464f9SBarry Smith   ierr = VecDestroy(&baij->lvec);CHKERRQ(ierr);
3076bf464f9SBarry Smith   ierr = VecScatterDestroy(&baij->Mvctx);CHKERRQ(ierr);
30801b2bd88SHong Zhang 
3096bf464f9SBarry Smith   ierr = VecDestroy(&baij->slvec0);CHKERRQ(ierr);
3106bf464f9SBarry Smith   ierr = VecDestroy(&baij->slvec0b);CHKERRQ(ierr);
3116bf464f9SBarry Smith   ierr = VecDestroy(&baij->slvec1);CHKERRQ(ierr);
3126bf464f9SBarry Smith   ierr = VecDestroy(&baij->slvec1a);CHKERRQ(ierr);
3136bf464f9SBarry Smith   ierr = VecDestroy(&baij->slvec1b);CHKERRQ(ierr);
31401b2bd88SHong Zhang 
315a30f8f8cSSatish Balay   if (baij->colmap) {
316a30f8f8cSSatish Balay #if defined(PETSC_USE_CTABLE)
3176bc0bbbfSBarry Smith     ierr = PetscTableDestroy(&baij->colmap);CHKERRQ(ierr);
318a30f8f8cSSatish Balay #else
319a30f8f8cSSatish Balay     ierr = PetscFree(baij->colmap);CHKERRQ(ierr);
3203bb1ff40SBarry Smith     ierr = PetscLogObjectMemory((PetscObject)A,-Bbaij->nbs*sizeof(PetscInt));CHKERRQ(ierr);
321a30f8f8cSSatish Balay #endif
322a30f8f8cSSatish Balay   }
323a30f8f8cSSatish Balay 
324a30f8f8cSSatish Balay   /* make sure that B is assembled so we can access its values */
325a30f8f8cSSatish Balay   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
326a30f8f8cSSatish Balay   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
327a30f8f8cSSatish Balay 
328a30f8f8cSSatish Balay   /* invent new B and copy stuff over */
329785e854fSJed Brown   ierr = PetscMalloc1(mbs,&nz);CHKERRQ(ierr);
330a30f8f8cSSatish Balay   for (i=0; i<mbs; i++) {
331a30f8f8cSSatish Balay     nz[i] = Bbaij->i[i+1]-Bbaij->i[i];
332a30f8f8cSSatish Balay   }
333f69a0ea3SMatthew Knepley   ierr = MatCreate(PETSC_COMM_SELF,&Bnew);CHKERRQ(ierr);
334f69a0ea3SMatthew Knepley   ierr = MatSetSizes(Bnew,m,n,m,n);CHKERRQ(ierr);
3357adad957SLisandro Dalcin   ierr = MatSetType(Bnew,((PetscObject)B)->type_name);CHKERRQ(ierr);
336d0f46423SBarry Smith   ierr = MatSeqBAIJSetPreallocation(Bnew,B->rmap->bs,0,nz);CHKERRQ(ierr);
33726fbe8dcSKarl Rupp 
3382576faa2SJed Brown   ((Mat_SeqSBAIJ*)Bnew->data)->nonew = Bbaij->nonew; /* Inherit insertion error options. */
33926fbe8dcSKarl Rupp 
340a30f8f8cSSatish Balay   ierr = PetscFree(nz);CHKERRQ(ierr);
341a30f8f8cSSatish Balay 
342785e854fSJed Brown   ierr = PetscMalloc1(bs,&rvals);CHKERRQ(ierr);
343a30f8f8cSSatish Balay   for (i=0; i<mbs; i++) {
344a30f8f8cSSatish Balay     rvals[0] = bs*i;
34526fbe8dcSKarl Rupp     for (j=1; j<bs; j++) rvals[j] = rvals[j-1] + 1;
346a30f8f8cSSatish Balay     for (j=Bbaij->i[i]; j<Bbaij->i[i+1]; j++) {
347a30f8f8cSSatish Balay       col = garray[Bbaij->j[j]]*bs;
348a30f8f8cSSatish Balay       for (k=0; k<bs; k++) {
349ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL_MAT_SINGLE)
350a30f8f8cSSatish Balay         for (l=0; l<bs; l++) atmp[l] = a[j*bs2+l];
351a30f8f8cSSatish Balay #else
35271730473SSatish Balay         atmp = a+j*bs2 + k*bs;
353a30f8f8cSSatish Balay #endif
354c8407628SSatish Balay         ierr = MatSetValues_SeqSBAIJ(Bnew,bs,rvals,1,&col,atmp,B->insertmode);CHKERRQ(ierr);
355a30f8f8cSSatish Balay         col++;
356a30f8f8cSSatish Balay       }
357a30f8f8cSSatish Balay     }
358a30f8f8cSSatish Balay   }
359ce63c4c1SBarry Smith #if defined(PETSC_USE_REAL_MAT_SINGLE)
360a30f8f8cSSatish Balay   ierr = PetscFree(atmp);CHKERRQ(ierr);
361a30f8f8cSSatish Balay #endif
362a30f8f8cSSatish Balay   ierr = PetscFree(baij->garray);CHKERRQ(ierr);
36326fbe8dcSKarl Rupp 
364a30f8f8cSSatish Balay   baij->garray = 0;
36526fbe8dcSKarl Rupp 
366a30f8f8cSSatish Balay   ierr = PetscFree(rvals);CHKERRQ(ierr);
3673bb1ff40SBarry Smith   ierr = PetscLogObjectMemory((PetscObject)A,-ec*sizeof(PetscInt));CHKERRQ(ierr);
3686bf464f9SBarry Smith   ierr = MatDestroy(&B);CHKERRQ(ierr);
3693bb1ff40SBarry Smith   ierr = PetscLogObjectParent((PetscObject)A,(PetscObject)Bnew);CHKERRQ(ierr);
37026fbe8dcSKarl Rupp 
371a30f8f8cSSatish Balay   baij->B = Bnew;
37226fbe8dcSKarl Rupp 
373a30f8f8cSSatish Balay   A->was_assembled = PETSC_FALSE;
374a30f8f8cSSatish Balay   PetscFunctionReturn(0);
375a30f8f8cSSatish Balay }
376a30f8f8cSSatish Balay 
377a30f8f8cSSatish Balay 
378a30f8f8cSSatish Balay 
379