Actual source code: mmbaij.c
1: /*
2: Support for the parallel BAIJ matrix vector multiply
3: */
4: #include src/mat/impls/baij/mpi/mpibaij.h
6: EXTERN PetscErrorCode MatSetValuesBlocked_SeqBAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
10: PetscErrorCode MatSetUpMultiply_MPIBAIJ(Mat mat)
11: {
12: Mat_MPIBAIJ *baij = (Mat_MPIBAIJ*)mat->data;
13: Mat_SeqBAIJ *B = (Mat_SeqBAIJ*)(baij->B->data);
15: PetscInt i,j,*aj = B->j,ec = 0,*garray;
16: PetscInt bs = mat->bs,*stmp;
17: IS from,to;
18: Vec gvec;
19: #if defined (PETSC_USE_CTABLE)
20: PetscTable gid1_lid1;
21: PetscTablePosition tpos;
22: PetscInt gid,lid;
23: #else
24: PetscInt Nbs = baij->Nbs,*indices;
25: #endif
29: #if defined (PETSC_USE_CTABLE)
30: /* use a table - Mark Adams */
31: PetscTableCreate(B->mbs,&gid1_lid1);
32: for (i=0; i<B->mbs; i++) {
33: for (j=0; j<B->ilen[i]; j++) {
34: PetscInt data,gid1 = aj[B->i[i]+j] + 1;
35: PetscTableFind(gid1_lid1,gid1,&data);
36: if (!data) {
37: /* one based table */
38: PetscTableAdd(gid1_lid1,gid1,++ec);
39: }
40: }
41: }
42: /* form array of columns we need */
43: PetscMalloc((ec+1)*sizeof(PetscInt),&garray);
44: PetscTableGetHeadPosition(gid1_lid1,&tpos);
45: while (tpos) {
46: PetscTableGetNext(gid1_lid1,&tpos,&gid,&lid);
47: gid--; lid--;
48: garray[lid] = gid;
49: }
50: PetscSortInt(ec,garray);
51: PetscTableRemoveAll(gid1_lid1);
52: for (i=0; i<ec; i++) {
53: PetscTableAdd(gid1_lid1,garray[i]+1,i+1);
54: }
55: /* compact out the extra columns in B */
56: for (i=0; i<B->mbs; i++) {
57: for (j=0; j<B->ilen[i]; j++) {
58: PetscInt gid1 = aj[B->i[i] + j] + 1;
59: PetscTableFind(gid1_lid1,gid1,&lid);
60: lid --;
61: aj[B->i[i]+j] = lid;
62: }
63: }
64: B->nbs = ec;
65: baij->B->n = ec*mat->bs;
66: PetscTableDelete(gid1_lid1);
67: /* Mark Adams */
68: #else
69: /* Make an array as long as the number of columns */
70: /* mark those columns that are in baij->B */
71: PetscMalloc((Nbs+1)*sizeof(PetscInt),&indices);
72: PetscMemzero(indices,Nbs*sizeof(PetscInt));
73: for (i=0; i<B->mbs; i++) {
74: for (j=0; j<B->ilen[i]; j++) {
75: if (!indices[aj[B->i[i] + j]]) ec++;
76: indices[aj[B->i[i] + j]] = 1;
77: }
78: }
80: /* form array of columns we need */
81: PetscMalloc((ec+1)*sizeof(PetscInt),&garray);
82: ec = 0;
83: for (i=0; i<Nbs; i++) {
84: if (indices[i]) {
85: garray[ec++] = i;
86: }
87: }
89: /* make indices now point into garray */
90: for (i=0; i<ec; i++) {
91: indices[garray[i]] = i;
92: }
94: /* compact out the extra columns in B */
95: for (i=0; i<B->mbs; i++) {
96: for (j=0; j<B->ilen[i]; j++) {
97: aj[B->i[i] + j] = indices[aj[B->i[i] + j]];
98: }
99: }
100: B->nbs = ec;
101: baij->B->n = ec*mat->bs;
102: PetscFree(indices);
103: #endif
105: /* create local vector that is used to scatter into */
106: VecCreateSeq(PETSC_COMM_SELF,ec*bs,&baij->lvec);
108: /* create two temporary index sets for building scatter-gather */
109: for (i=0; i<ec; i++) {
110: garray[i] = bs*garray[i];
111: }
112: ISCreateBlock(PETSC_COMM_SELF,bs,ec,garray,&from);
113: for (i=0; i<ec; i++) {
114: garray[i] = garray[i]/bs;
115: }
117: PetscMalloc((ec+1)*sizeof(PetscInt),&stmp);
118: for (i=0; i<ec; i++) { stmp[i] = bs*i; }
119: ISCreateBlock(PETSC_COMM_SELF,bs,ec,stmp,&to);
120: PetscFree(stmp);
122: /* create temporary global vector to generate scatter context */
123: /* this is inefficient, but otherwise we must do either
124: 1) save garray until the first actual scatter when the vector is known or
125: 2) have another way of generating a scatter context without a vector.*/
126: VecCreateMPI(mat->comm,mat->n,mat->N,&gvec);
128: /* gnerate the scatter context */
129: VecScatterCreate(gvec,from,baij->lvec,to,&baij->Mvctx);
131: /*
132: Post the receives for the first matrix vector product. We sync-chronize after
133: this on the chance that the user immediately calls MatMult() after assemblying
134: the matrix.
135: */
136: VecScatterPostRecvs(gvec,baij->lvec,INSERT_VALUES,SCATTER_FORWARD,baij->Mvctx);
137: MPI_Barrier(mat->comm);
139: PetscLogObjectParent(mat,baij->Mvctx);
140: PetscLogObjectParent(mat,baij->lvec);
141: PetscLogObjectParent(mat,from);
142: PetscLogObjectParent(mat,to);
143: baij->garray = garray;
144: PetscLogObjectMemory(mat,(ec+1)*sizeof(PetscInt));
145: ISDestroy(from);
146: ISDestroy(to);
147: VecDestroy(gvec);
148: return(0);
149: }
151: /*
152: Takes the local part of an already assembled MPIBAIJ matrix
153: and disassembles it. This is to allow new nonzeros into the matrix
154: that require more communication in the matrix vector multiply.
155: Thus certain data-structures must be rebuilt.
157: Kind of slow! But that's what application programmers get when
158: they are sloppy.
159: */
162: PetscErrorCode DisAssemble_MPIBAIJ(Mat A)
163: {
164: Mat_MPIBAIJ *baij = (Mat_MPIBAIJ*)A->data;
165: Mat B = baij->B,Bnew;
166: Mat_SeqBAIJ *Bbaij = (Mat_SeqBAIJ*)B->data;
168: PetscInt i,j,mbs=Bbaij->mbs,n = A->N,col,*garray=baij->garray;
169: PetscInt bs2 = baij->bs2,*nz,ec,m = A->m;
170: MatScalar *a = Bbaij->a;
171: PetscScalar *atmp;
172: #if defined(PETSC_USE_MAT_SINGLE)
173: PetscInt 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(PetscInt));
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(PetscInt),&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,B->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(PetscInt));
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 PetscInt *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: PetscErrorCode MatMPIBAIJDiagonalScaleLocalSetUp(Mat inA,Vec scale)
246: {
247: Mat_MPIBAIJ *ina = (Mat_MPIBAIJ*) inA->data; /*access private part of matrix */
248: Mat_SeqBAIJ *B = (Mat_SeqBAIJ*)ina->B->data;
250: PetscInt bs = inA->bs,i,n,nt,j,cstart,cend,no,*garray = ina->garray,*lindices;
251: PetscInt *r_rmapd,*r_rmapo;
252:
254: MatGetOwnershipRange(inA,&cstart,&cend);
255: MatGetSize(ina->A,PETSC_NULL,&n);
256: PetscMalloc((inA->bmapping->n+1)*sizeof(PetscInt),&r_rmapd);
257: PetscMemzero(r_rmapd,inA->bmapping->n*sizeof(PetscInt));
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(PETSC_ERR_PLIB,"Hmm nt*bs %D n %D",nt*bs,n);
266: PetscMalloc((n+1)*sizeof(PetscInt),&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(PetscInt),&lindices);
278: PetscMemzero(lindices,ina->Nbs*sizeof(PetscInt));
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(PetscInt),&r_rmapo);
284: PetscMemzero(r_rmapo,inA->bmapping->n*sizeof(PetscInt));
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(PETSC_ERR_PLIB,"Hmm nt %D no %D",nt,n);
293: PetscFree(lindices);
294: PetscMalloc((nt*bs+1)*sizeof(PetscInt),&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: PetscErrorCode MatMPIBAIJDiagonalScaleLocal(Mat A,Vec scale)
311: {
312: /* This routine should really be abandoned as it duplicates MatDiagonalScaleLocal */
313: PetscErrorCode 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: }
326: PetscErrorCode MatDiagonalScaleLocal_MPIBAIJ(Mat A,Vec scale)
327: {
328: Mat_MPIBAIJ *a = (Mat_MPIBAIJ*) A->data; /*access private part of matrix */
330: PetscInt n,i;
331: PetscScalar *d,*o,*s;
332:
334: if (!uglyrmapd) {
335: MatMPIBAIJDiagonalScaleLocalSetUp(A,scale);
336: }
338: VecGetArray(scale,&s);
339:
340: VecGetLocalSize(uglydd,&n);
341: VecGetArray(uglydd,&d);
342: for (i=0; i<n; i++) {
343: d[i] = s[uglyrmapd[i]]; /* copy "diagonal" (true local) portion of scale into dd vector */
344: }
345: VecRestoreArray(uglydd,&d);
346: /* column scale "diagonal" portion of local matrix */
347: MatDiagonalScale(a->A,PETSC_NULL,uglydd);
349: VecGetLocalSize(uglyoo,&n);
350: VecGetArray(uglyoo,&o);
351: for (i=0; i<n; i++) {
352: o[i] = s[uglyrmapo[i]]; /* copy "off-diagonal" portion of scale into oo vector */
353: }
354: VecRestoreArray(scale,&s);
355: VecRestoreArray(uglyoo,&o);
356: /* column scale "off-diagonal" portion of local matrix */
357: MatDiagonalScale(a->B,PETSC_NULL,uglyoo);
359: return(0);
360: }