Actual source code: mmsbaij.c

  1: /*$Id: mmsbaij.c,v 1.8 2001/03/23 23:22:26 balay Exp $*/

  3: /*
  4:    Support for the parallel SBAIJ matrix vector multiply
  5: */
  6: #include "src/mat/impls/baij/mpi/mpibaij.h"
  7: #include "src/vec/vecimpl.h"
  8: extern int MatSetValues_SeqSBAIJ(Mat,int,int*,int,int*,Scalar*,InsertMode);

 10: int MatSetUpMultiply_MPISBAIJ(Mat mat)
 11: {
 12:   Mat_MPIBAIJ        *baij = (Mat_MPIBAIJ*)mat->data;
 13:   Mat_SeqBAIJ        *B = (Mat_SeqBAIJ*)(baij->B->data);
 14:   int                Nbs = baij->Nbs,i,j,*indices,*aj = B->j,ierr,ec = 0,*garray;
 15:   int                col,bs = baij->bs,*tmp,*stmp;
 16:   IS                 from,to;
 17:   Vec                gvec;
 18: #if defined (PETSC_USE_CTABLE)
 19:   PetscTable         gid1_lid1;
 20:   PetscTablePosition tpos;
 21:   int                gid,lid;
 22: #endif  


 26: #if defined (PETSC_USE_CTABLE)
 27:   /* use a table - Mark Adams */
 28:   PetscTableCreate(B->mbs,&gid1_lid1);
 29:   for (i=0; i<B->mbs; i++) {
 30:     for (j=0; j<B->ilen[i]; j++) {
 31:       int 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(int),&garray);
 41:   PetscMalloc((ec*bs+1)*sizeof(int),&tmp);
 42:   PetscTableGetHeadPosition(gid1_lid1,&tpos);
 43:   while (tpos) {
 44:     PetscTableGetNext(gid1_lid1,&tpos,&gid,&lid);
 45:     gid--; lid--;
 46:     garray[lid] = gid;
 47:   }
 48:   PetscSortInt(ec,garray);
 49:   /* qsort(garray, ec, sizeof(int), intcomparcarc); */
 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<B->mbs; i++) {
 56:     for (j=0; j<B->ilen[i]; j++) {
 57:       int gid1 = aj[B->i[i] + j] + 1;
 58:       PetscTableFind(gid1_lid1,gid1,&lid);
 59:       lid --;
 60:       aj[B->i[i]+j] = lid;
 61:     }
 62:   }
 63:   B->nbs     = ec;
 64:   baij->B->n = ec*B->bs;
 65:   PetscTableDelete(gid1_lid1);
 66:   /* Mark Adams */
 67: #else
 68:   /* For the first stab we make an array as long as the number of columns */
 69:   /* mark those columns that are in baij->B */
 70:   PetscMalloc((Nbs+1)*sizeof(int),&indices);
 71:   PetscMemzero(indices,Nbs*sizeof(int));
 72:   for (i=0; i<B->mbs; i++) {
 73:     for (j=0; j<B->ilen[i]; j++) {
 74:       if (!indices[aj[B->i[i] + j]]) ec++;
 75:       indices[aj[B->i[i] + j] ] = 1;
 76:     }
 77:   }

 79:   /* form array of columns we need */
 80:   PetscMalloc((ec+1)*sizeof(int),&garray);
 81:   PetscMalloc((ec*bs+1)*sizeof(int),&tmp);
 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*B->bs;
102:   PetscFree(indices);
103: #endif  

105:   for (i=0,col=0; i<ec; i++) {
106:     for (j=0; j<bs; j++,col++) tmp[col] = garray[i]*bs+j;
107:   }
108:   /* create local vector that is used to scatter into */
109:   VecCreateSeq(PETSC_COMM_SELF,ec*bs,&baij->lvec);

111:   /* create two temporary index sets for building scatter-gather */

113:   /* ISCreateGeneral(PETSC_COMM_SELF,ec*bs,tmp,&from); */
114:   for (i=0; i<ec; i++) {
115:     garray[i] = bs*garray[i];
116:   }
117:   ISCreateBlock(PETSC_COMM_SELF,bs,ec,garray,&from);
118:   for (i=0,col=0; i<ec; i++) {
119:     garray[i] = garray[i]/bs;
120:   }

122:   PetscMalloc((ec+1)*sizeof(int),&stmp);
123:   for (i=0; i<ec; i++) { stmp[i] = bs*i; }
124:   ISCreateBlock(PETSC_COMM_SELF,bs,ec,stmp,&to);
125:   PetscFree(stmp);

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

133:   /* gnerate the scatter context */
134:   VecScatterCreate(gvec,from,baij->lvec,to,&baij->Mvctx);

136:   /*
137:       Post the receives for the first matrix vector product. We sync-chronize after
138:     this on the chance that the user immediately calls MatMult() after assemblying 
139:     the matrix.
140:   */
141:   VecScatterPostRecvs(gvec,baij->lvec,INSERT_VALUES,SCATTER_FORWARD,baij->Mvctx);
142:   MPI_Barrier(mat->comm);

144:   PetscLogObjectParent(mat,baij->Mvctx);
145:   PetscLogObjectParent(mat,baij->lvec);
146:   PetscLogObjectParent(mat,from);
147:   PetscLogObjectParent(mat,to);
148:   baij->garray = garray;
149:   PetscLogObjectMemory(mat,(ec+1)*sizeof(int));
150:   ISDestroy(from);
151:   ISDestroy(to);
152:   VecDestroy(gvec);
153:   PetscFree(tmp);
154:   return(0);
155: }


158: /*
159:      Takes the local part of an already assembled MPIBAIJ matrix
160:    and disassembles it. This is to allow new nonzeros into the matrix
161:    that require more communication in the matrix vector multiply. 
162:    Thus certain data-structures must be rebuilt.

164:    Kind of slow! But that's what application programmers get when 
165:    they are sloppy.
166: */
167: int DisAssemble_MPISBAIJ(Mat A)
168: {
169:   Mat_MPIBAIJ *baij = (Mat_MPIBAIJ*)A->data;
170:   Mat         B = baij->B,Bnew;
171:   Mat_SeqBAIJ *Bbaij = (Mat_SeqBAIJ*)B->data;
172:   int         ierr,i,j,mbs=Bbaij->mbs,n = A->N,col,*garray=baij->garray;
173:   int         k,bs=baij->bs,bs2=baij->bs2,*rvals,*nz,ec,m=A->m;
174:   MatScalar   *a = Bbaij->a;
175:   Scalar      *atmp;
176: #if defined(PETSC_USE_MAT_SINGLE)
177:   int         l;
178: #endif

181: #if defined(PETSC_USE_MAT_SINGLE)
182:   PetscMalloc(baij->bs*sizeof(Scalar),&atmp);
183: #endif
184:   /* free stuff related to matrix-vec multiply */
185:   VecGetSize(baij->lvec,&ec); /* needed for PetscLogObjectMemory below */
186:   VecDestroy(baij->lvec); baij->lvec = 0;
187:   VecScatterDestroy(baij->Mvctx); baij->Mvctx = 0;
188:   if (baij->colmap) {
189: #if defined (PETSC_USE_CTABLE)
190:     PetscTableDelete(baij->colmap); baij->colmap = 0;
191: #else
192:     PetscFree(baij->colmap);
193:     baij->colmap = 0;
194:     PetscLogObjectMemory(A,-Bbaij->nbs*sizeof(int));
195: #endif
196:   }

198:   /* make sure that B is assembled so we can access its values */
199:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
200:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);

202:   /* invent new B and copy stuff over */
203:   PetscMalloc(mbs*sizeof(int),&nz);
204:   for (i=0; i<mbs; i++) {
205:     nz[i] = Bbaij->i[i+1]-Bbaij->i[i];
206:   }
207:   MatCreateSeqBAIJ(PETSC_COMM_SELF,baij->bs,m,n,0,nz,&Bnew);
208:   PetscFree(nz);
209: 
210:   PetscMalloc(bs*sizeof(int),&rvals);
211:   for (i=0; i<mbs; i++) {
212:     rvals[0] = bs*i;
213:     for (j=1; j<bs; j++) { rvals[j] = rvals[j-1] + 1; }
214:     for (j=Bbaij->i[i]; j<Bbaij->i[i+1]; j++) {
215:       col = garray[Bbaij->j[j]]*bs;
216:       for (k=0; k<bs; k++) {
217: #if defined(PETSC_USE_MAT_SINGLE)
218:         for (l=0; l<bs; l++) atmp[l] = a[j*bs2+l];
219: #else
220:         atmp = a+j*bs2;
221: #endif
222:         MatSetValues_SeqSBAIJ(Bnew,bs,rvals,1,&col,atmp,B->insertmode);
223:         col++;
224:       }
225:     }
226:   }
227: #if defined(PETSC_USE_MAT_SINGLE)
228:   PetscFree(atmp);
229: #endif
230:   PetscFree(baij->garray);
231:   baij->garray = 0;
232:   PetscFree(rvals);
233:   PetscLogObjectMemory(A,-ec*sizeof(int));
234:   MatDestroy(B);
235:   PetscLogObjectParent(A,Bnew);
236:   baij->B = Bnew;
237:   A->was_assembled = PETSC_FALSE;
238:   return(0);
239: }