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