xref: /petsc/src/mat/impls/aij/seq/matmatmult.c (revision 26be04463c6db0f5711e8cc51d5b3a4aa4a79b1d)
1 /*
2   Defines matrix-matrix product routines for pairs of AIJ matrices
3           C = A * B
4 */
5 
6 #include "src/mat/impls/aij/seq/aij.h" /*I "petscmat.h" I*/
7 #include "src/mat/utils/freespace.h"
8 #include "src/mat/impls/aij/mpi/mpiaij.h"
9 
10 /*
11   Initialize a linked list
12   Input Parameters:
13     lnk_init  - the initial index value indicating the entry in the list is not set yet
14     nlnk      - max length of the list
15     lnk       - linked list(an integer array) that is allocated
16   output Parameters:
17     lnk       - the linked list with all values set as lnk_int
18 */
19 #define LNKLISTINITIALIZE(lnk_init,nlnk,lnk){\
20   int _i;\
21   for (_i=0; _i<nlnk; _i++) lnk[_i] = lnk_init;\
22 }
23 
24 /*
25   Add a index set into a sorted linked list
26   Input Parameters:
27     nidx      - number of input indices
28     indices   - interger array
29     lnk_head  - the header of the list
30     lnk_init  - the initial index value indicating the entry in the list is not set yet
31     lnk       - linked list(an integer array) that is created
32   output Parameters:
33     nlnk      - number of newly added indices
34     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from indices
35 */
36 #define LNKLISTADD(nidx,indices,lnk_head,lnk_init,nlnk,lnk){\
37   int _k,_entry,_lidx=lnk_head,_idx;\
38   nlnk = 0;\
39   _k=nidx;\
40   while (_k){/* assume indices are almost in increasing order, starting from its end saves computation */\
41     _entry = indices[--_k];\
42     if (lnk[_entry] == lnk_init) { /* new col */\
43       do {\
44         _idx = _lidx;\
45         _lidx  = lnk[_idx];\
46       } while (_entry > _lidx);\
47       lnk[_idx] = _entry;\
48       lnk[_entry] = _lidx;\
49       nlnk++;\
50     }\
51   }\
52 }
53 /*
54   Copy data on the list into an array, then initialize the list
55   Input Parameters:
56     lnk_head  - the header of the list
57     lnk_init  - the initial index value indicating the entry in the list is not set yet
58     nlnk      - number of data on the list to be copied
59     lnk       - linked list
60   output Parameters:
61     indices   - array that contains the copied data
62 */
63 #define LNKLISTCLEAR(lnk_head,lnk_init,nlnk,lnk,indices){\
64   int _j,_idx=lnk_head,_idx0;\
65   for (_j=0; _j<nlnk; _j++){\
66     _idx0 = _idx; _idx = lnk[_idx0];\
67     *(indices+_j) = _idx;\
68     lnk[_idx0] = lnk_init;\
69   }\
70   lnk[_idx] = lnk_init;\
71 }
72 
73 typedef struct { /* used by MatMatMult_MPIAIJ_MPIAIJ for reusing symbolic mat product */
74   IS        isrowa,isrowb,iscolb;
75   Mat       *aseq,*bseq,C_seq;
76 } Mat_MatMatMultMPI;
77 
78 static int logkey_matmatmult_symbolic = 0;
79 static int logkey_matmatmult_numeric  = 0;
80 
81 #undef __FUNCT__
82 #define __FUNCT__ "MatMatMult"
83 /*@
84    MatMatMult - Performs Matrix-Matrix Multiplication C=A*B.
85 
86    Collective on Mat
87 
88    Input Parameters:
89 +  A - the left matrix
90 .  B - the right matrix
91 .  scall - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX
92 -  fill - expected fill as ratio of nnz(C)/(nnz(A) + nnz(B))
93 
94    Output Parameters:
95 .  C - the product matrix
96 
97    Notes:
98    C will be created and must be destroyed by the user with MatDestroy().
99 
100    This routine is currently only implemented for pairs of SeqAIJ matrices and classes
101    which inherit from SeqAIJ.  C will be of type MATSEQAIJ.
102 
103    Level: intermediate
104 
105 .seealso: MatMatMultSymbolic(),MatMatMultNumeric()
106 @*/
107 int MatMatMult(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
108 {
109   int  ierr;
110 
111   PetscFunctionBegin;
112   PetscValidHeaderSpecific(A,MAT_COOKIE,1);
113   PetscValidType(A,1);
114   MatPreallocated(A);
115   if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
116   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
117   PetscValidHeaderSpecific(B,MAT_COOKIE,2);
118   PetscValidType(B,2);
119   MatPreallocated(B);
120   if (!B->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
121   if (B->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
122   PetscValidPointer(C,3);
123   if (B->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",B->M,A->N);
124 
125   if (fill <=0) SETERRQ1(PETSC_ERR_ARG_SIZ,"fill=%d must be > 0",fill);
126 
127   ierr = PetscLogEventBegin(MAT_MatMult,A,B,0,0);CHKERRQ(ierr);
128   ierr = (*A->ops->matmult)(A,B,scall,fill,C);CHKERRQ(ierr);
129   ierr = PetscLogEventEnd(MAT_MatMult,A,B,0,0);CHKERRQ(ierr);
130 
131   PetscFunctionReturn(0);
132 }
133 
134 #undef __FUNCT__
135 #define __FUNCT__ "MatMatMult_MPIAIJ_MPIAIJ"
136 int MatMatMult_MPIAIJ_MPIAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill, Mat *C)
137 {
138   int ierr;
139 
140   PetscFunctionBegin;
141   if (scall == MAT_INITIAL_MATRIX){
142     ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ(A,B,fill,C);CHKERRQ(ierr);
143   }
144   ierr = MatMatMultNumeric_MPIAIJ_MPIAIJ(A,B,*C);CHKERRQ(ierr);
145   PetscFunctionReturn(0);
146 }
147 
148 #undef __FUNCT__
149 #define __FUNCT__ "MatMatMult_SeqAIJ_SeqAIJ"
150 int MatMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) {
151   int ierr;
152 
153   PetscFunctionBegin;
154   if (scall == MAT_INITIAL_MATRIX){
155     ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
156   }
157   ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr);
158   PetscFunctionReturn(0);
159 }
160 
161 #undef __FUNCT__
162 #define __FUNCT__ "MatMatMultSymbolic"
163 /*@
164    MatMatMultSymbolic - Performs construction, preallocation, and computes the ij structure
165    of the matrix-matrix product C=A*B.  Call this routine before calling MatMatMultNumeric().
166 
167    Collective on Mat
168 
169    Input Parameters:
170 +  A - the left matrix
171 .  B - the right matrix
172 -  fill - expected fill as ratio of nnz(C)/(nnz(A) + nnz(B))
173 
174    Output Parameters:
175 .  C - the matrix containing the ij structure of product matrix
176 
177    Notes:
178    C will be created as a MATSEQAIJ matrix and must be destroyed by the user with MatDestroy().
179 
180    This routine is currently only implemented for SeqAIJ matrices and classes which inherit from SeqAIJ.
181 
182    Level: intermediate
183 
184 .seealso: MatMatMult(),MatMatMultNumeric()
185 @*/
186 int MatMatMultSymbolic(Mat A,Mat B,PetscReal fill,Mat *C) {
187   /* Perhaps this "interface" routine should be moved into the interface directory.*/
188   /* To facilitate implementations with varying types, QueryFunction is used.*/
189   /* It is assumed that implementations will be composed as "MatMatMultSymbolic_<type of A><type of B>". */
190   int  ierr;
191   char symfunct[80];
192   int  (*symbolic)(Mat,Mat,PetscReal,Mat *);
193 
194   PetscFunctionBegin;
195   PetscValidHeaderSpecific(A,MAT_COOKIE,1);
196   PetscValidType(A,1);
197   MatPreallocated(A);
198   if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
199   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
200 
201   PetscValidHeaderSpecific(B,MAT_COOKIE,2);
202   PetscValidType(B,2);
203   MatPreallocated(B);
204   if (!B->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
205   if (B->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
206   PetscValidPointer(C,3);
207 
208   if (B->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",B->M,A->N);
209   if (fill <=0) SETERRQ1(PETSC_ERR_ARG_SIZ,"fill=%d must be > 0",fill);
210 
211   /* Currently only _seqaijseqaij is implemented, so just query for it in A and B. */
212   /* When other implementations exist, attack the multiple dispatch problem. */
213   ierr = PetscStrcpy(symfunct,"MatMatMultSymbolic_seqaijseqaij");CHKERRQ(ierr);
214   ierr = PetscObjectQueryFunction((PetscObject)B,symfunct,(PetscVoidFunction)&symbolic);CHKERRQ(ierr);
215   if (!symbolic) SETERRQ1(PETSC_ERR_SUP,"C=A*B not implemented for B of type %s",B->type_name);
216   ierr = PetscObjectQueryFunction((PetscObject)A,symfunct,(PetscVoidFunction)&symbolic);CHKERRQ(ierr);
217   if (!symbolic) SETERRQ1(PETSC_ERR_SUP,"C=A*B not implemented for A of type %s",A->type_name);
218   ierr = (*symbolic)(A,B,fill,C);CHKERRQ(ierr);
219 
220   PetscFunctionReturn(0);
221 }
222 #undef __FUNCT__
223 #define __FUNCT__ "MatDestroy_MPIAIJ_MatMatMult"
224 int MatDestroy_MPIAIJ_MatMatMult(Mat A)
225 {
226   int               ierr;
227   Mat_MatMatMultMPI *mult = ( Mat_MatMatMultMPI*)A->spptr;
228 
229   PetscFunctionBegin;
230   /* printf(" ...MatDestroy_MPIAIJ_MatMatMult is called\n"); */
231   ierr = PetscFree(mult);CHKERRQ(ierr);
232   ierr = MatDestroy_MPIAIJ(A);CHKERRQ(ierr);
233 
234   PetscFunctionReturn(0);
235 }
236 
237 #undef __FUNCT__
238 #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ"
239 int MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
240 {
241   Mat           *aseq,*bseq,A_seq=PETSC_NULL,B_seq=PETSC_NULL,C_seq;
242   Mat_MPIAIJ    *a = (Mat_MPIAIJ*)A->data;
243   int           ierr,*idx,i,start,end,ncols,imark,nzA,nzB,*cmap;
244   IS            isrow,iscol;
245   Mat_MatMatMultMPI *mult;
246 
247   PetscFunctionBegin;
248   /* create a seq matrix B_seq = submatrix of B by taking rows of B that equal to nonzero col of A */
249   start = a->cstart;
250   cmap  = a->garray;
251   nzA   = a->A->n;
252   nzB   = a->B->n;
253   ierr  = PetscMalloc((nzA+nzB)*sizeof(int), &idx);CHKERRQ(ierr);
254   ncols = 0;
255   for (i=0; i<nzB; i++) {
256     if (cmap[i] < start) idx[ncols++] = cmap[i];
257     else break;
258   }
259   imark = i;
260   for (i=0; i<nzA; i++) idx[ncols++] = start + i;
261   for (i=imark; i<nzB; i++) idx[ncols++] = cmap[i];
262   ierr = ISCreateGeneral(PETSC_COMM_SELF,ncols,idx,&isrow);CHKERRQ(ierr); /* isrow of B = iscol of A! */
263   ierr = PetscFree(idx);CHKERRQ(ierr);
264   ierr = ISCreateStride(PETSC_COMM_SELF,B->N,0,1,&iscol);CHKERRQ(ierr);
265   ierr = MatGetSubMatrices(B,1,&isrow,&iscol,MAT_INITIAL_MATRIX,&bseq);CHKERRQ(ierr);
266   ierr = ISDestroy(iscol);CHKERRQ(ierr);
267   B_seq = bseq[0];
268 
269   /*  create a seq matrix A_seq = submatrix of A by taking all local rows of A */
270   start = a->rstart; end = a->rend;
271   ierr = ISCreateStride(PETSC_COMM_SELF,end-start,start,1,&iscol);CHKERRQ(ierr); /* isrow of A = iscol */
272   ierr = MatGetSubMatrices(A,1,&iscol,&isrow,MAT_INITIAL_MATRIX,&aseq);CHKERRQ(ierr);
273   ierr = ISDestroy(isrow);CHKERRQ(ierr);
274   ierr = ISDestroy(iscol);CHKERRQ(ierr);
275   A_seq = aseq[0];
276 
277   /* compute C_seq = A_seq * B_seq */
278   ierr = MatMatMult_SeqAIJ_SeqAIJ(A_seq,B_seq,MAT_INITIAL_MATRIX,fill,&C_seq);CHKERRQ(ierr);
279   /*
280   int rank;
281   ierr = MPI_Comm_rank(A->comm,&rank);CHKERRQ(ierr);
282   ierr = PetscPrintf(PETSC_COMM_SELF," [%d] A_seq: %d, %d; B_seq: %d, %d; C_seq: %d, %d;\n",rank,A_seq->m,A_seq->n,B_seq->m,B_seq->n,C_seq->m,C_seq->n);
283   */
284   ierr = MatDestroyMatrices(1,&aseq);CHKERRQ(ierr);
285   ierr = MatDestroyMatrices(1,&bseq);CHKERRQ(ierr);
286 
287   /* create mpi matrix C by concatinating C_seq */
288   ierr = MatMerge(A->comm,C_seq,C);CHKERRQ(ierr);
289 
290   /* attach the supporting struct to C for reuse of symbolic C */
291   ierr = PetscNew(Mat_MatMatMultMPI,&mult);CHKERRQ(ierr);
292   (*C)->spptr         = (void*)mult;
293   (*C)->ops->destroy  = MatDestroy_MPIAIJ_MatMatMult;
294 
295   PetscFunctionReturn(0);
296 }
297 
298 #undef __FUNCT__
299 #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ"
300 int MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
301 {
302   int            ierr;
303   FreeSpaceList  free_space=PETSC_NULL,current_space=PETSC_NULL;
304   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
305   int            *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj;
306   int            *ci,*cj,*lnk,idx0,idx;
307   int            am=A->M,bn=B->N,bm=B->M;
308   int            i,j,anzi,brow,bnzj,cnzi,nlnk,lnk_init=-1,nspacedouble=0;
309   MatScalar      *ca;
310 
311   PetscFunctionBegin;
312   /* Start timers */
313   ierr = PetscLogEventBegin(logkey_matmatmult_symbolic,A,B,0,0);CHKERRQ(ierr);
314   /* Set up */
315   /* Allocate ci array, arrays for fill computation and */
316   /* free space for accumulating nonzero column info */
317   ierr = PetscMalloc(((am+1)+1)*sizeof(int),&ci);CHKERRQ(ierr);
318   ci[0] = 0;
319 
320   ierr = PetscMalloc((bn+1)*sizeof(int),&lnk);CHKERRQ(ierr);
321   LNKLISTINITIALIZE(lnk_init,bn,lnk);
322 
323   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
324   ierr = GetMoreSpace(fill*(ai[am]+bi[bm]),&free_space);CHKERRQ(ierr);
325   current_space = free_space;
326 
327   /* Determine symbolic info for each row of the product: */
328   for (i=0;i<am;i++) {
329     anzi = ai[i+1] - ai[i];
330     cnzi = 0;
331     lnk[bn] = bn;
332     j       = anzi;
333     aj      = a->j + ai[i];
334     while (j){/* assume cols are almost in increasing order, starting from its end saves computation */
335       j--;
336       brow = *(aj + j);
337       bnzj = bi[brow+1] - bi[brow];
338       bjj  = bj + bi[brow];
339       /* add non-zero cols of B into the sorted linked list lnk */
340       LNKLISTADD(bnzj,bjj,bn,lnk_init,nlnk,lnk);
341       cnzi += nlnk;
342     }
343 
344     /* If free space is not available, make more free space */
345     /* Double the amount of total space in the list */
346     if (current_space->local_remaining<cnzi) {
347       /* printf("MatMatMultSymbolic_SeqAIJ_SeqAIJ()...%d -th row, double space ...\n",i);*/
348       ierr = GetMoreSpace(current_space->total_array_size,&current_space);CHKERRQ(ierr);
349       nspacedouble++;
350     }
351 
352     /* Copy data into free space, then initialize lnk */
353     LNKLISTCLEAR(bn,lnk_init,cnzi,lnk,current_space->array);
354     current_space->array += cnzi;
355 
356     current_space->local_used      += cnzi;
357     current_space->local_remaining -= cnzi;
358 
359     ci[i+1] = ci[i] + cnzi;
360   }
361 
362   /* Column indices are in the list of free space */
363   /* Allocate space for cj, initialize cj, and */
364   /* destroy list of free space and other temporary array(s) */
365   ierr = PetscMalloc((ci[am]+1)*sizeof(int),&cj);CHKERRQ(ierr);
366   ierr = MakeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
367   ierr = PetscFree(lnk);CHKERRQ(ierr);
368 
369   /* Allocate space for ca */
370   ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
371   ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr);
372 
373   /* put together the new symbolic matrix */
374   ierr = MatCreateSeqAIJWithArrays(A->comm,am,bn,ci,cj,ca,C);CHKERRQ(ierr);
375 
376   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
377   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
378   c = (Mat_SeqAIJ *)((*C)->data);
379   c->freedata = PETSC_TRUE;
380   c->nonew    = 0;
381 
382   ierr = PetscLogEventEnd(logkey_matmatmult_symbolic,A,B,0,0);CHKERRQ(ierr);
383   PetscLogInfo((PetscObject)(*C),"Number of calls to GetMoreSpace(): %d\n",nspacedouble);
384   PetscFunctionReturn(0);
385 }
386 
387 #undef __FUNCT__
388 #define __FUNCT__ "MatMatMultNumeric"
389 /*@
390    MatMatMultNumeric - Performs the numeric matrix-matrix product.
391    Call this routine after first calling MatMatMultSymbolic().
392 
393    Collective on Mat
394 
395    Input Parameters:
396 +  A - the left matrix
397 -  B - the right matrix
398 
399    Output Parameters:
400 .  C - the product matrix, whose ij structure was defined from MatMatMultSymbolic().
401 
402    Notes:
403    C must have been created with MatMatMultSymbolic.
404 
405    This routine is currently only implemented for SeqAIJ type matrices.
406 
407    Level: intermediate
408 
409 .seealso: MatMatMult(),MatMatMultSymbolic()
410 @*/
411 int MatMatMultNumeric(Mat A,Mat B,Mat C){
412   /* Perhaps this "interface" routine should be moved into the interface directory.*/
413   /* To facilitate implementations with varying types, QueryFunction is used.*/
414   /* It is assumed that implementations will be composed as "MatMatMultNumeric_<type of A><type of B>". */
415   int ierr;
416   char numfunct[80];
417   int (*numeric)(Mat,Mat,Mat);
418 
419   PetscFunctionBegin;
420 
421   PetscValidHeaderSpecific(A,MAT_COOKIE,1);
422   PetscValidType(A,1);
423   MatPreallocated(A);
424   if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
425   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
426 
427   PetscValidHeaderSpecific(B,MAT_COOKIE,2);
428   PetscValidType(B,2);
429   MatPreallocated(B);
430   if (!B->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
431   if (B->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
432 
433   PetscValidHeaderSpecific(C,MAT_COOKIE,3);
434   PetscValidType(C,3);
435   MatPreallocated(C);
436   if (!C->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
437   if (C->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
438 
439   if (B->N!=C->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",B->N,C->N);
440   if (B->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",B->M,A->N);
441   if (A->M!=C->M) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",A->M,C->M);
442 
443   /* Currently only _seqaijseqaij is implemented, so just query for it in A and B. */
444   /* When other implementations exist, attack the multiple dispatch problem. */
445   ierr = PetscStrcpy(numfunct,"MatMatMultNumeric_seqaijseqaij");CHKERRQ(ierr);
446   ierr = PetscObjectQueryFunction((PetscObject)A,numfunct,(PetscVoidFunction)&numeric);CHKERRQ(ierr);
447   if (!numeric) SETERRQ1(PETSC_ERR_SUP,"C=A*B not implemented for A of type %s",A->type_name);
448   ierr = PetscObjectQueryFunction((PetscObject)B,numfunct,(PetscVoidFunction)&numeric);CHKERRQ(ierr);
449   if (!numeric) SETERRQ1(PETSC_ERR_SUP,"C=A*B not implemented for B of type %s",B->type_name);
450 
451   ierr = (*numeric)(A,B,C);CHKERRQ(ierr);
452 
453   PetscFunctionReturn(0);
454 }
455 
456 #undef __FUNCT__
457 #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ"
458 int MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat A,Mat B,Mat C)
459 {
460   PetscFunctionBegin;
461   PetscFunctionReturn(0);
462 }
463 
464 #undef __FUNCT__
465 #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqAIJ"
466 int MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
467 {
468   int        ierr,flops=0;
469   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data;
470   Mat_SeqAIJ *b = (Mat_SeqAIJ *)B->data;
471   Mat_SeqAIJ *c = (Mat_SeqAIJ *)C->data;
472   int        *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j;
473   int        am=A->M,cn=C->N;
474   int        i,j,k,anzi,bnzi,cnzi,brow;
475   MatScalar  *aa=a->a,*ba=b->a,*baj,*ca=c->a,*temp;
476 
477   PetscFunctionBegin;
478 
479   /* Start timers */
480   ierr = PetscLogEventBegin(logkey_matmatmult_numeric,A,B,C,0);CHKERRQ(ierr);
481 
482   /* Allocate temp accumulation space to avoid searching for nonzero columns in C */
483   ierr = PetscMalloc((cn+1)*sizeof(MatScalar),&temp);CHKERRQ(ierr);
484   ierr = PetscMemzero(temp,cn*sizeof(MatScalar));CHKERRQ(ierr);
485   /* Traverse A row-wise. */
486   /* Build the ith row in C by summing over nonzero columns in A, */
487   /* the rows of B corresponding to nonzeros of A. */
488   for (i=0;i<am;i++) {
489     anzi = ai[i+1] - ai[i];
490     for (j=0;j<anzi;j++) {
491       brow = *aj++;
492       bnzi = bi[brow+1] - bi[brow];
493       bjj  = bj + bi[brow];
494       baj  = ba + bi[brow];
495       for (k=0;k<bnzi;k++) {
496         temp[bjj[k]] += (*aa)*baj[k];
497       }
498       flops += 2*bnzi;
499       aa++;
500     }
501     /* Store row back into C, and re-zero temp */
502     cnzi = ci[i+1] - ci[i];
503     for (j=0;j<cnzi;j++) {
504       ca[j] = temp[cj[j]];
505       temp[cj[j]] = 0.0;
506     }
507     ca += cnzi;
508     cj += cnzi;
509   }
510   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
511   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
512 
513   /* Free temp */
514   ierr = PetscFree(temp);CHKERRQ(ierr);
515   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
516   ierr = PetscLogEventEnd(logkey_matmatmult_numeric,A,B,C,0);CHKERRQ(ierr);
517   PetscFunctionReturn(0);
518 }
519 
520 #undef __FUNCT__
521 #define __FUNCT__ "RegisterMatMatMultRoutines_Private"
522 int RegisterMatMatMultRoutines_Private(Mat A) {
523   int ierr;
524 
525   PetscFunctionBegin;
526   if (!logkey_matmatmult_symbolic) {
527     ierr = PetscLogEventRegister(&logkey_matmatmult_symbolic,"MatMatMultSymbolic",MAT_COOKIE);CHKERRQ(ierr);
528   }
529   /*
530   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMatMultSymbolic_seqaijseqaij",
531                                            "MatMatMultSymbolic_SeqAIJ_SeqAIJ",
532                                            MatMatMultSymbolic_SeqAIJ_SeqAIJ);CHKERRQ(ierr);*/
533   if (!logkey_matmatmult_numeric) {
534     ierr = PetscLogEventRegister(&logkey_matmatmult_numeric,"MatMatMultNumeric",MAT_COOKIE);CHKERRQ(ierr);
535   }
536   /*
537   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMatMultNumeric_seqaijseqaij",
538                                            "MatMatMultNumeric_SeqAIJ_SeqAIJ",
539                                            MatMatMultNumeric_SeqAIJ_SeqAIJ);CHKERRQ(ierr);*/
540   PetscFunctionReturn(0);
541 }
542