Actual source code: mmsbaij.c
petsc-3.7.3 2016-08-01
2: /*
3: Support for the parallel SBAIJ matrix vector multiply
4: */
5: #include <../src/mat/impls/sbaij/mpi/mpisbaij.h>
7: extern PetscErrorCode MatSetValues_SeqSBAIJ(Mat,PetscInt,const PetscInt [],PetscInt,const PetscInt [],const PetscScalar [],InsertMode);
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,*ec_owner,k;
23: const PetscInt *sowners;
24: PetscScalar *ptr;
27: VecScatterDestroy(&sbaij->sMvctx);
29: /* For the first stab we make an array as long as the number of columns */
30: /* mark those columns that are in sbaij->B */
31: PetscCalloc1(Nbs,&indices);
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: PetscMalloc1(ec,&garray);
41: PetscMalloc2(2*ec,&sgarray,ec,&ec_owner);
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;
63: sbaij->B->cmap->n = sbaij->B->cmap->N = ec*mat->rmap->bs;
65: PetscLayoutSetUp((sbaij->B->cmap));
66: PetscFree(indices);
68: /* create local vector that is used to scatter into */
69: VecCreateSeq(PETSC_COMM_SELF,ec*bs,&sbaij->lvec);
71: /* create two temporary index sets for building scatter-gather */
72: PetscMalloc1(2*ec,&stmp);
73: ISCreateBlock(PETSC_COMM_SELF,bs,ec,garray,PETSC_COPY_VALUES,&from);
74: for (i=0; i<ec; i++) stmp[i] = i;
75: ISCreateBlock(PETSC_COMM_SELF,bs,ec,stmp,PETSC_COPY_VALUES,&to);
77: /* generate the scatter context
78: -- Mvctx and lvec are not used by MatMult_MPISBAIJ(), but usefule for some applications */
79: VecCreateMPIWithArray(PetscObjectComm((PetscObject)mat),1,mat->cmap->n,mat->cmap->N,NULL,&gvec);
80: VecScatterCreate(gvec,from,sbaij->lvec,to,&sbaij->Mvctx);
81: VecDestroy(&gvec);
83: sbaij->garray = garray;
85: PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->Mvctx);
86: PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->lvec);
87: PetscLogObjectParent((PetscObject)mat,(PetscObject)from);
88: PetscLogObjectParent((PetscObject)mat,(PetscObject)to);
90: ISDestroy(&from);
91: ISDestroy(&to);
93: /* create parallel vector that is used by SBAIJ matrix to scatter from/into */
94: lsize = (mbs + ec)*bs;
95: VecCreateMPI(PetscObjectComm((PetscObject)mat),lsize,PETSC_DETERMINE,&sbaij->slvec0);
96: VecDuplicate(sbaij->slvec0,&sbaij->slvec1);
97: VecGetSize(sbaij->slvec0,&vec_size);
99: VecGetOwnershipRanges(sbaij->slvec0,&sowners);
101: /* x index in the IS sfrom */
102: for (i=0; i<ec; i++) {
103: j = ec_owner[i];
104: sgarray[i] = garray[i] + (sowners[j]/bs - owners[j]);
105: }
106: /* b index in the IS sfrom */
107: k = sowners[rank]/bs + mbs;
108: for (i=ec,j=0; i< 2*ec; i++,j++) sgarray[i] = k + j;
109: ISCreateBlock(PETSC_COMM_SELF,bs,2*ec,sgarray,PETSC_COPY_VALUES,&from);
111: /* x index in the IS sto */
112: k = sowners[rank]/bs + mbs;
113: for (i=0; i<ec; i++) stmp[i] = (k + i);
114: /* b index in the IS sto */
115: for (i=ec; i<2*ec; i++) stmp[i] = sgarray[i-ec];
117: ISCreateBlock(PETSC_COMM_SELF,bs,2*ec,stmp,PETSC_COPY_VALUES,&to);
119: VecScatterCreate(sbaij->slvec0,from,sbaij->slvec1,to,&sbaij->sMvctx);
121: VecGetLocalSize(sbaij->slvec1,&nt);
122: VecGetArray(sbaij->slvec1,&ptr);
123: VecCreateSeqWithArray(PETSC_COMM_SELF,1,bs*mbs,ptr,&sbaij->slvec1a);
124: VecCreateSeqWithArray(PETSC_COMM_SELF,1,nt-bs*mbs,ptr+bs*mbs,&sbaij->slvec1b);
125: VecRestoreArray(sbaij->slvec1,&ptr);
127: VecGetArray(sbaij->slvec0,&ptr);
128: VecCreateSeqWithArray(PETSC_COMM_SELF,1,nt-bs*mbs,ptr+bs*mbs,&sbaij->slvec0b);
129: VecRestoreArray(sbaij->slvec0,&ptr);
131: PetscFree(stmp);
132: MPI_Barrier(PetscObjectComm((PetscObject)mat));
134: PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->sMvctx);
135: PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->slvec0);
136: PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->slvec1);
137: PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->slvec0b);
138: PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->slvec1a);
139: PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->slvec1b);
140: PetscLogObjectParent((PetscObject)mat,(PetscObject)from);
141: PetscLogObjectParent((PetscObject)mat,(PetscObject)to);
143: PetscLogObjectMemory((PetscObject)mat,(ec+1)*sizeof(PetscInt));
144: ISDestroy(&from);
145: ISDestroy(&to);
146: PetscFree2(sgarray,ec_owner);
147: return(0);
148: }
152: PetscErrorCode MatSetUpMultiply_MPISBAIJ_2comm(Mat mat)
153: {
154: Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)mat->data;
155: Mat_SeqBAIJ *B = (Mat_SeqBAIJ*)(baij->B->data);
157: PetscInt i,j,*aj = B->j,ec = 0,*garray;
158: PetscInt bs = mat->rmap->bs,*stmp;
159: IS from,to;
160: Vec gvec;
161: #if defined(PETSC_USE_CTABLE)
162: PetscTable gid1_lid1;
163: PetscTablePosition tpos;
164: PetscInt gid,lid;
165: #else
166: PetscInt Nbs = baij->Nbs,*indices;
167: #endif
170: #if defined(PETSC_USE_CTABLE)
171: /* use a table - Mark Adams */
172: PetscTableCreate(B->mbs,baij->Nbs+1,&gid1_lid1);
173: for (i=0; i<B->mbs; i++) {
174: for (j=0; j<B->ilen[i]; j++) {
175: PetscInt data,gid1 = aj[B->i[i]+j] + 1;
176: PetscTableFind(gid1_lid1,gid1,&data);
177: if (!data) {
178: /* one based table */
179: PetscTableAdd(gid1_lid1,gid1,++ec,INSERT_VALUES);
180: }
181: }
182: }
183: /* form array of columns we need */
184: PetscMalloc1(ec,&garray);
185: PetscTableGetHeadPosition(gid1_lid1,&tpos);
186: while (tpos) {
187: PetscTableGetNext(gid1_lid1,&tpos,&gid,&lid);
188: gid--; lid--;
189: garray[lid] = gid;
190: }
191: PetscSortInt(ec,garray);
192: PetscTableRemoveAll(gid1_lid1);
193: for (i=0; i<ec; i++) {
194: PetscTableAdd(gid1_lid1,garray[i]+1,i+1,INSERT_VALUES);
195: }
196: /* compact out the extra columns in B */
197: for (i=0; i<B->mbs; i++) {
198: for (j=0; j<B->ilen[i]; j++) {
199: PetscInt gid1 = aj[B->i[i] + j] + 1;
200: PetscTableFind(gid1_lid1,gid1,&lid);
201: lid--;
202: aj[B->i[i]+j] = lid;
203: }
204: }
205: B->nbs = ec;
207: baij->B->cmap->n = baij->B->cmap->N = ec*mat->rmap->bs;
209: PetscLayoutSetUp((baij->B->cmap));
210: PetscTableDestroy(&gid1_lid1);
211: #else
212: /* For the first stab we make an array as long as the number of columns */
213: /* mark those columns that are in baij->B */
214: PetscCalloc1(Nbs,&indices);
215: for (i=0; i<B->mbs; i++) {
216: for (j=0; j<B->ilen[i]; j++) {
217: if (!indices[aj[B->i[i] + j]]) ec++;
218: indices[aj[B->i[i] + j]] = 1;
219: }
220: }
222: /* form array of columns we need */
223: PetscMalloc1(ec,&garray);
224: ec = 0;
225: for (i=0; i<Nbs; i++) {
226: if (indices[i]) {
227: garray[ec++] = i;
228: }
229: }
231: /* make indices now point into garray */
232: for (i=0; i<ec; i++) indices[garray[i]] = i;
234: /* compact out the extra columns in B */
235: for (i=0; i<B->mbs; i++) {
236: for (j=0; j<B->ilen[i]; j++) {
237: aj[B->i[i] + j] = indices[aj[B->i[i] + j]];
238: }
239: }
240: B->nbs = ec;
242: baij->B->cmap->n = ec*mat->rmap->bs;
244: PetscFree(indices);
245: #endif
247: /* create local vector that is used to scatter into */
248: VecCreateSeq(PETSC_COMM_SELF,ec*bs,&baij->lvec);
250: /* create two temporary index sets for building scatter-gather */
251: ISCreateBlock(PETSC_COMM_SELF,bs,ec,garray,PETSC_COPY_VALUES,&from);
253: PetscMalloc1(ec,&stmp);
254: for (i=0; i<ec; i++) stmp[i] = i;
255: ISCreateBlock(PETSC_COMM_SELF,bs,ec,stmp,PETSC_OWN_POINTER,&to);
257: /* create temporary global vector to generate scatter context */
258: VecCreateMPIWithArray(PetscObjectComm((PetscObject)mat),1,mat->cmap->n,mat->cmap->N,NULL,&gvec);
259: VecScatterCreate(gvec,from,baij->lvec,to,&baij->Mvctx);
260: VecDestroy(&gvec);
262: PetscLogObjectParent((PetscObject)mat,(PetscObject)baij->Mvctx);
263: PetscLogObjectParent((PetscObject)mat,(PetscObject)baij->lvec);
264: PetscLogObjectParent((PetscObject)mat,(PetscObject)from);
265: PetscLogObjectParent((PetscObject)mat,(PetscObject)to);
267: baij->garray = garray;
269: PetscLogObjectMemory((PetscObject)mat,(ec+1)*sizeof(PetscInt));
270: ISDestroy(&from);
271: ISDestroy(&to);
272: return(0);
273: }
275: /*
276: Takes the local part of an already assembled MPISBAIJ matrix
277: and disassembles it. This is to allow new nonzeros into the matrix
278: that require more communication in the matrix vector multiply.
279: Thus certain data-structures must be rebuilt.
281: Kind of slow! But that's what application programmers get when
282: they are sloppy.
283: */
286: PetscErrorCode MatDisAssemble_MPISBAIJ(Mat A)
287: {
288: Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)A->data;
289: Mat B = baij->B,Bnew;
290: Mat_SeqBAIJ *Bbaij = (Mat_SeqBAIJ*)B->data;
292: PetscInt i,j,mbs=Bbaij->mbs,n = A->cmap->N,col,*garray=baij->garray;
293: PetscInt k,bs=A->rmap->bs,bs2=baij->bs2,*rvals,*nz,ec,m=A->rmap->n;
294: MatScalar *a = Bbaij->a;
295: PetscScalar *atmp;
296: #if defined(PETSC_USE_REAL_MAT_SINGLE)
297: PetscInt l;
298: #endif
301: #if defined(PETSC_USE_REAL_MAT_SINGLE)
302: PetscMalloc1(A->rmap->bs,&atmp);
303: #endif
304: /* free stuff related to matrix-vec multiply */
305: VecGetSize(baij->lvec,&ec); /* needed for PetscLogObjectMemory below */
306: VecDestroy(&baij->lvec);
307: VecScatterDestroy(&baij->Mvctx);
309: VecDestroy(&baij->slvec0);
310: VecDestroy(&baij->slvec0b);
311: VecDestroy(&baij->slvec1);
312: VecDestroy(&baij->slvec1a);
313: VecDestroy(&baij->slvec1b);
315: if (baij->colmap) {
316: #if defined(PETSC_USE_CTABLE)
317: PetscTableDestroy(&baij->colmap);
318: #else
319: PetscFree(baij->colmap);
320: PetscLogObjectMemory((PetscObject)A,-Bbaij->nbs*sizeof(PetscInt));
321: #endif
322: }
324: /* make sure that B is assembled so we can access its values */
325: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
326: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
328: /* invent new B and copy stuff over */
329: PetscMalloc1(mbs,&nz);
330: for (i=0; i<mbs; i++) {
331: nz[i] = Bbaij->i[i+1]-Bbaij->i[i];
332: }
333: MatCreate(PETSC_COMM_SELF,&Bnew);
334: MatSetSizes(Bnew,m,n,m,n);
335: MatSetType(Bnew,((PetscObject)B)->type_name);
336: MatSeqBAIJSetPreallocation(Bnew,B->rmap->bs,0,nz);
337: PetscFree(nz);
339: ((Mat_SeqSBAIJ*)Bnew->data)->nonew = Bbaij->nonew; /* Inherit insertion error options. */
340: /*
341: Ensure that B's nonzerostate is monotonically increasing.
342: Or should this follow the MatSetValues() loop to preserve B's nonzerstate across a MatDisAssemble() call?
343: */
344: Bnew->nonzerostate = B->nonzerostate;
345: PetscMalloc1(bs,&rvals);
346: for (i=0; i<mbs; i++) {
347: rvals[0] = bs*i;
348: for (j=1; j<bs; j++) rvals[j] = rvals[j-1] + 1;
349: for (j=Bbaij->i[i]; j<Bbaij->i[i+1]; j++) {
350: col = garray[Bbaij->j[j]]*bs;
351: for (k=0; k<bs; k++) {
352: #if defined(PETSC_USE_REAL_MAT_SINGLE)
353: for (l=0; l<bs; l++) atmp[l] = a[j*bs2+l];
354: #else
355: atmp = a+j*bs2 + k*bs;
356: #endif
357: MatSetValues_SeqSBAIJ(Bnew,bs,rvals,1,&col,atmp,B->insertmode);
358: col++;
359: }
360: }
361: }
362: #if defined(PETSC_USE_REAL_MAT_SINGLE)
363: PetscFree(atmp);
364: #endif
365: PetscFree(baij->garray);
367: baij->garray = 0;
369: PetscFree(rvals);
370: PetscLogObjectMemory((PetscObject)A,-ec*sizeof(PetscInt));
371: MatDestroy(&B);
372: PetscLogObjectParent((PetscObject)A,(PetscObject)Bnew);
374: baij->B = Bnew;
376: A->was_assembled = PETSC_FALSE;
377: return(0);
378: }