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