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