1 #ifndef lint 2 static char vcid[] = "$Id: baij2.c,v 1.3 1996/04/30 23:04:35 balay Exp balay $"; 3 #endif 4 5 #include "baij.h" 6 #include "petsc.h" 7 #include "src/inline/bitarray.h" 8 9 int MatIncreaseOverlap_SeqBAIJ(Mat A,int is_max,IS *is,int ov) 10 { 11 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 12 int row, i,j,k,l,m,n, *idx,ierr, *nidx, isz, val, ival; 13 int start, end, *ai, *aj,bs,*nidx2; 14 char *table; 15 16 m = a->mbs; 17 ai = a->i; 18 aj = a->j; 19 bs = a->bs; 20 21 if (ov < 0) SETERRQ(1,"MatIncreaseOverlap_SeqBAIJ: illegal overlap value used"); 22 23 table = (char *) PetscMalloc((m/BITSPERBYTE +1)*sizeof(char)); CHKPTRQ(table); 24 nidx = (int *) PetscMalloc((m+1)*sizeof(int)); CHKPTRQ(nidx); 25 nidx2 = (int *)PetscMalloc((a->m+1)*sizeof(int)); CHKPTRQ(nidx2); 26 27 for ( i=0; i<is_max; i++ ) { 28 /* Initialise the two local arrays */ 29 isz = 0; 30 PetscMemzero(table,(m/BITSPERBYTE +1)*sizeof(char)); 31 32 /* Extract the indices, assume there can be duplicate entries */ 33 ierr = ISGetIndices(is[i],&idx); CHKERRQ(ierr); 34 ierr = ISGetSize(is[i],&n); CHKERRQ(ierr); 35 36 /* Enter these into the temp arrays i.e mark table[row], enter row into new index */ 37 for ( j=0; j<n ; ++j){ 38 ival = idx[j]/bs; /* convert the indices into block indices */ 39 if (ival>m) SETERRQ(1,"MatIncreaseOverlap_SeqBAIJ: index greater than mat-dim"); 40 if(!BT_LOOKUP(table, ival)) { nidx[isz++] = ival;} 41 } 42 ierr = ISRestoreIndices(is[i],&idx); CHKERRQ(ierr); 43 ierr = ISDestroy(is[i]); CHKERRQ(ierr); 44 45 k = 0; 46 for ( j=0; j<ov; j++){ /* for each overlap*/ 47 n = isz; 48 for ( ; k<n ; k++){ /* do only those rows in nidx[k], which are not done yet */ 49 row = nidx[k]; 50 start = ai[row]; 51 end = ai[row+1]; 52 for ( l = start; l<end ; l++){ 53 val = aj[l]; 54 if (!BT_LOOKUP(table,val)) {nidx[isz++] = val;} 55 } 56 } 57 } 58 /* expand the Index Set */ 59 for (j=0; j<isz; j++ ) { 60 for (k=0; k<bs; k++ ) 61 nidx2[j*bs+k] = nidx[j]*bs+k; 62 } 63 ierr = ISCreateSeq(MPI_COMM_SELF, isz*bs, nidx2, (is+i)); CHKERRQ(ierr); 64 } 65 PetscFree(table); 66 PetscFree(nidx); 67 PetscFree(nidx2); 68 return 0; 69 } 70 int MatGetSubMatrix_SeqBAIJ_Private(Mat A,IS isrow,IS iscol,MatGetSubMatrixCall scall,Mat *B) 71 { 72 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data,*c; 73 int nznew, *smap, i, k, kstart, kend, ierr, oldcols = a->nbs,*lens; 74 int row,mat_i,*mat_j,tcol,*mat_ilen; 75 int *irow, *icol, nrows, ncols,*ssmap,bs=a->bs, bs2=a->bs2; 76 int *aj = a->j, *ai = a->i; 77 Scalar *mat_a; 78 Mat C; 79 80 ierr = ISSorted(iscol,(PetscTruth*)&i); 81 if (!i) SETERRQ(1,"MatGetSubmatrices_SeqBAIJ:IS is not sorted"); 82 83 ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr); 84 ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr); 85 ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr); 86 ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr); 87 88 smap = (int *) PetscMalloc((1+oldcols)*sizeof(int)); CHKPTRQ(smap); 89 ssmap = smap; 90 lens = (int *) PetscMalloc((1+nrows)*sizeof(int)); CHKPTRQ(lens); 91 PetscMemzero(smap,oldcols*sizeof(int)); 92 for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1; 93 /* determine lens of each row */ 94 for (i=0; i<nrows; i++) { 95 kstart = ai[irow[i]]; 96 kend = kstart + a->ilen[irow[i]]; 97 lens[i] = 0; 98 for ( k=kstart; k<kend; k++ ) { 99 if (ssmap[aj[k]]) { 100 lens[i]++; 101 } 102 } 103 } 104 /* Create and fill new matrix */ 105 if (scall == MAT_REUSE_MATRIX) { 106 c = (Mat_SeqBAIJ *)((*B)->data); 107 108 if (c->mbs!=nrows || c->nbs!=ncols || c->bs!=bs) 109 SETERRQ(1,"MatGetSubMatrix_SeqBAIJ:"); 110 if (PetscMemcmp(c->ilen,lens, c->mbs *sizeof(int))) { 111 SETERRQ(1,"MatGetSubmatrix_SeqBAIJ:Cannot reuse matrix. wrong no of nonzeros"); 112 } 113 PetscMemzero(c->ilen,c->mbs*sizeof(int)); 114 C = *B; 115 } 116 else { 117 ierr = MatCreateSeqBAIJ(A->comm,bs,nrows*bs,ncols*bs,0,lens,&C);CHKERRQ(ierr); 118 } 119 c = (Mat_SeqBAIJ *)(C->data); 120 for (i=0; i<nrows; i++) { 121 row = irow[i]; 122 nznew = 0; 123 kstart = ai[row]; 124 kend = kstart + a->ilen[row]; 125 mat_i = c->i[i]; 126 mat_j = c->j + mat_i; 127 mat_a = c->a + mat_i*bs2; 128 mat_ilen = c->ilen + i; 129 for ( k=kstart; k<kend; k++ ) { 130 if ((tcol=ssmap[a->j[k]])) { 131 *mat_j++ = tcol - 1; 132 PetscMemcpy(mat_a,a->a+k*bs2,bs2*sizeof(Scalar)); mat_a+=bs2; 133 (*mat_ilen)++; 134 135 } 136 } 137 } 138 139 /* Free work space */ 140 ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr); 141 PetscFree(smap); PetscFree(lens); 142 ierr = MatAssemblyBegin(C,FINAL_ASSEMBLY); CHKERRQ(ierr); 143 ierr = MatAssemblyEnd(C,FINAL_ASSEMBLY); CHKERRQ(ierr); 144 145 ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr); 146 *B = C; 147 return 0; 148 } 149 150 int MatGetSubMatrix_SeqBAIJ(Mat A,IS isrow,IS iscol,MatGetSubMatrixCall scall,Mat *B) 151 { 152 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 153 IS is1,is2; 154 int *vary,*iary,*irow,*icol,nrows,ncols,i,ierr,bs=a->bs,count; 155 156 ierr = ISGetIndices(isrow,&irow); CHKERRQ(ierr); 157 ierr = ISGetIndices(iscol,&icol); CHKERRQ(ierr); 158 ierr = ISGetSize(isrow,&nrows); CHKERRQ(ierr); 159 ierr = ISGetSize(iscol,&ncols); CHKERRQ(ierr); 160 161 /* Verify if the indices corespond to each elementin a block 162 and form the IS with compressed IS */ 163 vary = (int *) PetscMalloc(2*(a->mbs+1)*sizeof(int)); CHKPTRQ(vary); 164 iary = vary + a->mbs; 165 PetscMemzero(vary,(a->mbs)*sizeof(int)); 166 for ( i=0; i<nrows; i++) vary[irow[i]/bs]++; 167 count = 0; 168 for (i=0; i<a->mbs; i++) { 169 if (vary[i]!=0 && vary[i]!=bs) SETERRA(1,"MatGetSubmatrices_SeqBAIJ:"); 170 if (vary[i]==bs) iary[count++] = i; 171 } 172 ierr = ISCreateSeq(MPI_COMM_SELF, count, iary,&is1); CHKERRQ(ierr); 173 174 PetscMemzero(vary,(a->mbs)*sizeof(int)); 175 for ( i=0; i<ncols; i++) vary[icol[i]/bs]++; 176 count = 0; 177 for (i=0; i<a->mbs; i++) { 178 if (vary[i]!=0 && vary[i]!=bs) SETERRA(1,"MatGetSubmatrices_SeqBAIJ:"); 179 if (vary[i]==bs) iary[count++] = i; 180 } 181 ierr = ISCreateSeq(MPI_COMM_SELF, count, iary,&is2); CHKERRQ(ierr); 182 ierr = ISRestoreIndices(isrow,&irow); CHKERRQ(ierr); 183 ierr = ISRestoreIndices(iscol,&icol); CHKERRQ(ierr); 184 PetscFree(vary); 185 186 ierr = MatGetSubMatrix_SeqBAIJ_Private(A,is1,is2,scall,B); CHKERRQ(ierr); 187 ISDestroy(is1); 188 ISDestroy(is2); 189 return 0; 190 } 191 192 int MatGetSubMatrices_SeqBAIJ(Mat A,int n, IS *irow,IS *icol,MatGetSubMatrixCall scall, 193 Mat **B) 194 { 195 int ierr,i; 196 197 if (scall == MAT_INITIAL_MATRIX) { 198 *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) ); CHKPTRQ(*B); 199 } 200 201 for ( i=0; i<n; i++ ) { 202 ierr = MatGetSubMatrix(A,irow[i],icol[i],scall,&(*B)[i]); CHKERRQ(ierr); 203 } 204 return 0; 205 } 206 207 208 209 210 211 212 213