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