1 #ifdef PETSC_RCS_HEADER 2 static char vcid[] = "$Id: baij2.c,v 1.48 1999/06/30 23:51:46 balay Exp bsmith $"; 3 #endif 4 5 #include "sys.h" 6 #include "src/mat/impls/baij/seq/baij.h" 7 #include "src/vec/vecimpl.h" 8 #include "src/inline/spops.h" 9 #include "src/inline/ilu.h" 10 #include "bitarray.h" 11 12 #undef __FUNC__ 13 #define __FUNC__ "MatIncreaseOverlap_SeqBAIJ" 14 int MatIncreaseOverlap_SeqBAIJ(Mat A,int is_max,IS *is,int ov) 15 { 16 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 17 int row, i,j,k,l,m,n, *idx,ierr, *nidx, isz, val, ival; 18 int start, end, *ai, *aj,bs,*nidx2; 19 BT table; 20 21 PetscFunctionBegin; 22 m = a->mbs; 23 ai = a->i; 24 aj = a->j; 25 bs = a->bs; 26 27 if (ov < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Negative overlap specified"); 28 29 ierr = BTCreate(m,table);CHKERRQ(ierr); 30 nidx = (int *) PetscMalloc((m+1)*sizeof(int));CHKPTRQ(nidx); 31 nidx2 = (int *)PetscMalloc((a->m+1)*sizeof(int));CHKPTRQ(nidx2); 32 33 for ( i=0; i<is_max; i++ ) { 34 /* Initialise the two local arrays */ 35 isz = 0; 36 ierr = BTMemzero(m,table);CHKERRQ(ierr); 37 38 /* Extract the indices, assume there can be duplicate entries */ 39 ierr = ISGetIndices(is[i],&idx);CHKERRQ(ierr); 40 ierr = ISGetSize(is[i],&n);CHKERRQ(ierr); 41 42 /* Enter these into the temp arrays i.e mark table[row], enter row into new index */ 43 for ( j=0; j<n ; ++j){ 44 ival = idx[j]/bs; /* convert the indices into block indices */ 45 if (ival>m) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"index greater than mat-dim"); 46 if(!BTLookupSet(table, ival)) { nidx[isz++] = ival;} 47 } 48 ierr = ISRestoreIndices(is[i],&idx);CHKERRQ(ierr); 49 ierr = ISDestroy(is[i]);CHKERRQ(ierr); 50 51 k = 0; 52 for ( j=0; j<ov; j++){ /* for each overlap*/ 53 n = isz; 54 for ( ; k<n ; k++){ /* do only those rows in nidx[k], which are not done yet */ 55 row = nidx[k]; 56 start = ai[row]; 57 end = ai[row+1]; 58 for ( l = start; l<end ; l++){ 59 val = aj[l]; 60 if (!BTLookupSet(table,val)) {nidx[isz++] = val;} 61 } 62 } 63 } 64 /* expand the Index Set */ 65 for (j=0; j<isz; j++ ) { 66 for (k=0; k<bs; k++ ) 67 nidx2[j*bs+k] = nidx[j]*bs+k; 68 } 69 ierr = ISCreateGeneral(PETSC_COMM_SELF, isz*bs, nidx2, (is+i));CHKERRQ(ierr); 70 } 71 ierr = BTDestroy(table);CHKERRQ(ierr); 72 ierr = PetscFree(nidx);CHKERRQ(ierr); 73 ierr = PetscFree(nidx2);CHKERRQ(ierr); 74 PetscFunctionReturn(0); 75 } 76 77 #undef __FUNC__ 78 #define __FUNC__ "MatGetSubMatrix_SeqBAIJ_Private" 79 int MatGetSubMatrix_SeqBAIJ_Private(Mat A,IS isrow,IS iscol,int cs,MatReuse scall,Mat *B) 80 { 81 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data,*c; 82 int *smap, i, k, kstart, kend, ierr, oldcols = a->nbs,*lens; 83 int row,mat_i,*mat_j,tcol,*mat_ilen; 84 int *irow, *icol, nrows, ncols,*ssmap,bs=a->bs, bs2=a->bs2; 85 int *aj = a->j, *ai = a->i; 86 MatScalar *mat_a; 87 Mat C; 88 89 PetscFunctionBegin; 90 ierr = ISSorted(iscol,(PetscTruth*)&i);CHKERRQ(ierr); 91 if (!i) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"IS is not sorted"); 92 93 ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr); 94 ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr); 95 ierr = ISGetSize(isrow,&nrows);CHKERRQ(ierr); 96 ierr = ISGetSize(iscol,&ncols);CHKERRQ(ierr); 97 98 smap = (int *) PetscMalloc((1+oldcols)*sizeof(int));CHKPTRQ(smap); 99 ssmap = smap; 100 lens = (int *) PetscMalloc((1+nrows)*sizeof(int));CHKPTRQ(lens); 101 ierr = PetscMemzero(smap,oldcols*sizeof(int));CHKERRQ(ierr); 102 for ( i=0; i<ncols; i++ ) smap[icol[i]] = i+1; 103 /* determine lens of each row */ 104 for (i=0; i<nrows; i++) { 105 kstart = ai[irow[i]]; 106 kend = kstart + a->ilen[irow[i]]; 107 lens[i] = 0; 108 for ( k=kstart; k<kend; k++ ) { 109 if (ssmap[aj[k]]) { 110 lens[i]++; 111 } 112 } 113 } 114 /* Create and fill new matrix */ 115 if (scall == MAT_REUSE_MATRIX) { 116 c = (Mat_SeqBAIJ *)((*B)->data); 117 118 if (c->mbs!=nrows || c->nbs!=ncols || c->bs!=bs) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Submatrix wrong size"); 119 if (PetscMemcmp(c->ilen,lens, c->mbs *sizeof(int))) { 120 SETERRQ(PETSC_ERR_ARG_SIZ,0,"Cannot reuse matrix. wrong no of nonzeros"); 121 } 122 ierr = PetscMemzero(c->ilen,c->mbs*sizeof(int));CHKERRQ(ierr); 123 C = *B; 124 } else { 125 ierr = MatCreateSeqBAIJ(A->comm,bs,nrows*bs,ncols*bs,0,lens,&C);CHKERRQ(ierr); 126 } 127 c = (Mat_SeqBAIJ *)(C->data); 128 for (i=0; i<nrows; i++) { 129 row = irow[i]; 130 kstart = ai[row]; 131 kend = kstart + a->ilen[row]; 132 mat_i = c->i[i]; 133 mat_j = c->j + mat_i; 134 mat_a = c->a + mat_i*bs2; 135 mat_ilen = c->ilen + i; 136 for ( k=kstart; k<kend; k++ ) { 137 if ((tcol=ssmap[a->j[k]])) { 138 *mat_j++ = tcol - 1; 139 ierr = PetscMemcpy(mat_a,a->a+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 140 mat_a += bs2; 141 (*mat_ilen)++; 142 } 143 } 144 } 145 146 /* Free work space */ 147 ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr); 148 ierr = PetscFree(smap);CHKERRQ(ierr); 149 ierr = PetscFree(lens);CHKERRQ(ierr); 150 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 151 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 152 153 ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr); 154 *B = C; 155 PetscFunctionReturn(0); 156 } 157 158 #undef __FUNC__ 159 #define __FUNC__ "MatGetSubMatrix_SeqBAIJ" 160 int MatGetSubMatrix_SeqBAIJ(Mat A,IS isrow,IS iscol,int cs,MatReuse scall,Mat *B) 161 { 162 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 163 IS is1,is2; 164 int *vary,*iary,*irow,*icol,nrows,ncols,i,ierr,bs=a->bs,count; 165 166 PetscFunctionBegin; 167 ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr); 168 ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr); 169 ierr = ISGetSize(isrow,&nrows);CHKERRQ(ierr); 170 ierr = ISGetSize(iscol,&ncols);CHKERRQ(ierr); 171 172 /* Verify if the indices corespond to each element in a block 173 and form the IS with compressed IS */ 174 vary = (int *) PetscMalloc(2*(a->mbs+1)*sizeof(int));CHKPTRQ(vary); 175 iary = vary + a->mbs; 176 ierr = PetscMemzero(vary,(a->mbs)*sizeof(int));CHKERRQ(ierr); 177 for ( i=0; i<nrows; i++) vary[irow[i]/bs]++; 178 count = 0; 179 for (i=0; i<a->mbs; i++) { 180 if (vary[i]!=0 && vary[i]!=bs) SETERRA(1,0,"Index set does not match blocks"); 181 if (vary[i]==bs) iary[count++] = i; 182 } 183 ierr = ISCreateGeneral(PETSC_COMM_SELF, count, iary,&is1);CHKERRQ(ierr); 184 185 ierr = PetscMemzero(vary,(a->mbs)*sizeof(int));CHKERRQ(ierr); 186 for ( i=0; i<ncols; i++) vary[icol[i]/bs]++; 187 count = 0; 188 for (i=0; i<a->mbs; i++) { 189 if (vary[i]!=0 && vary[i]!=bs) SETERRA(1,0,"MatGetSubmatrices_SeqBAIJ:"); 190 if (vary[i]==bs) iary[count++] = i; 191 } 192 ierr = ISCreateGeneral(PETSC_COMM_SELF, count, iary,&is2);CHKERRQ(ierr); 193 ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr); 194 ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr); 195 ierr = PetscFree(vary);CHKERRQ(ierr); 196 197 ierr = MatGetSubMatrix_SeqBAIJ_Private(A,is1,is2,cs,scall,B);CHKERRQ(ierr); 198 ISDestroy(is1); 199 ISDestroy(is2); 200 PetscFunctionReturn(0); 201 } 202 203 #undef __FUNC__ 204 #define __FUNC__ "MatGetSubMatrices_SeqBAIJ" 205 int MatGetSubMatrices_SeqBAIJ(Mat A,int n, IS *irow,IS *icol,MatReuse scall,Mat **B) 206 { 207 int ierr,i; 208 209 PetscFunctionBegin; 210 if (scall == MAT_INITIAL_MATRIX) { 211 *B = (Mat *) PetscMalloc( (n+1)*sizeof(Mat) );CHKPTRQ(*B); 212 } 213 214 for ( i=0; i<n; i++ ) { 215 ierr = MatGetSubMatrix_SeqBAIJ(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr); 216 } 217 PetscFunctionReturn(0); 218 } 219 220 221 /* -------------------------------------------------------*/ 222 /* Should check that shapes of vectors and matrices match */ 223 /* -------------------------------------------------------*/ 224 #include "pinclude/blaslapack.h" 225 226 #undef __FUNC__ 227 #define __FUNC__ "MatMult_SeqBAIJ_1" 228 int MatMult_SeqBAIJ_1(Mat A,Vec xx,Vec zz) 229 { 230 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 231 Scalar *x,*z,sum; 232 MatScalar *v; 233 int mbs=a->mbs,i,*idx,*ii,n,ierr; 234 235 PetscFunctionBegin; 236 ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 237 ierr = VecGetArray(zz,&z);CHKERRQ(ierr); 238 239 idx = a->j; 240 v = a->a; 241 ii = a->i; 242 243 for ( i=0; i<mbs; i++ ) { 244 n = ii[1] - ii[0]; ii++; 245 sum = 0.0; 246 while (n--) sum += *v++ * x[*idx++]; 247 z[i] = sum; 248 } 249 ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 250 ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr); 251 PLogFlops(2*a->nz - a->m); 252 PetscFunctionReturn(0); 253 } 254 255 #undef __FUNC__ 256 #define __FUNC__ "MatMult_SeqBAIJ_2" 257 int MatMult_SeqBAIJ_2(Mat A,Vec xx,Vec zz) 258 { 259 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 260 Scalar *x,*z,*xb,sum1,sum2; 261 Scalar x1,x2; 262 MatScalar *v; 263 int ierr,mbs=a->mbs,i,*idx,*ii,j,n; 264 265 PetscFunctionBegin; 266 ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 267 ierr = VecGetArray(zz,&z);CHKERRQ(ierr); 268 269 idx = a->j; 270 v = a->a; 271 ii = a->i; 272 273 for ( i=0; i<mbs; i++ ) { 274 n = ii[1] - ii[0]; ii++; 275 sum1 = 0.0; sum2 = 0.0; 276 for ( j=0; j<n; j++ ) { 277 xb = x + 2*(*idx++); x1 = xb[0]; x2 = xb[1]; 278 sum1 += v[0]*x1 + v[2]*x2; 279 sum2 += v[1]*x1 + v[3]*x2; 280 v += 4; 281 } 282 z[0] = sum1; z[1] = sum2; 283 z += 2; 284 } 285 ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 286 ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr); 287 PLogFlops(8*a->nz - a->m); 288 PetscFunctionReturn(0); 289 } 290 291 #undef __FUNC__ 292 #define __FUNC__ "MatMult_SeqBAIJ_3" 293 int MatMult_SeqBAIJ_3(Mat A,Vec xx,Vec zz) 294 { 295 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 296 Scalar *x,*z,*xb,sum1,sum2,sum3,x1,x2,x3; 297 MatScalar *v; 298 int ierr,mbs=a->mbs,i,*idx,*ii,j,n; 299 300 #if defined(PETSC_HAVE_PRAGMA_DISJOINT) 301 #pragma disjoint(*v,*z,*xb) 302 #endif 303 304 PetscFunctionBegin; 305 ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 306 ierr = VecGetArray(zz,&z);CHKERRQ(ierr); 307 308 idx = a->j; 309 v = a->a; 310 ii = a->i; 311 312 for ( i=0; i<mbs; i++ ) { 313 n = ii[1] - ii[0]; ii++; 314 sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; 315 for ( j=0; j<n; j++ ) { 316 xb = x + 3*(*idx++); x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; 317 sum1 += v[0]*x1 + v[3]*x2 + v[6]*x3; 318 sum2 += v[1]*x1 + v[4]*x2 + v[7]*x3; 319 sum3 += v[2]*x1 + v[5]*x2 + v[8]*x3; 320 v += 9; 321 } 322 z[0] = sum1; z[1] = sum2; z[2] = sum3; 323 z += 3; 324 } 325 ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 326 ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr); 327 PLogFlops(18*a->nz - a->m); 328 PetscFunctionReturn(0); 329 } 330 331 #undef __FUNC__ 332 #define __FUNC__ "MatMult_SeqBAIJ_4" 333 int MatMult_SeqBAIJ_4(Mat A,Vec xx,Vec zz) 334 { 335 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 336 Scalar *x,*z,*xb,sum1,sum2,sum3,sum4,x1,x2,x3,x4; 337 MatScalar *v; 338 int ierr,mbs=a->mbs,i,*idx,*ii,j,n; 339 340 PetscFunctionBegin; 341 ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 342 ierr = VecGetArray(zz,&z);CHKERRQ(ierr); 343 344 idx = a->j; 345 v = a->a; 346 ii = a->i; 347 348 for ( i=0; i<mbs; i++ ) { 349 n = ii[1] - ii[0]; ii++; 350 sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; 351 for ( j=0; j<n; j++ ) { 352 xb = x + 4*(*idx++); 353 x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; 354 sum1 += v[0]*x1 + v[4]*x2 + v[8]*x3 + v[12]*x4; 355 sum2 += v[1]*x1 + v[5]*x2 + v[9]*x3 + v[13]*x4; 356 sum3 += v[2]*x1 + v[6]*x2 + v[10]*x3 + v[14]*x4; 357 sum4 += v[3]*x1 + v[7]*x2 + v[11]*x3 + v[15]*x4; 358 v += 16; 359 } 360 z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; 361 z += 4; 362 } 363 ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 364 ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr); 365 PLogFlops(32*a->nz - a->m); 366 PetscFunctionReturn(0); 367 } 368 369 #undef __FUNC__ 370 #define __FUNC__ "MatMult_SeqBAIJ_5" 371 int MatMult_SeqBAIJ_5(Mat A,Vec xx,Vec zz) 372 { 373 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 374 Scalar sum1,sum2,sum3,sum4,sum5,x1,x2,x3,x4,x5,*xb,*z,*x; 375 MatScalar *v; 376 int ierr,mbs=a->mbs,i,*idx,*ii,j,n; 377 378 PetscFunctionBegin; 379 ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 380 ierr = VecGetArray(zz,&z);CHKERRQ(ierr); 381 382 idx = a->j; 383 v = a->a; 384 ii = a->i; 385 386 for ( i=0; i<mbs; i++ ) { 387 n = ii[1] - ii[0]; ii++; 388 sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; 389 for ( j=0; j<n; j++ ) { 390 xb = x + 5*(*idx++); 391 x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; 392 sum1 += v[0]*x1 + v[5]*x2 + v[10]*x3 + v[15]*x4 + v[20]*x5; 393 sum2 += v[1]*x1 + v[6]*x2 + v[11]*x3 + v[16]*x4 + v[21]*x5; 394 sum3 += v[2]*x1 + v[7]*x2 + v[12]*x3 + v[17]*x4 + v[22]*x5; 395 sum4 += v[3]*x1 + v[8]*x2 + v[13]*x3 + v[18]*x4 + v[23]*x5; 396 sum5 += v[4]*x1 + v[9]*x2 + v[14]*x3 + v[19]*x4 + v[24]*x5; 397 v += 25; 398 } 399 z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; 400 z += 5; 401 } 402 ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 403 ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr); 404 PLogFlops(50*a->nz - a->m); 405 PetscFunctionReturn(0); 406 } 407 408 409 #undef __FUNC__ 410 #define __FUNC__ "MatMult_SeqBAIJ_6" 411 int MatMult_SeqBAIJ_6(Mat A,Vec xx,Vec zz) 412 { 413 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 414 Scalar *x,*z,*xb,sum1,sum2,sum3,sum4,sum5,sum6; 415 Scalar x1,x2,x3,x4,x5,x6; 416 MatScalar *v; 417 int ierr,mbs=a->mbs,i,*idx,*ii,j,n; 418 419 PetscFunctionBegin; 420 ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 421 ierr = VecGetArray(zz,&z);CHKERRQ(ierr); 422 423 idx = a->j; 424 v = a->a; 425 ii = a->i; 426 427 for ( i=0; i<mbs; i++ ) { 428 n = ii[1] - ii[0]; ii++; 429 sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0; 430 for ( j=0; j<n; j++ ) { 431 xb = x + 6*(*idx++); 432 x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; 433 sum1 += v[0]*x1 + v[6]*x2 + v[12]*x3 + v[18]*x4 + v[24]*x5 + v[30]*x6; 434 sum2 += v[1]*x1 + v[7]*x2 + v[13]*x3 + v[19]*x4 + v[25]*x5 + v[31]*x6; 435 sum3 += v[2]*x1 + v[8]*x2 + v[14]*x3 + v[20]*x4 + v[26]*x5 + v[32]*x6; 436 sum4 += v[3]*x1 + v[9]*x2 + v[15]*x3 + v[21]*x4 + v[27]*x5 + v[33]*x6; 437 sum5 += v[4]*x1 + v[10]*x2 + v[16]*x3 + v[22]*x4 + v[28]*x5 + v[34]*x6; 438 sum6 += v[5]*x1 + v[11]*x2 + v[17]*x3 + v[23]*x4 + v[29]*x5 + v[35]*x6; 439 v += 36; 440 } 441 z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; 442 z += 6; 443 } 444 445 ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 446 ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr); 447 PLogFlops(72*a->nz - a->m); 448 PetscFunctionReturn(0); 449 } 450 #undef __FUNC__ 451 #define __FUNC__ "MatMult_SeqBAIJ_7" 452 int MatMult_SeqBAIJ_7(Mat A,Vec xx,Vec zz) 453 { 454 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 455 Scalar *x,*z,*xb,sum1,sum2,sum3,sum4,sum5,sum6,sum7; 456 Scalar x1,x2,x3,x4,x5,x6,x7; 457 MatScalar *v; 458 int ierr,mbs=a->mbs,i,*idx,*ii,j,n; 459 460 PetscFunctionBegin; 461 ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 462 ierr = VecGetArray(zz,&z);CHKERRQ(ierr); 463 464 idx = a->j; 465 v = a->a; 466 ii = a->i; 467 468 for ( i=0; i<mbs; i++ ) { 469 n = ii[1] - ii[0]; ii++; 470 sum1 = 0.0; sum2 = 0.0; sum3 = 0.0; sum4 = 0.0; sum5 = 0.0; sum6 = 0.0; sum7 = 0.0; 471 for ( j=0; j<n; j++ ) { 472 xb = x + 7*(*idx++); 473 x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6]; 474 sum1 += v[0]*x1 + v[7]*x2 + v[14]*x3 + v[21]*x4 + v[28]*x5 + v[35]*x6 + v[42]*x7; 475 sum2 += v[1]*x1 + v[8]*x2 + v[15]*x3 + v[22]*x4 + v[29]*x5 + v[36]*x6 + v[43]*x7; 476 sum3 += v[2]*x1 + v[9]*x2 + v[16]*x3 + v[23]*x4 + v[30]*x5 + v[37]*x6 + v[44]*x7; 477 sum4 += v[3]*x1 + v[10]*x2 + v[17]*x3 + v[24]*x4 + v[31]*x5 + v[38]*x6 + v[45]*x7; 478 sum5 += v[4]*x1 + v[11]*x2 + v[18]*x3 + v[25]*x4 + v[32]*x5 + v[39]*x6 + v[46]*x7; 479 sum6 += v[5]*x1 + v[12]*x2 + v[19]*x3 + v[26]*x4 + v[33]*x5 + v[40]*x6 + v[47]*x7; 480 sum7 += v[6]*x1 + v[13]*x2 + v[20]*x3 + v[27]*x4 + v[34]*x5 + v[41]*x6 + v[48]*x7; 481 v += 49; 482 } 483 z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7; 484 z += 7; 485 } 486 487 ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 488 ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr); 489 PLogFlops(98*a->nz - a->m); 490 PetscFunctionReturn(0); 491 } 492 493 /* 494 This will not work with MatScalar == float because it calls the BLAS 495 */ 496 #undef __FUNC__ 497 #define __FUNC__ "MatMult_SeqBAIJ_N" 498 int MatMult_SeqBAIJ_N(Mat A,Vec xx,Vec zz) 499 { 500 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 501 Scalar *x,*z,*xb, *work,*workt; 502 MatScalar *v; 503 int ierr,mbs=a->mbs,i,*idx,*ii,bs=a->bs,j,n,bs2=a->bs2; 504 int ncols,k; 505 506 PetscFunctionBegin; 507 ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 508 ierr = VecGetArray(zz,&z);CHKERRQ(ierr); 509 510 idx = a->j; 511 v = a->a; 512 ii = a->i; 513 514 515 if (!a->mult_work) { 516 k = PetscMax(a->m,a->n); 517 a->mult_work = (Scalar *) PetscMalloc((k+1)*sizeof(Scalar));CHKPTRQ(a->mult_work); 518 } 519 work = a->mult_work; 520 for ( i=0; i<mbs; i++ ) { 521 n = ii[1] - ii[0]; ii++; 522 ncols = n*bs; 523 workt = work; 524 for ( j=0; j<n; j++ ) { 525 xb = x + bs*(*idx++); 526 for ( k=0; k<bs; k++ ) workt[k] = xb[k]; 527 workt += bs; 528 } 529 Kernel_w_gets_Ar_times_v(bs,ncols,work,v,z); 530 /* LAgemv_("N",&bs,&ncols,&_DOne,v,&bs,work,&_One,&_DZero,z,&_One); */ 531 v += n*bs2; 532 z += bs; 533 } 534 ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 535 ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr); 536 PLogFlops(2*a->nz*bs2 - a->m); 537 PetscFunctionReturn(0); 538 } 539 540 #undef __FUNC__ 541 #define __FUNC__ "MatMultAdd_SeqBAIJ_1" 542 int MatMultAdd_SeqBAIJ_1(Mat A,Vec xx,Vec yy,Vec zz) 543 { 544 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 545 Scalar *x,*y,*z,sum; 546 MatScalar *v; 547 int ierr,mbs=a->mbs,i,*idx,*ii,n; 548 549 PetscFunctionBegin; 550 ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 551 ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 552 if (zz != yy) { 553 ierr = VecGetArray(zz,&z);CHKERRQ(ierr); 554 } else { 555 z = y; 556 } 557 558 idx = a->j; 559 v = a->a; 560 ii = a->i; 561 562 for ( i=0; i<mbs; i++ ) { 563 n = ii[1] - ii[0]; ii++; 564 sum = y[i]; 565 while (n--) sum += *v++ * x[*idx++]; 566 z[i] = sum; 567 } 568 ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 569 ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 570 if (zz != yy) { 571 ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr); 572 } 573 PLogFlops(2*a->nz); 574 PetscFunctionReturn(0); 575 } 576 577 #undef __FUNC__ 578 #define __FUNC__ "MatMultAdd_SeqBAIJ_2" 579 int MatMultAdd_SeqBAIJ_2(Mat A,Vec xx,Vec yy,Vec zz) 580 { 581 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 582 Scalar *x,*y,*z,*xb,sum1,sum2; 583 Scalar x1,x2; 584 MatScalar *v; 585 int ierr,mbs=a->mbs,i,*idx,*ii,j,n; 586 587 PetscFunctionBegin; 588 ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 589 ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 590 if (zz != yy) { 591 ierr = VecGetArray(zz,&z);CHKERRQ(ierr); 592 } else { 593 z = y; 594 } 595 596 idx = a->j; 597 v = a->a; 598 ii = a->i; 599 600 for ( i=0; i<mbs; i++ ) { 601 n = ii[1] - ii[0]; ii++; 602 sum1 = y[0]; sum2 = y[1]; 603 for ( j=0; j<n; j++ ) { 604 xb = x + 2*(*idx++); x1 = xb[0]; x2 = xb[1]; 605 sum1 += v[0]*x1 + v[2]*x2; 606 sum2 += v[1]*x1 + v[3]*x2; 607 v += 4; 608 } 609 z[0] = sum1; z[1] = sum2; 610 z += 2; y += 2; 611 } 612 ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 613 ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 614 if (zz != yy) { 615 ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr); 616 } 617 PLogFlops(4*a->nz); 618 PetscFunctionReturn(0); 619 } 620 621 #undef __FUNC__ 622 #define __FUNC__ "MatMultAdd_SeqBAIJ_3" 623 int MatMultAdd_SeqBAIJ_3(Mat A,Vec xx,Vec yy,Vec zz) 624 { 625 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 626 Scalar *x,*y,*z,*xb,sum1,sum2,sum3,x1,x2,x3; 627 MatScalar *v; 628 int ierr,mbs=a->mbs,i,*idx,*ii,j,n; 629 630 PetscFunctionBegin; 631 ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 632 ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 633 if (zz != yy) { 634 ierr = VecGetArray(zz,&z);CHKERRQ(ierr); 635 } else { 636 z = y; 637 } 638 639 idx = a->j; 640 v = a->a; 641 ii = a->i; 642 643 for ( i=0; i<mbs; i++ ) { 644 n = ii[1] - ii[0]; ii++; 645 sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; 646 for ( j=0; j<n; j++ ) { 647 xb = x + 3*(*idx++); x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; 648 sum1 += v[0]*x1 + v[3]*x2 + v[6]*x3; 649 sum2 += v[1]*x1 + v[4]*x2 + v[7]*x3; 650 sum3 += v[2]*x1 + v[5]*x2 + v[8]*x3; 651 v += 9; 652 } 653 z[0] = sum1; z[1] = sum2; z[2] = sum3; 654 z += 3; y += 3; 655 } 656 ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 657 ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 658 if (zz != yy) { 659 ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr); 660 } 661 PLogFlops(18*a->nz); 662 PetscFunctionReturn(0); 663 } 664 665 #undef __FUNC__ 666 #define __FUNC__ "MatMultAdd_SeqBAIJ_4" 667 int MatMultAdd_SeqBAIJ_4(Mat A,Vec xx,Vec yy,Vec zz) 668 { 669 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 670 Scalar *x,*y,*z,*xb,sum1,sum2,sum3,sum4,x1,x2,x3,x4; 671 MatScalar *v; 672 int ierr,mbs=a->mbs,i,*idx,*ii; 673 int j,n; 674 675 PetscFunctionBegin; 676 ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 677 ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 678 if (zz != yy) { 679 ierr = VecGetArray(zz,&z);CHKERRQ(ierr); 680 } else { 681 z = y; 682 } 683 684 idx = a->j; 685 v = a->a; 686 ii = a->i; 687 688 for ( i=0; i<mbs; i++ ) { 689 n = ii[1] - ii[0]; ii++; 690 sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; 691 for ( j=0; j<n; j++ ) { 692 xb = x + 4*(*idx++); 693 x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; 694 sum1 += v[0]*x1 + v[4]*x2 + v[8]*x3 + v[12]*x4; 695 sum2 += v[1]*x1 + v[5]*x2 + v[9]*x3 + v[13]*x4; 696 sum3 += v[2]*x1 + v[6]*x2 + v[10]*x3 + v[14]*x4; 697 sum4 += v[3]*x1 + v[7]*x2 + v[11]*x3 + v[15]*x4; 698 v += 16; 699 } 700 z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; 701 z += 4; y += 4; 702 } 703 ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 704 ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 705 if (zz != yy) { 706 ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr); 707 } 708 PLogFlops(32*a->nz); 709 PetscFunctionReturn(0); 710 } 711 712 #undef __FUNC__ 713 #define __FUNC__ "MatMultAdd_SeqBAIJ_5" 714 int MatMultAdd_SeqBAIJ_5(Mat A,Vec xx,Vec yy,Vec zz) 715 { 716 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 717 Scalar *x,*y,*z,*xb,sum1,sum2,sum3,sum4,sum5,x1,x2,x3,x4,x5; 718 MatScalar *v; 719 int ierr,mbs=a->mbs,i,*idx,*ii,j,n; 720 721 PetscFunctionBegin; 722 ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 723 ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 724 if (zz != yy) { 725 ierr = VecGetArray(zz,&z);CHKERRQ(ierr); 726 } else { 727 z = y; 728 } 729 730 idx = a->j; 731 v = a->a; 732 ii = a->i; 733 734 for ( i=0; i<mbs; i++ ) { 735 n = ii[1] - ii[0]; ii++; 736 sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4]; 737 for ( j=0; j<n; j++ ) { 738 xb = x + 5*(*idx++); 739 x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; 740 sum1 += v[0]*x1 + v[5]*x2 + v[10]*x3 + v[15]*x4 + v[20]*x5; 741 sum2 += v[1]*x1 + v[6]*x2 + v[11]*x3 + v[16]*x4 + v[21]*x5; 742 sum3 += v[2]*x1 + v[7]*x2 + v[12]*x3 + v[17]*x4 + v[22]*x5; 743 sum4 += v[3]*x1 + v[8]*x2 + v[13]*x3 + v[18]*x4 + v[23]*x5; 744 sum5 += v[4]*x1 + v[9]*x2 + v[14]*x3 + v[19]*x4 + v[24]*x5; 745 v += 25; 746 } 747 z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; 748 z += 5; y += 5; 749 } 750 ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 751 ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 752 if (zz != yy) { 753 ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr); 754 } 755 PLogFlops(50*a->nz); 756 PetscFunctionReturn(0); 757 } 758 #undef __FUNC__ 759 #define __FUNC__ "MatMultAdd_SeqBAIJ_6" 760 int MatMultAdd_SeqBAIJ_6(Mat A,Vec xx,Vec yy,Vec zz) 761 { 762 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 763 Scalar *x,*y,*z,*xb,sum1,sum2,sum3,sum4,sum5,sum6; 764 Scalar x1,x2,x3,x4,x5,x6; 765 MatScalar *v; 766 int ierr,mbs=a->mbs,i,*idx,*ii,j,n; 767 768 PetscFunctionBegin; 769 ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 770 ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 771 if (zz != yy) { 772 ierr = VecGetArray(zz,&z);CHKERRQ(ierr); 773 } else { 774 z = y; 775 } 776 777 idx = a->j; 778 v = a->a; 779 ii = a->i; 780 781 for ( i=0; i<mbs; i++ ) { 782 n = ii[1] - ii[0]; ii++; 783 sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4]; sum6 = y[5]; 784 for ( j=0; j<n; j++ ) { 785 xb = x + 6*(*idx++); 786 x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; 787 sum1 += v[0]*x1 + v[6]*x2 + v[12]*x3 + v[18]*x4 + v[24]*x5 + v[30]*x6; 788 sum2 += v[1]*x1 + v[7]*x2 + v[13]*x3 + v[19]*x4 + v[25]*x5 + v[31]*x6; 789 sum3 += v[2]*x1 + v[8]*x2 + v[14]*x3 + v[20]*x4 + v[26]*x5 + v[32]*x6; 790 sum4 += v[3]*x1 + v[9]*x2 + v[15]*x3 + v[21]*x4 + v[27]*x5 + v[33]*x6; 791 sum5 += v[4]*x1 + v[10]*x2 + v[16]*x3 + v[22]*x4 + v[28]*x5 + v[34]*x6; 792 sum6 += v[5]*x1 + v[11]*x2 + v[17]*x3 + v[23]*x4 + v[29]*x5 + v[35]*x6; 793 v += 36; 794 } 795 z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; 796 z += 6; y += 6; 797 } 798 ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 799 ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 800 if (zz != yy) { 801 ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr); 802 } 803 PLogFlops(72*a->nz); 804 PetscFunctionReturn(0); 805 } 806 807 #undef __FUNC__ 808 #define __FUNC__ "MatMultAdd_SeqBAIJ_7" 809 int MatMultAdd_SeqBAIJ_7(Mat A,Vec xx,Vec yy,Vec zz) 810 { 811 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 812 Scalar *x,*y,*z,*xb,sum1,sum2,sum3,sum4,sum5,sum6,sum7; 813 Scalar x1,x2,x3,x4,x5,x6,x7; 814 MatScalar *v; 815 int ierr,mbs=a->mbs,i,*idx,*ii,j,n; 816 817 PetscFunctionBegin; 818 ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 819 ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 820 if (zz != yy) { 821 ierr = VecGetArray(zz,&z);CHKERRQ(ierr); 822 } else { 823 z = y; 824 } 825 826 idx = a->j; 827 v = a->a; 828 ii = a->i; 829 830 for ( i=0; i<mbs; i++ ) { 831 n = ii[1] - ii[0]; ii++; 832 sum1 = y[0]; sum2 = y[1]; sum3 = y[2]; sum4 = y[3]; sum5 = y[4]; sum6 = y[5]; sum7 = y[6]; 833 for ( j=0; j<n; j++ ) { 834 xb = x + 7*(*idx++); 835 x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; x5 = xb[4]; x6 = xb[5]; x7 = xb[6]; 836 sum1 += v[0]*x1 + v[7]*x2 + v[14]*x3 + v[21]*x4 + v[28]*x5 + v[35]*x6 + v[42]*x7; 837 sum2 += v[1]*x1 + v[8]*x2 + v[15]*x3 + v[22]*x4 + v[29]*x5 + v[36]*x6 + v[43]*x7; 838 sum3 += v[2]*x1 + v[9]*x2 + v[16]*x3 + v[23]*x4 + v[30]*x5 + v[37]*x6 + v[44]*x7; 839 sum4 += v[3]*x1 + v[10]*x2 + v[17]*x3 + v[24]*x4 + v[31]*x5 + v[38]*x6 + v[45]*x7; 840 sum5 += v[4]*x1 + v[11]*x2 + v[18]*x3 + v[25]*x4 + v[32]*x5 + v[39]*x6 + v[46]*x7; 841 sum6 += v[5]*x1 + v[12]*x2 + v[19]*x3 + v[26]*x4 + v[33]*x5 + v[40]*x6 + v[47]*x7; 842 sum7 += v[6]*x1 + v[13]*x2 + v[20]*x3 + v[27]*x4 + v[34]*x5 + v[41]*x6 + v[48]*x7; 843 v += 49; 844 } 845 z[0] = sum1; z[1] = sum2; z[2] = sum3; z[3] = sum4; z[4] = sum5; z[5] = sum6; z[6] = sum7; 846 z += 7; y += 7; 847 } 848 ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 849 ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 850 if (zz != yy) { 851 ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr); 852 } 853 PLogFlops(98*a->nz); 854 PetscFunctionReturn(0); 855 } 856 857 #undef __FUNC__ 858 #define __FUNC__ "MatMultAdd_SeqBAIJ_N" 859 int MatMultAdd_SeqBAIJ_N(Mat A,Vec xx,Vec yy,Vec zz) 860 { 861 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 862 Scalar *x,*z,*xb,*work,*workt; 863 MatScalar *v; 864 int mbs=a->mbs,i,*idx,*ii,bs=a->bs,j,n,bs2=a->bs2,ierr; 865 int ncols,k; 866 867 PetscFunctionBegin; 868 if ( xx != yy) { ierr = VecCopy(yy,zz);CHKERRQ(ierr); } 869 870 ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 871 ierr = VecGetArray(zz,&z);CHKERRQ(ierr); 872 873 idx = a->j; 874 v = a->a; 875 ii = a->i; 876 877 878 if (!a->mult_work) { 879 k = PetscMax(a->m,a->n); 880 a->mult_work = (Scalar *) PetscMalloc((k+1)*sizeof(Scalar));CHKPTRQ(a->mult_work); 881 } 882 work = a->mult_work; 883 for ( i=0; i<mbs; i++ ) { 884 n = ii[1] - ii[0]; ii++; 885 ncols = n*bs; 886 workt = work; 887 for ( j=0; j<n; j++ ) { 888 xb = x + bs*(*idx++); 889 for ( k=0; k<bs; k++ ) workt[k] = xb[k]; 890 workt += bs; 891 } 892 Kernel_w_gets_w_plus_Ar_times_v(bs,ncols,work,v,z); 893 /* LAgemv_("N",&bs,&ncols,&_DOne,v,&bs,work,&_One,&_DOne,z,&_One); */ 894 v += n*bs2; 895 z += bs; 896 } 897 ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 898 ierr = VecRestoreArray(zz,&z);CHKERRQ(ierr); 899 PLogFlops(2*a->nz*bs2 ); 900 PetscFunctionReturn(0); 901 } 902 903 #undef __FUNC__ 904 #define __FUNC__ "MatMultTrans_SeqBAIJ" 905 int MatMultTrans_SeqBAIJ(Mat A,Vec xx,Vec zz) 906 { 907 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 908 Scalar *xg,*zg,*zb; 909 Scalar *x,*z,*xb,x1,x2,x3,x4,x5; 910 MatScalar *v; 911 int mbs=a->mbs,i,*idx,*ii,*ai=a->i,rval,N=a->n; 912 int bs=a->bs,j,n,bs2=a->bs2,*ib,ierr; 913 914 915 PetscFunctionBegin; 916 ierr = VecGetArray(xx,&xg);CHKERRQ(ierr); x = xg; 917 ierr = VecGetArray(zz,&zg);CHKERRQ(ierr); z = zg; 918 ierr = PetscMemzero(z,N*sizeof(Scalar));CHKERRQ(ierr); 919 920 idx = a->j; 921 v = a->a; 922 ii = a->i; 923 924 switch (bs) { 925 case 1: 926 for ( i=0; i<mbs; i++ ) { 927 n = ii[1] - ii[0]; ii++; 928 xb = x + i; x1 = xb[0]; 929 ib = idx + ai[i]; 930 for ( j=0; j<n; j++ ) { 931 rval = ib[j]; 932 z[rval] += *v++ * x1; 933 } 934 } 935 break; 936 case 2: 937 for ( i=0; i<mbs; i++ ) { 938 n = ii[1] - ii[0]; ii++; 939 xb = x + 2*i; x1 = xb[0]; x2 = xb[1]; 940 ib = idx + ai[i]; 941 for ( j=0; j<n; j++ ) { 942 rval = ib[j]*2; 943 z[rval++] += v[0]*x1 + v[1]*x2; 944 z[rval++] += v[2]*x1 + v[3]*x2; 945 v += 4; 946 } 947 } 948 break; 949 case 3: 950 for ( i=0; i<mbs; i++ ) { 951 n = ii[1] - ii[0]; ii++; 952 xb = x + 3*i; x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; 953 ib = idx + ai[i]; 954 for ( j=0; j<n; j++ ) { 955 rval = ib[j]*3; 956 z[rval++] += v[0]*x1 + v[1]*x2 + v[2]*x3; 957 z[rval++] += v[3]*x1 + v[4]*x2 + v[5]*x3; 958 z[rval++] += v[6]*x1 + v[7]*x2 + v[8]*x3; 959 v += 9; 960 } 961 } 962 break; 963 case 4: 964 for ( i=0; i<mbs; i++ ) { 965 n = ii[1] - ii[0]; ii++; 966 xb = x + 4*i; x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; 967 ib = idx + ai[i]; 968 for ( j=0; j<n; j++ ) { 969 rval = ib[j]*4; 970 z[rval++] += v[0]*x1 + v[1]*x2 + v[2]*x3 + v[3]*x4; 971 z[rval++] += v[4]*x1 + v[5]*x2 + v[6]*x3 + v[7]*x4; 972 z[rval++] += v[8]*x1 + v[9]*x2 + v[10]*x3 + v[11]*x4; 973 z[rval++] += v[12]*x1 + v[13]*x2 + v[14]*x3 + v[15]*x4; 974 v += 16; 975 } 976 } 977 break; 978 case 5: 979 for ( i=0; i<mbs; i++ ) { 980 n = ii[1] - ii[0]; ii++; 981 xb = x + 5*i; x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; 982 x4 = xb[3]; x5 = xb[4]; 983 ib = idx + ai[i]; 984 for ( j=0; j<n; j++ ) { 985 rval = ib[j]*5; 986 z[rval++] += v[0]*x1 + v[1]*x2 + v[2]*x3 + v[3]*x4 + v[4]*x5; 987 z[rval++] += v[5]*x1 + v[6]*x2 + v[7]*x3 + v[8]*x4 + v[9]*x5; 988 z[rval++] += v[10]*x1 + v[11]*x2 + v[12]*x3 + v[13]*x4 + v[14]*x5; 989 z[rval++] += v[15]*x1 + v[16]*x2 + v[17]*x3 + v[18]*x4 + v[19]*x5; 990 z[rval++] += v[20]*x1 + v[21]*x2 + v[22]*x3 + v[23]*x4 + v[24]*x5; 991 v += 25; 992 } 993 } 994 break; 995 /* block sizes larger then 5 by 5 are handled by BLAS */ 996 default: { 997 int ncols,k; 998 Scalar *work,*workt; 999 1000 if (!a->mult_work) { 1001 k = PetscMax(a->m,a->n); 1002 a->mult_work = (Scalar *) PetscMalloc((k+1)*sizeof(Scalar));CHKPTRQ(a->mult_work); 1003 } 1004 work = a->mult_work; 1005 for ( i=0; i<mbs; i++ ) { 1006 n = ii[1] - ii[0]; ii++; 1007 ncols = n*bs; 1008 ierr = PetscMemzero(work,ncols*sizeof(Scalar));CHKERRQ(ierr); 1009 Kernel_w_gets_w_plus_trans_Ar_times_v(bs,ncols,x,v,work); 1010 /* LAgemv_("T",&bs,&ncols,&_DOne,v,&bs,x,&_One,&_DOne,work,&_One); */ 1011 v += n*bs2; 1012 x += bs; 1013 workt = work; 1014 for ( j=0; j<n; j++ ) { 1015 zb = z + bs*(*idx++); 1016 for ( k=0; k<bs; k++ ) zb[k] += workt[k] ; 1017 workt += bs; 1018 } 1019 } 1020 } 1021 } 1022 ierr = VecRestoreArray(xx,&xg);CHKERRQ(ierr); 1023 ierr = VecRestoreArray(zz,&zg);CHKERRQ(ierr); 1024 PLogFlops(2*a->nz*a->bs2 - a->n); 1025 PetscFunctionReturn(0); 1026 } 1027 1028 #undef __FUNC__ 1029 #define __FUNC__ "MatMultTransAdd_SeqBAIJ" 1030 int MatMultTransAdd_SeqBAIJ(Mat A,Vec xx,Vec yy,Vec zz) 1031 1032 { 1033 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 1034 Scalar *xg,*zg,*zb; 1035 Scalar *x,*z,*xb,x1,x2,x3,x4,x5; 1036 MatScalar *v; 1037 int mbs=a->mbs,i,*idx,*ii,*ai=a->i,rval,N=a->n; 1038 int bs=a->bs,j,n,bs2=a->bs2,*ib,ierr; 1039 1040 PetscFunctionBegin; 1041 ierr = VecGetArray(xx,&xg);CHKERRQ(ierr); x = xg; 1042 ierr = VecGetArray(zz,&zg);CHKERRQ(ierr); z = zg; 1043 1044 if ( yy != zz ) { ierr = VecCopy(yy,zz);CHKERRQ(ierr); } 1045 else { 1046 ierr = PetscMemzero(z,N*sizeof(Scalar));CHKERRQ(ierr); 1047 } 1048 1049 idx = a->j; 1050 v = a->a; 1051 ii = a->i; 1052 1053 switch (bs) { 1054 case 1: 1055 for ( i=0; i<mbs; i++ ) { 1056 n = ii[1] - ii[0]; ii++; 1057 xb = x + i; x1 = xb[0]; 1058 ib = idx + ai[i]; 1059 for ( j=0; j<n; j++ ) { 1060 rval = ib[j]; 1061 z[rval] += *v++ * x1; 1062 } 1063 } 1064 break; 1065 case 2: 1066 for ( i=0; i<mbs; i++ ) { 1067 n = ii[1] - ii[0]; ii++; 1068 xb = x + 2*i; x1 = xb[0]; x2 = xb[1]; 1069 ib = idx + ai[i]; 1070 for ( j=0; j<n; j++ ) { 1071 rval = ib[j]*2; 1072 z[rval++] += v[0]*x1 + v[1]*x2; 1073 z[rval++] += v[2]*x1 + v[3]*x2; 1074 v += 4; 1075 } 1076 } 1077 break; 1078 case 3: 1079 for ( i=0; i<mbs; i++ ) { 1080 n = ii[1] - ii[0]; ii++; 1081 xb = x + 3*i; x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; 1082 ib = idx + ai[i]; 1083 for ( j=0; j<n; j++ ) { 1084 rval = ib[j]*3; 1085 z[rval++] += v[0]*x1 + v[1]*x2 + v[2]*x3; 1086 z[rval++] += v[3]*x1 + v[4]*x2 + v[5]*x3; 1087 z[rval++] += v[6]*x1 + v[7]*x2 + v[8]*x3; 1088 v += 9; 1089 } 1090 } 1091 break; 1092 case 4: 1093 for ( i=0; i<mbs; i++ ) { 1094 n = ii[1] - ii[0]; ii++; 1095 xb = x + 4*i; x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; x4 = xb[3]; 1096 ib = idx + ai[i]; 1097 for ( j=0; j<n; j++ ) { 1098 rval = ib[j]*4; 1099 z[rval++] += v[0]*x1 + v[1]*x2 + v[2]*x3 + v[3]*x4; 1100 z[rval++] += v[4]*x1 + v[5]*x2 + v[6]*x3 + v[7]*x4; 1101 z[rval++] += v[8]*x1 + v[9]*x2 + v[10]*x3 + v[11]*x4; 1102 z[rval++] += v[12]*x1 + v[13]*x2 + v[14]*x3 + v[15]*x4; 1103 v += 16; 1104 } 1105 } 1106 break; 1107 case 5: 1108 for ( i=0; i<mbs; i++ ) { 1109 n = ii[1] - ii[0]; ii++; 1110 xb = x + 5*i; x1 = xb[0]; x2 = xb[1]; x3 = xb[2]; 1111 x4 = xb[3]; x5 = xb[4]; 1112 ib = idx + ai[i]; 1113 for ( j=0; j<n; j++ ) { 1114 rval = ib[j]*5; 1115 z[rval++] += v[0]*x1 + v[1]*x2 + v[2]*x3 + v[3]*x4 + v[4]*x5; 1116 z[rval++] += v[5]*x1 + v[6]*x2 + v[7]*x3 + v[8]*x4 + v[9]*x5; 1117 z[rval++] += v[10]*x1 + v[11]*x2 + v[12]*x3 + v[13]*x4 + v[14]*x5; 1118 z[rval++] += v[15]*x1 + v[16]*x2 + v[17]*x3 + v[18]*x4 + v[19]*x5; 1119 z[rval++] += v[20]*x1 + v[21]*x2 + v[22]*x3 + v[23]*x4 + v[24]*x5; 1120 v += 25; 1121 } 1122 } 1123 break; 1124 /* block sizes larger then 5 by 5 are handled by BLAS */ 1125 default: { 1126 int ncols,k; 1127 Scalar *work,*workt; 1128 1129 if (!a->mult_work) { 1130 k = PetscMax(a->m,a->n); 1131 a->mult_work = (Scalar *) PetscMalloc((k+1)*sizeof(Scalar));CHKPTRQ(a->mult_work); 1132 } 1133 work = a->mult_work; 1134 for ( i=0; i<mbs; i++ ) { 1135 n = ii[1] - ii[0]; ii++; 1136 ncols = n*bs; 1137 ierr = PetscMemzero(work,ncols*sizeof(Scalar));CHKERRQ(ierr); 1138 Kernel_w_gets_w_plus_trans_Ar_times_v(bs,ncols,x,v,work); 1139 /* LAgemv_("T",&bs,&ncols,&_DOne,v,&bs,x,&_One,&_DOne,work,&_One); */ 1140 v += n*bs2; 1141 x += bs; 1142 workt = work; 1143 for ( j=0; j<n; j++ ) { 1144 zb = z + bs*(*idx++); 1145 for ( k=0; k<bs; k++ ) zb[k] += workt[k] ; 1146 workt += bs; 1147 } 1148 } 1149 } 1150 } 1151 ierr = VecRestoreArray(xx,&xg);CHKERRQ(ierr); 1152 ierr = VecRestoreArray(zz,&zg);CHKERRQ(ierr); 1153 PLogFlops(2*a->nz*a->bs2); 1154 PetscFunctionReturn(0); 1155 } 1156 1157 #undef __FUNC__ 1158 #define __FUNC__ "MatScale_SeqBAIJ" 1159 int MatScale_SeqBAIJ(Scalar *alpha,Mat inA) 1160 { 1161 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) inA->data; 1162 int one = 1, totalnz = a->bs2*a->nz; 1163 1164 PetscFunctionBegin; 1165 BLscal_( &totalnz, alpha, a->a, &one ); 1166 PLogFlops(totalnz); 1167 PetscFunctionReturn(0); 1168 } 1169 1170 #undef __FUNC__ 1171 #define __FUNC__ "MatNorm_SeqBAIJ" 1172 int MatNorm_SeqBAIJ(Mat A,NormType type,double *norm) 1173 { 1174 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 1175 MatScalar *v = a->a; 1176 double sum = 0.0; 1177 int i,j,k,bs = a->bs, nz=a->nz,bs2=a->bs2; 1178 1179 PetscFunctionBegin; 1180 if (type == NORM_FROBENIUS) { 1181 for (i=0; i< bs2*nz; i++ ) { 1182 #if defined(PETSC_USE_COMPLEX) 1183 sum += PetscReal(PetscConj(*v)*(*v)); v++; 1184 #else 1185 sum += (*v)*(*v); v++; 1186 #endif 1187 } 1188 *norm = sqrt(sum); 1189 } else if (type == NORM_INFINITY) { /* maximum row sum */ 1190 *norm = 0.0; 1191 for ( k=0; k<bs; k++ ) { 1192 for ( j=0; j<a->m; j++ ) { 1193 v = a->a + bs2*a->i[j] + k; 1194 sum = 0.0; 1195 for ( i=0; i<a->i[j+1]-a->i[j]; i++ ) { 1196 sum += PetscAbsScalar(*v); 1197 v += bs; 1198 } 1199 if (sum > *norm) *norm = sum; 1200 } 1201 } 1202 } else { 1203 SETERRQ(PETSC_ERR_SUP,0,"No support for this norm yet"); 1204 } 1205 PetscFunctionReturn(0); 1206 } 1207 1208 1209 #undef __FUNC__ 1210 #define __FUNC__ "MatEqual_SeqBAIJ" 1211 int MatEqual_SeqBAIJ(Mat A,Mat B, PetscTruth* flg) 1212 { 1213 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)B->data; 1214 1215 PetscFunctionBegin; 1216 if (B->type !=MATSEQBAIJ) SETERRQ(PETSC_ERR_ARG_INCOMP,0,"Matrices must be same type"); 1217 1218 /* If the matrix/block dimensions are not equal, or no of nonzeros or shift */ 1219 if ((a->m != b->m) || (a->n !=b->n) || (a->bs != b->bs)|| (a->nz != b->nz)) { 1220 *flg = PETSC_FALSE; PetscFunctionReturn(0); 1221 } 1222 1223 /* if the a->i are the same */ 1224 if (PetscMemcmp(a->i,b->i, (a->mbs+1)*sizeof(int))) { 1225 *flg = PETSC_FALSE; PetscFunctionReturn(0); 1226 } 1227 1228 /* if a->j are the same */ 1229 if (PetscMemcmp(a->j,b->j,(a->nz)*sizeof(int))) { 1230 *flg = PETSC_FALSE; PetscFunctionReturn(0); 1231 } 1232 1233 /* if a->a are the same */ 1234 if (PetscMemcmp(a->a, b->a,(a->nz)*(a->bs)*(a->bs)*sizeof(Scalar))) { 1235 *flg = PETSC_FALSE; PetscFunctionReturn(0); 1236 } 1237 *flg = PETSC_TRUE; 1238 PetscFunctionReturn(0); 1239 1240 } 1241 1242 #undef __FUNC__ 1243 #define __FUNC__ "MatGetDiagonal_SeqBAIJ" 1244 int MatGetDiagonal_SeqBAIJ(Mat A,Vec v) 1245 { 1246 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 1247 int ierr,i,j,k,n,row,bs,*ai,*aj,ambs,bs2; 1248 Scalar *x, zero = 0.0; 1249 MatScalar *aa,*aa_j; 1250 1251 PetscFunctionBegin; 1252 if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Not for factored matrix"); 1253 bs = a->bs; 1254 aa = a->a; 1255 ai = a->i; 1256 aj = a->j; 1257 ambs = a->mbs; 1258 bs2 = a->bs2; 1259 1260 ierr = VecSet(&zero,v);CHKERRQ(ierr); 1261 ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1262 ierr = VecGetLocalSize(v,&n);CHKERRQ(ierr); 1263 if (n != a->m) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Nonconforming matrix and vector"); 1264 for ( i=0; i<ambs; i++ ) { 1265 for ( j=ai[i]; j<ai[i+1]; j++ ) { 1266 if (aj[j] == i) { 1267 row = i*bs; 1268 aa_j = aa+j*bs2; 1269 for (k=0; k<bs2; k+=(bs+1),row++) x[row] = aa_j[k]; 1270 break; 1271 } 1272 } 1273 } 1274 ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 1275 PetscFunctionReturn(0); 1276 } 1277 1278 #undef __FUNC__ 1279 #define __FUNC__ "MatDiagonalScale_SeqBAIJ" 1280 int MatDiagonalScale_SeqBAIJ(Mat A,Vec ll,Vec rr) 1281 { 1282 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 1283 Scalar *l,*r,x,*li,*ri; 1284 MatScalar *aa,*v; 1285 int ierr,i,j,k,lm,rn,M,m,n,*ai,*aj,mbs,tmp,bs,bs2; 1286 1287 PetscFunctionBegin; 1288 ai = a->i; 1289 aj = a->j; 1290 aa = a->a; 1291 m = a->m; 1292 n = a->n; 1293 bs = a->bs; 1294 mbs = a->mbs; 1295 bs2 = a->bs2; 1296 if (ll) { 1297 ierr = VecGetArray(ll,&l);CHKERRQ(ierr); 1298 ierr = VecGetLocalSize(ll,&lm);CHKERRQ(ierr); 1299 if (lm != m) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Left scaling vector wrong length"); 1300 for ( i=0; i<mbs; i++ ) { /* for each block row */ 1301 M = ai[i+1] - ai[i]; 1302 li = l + i*bs; 1303 v = aa + bs2*ai[i]; 1304 for ( j=0; j<M; j++ ) { /* for each block */ 1305 for ( k=0; k<bs2; k++ ) { 1306 (*v++) *= li[k%bs]; 1307 } 1308 } 1309 } 1310 ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr); 1311 PLogFlops(a->nz); 1312 } 1313 1314 if (rr) { 1315 ierr = VecGetArray(rr,&r);CHKERRQ(ierr); 1316 ierr = VecGetLocalSize(rr,&rn);CHKERRQ(ierr); 1317 if (rn != n) SETERRQ(PETSC_ERR_ARG_SIZ,0,"Right scaling vector wrong length"); 1318 for ( i=0; i<mbs; i++ ) { /* for each block row */ 1319 M = ai[i+1] - ai[i]; 1320 v = aa + bs2*ai[i]; 1321 for ( j=0; j<M; j++ ) { /* for each block */ 1322 ri = r + bs*aj[ai[i]+j]; 1323 for ( k=0; k<bs; k++ ) { 1324 x = ri[k]; 1325 for ( tmp=0; tmp<bs; tmp++ ) (*v++) *= x; 1326 } 1327 } 1328 } 1329 ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr); 1330 PLogFlops(a->nz); 1331 } 1332 PetscFunctionReturn(0); 1333 } 1334 1335 1336 #undef __FUNC__ 1337 #define __FUNC__ "MatGetInfo_SeqBAIJ" 1338 int MatGetInfo_SeqBAIJ(Mat A,MatInfoType flag,MatInfo *info) 1339 { 1340 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 1341 1342 PetscFunctionBegin; 1343 info->rows_global = (double)a->m; 1344 info->columns_global = (double)a->n; 1345 info->rows_local = (double)a->m; 1346 info->columns_local = (double)a->n; 1347 info->block_size = a->bs2; 1348 info->nz_allocated = a->maxnz; 1349 info->nz_used = a->bs2*a->nz; 1350 info->nz_unneeded = (double)(info->nz_allocated - info->nz_used); 1351 /* if (info->nz_unneeded != A->info.nz_unneeded) 1352 printf("space descrepancy: maxnz-nz = %d, nz_unneeded = %d\n",(int)info->nz_unneeded,(int)A->info.nz_unneeded); */ 1353 info->assemblies = A->num_ass; 1354 info->mallocs = a->reallocs; 1355 info->memory = A->mem; 1356 if (A->factor) { 1357 info->fill_ratio_given = A->info.fill_ratio_given; 1358 info->fill_ratio_needed = A->info.fill_ratio_needed; 1359 info->factor_mallocs = A->info.factor_mallocs; 1360 } else { 1361 info->fill_ratio_given = 0; 1362 info->fill_ratio_needed = 0; 1363 info->factor_mallocs = 0; 1364 } 1365 PetscFunctionReturn(0); 1366 } 1367 1368 1369 #undef __FUNC__ 1370 #define __FUNC__ "MatZeroEntries_SeqBAIJ" 1371 int MatZeroEntries_SeqBAIJ(Mat A) 1372 { 1373 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *) A->data; 1374 int ierr; 1375 1376 PetscFunctionBegin; 1377 ierr = PetscMemzero(a->a,a->bs2*a->i[a->mbs]*sizeof(MatScalar));CHKERRQ(ierr); 1378 PetscFunctionReturn(0); 1379 } 1380