1 2 /* 3 Defines matrix-matrix product routines for pairs of SeqAIJ matrices 4 C = A * B 5 */ 6 7 #include <../src/mat/impls/aij/seq/aij.h> /*I "petscmat.h" I*/ 8 #include <../src/mat/utils/freespace.h> 9 #include <petscbt.h> 10 #include <../src/mat/impls/dense/seq/dense.h> /*I "petscmat.h" I*/ 11 12 EXTERN_C_BEGIN 13 #undef __FUNCT__ 14 #define __FUNCT__ "MatMatMult_SeqAIJ_SeqAIJ" 15 PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 16 { 17 PetscErrorCode ierr; 18 19 PetscFunctionBegin; 20 if (scall == MAT_INITIAL_MATRIX){ 21 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr); 22 } 23 ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr); 24 PetscFunctionReturn(0); 25 } 26 EXTERN_C_END 27 28 #undef __FUNCT__ 29 #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ" 30 PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C) 31 { 32 PetscErrorCode ierr; 33 PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 34 Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c; 35 PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci,*cj; 36 PetscInt am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N; 37 PetscInt i,j,anzi,brow,bnzj,cnzi,nlnk,*lnk,nspacedouble=0; 38 MatScalar *ca; 39 PetscBT lnkbt; 40 PetscReal afill; 41 42 PetscFunctionBegin; 43 /* Set up */ 44 /* Allocate ci array, arrays for fill computation and */ 45 /* free space for accumulating nonzero column info */ 46 ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 47 ci[0] = 0; 48 49 /* create and initialize a linked list */ 50 nlnk = bn+1; 51 ierr = PetscLLCreate(bn,bn,nlnk,lnk,lnkbt);CHKERRQ(ierr); 52 53 /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */ 54 ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr); 55 current_space = free_space; 56 57 /* Determine symbolic info for each row of the product: */ 58 for (i=0;i<am;i++) { 59 anzi = ai[i+1] - ai[i]; 60 cnzi = 0; 61 j = anzi; 62 aj = a->j + ai[i]; 63 while (j){/* assume cols are almost in increasing order, starting from its end saves computation */ 64 j--; 65 brow = *(aj + j); 66 bnzj = bi[brow+1] - bi[brow]; 67 bjj = bj + bi[brow]; 68 /* add non-zero cols of B into the sorted linked list lnk */ 69 ierr = PetscLLAdd(bnzj,bjj,bn,nlnk,lnk,lnkbt);CHKERRQ(ierr); 70 cnzi += nlnk; 71 } 72 73 /* If free space is not available, make more free space */ 74 /* Double the amount of total space in the list */ 75 if (current_space->local_remaining<cnzi) { 76 ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 77 nspacedouble++; 78 } 79 80 /* Copy data into free space, then initialize lnk */ 81 ierr = PetscLLClean(bn,bn,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 82 current_space->array += cnzi; 83 current_space->local_used += cnzi; 84 current_space->local_remaining -= cnzi; 85 86 ci[i+1] = ci[i] + cnzi; 87 } 88 89 /* Column indices are in the list of free space */ 90 /* Allocate space for cj, initialize cj, and */ 91 /* destroy list of free space and other temporary array(s) */ 92 ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr); 93 ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 94 ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 95 96 /* Allocate space for ca */ 97 ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 98 ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr); 99 100 /* put together the new symbolic matrix */ 101 ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bn,ci,cj,ca,C);CHKERRQ(ierr); 102 103 /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 104 /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */ 105 c = (Mat_SeqAIJ *)((*C)->data); 106 c->free_a = PETSC_TRUE; 107 c->free_ij = PETSC_TRUE; 108 c->nonew = 0; 109 110 /* set MatInfo */ 111 afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5; 112 if (afill < 1.0) afill = 1.0; 113 c->maxnz = ci[am]; 114 c->nz = ci[am]; 115 (*C)->info.mallocs = nspacedouble; 116 (*C)->info.fill_ratio_given = fill; 117 (*C)->info.fill_ratio_needed = afill; 118 119 #if defined(PETSC_USE_INFO) 120 if (ci[am]) { 121 ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr); 122 ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr); 123 } else { 124 ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr); 125 } 126 #endif 127 PetscFunctionReturn(0); 128 } 129 130 131 #undef __FUNCT__ 132 #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqAIJ" 133 PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C) 134 { 135 PetscErrorCode ierr; 136 PetscLogDouble flops=0.0; 137 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 138 Mat_SeqAIJ *b = (Mat_SeqAIJ *)B->data; 139 Mat_SeqAIJ *c = (Mat_SeqAIJ *)C->data; 140 PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j; 141 PetscInt am=A->rmap->N,cm=C->rmap->N; 142 PetscInt i,j,k,anzi,bnzi,cnzi,brow,nextb; 143 MatScalar *aa=a->a,*ba=b->a,*baj,*ca=c->a; 144 145 PetscFunctionBegin; 146 /* clean old values in C */ 147 ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 148 /* Traverse A row-wise. */ 149 /* Build the ith row in C by summing over nonzero columns in A, */ 150 /* the rows of B corresponding to nonzeros of A. */ 151 for (i=0;i<am;i++) { 152 anzi = ai[i+1] - ai[i]; 153 for (j=0;j<anzi;j++) { 154 brow = *aj++; 155 bnzi = bi[brow+1] - bi[brow]; 156 bjj = bj + bi[brow]; 157 baj = ba + bi[brow]; 158 nextb = 0; 159 for (k=0; nextb<bnzi; k++) { 160 if (cj[k] == bjj[nextb]){ /* ccol == bcol */ 161 ca[k] += (*aa)*baj[nextb++]; 162 } 163 } 164 flops += 2*bnzi; 165 aa++; 166 } 167 cnzi = ci[i+1] - ci[i]; 168 ca += cnzi; 169 cj += cnzi; 170 } 171 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 172 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 173 174 ierr = PetscLogFlops(flops);CHKERRQ(ierr); 175 PetscFunctionReturn(0); 176 } 177 178 #undef __FUNCT__ 179 #define __FUNCT__ "MatMatMultTranspose_SeqAIJ_SeqAIJ" 180 PetscErrorCode MatMatMultTranspose_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 181 { 182 PetscErrorCode ierr; 183 184 PetscFunctionBegin; 185 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Not done yet"); 186 if (scall == MAT_INITIAL_MATRIX){ 187 ierr = MatMatMultTransposeSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr); 188 } 189 ierr = MatMatMultTransposeNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr); 190 PetscFunctionReturn(0); 191 } 192 193 #undef __FUNCT__ 194 #define __FUNCT__ "MatMatMultTransposeSymbolic_SeqAIJ_SeqAIJ" 195 PetscErrorCode MatMatMultTransposeSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C) 196 { 197 /* 198 PetscErrorCode ierr; 199 Mat At; 200 PetscInt *ati,*atj; 201 */ 202 PetscFunctionBegin; 203 204 PetscFunctionReturn(0); 205 } 206 207 #undef __FUNCT__ 208 #define __FUNCT__ "MatMatMultTransposeNumeric_SeqAIJ_SeqAIJ" 209 PetscErrorCode MatMatMultTransposeNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C) 210 { 211 /* 212 PetscErrorCode ierr; 213 Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data; 214 PetscInt am=A->rmap->n,anzi,*ai=a->i,*aj=a->j,*bi=b->i,*bj,bnzi,nextb; 215 PetscInt cm=C->rmap->n,*ci=c->i,*cj=c->j,crow,*cjj,i,j,k; 216 PetscLogDouble flops=0.0; 217 MatScalar *aa=a->a,*ba,*ca=c->a,*caj; 218 */ 219 220 PetscFunctionBegin; 221 222 PetscFunctionReturn(0); 223 } 224 225 #undef __FUNCT__ 226 #define __FUNCT__ "MatMatTransposeMult_SeqAIJ_SeqAIJ" 227 PetscErrorCode MatMatTransposeMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) { 228 PetscErrorCode ierr; 229 230 PetscFunctionBegin; 231 if (scall == MAT_INITIAL_MATRIX){ 232 ierr = MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr); 233 } 234 ierr = MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr); 235 PetscFunctionReturn(0); 236 } 237 238 #undef __FUNCT__ 239 #define __FUNCT__ "MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ" 240 PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C) 241 { 242 PetscErrorCode ierr; 243 Mat At; 244 PetscInt *ati,*atj; 245 246 PetscFunctionBegin; 247 /* create symbolic At */ 248 ierr = MatGetSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr); 249 ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,A->cmap->n,A->rmap->n,ati,atj,PETSC_NULL,&At);CHKERRQ(ierr); 250 251 /* get symbolic C=At*B */ 252 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(At,B,fill,C);CHKERRQ(ierr); 253 254 /* clean up */ 255 ierr = MatDestroy(&At);CHKERRQ(ierr); 256 ierr = MatRestoreSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr); 257 PetscFunctionReturn(0); 258 } 259 260 #undef __FUNCT__ 261 #define __FUNCT__ "MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ" 262 PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C) 263 { 264 PetscErrorCode ierr; 265 Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data; 266 PetscInt am=A->rmap->n,anzi,*ai=a->i,*aj=a->j,*bi=b->i,*bj,bnzi,nextb; 267 PetscInt cm=C->rmap->n,*ci=c->i,*cj=c->j,crow,*cjj,i,j,k; 268 PetscLogDouble flops=0.0; 269 MatScalar *aa=a->a,*ba,*ca=c->a,*caj; 270 271 PetscFunctionBegin; 272 /* clear old values in C */ 273 ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 274 275 /* compute A^T*B using outer product (A^T)[:,i]*B[i,:] */ 276 for (i=0;i<am;i++) { 277 bj = b->j + bi[i]; 278 ba = b->a + bi[i]; 279 bnzi = bi[i+1] - bi[i]; 280 anzi = ai[i+1] - ai[i]; 281 for (j=0; j<anzi; j++) { 282 nextb = 0; 283 crow = *aj++; 284 cjj = cj + ci[crow]; 285 caj = ca + ci[crow]; 286 /* perform sparse axpy operation. Note cjj includes bj. */ 287 for (k=0; nextb<bnzi; k++) { 288 if (cjj[k] == *(bj+nextb)) { /* ccol == bcol */ 289 caj[k] += (*aa)*(*(ba+nextb)); 290 nextb++; 291 } 292 } 293 flops += 2*bnzi; 294 aa++; 295 } 296 } 297 298 /* Assemble the final matrix and clean up */ 299 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 300 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 301 ierr = PetscLogFlops(flops);CHKERRQ(ierr); 302 PetscFunctionReturn(0); 303 } 304 305 EXTERN_C_BEGIN 306 #undef __FUNCT__ 307 #define __FUNCT__ "MatMatMult_SeqAIJ_SeqDense" 308 PetscErrorCode MatMatMult_SeqAIJ_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 309 { 310 PetscErrorCode ierr; 311 312 PetscFunctionBegin; 313 if (scall == MAT_INITIAL_MATRIX){ 314 ierr = MatMatMultSymbolic_SeqAIJ_SeqDense(A,B,fill,C);CHKERRQ(ierr); 315 } 316 ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,B,*C);CHKERRQ(ierr); 317 PetscFunctionReturn(0); 318 } 319 EXTERN_C_END 320 321 #undef __FUNCT__ 322 #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqDense" 323 PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C) 324 { 325 PetscErrorCode ierr; 326 327 PetscFunctionBegin; 328 ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,0.0,C);CHKERRQ(ierr); 329 PetscFunctionReturn(0); 330 } 331 332 #undef __FUNCT__ 333 #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqDense" 334 PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqDense(Mat A,Mat B,Mat C) 335 { 336 Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data; 337 PetscErrorCode ierr; 338 PetscScalar *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4; 339 MatScalar *aa; 340 PetscInt cm=C->rmap->n, cn=B->cmap->n, bm=B->rmap->n, col, i,j,n,*aj, am = A->rmap->n; 341 PetscInt am2 = 2*am, am3 = 3*am, bm4 = 4*bm,colam; 342 343 PetscFunctionBegin; 344 if (!cm || !cn) PetscFunctionReturn(0); 345 if (bm != A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Number columns in A %D not equal rows in B %D\n",A->cmap->n,bm); 346 if (A->rmap->n != C->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Number rows in C %D not equal rows in A %D\n",C->rmap->n,A->rmap->n); 347 if (B->cmap->n != C->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Number columns in B %D not equal columns in C %D\n",B->cmap->n,C->cmap->n); 348 ierr = MatGetArray(B,&b);CHKERRQ(ierr); 349 ierr = MatGetArray(C,&c);CHKERRQ(ierr); 350 b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm; 351 for (col=0; col<cn-4; col += 4){ /* over columns of C */ 352 colam = col*am; 353 for (i=0; i<am; i++) { /* over rows of C in those columns */ 354 r1 = r2 = r3 = r4 = 0.0; 355 n = a->i[i+1] - a->i[i]; 356 aj = a->j + a->i[i]; 357 aa = a->a + a->i[i]; 358 for (j=0; j<n; j++) { 359 r1 += (*aa)*b1[*aj]; 360 r2 += (*aa)*b2[*aj]; 361 r3 += (*aa)*b3[*aj]; 362 r4 += (*aa++)*b4[*aj++]; 363 } 364 c[colam + i] = r1; 365 c[colam + am + i] = r2; 366 c[colam + am2 + i] = r3; 367 c[colam + am3 + i] = r4; 368 } 369 b1 += bm4; 370 b2 += bm4; 371 b3 += bm4; 372 b4 += bm4; 373 } 374 for (;col<cn; col++){ /* over extra columns of C */ 375 for (i=0; i<am; i++) { /* over rows of C in those columns */ 376 r1 = 0.0; 377 n = a->i[i+1] - a->i[i]; 378 aj = a->j + a->i[i]; 379 aa = a->a + a->i[i]; 380 381 for (j=0; j<n; j++) { 382 r1 += (*aa++)*b1[*aj++]; 383 } 384 c[col*am + i] = r1; 385 } 386 b1 += bm; 387 } 388 ierr = PetscLogFlops(cn*(2.0*a->nz));CHKERRQ(ierr); 389 ierr = MatRestoreArray(B,&b);CHKERRQ(ierr); 390 ierr = MatRestoreArray(C,&c);CHKERRQ(ierr); 391 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 392 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 393 PetscFunctionReturn(0); 394 } 395 396 /* 397 Note very similar to MatMult_SeqAIJ(), should generate both codes from same base 398 */ 399 #undef __FUNCT__ 400 #define __FUNCT__ "MatMatMultNumericAdd_SeqAIJ_SeqDense" 401 PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat A,Mat B,Mat C) 402 { 403 Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data; 404 PetscErrorCode ierr; 405 PetscScalar *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4; 406 MatScalar *aa; 407 PetscInt cm=C->rmap->n, cn=B->cmap->n, bm=B->rmap->n, col, i,j,n,*aj, am = A->rmap->n,*ii,arm; 408 PetscInt am2 = 2*am, am3 = 3*am, bm4 = 4*bm,colam,*ridx; 409 410 PetscFunctionBegin; 411 if (!cm || !cn) PetscFunctionReturn(0); 412 ierr = MatGetArray(B,&b);CHKERRQ(ierr); 413 ierr = MatGetArray(C,&c);CHKERRQ(ierr); 414 b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm; 415 416 if (a->compressedrow.use){ /* use compressed row format */ 417 for (col=0; col<cn-4; col += 4){ /* over columns of C */ 418 colam = col*am; 419 arm = a->compressedrow.nrows; 420 ii = a->compressedrow.i; 421 ridx = a->compressedrow.rindex; 422 for (i=0; i<arm; i++) { /* over rows of C in those columns */ 423 r1 = r2 = r3 = r4 = 0.0; 424 n = ii[i+1] - ii[i]; 425 aj = a->j + ii[i]; 426 aa = a->a + ii[i]; 427 for (j=0; j<n; j++) { 428 r1 += (*aa)*b1[*aj]; 429 r2 += (*aa)*b2[*aj]; 430 r3 += (*aa)*b3[*aj]; 431 r4 += (*aa++)*b4[*aj++]; 432 } 433 c[colam + ridx[i]] += r1; 434 c[colam + am + ridx[i]] += r2; 435 c[colam + am2 + ridx[i]] += r3; 436 c[colam + am3 + ridx[i]] += r4; 437 } 438 b1 += bm4; 439 b2 += bm4; 440 b3 += bm4; 441 b4 += bm4; 442 } 443 for (;col<cn; col++){ /* over extra columns of C */ 444 colam = col*am; 445 arm = a->compressedrow.nrows; 446 ii = a->compressedrow.i; 447 ridx = a->compressedrow.rindex; 448 for (i=0; i<arm; i++) { /* over rows of C in those columns */ 449 r1 = 0.0; 450 n = ii[i+1] - ii[i]; 451 aj = a->j + ii[i]; 452 aa = a->a + ii[i]; 453 454 for (j=0; j<n; j++) { 455 r1 += (*aa++)*b1[*aj++]; 456 } 457 c[col*am + ridx[i]] += r1; 458 } 459 b1 += bm; 460 } 461 } else { 462 for (col=0; col<cn-4; col += 4){ /* over columns of C */ 463 colam = col*am; 464 for (i=0; i<am; i++) { /* over rows of C in those columns */ 465 r1 = r2 = r3 = r4 = 0.0; 466 n = a->i[i+1] - a->i[i]; 467 aj = a->j + a->i[i]; 468 aa = a->a + a->i[i]; 469 for (j=0; j<n; j++) { 470 r1 += (*aa)*b1[*aj]; 471 r2 += (*aa)*b2[*aj]; 472 r3 += (*aa)*b3[*aj]; 473 r4 += (*aa++)*b4[*aj++]; 474 } 475 c[colam + i] += r1; 476 c[colam + am + i] += r2; 477 c[colam + am2 + i] += r3; 478 c[colam + am3 + i] += r4; 479 } 480 b1 += bm4; 481 b2 += bm4; 482 b3 += bm4; 483 b4 += bm4; 484 } 485 for (;col<cn; col++){ /* over extra columns of C */ 486 for (i=0; i<am; i++) { /* over rows of C in those columns */ 487 r1 = 0.0; 488 n = a->i[i+1] - a->i[i]; 489 aj = a->j + a->i[i]; 490 aa = a->a + a->i[i]; 491 492 for (j=0; j<n; j++) { 493 r1 += (*aa++)*b1[*aj++]; 494 } 495 c[col*am + i] += r1; 496 } 497 b1 += bm; 498 } 499 } 500 ierr = PetscLogFlops(cn*2.0*a->nz);CHKERRQ(ierr); 501 ierr = MatRestoreArray(B,&b);CHKERRQ(ierr); 502 ierr = MatRestoreArray(C,&c);CHKERRQ(ierr); 503 PetscFunctionReturn(0); 504 } 505