Actual source code: mmsbaij.c

  1: /*$Id: mmsbaij.c,v 1.10 2001/08/07 03:03:05 balay Exp $*/

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

 10: int MatSetUpMultiply_MPISBAIJ(Mat mat)
 11: {
 12:   Mat_MPISBAIJ       *baij = (Mat_MPISBAIJ*)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                bs = baij->bs,*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:   PetscTableGetHeadPosition(gid1_lid1,&tpos);
 42:   while (tpos) {
 43:     PetscTableGetNext(gid1_lid1,&tpos,&gid,&lid);
 44:     gid--; lid--;
 45:     garray[lid] = gid;
 46:   }
 47:   PetscSortInt(ec,garray);
 48:   /* qsort(garray, ec, sizeof(int), intcomparcarc); */
 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<B->mbs; i++) {
 55:     for (j=0; j<B->ilen[i]; j++) {
 56:       int 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:   B->nbs     = ec;
 63:   baij->B->n = ec*B->bs;
 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 baij->B */
 69:   PetscMalloc((Nbs+1)*sizeof(int),&indices);
 70:   PetscMemzero(indices,Nbs*sizeof(int));
 71:   for (i=0; i<B->mbs; i++) {
 72:     for (j=0; j<B->ilen[i]; j++) {
 73:       if (!indices[aj[B->i[i] + j]]) ec++;
 74:       indices[aj[B->i[i] + j] ] = 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<Nbs; i++) {
 82:     if (indices[i]) {
 83:       garray[ec++] = i;
 84:     }
 85:   }

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

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

106:   /* create two temporary index sets for building scatter-gather */
107:   for (i=0; i<ec; i++) {
108:     garray[i] = bs*garray[i];
109:   }
110:   ISCreateBlock(PETSC_COMM_SELF,bs,ec,garray,&from);
111:   for (i=0; i<ec; i++) {
112:     garray[i] = garray[i]/bs;
113:   }

115:   PetscMalloc((ec+1)*sizeof(int),&stmp);
116:   for (i=0; i<ec; i++) { stmp[i] = bs*i; }
117:   ISCreateBlock(PETSC_COMM_SELF,bs,ec,stmp,&to);
118:   PetscFree(stmp);

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

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

129:   /*
130:       Post the receives for the first matrix vector product. We sync-chronize after
131:     this on the chance that the user immediately calls MatMult() after assemblying 
132:     the matrix.
133:   */
134:   VecScatterPostRecvs(gvec,baij->lvec,INSERT_VALUES,SCATTER_FORWARD,baij->Mvctx);
135:   MPI_Barrier(mat->comm);

137:   PetscLogObjectParent(mat,baij->Mvctx);
138:   PetscLogObjectParent(mat,baij->lvec);
139:   PetscLogObjectParent(mat,from);
140:   PetscLogObjectParent(mat,to);
141:   baij->garray = garray;
142:   PetscLogObjectMemory(mat,(ec+1)*sizeof(int));
143:   ISDestroy(from);
144:   ISDestroy(to);
145:   VecDestroy(gvec);
146:   return(0);
147: }


150: /*
151:      Takes the local part of an already assembled MPIBAIJ matrix
152:    and disassembles it. This is to allow new nonzeros into the matrix
153:    that require more communication in the matrix vector multiply. 
154:    Thus certain data-structures must be rebuilt.

156:    Kind of slow! But that's what application programmers get when 
157:    they are sloppy.
158: */
159: int DisAssemble_MPISBAIJ(Mat A)
160: {
161:   Mat_MPISBAIJ   *baij = (Mat_MPISBAIJ*)A->data;
162:   Mat           B = baij->B,Bnew;
163:   Mat_SeqBAIJ   *Bbaij = (Mat_SeqBAIJ*)B->data;
164:   int           ierr,i,j,mbs=Bbaij->mbs,n = A->N,col,*garray=baij->garray;
165:   int           k,bs=baij->bs,bs2=baij->bs2,*rvals,*nz,ec,m=A->m;
166:   MatScalar     *a = Bbaij->a;
167:   PetscScalar   *atmp;
168: #if defined(PETSC_USE_MAT_SINGLE)
169:   int           l;
170: #endif

173: #if defined(PETSC_USE_MAT_SINGLE)
174:   PetscMalloc(baij->bs*sizeof(PetscScalar),&atmp);
175: #endif
176:   /* free stuff related to matrix-vec multiply */
177:   VecGetSize(baij->lvec,&ec); /* needed for PetscLogObjectMemory below */
178:   VecDestroy(baij->lvec); baij->lvec = 0;
179:   VecScatterDestroy(baij->Mvctx); baij->Mvctx = 0;
180:   if (baij->colmap) {
181: #if defined (PETSC_USE_CTABLE)
182:     PetscTableDelete(baij->colmap); baij->colmap = 0;
183: #else
184:     PetscFree(baij->colmap);
185:     baij->colmap = 0;
186:     PetscLogObjectMemory(A,-Bbaij->nbs*sizeof(int));
187: #endif
188:   }

190:   /* make sure that B is assembled so we can access its values */
191:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
192:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);

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