Actual source code: mmaij.c

  1: /*
  2:    Support for the parallel AIJ matrix vector multiply
  3: */
 4:  #include src/mat/impls/aij/mpi/mpiaij.h

  8: PetscErrorCode MatSetUpMultiply_MPIAIJ(Mat mat)
  9: {
 10:   Mat_MPIAIJ         *aij = (Mat_MPIAIJ*)mat->data;
 11:   Mat_SeqAIJ         *B = (Mat_SeqAIJ*)(aij->B->data);
 12:   PetscErrorCode     ierr;
 13:   PetscInt           i,j,*aj = B->j,ec = 0,*garray;
 14:   IS                 from,to;
 15:   Vec                gvec;
 16: #if defined (PETSC_USE_CTABLE)
 17:   PetscTable         gid1_lid1;
 18:   PetscTablePosition tpos;
 19:   PetscInt           gid,lid;
 20: #else
 21:   PetscInt           N = mat->N,*indices;
 22: #endif


 26: #if defined (PETSC_USE_CTABLE)
 27:   /* use a table - Mark Adams */
 28:   PetscTableCreate(aij->B->m,&gid1_lid1);
 29:   for (i=0; i<aij->B->m; i++) {
 30:     for (j=0; j<B->ilen[i]; j++) {
 31:       PetscInt data,gid1 = aj[B->i[i] + j] + 1;
 32:       PetscTableFind(gid1_lid1,gid1,&data);
 33:       if (!data) {
 34:         /* one based table */
 35:         PetscTableAdd(gid1_lid1,gid1,++ec);
 36:       }
 37:     }
 38:   }
 39:   /* form array of columns we need */
 40:   PetscMalloc((ec+1)*sizeof(PetscInt),&garray);
 41:   PetscTableGetHeadPosition(gid1_lid1,&tpos);
 42:   while (tpos) {
 43:     PetscTableGetNext(gid1_lid1,&tpos,&gid,&lid);
 44:     gid--;
 45:     lid--;
 46:     garray[lid] = gid;
 47:   }
 48:   PetscSortInt(ec,garray); /* sort, and rebuild */
 49:   PetscTableRemoveAll(gid1_lid1);
 50:   for (i=0; i<ec; i++) {
 51:     PetscTableAdd(gid1_lid1,garray[i]+1,i+1);
 52:   }
 53:   /* compact out the extra columns in B */
 54:   for (i=0; i<aij->B->m; i++) {
 55:     for (j=0; j<B->ilen[i]; j++) {
 56:       PetscInt gid1 = aj[B->i[i] + j] + 1;
 57:       PetscTableFind(gid1_lid1,gid1,&lid);
 58:       lid --;
 59:       aj[B->i[i] + j]  = lid;
 60:     }
 61:   }
 62:   aij->B->n = aij->B->N = ec;
 63:   PetscTableDelete(gid1_lid1);
 64:   /* Mark Adams */
 65: #else
 66:   /* For the first stab we make an array as long as the number of columns */
 67:   /* mark those columns that are in aij->B */
 68:   PetscMalloc((N+1)*sizeof(PetscInt),&indices);
 69:   PetscMemzero(indices,N*sizeof(PetscInt));
 70:   for (i=0; i<aij->B->m; i++) {
 71:     for (j=0; j<B->ilen[i]; j++) {
 72:       if (!indices[aj[B->i[i] + j] ]) ec++;
 73:       indices[aj[B->i[i] + j] ] = 1;
 74:     }
 75:   }

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

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

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

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

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

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


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

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

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

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

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

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

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


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

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

256:   return(0);
257: }

261: PetscErrorCode MatMPIAIJDiagonalScaleLocal(Mat A,Vec scale)
262: {
263:   /* This routine should really be abandoned as it duplicates MatDiagonalScaleLocal */
264:   PetscErrorCode ierr,(*f)(Mat,Vec);

267:   PetscObjectQueryFunction((PetscObject)A,"MatDiagonalScaleLocal_C",(void (**)(void))&f);
268:   if (f) {
269:     (*f)(A,scale);
270:   }
271:   return(0);
272: }

277: PetscErrorCode MatDiagonalScaleLocal_MPIAIJ(Mat A,Vec scale)
278: {
279:   Mat_MPIAIJ     *a = (Mat_MPIAIJ*) A->data; /*access private part of matrix */
281:   PetscInt       n,i;
282:   PetscScalar    *d,*o,*s;
283: 
285:   if (!auglyrmapd) {
286:     MatMPIAIJDiagonalScaleLocalSetUp(A,scale);
287:   }

289:   VecGetArray(scale,&s);
290: 
291:   VecGetLocalSize(auglydd,&n);
292:   VecGetArray(auglydd,&d);
293:   for (i=0; i<n; i++) {
294:     d[i] = s[auglyrmapd[i]]; /* copy "diagonal" (true local) portion of scale into dd vector */
295:   }
296:   VecRestoreArray(auglydd,&d);
297:   /* column scale "diagonal" portion of local matrix */
298:   MatDiagonalScale(a->A,PETSC_NULL,auglydd);

300:   VecGetLocalSize(auglyoo,&n);
301:   VecGetArray(auglyoo,&o);
302:   for (i=0; i<n; i++) {
303:     o[i] = s[auglyrmapo[i]]; /* copy "off-diagonal" portion of scale into oo vector */
304:   }
305:   VecRestoreArray(scale,&s);
306:   VecRestoreArray(auglyoo,&o);
307:   /* column scale "off-diagonal" portion of local matrix */
308:   MatDiagonalScale(a->B,PETSC_NULL,auglyoo);

310:   return(0);
311: }