xref: /petsc/src/mat/impls/aij/seq/matmatmult.c (revision dfbe8321245cdf9338262824a41580dace4e4143)
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 PetscErrorCode MatMatMult(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
108 {
109   PetscErrorCode 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 PetscErrorCode MatMatMult_MPIAIJ_MPIAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill, Mat *C)
137 {
138   PetscErrorCode 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 PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) {
154   PetscErrorCode 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 PetscErrorCode 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   PetscErrorCode 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 PetscErrorCode MatDestroy_MPIAIJ(Mat);
227 #undef __FUNCT__
228 #define __FUNCT__ "MatDestroy_MPIAIJ_MatMatMult"
229 PetscErrorCode MatDestroy_MPIAIJ_MatMatMult(Mat A)
230 {
231   PetscErrorCode 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 PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
251 {
252   Mat_MPIAIJ        *a = (Mat_MPIAIJ*)A->data;
253   PetscErrorCode ierr;
254   int               *idx,i,start,end,ncols,imark,nzA,nzB,*cmap;
255   Mat_MatMatMultMPI *mult;
256 
257   PetscFunctionBegin;
258   ierr = PetscNew(Mat_MatMatMultMPI,&mult);CHKERRQ(ierr);
259 
260   /* create a seq matrix B_seq = submatrix of B by taking rows of B that equal to nonzero col of A */
261   start = a->cstart;
262   cmap  = a->garray;
263   nzA   = a->A->n;
264   nzB   = a->B->n;
265   ierr  = PetscMalloc((nzA+nzB)*sizeof(int), &idx);CHKERRQ(ierr);
266   ncols = 0;
267   for (i=0; i<nzB; i++) {
268     if (cmap[i] < start) idx[ncols++] = cmap[i];
269     else break;
270   }
271   imark = i;
272   for (i=0; i<nzA; i++) idx[ncols++] = start + i;
273   for (i=imark; i<nzB; i++) idx[ncols++] = cmap[i];
274   ierr = ISCreateGeneral(PETSC_COMM_SELF,ncols,idx,&mult->isrowb);CHKERRQ(ierr);
275   ierr = PetscFree(idx);CHKERRQ(ierr);
276   ierr = ISCreateStride(PETSC_COMM_SELF,B->N,0,1,&mult->iscolb);CHKERRQ(ierr);
277   ierr = MatGetSubMatrices(B,1,&mult->isrowb,&mult->iscolb,MAT_INITIAL_MATRIX,&mult->bseq);CHKERRQ(ierr)
278 
279   /*  create a seq matrix A_seq = submatrix of A by taking all local rows of A */
280   start = a->rstart; end = a->rend;
281   ierr = ISCreateStride(PETSC_COMM_SELF,end-start,start,1,&mult->isrowa);CHKERRQ(ierr);
282   ierr = MatGetSubMatrices(A,1,&mult->isrowa,&mult->isrowb,MAT_INITIAL_MATRIX,&mult->aseq);CHKERRQ(ierr);
283 
284   /* compute C_seq = A_seq * B_seq */
285   ierr = MatMatMult_SeqAIJ_SeqAIJ(mult->aseq[0],mult->bseq[0],MAT_INITIAL_MATRIX,fill,&mult->C_seq);CHKERRQ(ierr);
286 
287   /* create mpi matrix C by concatinating C_seq */
288   ierr = PetscObjectReference((PetscObject)mult->C_seq);CHKERRQ(ierr); /* prevent C_seq being destroyed by MatMerge() */
289   ierr = MatMerge(A->comm,mult->C_seq,MAT_INITIAL_MATRIX,C);CHKERRQ(ierr);
290 
291   /* attach the supporting struct to C for reuse of symbolic C */
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 PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
301 {
302   PetscErrorCode 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;
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((int)(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 PetscErrorCode 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   PetscErrorCode 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 /* This routine is called ONLY in the case of reusing previously computed symbolic C */
457 #undef __FUNCT__
458 #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ"
459 PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat A,Mat B,Mat C)
460 {
461   PetscErrorCode ierr;
462   Mat_MatMatMultMPI *mult=(Mat_MatMatMultMPI*)C->spptr;
463 
464   PetscFunctionBegin;
465   ierr = MatGetSubMatrices(B,1,&mult->isrowb,&mult->iscolb,MAT_REUSE_MATRIX,&mult->bseq);CHKERRQ(ierr)
466   ierr = MatGetSubMatrices(A,1,&mult->isrowa,&mult->isrowb,MAT_REUSE_MATRIX,&mult->aseq);CHKERRQ(ierr);
467   ierr = MatMatMult_SeqAIJ_SeqAIJ(mult->aseq[0],mult->bseq[0],MAT_REUSE_MATRIX,0.0,&mult->C_seq);CHKERRQ(ierr);
468 
469   ierr = PetscObjectReference((PetscObject)mult->C_seq);CHKERRQ(ierr);
470   ierr = MatMerge(A->comm,mult->C_seq,MAT_REUSE_MATRIX,&C);CHKERRQ(ierr);
471 
472   PetscFunctionReturn(0);
473 }
474 
475 #undef __FUNCT__
476 #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqAIJ"
477 PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
478 {
479   PetscErrorCode ierr;
480   int        flops=0;
481   Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data;
482   Mat_SeqAIJ *b = (Mat_SeqAIJ *)B->data;
483   Mat_SeqAIJ *c = (Mat_SeqAIJ *)C->data;
484   int        *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j;
485   int        am=A->M,cn=C->N;
486   int        i,j,k,anzi,bnzi,cnzi,brow;
487   MatScalar  *aa=a->a,*ba=b->a,*baj,*ca=c->a,*temp;
488 
489   PetscFunctionBegin;
490 
491   /* Start timers */
492   ierr = PetscLogEventBegin(logkey_matmatmult_numeric,A,B,C,0);CHKERRQ(ierr);
493 
494   /* Allocate temp accumulation space to avoid searching for nonzero columns in C */
495   ierr = PetscMalloc((cn+1)*sizeof(MatScalar),&temp);CHKERRQ(ierr);
496   ierr = PetscMemzero(temp,cn*sizeof(MatScalar));CHKERRQ(ierr);
497   /* Traverse A row-wise. */
498   /* Build the ith row in C by summing over nonzero columns in A, */
499   /* the rows of B corresponding to nonzeros of A. */
500   for (i=0;i<am;i++) {
501     anzi = ai[i+1] - ai[i];
502     for (j=0;j<anzi;j++) {
503       brow = *aj++;
504       bnzi = bi[brow+1] - bi[brow];
505       bjj  = bj + bi[brow];
506       baj  = ba + bi[brow];
507       for (k=0;k<bnzi;k++) {
508         temp[bjj[k]] += (*aa)*baj[k];
509       }
510       flops += 2*bnzi;
511       aa++;
512     }
513     /* Store row back into C, and re-zero temp */
514     cnzi = ci[i+1] - ci[i];
515     for (j=0;j<cnzi;j++) {
516       ca[j] = temp[cj[j]];
517       temp[cj[j]] = 0.0;
518     }
519     ca += cnzi;
520     cj += cnzi;
521   }
522   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
523   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
524 
525   /* Free temp */
526   ierr = PetscFree(temp);CHKERRQ(ierr);
527   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
528   ierr = PetscLogEventEnd(logkey_matmatmult_numeric,A,B,C,0);CHKERRQ(ierr);
529   PetscFunctionReturn(0);
530 }
531 
532 #undef __FUNCT__
533 #define __FUNCT__ "RegisterMatMatMultRoutines_Private"
534 PetscErrorCode RegisterMatMatMultRoutines_Private(Mat A)
535 {
536   PetscErrorCode ierr;
537 
538   PetscFunctionBegin;
539   if (!logkey_matmatmult_symbolic) {
540     ierr = PetscLogEventRegister(&logkey_matmatmult_symbolic,"MatMatMultSymbolic",MAT_COOKIE);CHKERRQ(ierr);
541   }
542   if (!logkey_matmatmult_numeric) {
543     ierr = PetscLogEventRegister(&logkey_matmatmult_numeric,"MatMatMultNumeric",MAT_COOKIE);CHKERRQ(ierr);
544   }
545 
546   PetscFunctionReturn(0);
547 }
548