Actual source code: mmbaij.c
1: /*$Id: mmbaij.c,v 1.46 2001/09/25 00:31:36 balay Exp $*/
3: /*
4: Support for the parallel BAIJ matrix vector multiply
5: */
6: #include src/mat/impls/baij/mpi/mpibaij.h
7: #include src/vec/vecimpl.h
9: EXTERN int MatSetValuesBlocked_SeqBAIJ(Mat,int,int*,int,int*,PetscScalar*,InsertMode);
11: int MatSetUpMultiply_MPIBAIJ(Mat mat)
12: {
13: Mat_MPIBAIJ *baij = (Mat_MPIBAIJ*)mat->data;
14: Mat_SeqBAIJ *B = (Mat_SeqBAIJ*)(baij->B->data);
15: int Nbs = baij->Nbs,i,j,*indices,*aj = B->j,ierr,ec = 0,*garray;
16: int bs = baij->bs,*stmp;
17: IS from,to;
18: Vec gvec;
19: #if defined (PETSC_USE_CTABLE)
20: PetscTable gid1_lid1;
21: PetscTablePosition tpos;
22: int gid,lid;
23: #endif
27: #if defined (PETSC_USE_CTABLE)
28: /* use a table - Mark Adams */
29: PetscTableCreate(B->mbs,&gid1_lid1);
30: for (i=0; i<B->mbs; i++) {
31: for (j=0; j<B->ilen[i]; j++) {
32: int data,gid1 = aj[B->i[i]+j] + 1;
33: PetscTableFind(gid1_lid1,gid1,&data) ;
34: if (!data) {
35: /* one based table */
36: PetscTableAdd(gid1_lid1,gid1,++ec);
37: }
38: }
39: }
40: /* form array of columns we need */
41: PetscMalloc((ec+1)*sizeof(int),&garray);
42: PetscTableGetHeadPosition(gid1_lid1,&tpos);
43: while (tpos) {
44: PetscTableGetNext(gid1_lid1,&tpos,&gid,&lid);
45: gid--; lid--;
46: garray[lid] = gid;
47: }
48: PetscSortInt(ec,garray);
49: /* qsort(garray, ec, sizeof(int), intcomparcarc); */
50: PetscTableRemoveAll(gid1_lid1);
51: for (i=0; i<ec; i++) {
52: PetscTableAdd(gid1_lid1,garray[i]+1,i+1);
53: }
54: /* compact out the extra columns in B */
55: for (i=0; i<B->mbs; i++) {
56: for (j=0; j<B->ilen[i]; j++) {
57: int gid1 = aj[B->i[i] + j] + 1;
58: PetscTableFind(gid1_lid1,gid1,&lid);
59: lid --;
60: aj[B->i[i]+j] = lid;
61: }
62: }
63: B->nbs = ec;
64: baij->B->n = ec*B->bs;
65: PetscTableDelete(gid1_lid1);
66: /* Mark Adams */
67: #else
68: /* Make an array as long as the number of columns */
69: /* mark those columns that are in baij->B */
70: PetscMalloc((Nbs+1)*sizeof(int),&indices);
71: PetscMemzero(indices,Nbs*sizeof(int));
72: for (i=0; i<B->mbs; i++) {
73: for (j=0; j<B->ilen[i]; j++) {
74: if (!indices[aj[B->i[i] + j]]) ec++;
75: indices[aj[B->i[i] + j]] = 1;
76: }
77: }
79: /* form array of columns we need */
80: PetscMalloc((ec+1)*sizeof(int),&garray);
81: ec = 0;
82: for (i=0; i<Nbs; i++) {
83: if (indices[i]) {
84: garray[ec++] = i;
85: }
86: }
88: /* make indices now point into garray */
89: for (i=0; i<ec; i++) {
90: indices[garray[i]] = i;
91: }
93: /* compact out the extra columns in B */
94: for (i=0; i<B->mbs; i++) {
95: for (j=0; j<B->ilen[i]; j++) {
96: aj[B->i[i] + j] = indices[aj[B->i[i] + j]];
97: }
98: }
99: B->nbs = ec;
100: baij->B->n = ec*B->bs;
101: PetscFree(indices);
102: #endif
104: /* create local vector that is used to scatter into */
105: VecCreateSeq(PETSC_COMM_SELF,ec*bs,&baij->lvec);
107: /* create two temporary index sets for building scatter-gather */
108: for (i=0; i<ec; i++) {
109: garray[i] = bs*garray[i];
110: }
111: ISCreateBlock(PETSC_COMM_SELF,bs,ec,garray,&from);
112: for (i=0; i<ec; i++) {
113: garray[i] = garray[i]/bs;
114: }
116: PetscMalloc((ec+1)*sizeof(int),&stmp);
117: for (i=0; i<ec; i++) { stmp[i] = bs*i; }
118: ISCreateBlock(PETSC_COMM_SELF,bs,ec,stmp,&to);
119: PetscFree(stmp);
121: /* create temporary global vector to generate scatter context */
122: /* this is inefficient, but otherwise we must do either
123: 1) save garray until the first actual scatter when the vector is known or
124: 2) have another way of generating a scatter context without a vector.*/
125: VecCreateMPI(mat->comm,mat->n,mat->N,&gvec);
127: /* gnerate the scatter context */
128: VecScatterCreate(gvec,from,baij->lvec,to,&baij->Mvctx);
130: /*
131: Post the receives for the first matrix vector product. We sync-chronize after
132: this on the chance that the user immediately calls MatMult() after assemblying
133: the matrix.
134: */
135: VecScatterPostRecvs(gvec,baij->lvec,INSERT_VALUES,SCATTER_FORWARD,baij->Mvctx);
136: MPI_Barrier(mat->comm);
138: PetscLogObjectParent(mat,baij->Mvctx);
139: PetscLogObjectParent(mat,baij->lvec);
140: PetscLogObjectParent(mat,from);
141: PetscLogObjectParent(mat,to);
142: baij->garray = garray;
143: PetscLogObjectMemory(mat,(ec+1)*sizeof(int));
144: ISDestroy(from);
145: ISDestroy(to);
146: VecDestroy(gvec);
147: return(0);
148: }
150: /*
151: Takes the local part of an already assembled MPIBAIJ matrix
152: and disassembles it. This is to allow new nonzeros into the matrix
153: that require more communication in the matrix vector multiply.
154: Thus certain data-structures must be rebuilt.
156: Kind of slow! But that's what application programmers get when
157: they are sloppy.
158: */
159: int DisAssemble_MPIBAIJ(Mat A)
160: {
161: Mat_MPIBAIJ *baij = (Mat_MPIBAIJ*)A->data;
162: Mat B = baij->B,Bnew;
163: Mat_SeqBAIJ *Bbaij = (Mat_SeqBAIJ*)B->data;
164: int ierr,i,j,mbs=Bbaij->mbs,n = A->N,col,*garray=baij->garray;
165: int bs2 = baij->bs2,*nz,ec,m = A->m;
166: MatScalar *a = Bbaij->a;
167: PetscScalar *atmp;
168: #if defined(PETSC_USE_MAT_SINGLE)
169: int k;
170: #endif
173: /* free stuff related to matrix-vec multiply */
174: VecGetSize(baij->lvec,&ec); /* needed for PetscLogObjectMemory below */
175: VecDestroy(baij->lvec); baij->lvec = 0;
176: VecScatterDestroy(baij->Mvctx); baij->Mvctx = 0;
177: if (baij->colmap) {
178: #if defined (PETSC_USE_CTABLE)
179: PetscTableDelete(baij->colmap); baij->colmap = 0;
180: #else
181: PetscFree(baij->colmap);
182: baij->colmap = 0;
183: PetscLogObjectMemory(A,-Bbaij->nbs*sizeof(int));
184: #endif
185: }
187: /* make sure that B is assembled so we can access its values */
188: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
189: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
191: /* invent new B and copy stuff over */
192: PetscMalloc(mbs*sizeof(int),&nz);
193: for (i=0; i<mbs; i++) {
194: nz[i] = Bbaij->i[i+1]-Bbaij->i[i];
195: }
196: MatCreateSeqBAIJ(PETSC_COMM_SELF,baij->bs,m,n,0,nz,&Bnew);
197: MatSetOption(Bnew,MAT_COLUMN_ORIENTED);
199: #if defined(PETSC_USE_MAT_SINGLE)
200: PetscMalloc(bs2*sizeof(PetscScalar),&atmp);
201: #endif
202: for (i=0; i<mbs; i++) {
203: for (j=Bbaij->i[i]; j<Bbaij->i[i+1]; j++) {
204: col = garray[Bbaij->j[j]];
205: #if defined(PETSC_USE_MAT_SINGLE)
206: for (k=0; k<bs2; k++) atmp[k] = a[j*bs2+k];
207: #else
208: atmp = a + j*bs2;
209: #endif
210: MatSetValuesBlocked_SeqBAIJ(Bnew,1,&i,1,&col,atmp,B->insertmode);
211: }
212: }
213: MatSetOption(Bnew,MAT_ROW_ORIENTED);
215: #if defined(PETSC_USE_MAT_SINGLE)
216: PetscFree(atmp);
217: #endif
219: PetscFree(nz);
220: PetscFree(baij->garray);
221: baij->garray = 0;
222: PetscLogObjectMemory(A,-ec*sizeof(int));
223: MatDestroy(B);
224: PetscLogObjectParent(A,Bnew);
225: baij->B = Bnew;
226: A->was_assembled = PETSC_FALSE;
227: return(0);
228: }
230: /* ugly stuff added for Glenn someday we should fix this up */
232: static int *uglyrmapd = 0,*uglyrmapo = 0; /* mapping from the local ordering to the "diagonal" and "off-diagonal"
233: parts of the local matrix */
234: static Vec uglydd = 0,uglyoo = 0; /* work vectors used to scale the two parts of the local matrix */
237: int MatMPIBAIJDiagonalScaleLocalSetUp(Mat inA,Vec scale)
238: {
239: Mat_MPIBAIJ *ina = (Mat_MPIBAIJ*) inA->data; /*access private part of matrix */
240: Mat_SeqBAIJ *A = (Mat_SeqBAIJ*)ina->A->data;
241: Mat_SeqBAIJ *B = (Mat_SeqBAIJ*)ina->B->data;
242: int ierr,bs = A->bs,i,n,nt,j,cstart,cend,no,*garray = ina->garray,*lindices;
243: int *r_rmapd,*r_rmapo;
244:
246: MatGetOwnershipRange(inA,&cstart,&cend);
247: MatGetSize(ina->A,PETSC_NULL,&n);
248: PetscMalloc((inA->bmapping->n+1)*sizeof(int),&r_rmapd);
249: PetscMemzero(r_rmapd,inA->bmapping->n*sizeof(int));
250: nt = 0;
251: for (i=0; i<inA->bmapping->n; i++) {
252: if (inA->bmapping->indices[i]*bs >= cstart && inA->bmapping->indices[i]*bs < cend) {
253: nt++;
254: r_rmapd[i] = inA->bmapping->indices[i] + 1;
255: }
256: }
257: if (nt*bs != n) SETERRQ2(1,"Hmm nt*bs %d n %d",nt*bs,n);
258: PetscMalloc((n+1)*sizeof(int),&uglyrmapd);
259: for (i=0; i<inA->bmapping->n; i++) {
260: if (r_rmapd[i]){
261: for (j=0; j<bs; j++) {
262: uglyrmapd[(r_rmapd[i]-1)*bs+j-cstart] = i*bs + j;
263: }
264: }
265: }
266: PetscFree(r_rmapd);
267: VecCreateSeq(PETSC_COMM_SELF,n,&uglydd);
269: PetscMalloc((ina->Nbs+1)*sizeof(int),&lindices);
270: PetscMemzero(lindices,ina->Nbs*sizeof(int));
271: for (i=0; i<B->nbs; i++) {
272: lindices[garray[i]] = i+1;
273: }
274: no = inA->bmapping->n - nt;
275: PetscMalloc((inA->bmapping->n+1)*sizeof(int),&r_rmapo);
276: PetscMemzero(r_rmapo,inA->bmapping->n*sizeof(int));
277: nt = 0;
278: for (i=0; i<inA->bmapping->n; i++) {
279: if (lindices[inA->bmapping->indices[i]]) {
280: nt++;
281: r_rmapo[i] = lindices[inA->bmapping->indices[i]];
282: }
283: }
284: if (nt > no) SETERRQ2(1,"Hmm nt %d no %d",nt,n);
285: PetscFree(lindices);
286: PetscMalloc((nt*bs+1)*sizeof(int),&uglyrmapo);
287: for (i=0; i<inA->bmapping->n; i++) {
288: if (r_rmapo[i]){
289: for (j=0; j<bs; j++) {
290: uglyrmapo[(r_rmapo[i]-1)*bs+j] = i*bs + j;
291: }
292: }
293: }
294: PetscFree(r_rmapo);
295: VecCreateSeq(PETSC_COMM_SELF,nt*bs,&uglyoo);
297: return(0);
298: }
300: int MatMPIBAIJDiagonalScaleLocal(Mat A,Vec scale)
301: {
302: Mat_MPIBAIJ *a = (Mat_MPIBAIJ*) A->data; /*access private part of matrix */
303: int ierr,n,i;
304: PetscScalar *d,*o,*s;
305:
307: if (!uglyrmapd) {
308: MatMPIBAIJDiagonalScaleLocalSetUp(A,scale);
309: }
311: VecGetArray(scale,&s);
312:
313: VecGetLocalSize(uglydd,&n);
314: VecGetArray(uglydd,&d);
315: for (i=0; i<n; i++) {
316: d[i] = s[uglyrmapd[i]]; /* copy "diagonal" (true local) portion of scale into dd vector */
317: }
318: VecRestoreArray(uglydd,&d);
319: /* column scale "diagonal" portion of local matrix */
320: MatDiagonalScale(a->A,PETSC_NULL,uglydd);
322: VecGetLocalSize(uglyoo,&n);
323: VecGetArray(uglyoo,&o);
324: for (i=0; i<n; i++) {
325: o[i] = s[uglyrmapo[i]]; /* copy "off-diagonal" portion of scale into oo vector */
326: }
327: VecRestoreArray(scale,&s);
328: VecRestoreArray(uglyoo,&o);
329: /* column scale "off-diagonal" portion of local matrix */
330: MatDiagonalScale(a->B,PETSC_NULL,uglyoo);
332: return(0);
333: }