1 2 #include "src/mat/impls/aij/seq/aij.h" 3 #include "src/mat/impls/baij/seq/baij.h" 4 #include "src/mat/impls/sbaij/seq/sbaij.h" 5 6 EXTERN_C_BEGIN 7 #undef __FUNCT__ 8 #define __FUNCT__ "MatConvert_SeqSBAI_SeqAIJ" 9 PetscErrorCode MatConvert_SeqSBAIJ_SeqAIJ(Mat A,const MatType newtype,Mat *newmat) 10 { 11 Mat B; 12 Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 13 Mat_SeqAIJ *b; 14 PetscErrorCode ierr; 15 int *ai=a->i,*aj=a->j,m=A->m,n=A->n,i,j,k,*bi,*bj, 16 *rowlengths,nz,*rowstart,itmp; 17 int bs=a->bs,bs2=bs*bs,mbs=A->m/bs; 18 PetscScalar *av,*bv; 19 20 PetscFunctionBegin; 21 22 /* compute rowlengths of newmat */ 23 ierr = PetscMalloc((2*m+1)*sizeof(int),&rowlengths);CHKERRQ(ierr); 24 rowstart = rowlengths + m; 25 26 for (i=0; i<mbs; i++) rowlengths[i*bs] = 0; 27 aj = a->j; 28 k = 0; 29 for (i=0; i<mbs; i++) { 30 nz = ai[i+1] - ai[i]; 31 aj++; /* skip diagonal */ 32 for (j=1; j<nz; j++) { /* no. of lower triangular blocks */ 33 rowlengths[(*aj)*bs]++; aj++; 34 } 35 rowlengths[k] += nz; /* no. of upper triangular blocks */ 36 rowlengths[k] *= bs; 37 for (j=1; j<bs; j++) { 38 rowlengths[k+j] = rowlengths[k]; 39 } 40 k += bs; 41 /* printf(" rowlengths[%d]: %d\n",i, rowlengths[i]); */ 42 } 43 44 ierr = MatCreate(A->comm,m,n,m,n,&B);CHKERRQ(ierr); 45 ierr = MatSetType(B,newtype);CHKERRQ(ierr); 46 ierr = MatSeqAIJSetPreallocation(B,0,rowlengths);CHKERRQ(ierr); 47 ierr = MatSetOption(B,MAT_COLUMN_ORIENTED);CHKERRQ(ierr); 48 ierr = MatSetOption(B,MAT_ROWS_SORTED);CHKERRQ(ierr); 49 ierr = MatSetOption(B,MAT_COLUMNS_SORTED);CHKERRQ(ierr); 50 51 b = (Mat_SeqAIJ*)(B->data); 52 bi = b->i; 53 bj = b->j; 54 bv = b->a; 55 56 /* set b->i */ 57 bi[0] = 0; rowstart[0] = 0; 58 for (i=0; i<mbs; i++){ 59 for (j=0; j<bs; j++){ 60 b->ilen[i*bs+j] = rowlengths[i*bs]; 61 rowstart[i*bs+j+1] = rowstart[i*bs+j] + rowlengths[i*bs]; 62 } 63 bi[i+1] = bi[i] + rowlengths[i*bs]/bs; 64 } 65 if (bi[mbs] != 2*a->nz - mbs) SETERRQ2(PETSC_ERR_PLIB,"bi[mbs]: %D != 2*a->nz-mbs: %D\n",bi[mbs],2*a->nz - mbs); 66 67 /* set b->j and b->a */ 68 aj = a->j; av = a->a; 69 for (i=0; i<mbs; i++) { 70 /* diagonal block */ 71 for (j=0; j<bs; j++){ /* row i*bs+j */ 72 itmp = i*bs+j; 73 for (k=0; k<bs; k++){ /* col i*bs+k */ 74 *(bj + rowstart[itmp]) = (*aj)*bs+k; 75 *(bv + rowstart[itmp]) = *(av+k*bs+j); 76 rowstart[itmp]++; 77 } 78 } 79 aj++; av += bs2; 80 81 nz = ai[i+1] - ai[i] -1; 82 while (nz--){ 83 /* lower triangular blocks */ 84 for (j=0; j<bs; j++){ /* row (*aj)*bs+j */ 85 itmp = (*aj)*bs+j; 86 for (k=0; k<bs; k++){ /* col i*bs+k */ 87 *(bj + rowstart[itmp]) = i*bs+k; 88 *(bv + rowstart[itmp]) = *(av+k*bs+j); 89 rowstart[itmp]++; 90 } 91 } 92 /* upper triangular blocks */ 93 for (j=0; j<bs; j++){ /* row i*bs+j */ 94 itmp = i*bs+j; 95 for (k=0; k<bs; k++){ /* col (*aj)*bs+k */ 96 *(bj + rowstart[itmp]) = (*aj)*bs+k; 97 *(bv + rowstart[itmp]) = *(av+k*bs+j); 98 rowstart[itmp]++; 99 } 100 } 101 aj++; av += bs2; 102 } 103 } 104 ierr = PetscFree(rowlengths);CHKERRQ(ierr); 105 ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 106 ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 107 108 /* Fake support for "inplace" convert. */ 109 if (*newmat == A) { 110 ierr = MatDestroy(A);CHKERRQ(ierr); 111 } 112 *newmat = B; 113 PetscFunctionReturn(0); 114 } 115 #undef __FUNCT__ 116 #define __FUNCT__ "MatConvert_SeqAIJ_SeqSBAIJ" 117 PetscErrorCode MatConvert_SeqAIJ_SeqSBAIJ(Mat A,const MatType newtype,Mat *newmat) { 118 Mat B; 119 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 120 Mat_SeqSBAIJ *b; 121 PetscErrorCode ierr; 122 int *ai=a->i,*aj,m=A->M,n=A->N,i,j, 123 *bi,*bj,*rowlengths; 124 PetscScalar *av,*bv; 125 126 PetscFunctionBegin; 127 if (n != m) SETERRQ(PETSC_ERR_ARG_WRONG,"Matrix must be square"); 128 ierr = MatMissingDiagonal_SeqAIJ(A);CHKERRQ(ierr); /* check for missing diagonals, then mark diag */ 129 130 ierr = PetscMalloc(m*sizeof(int),&rowlengths);CHKERRQ(ierr); 131 for (i=0; i<m; i++) { 132 rowlengths[i] = ai[i+1] - a->diag[i]; 133 } 134 ierr = MatCreate(A->comm,m,n,m,n,&B);CHKERRQ(ierr); 135 ierr = MatSetType(B,newtype);CHKERRQ(ierr); 136 ierr = MatSeqSBAIJSetPreallocation(B,1,0,rowlengths);CHKERRQ(ierr); 137 138 ierr = MatSetOption(B,MAT_ROW_ORIENTED);CHKERRQ(ierr); 139 ierr = MatSetOption(B,MAT_ROWS_SORTED);CHKERRQ(ierr); 140 ierr = MatSetOption(B,MAT_COLUMNS_SORTED);CHKERRQ(ierr); 141 142 b = (Mat_SeqSBAIJ*)(B->data); 143 bi = b->i; 144 bj = b->j; 145 bv = b->a; 146 147 bi[0] = 0; 148 for (i=0; i<m; i++) { 149 aj = a->j + a->diag[i]; 150 av = a->a + a->diag[i]; 151 for (j=0; j<rowlengths[i]; j++){ 152 *bj = *aj; bj++; aj++; 153 *bv = *av; bv++; av++; 154 } 155 bi[i+1] = bi[i] + rowlengths[i]; 156 b->ilen[i] = rowlengths[i]; 157 } 158 159 ierr = PetscFree(rowlengths);CHKERRQ(ierr); 160 ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 161 ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 162 163 /* Fake support for "inplace" convert. */ 164 if (*newmat == A) { 165 ierr = MatDestroy(A);CHKERRQ(ierr); 166 } 167 *newmat = B; 168 169 PetscFunctionReturn(0); 170 } 171 EXTERN_C_END 172 173 EXTERN_C_BEGIN 174 #undef __FUNCT__ 175 #define __FUNCT__ "MatConvert_SeqSBAI_SeqBAIJ" 176 PetscErrorCode MatConvert_SeqSBAIJ_SeqBAIJ(Mat A,const MatType newtype,Mat *newmat) 177 { 178 Mat B; 179 Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 180 Mat_SeqBAIJ *b; 181 PetscErrorCode ierr; 182 int *ai=a->i,*aj=a->j,m=A->m,n=A->n,i,k,*bi,*bj, 183 *browlengths,nz,*browstart,itmp; 184 int bs=a->bs,bs2=bs*bs,mbs=m/bs; 185 PetscScalar *av,*bv; 186 187 PetscFunctionBegin; 188 /* compute browlengths of newmat */ 189 ierr = PetscMalloc(2*mbs*sizeof(int),&browlengths);CHKERRQ(ierr); 190 browstart = browlengths + mbs; 191 for (i=0; i<mbs; i++) browlengths[i] = 0; 192 aj = a->j; 193 for (i=0; i<mbs; i++) { 194 nz = ai[i+1] - ai[i]; 195 aj++; /* skip diagonal */ 196 for (k=1; k<nz; k++) { /* no. of lower triangular blocks */ 197 browlengths[*aj]++; aj++; 198 } 199 browlengths[i] += nz; /* no. of upper triangular blocks */ 200 } 201 202 ierr = MatCreate(A->comm,m,n,m,n,&B);CHKERRQ(ierr); 203 ierr = MatSetType(B,newtype);CHKERRQ(ierr); 204 ierr = MatSeqBAIJSetPreallocation(B,bs,0,browlengths);CHKERRQ(ierr); 205 ierr = MatSetOption(B,MAT_ROW_ORIENTED);CHKERRQ(ierr); 206 ierr = MatSetOption(B,MAT_ROWS_SORTED);CHKERRQ(ierr); 207 ierr = MatSetOption(B,MAT_COLUMNS_SORTED);CHKERRQ(ierr); 208 209 b = (Mat_SeqBAIJ*)(B->data); 210 bi = b->i; 211 bj = b->j; 212 bv = b->a; 213 214 /* set b->i */ 215 bi[0] = 0; 216 for (i=0; i<mbs; i++){ 217 b->ilen[i] = browlengths[i]; 218 bi[i+1] = bi[i] + browlengths[i]; 219 browstart[i] = bi[i]; 220 } 221 if (bi[mbs] != 2*a->nz - mbs) SETERRQ2(PETSC_ERR_PLIB,"bi[mbs]: %D != 2*a->nz - mbs: %D\n",bi[mbs],2*a->nz - mbs); 222 223 /* set b->j and b->a */ 224 aj = a->j; av = a->a; 225 for (i=0; i<mbs; i++) { 226 /* diagonal block */ 227 *(bj + browstart[i]) = *aj; aj++; 228 itmp = bs2*browstart[i]; 229 for (k=0; k<bs2; k++){ 230 *(bv + itmp + k) = *av; av++; 231 } 232 browstart[i]++; 233 234 nz = ai[i+1] - ai[i] -1; 235 while (nz--){ 236 /* lower triangular blocks */ 237 *(bj + browstart[*aj]) = i; 238 itmp = bs2*browstart[*aj]; 239 for (k=0; k<bs2; k++){ 240 *(bv + itmp + k) = *(av + k); 241 } 242 browstart[*aj]++; 243 244 /* upper triangular blocks */ 245 *(bj + browstart[i]) = *aj; aj++; 246 itmp = bs2*browstart[i]; 247 for (k=0; k<bs2; k++){ 248 *(bv + itmp + k) = *av; av++; 249 } 250 browstart[i]++; 251 } 252 } 253 ierr = PetscFree(browlengths);CHKERRQ(ierr); 254 ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 255 ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 256 257 /* Fake support for "inplace" convert. */ 258 if (*newmat == A) { 259 ierr = MatDestroy(A);CHKERRQ(ierr); 260 } 261 *newmat = B; 262 PetscFunctionReturn(0); 263 } 264 #undef __FUNCT__ 265 #define __FUNCT__ "MatConvert_SeqBAIJ_SeqSBAIJ" 266 PetscErrorCode MatConvert_SeqBAIJ_SeqSBAIJ(Mat A,const MatType newtype,Mat *newmat) 267 { 268 Mat B; 269 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 270 Mat_SeqSBAIJ *b; 271 PetscErrorCode ierr; 272 int *ai=a->i,*aj,m=A->m,n=A->n,i,j,k, 273 *bi,*bj,*browlengths; 274 int bs=a->bs,bs2=bs*bs,mbs=m/bs; 275 PetscScalar *av,*bv; 276 277 PetscFunctionBegin; 278 if (n != m) SETERRQ(PETSC_ERR_ARG_WRONG,"Matrix must be square"); 279 ierr = MatMissingDiagonal_SeqBAIJ(A);CHKERRQ(ierr); /* check for missing diagonals, then mark diag */ 280 281 ierr = PetscMalloc(mbs*sizeof(int),&browlengths);CHKERRQ(ierr); 282 for (i=0; i<mbs; i++) { 283 browlengths[i] = ai[i+1] - a->diag[i]; 284 } 285 286 ierr = MatCreate(A->comm,m,n,m,n,&B);CHKERRQ(ierr); 287 ierr = MatSetType(B,newtype);CHKERRQ(ierr); 288 ierr = MatSeqSBAIJSetPreallocation(B,bs,0,browlengths);CHKERRQ(ierr); 289 ierr = MatSetOption(B,MAT_ROW_ORIENTED);CHKERRQ(ierr); 290 ierr = MatSetOption(B,MAT_ROWS_SORTED);CHKERRQ(ierr); 291 ierr = MatSetOption(B,MAT_COLUMNS_SORTED);CHKERRQ(ierr); 292 293 b = (Mat_SeqSBAIJ*)(B->data); 294 bi = b->i; 295 bj = b->j; 296 bv = b->a; 297 298 bi[0] = 0; 299 for (i=0; i<mbs; i++) { 300 aj = a->j + a->diag[i]; 301 av = a->a + (a->diag[i])*bs2; 302 for (j=0; j<browlengths[i]; j++){ 303 *bj = *aj; bj++; aj++; 304 for (k=0; k<bs2; k++){ 305 *bv = *av; bv++; av++; 306 } 307 } 308 bi[i+1] = bi[i] + browlengths[i]; 309 b->ilen[i] = browlengths[i]; 310 } 311 ierr = PetscFree(browlengths);CHKERRQ(ierr); 312 ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 313 ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 314 315 /* Fake support for "inplace" convert. */ 316 if (*newmat == A) { 317 ierr = MatDestroy(A);CHKERRQ(ierr); 318 } 319 *newmat = B; 320 321 PetscFunctionReturn(0); 322 } 323 EXTERN_C_END 324