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