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