xref: /petsc/src/mat/impls/aij/seq/matmatmult.c (revision 5d0c19d75c660d4fec594a5399ec8d8ba29c54a8)
1 #define PETSCMAT_DLL
2 
3 /*
4   Defines matrix-matrix product routines for pairs of SeqAIJ matrices
5           C = A * B
6 */
7 
8 #include "src/mat/impls/aij/seq/aij.h" /*I "petscmat.h" I*/
9 #include "src/mat/utils/freespace.h"
10 #include "petscbt.h"
11 #include "src/mat/impls/dense/seq/dense.h" /*I "petscmat.h" I*/
12 
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 
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 
41   PetscFunctionBegin;
42   /* Set up */
43   /* Allocate ci array, arrays for fill computation and */
44   /* free space for accumulating nonzero column info */
45   ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
46   ci[0] = 0;
47 
48   /* create and initialize a linked list */
49   nlnk = bn+1;
50   ierr = PetscLLCreate(bn,bn,nlnk,lnk,lnkbt);CHKERRQ(ierr);
51 
52   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
53   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr);
54   current_space = free_space;
55 
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       ierr = PetscLLAdd(bnzj,bjj,bn,nlnk,lnk,lnkbt);CHKERRQ(ierr);
69       cnzi += nlnk;
70     }
71 
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       ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,&current_space);CHKERRQ(ierr);
76       nspacedouble++;
77     }
78 
79     /* Copy data into free space, then initialize lnk */
80     ierr = PetscLLClean(bn,bn,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
81     current_space->array           += cnzi;
82     current_space->local_used      += cnzi;
83     current_space->local_remaining -= cnzi;
84 
85     ci[i+1] = ci[i] + cnzi;
86   }
87 
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   ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
92   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
93   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
94 
95   /* Allocate space for ca */
96   ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
97   ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr);
98 
99   /* put together the new symbolic matrix */
100   ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bn,ci,cj,ca,C);CHKERRQ(ierr);
101 
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;
108 
109 #if defined(PETSC_USE_INFO)
110   if (ci[am] != 0) {
111     PetscReal afill = ((PetscReal)ci[am])/(ai[am]+bi[bm]);
112     if (afill < 1.0) afill = 1.0;
113     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
114     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
115   } else {
116     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
117   }
118 #endif
119   PetscFunctionReturn(0);
120 }
121 
122 
123 #undef __FUNCT__
124 #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqAIJ"
125 PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
126 {
127   PetscErrorCode ierr;
128   PetscLogDouble flops=0.0;
129   Mat_SeqAIJ     *a = (Mat_SeqAIJ *)A->data;
130   Mat_SeqAIJ     *b = (Mat_SeqAIJ *)B->data;
131   Mat_SeqAIJ     *c = (Mat_SeqAIJ *)C->data;
132   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j;
133   PetscInt       am=A->rmap->N,cm=C->rmap->N;
134   PetscInt       i,j,k,anzi,bnzi,cnzi,brow,nextb;
135   MatScalar      *aa=a->a,*ba=b->a,*baj,*ca=c->a;
136 
137   PetscFunctionBegin;
138   /* clean old values in C */
139   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
140   /* Traverse A row-wise. */
141   /* Build the ith row in C by summing over nonzero columns in A, */
142   /* the rows of B corresponding to nonzeros of A. */
143   for (i=0;i<am;i++) {
144     anzi = ai[i+1] - ai[i];
145     for (j=0;j<anzi;j++) {
146       brow = *aj++;
147       bnzi = bi[brow+1] - bi[brow];
148       bjj  = bj + bi[brow];
149       baj  = ba + bi[brow];
150       nextb = 0;
151       for (k=0; nextb<bnzi; k++) {
152         if (cj[k] == bjj[nextb]){ /* ccol == bcol */
153           ca[k] += (*aa)*baj[nextb++];
154         }
155       }
156       flops += 2*bnzi;
157       aa++;
158     }
159     cnzi = ci[i+1] - ci[i];
160     ca += cnzi;
161     cj += cnzi;
162   }
163   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
164   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
165 
166   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
167   PetscFunctionReturn(0);
168 }
169 
170 
171 #undef __FUNCT__
172 #define __FUNCT__ "MatMatMultTranspose_SeqAIJ_SeqAIJ"
173 PetscErrorCode MatMatMultTranspose_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) {
174   PetscErrorCode ierr;
175 
176   PetscFunctionBegin;
177   if (scall == MAT_INITIAL_MATRIX){
178     ierr = MatMatMultTransposeSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
179   }
180   ierr = MatMatMultTransposeNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr);
181   PetscFunctionReturn(0);
182 }
183 
184 #undef __FUNCT__
185 #define __FUNCT__ "MatMatMultTransposeSymbolic_SeqAIJ_SeqAIJ"
186 PetscErrorCode MatMatMultTransposeSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
187 {
188   PetscErrorCode ierr;
189   Mat            At;
190   PetscInt       *ati,*atj;
191 
192   PetscFunctionBegin;
193   /* create symbolic At */
194   ierr = MatGetSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr);
195   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,A->cmap->n,A->rmap->n,ati,atj,PETSC_NULL,&At);CHKERRQ(ierr);
196 
197   /* get symbolic C=At*B */
198   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(At,B,fill,C);CHKERRQ(ierr);
199 
200   /* clean up */
201   ierr = MatDestroy(At);CHKERRQ(ierr);
202   ierr = MatRestoreSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr);
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   PetscErrorCode ierr;
212   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data;
213   PetscInt       am=A->rmap->n,anzi,*ai=a->i,*aj=a->j,*bi=b->i,*bj,bnzi,nextb;
214   PetscInt       cm=C->rmap->n,*ci=c->i,*cj=c->j,crow,*cjj,i,j,k;
215   PetscLogDouble flops=0.0;
216   MatScalar      *aa=a->a,*ba,*ca=c->a,*caj;
217 
218   PetscFunctionBegin;
219   /* clear old values in C */
220   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
221 
222   /* compute A^T*B using outer product (A^T)[:,i]*B[i,:] */
223   for (i=0;i<am;i++) {
224     bj   = b->j + bi[i];
225     ba   = b->a + bi[i];
226     bnzi = bi[i+1] - bi[i];
227     anzi = ai[i+1] - ai[i];
228     for (j=0; j<anzi; j++) {
229       nextb = 0;
230       crow  = *aj++;
231       cjj   = cj + ci[crow];
232       caj   = ca + ci[crow];
233       /* perform sparse axpy operation.  Note cjj includes bj. */
234       for (k=0; nextb<bnzi; k++) {
235         if (cjj[k] == *(bj+nextb)) { /* ccol == bcol */
236           caj[k] += (*aa)*(*(ba+nextb));
237           nextb++;
238         }
239       }
240       flops += 2*bnzi;
241       aa++;
242     }
243   }
244 
245   /* Assemble the final matrix and clean up */
246   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
247   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
248   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
249   PetscFunctionReturn(0);
250 }
251 
252 #undef __FUNCT__
253 #define __FUNCT__ "MatMatMult_SeqAIJ_SeqDense"
254 PetscErrorCode MatMatMult_SeqAIJ_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
255 {
256   PetscErrorCode ierr;
257 
258   PetscFunctionBegin;
259   if (scall == MAT_INITIAL_MATRIX){
260     ierr = MatMatMultSymbolic_SeqAIJ_SeqDense(A,B,fill,C);CHKERRQ(ierr);
261   }
262   ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,B,*C);CHKERRQ(ierr);
263   PetscFunctionReturn(0);
264 }
265 
266 #undef __FUNCT__
267 #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqDense"
268 PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
269 {
270   PetscErrorCode ierr;
271 
272   PetscFunctionBegin;
273   ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,0.0,C);
274   PetscFunctionReturn(0);
275 }
276 
277 #undef __FUNCT__
278 #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqDense"
279 PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
280 {
281   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
282   PetscErrorCode ierr;
283   PetscScalar    *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4;
284   MatScalar      *aa;
285   PetscInt       cm=C->rmap->n, cn=B->cmap->n, bm=B->rmap->n, col, i,j,n,*aj, am = A->rmap->n;
286   PetscInt       am2 = 2*am, am3 = 3*am,  bm4 = 4*bm,colam;
287 
288   PetscFunctionBegin;
289   if (!cm || !cn) PetscFunctionReturn(0);
290   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);
291   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);
292   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);
293   ierr = MatGetArray(B,&b);CHKERRQ(ierr);
294   ierr = MatGetArray(C,&c);CHKERRQ(ierr);
295   b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;
296   for (col=0; col<cn-4; col += 4){  /* over columns of C */
297     colam = col*am;
298     for (i=0; i<am; i++) {        /* over rows of C in those columns */
299       r1 = r2 = r3 = r4 = 0.0;
300       n   = a->i[i+1] - a->i[i];
301       aj  = a->j + a->i[i];
302       aa  = a->a + a->i[i];
303       for (j=0; j<n; j++) {
304         r1 += (*aa)*b1[*aj];
305         r2 += (*aa)*b2[*aj];
306         r3 += (*aa)*b3[*aj];
307         r4 += (*aa++)*b4[*aj++];
308       }
309       c[colam + i]       = r1;
310       c[colam + am + i]  = r2;
311       c[colam + am2 + i] = r3;
312       c[colam + am3 + i] = r4;
313     }
314     b1 += bm4;
315     b2 += bm4;
316     b3 += bm4;
317     b4 += bm4;
318   }
319   for (;col<cn; col++){     /* over extra columns of C */
320     for (i=0; i<am; i++) {  /* over rows of C in those columns */
321       r1 = 0.0;
322       n   = a->i[i+1] - a->i[i];
323       aj  = a->j + a->i[i];
324       aa  = a->a + a->i[i];
325 
326       for (j=0; j<n; j++) {
327         r1 += (*aa++)*b1[*aj++];
328       }
329       c[col*am + i]     = r1;
330     }
331     b1 += bm;
332   }
333   ierr = PetscLogFlops(cn*(2*a->nz));CHKERRQ(ierr);
334   ierr = MatRestoreArray(B,&b);CHKERRQ(ierr);
335   ierr = MatRestoreArray(C,&c);CHKERRQ(ierr);
336   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
337   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
338   PetscFunctionReturn(0);
339 }
340 
341 /*
342    Note very similar to MatMult_SeqAIJ(), should generate both codes from same base
343 */
344 #undef __FUNCT__
345 #define __FUNCT__ "MatMatMultNumericAdd_SeqAIJ_SeqDense"
346 PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
347 {
348   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
349   PetscErrorCode ierr;
350   PetscScalar    *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4;
351   MatScalar      *aa;
352   PetscInt       cm=C->rmap->n, cn=B->cmap->n, bm=B->rmap->n, col, i,j,n,*aj, am = A->rmap->n,*ii,arm;
353   PetscInt       am2 = 2*am, am3 = 3*am,  bm4 = 4*bm,colam,*ridx;
354 
355   PetscFunctionBegin;
356   if (!cm || !cn) PetscFunctionReturn(0);
357   ierr = MatGetArray(B,&b);CHKERRQ(ierr);
358   ierr = MatGetArray(C,&c);CHKERRQ(ierr);
359   b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;
360 
361   if (a->compressedrow.use){ /* use compressed row format */
362     for (col=0; col<cn-4; col += 4){  /* over columns of C */
363       colam = col*am;
364       arm   = a->compressedrow.nrows;
365       ii    = a->compressedrow.i;
366       ridx  = a->compressedrow.rindex;
367       for (i=0; i<arm; i++) {        /* over rows of C in those columns */
368 	r1 = r2 = r3 = r4 = 0.0;
369 	n   = ii[i+1] - ii[i];
370 	aj  = a->j + ii[i];
371 	aa  = a->a + ii[i];
372 	for (j=0; j<n; j++) {
373 	  r1 += (*aa)*b1[*aj];
374 	  r2 += (*aa)*b2[*aj];
375 	  r3 += (*aa)*b3[*aj];
376 	  r4 += (*aa++)*b4[*aj++];
377 	}
378 	c[colam       + ridx[i]] += r1;
379 	c[colam + am  + ridx[i]] += r2;
380 	c[colam + am2 + ridx[i]] += r3;
381 	c[colam + am3 + ridx[i]] += r4;
382       }
383       b1 += bm4;
384       b2 += bm4;
385       b3 += bm4;
386       b4 += bm4;
387     }
388     for (;col<cn; col++){     /* over extra columns of C */
389       colam = col*am;
390       arm   = a->compressedrow.nrows;
391       ii    = a->compressedrow.i;
392       ridx  = a->compressedrow.rindex;
393       for (i=0; i<arm; i++) {  /* over rows of C in those columns */
394 	r1 = 0.0;
395 	n   = ii[i+1] - ii[i];
396 	aj  = a->j + ii[i];
397 	aa  = a->a + ii[i];
398 
399 	for (j=0; j<n; j++) {
400 	  r1 += (*aa++)*b1[*aj++];
401 	}
402 	c[col*am + ridx[i]] += r1;
403       }
404       b1 += bm;
405     }
406   } else {
407     for (col=0; col<cn-4; col += 4){  /* over columns of C */
408       colam = col*am;
409       for (i=0; i<am; i++) {        /* over rows of C in those columns */
410 	r1 = r2 = r3 = r4 = 0.0;
411 	n   = a->i[i+1] - a->i[i];
412 	aj  = a->j + a->i[i];
413 	aa  = a->a + a->i[i];
414 	for (j=0; j<n; j++) {
415 	  r1 += (*aa)*b1[*aj];
416 	  r2 += (*aa)*b2[*aj];
417 	  r3 += (*aa)*b3[*aj];
418 	  r4 += (*aa++)*b4[*aj++];
419 	}
420 	c[colam + i]       += r1;
421 	c[colam + am + i]  += r2;
422 	c[colam + am2 + i] += r3;
423 	c[colam + am3 + i] += r4;
424       }
425       b1 += bm4;
426       b2 += bm4;
427       b3 += bm4;
428       b4 += bm4;
429     }
430     for (;col<cn; col++){     /* over extra columns of C */
431       for (i=0; i<am; i++) {  /* over rows of C in those columns */
432 	r1 = 0.0;
433 	n   = a->i[i+1] - a->i[i];
434 	aj  = a->j + a->i[i];
435 	aa  = a->a + a->i[i];
436 
437 	for (j=0; j<n; j++) {
438 	  r1 += (*aa++)*b1[*aj++];
439 	}
440 	c[col*am + i]     += r1;
441       }
442       b1 += bm;
443     }
444   }
445   ierr = PetscLogFlops(cn*2*a->nz);CHKERRQ(ierr);
446   ierr = MatRestoreArray(B,&b);CHKERRQ(ierr);
447   ierr = MatRestoreArray(C,&c);CHKERRQ(ierr);
448   PetscFunctionReturn(0);
449 }
450