Actual source code: mmaij.c

  1: /*$Id: mmaij.c,v 1.59 2001/08/07 03:02:49 balay Exp $*/

  3: /*
  4:    Support for the parallel AIJ matrix vector multiply
  5: */
 6:  #include src/mat/impls/aij/mpi/mpiaij.h
 7:  #include src/vec/vecimpl.h

  9: #undef __FUNCT__  
 11: int MatSetUpMultiply_MPIAIJ(Mat mat)
 12: {
 13:   Mat_MPIAIJ         *aij = (Mat_MPIAIJ*)mat->data;
 14:   Mat_SeqAIJ         *B = (Mat_SeqAIJ*)(aij->B->data);
 15:   int                N = mat->N,i,j,*indices,*aj = B->j,ierr,ec = 0,*garray;
 16:   int                shift = B->indexshift;
 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 (this has not been tested with "shift") */
 29:   PetscTableCreate(aij->B->m,&gid1_lid1);
 30:   for (i=0; i<aij->B->m; i++) {
 31:     for (j=0; j<B->ilen[i]; j++) {
 32:       int data,gid1 = aj[B->i[i] + shift + j] + 1 + shift;
 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--;
 46:     lid--;
 47:     garray[lid] = gid;
 48:   }
 49:   PetscSortInt(ec,garray); /* sort, and rebuild */
 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<aij->B->m; i++) {
 56:     for (j=0; j<B->ilen[i]; j++) {
 57:       int gid1 = aj[B->i[i] + shift + j] + 1 + shift;
 58:       PetscTableFind(gid1_lid1,gid1,&lid);
 59:       lid --;
 60:       aj[B->i[i] + shift + j]  = lid - shift;
 61:     }
 62:   }
 63:   aij->B->n = aij->B->N = ec;
 64:   PetscTableDelete(gid1_lid1);
 65:   /* Mark Adams */
 66: #else
 67:   /* For the first stab we make an array as long as the number of columns */
 68:   /* mark those columns that are in aij->B */
 69:   PetscMalloc((N+1)*sizeof(int),&indices);
 70:   PetscMemzero(indices,N*sizeof(int));
 71:   for (i=0; i<aij->B->m; i++) {
 72:     for (j=0; j<B->ilen[i]; j++) {
 73:       if (!indices[aj[B->i[i] +shift + j] + shift]) ec++;
 74:       indices[aj[B->i[i] + shift + j] + shift] = 1;
 75:     }
 76:   }

 78:   /* form array of columns we need */
 79:   PetscMalloc((ec+1)*sizeof(int),&garray);
 80:   ec = 0;
 81:   for (i=0; i<N; i++) {
 82:     if (indices[i]) garray[ec++] = i;
 83:   }

 85:   /* make indices now point into garray */
 86:   for (i=0; i<ec; i++) {
 87:     indices[garray[i]] = i-shift;
 88:   }

 90:   /* compact out the extra columns in B */
 91:   for (i=0; i<aij->B->m; i++) {
 92:     for (j=0; j<B->ilen[i]; j++) {
 93:       aj[B->i[i] + shift + j] = indices[aj[B->i[i] + shift + j]+shift];
 94:     }
 95:   }
 96:   aij->B->n = aij->B->N = ec;
 97:   PetscFree(indices);
 98: #endif  
 99:   /* create local vector that is used to scatter into */
100:   VecCreateSeq(PETSC_COMM_SELF,ec,&aij->lvec);

102:   /* create two temporary Index sets for build scatter gather */
103:   ISCreateGeneral(mat->comm,ec,garray,&from);
104:   ISCreateStride(PETSC_COMM_SELF,ec,0,1,&to);

106:   /* create temporary global vector to generate scatter context */
107:   /* this is inefficient, but otherwise we must do either 
108:      1) save garray until the first actual scatter when the vector is known or
109:      2) have another way of generating a scatter context without a vector.*/
110:   VecCreateMPI(mat->comm,mat->n,mat->N,&gvec);

112:   /* generate the scatter context */
113:   VecScatterCreate(gvec,from,aij->lvec,to,&aij->Mvctx);
114:   PetscLogObjectParent(mat,aij->Mvctx);
115:   PetscLogObjectParent(mat,aij->lvec);
116:   PetscLogObjectParent(mat,from);
117:   PetscLogObjectParent(mat,to);
118:   aij->garray = garray;
119:   PetscLogObjectMemory(mat,(ec+1)*sizeof(int));
120:   ISDestroy(from);
121:   ISDestroy(to);
122:   VecDestroy(gvec);
123:   return(0);
124: }


127: #undef __FUNCT__  
129: /*
130:      Takes the local part of an already assembled MPIAIJ matrix
131:    and disassembles it. This is to allow new nonzeros into the matrix
132:    that require more communication in the matrix vector multiply. 
133:    Thus certain data-structures must be rebuilt.

135:    Kind of slow! But that's what application programmers get when 
136:    they are sloppy.
137: */
138: int DisAssemble_MPIAIJ(Mat A)
139: {
140:   Mat_MPIAIJ   *aij = (Mat_MPIAIJ*)A->data;
141:   Mat          B = aij->B,Bnew;
142:   Mat_SeqAIJ   *Baij = (Mat_SeqAIJ*)B->data;
143:   int          ierr,i,j,m = B->m,n = A->N,col,ct = 0,*garray = aij->garray;
144:   int          *nz,ec,shift = Baij->indexshift;
145:   PetscScalar  v;

148:   /* free stuff related to matrix-vec multiply */
149:   VecGetSize(aij->lvec,&ec); /* needed for PetscLogObjectMemory below */
150:   VecDestroy(aij->lvec); aij->lvec = 0;
151:   VecScatterDestroy(aij->Mvctx); aij->Mvctx = 0;
152:   if (aij->colmap) {
153: #if defined (PETSC_USE_CTABLE)
154:     PetscTableDelete(aij->colmap);
155:     aij->colmap = 0;
156: #else
157:     PetscFree(aij->colmap);
158:     aij->colmap = 0;
159:     PetscLogObjectMemory(A,-aij->B->n*sizeof(int));
160: #endif
161:   }

163:   /* make sure that B is assembled so we can access its values */
164:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
165:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);

167:   /* invent new B and copy stuff over */
168:   PetscMalloc((m+1)*sizeof(int),&nz);
169:   for (i=0; i<m; i++) {
170:     nz[i] = Baij->i[i+1] - Baij->i[i];
171:   }
172:   MatCreateSeqAIJ(PETSC_COMM_SELF,m,n,0,nz,&Bnew);
173:   PetscFree(nz);
174:   for (i=0; i<m; i++) {
175:     for (j=Baij->i[i]+shift; j<Baij->i[i+1]+shift; j++) {
176:       col  = garray[Baij->j[ct]+shift];
177:       v    = Baij->a[ct++];
178:       MatSetValues(Bnew,1,&i,1,&col,&v,B->insertmode);
179:     }
180:   }
181:   PetscFree(aij->garray);
182:   aij->garray = 0;
183:   PetscLogObjectMemory(A,-ec*sizeof(int));
184:   MatDestroy(B);
185:   PetscLogObjectParent(A,Bnew);
186:   aij->B = Bnew;
187:   A->was_assembled = PETSC_FALSE;
188:   return(0);
189: }

191: /*      ugly stuff added for Glenn someday we should fix this up */

193: static int *auglyrmapd = 0,*auglyrmapo = 0;  /* mapping from the local ordering to the "diagonal" and "off-diagonal"
194:                                       parts of the local matrix */
195: static Vec auglydd = 0,auglyoo = 0;   /* work vectors used to scale the two parts of the local matrix */


198: #undef __FUNCT__
200: int MatMPIAIJDiagonalScaleLocalSetUp(Mat inA,Vec scale)
201: {
202:   Mat_MPIAIJ  *ina = (Mat_MPIAIJ*) inA->data; /*access private part of matrix */
203:   int          ierr,i,n,nt,cstart,cend,no,*garray = ina->garray,*lindices;
204:   int          *r_rmapd,*r_rmapo;
205: 
207:   MatGetOwnershipRange(inA,&cstart,&cend);
208:   MatGetSize(ina->A,PETSC_NULL,&n);
209:   PetscMalloc((inA->mapping->n+1)*sizeof(int),&r_rmapd);
210:   PetscMemzero(r_rmapd,inA->mapping->n*sizeof(int));
211:   nt   = 0;
212:   for (i=0; i<inA->mapping->n; i++) {
213:     if (inA->mapping->indices[i] >= cstart && inA->mapping->indices[i] < cend) {
214:       nt++;
215:       r_rmapd[i] = inA->mapping->indices[i] + 1;
216:     }
217:   }
218:   if (nt != n) SETERRQ2(1,"Hmm nt %d n %d",nt,n);
219:   PetscMalloc((n+1)*sizeof(int),&auglyrmapd);
220:   for (i=0; i<inA->mapping->n; i++) {
221:     if (r_rmapd[i]){
222:       auglyrmapd[(r_rmapd[i]-1)-cstart] = i;
223:     }
224:   }
225:   PetscFree(r_rmapd);
226:   VecCreateSeq(PETSC_COMM_SELF,n,&auglydd);

228:   PetscMalloc((inA->N+1)*sizeof(int),&lindices);
229:   PetscMemzero(lindices,inA->N*sizeof(int));
230:   for (i=0; i<ina->B->n; i++) {
231:     lindices[garray[i]] = i+1;
232:   }
233:   no   = inA->mapping->n - nt;
234:   PetscMalloc((inA->mapping->n+1)*sizeof(int),&r_rmapo);
235:   PetscMemzero(r_rmapo,inA->mapping->n*sizeof(int));
236:   nt   = 0;
237:   for (i=0; i<inA->mapping->n; i++) {
238:     if (lindices[inA->mapping->indices[i]]) {
239:       nt++;
240:       r_rmapo[i] = lindices[inA->mapping->indices[i]];
241:     }
242:   }
243:   if (nt > no) SETERRQ2(1,"Hmm nt %d no %d",nt,n);
244:   PetscFree(lindices);
245:   PetscMalloc((nt+1)*sizeof(int),&auglyrmapo);
246:   for (i=0; i<inA->mapping->n; i++) {
247:     if (r_rmapo[i]){
248:       auglyrmapo[(r_rmapo[i]-1)] = i;
249:     }
250:   }
251:   PetscFree(r_rmapo);
252:   VecCreateSeq(PETSC_COMM_SELF,nt,&auglyoo);

254:   return(0);
255: }

257: #undef __FUNCT__
259: int MatMPIAIJDiagonalScaleLocal(Mat A,Vec scale)
260: {
261:   Mat_MPIAIJ  *a = (Mat_MPIAIJ*) A->data; /*access private part of matrix */
262:   int         ierr,n,i;
263:   PetscScalar *d,*o,*s;
264: 
266:   if (!auglyrmapd) {
267:     MatMPIAIJDiagonalScaleLocalSetUp(A,scale);
268:   }

270:   VecGetArray(scale,&s);
271: 
272:   VecGetLocalSize(auglydd,&n);
273:   VecGetArray(auglydd,&d);
274:   for (i=0; i<n; i++) {
275:     d[i] = s[auglyrmapd[i]]; /* copy "diagonal" (true local) portion of scale into dd vector */
276:   }
277:   VecRestoreArray(auglydd,&d);
278:   /* column scale "diagonal" portion of local matrix */
279:   MatDiagonalScale(a->A,PETSC_NULL,auglydd);

281:   VecGetLocalSize(auglyoo,&n);
282:   VecGetArray(auglyoo,&o);
283:   for (i=0; i<n; i++) {
284:     o[i] = s[auglyrmapo[i]]; /* copy "off-diagonal" portion of scale into oo vector */
285:   }
286:   VecRestoreArray(scale,&s);
287:   VecRestoreArray(auglyoo,&o);
288:   /* column scale "off-diagonal" portion of local matrix */
289:   MatDiagonalScale(a->B,PETSC_NULL,auglyoo);

291:   return(0);
292: }