Actual source code: mmsbaij.c
2: /*
3: Support for the parallel SBAIJ matrix vector multiply
4: */
5: #include <../src/mat/impls/sbaij/mpi/mpisbaij.h>
12: PetscErrorCode MatSetUpMultiply_MPISBAIJ(Mat mat)
13: {
14: Mat_MPISBAIJ *sbaij = (Mat_MPISBAIJ*)mat->data;
15: Mat_SeqBAIJ *B = (Mat_SeqBAIJ*)(sbaij->B->data);
17: PetscInt Nbs = sbaij->Nbs,i,j,*indices,*aj = B->j,ec = 0,*garray,*sgarray;
18: PetscInt bs = mat->rmap->bs,*stmp,mbs=sbaij->mbs, vec_size,nt;
19: IS from,to;
20: Vec gvec;
21: PetscMPIInt rank=sbaij->rank,lsize,size=sbaij->size;
22: PetscInt *owners=sbaij->rangebs,*sowners,*ec_owner,k;
23: PetscScalar *ptr;
26: VecScatterDestroy(&sbaij->sMvctx);
27:
28: /* For the first stab we make an array as long as the number of columns */
29: /* mark those columns that are in sbaij->B */
30: PetscMalloc(Nbs*sizeof(PetscInt),&indices);
31: PetscMemzero(indices,Nbs*sizeof(PetscInt));
32: for (i=0; i<mbs; i++) {
33: for (j=0; j<B->ilen[i]; j++) {
34: if (!indices[aj[B->i[i] + j]]) ec++;
35: indices[aj[B->i[i] + j] ] = 1;
36: }
37: }
39: /* form arrays of columns we need */
40: PetscMalloc(ec*sizeof(PetscInt),&garray);
41: PetscMalloc2(2*ec,PetscInt,&sgarray,ec,PetscInt,&ec_owner);
42:
43: ec = 0;
44: for (j=0; j<size; j++){
45: for (i=owners[j]; i<owners[j+1]; i++){
46: if (indices[i]) {
47: garray[ec] = i;
48: ec_owner[ec] = j;
49: ec++;
50: }
51: }
52: }
54: /* make indices now point into garray */
55: for (i=0; i<ec; i++) indices[garray[i]] = i;
57: /* compact out the extra columns in B */
58: for (i=0; i<mbs; i++) {
59: for (j=0; j<B->ilen[i]; j++) aj[B->i[i] + j] = indices[aj[B->i[i] + j]];
60: }
61: B->nbs = ec;
62: sbaij->B->cmap->n = sbaij->B->cmap->N = ec*mat->rmap->bs;
63: PetscLayoutSetUp((sbaij->B->cmap));
64: PetscFree(indices);
66: /* create local vector that is used to scatter into */
67: VecCreateSeq(PETSC_COMM_SELF,ec*bs,&sbaij->lvec);
69: /* create two temporary index sets for building scatter-gather */
70: PetscMalloc(2*ec*sizeof(PetscInt),&stmp);
71: ISCreateBlock(PETSC_COMM_SELF,bs,ec,garray,PETSC_COPY_VALUES,&from);
72: for (i=0; i<ec; i++) { stmp[i] = i; }
73: ISCreateBlock(PETSC_COMM_SELF,bs,ec,stmp,PETSC_COPY_VALUES,&to);
75: /* generate the scatter context
76: -- Mvctx and lvec are not used by MatMult_MPISBAIJ(), but usefule for some applications */
77: VecCreateMPIWithArray(((PetscObject)mat)->comm,mat->cmap->n,mat->cmap->N,PETSC_NULL,&gvec);
78: VecScatterCreate(gvec,from,sbaij->lvec,to,&sbaij->Mvctx);
79: VecDestroy(&gvec);
81: sbaij->garray = garray;
82: PetscLogObjectParent(mat,sbaij->Mvctx);
83: PetscLogObjectParent(mat,sbaij->lvec);
84: PetscLogObjectParent(mat,from);
85: PetscLogObjectParent(mat,to);
87: ISDestroy(&from);
88: ISDestroy(&to);
90: /* create parallel vector that is used by SBAIJ matrix to scatter from/into */
91: lsize = (mbs + ec)*bs;
92: VecCreateMPI(((PetscObject)mat)->comm,lsize,PETSC_DETERMINE,&sbaij->slvec0);
93: VecDuplicate(sbaij->slvec0,&sbaij->slvec1);
94: VecGetSize(sbaij->slvec0,&vec_size);
96: sowners = sbaij->slvec0->map->range;
97:
98: /* x index in the IS sfrom */
99: for (i=0; i<ec; i++) {
100: j = ec_owner[i];
101: sgarray[i] = garray[i] + (sowners[j]/bs - owners[j]);
102: }
103: /* b index in the IS sfrom */
104: k = sowners[rank]/bs + mbs;
105: for (i=ec,j=0; i< 2*ec; i++,j++) sgarray[i] = k + j;
106: ISCreateBlock(PETSC_COMM_SELF,bs,2*ec,sgarray,PETSC_COPY_VALUES,&from);
107:
108: /* x index in the IS sto */
109: k = sowners[rank]/bs + mbs;
110: for (i=0; i<ec; i++) stmp[i] = (k + i);
111: /* b index in the IS sto */
112: for (i=ec; i<2*ec; i++) stmp[i] = sgarray[i-ec];
114: ISCreateBlock(PETSC_COMM_SELF,bs,2*ec,stmp,PETSC_COPY_VALUES,&to);
116: VecScatterCreate(sbaij->slvec0,from,sbaij->slvec1,to,&sbaij->sMvctx);
117:
118: VecGetLocalSize(sbaij->slvec1,&nt);
119: VecGetArray(sbaij->slvec1,&ptr);
120: VecCreateSeqWithArray(PETSC_COMM_SELF,bs*mbs,ptr,&sbaij->slvec1a);
121: VecCreateSeqWithArray(PETSC_COMM_SELF,nt-bs*mbs,ptr+bs*mbs,&sbaij->slvec1b);
122: VecRestoreArray(sbaij->slvec1,&ptr);
124: VecGetArray(sbaij->slvec0,&ptr);
125: VecCreateSeqWithArray(PETSC_COMM_SELF,nt-bs*mbs,ptr+bs*mbs,&sbaij->slvec0b);
126: VecRestoreArray(sbaij->slvec0,&ptr);
128: PetscFree(stmp);
129: MPI_Barrier(((PetscObject)mat)->comm);
130:
131: PetscLogObjectParent(mat,sbaij->sMvctx);
132: PetscLogObjectParent(mat,sbaij->slvec0);
133: PetscLogObjectParent(mat,sbaij->slvec1);
134: PetscLogObjectParent(mat,sbaij->slvec0b);
135: PetscLogObjectParent(mat,sbaij->slvec1a);
136: PetscLogObjectParent(mat,sbaij->slvec1b);
137: PetscLogObjectParent(mat,from);
138: PetscLogObjectParent(mat,to);
139:
140: PetscLogObjectMemory(mat,(ec+1)*sizeof(PetscInt));
141: ISDestroy(&from);
142: ISDestroy(&to);
143: PetscFree2(sgarray,ec_owner);
144: return(0);
145: }
149: PetscErrorCode MatSetUpMultiply_MPISBAIJ_2comm(Mat mat)
150: {
151: Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)mat->data;
152: Mat_SeqBAIJ *B = (Mat_SeqBAIJ*)(baij->B->data);
153: PetscErrorCode ierr;
154: PetscInt i,j,*aj = B->j,ec = 0,*garray;
155: PetscInt bs = mat->rmap->bs,*stmp;
156: IS from,to;
157: Vec gvec;
158: #if defined (PETSC_USE_CTABLE)
159: PetscTable gid1_lid1;
160: PetscTablePosition tpos;
161: PetscInt gid,lid;
162: #else
163: PetscInt Nbs = baij->Nbs,*indices;
164: #endif
167: #if defined (PETSC_USE_CTABLE)
168: /* use a table - Mark Adams */
169: PetscTableCreate(B->mbs,&gid1_lid1);
170: for (i=0; i<B->mbs; i++) {
171: for (j=0; j<B->ilen[i]; j++) {
172: PetscInt data,gid1 = aj[B->i[i]+j] + 1;
173: PetscTableFind(gid1_lid1,gid1,&data);
174: if (!data) {
175: /* one based table */
176: PetscTableAdd(gid1_lid1,gid1,++ec);
177: }
178: }
179: }
180: /* form array of columns we need */
181: PetscMalloc(ec*sizeof(PetscInt),&garray);
182: PetscTableGetHeadPosition(gid1_lid1,&tpos);
183: while (tpos) {
184: PetscTableGetNext(gid1_lid1,&tpos,&gid,&lid);
185: gid--; lid--;
186: garray[lid] = gid;
187: }
188: PetscSortInt(ec,garray);
189: PetscTableRemoveAll(gid1_lid1);
190: for (i=0; i<ec; i++) {
191: PetscTableAdd(gid1_lid1,garray[i]+1,i+1);
192: }
193: /* compact out the extra columns in B */
194: for (i=0; i<B->mbs; i++) {
195: for (j=0; j<B->ilen[i]; j++) {
196: PetscInt gid1 = aj[B->i[i] + j] + 1;
197: PetscTableFind(gid1_lid1,gid1,&lid);
198: lid --;
199: aj[B->i[i]+j] = lid;
200: }
201: }
202: B->nbs = ec;
203: baij->B->cmap->n = baij->B->cmap->N = ec*mat->rmap->bs;
204: PetscLayoutSetUp((baij->B->cmap));
205: PetscTableDestroy(&gid1_lid1);
206: #else
207: /* For the first stab we make an array as long as the number of columns */
208: /* mark those columns that are in baij->B */
209: PetscMalloc(Nbs*sizeof(PetscInt),&indices);
210: PetscMemzero(indices,Nbs*sizeof(PetscInt));
211: for (i=0; i<B->mbs; i++) {
212: for (j=0; j<B->ilen[i]; j++) {
213: if (!indices[aj[B->i[i] + j]]) ec++;
214: indices[aj[B->i[i] + j] ] = 1;
215: }
216: }
218: /* form array of columns we need */
219: PetscMalloc(ec*sizeof(PetscInt),&garray);
220: ec = 0;
221: for (i=0; i<Nbs; i++) {
222: if (indices[i]) {
223: garray[ec++] = i;
224: }
225: }
227: /* make indices now point into garray */
228: for (i=0; i<ec; i++) {
229: indices[garray[i]] = i;
230: }
232: /* compact out the extra columns in B */
233: for (i=0; i<B->mbs; i++) {
234: for (j=0; j<B->ilen[i]; j++) {
235: aj[B->i[i] + j] = indices[aj[B->i[i] + j]];
236: }
237: }
238: B->nbs = ec;
239: baij->B->cmap->n = ec*mat->rmap->bs;
240: PetscFree(indices);
241: #endif
242:
243: /* create local vector that is used to scatter into */
244: VecCreateSeq(PETSC_COMM_SELF,ec*bs,&baij->lvec);
246: /* create two temporary index sets for building scatter-gather */
247: ISCreateBlock(PETSC_COMM_SELF,bs,ec,garray,PETSC_COPY_VALUES,&from);
249: PetscMalloc(ec*sizeof(PetscInt),&stmp);
250: for (i=0; i<ec; i++) { stmp[i] = i; }
251: ISCreateBlock(PETSC_COMM_SELF,bs,ec,stmp,PETSC_OWN_POINTER,&to);
253: /* create temporary global vector to generate scatter context */
254: VecCreateMPIWithArray(((PetscObject)mat)->comm,mat->cmap->n,mat->cmap->N,PETSC_NULL,&gvec);
255: VecScatterCreate(gvec,from,baij->lvec,to,&baij->Mvctx);
256: VecDestroy(&gvec);
258: PetscLogObjectParent(mat,baij->Mvctx);
259: PetscLogObjectParent(mat,baij->lvec);
260: PetscLogObjectParent(mat,from);
261: PetscLogObjectParent(mat,to);
262: baij->garray = garray;
263: PetscLogObjectMemory(mat,(ec+1)*sizeof(PetscInt));
264: ISDestroy(&from);
265: ISDestroy(&to);
266: return(0);
267: }
269: /*
270: Takes the local part of an already assembled MPISBAIJ matrix
271: and disassembles it. This is to allow new nonzeros into the matrix
272: that require more communication in the matrix vector multiply.
273: Thus certain data-structures must be rebuilt.
275: Kind of slow! But that's what application programmers get when
276: they are sloppy.
277: */
280: PetscErrorCode DisAssemble_MPISBAIJ(Mat A)
281: {
282: Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)A->data;
283: Mat B = baij->B,Bnew;
284: Mat_SeqBAIJ *Bbaij = (Mat_SeqBAIJ*)B->data;
286: PetscInt i,j,mbs=Bbaij->mbs,n = A->cmap->N,col,*garray=baij->garray;
287: PetscInt k,bs=A->rmap->bs,bs2=baij->bs2,*rvals,*nz,ec,m=A->rmap->n;
288: MatScalar *a = Bbaij->a;
289: PetscScalar *atmp;
290: #if defined(PETSC_USE_REAL_MAT_SINGLE)
291: PetscInt l;
292: #endif
295: #if defined(PETSC_USE_REAL_MAT_SINGLE)
296: PetscMalloc(A->rmap->bs*sizeof(PetscScalar),&atmp);
297: #endif
298: /* free stuff related to matrix-vec multiply */
299: VecGetSize(baij->lvec,&ec); /* needed for PetscLogObjectMemory below */
300: VecDestroy(&baij->lvec);
301: VecScatterDestroy(&baij->Mvctx);
303: VecDestroy(&baij->slvec0);
304: VecDestroy(&baij->slvec0b);
305: VecDestroy(&baij->slvec1);
306: VecDestroy(&baij->slvec1a);
307: VecDestroy(&baij->slvec1b);
309: if (baij->colmap) {
310: #if defined (PETSC_USE_CTABLE)
311: PetscTableDestroy(&baij->colmap);
312: #else
313: PetscFree(baij->colmap);
314: PetscLogObjectMemory(A,-Bbaij->nbs*sizeof(PetscInt));
315: #endif
316: }
318: /* make sure that B is assembled so we can access its values */
319: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
320: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
322: /* invent new B and copy stuff over */
323: PetscMalloc(mbs*sizeof(PetscInt),&nz);
324: for (i=0; i<mbs; i++) {
325: nz[i] = Bbaij->i[i+1]-Bbaij->i[i];
326: }
327: MatCreate(PETSC_COMM_SELF,&Bnew);
328: MatSetSizes(Bnew,m,n,m,n);
329: MatSetType(Bnew,((PetscObject)B)->type_name);
330: MatSeqBAIJSetPreallocation(Bnew,B->rmap->bs,0,nz);
331: PetscFree(nz);
332:
333: PetscMalloc(bs*sizeof(PetscInt),&rvals);
334: for (i=0; i<mbs; i++) {
335: rvals[0] = bs*i;
336: for (j=1; j<bs; j++) { rvals[j] = rvals[j-1] + 1; }
337: for (j=Bbaij->i[i]; j<Bbaij->i[i+1]; j++) {
338: col = garray[Bbaij->j[j]]*bs;
339: for (k=0; k<bs; k++) {
340: #if defined(PETSC_USE_REAL_MAT_SINGLE)
341: for (l=0; l<bs; l++) atmp[l] = a[j*bs2+l];
342: #else
343: atmp = a+j*bs2 + k*bs;
344: #endif
345: MatSetValues_SeqSBAIJ(Bnew,bs,rvals,1,&col,atmp,B->insertmode);
346: col++;
347: }
348: }
349: }
350: #if defined(PETSC_USE_REAL_MAT_SINGLE)
351: PetscFree(atmp);
352: #endif
353: PetscFree(baij->garray);
354: baij->garray = 0;
355: PetscFree(rvals);
356: PetscLogObjectMemory(A,-ec*sizeof(PetscInt));
357: MatDestroy(&B);
358: PetscLogObjectParent(A,Bnew);
359: baij->B = Bnew;
360: A->was_assembled = PETSC_FALSE;
361: return(0);
362: }