Actual source code: mmsbaij.c
petsc-3.12.5 2020-03-29
2: /*
3: Support for the parallel SBAIJ matrix vector multiply
4: */
5: #include <../src/mat/impls/sbaij/mpi/mpisbaij.h>
7: PetscErrorCode MatSetUpMultiply_MPISBAIJ(Mat mat)
8: {
9: Mat_MPISBAIJ *sbaij = (Mat_MPISBAIJ*)mat->data;
10: Mat_SeqBAIJ *B = (Mat_SeqBAIJ*)(sbaij->B->data);
12: PetscInt Nbs = sbaij->Nbs,i,j,*indices,*aj = B->j,ec = 0,*garray,*sgarray;
13: PetscInt bs = mat->rmap->bs,*stmp,mbs=sbaij->mbs, vec_size,nt;
14: IS from,to;
15: Vec gvec;
16: PetscMPIInt rank =sbaij->rank,lsize,size=sbaij->size;
17: PetscInt *owners=sbaij->rangebs,*ec_owner,k;
18: const PetscInt *sowners;
19: PetscScalar *ptr;
22: VecScatterDestroy(&sbaij->sMvctx);
24: /* For the first stab we make an array as long as the number of columns */
25: /* mark those columns that are in sbaij->B */
26: PetscCalloc1(Nbs,&indices);
27: for (i=0; i<mbs; i++) {
28: for (j=0; j<B->ilen[i]; j++) {
29: if (!indices[aj[B->i[i] + j]]) ec++;
30: indices[aj[B->i[i] + j]] = 1;
31: }
32: }
34: /* form arrays of columns we need */
35: PetscMalloc1(ec,&garray);
36: PetscMalloc2(2*ec,&sgarray,ec,&ec_owner);
38: ec = 0;
39: for (j=0; j<size; j++) {
40: for (i=owners[j]; i<owners[j+1]; i++) {
41: if (indices[i]) {
42: garray[ec] = i;
43: ec_owner[ec] = j;
44: ec++;
45: }
46: }
47: }
49: /* make indices now point into garray */
50: for (i=0; i<ec; i++) indices[garray[i]] = i;
52: /* compact out the extra columns in B */
53: for (i=0; i<mbs; i++) {
54: for (j=0; j<B->ilen[i]; j++) aj[B->i[i] + j] = indices[aj[B->i[i] + j]];
55: }
56: B->nbs = ec;
57: PetscLayoutDestroy(&sbaij->B->cmap);
58: PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)sbaij->B),ec*mat->rmap->bs,ec*mat->rmap->bs,mat->rmap->bs,&sbaij->B->cmap);
59: PetscFree(indices);
61: /* create local vector that is used to scatter into */
62: VecCreateSeq(PETSC_COMM_SELF,ec*bs,&sbaij->lvec);
64: /* create two temporary index sets for building scatter-gather */
65: PetscMalloc1(2*ec,&stmp);
66: ISCreateBlock(PETSC_COMM_SELF,bs,ec,garray,PETSC_COPY_VALUES,&from);
67: for (i=0; i<ec; i++) stmp[i] = i;
68: ISCreateBlock(PETSC_COMM_SELF,bs,ec,stmp,PETSC_COPY_VALUES,&to);
70: /* generate the scatter context
71: -- Mvctx and lvec are not used by MatMult_MPISBAIJ(), but usefull for some applications */
72: VecCreateMPIWithArray(PetscObjectComm((PetscObject)mat),1,mat->cmap->n,mat->cmap->N,NULL,&gvec);
73: VecScatterCreate(gvec,from,sbaij->lvec,to,&sbaij->Mvctx);
74: VecDestroy(&gvec);
76: sbaij->garray = garray;
78: PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->Mvctx);
79: PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->lvec);
80: PetscLogObjectParent((PetscObject)mat,(PetscObject)from);
81: PetscLogObjectParent((PetscObject)mat,(PetscObject)to);
83: ISDestroy(&from);
84: ISDestroy(&to);
86: /* create parallel vector that is used by SBAIJ matrix to scatter from/into */
87: lsize = (mbs + ec)*bs;
88: VecCreateMPI(PetscObjectComm((PetscObject)mat),lsize,PETSC_DETERMINE,&sbaij->slvec0);
89: VecDuplicate(sbaij->slvec0,&sbaij->slvec1);
90: VecGetSize(sbaij->slvec0,&vec_size);
92: VecGetOwnershipRanges(sbaij->slvec0,&sowners);
94: /* x index in the IS sfrom */
95: for (i=0; i<ec; i++) {
96: j = ec_owner[i];
97: sgarray[i] = garray[i] + (sowners[j]/bs - owners[j]);
98: }
99: /* b index in the IS sfrom */
100: k = sowners[rank]/bs + mbs;
101: for (i=ec,j=0; i< 2*ec; i++,j++) sgarray[i] = k + j;
102: ISCreateBlock(PETSC_COMM_SELF,bs,2*ec,sgarray,PETSC_COPY_VALUES,&from);
104: /* x index in the IS sto */
105: k = sowners[rank]/bs + mbs;
106: for (i=0; i<ec; i++) stmp[i] = (k + i);
107: /* b index in the IS sto */
108: for (i=ec; i<2*ec; i++) stmp[i] = sgarray[i-ec];
110: ISCreateBlock(PETSC_COMM_SELF,bs,2*ec,stmp,PETSC_COPY_VALUES,&to);
112: VecScatterCreate(sbaij->slvec0,from,sbaij->slvec1,to,&sbaij->sMvctx);
114: VecGetLocalSize(sbaij->slvec1,&nt);
115: VecGetArray(sbaij->slvec1,&ptr);
116: VecCreateSeqWithArray(PETSC_COMM_SELF,1,bs*mbs,ptr,&sbaij->slvec1a);
117: VecCreateSeqWithArray(PETSC_COMM_SELF,1,nt-bs*mbs,ptr+bs*mbs,&sbaij->slvec1b);
118: VecRestoreArray(sbaij->slvec1,&ptr);
120: VecGetArray(sbaij->slvec0,&ptr);
121: VecCreateSeqWithArray(PETSC_COMM_SELF,1,nt-bs*mbs,ptr+bs*mbs,&sbaij->slvec0b);
122: VecRestoreArray(sbaij->slvec0,&ptr);
124: PetscFree(stmp);
125: MPI_Barrier(PetscObjectComm((PetscObject)mat));
127: PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->sMvctx);
128: PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->slvec0);
129: PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->slvec1);
130: PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->slvec0b);
131: PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->slvec1a);
132: PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->slvec1b);
133: PetscLogObjectParent((PetscObject)mat,(PetscObject)from);
134: PetscLogObjectParent((PetscObject)mat,(PetscObject)to);
136: PetscLogObjectMemory((PetscObject)mat,(ec+1)*sizeof(PetscInt));
137: ISDestroy(&from);
138: ISDestroy(&to);
139: PetscFree2(sgarray,ec_owner);
140: return(0);
141: }
143: PetscErrorCode MatSetUpMultiply_MPISBAIJ_2comm(Mat mat)
144: {
145: Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)mat->data;
146: Mat_SeqBAIJ *B = (Mat_SeqBAIJ*)(baij->B->data);
148: PetscInt i,j,*aj = B->j,ec = 0,*garray;
149: PetscInt bs = mat->rmap->bs,*stmp;
150: IS from,to;
151: Vec gvec;
152: #if defined(PETSC_USE_CTABLE)
153: PetscTable gid1_lid1;
154: PetscTablePosition tpos;
155: PetscInt gid,lid;
156: #else
157: PetscInt Nbs = baij->Nbs,*indices;
158: #endif
161: #if defined(PETSC_USE_CTABLE)
162: /* use a table - Mark Adams */
163: PetscTableCreate(B->mbs,baij->Nbs+1,&gid1_lid1);
164: for (i=0; i<B->mbs; i++) {
165: for (j=0; j<B->ilen[i]; j++) {
166: PetscInt data,gid1 = aj[B->i[i]+j] + 1;
167: PetscTableFind(gid1_lid1,gid1,&data);
168: if (!data) {
169: /* one based table */
170: PetscTableAdd(gid1_lid1,gid1,++ec,INSERT_VALUES);
171: }
172: }
173: }
174: /* form array of columns we need */
175: PetscMalloc1(ec,&garray);
176: PetscTableGetHeadPosition(gid1_lid1,&tpos);
177: while (tpos) {
178: PetscTableGetNext(gid1_lid1,&tpos,&gid,&lid);
179: gid--; lid--;
180: garray[lid] = gid;
181: }
182: PetscSortInt(ec,garray);
183: PetscTableRemoveAll(gid1_lid1);
184: for (i=0; i<ec; i++) {
185: PetscTableAdd(gid1_lid1,garray[i]+1,i+1,INSERT_VALUES);
186: }
187: /* compact out the extra columns in B */
188: for (i=0; i<B->mbs; i++) {
189: for (j=0; j<B->ilen[i]; j++) {
190: PetscInt gid1 = aj[B->i[i] + j] + 1;
191: PetscTableFind(gid1_lid1,gid1,&lid);
192: lid--;
193: aj[B->i[i]+j] = lid;
194: }
195: }
196: B->nbs = ec;
197: PetscLayoutDestroy(&baij->B->cmap);
198: PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)baij->B),ec*mat->rmap->bs,ec*mat->rmap->bs,mat->rmap->bs,&baij->B->cmap);
199: PetscTableDestroy(&gid1_lid1);
200: #else
201: /* For the first stab we make an array as long as the number of columns */
202: /* mark those columns that are in baij->B */
203: PetscCalloc1(Nbs,&indices);
204: for (i=0; i<B->mbs; i++) {
205: for (j=0; j<B->ilen[i]; j++) {
206: if (!indices[aj[B->i[i] + j]]) ec++;
207: indices[aj[B->i[i] + j]] = 1;
208: }
209: }
211: /* form array of columns we need */
212: PetscMalloc1(ec,&garray);
213: ec = 0;
214: for (i=0; i<Nbs; i++) {
215: if (indices[i]) {
216: garray[ec++] = i;
217: }
218: }
220: /* make indices now point into garray */
221: for (i=0; i<ec; i++) indices[garray[i]] = i;
223: /* compact out the extra columns in B */
224: for (i=0; i<B->mbs; i++) {
225: for (j=0; j<B->ilen[i]; j++) {
226: aj[B->i[i] + j] = indices[aj[B->i[i] + j]];
227: }
228: }
229: B->nbs = ec;
231: baij->B->cmap->n = ec*mat->rmap->bs;
233: PetscFree(indices);
234: #endif
236: /* create local vector that is used to scatter into */
237: VecCreateSeq(PETSC_COMM_SELF,ec*bs,&baij->lvec);
239: /* create two temporary index sets for building scatter-gather */
240: ISCreateBlock(PETSC_COMM_SELF,bs,ec,garray,PETSC_COPY_VALUES,&from);
242: PetscMalloc1(ec,&stmp);
243: for (i=0; i<ec; i++) stmp[i] = i;
244: ISCreateBlock(PETSC_COMM_SELF,bs,ec,stmp,PETSC_OWN_POINTER,&to);
246: /* create temporary global vector to generate scatter context */
247: VecCreateMPIWithArray(PetscObjectComm((PetscObject)mat),1,mat->cmap->n,mat->cmap->N,NULL,&gvec);
248: VecScatterCreate(gvec,from,baij->lvec,to,&baij->Mvctx);
249: VecDestroy(&gvec);
251: PetscLogObjectParent((PetscObject)mat,(PetscObject)baij->Mvctx);
252: PetscLogObjectParent((PetscObject)mat,(PetscObject)baij->lvec);
253: PetscLogObjectParent((PetscObject)mat,(PetscObject)from);
254: PetscLogObjectParent((PetscObject)mat,(PetscObject)to);
256: baij->garray = garray;
258: PetscLogObjectMemory((PetscObject)mat,(ec+1)*sizeof(PetscInt));
259: ISDestroy(&from);
260: ISDestroy(&to);
261: return(0);
262: }
264: /*
265: Takes the local part of an already assembled MPISBAIJ matrix
266: and disassembles it. This is to allow new nonzeros into the matrix
267: that require more communication in the matrix vector multiply.
268: Thus certain data-structures must be rebuilt.
270: Kind of slow! But that's what application programmers get when
271: they are sloppy.
272: */
273: PetscErrorCode MatDisAssemble_MPISBAIJ(Mat A)
274: {
275: Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)A->data;
276: Mat B = baij->B,Bnew;
277: Mat_SeqBAIJ *Bbaij = (Mat_SeqBAIJ*)B->data;
279: PetscInt i,j,mbs=Bbaij->mbs,n = A->cmap->N,col,*garray=baij->garray;
280: PetscInt k,bs=A->rmap->bs,bs2=baij->bs2,*rvals,*nz,ec,m=A->rmap->n;
281: MatScalar *a = Bbaij->a;
282: PetscScalar *atmp;
283: #if defined(PETSC_USE_REAL_MAT_SINGLE)
284: PetscInt l;
285: #endif
288: #if defined(PETSC_USE_REAL_MAT_SINGLE)
289: PetscMalloc1(A->rmap->bs,&atmp);
290: #endif
291: /* free stuff related to matrix-vec multiply */
292: VecGetSize(baij->lvec,&ec); /* needed for PetscLogObjectMemory below */
293: VecDestroy(&baij->lvec);
294: VecScatterDestroy(&baij->Mvctx);
296: VecDestroy(&baij->slvec0);
297: VecDestroy(&baij->slvec0b);
298: VecDestroy(&baij->slvec1);
299: VecDestroy(&baij->slvec1a);
300: VecDestroy(&baij->slvec1b);
302: if (baij->colmap) {
303: #if defined(PETSC_USE_CTABLE)
304: PetscTableDestroy(&baij->colmap);
305: #else
306: PetscFree(baij->colmap);
307: PetscLogObjectMemory((PetscObject)A,-Bbaij->nbs*sizeof(PetscInt));
308: #endif
309: }
311: /* make sure that B is assembled so we can access its values */
312: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
313: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
315: /* invent new B and copy stuff over */
316: PetscMalloc1(mbs,&nz);
317: for (i=0; i<mbs; i++) {
318: nz[i] = Bbaij->i[i+1]-Bbaij->i[i];
319: }
320: MatCreate(PETSC_COMM_SELF,&Bnew);
321: MatSetSizes(Bnew,m,n,m,n);
322: MatSetType(Bnew,((PetscObject)B)->type_name);
323: MatSeqBAIJSetPreallocation(Bnew,B->rmap->bs,0,nz);
324: PetscFree(nz);
326: if (Bbaij->nonew >= 0) { /* Inherit insertion error options (if positive). */
327: ((Mat_SeqSBAIJ*)Bnew->data)->nonew = Bbaij->nonew;
328: }
330: /*
331: Ensure that B's nonzerostate is monotonically increasing.
332: Or should this follow the MatSetValues() loop to preserve B's nonzerstate across a MatDisAssemble() call?
333: */
334: Bnew->nonzerostate = B->nonzerostate;
335: PetscMalloc1(bs,&rvals);
336: for (i=0; i<mbs; i++) {
337: rvals[0] = bs*i;
338: for (j=1; j<bs; j++) rvals[j] = rvals[j-1] + 1;
339: for (j=Bbaij->i[i]; j<Bbaij->i[i+1]; j++) {
340: col = garray[Bbaij->j[j]]*bs;
341: for (k=0; k<bs; k++) {
342: #if defined(PETSC_USE_REAL_MAT_SINGLE)
343: for (l=0; l<bs; l++) atmp[l] = a[j*bs2+l];
344: #else
345: atmp = a+j*bs2 + k*bs;
346: #endif
347: MatSetValues_SeqSBAIJ(Bnew,bs,rvals,1,&col,atmp,B->insertmode);
348: col++;
349: }
350: }
351: }
352: #if defined(PETSC_USE_REAL_MAT_SINGLE)
353: PetscFree(atmp);
354: #endif
355: PetscFree(baij->garray);
357: baij->garray = 0;
359: PetscFree(rvals);
360: PetscLogObjectMemory((PetscObject)A,-ec*sizeof(PetscInt));
361: MatDestroy(&B);
362: PetscLogObjectParent((PetscObject)A,(PetscObject)Bnew);
364: baij->B = Bnew;
366: A->was_assembled = PETSC_FALSE;
367: return(0);
368: }