xref: /petsc/src/mat/impls/aij/seq/matmatmult.c (revision 9af31e4ad595286b4e2df8194fee047feeccfe42)
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.0) SETERRQ1(PETSC_ERR_ARG_SIZ,"fill=%g must be > 0.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);/* numeric product is computed as well */
143   } else if (scall == MAT_REUSE_MATRIX){
144     ierr = MatMatMultNumeric_MPIAIJ_MPIAIJ(A,B,*C);CHKERRQ(ierr);
145   } else {
146     SETERRQ1(PETSC_ERR_ARG_WRONG,"Invalid MatReuse %d",scall);
147   }
148   PetscFunctionReturn(0);
149 }
150 
151 #undef __FUNCT__
152 #define __FUNCT__ "MatMatMult_SeqAIJ_SeqAIJ"
153 int MatMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) {
154   int ierr;
155 
156   PetscFunctionBegin;
157   if (scall == MAT_INITIAL_MATRIX){
158     ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
159   }
160   ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr);
161   PetscFunctionReturn(0);
162 }
163 
164 #undef __FUNCT__
165 #define __FUNCT__ "MatMatMultSymbolic"
166 /*@
167    MatMatMultSymbolic - Performs construction, preallocation, and computes the ij structure
168    of the matrix-matrix product C=A*B.  Call this routine before calling MatMatMultNumeric().
169 
170    Collective on Mat
171 
172    Input Parameters:
173 +  A - the left matrix
174 .  B - the right matrix
175 -  fill - expected fill as ratio of nnz(C)/(nnz(A) + nnz(B))
176 
177    Output Parameters:
178 .  C - the matrix containing the ij structure of product matrix
179 
180    Notes:
181    C will be created as a MATSEQAIJ matrix and must be destroyed by the user with MatDestroy().
182 
183    This routine is currently only implemented for SeqAIJ matrices and classes which inherit from SeqAIJ.
184 
185    Level: intermediate
186 
187 .seealso: MatMatMult(),MatMatMultNumeric()
188 @*/
189 int MatMatMultSymbolic(Mat A,Mat B,PetscReal fill,Mat *C) {
190   /* Perhaps this "interface" routine should be moved into the interface directory.*/
191   /* To facilitate implementations with varying types, QueryFunction is used.*/
192   /* It is assumed that implementations will be composed as "MatMatMultSymbolic_<type of A><type of B>". */
193   int  ierr;
194   char symfunct[80];
195   int  (*symbolic)(Mat,Mat,PetscReal,Mat *);
196 
197   PetscFunctionBegin;
198   PetscValidHeaderSpecific(A,MAT_COOKIE,1);
199   PetscValidType(A,1);
200   MatPreallocated(A);
201   if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
202   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
203 
204   PetscValidHeaderSpecific(B,MAT_COOKIE,2);
205   PetscValidType(B,2);
206   MatPreallocated(B);
207   if (!B->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
208   if (B->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
209   PetscValidPointer(C,3);
210 
211   if (B->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",B->M,A->N);
212   if (fill <=0.0) SETERRQ1(PETSC_ERR_ARG_SIZ,"fill=%g must be > 0.0",fill);
213 
214   /* Currently only _seqaijseqaij is implemented, so just query for it in A and B. */
215   /* When other implementations exist, attack the multiple dispatch problem. */
216   ierr = PetscStrcpy(symfunct,"MatMatMultSymbolic_seqaijseqaij");CHKERRQ(ierr);
217   ierr = PetscObjectQueryFunction((PetscObject)B,symfunct,(PetscVoidFunction)&symbolic);CHKERRQ(ierr);
218   if (!symbolic) SETERRQ1(PETSC_ERR_SUP,"C=A*B not implemented for B of type %s",B->type_name);
219   ierr = PetscObjectQueryFunction((PetscObject)A,symfunct,(PetscVoidFunction)&symbolic);CHKERRQ(ierr);
220   if (!symbolic) SETERRQ1(PETSC_ERR_SUP,"C=A*B not implemented for A of type %s",A->type_name);
221   ierr = (*symbolic)(A,B,fill,C);CHKERRQ(ierr);
222 
223   PetscFunctionReturn(0);
224 }
225 
226 EXTERN int MatDestroy_MPIAIJ(Mat);
227 #undef __FUNCT__
228 #define __FUNCT__ "MatDestroy_MPIAIJ_MatMatMult"
229 int MatDestroy_MPIAIJ_MatMatMult(Mat A)
230 {
231   int               ierr;
232   Mat_MatMatMultMPI *mult=(Mat_MatMatMultMPI*)A->spptr;
233 
234   PetscFunctionBegin;
235   ierr = ISDestroy(mult->isrowb);CHKERRQ(ierr);
236   ierr = ISDestroy(mult->iscolb);CHKERRQ(ierr);
237   ierr = ISDestroy(mult->isrowa);CHKERRQ(ierr);
238   ierr = MatDestroyMatrices(1,&mult->aseq);CHKERRQ(ierr);
239   ierr = MatDestroyMatrices(1,&mult->bseq);CHKERRQ(ierr);
240   ierr = MatDestroy(mult->C_seq);CHKERRQ(ierr);
241   ierr = PetscFree(mult);CHKERRQ(ierr);
242 
243   ierr = MatDestroy_MPIAIJ(A);CHKERRQ(ierr);
244 
245   PetscFunctionReturn(0);
246 }
247 
248 #undef __FUNCT__
249 #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ"
250 int MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
251 {
252   Mat_MPIAIJ        *a = (Mat_MPIAIJ*)A->data;
253   int               ierr,*idx,i,start,end,ncols,imark,nzA,nzB,*cmap;
254   Mat_MatMatMultMPI *mult;
255 
256   PetscFunctionBegin;
257   ierr = PetscNew(Mat_MatMatMultMPI,&mult);CHKERRQ(ierr);
258 
259   /* create a seq matrix B_seq = submatrix of B by taking rows of B that equal to nonzero col of A */
260   start = a->cstart;
261   cmap  = a->garray;
262   nzA   = a->A->n;
263   nzB   = a->B->n;
264   ierr  = PetscMalloc((nzA+nzB)*sizeof(int), &idx);CHKERRQ(ierr);
265   ncols = 0;
266   for (i=0; i<nzB; i++) {
267     if (cmap[i] < start) idx[ncols++] = cmap[i];
268     else break;
269   }
270   imark = i;
271   for (i=0; i<nzA; i++) idx[ncols++] = start + i;
272   for (i=imark; i<nzB; i++) idx[ncols++] = cmap[i];
273   ierr = ISCreateGeneral(PETSC_COMM_SELF,ncols,idx,&mult->isrowb);CHKERRQ(ierr);
274   ierr = PetscFree(idx);CHKERRQ(ierr);
275   ierr = ISCreateStride(PETSC_COMM_SELF,B->N,0,1,&mult->iscolb);CHKERRQ(ierr);
276   ierr = MatGetSubMatrices(B,1,&mult->isrowb,&mult->iscolb,MAT_INITIAL_MATRIX,&mult->bseq);CHKERRQ(ierr)
277 
278   /*  create a seq matrix A_seq = submatrix of A by taking all local rows of A */
279   start = a->rstart; end = a->rend;
280   ierr = ISCreateStride(PETSC_COMM_SELF,end-start,start,1,&mult->isrowa);CHKERRQ(ierr);
281   ierr = MatGetSubMatrices(A,1,&mult->isrowa,&mult->isrowb,MAT_INITIAL_MATRIX,&mult->aseq);CHKERRQ(ierr);
282 
283   /* compute C_seq = A_seq * B_seq */
284   ierr = MatMatMult_SeqAIJ_SeqAIJ(mult->aseq[0],mult->bseq[0],MAT_INITIAL_MATRIX,fill,&mult->C_seq);CHKERRQ(ierr);
285 
286   /* create mpi matrix C by concatinating C_seq */
287   ierr = PetscObjectReference((PetscObject)mult->C_seq);CHKERRQ(ierr); /* prevent C_seq being destroyed by MatMerge() */
288   ierr = MatMerge(A->comm,mult->C_seq,MAT_INITIAL_MATRIX,C);CHKERRQ(ierr);
289 
290   /* attach the supporting struct to C for reuse of symbolic C */
291   (*C)->spptr         = (void*)mult;
292   (*C)->ops->destroy  = MatDestroy_MPIAIJ_MatMatMult;
293 
294   PetscFunctionReturn(0);
295 }
296 
297 #undef __FUNCT__
298 #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ"
299 int MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
300 {
301   int            ierr;
302   FreeSpaceList  free_space=PETSC_NULL,current_space=PETSC_NULL;
303   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
304   int            *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj;
305   int            *ci,*cj,*lnk;
306   int            am=A->M,bn=B->N,bm=B->M;
307   int            i,j,anzi,brow,bnzj,cnzi,nlnk,lnk_init=-1,nspacedouble=0;
308   MatScalar      *ca;
309 
310   PetscFunctionBegin;
311   /* Start timers */
312   ierr = PetscLogEventBegin(logkey_matmatmult_symbolic,A,B,0,0);CHKERRQ(ierr);
313   /* Set up */
314   /* Allocate ci array, arrays for fill computation and */
315   /* free space for accumulating nonzero column info */
316   ierr = PetscMalloc(((am+1)+1)*sizeof(int),&ci);CHKERRQ(ierr);
317   ci[0] = 0;
318 
319   ierr = PetscMalloc((bn+1)*sizeof(int),&lnk);CHKERRQ(ierr);
320   LNKLISTINITIALIZE(lnk_init,bn,lnk);
321 
322   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
323   ierr = GetMoreSpace((int)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr);
324   current_space = free_space;
325 
326   /* Determine symbolic info for each row of the product: */
327   for (i=0;i<am;i++) {
328     anzi = ai[i+1] - ai[i];
329     cnzi = 0;
330     lnk[bn] = bn;
331     j       = anzi;
332     aj      = a->j + ai[i];
333     while (j){/* assume cols are almost in increasing order, starting from its end saves computation */
334       j--;
335       brow = *(aj + j);
336       bnzj = bi[brow+1] - bi[brow];
337       bjj  = bj + bi[brow];
338       /* add non-zero cols of B into the sorted linked list lnk */
339       LNKLISTADD(bnzj,bjj,bn,lnk_init,nlnk,lnk);
340       cnzi += nlnk;
341     }
342 
343     /* If free space is not available, make more free space */
344     /* Double the amount of total space in the list */
345     if (current_space->local_remaining<cnzi) {
346       /* printf("MatMatMultSymbolic_SeqAIJ_SeqAIJ()...%d -th row, double space ...\n",i);*/
347       ierr = GetMoreSpace(current_space->total_array_size,&current_space);CHKERRQ(ierr);
348       nspacedouble++;
349     }
350 
351     /* Copy data into free space, then initialize lnk */
352     LNKLISTCLEAR(bn,lnk_init,cnzi,lnk,current_space->array);
353     current_space->array += cnzi;
354 
355     current_space->local_used      += cnzi;
356     current_space->local_remaining -= cnzi;
357 
358     ci[i+1] = ci[i] + cnzi;
359   }
360 
361   /* Column indices are in the list of free space */
362   /* Allocate space for cj, initialize cj, and */
363   /* destroy list of free space and other temporary array(s) */
364   ierr = PetscMalloc((ci[am]+1)*sizeof(int),&cj);CHKERRQ(ierr);
365   ierr = MakeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
366   ierr = PetscFree(lnk);CHKERRQ(ierr);
367 
368   /* Allocate space for ca */
369   ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
370   ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr);
371 
372   /* put together the new symbolic matrix */
373   ierr = MatCreateSeqAIJWithArrays(A->comm,am,bn,ci,cj,ca,C);CHKERRQ(ierr);
374 
375   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
376   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
377   c = (Mat_SeqAIJ *)((*C)->data);
378   c->freedata = PETSC_TRUE;
379   c->nonew    = 0;
380 
381   ierr = PetscLogEventEnd(logkey_matmatmult_symbolic,A,B,0,0);CHKERRQ(ierr);
382   PetscLogInfo((PetscObject)(*C),"Number of calls to GetMoreSpace(): %d\n",nspacedouble);
383   PetscFunctionReturn(0);
384 }
385 
386 #undef __FUNCT__
387 #define __FUNCT__ "MatMatMultNumeric"
388 /*@
389    MatMatMultNumeric - Performs the numeric matrix-matrix product.
390    Call this routine after first calling MatMatMultSymbolic().
391 
392    Collective on Mat
393 
394    Input Parameters:
395 +  A - the left matrix
396 -  B - the right matrix
397 
398    Output Parameters:
399 .  C - the product matrix, whose ij structure was defined from MatMatMultSymbolic().
400 
401    Notes:
402    C must have been created with MatMatMultSymbolic.
403 
404    This routine is currently only implemented for SeqAIJ type matrices.
405 
406    Level: intermediate
407 
408 .seealso: MatMatMult(),MatMatMultSymbolic()
409 @*/
410 int MatMatMultNumeric(Mat A,Mat B,Mat C){
411   /* Perhaps this "interface" routine should be moved into the interface directory.*/
412   /* To facilitate implementations with varying types, QueryFunction is used.*/
413   /* It is assumed that implementations will be composed as "MatMatMultNumeric_<type of A><type of B>". */
414   int ierr;
415   char numfunct[80];
416   int (*numeric)(Mat,Mat,Mat);
417 
418   PetscFunctionBegin;
419 
420   PetscValidHeaderSpecific(A,MAT_COOKIE,1);
421   PetscValidType(A,1);
422   MatPreallocated(A);
423   if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
424   if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
425 
426   PetscValidHeaderSpecific(B,MAT_COOKIE,2);
427   PetscValidType(B,2);
428   MatPreallocated(B);
429   if (!B->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
430   if (B->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
431 
432   PetscValidHeaderSpecific(C,MAT_COOKIE,3);
433   PetscValidType(C,3);
434   MatPreallocated(C);
435   if (!C->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
436   if (C->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
437 
438   if (B->N!=C->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",B->N,C->N);
439   if (B->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",B->M,A->N);
440   if (A->M!=C->M) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",A->M,C->M);
441 
442   /* Currently only _seqaijseqaij is implemented, so just query for it in A and B. */
443   /* When other implementations exist, attack the multiple dispatch problem. */
444   ierr = PetscStrcpy(numfunct,"MatMatMultNumeric_seqaijseqaij");CHKERRQ(ierr);
445   ierr = PetscObjectQueryFunction((PetscObject)A,numfunct,(PetscVoidFunction)&numeric);CHKERRQ(ierr);
446   if (!numeric) SETERRQ1(PETSC_ERR_SUP,"C=A*B not implemented for A of type %s",A->type_name);
447   ierr = PetscObjectQueryFunction((PetscObject)B,numfunct,(PetscVoidFunction)&numeric);CHKERRQ(ierr);
448   if (!numeric) SETERRQ1(PETSC_ERR_SUP,"C=A*B not implemented for B of type %s",B->type_name);
449 
450   ierr = (*numeric)(A,B,C);CHKERRQ(ierr);
451 
452   PetscFunctionReturn(0);
453 }
454 
455 /* This routine is called ONLY in the case of reusing previously computed symbolic C */
456 #undef __FUNCT__
457 #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ"
458 int MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat A,Mat B,Mat C)
459 {
460   int               ierr;
461   Mat_MatMatMultMPI *mult=(Mat_MatMatMultMPI*)C->spptr;
462 
463   PetscFunctionBegin;
464   ierr = MatGetSubMatrices(B,1,&mult->isrowb,&mult->iscolb,MAT_REUSE_MATRIX,&mult->bseq);CHKERRQ(ierr)
465   ierr = MatGetSubMatrices(A,1,&mult->isrowa,&mult->isrowb,MAT_REUSE_MATRIX,&mult->aseq);CHKERRQ(ierr);
466   ierr = MatMatMult_SeqAIJ_SeqAIJ(mult->aseq[0],mult->bseq[0],MAT_REUSE_MATRIX,0.0,&mult->C_seq);CHKERRQ(ierr);
467 
468   ierr = PetscObjectReference((PetscObject)mult->C_seq);CHKERRQ(ierr);
469   ierr = MatMerge(A->comm,mult->C_seq,MAT_REUSE_MATRIX,&C);CHKERRQ(ierr);
470 
471   PetscFunctionReturn(0);
472 }
473 
474 #undef __FUNCT__
475 #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqAIJ"
476 int MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
477 {
478   int        ierr,flops=0;
479   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data;
480   Mat_SeqAIJ *b = (Mat_SeqAIJ *)B->data;
481   Mat_SeqAIJ *c = (Mat_SeqAIJ *)C->data;
482   int        *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j;
483   int        am=A->M,cn=C->N;
484   int        i,j,k,anzi,bnzi,cnzi,brow;
485   MatScalar  *aa=a->a,*ba=b->a,*baj,*ca=c->a,*temp;
486 
487   PetscFunctionBegin;
488 
489   /* Start timers */
490   ierr = PetscLogEventBegin(logkey_matmatmult_numeric,A,B,C,0);CHKERRQ(ierr);
491 
492   /* Allocate temp accumulation space to avoid searching for nonzero columns in C */
493   ierr = PetscMalloc((cn+1)*sizeof(MatScalar),&temp);CHKERRQ(ierr);
494   ierr = PetscMemzero(temp,cn*sizeof(MatScalar));CHKERRQ(ierr);
495   /* Traverse A row-wise. */
496   /* Build the ith row in C by summing over nonzero columns in A, */
497   /* the rows of B corresponding to nonzeros of A. */
498   for (i=0;i<am;i++) {
499     anzi = ai[i+1] - ai[i];
500     for (j=0;j<anzi;j++) {
501       brow = *aj++;
502       bnzi = bi[brow+1] - bi[brow];
503       bjj  = bj + bi[brow];
504       baj  = ba + bi[brow];
505       for (k=0;k<bnzi;k++) {
506         temp[bjj[k]] += (*aa)*baj[k];
507       }
508       flops += 2*bnzi;
509       aa++;
510     }
511     /* Store row back into C, and re-zero temp */
512     cnzi = ci[i+1] - ci[i];
513     for (j=0;j<cnzi;j++) {
514       ca[j] = temp[cj[j]];
515       temp[cj[j]] = 0.0;
516     }
517     ca += cnzi;
518     cj += cnzi;
519   }
520   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
521   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
522 
523   /* Free temp */
524   ierr = PetscFree(temp);CHKERRQ(ierr);
525   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
526   ierr = PetscLogEventEnd(logkey_matmatmult_numeric,A,B,C,0);CHKERRQ(ierr);
527   PetscFunctionReturn(0);
528 }
529 
530 #undef __FUNCT__
531 #define __FUNCT__ "RegisterMatMatMultRoutines_Private"
532 int RegisterMatMatMultRoutines_Private(Mat A) {
533   int ierr;
534 
535   PetscFunctionBegin;
536   if (!logkey_matmatmult_symbolic) {
537     ierr = PetscLogEventRegister(&logkey_matmatmult_symbolic,"MatMatMultSymbolic",MAT_COOKIE);CHKERRQ(ierr);
538   }
539   if (!logkey_matmatmult_numeric) {
540     ierr = PetscLogEventRegister(&logkey_matmatmult_numeric,"MatMatMultNumeric",MAT_COOKIE);CHKERRQ(ierr);
541   }
542 
543   PetscFunctionReturn(0);
544 }
545