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
8: EXTERN int MatSetValuesBlocked_SeqBAIJ(Mat,int,const int[],int,const int[],const PetscScalar[],InsertMode);
12: int MatSetUpMultiply_MPIBAIJ(Mat mat)
13: {
14: Mat_MPIBAIJ *baij = (Mat_MPIBAIJ*)mat->data;
15: Mat_SeqBAIJ *B = (Mat_SeqBAIJ*)(baij->B->data);
16: int i,j,*aj = B->j,ierr,ec = 0,*garray;
17: int bs = baij->bs,*stmp;
18: IS from,to;
19: Vec gvec;
20: #if defined (PETSC_USE_CTABLE)
21: PetscTable gid1_lid1;
22: PetscTablePosition tpos;
23: int gid,lid;
24: #else
25: int Nbs = baij->Nbs,*indices;
26: #endif
30: #if defined (PETSC_USE_CTABLE)
31: /* use a table - Mark Adams */
32: PetscTableCreate(B->mbs,&gid1_lid1);
33: for (i=0; i<B->mbs; i++) {
34: for (j=0; j<B->ilen[i]; j++) {
35: int data,gid1 = aj[B->i[i]+j] + 1;
36: PetscTableFind(gid1_lid1,gid1,&data) ;
37: if (!data) {
38: /* one based table */
39: PetscTableAdd(gid1_lid1,gid1,++ec);
40: }
41: }
42: }
43: /* form array of columns we need */
44: PetscMalloc((ec+1)*sizeof(int),&garray);
45: PetscTableGetHeadPosition(gid1_lid1,&tpos);
46: while (tpos) {
47: PetscTableGetNext(gid1_lid1,&tpos,&gid,&lid);
48: gid--; lid--;
49: garray[lid] = gid;
50: }
51: PetscSortInt(ec,garray);
52: PetscTableRemoveAll(gid1_lid1);
53: for (i=0; i<ec; i++) {
54: PetscTableAdd(gid1_lid1,garray[i]+1,i+1);
55: }
56: /* compact out the extra columns in B */
57: for (i=0; i<B->mbs; i++) {
58: for (j=0; j<B->ilen[i]; j++) {
59: int gid1 = aj[B->i[i] + j] + 1;
60: PetscTableFind(gid1_lid1,gid1,&lid);
61: lid --;
62: aj[B->i[i]+j] = lid;
63: }
64: }
65: B->nbs = ec;
66: baij->B->n = ec*B->bs;
67: PetscTableDelete(gid1_lid1);
68: /* Mark Adams */
69: #else
70: /* Make an array as long as the number of columns */
71: /* mark those columns that are in baij->B */
72: PetscMalloc((Nbs+1)*sizeof(int),&indices);
73: PetscMemzero(indices,Nbs*sizeof(int));
74: for (i=0; i<B->mbs; i++) {
75: for (j=0; j<B->ilen[i]; j++) {
76: if (!indices[aj[B->i[i] + j]]) ec++;
77: indices[aj[B->i[i] + j]] = 1;
78: }
79: }
81: /* form array of columns we need */
82: PetscMalloc((ec+1)*sizeof(int),&garray);
83: ec = 0;
84: for (i=0; i<Nbs; i++) {
85: if (indices[i]) {
86: garray[ec++] = i;
87: }
88: }
90: /* make indices now point into garray */
91: for (i=0; i<ec; i++) {
92: indices[garray[i]] = i;
93: }
95: /* compact out the extra columns in B */
96: for (i=0; i<B->mbs; i++) {
97: for (j=0; j<B->ilen[i]; j++) {
98: aj[B->i[i] + j] = indices[aj[B->i[i] + j]];
99: }
100: }
101: B->nbs = ec;
102: baij->B->n = ec*B->bs;
103: PetscFree(indices);
104: #endif
106: /* create local vector that is used to scatter into */
107: VecCreateSeq(PETSC_COMM_SELF,ec*bs,&baij->lvec);
109: /* create two temporary index sets for building scatter-gather */
110: for (i=0; i<ec; i++) {
111: garray[i] = bs*garray[i];
112: }
113: ISCreateBlock(PETSC_COMM_SELF,bs,ec,garray,&from);
114: for (i=0; i<ec; i++) {
115: garray[i] = garray[i]/bs;
116: }
118: PetscMalloc((ec+1)*sizeof(int),&stmp);
119: for (i=0; i<ec; i++) { stmp[i] = bs*i; }
120: ISCreateBlock(PETSC_COMM_SELF,bs,ec,stmp,&to);
121: PetscFree(stmp);
123: /* create temporary global vector to generate scatter context */
124: /* this is inefficient, but otherwise we must do either
125: 1) save garray until the first actual scatter when the vector is known or
126: 2) have another way of generating a scatter context without a vector.*/
127: VecCreateMPI(mat->comm,mat->n,mat->N,&gvec);
129: /* gnerate the scatter context */
130: VecScatterCreate(gvec,from,baij->lvec,to,&baij->Mvctx);
132: /*
133: Post the receives for the first matrix vector product. We sync-chronize after
134: this on the chance that the user immediately calls MatMult() after assemblying
135: the matrix.
136: */
137: VecScatterPostRecvs(gvec,baij->lvec,INSERT_VALUES,SCATTER_FORWARD,baij->Mvctx);
138: MPI_Barrier(mat->comm);
140: PetscLogObjectParent(mat,baij->Mvctx);
141: PetscLogObjectParent(mat,baij->lvec);
142: PetscLogObjectParent(mat,from);
143: PetscLogObjectParent(mat,to);
144: baij->garray = garray;
145: PetscLogObjectMemory(mat,(ec+1)*sizeof(int));
146: ISDestroy(from);
147: ISDestroy(to);
148: VecDestroy(gvec);
149: return(0);
150: }
152: /*
153: Takes the local part of an already assembled MPIBAIJ matrix
154: and disassembles it. This is to allow new nonzeros into the matrix
155: that require more communication in the matrix vector multiply.
156: Thus certain data-structures must be rebuilt.
158: Kind of slow! But that's what application programmers get when
159: they are sloppy.
160: */
163: int DisAssemble_MPIBAIJ(Mat A)
164: {
165: Mat_MPIBAIJ *baij = (Mat_MPIBAIJ*)A->data;
166: Mat B = baij->B,Bnew;
167: Mat_SeqBAIJ *Bbaij = (Mat_SeqBAIJ*)B->data;
168: int ierr,i,j,mbs=Bbaij->mbs,n = A->N,col,*garray=baij->garray;
169: int bs2 = baij->bs2,*nz,ec,m = A->m;
170: MatScalar *a = Bbaij->a;
171: PetscScalar *atmp;
172: #if defined(PETSC_USE_MAT_SINGLE)
173: int k;
174: #endif
177: /* free stuff related to matrix-vec multiply */
178: VecGetSize(baij->lvec,&ec); /* needed for PetscLogObjectMemory below */
179: VecDestroy(baij->lvec); baij->lvec = 0;
180: VecScatterDestroy(baij->Mvctx); baij->Mvctx = 0;
181: if (baij->colmap) {
182: #if defined (PETSC_USE_CTABLE)
183: PetscTableDelete(baij->colmap); baij->colmap = 0;
184: #else
185: PetscFree(baij->colmap);
186: baij->colmap = 0;
187: PetscLogObjectMemory(A,-Bbaij->nbs*sizeof(int));
188: #endif
189: }
191: /* make sure that B is assembled so we can access its values */
192: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
193: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
195: /* invent new B and copy stuff over */
196: PetscMalloc(mbs*sizeof(int),&nz);
197: for (i=0; i<mbs; i++) {
198: nz[i] = Bbaij->i[i+1]-Bbaij->i[i];
199: }
200: MatCreate(B->comm,m,n,m,n,&Bnew);
201: MatSetType(Bnew,B->type_name);
202: MatSeqBAIJSetPreallocation(Bnew,baij->bs,0,nz);
203: MatSetOption(Bnew,MAT_COLUMN_ORIENTED);
205: #if defined(PETSC_USE_MAT_SINGLE)
206: PetscMalloc(bs2*sizeof(PetscScalar),&atmp);
207: #endif
208: for (i=0; i<mbs; i++) {
209: for (j=Bbaij->i[i]; j<Bbaij->i[i+1]; j++) {
210: col = garray[Bbaij->j[j]];
211: #if defined(PETSC_USE_MAT_SINGLE)
212: for (k=0; k<bs2; k++) atmp[k] = a[j*bs2+k];
213: #else
214: atmp = a + j*bs2;
215: #endif
216: MatSetValuesBlocked_SeqBAIJ(Bnew,1,&i,1,&col,atmp,B->insertmode);
217: }
218: }
219: MatSetOption(Bnew,MAT_ROW_ORIENTED);
221: #if defined(PETSC_USE_MAT_SINGLE)
222: PetscFree(atmp);
223: #endif
225: PetscFree(nz);
226: PetscFree(baij->garray);
227: baij->garray = 0;
228: PetscLogObjectMemory(A,-ec*sizeof(int));
229: MatDestroy(B);
230: PetscLogObjectParent(A,Bnew);
231: baij->B = Bnew;
232: A->was_assembled = PETSC_FALSE;
233: return(0);
234: }
236: /* ugly stuff added for Glenn someday we should fix this up */
238: static int *uglyrmapd = 0,*uglyrmapo = 0; /* mapping from the local ordering to the "diagonal" and "off-diagonal"
239: parts of the local matrix */
240: static Vec uglydd = 0,uglyoo = 0; /* work vectors used to scale the two parts of the local matrix */
245: int MatMPIBAIJDiagonalScaleLocalSetUp(Mat inA,Vec scale)
246: {
247: Mat_MPIBAIJ *ina = (Mat_MPIBAIJ*) inA->data; /*access private part of matrix */
248: Mat_SeqBAIJ *A = (Mat_SeqBAIJ*)ina->A->data;
249: Mat_SeqBAIJ *B = (Mat_SeqBAIJ*)ina->B->data;
250: int ierr,bs = A->bs,i,n,nt,j,cstart,cend,no,*garray = ina->garray,*lindices;
251: int *r_rmapd,*r_rmapo;
252:
254: MatGetOwnershipRange(inA,&cstart,&cend);
255: MatGetSize(ina->A,PETSC_NULL,&n);
256: PetscMalloc((inA->bmapping->n+1)*sizeof(int),&r_rmapd);
257: PetscMemzero(r_rmapd,inA->bmapping->n*sizeof(int));
258: nt = 0;
259: for (i=0; i<inA->bmapping->n; i++) {
260: if (inA->bmapping->indices[i]*bs >= cstart && inA->bmapping->indices[i]*bs < cend) {
261: nt++;
262: r_rmapd[i] = inA->bmapping->indices[i] + 1;
263: }
264: }
265: if (nt*bs != n) SETERRQ2(1,"Hmm nt*bs %d n %d",nt*bs,n);
266: PetscMalloc((n+1)*sizeof(int),&uglyrmapd);
267: for (i=0; i<inA->bmapping->n; i++) {
268: if (r_rmapd[i]){
269: for (j=0; j<bs; j++) {
270: uglyrmapd[(r_rmapd[i]-1)*bs+j-cstart] = i*bs + j;
271: }
272: }
273: }
274: PetscFree(r_rmapd);
275: VecCreateSeq(PETSC_COMM_SELF,n,&uglydd);
277: PetscMalloc((ina->Nbs+1)*sizeof(int),&lindices);
278: PetscMemzero(lindices,ina->Nbs*sizeof(int));
279: for (i=0; i<B->nbs; i++) {
280: lindices[garray[i]] = i+1;
281: }
282: no = inA->bmapping->n - nt;
283: PetscMalloc((inA->bmapping->n+1)*sizeof(int),&r_rmapo);
284: PetscMemzero(r_rmapo,inA->bmapping->n*sizeof(int));
285: nt = 0;
286: for (i=0; i<inA->bmapping->n; i++) {
287: if (lindices[inA->bmapping->indices[i]]) {
288: nt++;
289: r_rmapo[i] = lindices[inA->bmapping->indices[i]];
290: }
291: }
292: if (nt > no) SETERRQ2(1,"Hmm nt %d no %d",nt,n);
293: PetscFree(lindices);
294: PetscMalloc((nt*bs+1)*sizeof(int),&uglyrmapo);
295: for (i=0; i<inA->bmapping->n; i++) {
296: if (r_rmapo[i]){
297: for (j=0; j<bs; j++) {
298: uglyrmapo[(r_rmapo[i]-1)*bs+j] = i*bs + j;
299: }
300: }
301: }
302: PetscFree(r_rmapo);
303: VecCreateSeq(PETSC_COMM_SELF,nt*bs,&uglyoo);
305: return(0);
306: }
310: int MatMPIBAIJDiagonalScaleLocal(Mat A,Vec scale)
311: {
312: /* This routine should really be abandoned as it duplicates MatDiagonalScaleLocal */
313: int ierr,(*f)(Mat,Vec);
316: PetscObjectQueryFunction((PetscObject)A,"MatDiagonalScaleLocal_C",(void (**)(void))&f);
317: if (f) {
318: (*f)(A,scale);
319: }
320: return(0);
321: }
323: EXTERN_C_BEGIN
326: int MatDiagonalScaleLocal_MPIBAIJ(Mat A,Vec scale)
327: {
328: Mat_MPIBAIJ *a = (Mat_MPIBAIJ*) A->data; /*access private part of matrix */
329: int ierr,n,i;
330: PetscScalar *d,*o,*s;
331:
333: if (!uglyrmapd) {
334: MatMPIBAIJDiagonalScaleLocalSetUp(A,scale);
335: }
337: VecGetArray(scale,&s);
338:
339: VecGetLocalSize(uglydd,&n);
340: VecGetArray(uglydd,&d);
341: for (i=0; i<n; i++) {
342: d[i] = s[uglyrmapd[i]]; /* copy "diagonal" (true local) portion of scale into dd vector */
343: }
344: VecRestoreArray(uglydd,&d);
345: /* column scale "diagonal" portion of local matrix */
346: MatDiagonalScale(a->A,PETSC_NULL,uglydd);
348: VecGetLocalSize(uglyoo,&n);
349: VecGetArray(uglyoo,&o);
350: for (i=0; i<n; i++) {
351: o[i] = s[uglyrmapo[i]]; /* copy "off-diagonal" portion of scale into oo vector */
352: }
353: VecRestoreArray(scale,&s);
354: VecRestoreArray(uglyoo,&o);
355: /* column scale "off-diagonal" portion of local matrix */
356: MatDiagonalScale(a->B,PETSC_NULL,uglyoo);
358: return(0);
359: }
360: EXTERN_C_END