Actual source code: mmbaij.c
petsc-3.3-p5 2012-12-01
2: /*
3: Support for the parallel BAIJ matrix vector multiply
4: */
5: #include <../src/mat/impls/baij/mpi/mpibaij.h>
7: extern PetscErrorCode MatSetValuesBlocked_SeqBAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
11: PetscErrorCode MatSetUpMultiply_MPIBAIJ(Mat mat)
12: {
13: Mat_MPIBAIJ *baij = (Mat_MPIBAIJ*)mat->data;
14: Mat_SeqBAIJ *B = (Mat_SeqBAIJ*)(baij->B->data);
16: PetscInt i,j,*aj = B->j,ec = 0,*garray;
17: PetscInt bs = mat->rmap->bs,*stmp;
18: IS from,to;
19: Vec gvec;
20: #if defined (PETSC_USE_CTABLE)
21: PetscTable gid1_lid1;
22: PetscTablePosition tpos;
23: PetscInt gid,lid;
24: #else
25: PetscInt Nbs = baij->Nbs,*indices;
26: #endif
30: #if defined (PETSC_USE_CTABLE)
31: /* use a table - Mark Adams */
32: PetscTableCreate(B->mbs,baij->Nbs+1,&gid1_lid1);
33: for (i=0; i<B->mbs; i++) {
34: for (j=0; j<B->ilen[i]; j++) {
35: PetscInt 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,INSERT_VALUES);
40: }
41: }
42: }
43: /* form array of columns we need */
44: PetscMalloc((ec+1)*sizeof(PetscInt),&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,INSERT_VALUES);
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: PetscInt 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->cmap->n = baij->B->cmap->N = ec*mat->rmap->bs;
67: PetscLayoutSetUp((baij->B->cmap));
68: PetscTableDestroy(&gid1_lid1);
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(PetscInt),&indices);
73: PetscMemzero(indices,Nbs*sizeof(PetscInt));
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(PetscInt),&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->cmap->n =baij->B->cmap->N = ec*mat->rmap->bs;
103: PetscLayoutSetUp((baij->B->cmap));
104: PetscFree(indices);
105: #endif
107: /* create local vector that is used to scatter into */
108: VecCreateSeq(PETSC_COMM_SELF,ec*bs,&baij->lvec);
110: /* create two temporary index sets for building scatter-gather */
111: ISCreateBlock(PETSC_COMM_SELF,bs,ec,garray,PETSC_COPY_VALUES,&from);
113: PetscMalloc((ec+1)*sizeof(PetscInt),&stmp);
114: for (i=0; i<ec; i++) { stmp[i] = i; }
115: ISCreateBlock(PETSC_COMM_SELF,bs,ec,stmp,PETSC_OWN_POINTER,&to);
117: /* create temporary global vector to generate scatter context */
118: VecCreateMPIWithArray(((PetscObject)mat)->comm,1,mat->cmap->n,mat->cmap->N,PETSC_NULL,&gvec);
120: VecScatterCreate(gvec,from,baij->lvec,to,&baij->Mvctx);
122: PetscLogObjectParent(mat,baij->Mvctx);
123: PetscLogObjectParent(mat,baij->lvec);
124: PetscLogObjectParent(mat,from);
125: PetscLogObjectParent(mat,to);
126: baij->garray = garray;
127: PetscLogObjectMemory(mat,(ec+1)*sizeof(PetscInt));
128: ISDestroy(&from);
129: ISDestroy(&to);
130: VecDestroy(&gvec);
131: return(0);
132: }
134: /*
135: Takes the local part of an already assembled MPIBAIJ matrix
136: and disassembles it. This is to allow new nonzeros into the matrix
137: that require more communication in the matrix vector multiply.
138: Thus certain data-structures must be rebuilt.
140: Kind of slow! But that's what application programmers get when
141: they are sloppy.
142: */
145: PetscErrorCode MatDisAssemble_MPIBAIJ(Mat A)
146: {
147: Mat_MPIBAIJ *baij = (Mat_MPIBAIJ*)A->data;
148: Mat B = baij->B,Bnew;
149: Mat_SeqBAIJ *Bbaij = (Mat_SeqBAIJ*)B->data;
151: PetscInt i,j,mbs=Bbaij->mbs,n = A->cmap->N,col,*garray=baij->garray;
152: PetscInt bs2 = baij->bs2,*nz,ec,m = A->rmap->n;
153: MatScalar *a = Bbaij->a;
154: MatScalar *atmp;
158: /* free stuff related to matrix-vec multiply */
159: VecGetSize(baij->lvec,&ec); /* needed for PetscLogObjectMemory below */
160: VecDestroy(&baij->lvec); baij->lvec = 0;
161: VecScatterDestroy(&baij->Mvctx); baij->Mvctx = 0;
162: if (baij->colmap) {
163: #if defined (PETSC_USE_CTABLE)
164: PetscTableDestroy(&baij->colmap);
165: #else
166: PetscFree(baij->colmap);
167: PetscLogObjectMemory(A,-Bbaij->nbs*sizeof(PetscInt));
168: #endif
169: }
171: /* make sure that B is assembled so we can access its values */
172: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
173: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
175: /* invent new B and copy stuff over */
176: PetscMalloc(mbs*sizeof(PetscInt),&nz);
177: for (i=0; i<mbs; i++) {
178: nz[i] = Bbaij->i[i+1]-Bbaij->i[i];
179: }
180: MatCreate(((PetscObject)B)->comm,&Bnew);
181: MatSetSizes(Bnew,m,n,m,n);
182: MatSetType(Bnew,((PetscObject)B)->type_name);
183: MatSeqBAIJSetPreallocation(Bnew,B->rmap->bs,0,nz);
184: ((Mat_SeqBAIJ*)Bnew->data)->nonew = Bbaij->nonew; /* Inherit insertion error options. */
185: MatSetOption(Bnew,MAT_ROW_ORIENTED,PETSC_FALSE);
187: for (i=0; i<mbs; i++) {
188: for (j=Bbaij->i[i]; j<Bbaij->i[i+1]; j++) {
189: col = garray[Bbaij->j[j]];
190: atmp = a + j*bs2;
191: MatSetValuesBlocked_SeqBAIJ(Bnew,1,&i,1,&col,atmp,B->insertmode);
192: }
193: }
194: MatSetOption(Bnew,MAT_ROW_ORIENTED,PETSC_TRUE);
196: PetscFree(nz);
197: PetscFree(baij->garray);
198: PetscLogObjectMemory(A,-ec*sizeof(PetscInt));
199: MatDestroy(&B);
200: PetscLogObjectParent(A,Bnew);
201: baij->B = Bnew;
202: A->was_assembled = PETSC_FALSE;
203: A->assembled = PETSC_FALSE;
204: return(0);
205: }
207: /* ugly stuff added for Glenn someday we should fix this up */
209: static PetscInt *uglyrmapd = 0,*uglyrmapo = 0; /* mapping from the local ordering to the "diagonal" and "off-diagonal"
210: parts of the local matrix */
211: static Vec uglydd = 0,uglyoo = 0; /* work vectors used to scale the two parts of the local matrix */
216: PetscErrorCode MatMPIBAIJDiagonalScaleLocalSetUp(Mat inA,Vec scale)
217: {
218: Mat_MPIBAIJ *ina = (Mat_MPIBAIJ*) inA->data; /*access private part of matrix */
219: Mat_SeqBAIJ *B = (Mat_SeqBAIJ*)ina->B->data;
221: PetscInt bs = inA->rmap->bs,i,n,nt,j,cstart,cend,no,*garray = ina->garray,*lindices;
222: PetscInt *r_rmapd,*r_rmapo;
223:
225: MatGetOwnershipRange(inA,&cstart,&cend);
226: MatGetSize(ina->A,PETSC_NULL,&n);
227: PetscMalloc((inA->rmap->bmapping->n+1)*sizeof(PetscInt),&r_rmapd);
228: PetscMemzero(r_rmapd,inA->rmap->bmapping->n*sizeof(PetscInt));
229: nt = 0;
230: for (i=0; i<inA->rmap->bmapping->n; i++) {
231: if (inA->rmap->bmapping->indices[i]*bs >= cstart && inA->rmap->bmapping->indices[i]*bs < cend) {
232: nt++;
233: r_rmapd[i] = inA->rmap->bmapping->indices[i] + 1;
234: }
235: }
236: if (nt*bs != n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hmm nt*bs %D n %D",nt*bs,n);
237: PetscMalloc((n+1)*sizeof(PetscInt),&uglyrmapd);
238: for (i=0; i<inA->rmap->bmapping->n; i++) {
239: if (r_rmapd[i]){
240: for (j=0; j<bs; j++) {
241: uglyrmapd[(r_rmapd[i]-1)*bs+j-cstart] = i*bs + j;
242: }
243: }
244: }
245: PetscFree(r_rmapd);
246: VecCreateSeq(PETSC_COMM_SELF,n,&uglydd);
248: PetscMalloc((ina->Nbs+1)*sizeof(PetscInt),&lindices);
249: PetscMemzero(lindices,ina->Nbs*sizeof(PetscInt));
250: for (i=0; i<B->nbs; i++) {
251: lindices[garray[i]] = i+1;
252: }
253: no = inA->rmap->bmapping->n - nt;
254: PetscMalloc((inA->rmap->bmapping->n+1)*sizeof(PetscInt),&r_rmapo);
255: PetscMemzero(r_rmapo,inA->rmap->bmapping->n*sizeof(PetscInt));
256: nt = 0;
257: for (i=0; i<inA->rmap->bmapping->n; i++) {
258: if (lindices[inA->rmap->bmapping->indices[i]]) {
259: nt++;
260: r_rmapo[i] = lindices[inA->rmap->bmapping->indices[i]];
261: }
262: }
263: if (nt > no) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hmm nt %D no %D",nt,n);
264: PetscFree(lindices);
265: PetscMalloc((nt*bs+1)*sizeof(PetscInt),&uglyrmapo);
266: for (i=0; i<inA->rmap->bmapping->n; i++) {
267: if (r_rmapo[i]){
268: for (j=0; j<bs; j++) {
269: uglyrmapo[(r_rmapo[i]-1)*bs+j] = i*bs + j;
270: }
271: }
272: }
273: PetscFree(r_rmapo);
274: VecCreateSeq(PETSC_COMM_SELF,nt*bs,&uglyoo);
276: return(0);
277: }
281: PetscErrorCode MatMPIBAIJDiagonalScaleLocal(Mat A,Vec scale)
282: {
283: /* This routine should really be abandoned as it duplicates MatDiagonalScaleLocal */
287: PetscTryMethod(A,"MatDiagonalScaleLocal_C",(Mat,Vec),(A,scale));
288: return(0);
289: }
291: EXTERN_C_BEGIN
294: PetscErrorCode MatDiagonalScaleLocal_MPIBAIJ(Mat A,Vec scale)
295: {
296: Mat_MPIBAIJ *a = (Mat_MPIBAIJ*) A->data; /*access private part of matrix */
298: PetscInt n,i;
299: PetscScalar *d,*o,*s;
300:
302: if (!uglyrmapd) {
303: MatMPIBAIJDiagonalScaleLocalSetUp(A,scale);
304: }
306: VecGetArray(scale,&s);
307:
308: VecGetLocalSize(uglydd,&n);
309: VecGetArray(uglydd,&d);
310: for (i=0; i<n; i++) {
311: d[i] = s[uglyrmapd[i]]; /* copy "diagonal" (true local) portion of scale into dd vector */
312: }
313: VecRestoreArray(uglydd,&d);
314: /* column scale "diagonal" portion of local matrix */
315: MatDiagonalScale(a->A,PETSC_NULL,uglydd);
317: VecGetLocalSize(uglyoo,&n);
318: VecGetArray(uglyoo,&o);
319: for (i=0; i<n; i++) {
320: o[i] = s[uglyrmapo[i]]; /* copy "off-diagonal" portion of scale into oo vector */
321: }
322: VecRestoreArray(scale,&s);
323: VecRestoreArray(uglyoo,&o);
324: /* column scale "off-diagonal" portion of local matrix */
325: MatDiagonalScale(a->B,PETSC_NULL,uglyoo);
327: return(0);
328: }
329: EXTERN_C_END