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