Actual source code: matmatmult.c

  1: #define PETSCMAT_DLL

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

 8:  #include src/mat/impls/aij/seq/aij.h
 9:  #include src/mat/utils/freespace.h
 10:  #include petscbt.h
 11:  #include src/mat/impls/dense/seq/dense.h

 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;

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

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

 56:   /* Determine symbolic info for each row of the product: */
 57:   for (i=0;i<am;i++) {
 58:     anzi = ai[i+1] - ai[i];
 59:     cnzi = 0;
 60:     j    = anzi;
 61:     aj   = a->j + ai[i];
 62:     while (j){/* assume cols are almost in increasing order, starting from its end saves computation */
 63:       j--;
 64:       brow = *(aj + j);
 65:       bnzj = bi[brow+1] - bi[brow];
 66:       bjj  = bj + bi[brow];
 67:       /* add non-zero cols of B into the sorted linked list lnk */
 68:       PetscLLAdd(bnzj,bjj,bn,nlnk,lnk,lnkbt);
 69:       cnzi += nlnk;
 70:     }

 72:     /* If free space is not available, make more free space */
 73:     /* Double the amount of total space in the list */
 74:     if (current_space->local_remaining<cnzi) {
 75:       PetscFreeSpaceGet(current_space->total_array_size,&current_space);
 76:       nspacedouble++;
 77:     }

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

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

 88:   /* Column indices are in the list of free space */
 89:   /* Allocate space for cj, initialize cj, and */
 90:   /* destroy list of free space and other temporary array(s) */
 91:   PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);
 92:   PetscFreeSpaceContiguous(&free_space,cj);
 93:   PetscLLDestroy(lnk,lnkbt);
 94: 
 95:   /* Allocate space for ca */
 96:   PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);
 97:   PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));
 98: 
 99:   /* put together the new symbolic matrix */
100:   MatCreateSeqAIJWithArrays(A->comm,am,bn,ci,cj,ca,C);

102:   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
103:   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
104:   c = (Mat_SeqAIJ *)((*C)->data);
105:   c->free_a   = PETSC_TRUE;
106:   c->free_ij  = PETSC_TRUE;
107:   c->nonew    = 0;

109:   if (nspacedouble){
110:     PetscInfo5((*C),"nspacedouble:%D, nnz(A):%D, nnz(B):%D, fill:%G, nnz(C):%D\n",nspacedouble,ai[am],bi[bm],fill,ci[am]);
111:   }
112:   return(0);
113: }


118: PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
119: {
121:   PetscInt       flops=0;
122:   Mat_SeqAIJ     *a = (Mat_SeqAIJ *)A->data;
123:   Mat_SeqAIJ     *b = (Mat_SeqAIJ *)B->data;
124:   Mat_SeqAIJ     *c = (Mat_SeqAIJ *)C->data;
125:   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j;
126:   PetscInt       am=A->rmap.N,cm=C->rmap.N;
127:   PetscInt       i,j,k,anzi,bnzi,cnzi,brow,nextb;
128:   MatScalar      *aa=a->a,*ba=b->a,*baj,*ca=c->a;

131:   /* clean old values in C */
132:   PetscMemzero(ca,ci[cm]*sizeof(MatScalar));
133:   /* Traverse A row-wise. */
134:   /* Build the ith row in C by summing over nonzero columns in A, */
135:   /* the rows of B corresponding to nonzeros of A. */
136:   for (i=0;i<am;i++) {
137:     anzi = ai[i+1] - ai[i];
138:     for (j=0;j<anzi;j++) {
139:       brow = *aj++;
140:       bnzi = bi[brow+1] - bi[brow];
141:       bjj  = bj + bi[brow];
142:       baj  = ba + bi[brow];
143:       nextb = 0;
144:       for (k=0; nextb<bnzi; k++) {
145:         if (cj[k] == bjj[nextb]){ /* ccol == bcol */
146:           ca[k] += (*aa)*baj[nextb++];
147:         }
148:       }
149:       flops += 2*bnzi;
150:       aa++;
151:     }
152:     cnzi = ci[i+1] - ci[i];
153:     ca += cnzi;
154:     cj += cnzi;
155:   }
156:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
157:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

159:   PetscLogFlops(flops);
160:   return(0);
161: }


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

170:   if (scall == MAT_INITIAL_MATRIX){
171:     MatMatMultTransposeSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);
172:   }
173:   MatMatMultTransposeNumeric_SeqAIJ_SeqAIJ(A,B,*C);
174:   return(0);
175: }

179: PetscErrorCode MatMatMultTransposeSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
180: {
182:   Mat            At;
183:   PetscInt       *ati,*atj;

186:   /* create symbolic At */
187:   MatGetSymbolicTranspose_SeqAIJ(A,&ati,&atj);
188:   MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,A->cmap.n,A->rmap.n,ati,atj,PETSC_NULL,&At);

190:   /* get symbolic C=At*B */
191:   MatMatMultSymbolic_SeqAIJ_SeqAIJ(At,B,fill,C);

193:   /* clean up */
194:   MatDestroy(At);
195:   MatRestoreSymbolicTranspose_SeqAIJ(A,&ati,&atj);
196: 
197:   return(0);
198: }

202: PetscErrorCode MatMatMultTransposeNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
203: {
205:   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data;
206:   PetscInt       am=A->rmap.n,anzi,*ai=a->i,*aj=a->j,*bi=b->i,*bj,bnzi,nextb;
207:   PetscInt       cm=C->rmap.n,*ci=c->i,*cj=c->j,crow,*cjj,i,j,k,flops=0;
208:   MatScalar      *aa=a->a,*ba,*ca=c->a,*caj;
209: 
211:   /* clear old values in C */
212:   PetscMemzero(ca,ci[cm]*sizeof(MatScalar));

214:   /* compute A^T*B using outer product (A^T)[:,i]*B[i,:] */
215:   for (i=0;i<am;i++) {
216:     bj   = b->j + bi[i];
217:     ba   = b->a + bi[i];
218:     bnzi = bi[i+1] - bi[i];
219:     anzi = ai[i+1] - ai[i];
220:     for (j=0; j<anzi; j++) {
221:       nextb = 0;
222:       crow  = *aj++;
223:       cjj   = cj + ci[crow];
224:       caj   = ca + ci[crow];
225:       /* perform sparse axpy operation.  Note cjj includes bj. */
226:       for (k=0; nextb<bnzi; k++) {
227:         if (cjj[k] == *(bj+nextb)) { /* ccol == bcol */
228:           caj[k] += (*aa)*(*(ba+nextb));
229:           nextb++;
230:         }
231:       }
232:       flops += 2*bnzi;
233:       aa++;
234:     }
235:   }

237:   /* Assemble the final matrix and clean up */
238:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
239:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
240:   PetscLogFlops(flops);
241:   return(0);
242: }

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

251:   if (scall == MAT_INITIAL_MATRIX){
252:     MatMatMultSymbolic_SeqAIJ_SeqDense(A,B,fill,C);
253:   }
254:   MatMatMultNumeric_SeqAIJ_SeqDense(A,B,*C);
255:   return(0);
256: }

260: PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
261: {

265:   MatMatMultSymbolic_SeqDense_SeqDense(A,B,0.0,C);
266:   return(0);
267: }

271: PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
272: {
273:   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
275:   PetscScalar    *b,*c,r1,r2,r3,r4,*aa,*b1,*b2,*b3,*b4;
276:   PetscInt       cm=C->rmap.n, cn=B->cmap.n, bm=B->rmap.n, col, i,j,n,*aj, am = A->rmap.n;
277:   PetscInt       am2 = 2*am, am3 = 3*am,  bm4 = 4*bm,colam;

280:   if (!cm || !cn) return(0);
281:   if (bm != A->cmap.n) SETERRQ2(PETSC_ERR_ARG_SIZ,"Number columns in A %D not equal rows in B %D\n",A->cmap.n,bm);
282:   if (A->rmap.n != C->rmap.n) SETERRQ2(PETSC_ERR_ARG_SIZ,"Number rows in C %D not equal rows in A %D\n",C->rmap.n,A->rmap.n);
283:   if (B->cmap.n != C->cmap.n) SETERRQ2(PETSC_ERR_ARG_SIZ,"Number columns in B %D not equal columns in C %D\n",B->cmap.n,C->cmap.n);
284:   MatGetArray(B,&b);
285:   MatGetArray(C,&c);
286:   b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;
287:   for (col=0; col<cn-4; col += 4){  /* over columns of C */
288:     colam = col*am;
289:     for (i=0; i<am; i++) {        /* over rows of C in those columns */
290:       r1 = r2 = r3 = r4 = 0.0;
291:       n   = a->i[i+1] - a->i[i];
292:       aj  = a->j + a->i[i];
293:       aa  = a->a + a->i[i];
294:       for (j=0; j<n; j++) {
295:         r1 += (*aa)*b1[*aj];
296:         r2 += (*aa)*b2[*aj];
297:         r3 += (*aa)*b3[*aj];
298:         r4 += (*aa++)*b4[*aj++];
299:       }
300:       c[colam + i]       = r1;
301:       c[colam + am + i]  = r2;
302:       c[colam + am2 + i] = r3;
303:       c[colam + am3 + i] = r4;
304:     }
305:     b1 += bm4;
306:     b2 += bm4;
307:     b3 += bm4;
308:     b4 += bm4;
309:   }
310:   for (;col<cn; col++){     /* over extra columns of C */
311:     for (i=0; i<am; i++) {  /* over rows of C in those columns */
312:       r1 = 0.0;
313:       n   = a->i[i+1] - a->i[i];
314:       aj  = a->j + a->i[i];
315:       aa  = a->a + a->i[i];

317:       for (j=0; j<n; j++) {
318:         r1 += (*aa++)*b1[*aj++];
319:       }
320:       c[col*am + i]     = r1;
321:     }
322:     b1 += bm;
323:   }
324:   PetscLogFlops(cn*(2*a->nz - A->rmap.n));
325:   MatRestoreArray(B,&b);
326:   MatRestoreArray(C,&c);
327:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
328:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
329:   return(0);
330: }

332: /*
333:    Note very similar to MatMult_SeqAIJ(), should generate both codes from same base
334: */
337: PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
338: {
339:   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
341:   PetscScalar    *b,*c,r1,r2,r3,r4,*aa,*b1,*b2,*b3,*b4;
342:   PetscInt       cm=C->rmap.n, cn=B->cmap.n, bm=B->rmap.n, col, i,j,n,*aj, am = A->rmap.n,*ii,arm;
343:   PetscInt       am2 = 2*am, am3 = 3*am,  bm4 = 4*bm,colam,*ridx;

346:   if (!cm || !cn) return(0);
347:   MatGetArray(B,&b);
348:   MatGetArray(C,&c);
349:   b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;

351:   if (a->compressedrow.use){ /* use compressed row format */
352:     for (col=0; col<cn-4; col += 4){  /* over columns of C */
353:       colam = col*am;
354:       arm   = a->compressedrow.nrows;
355:       ii    = a->compressedrow.i;
356:       ridx  = a->compressedrow.rindex;
357:       for (i=0; i<arm; i++) {        /* over rows of C in those columns */
358:         r1 = r2 = r3 = r4 = 0.0;
359:         n   = ii[i+1] - ii[i];
360:         aj  = a->j + ii[i];
361:         aa  = a->a + ii[i];
362:         for (j=0; j<n; j++) {
363:           r1 += (*aa)*b1[*aj];
364:           r2 += (*aa)*b2[*aj];
365:           r3 += (*aa)*b3[*aj];
366:           r4 += (*aa++)*b4[*aj++];
367:         }
368:         c[colam       + ridx[i]] += r1;
369:         c[colam + am  + ridx[i]] += r2;
370:         c[colam + am2 + ridx[i]] += r3;
371:         c[colam + am3 + ridx[i]] += r4;
372:       }
373:       b1 += bm4;
374:       b2 += bm4;
375:       b3 += bm4;
376:       b4 += bm4;
377:     }
378:     for (;col<cn; col++){     /* over extra columns of C */
379:       colam = col*am;
380:       arm   = a->compressedrow.nrows;
381:       ii    = a->compressedrow.i;
382:       ridx  = a->compressedrow.rindex;
383:       for (i=0; i<arm; i++) {  /* over rows of C in those columns */
384:         r1 = 0.0;
385:         n   = ii[i+1] - ii[i];
386:         aj  = a->j + ii[i];
387:         aa  = a->a + ii[i];

389:         for (j=0; j<n; j++) {
390:           r1 += (*aa++)*b1[*aj++];
391:         }
392:         c[col*am + ridx[i]] += r1;
393:       }
394:       b1 += bm;
395:     }
396:   } else {
397:     for (col=0; col<cn-4; col += 4){  /* over columns of C */
398:       colam = col*am;
399:       for (i=0; i<am; i++) {        /* over rows of C in those columns */
400:         r1 = r2 = r3 = r4 = 0.0;
401:         n   = a->i[i+1] - a->i[i];
402:         aj  = a->j + a->i[i];
403:         aa  = a->a + a->i[i];
404:         for (j=0; j<n; j++) {
405:           r1 += (*aa)*b1[*aj];
406:           r2 += (*aa)*b2[*aj];
407:           r3 += (*aa)*b3[*aj];
408:           r4 += (*aa++)*b4[*aj++];
409:         }
410:         c[colam + i]       += r1;
411:         c[colam + am + i]  += r2;
412:         c[colam + am2 + i] += r3;
413:         c[colam + am3 + i] += r4;
414:       }
415:       b1 += bm4;
416:       b2 += bm4;
417:       b3 += bm4;
418:       b4 += bm4;
419:     }
420:     for (;col<cn; col++){     /* over extra columns of C */
421:       for (i=0; i<am; i++) {  /* over rows of C in those columns */
422:         r1 = 0.0;
423:         n   = a->i[i+1] - a->i[i];
424:         aj  = a->j + a->i[i];
425:         aa  = a->a + a->i[i];

427:         for (j=0; j<n; j++) {
428:           r1 += (*aa++)*b1[*aj++];
429:         }
430:         c[col*am + i]     += r1;
431:       }
432:       b1 += bm;
433:     }
434:   }
435:   PetscLogFlops(cn*2*a->nz);
436:   MatRestoreArray(B,&b);
437:   MatRestoreArray(C,&c);
438:   return(0);
439: }