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>
7: PetscErrorCode MatSetUpMultiply_MPISBAIJ(Mat mat)
8: {
9: Mat_MPISBAIJ *sbaij = (Mat_MPISBAIJ*)mat->data;
10: Mat_SeqBAIJ *B = (Mat_SeqBAIJ*)(sbaij->B->data);
11: PetscInt Nbs = sbaij->Nbs,i,j,*aj = B->j,ec = 0,*garray,*sgarray;
12: PetscInt bs = mat->rmap->bs,*stmp,mbs=sbaij->mbs, vec_size,nt;
13: IS from,to;
14: Vec gvec;
15: PetscMPIInt rank = sbaij->rank,lsize;
16: PetscInt *owners = sbaij->rangebs,*ec_owner,k;
17: const PetscInt *sowners;
18: PetscScalar *ptr;
19: #if defined(PETSC_USE_CTABLE)
20: PetscTable gid1_lid1; /* one-based gid to lid table */
21: PetscTablePosition tpos;
22: PetscInt gid,lid;
23: #else
24: PetscInt *indices;
25: #endif
27: #if defined(PETSC_USE_CTABLE)
28: PetscTableCreate(mbs,Nbs+1,&gid1_lid1);
29: for (i=0; i<B->mbs; i++) {
30: for (j=0; j<B->ilen[i]; j++) {
31: PetscInt data,gid1 = aj[B->i[i]+j] + 1;
32: PetscTableFind(gid1_lid1,gid1,&data);
33: if (!data) PetscTableAdd(gid1_lid1,gid1,++ec,INSERT_VALUES);
34: }
35: }
36: /* form array of columns we need */
37: PetscMalloc1(ec,&garray);
38: PetscTableGetHeadPosition(gid1_lid1,&tpos);
39: while (tpos) {
40: PetscTableGetNext(gid1_lid1,&tpos,&gid,&lid);
41: gid--; lid--;
42: garray[lid] = gid;
43: }
44: PetscSortInt(ec,garray);
45: PetscTableRemoveAll(gid1_lid1);
46: for (i=0; i<ec; i++) PetscTableAdd(gid1_lid1,garray[i]+1,i+1,INSERT_VALUES);
47: /* compact out the extra columns in B */
48: for (i=0; i<B->mbs; i++) {
49: for (j=0; j<B->ilen[i]; j++) {
50: PetscInt gid1 = aj[B->i[i] + j] + 1;
51: PetscTableFind(gid1_lid1,gid1,&lid);
52: lid--;
53: aj[B->i[i]+j] = lid;
54: }
55: }
56: PetscTableDestroy(&gid1_lid1);
57: PetscMalloc2(2*ec,&sgarray,ec,&ec_owner);
58: for (i=j=0; i<ec; i++) {
59: while (garray[i]>=owners[j+1]) j++;
60: ec_owner[i] = j;
61: }
62: #else
63: /* For the first stab we make an array as long as the number of columns */
64: /* mark those columns that are in sbaij->B */
65: PetscCalloc1(Nbs,&indices);
66: for (i=0; i<mbs; i++) {
67: for (j=0; j<B->ilen[i]; j++) {
68: if (!indices[aj[B->i[i] + j]]) ec++;
69: indices[aj[B->i[i] + j]] = 1;
70: }
71: }
73: /* form arrays of columns we need */
74: PetscMalloc1(ec,&garray);
75: PetscMalloc2(2*ec,&sgarray,ec,&ec_owner);
77: ec = 0;
78: for (j=0; j<sbaij->size; j++) {
79: for (i=owners[j]; i<owners[j+1]; i++) {
80: if (indices[i]) {
81: garray[ec] = i;
82: ec_owner[ec] = j;
83: ec++;
84: }
85: }
86: }
88: /* make indices now point into garray */
89: for (i=0; i<ec; i++) indices[garray[i]] = i;
91: /* compact out the extra columns in B */
92: for (i=0; i<mbs; i++) {
93: for (j=0; j<B->ilen[i]; j++) aj[B->i[i] + j] = indices[aj[B->i[i] + j]];
94: }
95: PetscFree(indices);
96: #endif
97: B->nbs = ec;
98: PetscLayoutDestroy(&sbaij->B->cmap);
99: PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)sbaij->B),ec*mat->rmap->bs,ec*mat->rmap->bs,mat->rmap->bs,&sbaij->B->cmap);
101: VecScatterDestroy(&sbaij->sMvctx);
102: /* create local vector that is used to scatter into */
103: VecCreateSeq(PETSC_COMM_SELF,ec*bs,&sbaij->lvec);
105: /* create two temporary index sets for building scatter-gather */
106: PetscMalloc1(2*ec,&stmp);
107: ISCreateBlock(PETSC_COMM_SELF,bs,ec,garray,PETSC_COPY_VALUES,&from);
108: for (i=0; i<ec; i++) stmp[i] = i;
109: ISCreateBlock(PETSC_COMM_SELF,bs,ec,stmp,PETSC_COPY_VALUES,&to);
111: /* generate the scatter context
112: -- Mvctx and lvec are not used by MatMult_MPISBAIJ(), but have other uses, such as in MatDiagonalScale_MPISBAIJ */
113: VecCreateMPIWithArray(PetscObjectComm((PetscObject)mat),1,mat->cmap->n,mat->cmap->N,NULL,&gvec);
114: VecScatterCreate(gvec,from,sbaij->lvec,to,&sbaij->Mvctx);
115: VecDestroy(&gvec);
117: sbaij->garray = garray;
119: PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->Mvctx);
120: PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->lvec);
121: PetscLogObjectParent((PetscObject)mat,(PetscObject)from);
122: PetscLogObjectParent((PetscObject)mat,(PetscObject)to);
124: ISDestroy(&from);
125: ISDestroy(&to);
127: /* create parallel vector that is used by SBAIJ matrix to scatter from/into */
128: lsize = (mbs + ec)*bs;
129: VecCreateMPI(PetscObjectComm((PetscObject)mat),lsize,PETSC_DETERMINE,&sbaij->slvec0);
130: VecDuplicate(sbaij->slvec0,&sbaij->slvec1);
131: VecGetSize(sbaij->slvec0,&vec_size);
133: VecGetOwnershipRanges(sbaij->slvec0,&sowners);
135: /* x index in the IS sfrom */
136: for (i=0; i<ec; i++) {
137: j = ec_owner[i];
138: sgarray[i] = garray[i] + (sowners[j]/bs - owners[j]);
139: }
140: /* b index in the IS sfrom */
141: k = sowners[rank]/bs + mbs;
142: for (i=ec,j=0; i< 2*ec; i++,j++) sgarray[i] = k + j;
143: ISCreateBlock(PETSC_COMM_SELF,bs,2*ec,sgarray,PETSC_COPY_VALUES,&from);
145: /* x index in the IS sto */
146: k = sowners[rank]/bs + mbs;
147: for (i=0; i<ec; i++) stmp[i] = (k + i);
148: /* b index in the IS sto */
149: for (i=ec; i<2*ec; i++) stmp[i] = sgarray[i-ec];
151: ISCreateBlock(PETSC_COMM_SELF,bs,2*ec,stmp,PETSC_COPY_VALUES,&to);
153: VecScatterCreate(sbaij->slvec0,from,sbaij->slvec1,to,&sbaij->sMvctx);
154: VecScatterViewFromOptions(sbaij->sMvctx,(PetscObject)mat,"-matmult_vecscatter_view");
156: VecGetLocalSize(sbaij->slvec1,&nt);
157: VecGetArray(sbaij->slvec1,&ptr);
158: VecCreateSeqWithArray(PETSC_COMM_SELF,1,bs*mbs,ptr,&sbaij->slvec1a);
159: VecCreateSeqWithArray(PETSC_COMM_SELF,1,nt-bs*mbs,ptr+bs*mbs,&sbaij->slvec1b);
160: VecRestoreArray(sbaij->slvec1,&ptr);
162: VecGetArray(sbaij->slvec0,&ptr);
163: VecCreateSeqWithArray(PETSC_COMM_SELF,1,nt-bs*mbs,ptr+bs*mbs,&sbaij->slvec0b);
164: VecRestoreArray(sbaij->slvec0,&ptr);
166: PetscFree(stmp);
168: PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->sMvctx);
169: PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->slvec0);
170: PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->slvec1);
171: PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->slvec0b);
172: PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->slvec1a);
173: PetscLogObjectParent((PetscObject)mat,(PetscObject)sbaij->slvec1b);
174: PetscLogObjectParent((PetscObject)mat,(PetscObject)from);
175: PetscLogObjectParent((PetscObject)mat,(PetscObject)to);
177: PetscLogObjectMemory((PetscObject)mat,ec*sizeof(PetscInt));
178: ISDestroy(&from);
179: ISDestroy(&to);
180: PetscFree2(sgarray,ec_owner);
181: return 0;
182: }
184: /*
185: Takes the local part of an already assembled MPISBAIJ matrix
186: and disassembles it. This is to allow new nonzeros into the matrix
187: that require more communication in the matrix vector multiply.
188: Thus certain data-structures must be rebuilt.
190: Kind of slow! But that's what application programmers get when
191: they are sloppy.
192: */
193: PetscErrorCode MatDisAssemble_MPISBAIJ(Mat A)
194: {
195: Mat_MPISBAIJ *baij = (Mat_MPISBAIJ*)A->data;
196: Mat B = baij->B,Bnew;
197: Mat_SeqBAIJ *Bbaij = (Mat_SeqBAIJ*)B->data;
198: PetscInt i,j,mbs=Bbaij->mbs,n = A->cmap->N,col,*garray=baij->garray;
199: PetscInt k,bs=A->rmap->bs,bs2=baij->bs2,*rvals,*nz,ec,m=A->rmap->n;
200: MatScalar *a = Bbaij->a;
201: PetscScalar *atmp;
202: #if defined(PETSC_USE_REAL_MAT_SINGLE)
203: PetscInt l;
204: #endif
206: #if defined(PETSC_USE_REAL_MAT_SINGLE)
207: PetscMalloc1(A->rmap->bs,&atmp);
208: #endif
209: /* free stuff related to matrix-vec multiply */
210: VecGetSize(baij->lvec,&ec); /* needed for PetscLogObjectMemory below */
211: VecDestroy(&baij->lvec);
212: VecScatterDestroy(&baij->Mvctx);
214: VecDestroy(&baij->slvec0);
215: VecDestroy(&baij->slvec0b);
216: VecDestroy(&baij->slvec1);
217: VecDestroy(&baij->slvec1a);
218: VecDestroy(&baij->slvec1b);
220: if (baij->colmap) {
221: #if defined(PETSC_USE_CTABLE)
222: PetscTableDestroy(&baij->colmap);
223: #else
224: PetscFree(baij->colmap);
225: PetscLogObjectMemory((PetscObject)A,-Bbaij->nbs*sizeof(PetscInt));
226: #endif
227: }
229: /* make sure that B is assembled so we can access its values */
230: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
231: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
233: /* invent new B and copy stuff over */
234: PetscMalloc1(mbs,&nz);
235: for (i=0; i<mbs; i++) {
236: nz[i] = Bbaij->i[i+1]-Bbaij->i[i];
237: }
238: MatCreate(PETSC_COMM_SELF,&Bnew);
239: MatSetSizes(Bnew,m,n,m,n);
240: MatSetType(Bnew,((PetscObject)B)->type_name);
241: MatSeqBAIJSetPreallocation(Bnew,B->rmap->bs,0,nz);
242: PetscFree(nz);
244: if (Bbaij->nonew >= 0) { /* Inherit insertion error options (if positive). */
245: ((Mat_SeqSBAIJ*)Bnew->data)->nonew = Bbaij->nonew;
246: }
248: /*
249: Ensure that B's nonzerostate is monotonically increasing.
250: Or should this follow the MatSetValues() loop to preserve B's nonzerstate across a MatDisAssemble() call?
251: */
252: Bnew->nonzerostate = B->nonzerostate;
253: PetscMalloc1(bs,&rvals);
254: for (i=0; i<mbs; i++) {
255: rvals[0] = bs*i;
256: for (j=1; j<bs; j++) rvals[j] = rvals[j-1] + 1;
257: for (j=Bbaij->i[i]; j<Bbaij->i[i+1]; j++) {
258: col = garray[Bbaij->j[j]]*bs;
259: for (k=0; k<bs; k++) {
260: #if defined(PETSC_USE_REAL_MAT_SINGLE)
261: for (l=0; l<bs; l++) atmp[l] = a[j*bs2+l];
262: #else
263: atmp = a+j*bs2 + k*bs;
264: #endif
265: MatSetValues_SeqSBAIJ(Bnew,bs,rvals,1,&col,atmp,B->insertmode);
266: col++;
267: }
268: }
269: }
270: #if defined(PETSC_USE_REAL_MAT_SINGLE)
271: PetscFree(atmp);
272: #endif
273: PetscFree(baij->garray);
275: baij->garray = NULL;
277: PetscFree(rvals);
278: PetscLogObjectMemory((PetscObject)A,-ec*sizeof(PetscInt));
279: MatDestroy(&B);
280: PetscLogObjectParent((PetscObject)A,(PetscObject)Bnew);
282: baij->B = Bnew;
284: A->was_assembled = PETSC_FALSE;
285: return 0;
286: }