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