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