Actual source code: matmatmult.c

  2: /*
  3:   Defines matrix-matrix product routines for pairs of SeqAIJ matrices
  4:           C = A * B
  5: */

  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*/

 15: PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
 16: {

 20:   if (scall == MAT_INITIAL_MATRIX){
 21:     MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);
 22:   }
 23:   MatMatMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);
 24:   return(0);
 25: }

 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;

 43:   /* Set up */
 44:   /* Allocate ci array, arrays for fill computation and */
 45:   /* free space for accumulating nonzero column info */
 46:   PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);
 47:   ci[0] = 0;
 48: 
 49:   /* create and initialize a linked list */
 50:   nlnk = bn+1;
 51:   PetscLLCreate(bn,bn,nlnk,lnk,lnkbt);

 53:   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
 54:   PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);
 55:   current_space = free_space;

 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:       PetscLLAdd(bnzj,bjj,bn,nlnk,lnk,lnkbt);
 70:       cnzi += nlnk;
 71:     }

 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:       PetscFreeSpaceGet(cnzi+current_space->total_array_size,&current_space);
 77:       nspacedouble++;
 78:     }

 80:     /* Copy data into free space, then initialize lnk */
 81:     PetscLLClean(bn,bn,cnzi,lnk,current_space->array,lnkbt);
 82:     current_space->array           += cnzi;
 83:     current_space->local_used      += cnzi;
 84:     current_space->local_remaining -= cnzi;

 86:     ci[i+1] = ci[i] + cnzi;
 87:   }

 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:   PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);
 93:   PetscFreeSpaceContiguous(&free_space,cj);
 94:   PetscLLDestroy(lnk,lnkbt);
 95: 
 96:   /* Allocate space for ca */
 97:   PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);
 98:   PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));
 99: 
100:   /* put together the new symbolic matrix */
101:   MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bn,ci,cj,ca,C);

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;

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;

119: #if defined(PETSC_USE_INFO)
120:   if (ci[am]) {
121:     PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);
122:     PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);
123:   } else {
124:     PetscInfo((*C),"Empty matrix product\n");
125:   }
126: #endif
127:   return(0);
128: }


133: PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
134: {
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;

146:   /* clean old values in C */
147:   PetscMemzero(ca,ci[cm]*sizeof(MatScalar));
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:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
172:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

174:   PetscLogFlops(flops);
175:   return(0);
176: }


181: PetscErrorCode MatMatMultTranspose_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) {

185:   if (scall == MAT_INITIAL_MATRIX){
186:     MatMatMultTransposeSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);
187:   }
188:   MatMatMultTransposeNumeric_SeqAIJ_SeqAIJ(A,B,*C);
189:   return(0);
190: }

194: PetscErrorCode MatMatMultTransposeSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
195: {
197:   Mat            At;
198:   PetscInt       *ati,*atj;

201:   /* create symbolic At */
202:   MatGetSymbolicTranspose_SeqAIJ(A,&ati,&atj);
203:   MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,A->cmap->n,A->rmap->n,ati,atj,PETSC_NULL,&At);

205:   /* get symbolic C=At*B */
206:   MatMatMultSymbolic_SeqAIJ_SeqAIJ(At,B,fill,C);

208:   /* clean up */
209:   MatDestroy(&At);
210:   MatRestoreSymbolicTranspose_SeqAIJ(A,&ati,&atj);
211:   return(0);
212: }

216: PetscErrorCode MatMatMultTransposeNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
217: {
219:   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data;
220:   PetscInt       am=A->rmap->n,anzi,*ai=a->i,*aj=a->j,*bi=b->i,*bj,bnzi,nextb;
221:   PetscInt       cm=C->rmap->n,*ci=c->i,*cj=c->j,crow,*cjj,i,j,k;
222:   PetscLogDouble flops=0.0;
223:   MatScalar      *aa=a->a,*ba,*ca=c->a,*caj;
224: 
226:   /* clear old values in C */
227:   PetscMemzero(ca,ci[cm]*sizeof(MatScalar));

229:   /* compute A^T*B using outer product (A^T)[:,i]*B[i,:] */
230:   for (i=0;i<am;i++) {
231:     bj   = b->j + bi[i];
232:     ba   = b->a + bi[i];
233:     bnzi = bi[i+1] - bi[i];
234:     anzi = ai[i+1] - ai[i];
235:     for (j=0; j<anzi; j++) {
236:       nextb = 0;
237:       crow  = *aj++;
238:       cjj   = cj + ci[crow];
239:       caj   = ca + ci[crow];
240:       /* perform sparse axpy operation.  Note cjj includes bj. */
241:       for (k=0; nextb<bnzi; k++) {
242:         if (cjj[k] == *(bj+nextb)) { /* ccol == bcol */
243:           caj[k] += (*aa)*(*(ba+nextb));
244:           nextb++;
245:         }
246:       }
247:       flops += 2*bnzi;
248:       aa++;
249:     }
250:   }

252:   /* Assemble the final matrix and clean up */
253:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
254:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
255:   PetscLogFlops(flops);
256:   return(0);
257: }

262: PetscErrorCode MatMatMult_SeqAIJ_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
263: {

267:   if (scall == MAT_INITIAL_MATRIX){
268:     MatMatMultSymbolic_SeqAIJ_SeqDense(A,B,fill,C);
269:   }
270:   MatMatMultNumeric_SeqAIJ_SeqDense(A,B,*C);
271:   return(0);
272: }

277: PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
278: {

282:   MatMatMultSymbolic_SeqDense_SeqDense(A,B,0.0,C);
283:   return(0);
284: }

288: PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
289: {
290:   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
292:   PetscScalar    *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4;
293:   MatScalar      *aa;
294:   PetscInt       cm=C->rmap->n, cn=B->cmap->n, bm=B->rmap->n, col, i,j,n,*aj, am = A->rmap->n;
295:   PetscInt       am2 = 2*am, am3 = 3*am,  bm4 = 4*bm,colam;

298:   if (!cm || !cn) return(0);
299:   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);
300:   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);
301:   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);
302:   MatGetArray(B,&b);
303:   MatGetArray(C,&c);
304:   b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;
305:   for (col=0; col<cn-4; col += 4){  /* over columns of C */
306:     colam = col*am;
307:     for (i=0; i<am; i++) {        /* over rows of C in those columns */
308:       r1 = r2 = r3 = r4 = 0.0;
309:       n   = a->i[i+1] - a->i[i];
310:       aj  = a->j + a->i[i];
311:       aa  = a->a + a->i[i];
312:       for (j=0; j<n; j++) {
313:         r1 += (*aa)*b1[*aj];
314:         r2 += (*aa)*b2[*aj];
315:         r3 += (*aa)*b3[*aj];
316:         r4 += (*aa++)*b4[*aj++];
317:       }
318:       c[colam + i]       = r1;
319:       c[colam + am + i]  = r2;
320:       c[colam + am2 + i] = r3;
321:       c[colam + am3 + i] = r4;
322:     }
323:     b1 += bm4;
324:     b2 += bm4;
325:     b3 += bm4;
326:     b4 += bm4;
327:   }
328:   for (;col<cn; col++){     /* over extra columns of C */
329:     for (i=0; i<am; i++) {  /* over rows of C in those columns */
330:       r1 = 0.0;
331:       n   = a->i[i+1] - a->i[i];
332:       aj  = a->j + a->i[i];
333:       aa  = a->a + a->i[i];

335:       for (j=0; j<n; j++) {
336:         r1 += (*aa++)*b1[*aj++];
337:       }
338:       c[col*am + i]     = r1;
339:     }
340:     b1 += bm;
341:   }
342:   PetscLogFlops(cn*(2.0*a->nz));
343:   MatRestoreArray(B,&b);
344:   MatRestoreArray(C,&c);
345:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
346:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
347:   return(0);
348: }

350: /*
351:    Note very similar to MatMult_SeqAIJ(), should generate both codes from same base
352: */
355: PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
356: {
357:   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
359:   PetscScalar    *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4;
360:   MatScalar      *aa;
361:   PetscInt       cm=C->rmap->n, cn=B->cmap->n, bm=B->rmap->n, col, i,j,n,*aj, am = A->rmap->n,*ii,arm;
362:   PetscInt       am2 = 2*am, am3 = 3*am,  bm4 = 4*bm,colam,*ridx;

365:   if (!cm || !cn) return(0);
366:   MatGetArray(B,&b);
367:   MatGetArray(C,&c);
368:   b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;

370:   if (a->compressedrow.use){ /* use compressed row format */
371:     for (col=0; col<cn-4; col += 4){  /* over columns of C */
372:       colam = col*am;
373:       arm   = a->compressedrow.nrows;
374:       ii    = a->compressedrow.i;
375:       ridx  = a->compressedrow.rindex;
376:       for (i=0; i<arm; i++) {        /* over rows of C in those columns */
377:         r1 = r2 = r3 = r4 = 0.0;
378:         n   = ii[i+1] - ii[i];
379:         aj  = a->j + ii[i];
380:         aa  = a->a + ii[i];
381:         for (j=0; j<n; j++) {
382:           r1 += (*aa)*b1[*aj];
383:           r2 += (*aa)*b2[*aj];
384:           r3 += (*aa)*b3[*aj];
385:           r4 += (*aa++)*b4[*aj++];
386:         }
387:         c[colam       + ridx[i]] += r1;
388:         c[colam + am  + ridx[i]] += r2;
389:         c[colam + am2 + ridx[i]] += r3;
390:         c[colam + am3 + ridx[i]] += r4;
391:       }
392:       b1 += bm4;
393:       b2 += bm4;
394:       b3 += bm4;
395:       b4 += bm4;
396:     }
397:     for (;col<cn; col++){     /* over extra columns of C */
398:       colam = col*am;
399:       arm   = a->compressedrow.nrows;
400:       ii    = a->compressedrow.i;
401:       ridx  = a->compressedrow.rindex;
402:       for (i=0; i<arm; i++) {  /* over rows of C in those columns */
403:         r1 = 0.0;
404:         n   = ii[i+1] - ii[i];
405:         aj  = a->j + ii[i];
406:         aa  = a->a + ii[i];

408:         for (j=0; j<n; j++) {
409:           r1 += (*aa++)*b1[*aj++];
410:         }
411:         c[col*am + ridx[i]] += r1;
412:       }
413:       b1 += bm;
414:     }
415:   } else {
416:     for (col=0; col<cn-4; col += 4){  /* over columns of C */
417:       colam = col*am;
418:       for (i=0; i<am; i++) {        /* over rows of C in those columns */
419:         r1 = r2 = r3 = r4 = 0.0;
420:         n   = a->i[i+1] - a->i[i];
421:         aj  = a->j + a->i[i];
422:         aa  = a->a + a->i[i];
423:         for (j=0; j<n; j++) {
424:           r1 += (*aa)*b1[*aj];
425:           r2 += (*aa)*b2[*aj];
426:           r3 += (*aa)*b3[*aj];
427:           r4 += (*aa++)*b4[*aj++];
428:         }
429:         c[colam + i]       += r1;
430:         c[colam + am + i]  += r2;
431:         c[colam + am2 + i] += r3;
432:         c[colam + am3 + i] += r4;
433:       }
434:       b1 += bm4;
435:       b2 += bm4;
436:       b3 += bm4;
437:       b4 += bm4;
438:     }
439:     for (;col<cn; col++){     /* over extra columns of C */
440:       for (i=0; i<am; i++) {  /* over rows of C in those columns */
441:         r1 = 0.0;
442:         n   = a->i[i+1] - a->i[i];
443:         aj  = a->j + a->i[i];
444:         aa  = a->a + a->i[i];

446:         for (j=0; j<n; j++) {
447:           r1 += (*aa++)*b1[*aj++];
448:         }
449:         c[col*am + i]     += r1;
450:       }
451:       b1 += bm;
452:     }
453:   }
454:   PetscLogFlops(cn*2.0*a->nz);
455:   MatRestoreArray(B,&b);
456:   MatRestoreArray(C,&c);
457:   return(0);
458: }