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
16: PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
17: {
21: if (scall == MAT_INITIAL_MATRIX){
22: MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);
23: }
24: MatMatMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);
25: return(0);
26: }
31: PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
32: {
33: PetscErrorCode ierr;
34: PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
35: Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
36: PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci,*cj;
37: PetscInt am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
38: PetscInt i,j,anzi,brow,bnzj,cnzi,nlnk,*lnk,nspacedouble=0;
39: MatScalar *ca;
40: PetscBT lnkbt;
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,¤t_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: #if defined(PETSC_USE_INFO)
111: if (ci[am] != 0) {
112: PetscReal afill = ((PetscReal)ci[am])/(ai[am]+bi[bm]);
113: if (afill < 1.0) afill = 1.0;
114: PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);
115: PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);
116: } else {
117: PetscInfo((*C),"Empty matrix product\n");
118: }
119: #endif
120: return(0);
121: }
126: PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
127: {
129: PetscLogDouble flops=0.0;
130: Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data;
131: Mat_SeqAIJ *b = (Mat_SeqAIJ *)B->data;
132: Mat_SeqAIJ *c = (Mat_SeqAIJ *)C->data;
133: PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j;
134: PetscInt am=A->rmap->N,cm=C->rmap->N;
135: PetscInt i,j,k,anzi,bnzi,cnzi,brow,nextb;
136: MatScalar *aa=a->a,*ba=b->a,*baj,*ca=c->a;
139: /* clean old values in C */
140: PetscMemzero(ca,ci[cm]*sizeof(MatScalar));
141: /* Traverse A row-wise. */
142: /* Build the ith row in C by summing over nonzero columns in A, */
143: /* the rows of B corresponding to nonzeros of A. */
144: for (i=0;i<am;i++) {
145: anzi = ai[i+1] - ai[i];
146: for (j=0;j<anzi;j++) {
147: brow = *aj++;
148: bnzi = bi[brow+1] - bi[brow];
149: bjj = bj + bi[brow];
150: baj = ba + bi[brow];
151: nextb = 0;
152: for (k=0; nextb<bnzi; k++) {
153: if (cj[k] == bjj[nextb]){ /* ccol == bcol */
154: ca[k] += (*aa)*baj[nextb++];
155: }
156: }
157: flops += 2*bnzi;
158: aa++;
159: }
160: cnzi = ci[i+1] - ci[i];
161: ca += cnzi;
162: cj += cnzi;
163: }
164: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
165: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
167: PetscLogFlops(flops);
168: return(0);
169: }
174: PetscErrorCode MatMatMultTranspose_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) {
178: if (scall == MAT_INITIAL_MATRIX){
179: MatMatMultTransposeSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);
180: }
181: MatMatMultTransposeNumeric_SeqAIJ_SeqAIJ(A,B,*C);
182: return(0);
183: }
187: PetscErrorCode MatMatMultTransposeSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
188: {
190: Mat At;
191: PetscInt *ati,*atj;
194: /* create symbolic At */
195: MatGetSymbolicTranspose_SeqAIJ(A,&ati,&atj);
196: MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,A->cmap->n,A->rmap->n,ati,atj,PETSC_NULL,&At);
198: /* get symbolic C=At*B */
199: MatMatMultSymbolic_SeqAIJ_SeqAIJ(At,B,fill,C);
201: /* clean up */
202: MatDestroy(At);
203: MatRestoreSymbolicTranspose_SeqAIJ(A,&ati,&atj);
204:
205: return(0);
206: }
210: PetscErrorCode MatMatMultTransposeNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
211: {
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:
220: /* clear old values in C */
221: PetscMemzero(ca,ci[cm]*sizeof(MatScalar));
223: /* compute A^T*B using outer product (A^T)[:,i]*B[i,:] */
224: for (i=0;i<am;i++) {
225: bj = b->j + bi[i];
226: ba = b->a + bi[i];
227: bnzi = bi[i+1] - bi[i];
228: anzi = ai[i+1] - ai[i];
229: for (j=0; j<anzi; j++) {
230: nextb = 0;
231: crow = *aj++;
232: cjj = cj + ci[crow];
233: caj = ca + ci[crow];
234: /* perform sparse axpy operation. Note cjj includes bj. */
235: for (k=0; nextb<bnzi; k++) {
236: if (cjj[k] == *(bj+nextb)) { /* ccol == bcol */
237: caj[k] += (*aa)*(*(ba+nextb));
238: nextb++;
239: }
240: }
241: flops += 2*bnzi;
242: aa++;
243: }
244: }
246: /* Assemble the final matrix and clean up */
247: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
248: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
249: PetscLogFlops(flops);
250: return(0);
251: }
256: PetscErrorCode MatMatMult_SeqAIJ_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
257: {
261: if (scall == MAT_INITIAL_MATRIX){
262: MatMatMultSymbolic_SeqAIJ_SeqDense(A,B,fill,C);
263: }
264: MatMatMultNumeric_SeqAIJ_SeqDense(A,B,*C);
265: return(0);
266: }
271: PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
272: {
276: MatMatMultSymbolic_SeqDense_SeqDense(A,B,0.0,C);
277: return(0);
278: }
282: PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
283: {
284: Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data;
286: PetscScalar *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4;
287: MatScalar *aa;
288: PetscInt cm=C->rmap->n, cn=B->cmap->n, bm=B->rmap->n, col, i,j,n,*aj, am = A->rmap->n;
289: PetscInt am2 = 2*am, am3 = 3*am, bm4 = 4*bm,colam;
292: if (!cm || !cn) return(0);
293: 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);
294: 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);
295: 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);
296: MatGetArray(B,&b);
297: MatGetArray(C,&c);
298: b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;
299: for (col=0; col<cn-4; col += 4){ /* over columns of C */
300: colam = col*am;
301: for (i=0; i<am; i++) { /* over rows of C in those columns */
302: r1 = r2 = r3 = r4 = 0.0;
303: n = a->i[i+1] - a->i[i];
304: aj = a->j + a->i[i];
305: aa = a->a + a->i[i];
306: for (j=0; j<n; j++) {
307: r1 += (*aa)*b1[*aj];
308: r2 += (*aa)*b2[*aj];
309: r3 += (*aa)*b3[*aj];
310: r4 += (*aa++)*b4[*aj++];
311: }
312: c[colam + i] = r1;
313: c[colam + am + i] = r2;
314: c[colam + am2 + i] = r3;
315: c[colam + am3 + i] = r4;
316: }
317: b1 += bm4;
318: b2 += bm4;
319: b3 += bm4;
320: b4 += bm4;
321: }
322: for (;col<cn; col++){ /* over extra columns of C */
323: for (i=0; i<am; i++) { /* over rows of C in those columns */
324: r1 = 0.0;
325: n = a->i[i+1] - a->i[i];
326: aj = a->j + a->i[i];
327: aa = a->a + a->i[i];
329: for (j=0; j<n; j++) {
330: r1 += (*aa++)*b1[*aj++];
331: }
332: c[col*am + i] = r1;
333: }
334: b1 += bm;
335: }
336: PetscLogFlops(cn*(2.0*a->nz));
337: MatRestoreArray(B,&b);
338: MatRestoreArray(C,&c);
339: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
340: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
341: return(0);
342: }
344: /*
345: Note very similar to MatMult_SeqAIJ(), should generate both codes from same base
346: */
349: PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
350: {
351: Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data;
353: PetscScalar *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4;
354: MatScalar *aa;
355: PetscInt cm=C->rmap->n, cn=B->cmap->n, bm=B->rmap->n, col, i,j,n,*aj, am = A->rmap->n,*ii,arm;
356: PetscInt am2 = 2*am, am3 = 3*am, bm4 = 4*bm,colam,*ridx;
359: if (!cm || !cn) return(0);
360: MatGetArray(B,&b);
361: MatGetArray(C,&c);
362: b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;
364: if (a->compressedrow.use){ /* use compressed row format */
365: for (col=0; col<cn-4; col += 4){ /* over columns of C */
366: colam = col*am;
367: arm = a->compressedrow.nrows;
368: ii = a->compressedrow.i;
369: ridx = a->compressedrow.rindex;
370: for (i=0; i<arm; i++) { /* over rows of C in those columns */
371: r1 = r2 = r3 = r4 = 0.0;
372: n = ii[i+1] - ii[i];
373: aj = a->j + ii[i];
374: aa = a->a + ii[i];
375: for (j=0; j<n; j++) {
376: r1 += (*aa)*b1[*aj];
377: r2 += (*aa)*b2[*aj];
378: r3 += (*aa)*b3[*aj];
379: r4 += (*aa++)*b4[*aj++];
380: }
381: c[colam + ridx[i]] += r1;
382: c[colam + am + ridx[i]] += r2;
383: c[colam + am2 + ridx[i]] += r3;
384: c[colam + am3 + ridx[i]] += r4;
385: }
386: b1 += bm4;
387: b2 += bm4;
388: b3 += bm4;
389: b4 += bm4;
390: }
391: for (;col<cn; col++){ /* over extra columns of C */
392: colam = col*am;
393: arm = a->compressedrow.nrows;
394: ii = a->compressedrow.i;
395: ridx = a->compressedrow.rindex;
396: for (i=0; i<arm; i++) { /* over rows of C in those columns */
397: r1 = 0.0;
398: n = ii[i+1] - ii[i];
399: aj = a->j + ii[i];
400: aa = a->a + ii[i];
402: for (j=0; j<n; j++) {
403: r1 += (*aa++)*b1[*aj++];
404: }
405: c[col*am + ridx[i]] += r1;
406: }
407: b1 += bm;
408: }
409: } else {
410: for (col=0; col<cn-4; col += 4){ /* over columns of C */
411: colam = col*am;
412: for (i=0; i<am; i++) { /* over rows of C in those columns */
413: r1 = r2 = r3 = r4 = 0.0;
414: n = a->i[i+1] - a->i[i];
415: aj = a->j + a->i[i];
416: aa = a->a + a->i[i];
417: for (j=0; j<n; j++) {
418: r1 += (*aa)*b1[*aj];
419: r2 += (*aa)*b2[*aj];
420: r3 += (*aa)*b3[*aj];
421: r4 += (*aa++)*b4[*aj++];
422: }
423: c[colam + i] += r1;
424: c[colam + am + i] += r2;
425: c[colam + am2 + i] += r3;
426: c[colam + am3 + i] += r4;
427: }
428: b1 += bm4;
429: b2 += bm4;
430: b3 += bm4;
431: b4 += bm4;
432: }
433: for (;col<cn; col++){ /* over extra columns of C */
434: for (i=0; i<am; i++) { /* over rows of C in those columns */
435: r1 = 0.0;
436: n = a->i[i+1] - a->i[i];
437: aj = a->j + a->i[i];
438: aa = a->a + a->i[i];
440: for (j=0; j<n; j++) {
441: r1 += (*aa++)*b1[*aj++];
442: }
443: c[col*am + i] += r1;
444: }
445: b1 += bm;
446: }
447: }
448: PetscLogFlops(cn*2.0*a->nz);
449: MatRestoreArray(B,&b);
450: MatRestoreArray(C,&c);
451: return(0);
452: }