xref: /petsc/src/mat/impls/aij/seq/matmatmult.c (revision 5dc64a9727a0cc1147601b21df2b9dcd374fa7c9)
1 
2  /*
3    Defines matrix-matrix product routines for pairs of SeqAIJ matrices
4            C = A * B
5  */
6 
7  #include <../src/mat/impls/aij/seq/aij.h> /*I "petscmat.h" I*/
8  #include <../src/mat/utils/freespace.h>
9  #include <petscbt.h>
10  #include <petsc/private/isimpl.h>
11  #include <../src/mat/impls/dense/seq/dense.h>
12 
13 
14  PETSC_INTERN PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
15  {
16    PetscErrorCode ierr;
17 
18    PetscFunctionBegin;
19    if (scall == MAT_INITIAL_MATRIX) {
20      ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
21      ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
22      ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
23    }
24 
25    ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
26    ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr);
27    ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
28    PetscFunctionReturn(0);
29  }
30 
31  PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
32  {
33    PetscErrorCode ierr;
34 
35    PetscFunctionBegin;
36    if (C->ops->matmultnumeric) {
37      ierr = (*C->ops->matmultnumeric)(A,B,C);CHKERRQ(ierr);
38    } else {
39      ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ_Sorted(A,B,C);CHKERRQ(ierr);
40    }
41    PetscFunctionReturn(0);
42  }
43 
44  PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
45  {
46    PetscErrorCode ierr;
47  #if !defined(PETSC_HAVE_HYPRE)
48    const char     *algTypes[8] = {"sorted","scalable","scalable_fast","heap","btheap","llcondensed","combined","rowmerge"};
49    PetscInt       nalg = 8;
50  #else
51    const char     *algTypes[9] = {"sorted","scalable","scalable_fast","heap","btheap","llcondensed","combined","rowmerge","hypre"};
52    PetscInt       nalg = 9;
53  #endif
54    PetscInt       alg = 0; /* set default algorithm */
55 
56    PetscFunctionBegin;
57    ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"MatMatMult","Mat");CHKERRQ(ierr);
58    ierr = PetscOptionsEList("-matmatmult_via","Algorithmic approach","MatMatMult",algTypes,nalg,algTypes[0],&alg,NULL);CHKERRQ(ierr);
59    ierr = PetscOptionsEnd();CHKERRQ(ierr);
60    switch (alg) {
61    case 1:
62      ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(A,B,fill,C);CHKERRQ(ierr);
63      break;
64    case 2:
65      ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(A,B,fill,C);CHKERRQ(ierr);
66      break;
67    case 3:
68      ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(A,B,fill,C);CHKERRQ(ierr);
69      break;
70    case 4:
71      ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(A,B,fill,C);CHKERRQ(ierr);
72      break;
73    case 5:
74      ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_LLCondensed(A,B,fill,C);CHKERRQ(ierr);
75      break;
76    case 6:
77      ierr = MatMatMult_SeqAIJ_SeqAIJ_Combined(A,B,fill,C);CHKERRQ(ierr);
78      break;
79    case 7:
80      ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_RowMerge(A,B,fill,C);CHKERRQ(ierr);
81      break;
82  #if defined(PETSC_HAVE_HYPRE)
83    case 8:
84      ierr = MatMatMultSymbolic_AIJ_AIJ_wHYPRE(A,B,fill,C);CHKERRQ(ierr);
85      break;
86  #endif
87    default:
88      ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Sorted(A,B,fill,C);CHKERRQ(ierr);
89      break;
90    }
91    PetscFunctionReturn(0);
92  }
93 
94  PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_LLCondensed(Mat A,Mat B,PetscReal fill,Mat *C)
95  {
96    PetscErrorCode     ierr;
97    Mat_SeqAIJ         *a =(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
98    PetscInt           *ai=a->i,*bi=b->i,*ci,*cj;
99    PetscInt           am =A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
100    PetscReal          afill;
101    PetscInt           i,j,anzi,brow,bnzj,cnzi,*bj,*aj,*lnk,ndouble=0,Crmax;
102    PetscTable         ta;
103    PetscBT            lnkbt;
104    PetscFreeSpaceList free_space=NULL,current_space=NULL;
105 
106    PetscFunctionBegin;
107    /* Get ci and cj */
108    /*---------------*/
109    /* Allocate ci array, arrays for fill computation and */
110    /* free space for accumulating nonzero column info */
111    ierr  = PetscMalloc1(am+2,&ci);CHKERRQ(ierr);
112    ci[0] = 0;
113 
114    /* create and initialize a linked list */
115    ierr = PetscTableCreate(bn,bn,&ta);CHKERRQ(ierr);
116    MatRowMergeMax_SeqAIJ(b,bm,ta);
117    ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr);
118    ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
119 
120    ierr = PetscLLCondensedCreate(Crmax,bn,&lnk,&lnkbt);CHKERRQ(ierr);
121 
122    /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
123    ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am],bi[bm])),&free_space);CHKERRQ(ierr);
124 
125    current_space = free_space;
126 
127    /* Determine ci and cj */
128    for (i=0; i<am; i++) {
129      anzi = ai[i+1] - ai[i];
130      aj   = a->j + ai[i];
131      for (j=0; j<anzi; j++) {
132        brow = aj[j];
133        bnzj = bi[brow+1] - bi[brow];
134        bj   = b->j + bi[brow];
135        /* add non-zero cols of B into the sorted linked list lnk */
136        ierr = PetscLLCondensedAddSorted(bnzj,bj,lnk,lnkbt);CHKERRQ(ierr);
137      }
138      cnzi = lnk[0];
139 
140      /* If free space is not available, make more free space */
141      /* Double the amount of total space in the list */
142      if (current_space->local_remaining<cnzi) {
143        ierr = PetscFreeSpaceGet(PetscIntSumTruncate(cnzi,current_space->total_array_size),&current_space);CHKERRQ(ierr);
144        ndouble++;
145      }
146 
147      /* Copy data into free space, then initialize lnk */
148      ierr = PetscLLCondensedClean(bn,cnzi,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
149 
150      current_space->array           += cnzi;
151      current_space->local_used      += cnzi;
152      current_space->local_remaining -= cnzi;
153 
154      ci[i+1] = ci[i] + cnzi;
155    }
156 
157    /* Column indices are in the list of free space */
158    /* Allocate space for cj, initialize cj, and */
159    /* destroy list of free space and other temporary array(s) */
160    ierr = PetscMalloc1(ci[am]+1,&cj);CHKERRQ(ierr);
161    ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
162    ierr = PetscLLCondensedDestroy(lnk,lnkbt);CHKERRQ(ierr);
163 
164    /* put together the new symbolic matrix */
165    ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,C);CHKERRQ(ierr);
166    ierr = MatSetBlockSizesFromMats(*C,A,B);CHKERRQ(ierr);
167    ierr = MatSetType(*C,((PetscObject)A)->type_name);CHKERRQ(ierr);
168 
169   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
170   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
171   c                         = (Mat_SeqAIJ*)((*C)->data);
172   c->free_a                 = PETSC_FALSE;
173   c->free_ij                = PETSC_TRUE;
174   c->nonew                  = 0;
175   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Sorted; /* fast, needs non-scalable O(bn) array 'abdense' */
176 
177   /* set MatInfo */
178   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
179   if (afill < 1.0) afill = 1.0;
180   c->maxnz                     = ci[am];
181   c->nz                        = ci[am];
182   (*C)->info.mallocs           = ndouble;
183   (*C)->info.fill_ratio_given  = fill;
184   (*C)->info.fill_ratio_needed = afill;
185 
186 #if defined(PETSC_USE_INFO)
187   if (ci[am]) {
188     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);CHKERRQ(ierr);
189     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
190   } else {
191     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
192   }
193 #endif
194   PetscFunctionReturn(0);
195 }
196 
197 PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Sorted(Mat A,Mat B,Mat C)
198 {
199   PetscErrorCode ierr;
200   PetscLogDouble flops=0.0;
201   Mat_SeqAIJ     *a   = (Mat_SeqAIJ*)A->data;
202   Mat_SeqAIJ     *b   = (Mat_SeqAIJ*)B->data;
203   Mat_SeqAIJ     *c   = (Mat_SeqAIJ*)C->data;
204   PetscInt       *ai  =a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j;
205   PetscInt       am   =A->rmap->n,cm=C->rmap->n;
206   PetscInt       i,j,k,anzi,bnzi,cnzi,brow;
207   PetscScalar    *aa=a->a,*ba=b->a,*baj,*ca,valtmp;
208   PetscScalar    *ab_dense;
209 
210   PetscFunctionBegin;
211   if (!c->a) { /* first call of MatMatMultNumeric_SeqAIJ_SeqAIJ, allocate ca and matmult_abdense */
212     ierr      = PetscMalloc1(ci[cm]+1,&ca);CHKERRQ(ierr);
213     c->a      = ca;
214     c->free_a = PETSC_TRUE;
215   } else {
216     ca        = c->a;
217   }
218   if (!c->matmult_abdense) {
219     ierr = PetscCalloc1(B->cmap->N,&ab_dense);CHKERRQ(ierr);
220     c->matmult_abdense = ab_dense;
221   } else {
222     ab_dense = c->matmult_abdense;
223   }
224 
225   /* clean old values in C */
226   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
227   /* Traverse A row-wise. */
228   /* Build the ith row in C by summing over nonzero columns in A, */
229   /* the rows of B corresponding to nonzeros of A. */
230   for (i=0; i<am; i++) {
231     anzi = ai[i+1] - ai[i];
232     for (j=0; j<anzi; j++) {
233       brow = aj[j];
234       bnzi = bi[brow+1] - bi[brow];
235       bjj  = bj + bi[brow];
236       baj  = ba + bi[brow];
237       /* perform dense axpy */
238       valtmp = aa[j];
239       for (k=0; k<bnzi; k++) {
240         ab_dense[bjj[k]] += valtmp*baj[k];
241       }
242       flops += 2*bnzi;
243     }
244     aj += anzi; aa += anzi;
245 
246     cnzi = ci[i+1] - ci[i];
247     for (k=0; k<cnzi; k++) {
248       ca[k]          += ab_dense[cj[k]];
249       ab_dense[cj[k]] = 0.0; /* zero ab_dense */
250     }
251     flops += cnzi;
252     cj    += cnzi; ca += cnzi;
253   }
254   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
255   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
256   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
257   PetscFunctionReturn(0);
258 }
259 
260 PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(Mat A,Mat B,Mat C)
261 {
262   PetscErrorCode ierr;
263   PetscLogDouble flops=0.0;
264   Mat_SeqAIJ     *a   = (Mat_SeqAIJ*)A->data;
265   Mat_SeqAIJ     *b   = (Mat_SeqAIJ*)B->data;
266   Mat_SeqAIJ     *c   = (Mat_SeqAIJ*)C->data;
267   PetscInt       *ai  = a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j;
268   PetscInt       am   = A->rmap->N,cm=C->rmap->N;
269   PetscInt       i,j,k,anzi,bnzi,cnzi,brow;
270   PetscScalar    *aa=a->a,*ba=b->a,*baj,*ca=c->a,valtmp;
271   PetscInt       nextb;
272 
273   PetscFunctionBegin;
274   if (!ca) { /* first call of MatMatMultNumeric_SeqAIJ_SeqAIJ, allocate ca and matmult_abdense */
275     ierr      = PetscMalloc1(ci[cm]+1,&ca);CHKERRQ(ierr);
276     c->a      = ca;
277     c->free_a = PETSC_TRUE;
278   }
279 
280   /* clean old values in C */
281   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
282   /* Traverse A row-wise. */
283   /* Build the ith row in C by summing over nonzero columns in A, */
284   /* the rows of B corresponding to nonzeros of A. */
285   for (i=0; i<am; i++) {
286     anzi = ai[i+1] - ai[i];
287     cnzi = ci[i+1] - ci[i];
288     for (j=0; j<anzi; j++) {
289       brow = aj[j];
290       bnzi = bi[brow+1] - bi[brow];
291       bjj  = bj + bi[brow];
292       baj  = ba + bi[brow];
293       /* perform sparse axpy */
294       valtmp = aa[j];
295       nextb  = 0;
296       for (k=0; nextb<bnzi; k++) {
297         if (cj[k] == bjj[nextb]) { /* ccol == bcol */
298           ca[k] += valtmp*baj[nextb++];
299         }
300       }
301       flops += 2*bnzi;
302     }
303     aj += anzi; aa += anzi;
304     cj += cnzi; ca += cnzi;
305   }
306 
307   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
308   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
309   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
310   PetscFunctionReturn(0);
311 }
312 
313 PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(Mat A,Mat B,PetscReal fill,Mat *C)
314 {
315   PetscErrorCode     ierr;
316   Mat_SeqAIJ         *a  = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
317   PetscInt           *ai = a->i,*bi=b->i,*ci,*cj;
318   PetscInt           am  = A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
319   MatScalar          *ca;
320   PetscReal          afill;
321   PetscInt           i,j,anzi,brow,bnzj,cnzi,*bj,*aj,*lnk,ndouble=0,Crmax;
322   PetscTable         ta;
323   PetscFreeSpaceList free_space=NULL,current_space=NULL;
324 
325   PetscFunctionBegin;
326   /* Get ci and cj - same as MatMatMultSymbolic_SeqAIJ_SeqAIJ except using PetscLLxxx_fast() */
327   /*-----------------------------------------------------------------------------------------*/
328   /* Allocate arrays for fill computation and free space for accumulating nonzero column */
329   ierr  = PetscMalloc1(am+2,&ci);CHKERRQ(ierr);
330   ci[0] = 0;
331 
332   /* create and initialize a linked list */
333   ierr = PetscTableCreate(bn,bn,&ta);CHKERRQ(ierr);
334   MatRowMergeMax_SeqAIJ(b,bm,ta);
335   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr);
336   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
337 
338   ierr = PetscLLCondensedCreate_fast(Crmax,&lnk);CHKERRQ(ierr);
339 
340   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
341   ierr          = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am],bi[bm])),&free_space);CHKERRQ(ierr);
342   current_space = free_space;
343 
344   /* Determine ci and cj */
345   for (i=0; i<am; i++) {
346     anzi = ai[i+1] - ai[i];
347     aj   = a->j + ai[i];
348     for (j=0; j<anzi; j++) {
349       brow = aj[j];
350       bnzj = bi[brow+1] - bi[brow];
351       bj   = b->j + bi[brow];
352       /* add non-zero cols of B into the sorted linked list lnk */
353       ierr = PetscLLCondensedAddSorted_fast(bnzj,bj,lnk);CHKERRQ(ierr);
354     }
355     cnzi = lnk[1];
356 
357     /* If free space is not available, make more free space */
358     /* Double the amount of total space in the list */
359     if (current_space->local_remaining<cnzi) {
360       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(cnzi,current_space->total_array_size),&current_space);CHKERRQ(ierr);
361       ndouble++;
362     }
363 
364     /* Copy data into free space, then initialize lnk */
365     ierr = PetscLLCondensedClean_fast(cnzi,current_space->array,lnk);CHKERRQ(ierr);
366 
367     current_space->array           += cnzi;
368     current_space->local_used      += cnzi;
369     current_space->local_remaining -= cnzi;
370 
371     ci[i+1] = ci[i] + cnzi;
372   }
373 
374   /* Column indices are in the list of free space */
375   /* Allocate space for cj, initialize cj, and */
376   /* destroy list of free space and other temporary array(s) */
377   ierr = PetscMalloc1(ci[am]+1,&cj);CHKERRQ(ierr);
378   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
379   ierr = PetscLLCondensedDestroy_fast(lnk);CHKERRQ(ierr);
380 
381   /* Allocate space for ca */
382   ierr = PetscMalloc1(ci[am]+1,&ca);CHKERRQ(ierr);
383   ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr);
384 
385   /* put together the new symbolic matrix */
386   ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,ca,C);CHKERRQ(ierr);
387   ierr = MatSetBlockSizesFromMats(*C,A,B);CHKERRQ(ierr);
388   ierr = MatSetType(*C,((PetscObject)A)->type_name);CHKERRQ(ierr);
389 
390   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
391   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
392   c          = (Mat_SeqAIJ*)((*C)->data);
393   c->free_a  = PETSC_TRUE;
394   c->free_ij = PETSC_TRUE;
395   c->nonew   = 0;
396 
397   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable; /* slower, less memory */
398 
399   /* set MatInfo */
400   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
401   if (afill < 1.0) afill = 1.0;
402   c->maxnz                     = ci[am];
403   c->nz                        = ci[am];
404   (*C)->info.mallocs           = ndouble;
405   (*C)->info.fill_ratio_given  = fill;
406   (*C)->info.fill_ratio_needed = afill;
407 
408 #if defined(PETSC_USE_INFO)
409   if (ci[am]) {
410     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);CHKERRQ(ierr);
411     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
412   } else {
413     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
414   }
415 #endif
416   PetscFunctionReturn(0);
417 }
418 
419 
420 PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(Mat A,Mat B,PetscReal fill,Mat *C)
421 {
422   PetscErrorCode     ierr;
423   Mat_SeqAIJ         *a  = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
424   PetscInt           *ai = a->i,*bi=b->i,*ci,*cj;
425   PetscInt           am  = A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
426   MatScalar          *ca;
427   PetscReal          afill;
428   PetscInt           i,j,anzi,brow,bnzj,cnzi,*bj,*aj,*lnk,ndouble=0,Crmax;
429   PetscTable         ta;
430   PetscFreeSpaceList free_space=NULL,current_space=NULL;
431 
432   PetscFunctionBegin;
433   /* Get ci and cj - same as MatMatMultSymbolic_SeqAIJ_SeqAIJ except using PetscLLxxx_Scalalbe() */
434   /*---------------------------------------------------------------------------------------------*/
435   /* Allocate arrays for fill computation and free space for accumulating nonzero column */
436   ierr  = PetscMalloc1(am+2,&ci);CHKERRQ(ierr);
437   ci[0] = 0;
438 
439   /* create and initialize a linked list */
440   ierr = PetscTableCreate(bn,bn,&ta);CHKERRQ(ierr);
441   MatRowMergeMax_SeqAIJ(b,bm,ta);
442   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr);
443   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
444   ierr = PetscLLCondensedCreate_Scalable(Crmax,&lnk);CHKERRQ(ierr);
445 
446   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
447   ierr          = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am],bi[bm])),&free_space);CHKERRQ(ierr);
448   current_space = free_space;
449 
450   /* Determine ci and cj */
451   for (i=0; i<am; i++) {
452     anzi = ai[i+1] - ai[i];
453     aj   = a->j + ai[i];
454     for (j=0; j<anzi; j++) {
455       brow = aj[j];
456       bnzj = bi[brow+1] - bi[brow];
457       bj   = b->j + bi[brow];
458       /* add non-zero cols of B into the sorted linked list lnk */
459       ierr = PetscLLCondensedAddSorted_Scalable(bnzj,bj,lnk);CHKERRQ(ierr);
460     }
461     cnzi = lnk[0];
462 
463     /* If free space is not available, make more free space */
464     /* Double the amount of total space in the list */
465     if (current_space->local_remaining<cnzi) {
466       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(cnzi,current_space->total_array_size),&current_space);CHKERRQ(ierr);
467       ndouble++;
468     }
469 
470     /* Copy data into free space, then initialize lnk */
471     ierr = PetscLLCondensedClean_Scalable(cnzi,current_space->array,lnk);CHKERRQ(ierr);
472 
473     current_space->array           += cnzi;
474     current_space->local_used      += cnzi;
475     current_space->local_remaining -= cnzi;
476 
477     ci[i+1] = ci[i] + cnzi;
478   }
479 
480   /* Column indices are in the list of free space */
481   /* Allocate space for cj, initialize cj, and */
482   /* destroy list of free space and other temporary array(s) */
483   ierr = PetscMalloc1(ci[am]+1,&cj);CHKERRQ(ierr);
484   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
485   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
486 
487   /* Allocate space for ca */
488   /*-----------------------*/
489   ierr = PetscMalloc1(ci[am]+1,&ca);CHKERRQ(ierr);
490   ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr);
491 
492   /* put together the new symbolic matrix */
493   ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,ca,C);CHKERRQ(ierr);
494   ierr = MatSetBlockSizesFromMats(*C,A,B);CHKERRQ(ierr);
495   ierr = MatSetType(*C,((PetscObject)A)->type_name);CHKERRQ(ierr);
496 
497   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
498   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
499   c          = (Mat_SeqAIJ*)((*C)->data);
500   c->free_a  = PETSC_TRUE;
501   c->free_ij = PETSC_TRUE;
502   c->nonew   = 0;
503 
504   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable; /* slower, less memory */
505 
506   /* set MatInfo */
507   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
508   if (afill < 1.0) afill = 1.0;
509   c->maxnz                     = ci[am];
510   c->nz                        = ci[am];
511   (*C)->info.mallocs           = ndouble;
512   (*C)->info.fill_ratio_given  = fill;
513   (*C)->info.fill_ratio_needed = afill;
514 
515 #if defined(PETSC_USE_INFO)
516   if (ci[am]) {
517     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);CHKERRQ(ierr);
518     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
519   } else {
520     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
521   }
522 #endif
523   PetscFunctionReturn(0);
524 }
525 
526 PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(Mat A,Mat B,PetscReal fill,Mat *C)
527 {
528   PetscErrorCode     ierr;
529   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
530   const PetscInt     *ai=a->i,*bi=b->i,*aj=a->j,*bj=b->j;
531   PetscInt           *ci,*cj,*bb;
532   PetscInt           am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
533   PetscReal          afill;
534   PetscInt           i,j,col,ndouble = 0;
535   PetscFreeSpaceList free_space=NULL,current_space=NULL;
536   PetscHeap          h;
537 
538   PetscFunctionBegin;
539   /* Get ci and cj - by merging sorted rows using a heap */
540   /*---------------------------------------------------------------------------------------------*/
541   /* Allocate arrays for fill computation and free space for accumulating nonzero column */
542   ierr  = PetscMalloc1(am+2,&ci);CHKERRQ(ierr);
543   ci[0] = 0;
544 
545   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
546   ierr          = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am],bi[bm])),&free_space);CHKERRQ(ierr);
547   current_space = free_space;
548 
549   ierr = PetscHeapCreate(a->rmax,&h);CHKERRQ(ierr);
550   ierr = PetscMalloc1(a->rmax,&bb);CHKERRQ(ierr);
551 
552   /* Determine ci and cj */
553   for (i=0; i<am; i++) {
554     const PetscInt anzi  = ai[i+1] - ai[i]; /* number of nonzeros in this row of A, this is the number of rows of B that we merge */
555     const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */
556     ci[i+1] = ci[i];
557     /* Populate the min heap */
558     for (j=0; j<anzi; j++) {
559       bb[j] = bi[acol[j]];         /* bb points at the start of the row */
560       if (bb[j] < bi[acol[j]+1]) { /* Add if row is nonempty */
561         ierr = PetscHeapAdd(h,j,bj[bb[j]++]);CHKERRQ(ierr);
562       }
563     }
564     /* Pick off the min element, adding it to free space */
565     ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr);
566     while (j >= 0) {
567       if (current_space->local_remaining < 1) { /* double the size, but don't exceed 16 MiB */
568         ierr = PetscFreeSpaceGet(PetscMin(PetscIntMultTruncate(2,current_space->total_array_size),16 << 20),&current_space);CHKERRQ(ierr);
569         ndouble++;
570       }
571       *(current_space->array++) = col;
572       current_space->local_used++;
573       current_space->local_remaining--;
574       ci[i+1]++;
575 
576       /* stash if anything else remains in this row of B */
577       if (bb[j] < bi[acol[j]+1]) {ierr = PetscHeapStash(h,j,bj[bb[j]++]);CHKERRQ(ierr);}
578       while (1) {               /* pop and stash any other rows of B that also had an entry in this column */
579         PetscInt j2,col2;
580         ierr = PetscHeapPeek(h,&j2,&col2);CHKERRQ(ierr);
581         if (col2 != col) break;
582         ierr = PetscHeapPop(h,&j2,&col2);CHKERRQ(ierr);
583         if (bb[j2] < bi[acol[j2]+1]) {ierr = PetscHeapStash(h,j2,bj[bb[j2]++]);CHKERRQ(ierr);}
584       }
585       /* Put any stashed elements back into the min heap */
586       ierr = PetscHeapUnstash(h);CHKERRQ(ierr);
587       ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr);
588     }
589   }
590   ierr = PetscFree(bb);CHKERRQ(ierr);
591   ierr = PetscHeapDestroy(&h);CHKERRQ(ierr);
592 
593   /* Column indices are in the list of free space */
594   /* Allocate space for cj, initialize cj, and */
595   /* destroy list of free space and other temporary array(s) */
596   ierr = PetscMalloc1(ci[am],&cj);CHKERRQ(ierr);
597   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
598 
599   /* put together the new symbolic matrix */
600   ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,C);CHKERRQ(ierr);
601   ierr = MatSetBlockSizesFromMats(*C,A,B);CHKERRQ(ierr);
602   ierr = MatSetType(*C,((PetscObject)A)->type_name);CHKERRQ(ierr);
603 
604   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
605   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
606   c          = (Mat_SeqAIJ*)((*C)->data);
607   c->free_a  = PETSC_TRUE;
608   c->free_ij = PETSC_TRUE;
609   c->nonew   = 0;
610 
611   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Sorted;
612 
613   /* set MatInfo */
614   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
615   if (afill < 1.0) afill = 1.0;
616   c->maxnz                     = ci[am];
617   c->nz                        = ci[am];
618   (*C)->info.mallocs           = ndouble;
619   (*C)->info.fill_ratio_given  = fill;
620   (*C)->info.fill_ratio_needed = afill;
621 
622 #if defined(PETSC_USE_INFO)
623   if (ci[am]) {
624     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);CHKERRQ(ierr);
625     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
626   } else {
627     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
628   }
629 #endif
630   PetscFunctionReturn(0);
631 }
632 
633 PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(Mat A,Mat B,PetscReal fill,Mat *C)
634 {
635   PetscErrorCode     ierr;
636   Mat_SeqAIJ         *a  = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
637   const PetscInt     *ai = a->i,*bi=b->i,*aj=a->j,*bj=b->j;
638   PetscInt           *ci,*cj,*bb;
639   PetscInt           am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
640   PetscReal          afill;
641   PetscInt           i,j,col,ndouble = 0;
642   PetscFreeSpaceList free_space=NULL,current_space=NULL;
643   PetscHeap          h;
644   PetscBT            bt;
645 
646   PetscFunctionBegin;
647   /* Get ci and cj - using a heap for the sorted rows, but use BT so that each index is only added once */
648   /*---------------------------------------------------------------------------------------------*/
649   /* Allocate arrays for fill computation and free space for accumulating nonzero column */
650   ierr  = PetscMalloc1(am+2,&ci);CHKERRQ(ierr);
651   ci[0] = 0;
652 
653   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
654   ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am],bi[bm])),&free_space);CHKERRQ(ierr);
655 
656   current_space = free_space;
657 
658   ierr = PetscHeapCreate(a->rmax,&h);CHKERRQ(ierr);
659   ierr = PetscMalloc1(a->rmax,&bb);CHKERRQ(ierr);
660   ierr = PetscBTCreate(bn,&bt);CHKERRQ(ierr);
661 
662   /* Determine ci and cj */
663   for (i=0; i<am; i++) {
664     const PetscInt anzi  = ai[i+1] - ai[i]; /* number of nonzeros in this row of A, this is the number of rows of B that we merge */
665     const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */
666     const PetscInt *fptr = current_space->array; /* Save beginning of the row so we can clear the BT later */
667     ci[i+1] = ci[i];
668     /* Populate the min heap */
669     for (j=0; j<anzi; j++) {
670       PetscInt brow = acol[j];
671       for (bb[j] = bi[brow]; bb[j] < bi[brow+1]; bb[j]++) {
672         PetscInt bcol = bj[bb[j]];
673         if (!PetscBTLookupSet(bt,bcol)) { /* new entry */
674           ierr = PetscHeapAdd(h,j,bcol);CHKERRQ(ierr);
675           bb[j]++;
676           break;
677         }
678       }
679     }
680     /* Pick off the min element, adding it to free space */
681     ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr);
682     while (j >= 0) {
683       if (current_space->local_remaining < 1) { /* double the size, but don't exceed 16 MiB */
684         fptr = NULL;                      /* need PetscBTMemzero */
685         ierr = PetscFreeSpaceGet(PetscMin(PetscIntMultTruncate(2,current_space->total_array_size),16 << 20),&current_space);CHKERRQ(ierr);
686         ndouble++;
687       }
688       *(current_space->array++) = col;
689       current_space->local_used++;
690       current_space->local_remaining--;
691       ci[i+1]++;
692 
693       /* stash if anything else remains in this row of B */
694       for (; bb[j] < bi[acol[j]+1]; bb[j]++) {
695         PetscInt bcol = bj[bb[j]];
696         if (!PetscBTLookupSet(bt,bcol)) { /* new entry */
697           ierr = PetscHeapAdd(h,j,bcol);CHKERRQ(ierr);
698           bb[j]++;
699           break;
700         }
701       }
702       ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr);
703     }
704     if (fptr) {                 /* Clear the bits for this row */
705       for (; fptr<current_space->array; fptr++) {ierr = PetscBTClear(bt,*fptr);CHKERRQ(ierr);}
706     } else {                    /* We reallocated so we don't remember (easily) how to clear only the bits we changed */
707       ierr = PetscBTMemzero(bn,bt);CHKERRQ(ierr);
708     }
709   }
710   ierr = PetscFree(bb);CHKERRQ(ierr);
711   ierr = PetscHeapDestroy(&h);CHKERRQ(ierr);
712   ierr = PetscBTDestroy(&bt);CHKERRQ(ierr);
713 
714   /* Column indices are in the list of free space */
715   /* Allocate space for cj, initialize cj, and */
716   /* destroy list of free space and other temporary array(s) */
717   ierr = PetscMalloc1(ci[am],&cj);CHKERRQ(ierr);
718   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
719 
720   /* put together the new symbolic matrix */
721   ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,C);CHKERRQ(ierr);
722   ierr = MatSetBlockSizesFromMats(*C,A,B);CHKERRQ(ierr);
723   ierr = MatSetType(*C,((PetscObject)A)->type_name);CHKERRQ(ierr);
724 
725   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
726   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
727   c          = (Mat_SeqAIJ*)((*C)->data);
728   c->free_a  = PETSC_TRUE;
729   c->free_ij = PETSC_TRUE;
730   c->nonew   = 0;
731 
732   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Sorted;
733 
734   /* set MatInfo */
735   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
736   if (afill < 1.0) afill = 1.0;
737   c->maxnz                     = ci[am];
738   c->nz                        = ci[am];
739   (*C)->info.mallocs           = ndouble;
740   (*C)->info.fill_ratio_given  = fill;
741   (*C)->info.fill_ratio_needed = afill;
742 
743 #if defined(PETSC_USE_INFO)
744   if (ci[am]) {
745     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);CHKERRQ(ierr);
746     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
747   } else {
748     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
749   }
750 #endif
751   PetscFunctionReturn(0);
752 }
753 
754 
755 PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_RowMerge(Mat A,Mat B,PetscReal fill,Mat *C)
756 {
757   PetscErrorCode     ierr;
758   Mat_SeqAIJ         *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
759   const PetscInt     *ai=a->i,*bi=b->i,*aj=a->j,*bj=b->j,*inputi,*inputj,*inputcol,*inputcol_L1;
760   PetscInt           *ci,*cj,*outputj,worki_L1[9],worki_L2[9];
761   PetscInt           c_maxmem,a_maxrownnz=0,a_rownnz;
762   const PetscInt     workcol[8]={0,1,2,3,4,5,6,7};
763   const PetscInt     am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
764   const PetscInt     *brow_ptr[8],*brow_end[8];
765   PetscInt           window[8];
766   PetscInt           window_min,old_window_min,ci_nnz,outputi_nnz=0,L1_nrows,L2_nrows;
767   PetscInt           i,k,ndouble=0,L1_rowsleft,rowsleft;
768   PetscReal          afill;
769   PetscInt           *workj_L1,*workj_L2,*workj_L3;
770   PetscInt           L1_nnz,L2_nnz;
771 
772   /* Step 1: Get upper bound on memory required for allocation.
773              Because of the way virtual memory works,
774              only the memory pages that are actually needed will be physically allocated. */
775   PetscFunctionBegin;
776   ierr = PetscMalloc1(am+1,&ci);CHKERRQ(ierr);
777   for (i=0; i<am; i++) {
778     const PetscInt anzi  = ai[i+1] - ai[i]; /* number of nonzeros in this row of A, this is the number of rows of B that we merge */
779     const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */
780     a_rownnz = 0;
781     for (k=0; k<anzi; ++k) {
782       a_rownnz += bi[acol[k]+1] - bi[acol[k]];
783       if (a_rownnz > bn) {
784         a_rownnz = bn;
785         break;
786       }
787     }
788     a_maxrownnz = PetscMax(a_maxrownnz, a_rownnz);
789   }
790   /* temporary work areas for merging rows */
791   ierr = PetscMalloc1(a_maxrownnz*8,&workj_L1);CHKERRQ(ierr);
792   ierr = PetscMalloc1(a_maxrownnz*8,&workj_L2);CHKERRQ(ierr);
793   ierr = PetscMalloc1(a_maxrownnz,&workj_L3);CHKERRQ(ierr);
794 
795   /* This should be enough for almost all matrices. If not, memory is reallocated later. */
796   c_maxmem = 8*(ai[am]+bi[bm]);
797   /* Step 2: Populate pattern for C */
798   ierr  = PetscMalloc1(c_maxmem,&cj);CHKERRQ(ierr);
799 
800   ci_nnz       = 0;
801   ci[0]        = 0;
802   worki_L1[0]  = 0;
803   worki_L2[0]  = 0;
804   for (i=0; i<am; i++) {
805     const PetscInt anzi  = ai[i+1] - ai[i]; /* number of nonzeros in this row of A, this is the number of rows of B that we merge */
806     const PetscInt *acol = aj + ai[i];      /* column indices of nonzero entries in this row */
807     rowsleft             = anzi;
808     inputcol_L1          = acol;
809     L2_nnz               = 0;
810     L2_nrows             = 1;  /* Number of rows to be merged on Level 3. output of L3 already exists -> initial value 1   */
811     worki_L2[1]          = 0;
812     outputi_nnz          = 0;
813 
814     /* If the number of indices in C so far + the max number of columns in the next row > c_maxmem  -> allocate more memory */
815     while (ci_nnz+a_maxrownnz > c_maxmem) {
816       c_maxmem *= 2;
817       ndouble++;
818       ierr = PetscRealloc(sizeof(PetscInt)*c_maxmem,&cj);CHKERRQ(ierr);
819     }
820 
821     while (rowsleft) {
822       L1_rowsleft = PetscMin(64, rowsleft); /* In the inner loop max 64 rows of B can be merged */
823       L1_nrows    = 0;
824       L1_nnz      = 0;
825       inputcol    = inputcol_L1;
826       inputi      = bi;
827       inputj      = bj;
828 
829       /* The following macro is used to specialize for small rows in A.
830          This helps with compiler unrolling, improving performance substantially.
831           Input:  inputj   inputi  inputcol  bn
832           Output: outputj  outputi_nnz                       */
833        #define MatMatMultSymbolic_RowMergeMacro(ANNZ)                        \
834          window_min  = bn;                                                   \
835          outputi_nnz = 0;                                                    \
836          for (k=0; k<ANNZ; ++k) {                                            \
837            brow_ptr[k] = inputj + inputi[inputcol[k]];                       \
838            brow_end[k] = inputj + inputi[inputcol[k]+1];                     \
839            window[k]   = (brow_ptr[k] != brow_end[k]) ? *brow_ptr[k] : bn;   \
840            window_min  = PetscMin(window[k], window_min);                    \
841          }                                                                   \
842          while (window_min < bn) {                                           \
843            outputj[outputi_nnz++] = window_min;                              \
844            /* advance front and compute new minimum */                       \
845            old_window_min = window_min;                                      \
846            window_min = bn;                                                  \
847            for (k=0; k<ANNZ; ++k) {                                          \
848              if (window[k] == old_window_min) {                              \
849                brow_ptr[k]++;                                                \
850                window[k] = (brow_ptr[k] != brow_end[k]) ? *brow_ptr[k] : bn; \
851              }                                                               \
852              window_min = PetscMin(window[k], window_min);                   \
853            }                                                                 \
854          }
855 
856       /************** L E V E L  1 ***************/
857       /* Merge up to 8 rows of B to L1 work array*/
858       while (L1_rowsleft) {
859         outputi_nnz = 0;
860         if (anzi > 8)  outputj = workj_L1 + L1_nnz;     /* Level 1 rowmerge*/
861         else           outputj = cj + ci_nnz;           /* Merge directly to C */
862 
863         switch (L1_rowsleft) {
864         case 1:  brow_ptr[0] = inputj + inputi[inputcol[0]];
865                  brow_end[0] = inputj + inputi[inputcol[0]+1];
866                  for (; brow_ptr[0] != brow_end[0]; ++brow_ptr[0]) outputj[outputi_nnz++] = *brow_ptr[0]; /* copy row in b over */
867                  inputcol    += L1_rowsleft;
868                  rowsleft    -= L1_rowsleft;
869                  L1_rowsleft  = 0;
870                  break;
871         case 2:  MatMatMultSymbolic_RowMergeMacro(2);
872                  inputcol    += L1_rowsleft;
873                  rowsleft    -= L1_rowsleft;
874                  L1_rowsleft  = 0;
875                  break;
876         case 3: MatMatMultSymbolic_RowMergeMacro(3);
877                  inputcol    += L1_rowsleft;
878                  rowsleft    -= L1_rowsleft;
879                  L1_rowsleft  = 0;
880                  break;
881         case 4:  MatMatMultSymbolic_RowMergeMacro(4);
882                  inputcol    += L1_rowsleft;
883                  rowsleft    -= L1_rowsleft;
884                  L1_rowsleft  = 0;
885                  break;
886         case 5:  MatMatMultSymbolic_RowMergeMacro(5);
887                  inputcol    += L1_rowsleft;
888                  rowsleft    -= L1_rowsleft;
889                  L1_rowsleft  = 0;
890                  break;
891         case 6:  MatMatMultSymbolic_RowMergeMacro(6);
892                  inputcol    += L1_rowsleft;
893                  rowsleft    -= L1_rowsleft;
894                  L1_rowsleft  = 0;
895                  break;
896         case 7:  MatMatMultSymbolic_RowMergeMacro(7);
897                  inputcol    += L1_rowsleft;
898                  rowsleft    -= L1_rowsleft;
899                  L1_rowsleft  = 0;
900                  break;
901         default: MatMatMultSymbolic_RowMergeMacro(8);
902                  inputcol    += 8;
903                  rowsleft    -= 8;
904                  L1_rowsleft -= 8;
905                  break;
906         }
907         inputcol_L1           = inputcol;
908         L1_nnz               += outputi_nnz;
909         worki_L1[++L1_nrows]  = L1_nnz;
910       }
911 
912       /********************** L E V E L  2 ************************/
913       /* Merge from L1 work array to either C or to L2 work array */
914       if (anzi > 8) {
915         inputi      = worki_L1;
916         inputj      = workj_L1;
917         inputcol    = workcol;
918         outputi_nnz = 0;
919 
920         if (anzi <= 64) outputj = cj + ci_nnz;        /* Merge from L1 work array to C */
921         else            outputj = workj_L2 + L2_nnz;  /* Merge from L1 work array to L2 work array */
922 
923         switch (L1_nrows) {
924         case 1:  brow_ptr[0] = inputj + inputi[inputcol[0]];
925                  brow_end[0] = inputj + inputi[inputcol[0]+1];
926                  for (; brow_ptr[0] != brow_end[0]; ++brow_ptr[0]) outputj[outputi_nnz++] = *brow_ptr[0]; /* copy row in b over */
927                  break;
928         case 2:  MatMatMultSymbolic_RowMergeMacro(2); break;
929         case 3:  MatMatMultSymbolic_RowMergeMacro(3); break;
930         case 4:  MatMatMultSymbolic_RowMergeMacro(4); break;
931         case 5:  MatMatMultSymbolic_RowMergeMacro(5); break;
932         case 6:  MatMatMultSymbolic_RowMergeMacro(6); break;
933         case 7:  MatMatMultSymbolic_RowMergeMacro(7); break;
934         case 8:  MatMatMultSymbolic_RowMergeMacro(8); break;
935         default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatMatMult logic error: Not merging 1-8 rows from L1 work array!");
936         }
937         L2_nnz               += outputi_nnz;
938         worki_L2[++L2_nrows]  = L2_nnz;
939 
940         /************************ L E V E L  3 **********************/
941         /* Merge from L2 work array to either C or to L2 work array */
942         if (anzi > 64 && (L2_nrows == 8 || rowsleft == 0)) {
943           inputi      = worki_L2;
944           inputj      = workj_L2;
945           inputcol    = workcol;
946           outputi_nnz = 0;
947           if (rowsleft) outputj = workj_L3;
948           else          outputj = cj + ci_nnz;
949           switch (L2_nrows) {
950           case 1:  brow_ptr[0] = inputj + inputi[inputcol[0]];
951                    brow_end[0] = inputj + inputi[inputcol[0]+1];
952                    for (; brow_ptr[0] != brow_end[0]; ++brow_ptr[0]) outputj[outputi_nnz++] = *brow_ptr[0]; /* copy row in b over */
953                    break;
954           case 2:  MatMatMultSymbolic_RowMergeMacro(2); break;
955           case 3:  MatMatMultSymbolic_RowMergeMacro(3); break;
956           case 4:  MatMatMultSymbolic_RowMergeMacro(4); break;
957           case 5:  MatMatMultSymbolic_RowMergeMacro(5); break;
958           case 6:  MatMatMultSymbolic_RowMergeMacro(6); break;
959           case 7:  MatMatMultSymbolic_RowMergeMacro(7); break;
960           case 8:  MatMatMultSymbolic_RowMergeMacro(8); break;
961           default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatMatMult logic error: Not merging 1-8 rows from L2 work array!");
962           }
963           L2_nrows    = 1;
964           L2_nnz      = outputi_nnz;
965           worki_L2[1] = outputi_nnz;
966           /* Copy to workj_L2 */
967           if (rowsleft) {
968             for (k=0; k<outputi_nnz; ++k)  workj_L2[k] = outputj[k];
969           }
970         }
971       }
972     }  /* while (rowsleft) */
973 #undef MatMatMultSymbolic_RowMergeMacro
974 
975     /* terminate current row */
976     ci_nnz += outputi_nnz;
977     ci[i+1] = ci_nnz;
978   }
979 
980   /* Step 3: Create the new symbolic matrix */
981   ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,C);CHKERRQ(ierr);
982   ierr = MatSetBlockSizesFromMats(*C,A,B);CHKERRQ(ierr);
983   ierr = MatSetType(*C,((PetscObject)A)->type_name);CHKERRQ(ierr);
984 
985   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
986   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
987   c          = (Mat_SeqAIJ*)((*C)->data);
988   c->free_a  = PETSC_TRUE;
989   c->free_ij = PETSC_TRUE;
990   c->nonew   = 0;
991 
992   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Sorted;
993 
994   /* set MatInfo */
995   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
996   if (afill < 1.0) afill = 1.0;
997   c->maxnz                     = ci[am];
998   c->nz                        = ci[am];
999   (*C)->info.mallocs           = ndouble;
1000   (*C)->info.fill_ratio_given  = fill;
1001   (*C)->info.fill_ratio_needed = afill;
1002 
1003 #if defined(PETSC_USE_INFO)
1004   if (ci[am]) {
1005     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);CHKERRQ(ierr);
1006     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
1007   } else {
1008     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
1009   }
1010 #endif
1011 
1012   /* Step 4: Free temporary work areas */
1013   ierr = PetscFree(workj_L1);CHKERRQ(ierr);
1014   ierr = PetscFree(workj_L2);CHKERRQ(ierr);
1015   ierr = PetscFree(workj_L3);CHKERRQ(ierr);
1016   PetscFunctionReturn(0);
1017 }
1018 
1019 /* concatenate unique entries and then sort */
1020 PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Sorted(Mat A,Mat B,PetscReal fill,Mat *C)
1021 {
1022   PetscErrorCode     ierr;
1023   Mat_SeqAIJ         *a  = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
1024   const PetscInt     *ai = a->i,*bi=b->i,*aj=a->j,*bj=b->j;
1025   PetscInt           *ci,*cj;
1026   PetscInt           am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
1027   PetscReal          afill;
1028   PetscInt           i,j,ndouble = 0;
1029   PetscSegBuffer     seg,segrow;
1030   char               *seen;
1031 
1032   PetscFunctionBegin;
1033   ierr  = PetscMalloc1(am+1,&ci);CHKERRQ(ierr);
1034   ci[0] = 0;
1035 
1036   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
1037   ierr = PetscSegBufferCreate(sizeof(PetscInt),(PetscInt)(fill*(ai[am]+bi[bm])),&seg);CHKERRQ(ierr);
1038   ierr = PetscSegBufferCreate(sizeof(PetscInt),100,&segrow);CHKERRQ(ierr);
1039   ierr = PetscMalloc1(bn,&seen);CHKERRQ(ierr);
1040   ierr = PetscMemzero(seen,bn*sizeof(char));CHKERRQ(ierr);
1041 
1042   /* Determine ci and cj */
1043   for (i=0; i<am; i++) {
1044     const PetscInt anzi  = ai[i+1] - ai[i]; /* number of nonzeros in this row of A, this is the number of rows of B that we merge */
1045     const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */
1046     PetscInt packlen = 0,*PETSC_RESTRICT crow;
1047     /* Pack segrow */
1048     for (j=0; j<anzi; j++) {
1049       PetscInt brow = acol[j],bjstart = bi[brow],bjend = bi[brow+1],k;
1050       for (k=bjstart; k<bjend; k++) {
1051         PetscInt bcol = bj[k];
1052         if (!seen[bcol]) { /* new entry */
1053           PetscInt *PETSC_RESTRICT slot;
1054           ierr = PetscSegBufferGetInts(segrow,1,&slot);CHKERRQ(ierr);
1055           *slot = bcol;
1056           seen[bcol] = 1;
1057           packlen++;
1058         }
1059       }
1060     }
1061     ierr = PetscSegBufferGetInts(seg,packlen,&crow);CHKERRQ(ierr);
1062     ierr = PetscSegBufferExtractTo(segrow,crow);CHKERRQ(ierr);
1063     ierr = PetscSortInt(packlen,crow);CHKERRQ(ierr);
1064     ci[i+1] = ci[i] + packlen;
1065     for (j=0; j<packlen; j++) seen[crow[j]] = 0;
1066   }
1067   ierr = PetscSegBufferDestroy(&segrow);CHKERRQ(ierr);
1068   ierr = PetscFree(seen);CHKERRQ(ierr);
1069 
1070   /* Column indices are in the segmented buffer */
1071   ierr = PetscSegBufferExtractAlloc(seg,&cj);CHKERRQ(ierr);
1072   ierr = PetscSegBufferDestroy(&seg);CHKERRQ(ierr);
1073 
1074   /* put together the new symbolic matrix */
1075   ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,C);CHKERRQ(ierr);
1076   ierr = MatSetBlockSizesFromMats(*C,A,B);CHKERRQ(ierr);
1077   ierr = MatSetType(*C,((PetscObject)A)->type_name);CHKERRQ(ierr);
1078 
1079   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
1080   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
1081   c          = (Mat_SeqAIJ*)((*C)->data);
1082   c->free_a  = PETSC_TRUE;
1083   c->free_ij = PETSC_TRUE;
1084   c->nonew   = 0;
1085 
1086   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Sorted;
1087 
1088   /* set MatInfo */
1089   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
1090   if (afill < 1.0) afill = 1.0;
1091   c->maxnz                     = ci[am];
1092   c->nz                        = ci[am];
1093   (*C)->info.mallocs           = ndouble;
1094   (*C)->info.fill_ratio_given  = fill;
1095   (*C)->info.fill_ratio_needed = afill;
1096 
1097 #if defined(PETSC_USE_INFO)
1098   if (ci[am]) {
1099     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);CHKERRQ(ierr);
1100     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
1101   } else {
1102     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
1103   }
1104 #endif
1105   PetscFunctionReturn(0);
1106 }
1107 
1108 /* This routine is not used. Should be removed! */
1109 PetscErrorCode MatMatTransposeMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
1110 {
1111   PetscErrorCode ierr;
1112 
1113   PetscFunctionBegin;
1114   if (scall == MAT_INITIAL_MATRIX) {
1115     ierr = PetscLogEventBegin(MAT_MatTransposeMultSymbolic,A,B,0,0);CHKERRQ(ierr);
1116     ierr = MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
1117     ierr = PetscLogEventEnd(MAT_MatTransposeMultSymbolic,A,B,0,0);CHKERRQ(ierr);
1118   }
1119   ierr = PetscLogEventBegin(MAT_MatTransposeMultNumeric,A,B,0,0);CHKERRQ(ierr);
1120   ierr = MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr);
1121   ierr = PetscLogEventEnd(MAT_MatTransposeMultNumeric,A,B,0,0);CHKERRQ(ierr);
1122   PetscFunctionReturn(0);
1123 }
1124 
1125 PetscErrorCode MatDestroy_SeqAIJ_MatMatMultTrans(Mat A)
1126 {
1127   PetscErrorCode      ierr;
1128   Mat_SeqAIJ          *a=(Mat_SeqAIJ*)A->data;
1129   Mat_MatMatTransMult *abt=a->abt;
1130 
1131   PetscFunctionBegin;
1132   ierr = (abt->destroy)(A);CHKERRQ(ierr);
1133   ierr = MatTransposeColoringDestroy(&abt->matcoloring);CHKERRQ(ierr);
1134   ierr = MatDestroy(&abt->Bt_den);CHKERRQ(ierr);
1135   ierr = MatDestroy(&abt->ABt_den);CHKERRQ(ierr);
1136   ierr = PetscFree(abt);CHKERRQ(ierr);
1137   PetscFunctionReturn(0);
1138 }
1139 
1140 PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
1141 {
1142   PetscErrorCode      ierr;
1143   Mat                 Bt;
1144   PetscInt            *bti,*btj;
1145   Mat_MatMatTransMult *abt;
1146   Mat_SeqAIJ          *c;
1147 
1148   PetscFunctionBegin;
1149   /* create symbolic Bt */
1150   ierr = MatGetSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr);
1151   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,B->cmap->n,B->rmap->n,bti,btj,NULL,&Bt);CHKERRQ(ierr);
1152   ierr = MatSetBlockSizes(Bt,PetscAbs(A->cmap->bs),PetscAbs(B->cmap->bs));CHKERRQ(ierr);
1153   ierr = MatSetType(Bt,((PetscObject)A)->type_name);CHKERRQ(ierr);
1154 
1155   /* get symbolic C=A*Bt */
1156   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,Bt,fill,C);CHKERRQ(ierr);
1157 
1158   /* create a supporting struct for reuse intermidiate dense matrices with matcoloring */
1159   ierr   = PetscNew(&abt);CHKERRQ(ierr);
1160   c      = (Mat_SeqAIJ*)(*C)->data;
1161   c->abt = abt;
1162 
1163   abt->usecoloring = PETSC_FALSE;
1164   abt->destroy     = (*C)->ops->destroy;
1165   (*C)->ops->destroy     = MatDestroy_SeqAIJ_MatMatMultTrans;
1166 
1167   ierr = PetscOptionsGetBool(((PetscObject)A)->options,NULL,"-matmattransmult_color",&abt->usecoloring,NULL);CHKERRQ(ierr);
1168   if (abt->usecoloring) {
1169     /* Create MatTransposeColoring from symbolic C=A*B^T */
1170     MatTransposeColoring matcoloring;
1171     MatColoring          coloring;
1172     ISColoring           iscoloring;
1173     Mat                  Bt_dense,C_dense;
1174     Mat_SeqAIJ           *c=(Mat_SeqAIJ*)(*C)->data;
1175     /* inode causes memory problem, don't know why */
1176     if (c->inode.use) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MAT_USE_INODES is not supported. Use '-mat_no_inode'");
1177 
1178     ierr = MatColoringCreate(*C,&coloring);CHKERRQ(ierr);
1179     ierr = MatColoringSetDistance(coloring,2);CHKERRQ(ierr);
1180     ierr = MatColoringSetType(coloring,MATCOLORINGSL);CHKERRQ(ierr);
1181     ierr = MatColoringSetFromOptions(coloring);CHKERRQ(ierr);
1182     ierr = MatColoringApply(coloring,&iscoloring);CHKERRQ(ierr);
1183     ierr = MatColoringDestroy(&coloring);CHKERRQ(ierr);
1184     ierr = MatTransposeColoringCreate(*C,iscoloring,&matcoloring);CHKERRQ(ierr);
1185 
1186     abt->matcoloring = matcoloring;
1187 
1188     ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr);
1189 
1190     /* Create Bt_dense and C_dense = A*Bt_dense */
1191     ierr = MatCreate(PETSC_COMM_SELF,&Bt_dense);CHKERRQ(ierr);
1192     ierr = MatSetSizes(Bt_dense,A->cmap->n,matcoloring->ncolors,A->cmap->n,matcoloring->ncolors);CHKERRQ(ierr);
1193     ierr = MatSetType(Bt_dense,MATSEQDENSE);CHKERRQ(ierr);
1194     ierr = MatSeqDenseSetPreallocation(Bt_dense,NULL);CHKERRQ(ierr);
1195 
1196     Bt_dense->assembled = PETSC_TRUE;
1197     abt->Bt_den   = Bt_dense;
1198 
1199     ierr = MatCreate(PETSC_COMM_SELF,&C_dense);CHKERRQ(ierr);
1200     ierr = MatSetSizes(C_dense,A->rmap->n,matcoloring->ncolors,A->rmap->n,matcoloring->ncolors);CHKERRQ(ierr);
1201     ierr = MatSetType(C_dense,MATSEQDENSE);CHKERRQ(ierr);
1202     ierr = MatSeqDenseSetPreallocation(C_dense,NULL);CHKERRQ(ierr);
1203 
1204     Bt_dense->assembled = PETSC_TRUE;
1205     abt->ABt_den  = C_dense;
1206 
1207 #if defined(PETSC_USE_INFO)
1208     {
1209       Mat_SeqAIJ *c = (Mat_SeqAIJ*)(*C)->data;
1210       ierr = PetscInfo7(*C,"Use coloring of C=A*B^T; B^T: %D %D, Bt_dense: %D,%D; Cnz %D / (cm*ncolors %D) = %g\n",B->cmap->n,B->rmap->n,Bt_dense->rmap->n,Bt_dense->cmap->n,c->nz,A->rmap->n*matcoloring->ncolors,(PetscReal)(c->nz)/(A->rmap->n*matcoloring->ncolors));CHKERRQ(ierr);
1211     }
1212 #endif
1213   }
1214   /* clean up */
1215   ierr = MatDestroy(&Bt);CHKERRQ(ierr);
1216   ierr = MatRestoreSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr);
1217   PetscFunctionReturn(0);
1218 }
1219 
1220 PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
1221 {
1222   PetscErrorCode      ierr;
1223   Mat_SeqAIJ          *a   =(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data;
1224   PetscInt            *ai  =a->i,*aj=a->j,*bi=b->i,*bj=b->j,anzi,bnzj,nexta,nextb,*acol,*bcol,brow;
1225   PetscInt            cm   =C->rmap->n,*ci=c->i,*cj=c->j,i,j,cnzi,*ccol;
1226   PetscLogDouble      flops=0.0;
1227   MatScalar           *aa  =a->a,*aval,*ba=b->a,*bval,*ca,*cval;
1228   Mat_MatMatTransMult *abt = c->abt;
1229 
1230   PetscFunctionBegin;
1231   /* clear old values in C */
1232   if (!c->a) {
1233     ierr      = PetscMalloc1(ci[cm]+1,&ca);CHKERRQ(ierr);
1234     c->a      = ca;
1235     c->free_a = PETSC_TRUE;
1236   } else {
1237     ca =  c->a;
1238   }
1239   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
1240 
1241   if (abt->usecoloring) {
1242     MatTransposeColoring matcoloring = abt->matcoloring;
1243     Mat                  Bt_dense,C_dense = abt->ABt_den;
1244 
1245     /* Get Bt_dense by Apply MatTransposeColoring to B */
1246     Bt_dense = abt->Bt_den;
1247     ierr = MatTransColoringApplySpToDen(matcoloring,B,Bt_dense);CHKERRQ(ierr);
1248 
1249     /* C_dense = A*Bt_dense */
1250     ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,Bt_dense,C_dense);CHKERRQ(ierr);
1251 
1252     /* Recover C from C_dense */
1253     ierr = MatTransColoringApplyDenToSp(matcoloring,C_dense,C);CHKERRQ(ierr);
1254     PetscFunctionReturn(0);
1255   }
1256 
1257   for (i=0; i<cm; i++) {
1258     anzi = ai[i+1] - ai[i];
1259     acol = aj + ai[i];
1260     aval = aa + ai[i];
1261     cnzi = ci[i+1] - ci[i];
1262     ccol = cj + ci[i];
1263     cval = ca + ci[i];
1264     for (j=0; j<cnzi; j++) {
1265       brow = ccol[j];
1266       bnzj = bi[brow+1] - bi[brow];
1267       bcol = bj + bi[brow];
1268       bval = ba + bi[brow];
1269 
1270       /* perform sparse inner-product c(i,j)=A[i,:]*B[j,:]^T */
1271       nexta = 0; nextb = 0;
1272       while (nexta<anzi && nextb<bnzj) {
1273         while (nexta < anzi && acol[nexta] < bcol[nextb]) nexta++;
1274         if (nexta == anzi) break;
1275         while (nextb < bnzj && acol[nexta] > bcol[nextb]) nextb++;
1276         if (nextb == bnzj) break;
1277         if (acol[nexta] == bcol[nextb]) {
1278           cval[j] += aval[nexta]*bval[nextb];
1279           nexta++; nextb++;
1280           flops += 2;
1281         }
1282       }
1283     }
1284   }
1285   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1286   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1287   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
1288   PetscFunctionReturn(0);
1289 }
1290 
1291 PetscErrorCode MatDestroy_SeqAIJ_MatTransMatMult(Mat A)
1292 {
1293   PetscErrorCode      ierr;
1294   Mat_SeqAIJ          *a = (Mat_SeqAIJ*)A->data;
1295   Mat_MatTransMatMult *atb = a->atb;
1296 
1297   PetscFunctionBegin;
1298   if (atb) {
1299     ierr = MatDestroy(&atb->At);CHKERRQ(ierr);
1300     ierr = (*atb->destroy)(A);CHKERRQ(ierr);
1301   }
1302   ierr = PetscFree(atb);CHKERRQ(ierr);
1303   PetscFunctionReturn(0);
1304 }
1305 
1306 PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
1307 {
1308   PetscErrorCode      ierr;
1309   const char          *algTypes[2] = {"matmatmult","outerproduct"};
1310   PetscInt            alg=0; /* set default algorithm */
1311   Mat                 At;
1312   Mat_MatTransMatMult *atb;
1313   Mat_SeqAIJ          *c;
1314 
1315   PetscFunctionBegin;
1316   if (scall == MAT_INITIAL_MATRIX) {
1317     ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"MatTransposeMatMult","Mat");CHKERRQ(ierr);
1318     ierr = PetscOptionsEList("-mattransposematmult_via","Algorithmic approach","MatTransposeMatMult",algTypes,2,algTypes[0],&alg,NULL);CHKERRQ(ierr);
1319     ierr = PetscOptionsEnd();CHKERRQ(ierr);
1320 
1321     switch (alg) {
1322     case 1:
1323       ierr = MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
1324       break;
1325     default:
1326       ierr = PetscNew(&atb);CHKERRQ(ierr);
1327       ierr = MatTranspose_SeqAIJ(A,MAT_INITIAL_MATRIX,&At);CHKERRQ(ierr);
1328       ierr = MatMatMult_SeqAIJ_SeqAIJ(At,B,MAT_INITIAL_MATRIX,fill,C);CHKERRQ(ierr);
1329 
1330       c                  = (Mat_SeqAIJ*)(*C)->data;
1331       c->atb             = atb;
1332       atb->At            = At;
1333       atb->destroy       = (*C)->ops->destroy;
1334       (*C)->ops->destroy = MatDestroy_SeqAIJ_MatTransMatMult;
1335 
1336       break;
1337     }
1338   }
1339   if (alg) {
1340     ierr = (*(*C)->ops->mattransposemultnumeric)(A,B,*C);CHKERRQ(ierr);
1341   } else if (!alg && scall == MAT_REUSE_MATRIX) {
1342     c   = (Mat_SeqAIJ*)(*C)->data;
1343     atb = c->atb;
1344     At  = atb->At;
1345     ierr = MatTranspose_SeqAIJ(A,MAT_REUSE_MATRIX,&At);CHKERRQ(ierr);
1346     ierr = MatMatMult_SeqAIJ_SeqAIJ(At,B,MAT_REUSE_MATRIX,fill,C);CHKERRQ(ierr);
1347   }
1348   PetscFunctionReturn(0);
1349 }
1350 
1351 PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
1352 {
1353   PetscErrorCode ierr;
1354   Mat            At;
1355   PetscInt       *ati,*atj;
1356 
1357   PetscFunctionBegin;
1358   /* create symbolic At */
1359   ierr = MatGetSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr);
1360   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,A->cmap->n,A->rmap->n,ati,atj,NULL,&At);CHKERRQ(ierr);
1361   ierr = MatSetBlockSizes(At,PetscAbs(A->cmap->bs),PetscAbs(B->cmap->bs));CHKERRQ(ierr);
1362   ierr = MatSetType(At,((PetscObject)A)->type_name);CHKERRQ(ierr);
1363 
1364   /* get symbolic C=At*B */
1365   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(At,B,fill,C);CHKERRQ(ierr);
1366 
1367   /* clean up */
1368   ierr = MatDestroy(&At);CHKERRQ(ierr);
1369   ierr = MatRestoreSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr);
1370 
1371   (*C)->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ;
1372   PetscFunctionReturn(0);
1373 }
1374 
1375 PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
1376 {
1377   PetscErrorCode ierr;
1378   Mat_SeqAIJ     *a   =(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data;
1379   PetscInt       am   =A->rmap->n,anzi,*ai=a->i,*aj=a->j,*bi=b->i,*bj,bnzi,nextb;
1380   PetscInt       cm   =C->rmap->n,*ci=c->i,*cj=c->j,crow,*cjj,i,j,k;
1381   PetscLogDouble flops=0.0;
1382   MatScalar      *aa  =a->a,*ba,*ca,*caj;
1383 
1384   PetscFunctionBegin;
1385   if (!c->a) {
1386     ierr = PetscMalloc1(ci[cm]+1,&ca);CHKERRQ(ierr);
1387 
1388     c->a      = ca;
1389     c->free_a = PETSC_TRUE;
1390   } else {
1391     ca = c->a;
1392   }
1393   /* clear old values in C */
1394   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
1395 
1396   /* compute A^T*B using outer product (A^T)[:,i]*B[i,:] */
1397   for (i=0; i<am; i++) {
1398     bj   = b->j + bi[i];
1399     ba   = b->a + bi[i];
1400     bnzi = bi[i+1] - bi[i];
1401     anzi = ai[i+1] - ai[i];
1402     for (j=0; j<anzi; j++) {
1403       nextb = 0;
1404       crow  = *aj++;
1405       cjj   = cj + ci[crow];
1406       caj   = ca + ci[crow];
1407       /* perform sparse axpy operation.  Note cjj includes bj. */
1408       for (k=0; nextb<bnzi; k++) {
1409         if (cjj[k] == *(bj+nextb)) { /* ccol == bcol */
1410           caj[k] += (*aa)*(*(ba+nextb));
1411           nextb++;
1412         }
1413       }
1414       flops += 2*bnzi;
1415       aa++;
1416     }
1417   }
1418 
1419   /* Assemble the final matrix and clean up */
1420   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1421   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1422   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
1423   PetscFunctionReturn(0);
1424 }
1425 
1426 PETSC_INTERN PetscErrorCode MatMatMult_SeqAIJ_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
1427 {
1428   PetscErrorCode ierr;
1429 
1430   PetscFunctionBegin;
1431   if (scall == MAT_INITIAL_MATRIX) {
1432     ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
1433     ierr = MatMatMultSymbolic_SeqAIJ_SeqDense(A,B,fill,C);CHKERRQ(ierr);
1434     ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
1435   }
1436   ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
1437   ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,B,*C);CHKERRQ(ierr);
1438   ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
1439   PetscFunctionReturn(0);
1440 }
1441 
1442 PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
1443 {
1444   PetscErrorCode ierr;
1445 
1446   PetscFunctionBegin;
1447   ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,0.0,C);CHKERRQ(ierr);
1448 
1449   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqDense;
1450   PetscFunctionReturn(0);
1451 }
1452 
1453 PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
1454 {
1455   Mat_SeqAIJ        *a=(Mat_SeqAIJ*)A->data;
1456   Mat_SeqDense      *bd = (Mat_SeqDense*)B->data;
1457   PetscErrorCode    ierr;
1458   PetscScalar       *c,*b,r1,r2,r3,r4,*c1,*c2,*c3,*c4,aatmp;
1459   const PetscScalar *aa,*b1,*b2,*b3,*b4;
1460   const PetscInt    *aj;
1461   PetscInt          cm=C->rmap->n,cn=B->cmap->n,bm=bd->lda,am=A->rmap->n;
1462   PetscInt          am4=4*am,bm4=4*bm,col,i,j,n,ajtmp;
1463 
1464   PetscFunctionBegin;
1465   if (!cm || !cn) PetscFunctionReturn(0);
1466   if (B->rmap->n != A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Number columns in A %D not equal rows in B %D\n",A->cmap->n,B->rmap->n);
1467   if (A->rmap->n != C->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Number rows in C %D not equal rows in A %D\n",C->rmap->n,A->rmap->n);
1468   if (B->cmap->n != C->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Number columns in B %D not equal columns in C %D\n",B->cmap->n,C->cmap->n);
1469   b = bd->v;
1470   ierr = MatDenseGetArray(C,&c);CHKERRQ(ierr);
1471   b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;
1472   c1 = c; c2 = c1 + am; c3 = c2 + am; c4 = c3 + am;
1473   for (col=0; col<cn-4; col += 4) {  /* over columns of C */
1474     for (i=0; i<am; i++) {        /* over rows of C in those columns */
1475       r1 = r2 = r3 = r4 = 0.0;
1476       n  = a->i[i+1] - a->i[i];
1477       aj = a->j + a->i[i];
1478       aa = a->a + a->i[i];
1479       for (j=0; j<n; j++) {
1480         aatmp = aa[j]; ajtmp = aj[j];
1481         r1 += aatmp*b1[ajtmp];
1482         r2 += aatmp*b2[ajtmp];
1483         r3 += aatmp*b3[ajtmp];
1484         r4 += aatmp*b4[ajtmp];
1485       }
1486       c1[i] = r1;
1487       c2[i] = r2;
1488       c3[i] = r3;
1489       c4[i] = r4;
1490     }
1491     b1 += bm4; b2 += bm4; b3 += bm4; b4 += bm4;
1492     c1 += am4; c2 += am4; c3 += am4; c4 += am4;
1493   }
1494   for (; col<cn; col++) {   /* over extra columns of C */
1495     for (i=0; i<am; i++) {  /* over rows of C in those columns */
1496       r1 = 0.0;
1497       n  = a->i[i+1] - a->i[i];
1498       aj = a->j + a->i[i];
1499       aa = a->a + a->i[i];
1500       for (j=0; j<n; j++) {
1501         r1 += aa[j]*b1[aj[j]];
1502       }
1503       c1[i] = r1;
1504     }
1505     b1 += bm;
1506     c1 += am;
1507   }
1508   ierr = PetscLogFlops(cn*(2.0*a->nz));CHKERRQ(ierr);
1509   ierr = MatDenseRestoreArray(C,&c);CHKERRQ(ierr);
1510   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1511   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1512   PetscFunctionReturn(0);
1513 }
1514 
1515 /*
1516    Note very similar to MatMult_SeqAIJ(), should generate both codes from same base
1517 */
1518 PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
1519 {
1520   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
1521   Mat_SeqDense   *bd = (Mat_SeqDense*)B->data;
1522   PetscErrorCode ierr;
1523   PetscScalar    *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4;
1524   MatScalar      *aa;
1525   PetscInt       cm  = C->rmap->n, cn=B->cmap->n, bm=bd->lda, col, i,j,n,*aj, am = A->rmap->n,*ii,arm;
1526   PetscInt       am2 = 2*am, am3 = 3*am,  bm4 = 4*bm,colam,*ridx;
1527 
1528   PetscFunctionBegin;
1529   if (!cm || !cn) PetscFunctionReturn(0);
1530   b = bd->v;
1531   ierr = MatDenseGetArray(C,&c);CHKERRQ(ierr);
1532   b1   = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;
1533 
1534   if (a->compressedrow.use) { /* use compressed row format */
1535     for (col=0; col<cn-4; col += 4) {  /* over columns of C */
1536       colam = col*am;
1537       arm   = a->compressedrow.nrows;
1538       ii    = a->compressedrow.i;
1539       ridx  = a->compressedrow.rindex;
1540       for (i=0; i<arm; i++) {        /* over rows of C in those columns */
1541         r1 = r2 = r3 = r4 = 0.0;
1542         n  = ii[i+1] - ii[i];
1543         aj = a->j + ii[i];
1544         aa = a->a + ii[i];
1545         for (j=0; j<n; j++) {
1546           r1 += (*aa)*b1[*aj];
1547           r2 += (*aa)*b2[*aj];
1548           r3 += (*aa)*b3[*aj];
1549           r4 += (*aa++)*b4[*aj++];
1550         }
1551         c[colam       + ridx[i]] += r1;
1552         c[colam + am  + ridx[i]] += r2;
1553         c[colam + am2 + ridx[i]] += r3;
1554         c[colam + am3 + ridx[i]] += r4;
1555       }
1556       b1 += bm4;
1557       b2 += bm4;
1558       b3 += bm4;
1559       b4 += bm4;
1560     }
1561     for (; col<cn; col++) {     /* over extra columns of C */
1562       colam = col*am;
1563       arm   = a->compressedrow.nrows;
1564       ii    = a->compressedrow.i;
1565       ridx  = a->compressedrow.rindex;
1566       for (i=0; i<arm; i++) {  /* over rows of C in those columns */
1567         r1 = 0.0;
1568         n  = ii[i+1] - ii[i];
1569         aj = a->j + ii[i];
1570         aa = a->a + ii[i];
1571 
1572         for (j=0; j<n; j++) {
1573           r1 += (*aa++)*b1[*aj++];
1574         }
1575         c[colam + ridx[i]] += r1;
1576       }
1577       b1 += bm;
1578     }
1579   } else {
1580     for (col=0; col<cn-4; col += 4) {  /* over columns of C */
1581       colam = col*am;
1582       for (i=0; i<am; i++) {        /* over rows of C in those columns */
1583         r1 = r2 = r3 = r4 = 0.0;
1584         n  = a->i[i+1] - a->i[i];
1585         aj = a->j + a->i[i];
1586         aa = a->a + a->i[i];
1587         for (j=0; j<n; j++) {
1588           r1 += (*aa)*b1[*aj];
1589           r2 += (*aa)*b2[*aj];
1590           r3 += (*aa)*b3[*aj];
1591           r4 += (*aa++)*b4[*aj++];
1592         }
1593         c[colam + i]       += r1;
1594         c[colam + am + i]  += r2;
1595         c[colam + am2 + i] += r3;
1596         c[colam + am3 + i] += r4;
1597       }
1598       b1 += bm4;
1599       b2 += bm4;
1600       b3 += bm4;
1601       b4 += bm4;
1602     }
1603     for (; col<cn; col++) {     /* over extra columns of C */
1604       colam = col*am;
1605       for (i=0; i<am; i++) {  /* over rows of C in those columns */
1606         r1 = 0.0;
1607         n  = a->i[i+1] - a->i[i];
1608         aj = a->j + a->i[i];
1609         aa = a->a + a->i[i];
1610 
1611         for (j=0; j<n; j++) {
1612           r1 += (*aa++)*b1[*aj++];
1613         }
1614         c[colam + i] += r1;
1615       }
1616       b1 += bm;
1617     }
1618   }
1619   ierr = PetscLogFlops(cn*2.0*a->nz);CHKERRQ(ierr);
1620   ierr = MatDenseRestoreArray(C,&c);CHKERRQ(ierr);
1621   PetscFunctionReturn(0);
1622 }
1623 
1624 PetscErrorCode  MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring coloring,Mat B,Mat Btdense)
1625 {
1626   PetscErrorCode ierr;
1627   Mat_SeqAIJ     *b       = (Mat_SeqAIJ*)B->data;
1628   Mat_SeqDense   *btdense = (Mat_SeqDense*)Btdense->data;
1629   PetscInt       *bi      = b->i,*bj=b->j;
1630   PetscInt       m        = Btdense->rmap->n,n=Btdense->cmap->n,j,k,l,col,anz,*btcol,brow,ncolumns;
1631   MatScalar      *btval,*btval_den,*ba=b->a;
1632   PetscInt       *columns=coloring->columns,*colorforcol=coloring->colorforcol,ncolors=coloring->ncolors;
1633 
1634   PetscFunctionBegin;
1635   btval_den=btdense->v;
1636   ierr     = PetscMemzero(btval_den,(m*n)*sizeof(MatScalar));CHKERRQ(ierr);
1637   for (k=0; k<ncolors; k++) {
1638     ncolumns = coloring->ncolumns[k];
1639     for (l=0; l<ncolumns; l++) { /* insert a row of B to a column of Btdense */
1640       col   = *(columns + colorforcol[k] + l);
1641       btcol = bj + bi[col];
1642       btval = ba + bi[col];
1643       anz   = bi[col+1] - bi[col];
1644       for (j=0; j<anz; j++) {
1645         brow            = btcol[j];
1646         btval_den[brow] = btval[j];
1647       }
1648     }
1649     btval_den += m;
1650   }
1651   PetscFunctionReturn(0);
1652 }
1653 
1654 PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring matcoloring,Mat Cden,Mat Csp)
1655 {
1656   PetscErrorCode    ierr;
1657   Mat_SeqAIJ        *csp = (Mat_SeqAIJ*)Csp->data;
1658   const PetscScalar *ca_den,*ca_den_ptr;
1659   PetscScalar       *ca=csp->a;
1660   PetscInt          k,l,m=Cden->rmap->n,ncolors=matcoloring->ncolors;
1661   PetscInt          brows=matcoloring->brows,*den2sp=matcoloring->den2sp;
1662   PetscInt          nrows,*row,*idx;
1663   PetscInt          *rows=matcoloring->rows,*colorforrow=matcoloring->colorforrow;
1664 
1665   PetscFunctionBegin;
1666   ierr   = MatDenseGetArrayRead(Cden,&ca_den);CHKERRQ(ierr);
1667 
1668   if (brows > 0) {
1669     PetscInt *lstart,row_end,row_start;
1670     lstart = matcoloring->lstart;
1671     ierr = PetscMemzero(lstart,ncolors*sizeof(PetscInt));CHKERRQ(ierr);
1672 
1673     row_end = brows;
1674     if (row_end > m) row_end = m;
1675     for (row_start=0; row_start<m; row_start+=brows) { /* loop over row blocks of Csp */
1676       ca_den_ptr = ca_den;
1677       for (k=0; k<ncolors; k++) { /* loop over colors (columns of Cden) */
1678         nrows = matcoloring->nrows[k];
1679         row   = rows  + colorforrow[k];
1680         idx   = den2sp + colorforrow[k];
1681         for (l=lstart[k]; l<nrows; l++) {
1682           if (row[l] >= row_end) {
1683             lstart[k] = l;
1684             break;
1685           } else {
1686             ca[idx[l]] = ca_den_ptr[row[l]];
1687           }
1688         }
1689         ca_den_ptr += m;
1690       }
1691       row_end += brows;
1692       if (row_end > m) row_end = m;
1693     }
1694   } else { /* non-blocked impl: loop over columns of Csp - slow if Csp is large */
1695     ca_den_ptr = ca_den;
1696     for (k=0; k<ncolors; k++) {
1697       nrows = matcoloring->nrows[k];
1698       row   = rows  + colorforrow[k];
1699       idx   = den2sp + colorforrow[k];
1700       for (l=0; l<nrows; l++) {
1701         ca[idx[l]] = ca_den_ptr[row[l]];
1702       }
1703       ca_den_ptr += m;
1704     }
1705   }
1706 
1707   ierr = MatDenseRestoreArrayRead(Cden,&ca_den);CHKERRQ(ierr);
1708 #if defined(PETSC_USE_INFO)
1709   if (matcoloring->brows > 0) {
1710     ierr = PetscInfo1(Csp,"Loop over %D row blocks for den2sp\n",brows);CHKERRQ(ierr);
1711   } else {
1712     ierr = PetscInfo(Csp,"Loop over colors/columns of Cden, inefficient for large sparse matrix product \n");CHKERRQ(ierr);
1713   }
1714 #endif
1715   PetscFunctionReturn(0);
1716 }
1717 
1718 PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat mat,ISColoring iscoloring,MatTransposeColoring c)
1719 {
1720   PetscErrorCode ierr;
1721   PetscInt       i,n,nrows,Nbs,j,k,m,ncols,col,cm;
1722   const PetscInt *is,*ci,*cj,*row_idx;
1723   PetscInt       nis = iscoloring->n,*rowhit,bs = 1;
1724   IS             *isa;
1725   Mat_SeqAIJ     *csp = (Mat_SeqAIJ*)mat->data;
1726   PetscInt       *colorforrow,*rows,*rows_i,*idxhit,*spidx,*den2sp,*den2sp_i;
1727   PetscInt       *colorforcol,*columns,*columns_i,brows;
1728   PetscBool      flg;
1729 
1730   PetscFunctionBegin;
1731   ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr);
1732 
1733   /* bs >1 is not being tested yet! */
1734   Nbs       = mat->cmap->N/bs;
1735   c->M      = mat->rmap->N/bs;  /* set total rows, columns and local rows */
1736   c->N      = Nbs;
1737   c->m      = c->M;
1738   c->rstart = 0;
1739   c->brows  = 100;
1740 
1741   c->ncolors = nis;
1742   ierr = PetscMalloc3(nis,&c->ncolumns,nis,&c->nrows,nis+1,&colorforrow);CHKERRQ(ierr);
1743   ierr = PetscMalloc1(csp->nz+1,&rows);CHKERRQ(ierr);
1744   ierr = PetscMalloc1(csp->nz+1,&den2sp);CHKERRQ(ierr);
1745 
1746   brows = c->brows;
1747   ierr = PetscOptionsGetInt(NULL,NULL,"-matden2sp_brows",&brows,&flg);CHKERRQ(ierr);
1748   if (flg) c->brows = brows;
1749   if (brows > 0) {
1750     ierr = PetscMalloc1(nis+1,&c->lstart);CHKERRQ(ierr);
1751   }
1752 
1753   colorforrow[0] = 0;
1754   rows_i         = rows;
1755   den2sp_i       = den2sp;
1756 
1757   ierr = PetscMalloc1(nis+1,&colorforcol);CHKERRQ(ierr);
1758   ierr = PetscMalloc1(Nbs+1,&columns);CHKERRQ(ierr);
1759 
1760   colorforcol[0] = 0;
1761   columns_i      = columns;
1762 
1763   /* get column-wise storage of mat */
1764   ierr = MatGetColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr);
1765 
1766   cm   = c->m;
1767   ierr = PetscMalloc1(cm+1,&rowhit);CHKERRQ(ierr);
1768   ierr = PetscMalloc1(cm+1,&idxhit);CHKERRQ(ierr);
1769   for (i=0; i<nis; i++) { /* loop over color */
1770     ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr);
1771     ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr);
1772 
1773     c->ncolumns[i] = n;
1774     if (n) {
1775       ierr = PetscMemcpy(columns_i,is,n*sizeof(PetscInt));CHKERRQ(ierr);
1776     }
1777     colorforcol[i+1] = colorforcol[i] + n;
1778     columns_i       += n;
1779 
1780     /* fast, crude version requires O(N*N) work */
1781     ierr = PetscMemzero(rowhit,cm*sizeof(PetscInt));CHKERRQ(ierr);
1782 
1783     for (j=0; j<n; j++) { /* loop over columns*/
1784       col     = is[j];
1785       row_idx = cj + ci[col];
1786       m       = ci[col+1] - ci[col];
1787       for (k=0; k<m; k++) { /* loop over columns marking them in rowhit */
1788         idxhit[*row_idx]   = spidx[ci[col] + k];
1789         rowhit[*row_idx++] = col + 1;
1790       }
1791     }
1792     /* count the number of hits */
1793     nrows = 0;
1794     for (j=0; j<cm; j++) {
1795       if (rowhit[j]) nrows++;
1796     }
1797     c->nrows[i]      = nrows;
1798     colorforrow[i+1] = colorforrow[i] + nrows;
1799 
1800     nrows = 0;
1801     for (j=0; j<cm; j++) { /* loop over rows */
1802       if (rowhit[j]) {
1803         rows_i[nrows]   = j;
1804         den2sp_i[nrows] = idxhit[j];
1805         nrows++;
1806       }
1807     }
1808     den2sp_i += nrows;
1809 
1810     ierr    = ISRestoreIndices(isa[i],&is);CHKERRQ(ierr);
1811     rows_i += nrows;
1812   }
1813   ierr = MatRestoreColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr);
1814   ierr = PetscFree(rowhit);CHKERRQ(ierr);
1815   ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr);
1816   if (csp->nz != colorforrow[nis]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"csp->nz %d != colorforrow[nis] %d",csp->nz,colorforrow[nis]);
1817 
1818   c->colorforrow = colorforrow;
1819   c->rows        = rows;
1820   c->den2sp      = den2sp;
1821   c->colorforcol = colorforcol;
1822   c->columns     = columns;
1823 
1824   ierr = PetscFree(idxhit);CHKERRQ(ierr);
1825   PetscFunctionReturn(0);
1826 }
1827 
1828 /* The combine method has done the symbolic and numeric in the first phase, and so we just return here */
1829 PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Combined(Mat A,Mat B,Mat C)
1830 {
1831   PetscFunctionBegin;
1832   /* Empty function */
1833   PetscFunctionReturn(0);
1834 }
1835 
1836 /* This algorithm combines the symbolic and numeric phase of matrix-matrix multiplication. */
1837 PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ_Combined(Mat A,Mat B,PetscReal fill,Mat *C)
1838 {
1839   PetscErrorCode     ierr;
1840   PetscLogDouble     flops=0.0;
1841   Mat_SeqAIJ         *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
1842   const PetscInt     *ai=a->i,*bi=b->i;
1843   PetscInt           *ci,*cj,*cj_i;
1844   PetscScalar        *ca,*ca_i;
1845   PetscInt           b_maxmemrow,c_maxmem,a_col;
1846   PetscInt           am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
1847   PetscInt           i,k,ndouble=0;
1848   PetscReal          afill;
1849   PetscScalar        *c_row_val_dense;
1850   PetscBool          *c_row_idx_flags;
1851   PetscInt           *aj_i=a->j;
1852   PetscScalar        *aa_i=a->a;
1853 
1854   PetscFunctionBegin;
1855 
1856   /* Step 1: Determine upper bounds on memory for C and allocate memory */
1857   /* This should be enough for almost all matrices. If still more memory is needed, it is reallocated later. */
1858   c_maxmem    = 8*(ai[am]+bi[bm]);
1859   b_maxmemrow = PetscMin(bi[bm],bn);
1860   ierr  = PetscMalloc1(am+1,&ci);CHKERRQ(ierr);
1861   ierr  = PetscMalloc1(bn,&c_row_val_dense);CHKERRQ(ierr);
1862   ierr  = PetscMalloc1(bn,&c_row_idx_flags);CHKERRQ(ierr);
1863   ierr  = PetscMalloc1(c_maxmem,&cj);CHKERRQ(ierr);
1864   ierr  = PetscMalloc1(c_maxmem,&ca);CHKERRQ(ierr);
1865   ca_i  = ca;
1866   cj_i  = cj;
1867   ci[0] = 0;
1868   ierr  = PetscMemzero(c_row_val_dense,bn*sizeof(PetscScalar));CHKERRQ(ierr);
1869   ierr  = PetscMemzero(c_row_idx_flags,bn*sizeof(PetscBool));CHKERRQ(ierr);
1870   for (i=0; i<am; i++) {
1871     /* Step 2: Initialize the dense row vector for C  */
1872     const PetscInt anzi = ai[i+1] - ai[i]; /* number of nonzeros in this row of A, this is the number of rows of B that we merge */
1873     PetscInt       cnzi = 0;
1874     PetscInt       *bj_i;
1875     PetscScalar    *ba_i;
1876     /* If the number of indices in C so far + the max number of columns in the next row > c_maxmem  -> allocate more memory
1877        Usually, there is enough memory in the first place, so this is not executed. */
1878     while (ci[i] + b_maxmemrow > c_maxmem) {
1879       c_maxmem *= 2;
1880       ndouble++;
1881       ierr = PetscRealloc(sizeof(PetscInt)*c_maxmem,&cj);CHKERRQ(ierr);
1882       ierr = PetscRealloc(sizeof(PetscScalar)*c_maxmem,&ca);CHKERRQ(ierr);
1883     }
1884 
1885     /* Step 3: Do the numerical calculations */
1886     for (a_col=0; a_col<anzi; a_col++) {          /* iterate over all non zero values in a row of A */
1887       PetscInt       a_col_index = aj_i[a_col];
1888       const PetscInt bnzi        = bi[a_col_index+1] - bi[a_col_index];
1889       flops += 2*bnzi;
1890       bj_i   = b->j + bi[a_col_index];   /* points to the current row in bj */
1891       ba_i   = b->a + bi[a_col_index];   /* points to the current row in ba */
1892       for (k=0; k<bnzi; ++k) { /* iterate over all non zeros of this row in B */
1893         if (c_row_idx_flags[bj_i[k]] == PETSC_FALSE) {
1894           cj_i[cnzi++]             = bj_i[k];
1895           c_row_idx_flags[bj_i[k]] = PETSC_TRUE;
1896         }
1897         c_row_val_dense[bj_i[k]] += aa_i[a_col] * ba_i[k];
1898       }
1899     }
1900 
1901     /* Sort array */
1902     ierr = PetscSortInt(cnzi,cj_i);CHKERRQ(ierr);
1903     /* Step 4 */
1904     for (k=0; k<cnzi; k++) {
1905       ca_i[k]                  = c_row_val_dense[cj_i[k]];
1906       c_row_val_dense[cj_i[k]] = 0.;
1907       c_row_idx_flags[cj_i[k]] = PETSC_FALSE;
1908     }
1909     /* terminate current row */
1910     aa_i   += anzi;
1911     aj_i   += anzi;
1912     ca_i   += cnzi;
1913     cj_i   += cnzi;
1914     ci[i+1] = ci[i] + cnzi;
1915     flops  += cnzi;
1916   }
1917 
1918   /* Step 5 */
1919   /* Create the new matrix */
1920   ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,C);CHKERRQ(ierr);
1921   ierr = MatSetBlockSizesFromMats(*C,A,B);CHKERRQ(ierr);
1922   ierr = MatSetType(*C,((PetscObject)A)->type_name);CHKERRQ(ierr);
1923 
1924   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
1925   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
1926   c          = (Mat_SeqAIJ*)((*C)->data);
1927   c->a       = ca;
1928   c->free_a  = PETSC_TRUE;
1929   c->free_ij = PETSC_TRUE;
1930   c->nonew   = 0;
1931 
1932   /* set MatInfo */
1933   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
1934   if (afill < 1.0) afill = 1.0;
1935   c->maxnz                     = ci[am];
1936   c->nz                        = ci[am];
1937   (*C)->info.mallocs           = ndouble;
1938   (*C)->info.fill_ratio_given  = fill;
1939   (*C)->info.fill_ratio_needed = afill;
1940   (*C)->ops->matmultnumeric    = MatMatMultNumeric_SeqAIJ_SeqAIJ_Combined;
1941 
1942   ierr = PetscFree(c_row_val_dense);CHKERRQ(ierr);
1943   ierr = PetscFree(c_row_idx_flags);CHKERRQ(ierr);
1944 
1945   ierr = MatAssemblyBegin(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1946   ierr = MatAssemblyEnd(*C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1947   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
1948   PetscFunctionReturn(0);
1949 }
1950