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