Actual source code: mmbaij.c

  1: /*$Id: mmbaij.c,v 1.46 2001/09/25 00:31:36 balay Exp $*/

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

  9: EXTERN int MatSetValuesBlocked_SeqBAIJ(Mat,int,int*,int,int*,PetscScalar*,InsertMode);

 11: int MatSetUpMultiply_MPIBAIJ(Mat mat)
 12: {
 13:   Mat_MPIBAIJ        *baij = (Mat_MPIBAIJ*)mat->data;
 14:   Mat_SeqBAIJ        *B = (Mat_SeqBAIJ*)(baij->B->data);
 15:   int                Nbs = baij->Nbs,i,j,*indices,*aj = B->j,ierr,ec = 0,*garray;
 16:   int                bs = baij->bs,*stmp;
 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 */
 29:   PetscTableCreate(B->mbs,&gid1_lid1);
 30:   for (i=0; i<B->mbs; i++) {
 31:     for (j=0; j<B->ilen[i]; j++) {
 32:       int data,gid1 = aj[B->i[i]+j] + 1;
 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--; 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:   /* 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:   ec = 0;
 82:   for (i=0; i<Nbs; i++) {
 83:     if (indices[i]) {
 84:       garray[ec++] = i;
 85:     }
 86:   }

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

 93:   /* compact out the extra columns in B */
 94:   for (i=0; i<B->mbs; i++) {
 95:     for (j=0; j<B->ilen[i]; j++) {
 96:       aj[B->i[i] + j] = indices[aj[B->i[i] + j]];
 97:     }
 98:   }
 99:   B->nbs       = ec;
100:   baij->B->n   = ec*B->bs;
101:   PetscFree(indices);
102: #endif  

104:   /* create local vector that is used to scatter into */
105:   VecCreateSeq(PETSC_COMM_SELF,ec*bs,&baij->lvec);

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

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

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

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

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

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

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_MPIBAIJ(Mat A)
160: {
161:   Mat_MPIBAIJ  *baij = (Mat_MPIBAIJ*)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          bs2=baij->bs2,*nz,ec,m = A->m;
166:   MatScalar    *a = Bbaij->a;
167:   PetscScalar  *atmp;

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

184:   /* make sure that B is assembled so we can access its values */
185:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
186:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);

188:   /* invent new B and copy stuff over */
189:   PetscMalloc(mbs*sizeof(int),&nz);
190:   for (i=0; i<mbs; i++) {
191:     nz[i] = Bbaij->i[i+1]-Bbaij->i[i];
192:   }
193:   MatCreateSeqBAIJ(PETSC_COMM_SELF,baij->bs,m,n,0,nz,&Bnew);
194:   MatSetOption(Bnew,MAT_COLUMN_ORIENTED);

196: #if defined(PETSC_USE_MAT_SINGLE)
197:   PetscMalloc(bs2*sizeof(PetscScalar),&atmp);
198: #endif
199:     for (i=0; i<mbs; i++) {
200:       for (j=Bbaij->i[i]; j<Bbaij->i[i+1]; j++) {
201:         col  = garray[Bbaij->j[j]];
202: #if defined(PETSC_USE_MAT_SINGLE)
203:         for (k=0; k<bs2; k++) atmp[k] = a[j*bs2+k];
204: #else
205:         atmp = a + j*bs2;
206: #endif
207:         MatSetValuesBlocked_SeqBAIJ(Bnew,1,&i,1,&col,atmp,B->insertmode);
208:       }
209:     }
210:   MatSetOption(Bnew,MAT_ROW_ORIENTED);

212: #if defined(PETSC_USE_MAT_SINGLE)
213:   PetscFree(atmp);
214: #endif

216:   PetscFree(nz);
217:   PetscFree(baij->garray);
218:   baij->garray = 0;
219:   PetscLogObjectMemory(A,-ec*sizeof(int));
220:   MatDestroy(B);
221:   PetscLogObjectParent(A,Bnew);
222:   baij->B = Bnew;
223:   A->was_assembled = PETSC_FALSE;
224:   return(0);
225: }