xref: /petsc/src/mat/impls/aij/seq/matmatmult.c (revision 95a2cb335deee435f0b06953e4461b2237b5f64e)
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 PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
14 {
15   PetscErrorCode ierr;
16 
17   PetscFunctionBegin;
18   if (C->ops->matmultnumeric) {
19     if (C->ops->matmultnumeric == MatMatMultNumeric_SeqAIJ_SeqAIJ) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Recursive call");
20     ierr = (*C->ops->matmultnumeric)(A,B,C);CHKERRQ(ierr);
21   } else {
22     ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ_Sorted(A,B,C);CHKERRQ(ierr);
23   }
24   PetscFunctionReturn(0);
25 }
26 
27 /* Modified from MatCreateSeqAIJWithArrays() */
28 PETSC_INTERN PetscErrorCode MatSetSeqAIJWithArrays_private(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt i[],PetscInt j[],PetscScalar a[],MatType mtype,Mat mat)
29 {
30   PetscErrorCode ierr;
31   PetscInt       ii;
32   Mat_SeqAIJ     *aij;
33   PetscBool      isseqaij;
34 
35   PetscFunctionBegin;
36   if (m > 0 && i[0]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"i (row indices) must start with 0");
37   ierr = MatSetSizes(mat,m,n,m,n);CHKERRQ(ierr);
38 
39   if (!mtype) {
40     ierr = PetscObjectBaseTypeCompare((PetscObject)mat,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
41     if (!isseqaij) { ierr = MatSetType(mat,MATSEQAIJ);CHKERRQ(ierr); }
42   } else {
43     ierr = MatSetType(mat,mtype);CHKERRQ(ierr);
44   }
45   ierr = MatSeqAIJSetPreallocation_SeqAIJ(mat,MAT_SKIP_ALLOCATION,NULL);CHKERRQ(ierr);
46   aij  = (Mat_SeqAIJ*)(mat)->data;
47   ierr = PetscMalloc1(m,&aij->imax);CHKERRQ(ierr);
48   ierr = PetscMalloc1(m,&aij->ilen);CHKERRQ(ierr);
49 
50   aij->i            = i;
51   aij->j            = j;
52   aij->a            = a;
53   aij->singlemalloc = PETSC_FALSE;
54   aij->nonew        = -1; /*this indicates that inserting a new value in the matrix that generates a new nonzero is an error*/
55   aij->free_a       = PETSC_FALSE;
56   aij->free_ij      = PETSC_FALSE;
57 
58   for (ii=0; ii<m; ii++) {
59     aij->ilen[ii] = aij->imax[ii] = i[ii+1] - i[ii];
60   }
61 
62   PetscFunctionReturn(0);
63 }
64 
65 PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat C)
66 {
67   PetscErrorCode      ierr;
68   Mat_Product         *product = C->product;
69   MatProductAlgorithm alg;
70   PetscBool           flg;
71 
72   PetscFunctionBegin;
73   if (product) {
74     alg = product->alg;
75   } else {
76     alg = "sorted";
77   }
78   /* sorted */
79   ierr = PetscStrcmp(alg,"sorted",&flg);CHKERRQ(ierr);
80   if (flg) {
81     ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Sorted(A,B,fill,C);CHKERRQ(ierr);
82     PetscFunctionReturn(0);
83   }
84 
85   /* scalable */
86   ierr = PetscStrcmp(alg,"scalable",&flg);CHKERRQ(ierr);
87   if (flg) {
88     ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(A,B,fill,C);CHKERRQ(ierr);
89     PetscFunctionReturn(0);
90   }
91 
92   /* scalable_fast */
93   ierr = PetscStrcmp(alg,"scalable_fast",&flg);CHKERRQ(ierr);
94   if (flg) {
95     ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(A,B,fill,C);CHKERRQ(ierr);
96     PetscFunctionReturn(0);
97   }
98 
99   /* heap */
100   ierr = PetscStrcmp(alg,"heap",&flg);CHKERRQ(ierr);
101   if (flg) {
102     ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(A,B,fill,C);CHKERRQ(ierr);
103     PetscFunctionReturn(0);
104   }
105 
106   /* btheap */
107   ierr = PetscStrcmp(alg,"btheap",&flg);CHKERRQ(ierr);
108   if (flg) {
109     ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(A,B,fill,C);CHKERRQ(ierr);
110     PetscFunctionReturn(0);
111   }
112 
113   /* llcondensed */
114   ierr = PetscStrcmp(alg,"llcondensed",&flg);CHKERRQ(ierr);
115   if (flg) {
116     ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_LLCondensed(A,B,fill,C);CHKERRQ(ierr);
117     PetscFunctionReturn(0);
118   }
119 
120   /* rowmerge */
121   ierr = PetscStrcmp(alg,"rowmerge",&flg);CHKERRQ(ierr);
122   if (flg) {
123     ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_RowMerge(A,B,fill,C);CHKERRQ(ierr);
124     PetscFunctionReturn(0);
125   }
126 
127 #if defined(PETSC_HAVE_HYPRE)
128   ierr = PetscStrcmp(alg,"hypre",&flg);CHKERRQ(ierr);
129   if (flg) {
130     ierr = MatMatMultSymbolic_AIJ_AIJ_wHYPRE(A,B,fill,C);CHKERRQ(ierr);
131     PetscFunctionReturn(0);
132   }
133 #endif
134 
135   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Mat Product Algorithm is not supported");
136 }
137 
138 PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_LLCondensed(Mat A,Mat B,PetscReal fill,Mat C)
139 {
140   PetscErrorCode     ierr;
141   Mat_SeqAIJ         *a =(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
142   PetscInt           *ai=a->i,*bi=b->i,*ci,*cj;
143   PetscInt           am =A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
144   PetscReal          afill;
145   PetscInt           i,j,anzi,brow,bnzj,cnzi,*bj,*aj,*lnk,ndouble=0,Crmax;
146   PetscTable         ta;
147   PetscBT            lnkbt;
148   PetscFreeSpaceList free_space=NULL,current_space=NULL;
149 
150   PetscFunctionBegin;
151   /* Get ci and cj */
152   /*---------------*/
153   /* Allocate ci array, arrays for fill computation and */
154   /* free space for accumulating nonzero column info */
155   ierr  = PetscMalloc1(am+2,&ci);CHKERRQ(ierr);
156   ci[0] = 0;
157 
158   /* create and initialize a linked list */
159   ierr = PetscTableCreate(bn,bn,&ta);CHKERRQ(ierr);
160   MatRowMergeMax_SeqAIJ(b,bm,ta);
161   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr);
162   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
163 
164   ierr = PetscLLCondensedCreate(Crmax,bn,&lnk,&lnkbt);CHKERRQ(ierr);
165 
166   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
167   ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am],bi[bm])),&free_space);CHKERRQ(ierr);
168 
169   current_space = free_space;
170 
171   /* Determine ci and cj */
172   for (i=0; i<am; i++) {
173     anzi = ai[i+1] - ai[i];
174     aj   = a->j + ai[i];
175     for (j=0; j<anzi; j++) {
176       brow = aj[j];
177       bnzj = bi[brow+1] - bi[brow];
178       bj   = b->j + bi[brow];
179       /* add non-zero cols of B into the sorted linked list lnk */
180       ierr = PetscLLCondensedAddSorted(bnzj,bj,lnk,lnkbt);CHKERRQ(ierr);
181     }
182     /* add possible missing diagonal entry */
183     if (C->force_diagonals) {
184       ierr = PetscLLCondensedAddSorted(1,&i,lnk,lnkbt);CHKERRQ(ierr);
185     }
186     cnzi = lnk[0];
187 
188     /* If free space is not available, make more free space */
189     /* Double the amount of total space in the list */
190     if (current_space->local_remaining<cnzi) {
191       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(cnzi,current_space->total_array_size),&current_space);CHKERRQ(ierr);
192       ndouble++;
193     }
194 
195     /* Copy data into free space, then initialize lnk */
196     ierr = PetscLLCondensedClean(bn,cnzi,current_space->array,lnk,lnkbt);CHKERRQ(ierr);
197 
198     current_space->array           += cnzi;
199     current_space->local_used      += cnzi;
200     current_space->local_remaining -= cnzi;
201 
202     ci[i+1] = ci[i] + cnzi;
203   }
204 
205   /* Column indices are in the list of free space */
206   /* Allocate space for cj, initialize cj, and */
207   /* destroy list of free space and other temporary array(s) */
208   ierr = PetscMalloc1(ci[am]+1,&cj);CHKERRQ(ierr);
209   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
210   ierr = PetscLLCondensedDestroy(lnk,lnkbt);CHKERRQ(ierr);
211 
212   /* put together the new symbolic matrix */
213   ierr = MatSetSeqAIJWithArrays_private(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,((PetscObject)A)->type_name,C);CHKERRQ(ierr);
214   ierr = MatSetBlockSizesFromMats(C,A,B);CHKERRQ(ierr);
215 
216   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
217   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
218   c          = (Mat_SeqAIJ*)(C->data);
219   c->free_a  = PETSC_FALSE;
220   c->free_ij = PETSC_TRUE;
221   c->nonew   = 0;
222 
223   /* fast, needs non-scalable O(bn) array 'abdense' */
224   C->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Sorted;
225 
226   /* set MatInfo */
227   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
228   if (afill < 1.0) afill = 1.0;
229   c->maxnz                  = ci[am];
230   c->nz                     = ci[am];
231   C->info.mallocs           = ndouble;
232   C->info.fill_ratio_given  = fill;
233   C->info.fill_ratio_needed = afill;
234 
235 #if defined(PETSC_USE_INFO)
236   if (ci[am]) {
237     ierr = PetscInfo3(C,"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);CHKERRQ(ierr);
238     ierr = PetscInfo1(C,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
239   } else {
240     ierr = PetscInfo(C,"Empty matrix product\n");CHKERRQ(ierr);
241   }
242 #endif
243   PetscFunctionReturn(0);
244 }
245 
246 PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Sorted(Mat A,Mat B,Mat C)
247 {
248   PetscErrorCode ierr;
249   PetscLogDouble flops=0.0;
250   Mat_SeqAIJ     *a   = (Mat_SeqAIJ*)A->data;
251   Mat_SeqAIJ     *b   = (Mat_SeqAIJ*)B->data;
252   Mat_SeqAIJ     *c   = (Mat_SeqAIJ*)C->data;
253   PetscInt       *ai  =a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j;
254   PetscInt       am   =A->rmap->n,cm=C->rmap->n;
255   PetscInt       i,j,k,anzi,bnzi,cnzi,brow;
256   PetscScalar    *aa=a->a,*ba=b->a,*baj,*ca,valtmp;
257   PetscScalar    *ab_dense;
258   PetscContainer cab_dense;
259 
260   PetscFunctionBegin;
261   if (!c->a) { /* first call of MatMatMultNumeric_SeqAIJ_SeqAIJ, allocate ca and matmult_abdense */
262     ierr      = PetscMalloc1(ci[cm]+1,&ca);CHKERRQ(ierr);
263     c->a      = ca;
264     c->free_a = PETSC_TRUE;
265   } else ca = c->a;
266 
267   /* TODO this should be done in the symbolic phase */
268   /* However, this function is so heavily used (sometimes in an hidden way through multnumeric function pointers
269      that is hard to eradicate) */
270   ierr = PetscObjectQuery((PetscObject)C,"__PETSc__ab_dense",(PetscObject*)&cab_dense);CHKERRQ(ierr);
271   if (!cab_dense) {
272     ierr = PetscMalloc1(B->cmap->N,&ab_dense);CHKERRQ(ierr);
273     ierr = PetscContainerCreate(PETSC_COMM_SELF,&cab_dense);CHKERRQ(ierr);
274     ierr = PetscContainerSetPointer(cab_dense,ab_dense);CHKERRQ(ierr);
275     ierr = PetscContainerSetUserDestroy(cab_dense,PetscContainerUserDestroyDefault);CHKERRQ(ierr);
276     ierr = PetscObjectCompose((PetscObject)C,"__PETSc__ab_dense",(PetscObject)cab_dense);CHKERRQ(ierr);
277     ierr = PetscObjectDereference((PetscObject)cab_dense);CHKERRQ(ierr);
278   }
279   ierr = PetscContainerGetPointer(cab_dense,(void**)&ab_dense);CHKERRQ(ierr);
280   ierr = PetscArrayzero(ab_dense,B->cmap->N);CHKERRQ(ierr);
281 
282   /* clean old values in C */
283   ierr = PetscArrayzero(ca,ci[cm]);CHKERRQ(ierr);
284   /* Traverse A row-wise. */
285   /* Build the ith row in C by summing over nonzero columns in A, */
286   /* the rows of B corresponding to nonzeros of A. */
287   for (i=0; i<am; i++) {
288     anzi = ai[i+1] - ai[i];
289     for (j=0; j<anzi; j++) {
290       brow = aj[j];
291       bnzi = bi[brow+1] - bi[brow];
292       bjj  = bj + bi[brow];
293       baj  = ba + bi[brow];
294       /* perform dense axpy */
295       valtmp = aa[j];
296       for (k=0; k<bnzi; k++) {
297         ab_dense[bjj[k]] += valtmp*baj[k];
298       }
299       flops += 2*bnzi;
300     }
301     aj += anzi; aa += anzi;
302 
303     cnzi = ci[i+1] - ci[i];
304     for (k=0; k<cnzi; k++) {
305       ca[k]          += ab_dense[cj[k]];
306       ab_dense[cj[k]] = 0.0; /* zero ab_dense */
307     }
308     flops += cnzi;
309     cj    += cnzi; ca += cnzi;
310   }
311   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
312   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
313   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
314   PetscFunctionReturn(0);
315 }
316 
317 PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(Mat A,Mat B,Mat C)
318 {
319   PetscErrorCode ierr;
320   PetscLogDouble flops=0.0;
321   Mat_SeqAIJ     *a   = (Mat_SeqAIJ*)A->data;
322   Mat_SeqAIJ     *b   = (Mat_SeqAIJ*)B->data;
323   Mat_SeqAIJ     *c   = (Mat_SeqAIJ*)C->data;
324   PetscInt       *ai  = a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j;
325   PetscInt       am   = A->rmap->N,cm=C->rmap->N;
326   PetscInt       i,j,k,anzi,bnzi,cnzi,brow;
327   PetscScalar    *aa=a->a,*ba=b->a,*baj,*ca=c->a,valtmp;
328   PetscInt       nextb;
329 
330   PetscFunctionBegin;
331   if (!ca) { /* first call of MatMatMultNumeric_SeqAIJ_SeqAIJ, allocate ca and matmult_abdense */
332     ierr      = PetscMalloc1(ci[cm]+1,&ca);CHKERRQ(ierr);
333     c->a      = ca;
334     c->free_a = PETSC_TRUE;
335   }
336 
337   /* clean old values in C */
338   ierr = PetscArrayzero(ca,ci[cm]);CHKERRQ(ierr);
339   /* Traverse A row-wise. */
340   /* Build the ith row in C by summing over nonzero columns in A, */
341   /* the rows of B corresponding to nonzeros of A. */
342   for (i=0; i<am; i++) {
343     anzi = ai[i+1] - ai[i];
344     cnzi = ci[i+1] - ci[i];
345     for (j=0; j<anzi; j++) {
346       brow = aj[j];
347       bnzi = bi[brow+1] - bi[brow];
348       bjj  = bj + bi[brow];
349       baj  = ba + bi[brow];
350       /* perform sparse axpy */
351       valtmp = aa[j];
352       nextb  = 0;
353       for (k=0; nextb<bnzi; k++) {
354         if (cj[k] == bjj[nextb]) { /* ccol == bcol */
355           ca[k] += valtmp*baj[nextb++];
356         }
357       }
358       flops += 2*bnzi;
359     }
360     aj += anzi; aa += anzi;
361     cj += cnzi; ca += cnzi;
362   }
363   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
364   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
365   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
366   PetscFunctionReturn(0);
367 }
368 
369 PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(Mat A,Mat B,PetscReal fill,Mat C)
370 {
371   PetscErrorCode     ierr;
372   Mat_SeqAIJ         *a  = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
373   PetscInt           *ai = a->i,*bi=b->i,*ci,*cj;
374   PetscInt           am  = A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
375   MatScalar          *ca;
376   PetscReal          afill;
377   PetscInt           i,j,anzi,brow,bnzj,cnzi,*bj,*aj,*lnk,ndouble=0,Crmax;
378   PetscTable         ta;
379   PetscFreeSpaceList free_space=NULL,current_space=NULL;
380 
381   PetscFunctionBegin;
382   /* Get ci and cj - same as MatMatMultSymbolic_SeqAIJ_SeqAIJ except using PetscLLxxx_fast() */
383   /*-----------------------------------------------------------------------------------------*/
384   /* Allocate arrays for fill computation and free space for accumulating nonzero column */
385   ierr  = PetscMalloc1(am+2,&ci);CHKERRQ(ierr);
386   ci[0] = 0;
387 
388   /* create and initialize a linked list */
389   ierr = PetscTableCreate(bn,bn,&ta);CHKERRQ(ierr);
390   MatRowMergeMax_SeqAIJ(b,bm,ta);
391   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr);
392   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
393 
394   ierr = PetscLLCondensedCreate_fast(Crmax,&lnk);CHKERRQ(ierr);
395 
396   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
397   ierr          = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am],bi[bm])),&free_space);CHKERRQ(ierr);
398   current_space = free_space;
399 
400   /* Determine ci and cj */
401   for (i=0; i<am; i++) {
402     anzi = ai[i+1] - ai[i];
403     aj   = a->j + ai[i];
404     for (j=0; j<anzi; j++) {
405       brow = aj[j];
406       bnzj = bi[brow+1] - bi[brow];
407       bj   = b->j + bi[brow];
408       /* add non-zero cols of B into the sorted linked list lnk */
409       ierr = PetscLLCondensedAddSorted_fast(bnzj,bj,lnk);CHKERRQ(ierr);
410     }
411     /* add possible missing diagonal entry */
412     if (C->force_diagonals) {
413       ierr = PetscLLCondensedAddSorted_fast(1,&i,lnk);CHKERRQ(ierr);
414     }
415     cnzi = lnk[1];
416 
417     /* If free space is not available, make more free space */
418     /* Double the amount of total space in the list */
419     if (current_space->local_remaining<cnzi) {
420       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(cnzi,current_space->total_array_size),&current_space);CHKERRQ(ierr);
421       ndouble++;
422     }
423 
424     /* Copy data into free space, then initialize lnk */
425     ierr = PetscLLCondensedClean_fast(cnzi,current_space->array,lnk);CHKERRQ(ierr);
426 
427     current_space->array           += cnzi;
428     current_space->local_used      += cnzi;
429     current_space->local_remaining -= cnzi;
430 
431     ci[i+1] = ci[i] + cnzi;
432   }
433 
434   /* Column indices are in the list of free space */
435   /* Allocate space for cj, initialize cj, and */
436   /* destroy list of free space and other temporary array(s) */
437   ierr = PetscMalloc1(ci[am]+1,&cj);CHKERRQ(ierr);
438   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
439   ierr = PetscLLCondensedDestroy_fast(lnk);CHKERRQ(ierr);
440 
441   /* Allocate space for ca */
442   ierr = PetscCalloc1(ci[am]+1,&ca);CHKERRQ(ierr);
443 
444   /* put together the new symbolic matrix */
445   ierr = MatSetSeqAIJWithArrays_private(PetscObjectComm((PetscObject)A),am,bn,ci,cj,ca,((PetscObject)A)->type_name,C);CHKERRQ(ierr);
446   ierr = MatSetBlockSizesFromMats(C,A,B);CHKERRQ(ierr);
447 
448   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
449   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
450   c          = (Mat_SeqAIJ*)(C->data);
451   c->free_a  = PETSC_TRUE;
452   c->free_ij = PETSC_TRUE;
453   c->nonew   = 0;
454 
455   /* slower, less memory */
456   C->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable;
457 
458   /* set MatInfo */
459   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
460   if (afill < 1.0) afill = 1.0;
461   c->maxnz                     = ci[am];
462   c->nz                        = ci[am];
463   C->info.mallocs           = ndouble;
464   C->info.fill_ratio_given  = fill;
465   C->info.fill_ratio_needed = afill;
466 
467 #if defined(PETSC_USE_INFO)
468   if (ci[am]) {
469     ierr = PetscInfo3(C,"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);CHKERRQ(ierr);
470     ierr = PetscInfo1(C,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
471   } else {
472     ierr = PetscInfo(C,"Empty matrix product\n");CHKERRQ(ierr);
473   }
474 #endif
475   PetscFunctionReturn(0);
476 }
477 
478 PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(Mat A,Mat B,PetscReal fill,Mat C)
479 {
480   PetscErrorCode     ierr;
481   Mat_SeqAIJ         *a  = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
482   PetscInt           *ai = a->i,*bi=b->i,*ci,*cj;
483   PetscInt           am  = A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
484   MatScalar          *ca;
485   PetscReal          afill;
486   PetscInt           i,j,anzi,brow,bnzj,cnzi,*bj,*aj,*lnk,ndouble=0,Crmax;
487   PetscTable         ta;
488   PetscFreeSpaceList free_space=NULL,current_space=NULL;
489 
490   PetscFunctionBegin;
491   /* Get ci and cj - same as MatMatMultSymbolic_SeqAIJ_SeqAIJ except using PetscLLxxx_Scalalbe() */
492   /*---------------------------------------------------------------------------------------------*/
493   /* Allocate arrays for fill computation and free space for accumulating nonzero column */
494   ierr  = PetscMalloc1(am+2,&ci);CHKERRQ(ierr);
495   ci[0] = 0;
496 
497   /* create and initialize a linked list */
498   ierr = PetscTableCreate(bn,bn,&ta);CHKERRQ(ierr);
499   MatRowMergeMax_SeqAIJ(b,bm,ta);
500   ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr);
501   ierr = PetscTableDestroy(&ta);CHKERRQ(ierr);
502   ierr = PetscLLCondensedCreate_Scalable(Crmax,&lnk);CHKERRQ(ierr);
503 
504   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
505   ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am],bi[bm])),&free_space);CHKERRQ(ierr);
506   current_space = free_space;
507 
508   /* Determine ci and cj */
509   for (i=0; i<am; i++) {
510     anzi = ai[i+1] - ai[i];
511     aj   = a->j + ai[i];
512     for (j=0; j<anzi; j++) {
513       brow = aj[j];
514       bnzj = bi[brow+1] - bi[brow];
515       bj   = b->j + bi[brow];
516       /* add non-zero cols of B into the sorted linked list lnk */
517       ierr = PetscLLCondensedAddSorted_Scalable(bnzj,bj,lnk);CHKERRQ(ierr);
518     }
519     /* add possible missing diagonal entry */
520     if (C->force_diagonals) {
521       ierr = PetscLLCondensedAddSorted_Scalable(1,&i,lnk);CHKERRQ(ierr);
522     }
523 
524     cnzi = lnk[0];
525 
526     /* If free space is not available, make more free space */
527     /* Double the amount of total space in the list */
528     if (current_space->local_remaining<cnzi) {
529       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(cnzi,current_space->total_array_size),&current_space);CHKERRQ(ierr);
530       ndouble++;
531     }
532 
533     /* Copy data into free space, then initialize lnk */
534     ierr = PetscLLCondensedClean_Scalable(cnzi,current_space->array,lnk);CHKERRQ(ierr);
535 
536     current_space->array           += cnzi;
537     current_space->local_used      += cnzi;
538     current_space->local_remaining -= cnzi;
539 
540     ci[i+1] = ci[i] + cnzi;
541   }
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 = PetscMalloc1(ci[am]+1,&cj);CHKERRQ(ierr);
547   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
548   ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr);
549 
550   /* Allocate space for ca */
551   /*-----------------------*/
552   ierr = PetscCalloc1(ci[am]+1,&ca);CHKERRQ(ierr);
553 
554   /* put together the new symbolic matrix */
555   ierr = MatSetSeqAIJWithArrays_private(PetscObjectComm((PetscObject)A),am,bn,ci,cj,ca,((PetscObject)A)->type_name,C);CHKERRQ(ierr);
556   ierr = MatSetBlockSizesFromMats(C,A,B);CHKERRQ(ierr);
557 
558   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
559   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
560   c          = (Mat_SeqAIJ*)(C->data);
561   c->free_a  = PETSC_TRUE;
562   c->free_ij = PETSC_TRUE;
563   c->nonew   = 0;
564 
565   /* slower, less memory */
566   C->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable;
567 
568   /* set MatInfo */
569   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
570   if (afill < 1.0) afill = 1.0;
571   c->maxnz                     = ci[am];
572   c->nz                        = ci[am];
573   C->info.mallocs           = ndouble;
574   C->info.fill_ratio_given  = fill;
575   C->info.fill_ratio_needed = afill;
576 
577 #if defined(PETSC_USE_INFO)
578   if (ci[am]) {
579     ierr = PetscInfo3(C,"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);CHKERRQ(ierr);
580     ierr = PetscInfo1(C,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
581   } else {
582     ierr = PetscInfo(C,"Empty matrix product\n");CHKERRQ(ierr);
583   }
584 #endif
585   PetscFunctionReturn(0);
586 }
587 
588 PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(Mat A,Mat B,PetscReal fill,Mat C)
589 {
590   PetscErrorCode     ierr;
591   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
592   const PetscInt     *ai=a->i,*bi=b->i,*aj=a->j,*bj=b->j;
593   PetscInt           *ci,*cj,*bb;
594   PetscInt           am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
595   PetscReal          afill;
596   PetscInt           i,j,col,ndouble = 0;
597   PetscFreeSpaceList free_space=NULL,current_space=NULL;
598   PetscHeap          h;
599 
600   PetscFunctionBegin;
601   /* Get ci and cj - by merging sorted rows using a heap */
602   /*---------------------------------------------------------------------------------------------*/
603   /* Allocate arrays for fill computation and free space for accumulating nonzero column */
604   ierr  = PetscMalloc1(am+2,&ci);CHKERRQ(ierr);
605   ci[0] = 0;
606 
607   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
608   ierr          = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am],bi[bm])),&free_space);CHKERRQ(ierr);
609   current_space = free_space;
610 
611   ierr = PetscHeapCreate(a->rmax,&h);CHKERRQ(ierr);
612   ierr = PetscMalloc1(a->rmax,&bb);CHKERRQ(ierr);
613 
614   /* Determine ci and cj */
615   for (i=0; i<am; i++) {
616     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 */
617     const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */
618     ci[i+1] = ci[i];
619     /* Populate the min heap */
620     for (j=0; j<anzi; j++) {
621       bb[j] = bi[acol[j]];         /* bb points at the start of the row */
622       if (bb[j] < bi[acol[j]+1]) { /* Add if row is nonempty */
623         ierr = PetscHeapAdd(h,j,bj[bb[j]++]);CHKERRQ(ierr);
624       }
625     }
626     /* Pick off the min element, adding it to free space */
627     ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr);
628     while (j >= 0) {
629       if (current_space->local_remaining < 1) { /* double the size, but don't exceed 16 MiB */
630         ierr = PetscFreeSpaceGet(PetscMin(PetscIntMultTruncate(2,current_space->total_array_size),16 << 20),&current_space);CHKERRQ(ierr);
631         ndouble++;
632       }
633       *(current_space->array++) = col;
634       current_space->local_used++;
635       current_space->local_remaining--;
636       ci[i+1]++;
637 
638       /* stash if anything else remains in this row of B */
639       if (bb[j] < bi[acol[j]+1]) {ierr = PetscHeapStash(h,j,bj[bb[j]++]);CHKERRQ(ierr);}
640       while (1) {               /* pop and stash any other rows of B that also had an entry in this column */
641         PetscInt j2,col2;
642         ierr = PetscHeapPeek(h,&j2,&col2);CHKERRQ(ierr);
643         if (col2 != col) break;
644         ierr = PetscHeapPop(h,&j2,&col2);CHKERRQ(ierr);
645         if (bb[j2] < bi[acol[j2]+1]) {ierr = PetscHeapStash(h,j2,bj[bb[j2]++]);CHKERRQ(ierr);}
646       }
647       /* Put any stashed elements back into the min heap */
648       ierr = PetscHeapUnstash(h);CHKERRQ(ierr);
649       ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr);
650     }
651   }
652   ierr = PetscFree(bb);CHKERRQ(ierr);
653   ierr = PetscHeapDestroy(&h);CHKERRQ(ierr);
654 
655   /* Column indices are in the list of free space */
656   /* Allocate space for cj, initialize cj, and */
657   /* destroy list of free space and other temporary array(s) */
658   ierr = PetscMalloc1(ci[am],&cj);CHKERRQ(ierr);
659   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
660 
661   /* put together the new symbolic matrix */
662   ierr = MatSetSeqAIJWithArrays_private(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,((PetscObject)A)->type_name,C);CHKERRQ(ierr);
663   ierr = MatSetBlockSizesFromMats(C,A,B);CHKERRQ(ierr);
664 
665   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
666   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
667   c          = (Mat_SeqAIJ*)(C->data);
668   c->free_a  = PETSC_TRUE;
669   c->free_ij = PETSC_TRUE;
670   c->nonew   = 0;
671 
672   C->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Sorted;
673 
674   /* set MatInfo */
675   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
676   if (afill < 1.0) afill = 1.0;
677   c->maxnz                     = ci[am];
678   c->nz                        = ci[am];
679   C->info.mallocs           = ndouble;
680   C->info.fill_ratio_given  = fill;
681   C->info.fill_ratio_needed = afill;
682 
683 #if defined(PETSC_USE_INFO)
684   if (ci[am]) {
685     ierr = PetscInfo3(C,"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);CHKERRQ(ierr);
686     ierr = PetscInfo1(C,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
687   } else {
688     ierr = PetscInfo(C,"Empty matrix product\n");CHKERRQ(ierr);
689   }
690 #endif
691   PetscFunctionReturn(0);
692 }
693 
694 PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(Mat A,Mat B,PetscReal fill,Mat C)
695 {
696   PetscErrorCode     ierr;
697   Mat_SeqAIJ         *a  = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
698   const PetscInt     *ai = a->i,*bi=b->i,*aj=a->j,*bj=b->j;
699   PetscInt           *ci,*cj,*bb;
700   PetscInt           am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
701   PetscReal          afill;
702   PetscInt           i,j,col,ndouble = 0;
703   PetscFreeSpaceList free_space=NULL,current_space=NULL;
704   PetscHeap          h;
705   PetscBT            bt;
706 
707   PetscFunctionBegin;
708   /* Get ci and cj - using a heap for the sorted rows, but use BT so that each index is only added once */
709   /*---------------------------------------------------------------------------------------------*/
710   /* Allocate arrays for fill computation and free space for accumulating nonzero column */
711   ierr  = PetscMalloc1(am+2,&ci);CHKERRQ(ierr);
712   ci[0] = 0;
713 
714   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
715   ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am],bi[bm])),&free_space);CHKERRQ(ierr);
716 
717   current_space = free_space;
718 
719   ierr = PetscHeapCreate(a->rmax,&h);CHKERRQ(ierr);
720   ierr = PetscMalloc1(a->rmax,&bb);CHKERRQ(ierr);
721   ierr = PetscBTCreate(bn,&bt);CHKERRQ(ierr);
722 
723   /* Determine ci and cj */
724   for (i=0; i<am; i++) {
725     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 */
726     const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */
727     const PetscInt *fptr = current_space->array; /* Save beginning of the row so we can clear the BT later */
728     ci[i+1] = ci[i];
729     /* Populate the min heap */
730     for (j=0; j<anzi; j++) {
731       PetscInt brow = acol[j];
732       for (bb[j] = bi[brow]; bb[j] < bi[brow+1]; bb[j]++) {
733         PetscInt bcol = bj[bb[j]];
734         if (!PetscBTLookupSet(bt,bcol)) { /* new entry */
735           ierr = PetscHeapAdd(h,j,bcol);CHKERRQ(ierr);
736           bb[j]++;
737           break;
738         }
739       }
740     }
741     /* Pick off the min element, adding it to free space */
742     ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr);
743     while (j >= 0) {
744       if (current_space->local_remaining < 1) { /* double the size, but don't exceed 16 MiB */
745         fptr = NULL;                      /* need PetscBTMemzero */
746         ierr = PetscFreeSpaceGet(PetscMin(PetscIntMultTruncate(2,current_space->total_array_size),16 << 20),&current_space);CHKERRQ(ierr);
747         ndouble++;
748       }
749       *(current_space->array++) = col;
750       current_space->local_used++;
751       current_space->local_remaining--;
752       ci[i+1]++;
753 
754       /* stash if anything else remains in this row of B */
755       for (; bb[j] < bi[acol[j]+1]; bb[j]++) {
756         PetscInt bcol = bj[bb[j]];
757         if (!PetscBTLookupSet(bt,bcol)) { /* new entry */
758           ierr = PetscHeapAdd(h,j,bcol);CHKERRQ(ierr);
759           bb[j]++;
760           break;
761         }
762       }
763       ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr);
764     }
765     if (fptr) {                 /* Clear the bits for this row */
766       for (; fptr<current_space->array; fptr++) {ierr = PetscBTClear(bt,*fptr);CHKERRQ(ierr);}
767     } else {                    /* We reallocated so we don't remember (easily) how to clear only the bits we changed */
768       ierr = PetscBTMemzero(bn,bt);CHKERRQ(ierr);
769     }
770   }
771   ierr = PetscFree(bb);CHKERRQ(ierr);
772   ierr = PetscHeapDestroy(&h);CHKERRQ(ierr);
773   ierr = PetscBTDestroy(&bt);CHKERRQ(ierr);
774 
775   /* Column indices are in the list of free space */
776   /* Allocate space for cj, initialize cj, and */
777   /* destroy list of free space and other temporary array(s) */
778   ierr = PetscMalloc1(ci[am],&cj);CHKERRQ(ierr);
779   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
780 
781   /* put together the new symbolic matrix */
782   ierr = MatSetSeqAIJWithArrays_private(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,((PetscObject)A)->type_name,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_Sorted;
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 
815 PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_RowMerge(Mat A,Mat B,PetscReal fill,Mat C)
816 {
817   PetscErrorCode     ierr;
818   Mat_SeqAIJ         *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
819   const PetscInt     *ai=a->i,*bi=b->i,*aj=a->j,*bj=b->j,*inputi,*inputj,*inputcol,*inputcol_L1;
820   PetscInt           *ci,*cj,*outputj,worki_L1[9],worki_L2[9];
821   PetscInt           c_maxmem,a_maxrownnz=0,a_rownnz;
822   const PetscInt     workcol[8]={0,1,2,3,4,5,6,7};
823   const PetscInt     am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
824   const PetscInt     *brow_ptr[8],*brow_end[8];
825   PetscInt           window[8];
826   PetscInt           window_min,old_window_min,ci_nnz,outputi_nnz=0,L1_nrows,L2_nrows;
827   PetscInt           i,k,ndouble=0,L1_rowsleft,rowsleft;
828   PetscReal          afill;
829   PetscInt           *workj_L1,*workj_L2,*workj_L3;
830   PetscInt           L1_nnz,L2_nnz;
831 
832   /* Step 1: Get upper bound on memory required for allocation.
833              Because of the way virtual memory works,
834              only the memory pages that are actually needed will be physically allocated. */
835   PetscFunctionBegin;
836   ierr = PetscMalloc1(am+1,&ci);CHKERRQ(ierr);
837   for (i=0; i<am; i++) {
838     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 */
839     const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */
840     a_rownnz = 0;
841     for (k=0; k<anzi; ++k) {
842       a_rownnz += bi[acol[k]+1] - bi[acol[k]];
843       if (a_rownnz > bn) {
844         a_rownnz = bn;
845         break;
846       }
847     }
848     a_maxrownnz = PetscMax(a_maxrownnz, a_rownnz);
849   }
850   /* temporary work areas for merging rows */
851   ierr = PetscMalloc1(a_maxrownnz*8,&workj_L1);CHKERRQ(ierr);
852   ierr = PetscMalloc1(a_maxrownnz*8,&workj_L2);CHKERRQ(ierr);
853   ierr = PetscMalloc1(a_maxrownnz,&workj_L3);CHKERRQ(ierr);
854 
855   /* This should be enough for almost all matrices. If not, memory is reallocated later. */
856   c_maxmem = 8*(ai[am]+bi[bm]);
857   /* Step 2: Populate pattern for C */
858   ierr  = PetscMalloc1(c_maxmem,&cj);CHKERRQ(ierr);
859 
860   ci_nnz       = 0;
861   ci[0]        = 0;
862   worki_L1[0]  = 0;
863   worki_L2[0]  = 0;
864   for (i=0; i<am; i++) {
865     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 */
866     const PetscInt *acol = aj + ai[i];      /* column indices of nonzero entries in this row */
867     rowsleft             = anzi;
868     inputcol_L1          = acol;
869     L2_nnz               = 0;
870     L2_nrows             = 1;  /* Number of rows to be merged on Level 3. output of L3 already exists -> initial value 1   */
871     worki_L2[1]          = 0;
872     outputi_nnz          = 0;
873 
874     /* If the number of indices in C so far + the max number of columns in the next row > c_maxmem  -> allocate more memory */
875     while (ci_nnz+a_maxrownnz > c_maxmem) {
876       c_maxmem *= 2;
877       ndouble++;
878       ierr = PetscRealloc(sizeof(PetscInt)*c_maxmem,&cj);CHKERRQ(ierr);
879     }
880 
881     while (rowsleft) {
882       L1_rowsleft = PetscMin(64, rowsleft); /* In the inner loop max 64 rows of B can be merged */
883       L1_nrows    = 0;
884       L1_nnz      = 0;
885       inputcol    = inputcol_L1;
886       inputi      = bi;
887       inputj      = bj;
888 
889       /* The following macro is used to specialize for small rows in A.
890          This helps with compiler unrolling, improving performance substantially.
891           Input:  inputj   inputi  inputcol  bn
892           Output: outputj  outputi_nnz                       */
893        #define MatMatMultSymbolic_RowMergeMacro(ANNZ)                        \
894          window_min  = bn;                                                   \
895          outputi_nnz = 0;                                                    \
896          for (k=0; k<ANNZ; ++k) {                                            \
897            brow_ptr[k] = inputj + inputi[inputcol[k]];                       \
898            brow_end[k] = inputj + inputi[inputcol[k]+1];                     \
899            window[k]   = (brow_ptr[k] != brow_end[k]) ? *brow_ptr[k] : bn;   \
900            window_min  = PetscMin(window[k], window_min);                    \
901          }                                                                   \
902          while (window_min < bn) {                                           \
903            outputj[outputi_nnz++] = window_min;                              \
904            /* advance front and compute new minimum */                       \
905            old_window_min = window_min;                                      \
906            window_min = bn;                                                  \
907            for (k=0; k<ANNZ; ++k) {                                          \
908              if (window[k] == old_window_min) {                              \
909                brow_ptr[k]++;                                                \
910                window[k] = (brow_ptr[k] != brow_end[k]) ? *brow_ptr[k] : bn; \
911              }                                                               \
912              window_min = PetscMin(window[k], window_min);                   \
913            }                                                                 \
914          }
915 
916       /************** L E V E L  1 ***************/
917       /* Merge up to 8 rows of B to L1 work array*/
918       while (L1_rowsleft) {
919         outputi_nnz = 0;
920         if (anzi > 8)  outputj = workj_L1 + L1_nnz;     /* Level 1 rowmerge*/
921         else           outputj = cj + ci_nnz;           /* Merge directly to C */
922 
923         switch (L1_rowsleft) {
924         case 1:  brow_ptr[0] = inputj + inputi[inputcol[0]];
925                  brow_end[0] = inputj + inputi[inputcol[0]+1];
926                  for (; brow_ptr[0] != brow_end[0]; ++brow_ptr[0]) outputj[outputi_nnz++] = *brow_ptr[0]; /* copy row in b over */
927                  inputcol    += L1_rowsleft;
928                  rowsleft    -= L1_rowsleft;
929                  L1_rowsleft  = 0;
930                  break;
931         case 2:  MatMatMultSymbolic_RowMergeMacro(2);
932                  inputcol    += L1_rowsleft;
933                  rowsleft    -= L1_rowsleft;
934                  L1_rowsleft  = 0;
935                  break;
936         case 3: MatMatMultSymbolic_RowMergeMacro(3);
937                  inputcol    += L1_rowsleft;
938                  rowsleft    -= L1_rowsleft;
939                  L1_rowsleft  = 0;
940                  break;
941         case 4:  MatMatMultSymbolic_RowMergeMacro(4);
942                  inputcol    += L1_rowsleft;
943                  rowsleft    -= L1_rowsleft;
944                  L1_rowsleft  = 0;
945                  break;
946         case 5:  MatMatMultSymbolic_RowMergeMacro(5);
947                  inputcol    += L1_rowsleft;
948                  rowsleft    -= L1_rowsleft;
949                  L1_rowsleft  = 0;
950                  break;
951         case 6:  MatMatMultSymbolic_RowMergeMacro(6);
952                  inputcol    += L1_rowsleft;
953                  rowsleft    -= L1_rowsleft;
954                  L1_rowsleft  = 0;
955                  break;
956         case 7:  MatMatMultSymbolic_RowMergeMacro(7);
957                  inputcol    += L1_rowsleft;
958                  rowsleft    -= L1_rowsleft;
959                  L1_rowsleft  = 0;
960                  break;
961         default: MatMatMultSymbolic_RowMergeMacro(8);
962                  inputcol    += 8;
963                  rowsleft    -= 8;
964                  L1_rowsleft -= 8;
965                  break;
966         }
967         inputcol_L1           = inputcol;
968         L1_nnz               += outputi_nnz;
969         worki_L1[++L1_nrows]  = L1_nnz;
970       }
971 
972       /********************** L E V E L  2 ************************/
973       /* Merge from L1 work array to either C or to L2 work array */
974       if (anzi > 8) {
975         inputi      = worki_L1;
976         inputj      = workj_L1;
977         inputcol    = workcol;
978         outputi_nnz = 0;
979 
980         if (anzi <= 64) outputj = cj + ci_nnz;        /* Merge from L1 work array to C */
981         else            outputj = workj_L2 + L2_nnz;  /* Merge from L1 work array to L2 work array */
982 
983         switch (L1_nrows) {
984         case 1:  brow_ptr[0] = inputj + inputi[inputcol[0]];
985                  brow_end[0] = inputj + inputi[inputcol[0]+1];
986                  for (; brow_ptr[0] != brow_end[0]; ++brow_ptr[0]) outputj[outputi_nnz++] = *brow_ptr[0]; /* copy row in b over */
987                  break;
988         case 2:  MatMatMultSymbolic_RowMergeMacro(2); break;
989         case 3:  MatMatMultSymbolic_RowMergeMacro(3); break;
990         case 4:  MatMatMultSymbolic_RowMergeMacro(4); break;
991         case 5:  MatMatMultSymbolic_RowMergeMacro(5); break;
992         case 6:  MatMatMultSymbolic_RowMergeMacro(6); break;
993         case 7:  MatMatMultSymbolic_RowMergeMacro(7); break;
994         case 8:  MatMatMultSymbolic_RowMergeMacro(8); break;
995         default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatMatMult logic error: Not merging 1-8 rows from L1 work array!");
996         }
997         L2_nnz               += outputi_nnz;
998         worki_L2[++L2_nrows]  = L2_nnz;
999 
1000         /************************ L E V E L  3 **********************/
1001         /* Merge from L2 work array to either C or to L2 work array */
1002         if (anzi > 64 && (L2_nrows == 8 || rowsleft == 0)) {
1003           inputi      = worki_L2;
1004           inputj      = workj_L2;
1005           inputcol    = workcol;
1006           outputi_nnz = 0;
1007           if (rowsleft) outputj = workj_L3;
1008           else          outputj = cj + ci_nnz;
1009           switch (L2_nrows) {
1010           case 1:  brow_ptr[0] = inputj + inputi[inputcol[0]];
1011                    brow_end[0] = inputj + inputi[inputcol[0]+1];
1012                    for (; brow_ptr[0] != brow_end[0]; ++brow_ptr[0]) outputj[outputi_nnz++] = *brow_ptr[0]; /* copy row in b over */
1013                    break;
1014           case 2:  MatMatMultSymbolic_RowMergeMacro(2); break;
1015           case 3:  MatMatMultSymbolic_RowMergeMacro(3); break;
1016           case 4:  MatMatMultSymbolic_RowMergeMacro(4); break;
1017           case 5:  MatMatMultSymbolic_RowMergeMacro(5); break;
1018           case 6:  MatMatMultSymbolic_RowMergeMacro(6); break;
1019           case 7:  MatMatMultSymbolic_RowMergeMacro(7); break;
1020           case 8:  MatMatMultSymbolic_RowMergeMacro(8); break;
1021           default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatMatMult logic error: Not merging 1-8 rows from L2 work array!");
1022           }
1023           L2_nrows    = 1;
1024           L2_nnz      = outputi_nnz;
1025           worki_L2[1] = outputi_nnz;
1026           /* Copy to workj_L2 */
1027           if (rowsleft) {
1028             for (k=0; k<outputi_nnz; ++k)  workj_L2[k] = outputj[k];
1029           }
1030         }
1031       }
1032     }  /* while (rowsleft) */
1033 #undef MatMatMultSymbolic_RowMergeMacro
1034 
1035     /* terminate current row */
1036     ci_nnz += outputi_nnz;
1037     ci[i+1] = ci_nnz;
1038   }
1039 
1040   /* Step 3: Create the new symbolic matrix */
1041   ierr = MatSetSeqAIJWithArrays_private(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,((PetscObject)A)->type_name,C);CHKERRQ(ierr);
1042   ierr = MatSetBlockSizesFromMats(C,A,B);CHKERRQ(ierr);
1043 
1044   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
1045   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
1046   c          = (Mat_SeqAIJ*)(C->data);
1047   c->free_a  = PETSC_TRUE;
1048   c->free_ij = PETSC_TRUE;
1049   c->nonew   = 0;
1050 
1051   C->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Sorted;
1052 
1053   /* set MatInfo */
1054   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
1055   if (afill < 1.0) afill = 1.0;
1056   c->maxnz                     = ci[am];
1057   c->nz                        = ci[am];
1058   C->info.mallocs           = ndouble;
1059   C->info.fill_ratio_given  = fill;
1060   C->info.fill_ratio_needed = afill;
1061 
1062 #if defined(PETSC_USE_INFO)
1063   if (ci[am]) {
1064     ierr = PetscInfo3(C,"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);CHKERRQ(ierr);
1065     ierr = PetscInfo1(C,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
1066   } else {
1067     ierr = PetscInfo(C,"Empty matrix product\n");CHKERRQ(ierr);
1068   }
1069 #endif
1070 
1071   /* Step 4: Free temporary work areas */
1072   ierr = PetscFree(workj_L1);CHKERRQ(ierr);
1073   ierr = PetscFree(workj_L2);CHKERRQ(ierr);
1074   ierr = PetscFree(workj_L3);CHKERRQ(ierr);
1075   PetscFunctionReturn(0);
1076 }
1077 
1078 /* concatenate unique entries and then sort */
1079 PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Sorted(Mat A,Mat B,PetscReal fill,Mat C)
1080 {
1081   PetscErrorCode ierr;
1082   Mat_SeqAIJ     *a  = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
1083   const PetscInt *ai = a->i,*bi=b->i,*aj=a->j,*bj=b->j;
1084   PetscInt       *ci,*cj,bcol;
1085   PetscInt       am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
1086   PetscReal      afill;
1087   PetscInt       i,j,ndouble = 0;
1088   PetscSegBuffer seg,segrow;
1089   char           *seen;
1090 
1091   PetscFunctionBegin;
1092   ierr  = PetscMalloc1(am+1,&ci);CHKERRQ(ierr);
1093   ci[0] = 0;
1094 
1095   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
1096   ierr = PetscSegBufferCreate(sizeof(PetscInt),(PetscInt)(fill*(ai[am]+bi[bm])),&seg);CHKERRQ(ierr);
1097   ierr = PetscSegBufferCreate(sizeof(PetscInt),100,&segrow);CHKERRQ(ierr);
1098   ierr = PetscCalloc1(bn,&seen);CHKERRQ(ierr);
1099 
1100   /* Determine ci and cj */
1101   for (i=0; i<am; i++) {
1102     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 */
1103     const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */
1104     PetscInt packlen = 0,*PETSC_RESTRICT crow;
1105 
1106     /* Pack segrow */
1107     for (j=0; j<anzi; j++) {
1108       PetscInt brow = acol[j],bjstart = bi[brow],bjend = bi[brow+1],k;
1109       for (k=bjstart; k<bjend; k++) {
1110         bcol = bj[k];
1111         if (!seen[bcol]) { /* new entry */
1112           PetscInt *PETSC_RESTRICT slot;
1113           ierr = PetscSegBufferGetInts(segrow,1,&slot);CHKERRQ(ierr);
1114           *slot = bcol;
1115           seen[bcol] = 1;
1116           packlen++;
1117         }
1118       }
1119     }
1120 
1121     /* Check i-th diagonal entry */
1122     if (C->force_diagonals && !seen[i]) {
1123       PetscInt *PETSC_RESTRICT slot;
1124       ierr = PetscSegBufferGetInts(segrow,1,&slot);CHKERRQ(ierr);
1125       *slot   = i;
1126       seen[i] = 1;
1127       packlen++;
1128     }
1129 
1130     ierr = PetscSegBufferGetInts(seg,packlen,&crow);CHKERRQ(ierr);
1131     ierr = PetscSegBufferExtractTo(segrow,crow);CHKERRQ(ierr);
1132     ierr = PetscSortInt(packlen,crow);CHKERRQ(ierr);
1133     ci[i+1] = ci[i] + packlen;
1134     for (j=0; j<packlen; j++) seen[crow[j]] = 0;
1135   }
1136   ierr = PetscSegBufferDestroy(&segrow);CHKERRQ(ierr);
1137   ierr = PetscFree(seen);CHKERRQ(ierr);
1138 
1139   /* Column indices are in the segmented buffer */
1140   ierr = PetscSegBufferExtractAlloc(seg,&cj);CHKERRQ(ierr);
1141   ierr = PetscSegBufferDestroy(&seg);CHKERRQ(ierr);
1142 
1143   /* put together the new symbolic matrix */
1144   ierr = MatSetSeqAIJWithArrays_private(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,((PetscObject)A)->type_name,C);CHKERRQ(ierr);
1145   ierr = MatSetBlockSizesFromMats(C,A,B);CHKERRQ(ierr);
1146 
1147   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
1148   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
1149   c          = (Mat_SeqAIJ*)(C->data);
1150   c->free_a  = PETSC_TRUE;
1151   c->free_ij = PETSC_TRUE;
1152   c->nonew   = 0;
1153 
1154   C->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Sorted;
1155 
1156   /* set MatInfo */
1157   afill = (PetscReal)ci[am]/PetscMax(ai[am]+bi[bm],1) + 1.e-5;
1158   if (afill < 1.0) afill = 1.0;
1159   c->maxnz                  = ci[am];
1160   c->nz                     = ci[am];
1161   C->info.mallocs           = ndouble;
1162   C->info.fill_ratio_given  = fill;
1163   C->info.fill_ratio_needed = afill;
1164 
1165 #if defined(PETSC_USE_INFO)
1166   if (ci[am]) {
1167     ierr = PetscInfo3(C,"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);CHKERRQ(ierr);
1168     ierr = PetscInfo1(C,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr);
1169   } else {
1170     ierr = PetscInfo(C,"Empty matrix product\n");CHKERRQ(ierr);
1171   }
1172 #endif
1173   PetscFunctionReturn(0);
1174 }
1175 
1176 PetscErrorCode MatDestroy_SeqAIJ_MatMatMultTrans(void *data)
1177 {
1178   PetscErrorCode      ierr;
1179   Mat_MatMatTransMult *abt=(Mat_MatMatTransMult *)data;
1180 
1181   PetscFunctionBegin;
1182   ierr = MatTransposeColoringDestroy(&abt->matcoloring);CHKERRQ(ierr);
1183   ierr = MatDestroy(&abt->Bt_den);CHKERRQ(ierr);
1184   ierr = MatDestroy(&abt->ABt_den);CHKERRQ(ierr);
1185   ierr = PetscFree(abt);CHKERRQ(ierr);
1186   PetscFunctionReturn(0);
1187 }
1188 
1189 PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat C)
1190 {
1191   PetscErrorCode      ierr;
1192   Mat                 Bt;
1193   PetscInt            *bti,*btj;
1194   Mat_MatMatTransMult *abt;
1195   Mat_Product         *product = C->product;
1196   char                *alg;
1197 
1198   PetscFunctionBegin;
1199   if (!product) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing product struct");
1200   if (product->data) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Extra product struct not empty");
1201 
1202   /* create symbolic Bt */
1203   ierr = MatGetSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr);
1204   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,B->cmap->n,B->rmap->n,bti,btj,NULL,&Bt);CHKERRQ(ierr);
1205   ierr = MatSetBlockSizes(Bt,PetscAbs(A->cmap->bs),PetscAbs(B->cmap->bs));CHKERRQ(ierr);
1206   ierr = MatSetType(Bt,((PetscObject)A)->type_name);CHKERRQ(ierr);
1207 
1208   /* get symbolic C=A*Bt */
1209   ierr = PetscStrallocpy(product->alg,&alg);CHKERRQ(ierr);
1210   ierr = MatProductSetAlgorithm(C,"sorted");CHKERRQ(ierr); /* set algorithm for C = A*Bt */
1211   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,Bt,fill,C);CHKERRQ(ierr);
1212   ierr = MatProductSetAlgorithm(C,alg);CHKERRQ(ierr); /* resume original algorithm for ABt product */
1213   ierr = PetscFree(alg);CHKERRQ(ierr);
1214 
1215   /* create a supporting struct for reuse intermidiate dense matrices with matcoloring */
1216   ierr = PetscNew(&abt);CHKERRQ(ierr);
1217 
1218   product->data    = abt;
1219   product->destroy = MatDestroy_SeqAIJ_MatMatMultTrans;
1220 
1221   C->ops->mattransposemultnumeric = MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ;
1222 
1223   abt->usecoloring = PETSC_FALSE;
1224   ierr = PetscStrcmp(product->alg,"color",&abt->usecoloring);CHKERRQ(ierr);
1225   if (abt->usecoloring) {
1226     /* Create MatTransposeColoring from symbolic C=A*B^T */
1227     MatTransposeColoring matcoloring;
1228     MatColoring          coloring;
1229     ISColoring           iscoloring;
1230     Mat                  Bt_dense,C_dense;
1231 
1232     /* inode causes memory problem */
1233     ierr = MatSetOption(C,MAT_USE_INODES,PETSC_FALSE);CHKERRQ(ierr);
1234 
1235     ierr = MatColoringCreate(C,&coloring);CHKERRQ(ierr);
1236     ierr = MatColoringSetDistance(coloring,2);CHKERRQ(ierr);
1237     ierr = MatColoringSetType(coloring,MATCOLORINGSL);CHKERRQ(ierr);
1238     ierr = MatColoringSetFromOptions(coloring);CHKERRQ(ierr);
1239     ierr = MatColoringApply(coloring,&iscoloring);CHKERRQ(ierr);
1240     ierr = MatColoringDestroy(&coloring);CHKERRQ(ierr);
1241     ierr = MatTransposeColoringCreate(C,iscoloring,&matcoloring);CHKERRQ(ierr);
1242 
1243     abt->matcoloring = matcoloring;
1244 
1245     ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr);
1246 
1247     /* Create Bt_dense and C_dense = A*Bt_dense */
1248     ierr = MatCreate(PETSC_COMM_SELF,&Bt_dense);CHKERRQ(ierr);
1249     ierr = MatSetSizes(Bt_dense,A->cmap->n,matcoloring->ncolors,A->cmap->n,matcoloring->ncolors);CHKERRQ(ierr);
1250     ierr = MatSetType(Bt_dense,MATSEQDENSE);CHKERRQ(ierr);
1251     ierr = MatSeqDenseSetPreallocation(Bt_dense,NULL);CHKERRQ(ierr);
1252 
1253     Bt_dense->assembled = PETSC_TRUE;
1254     abt->Bt_den         = Bt_dense;
1255 
1256     ierr = MatCreate(PETSC_COMM_SELF,&C_dense);CHKERRQ(ierr);
1257     ierr = MatSetSizes(C_dense,A->rmap->n,matcoloring->ncolors,A->rmap->n,matcoloring->ncolors);CHKERRQ(ierr);
1258     ierr = MatSetType(C_dense,MATSEQDENSE);CHKERRQ(ierr);
1259     ierr = MatSeqDenseSetPreallocation(C_dense,NULL);CHKERRQ(ierr);
1260 
1261     Bt_dense->assembled = PETSC_TRUE;
1262     abt->ABt_den  = C_dense;
1263 
1264 #if defined(PETSC_USE_INFO)
1265     {
1266       Mat_SeqAIJ *c = (Mat_SeqAIJ*)C->data;
1267       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);
1268     }
1269 #endif
1270   }
1271   /* clean up */
1272   ierr = MatDestroy(&Bt);CHKERRQ(ierr);
1273   ierr = MatRestoreSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr);
1274   PetscFunctionReturn(0);
1275 }
1276 
1277 PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
1278 {
1279   PetscErrorCode      ierr;
1280   Mat_SeqAIJ          *a   =(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data;
1281   PetscInt            *ai  =a->i,*aj=a->j,*bi=b->i,*bj=b->j,anzi,bnzj,nexta,nextb,*acol,*bcol,brow;
1282   PetscInt            cm   =C->rmap->n,*ci=c->i,*cj=c->j,i,j,cnzi,*ccol;
1283   PetscLogDouble      flops=0.0;
1284   MatScalar           *aa  =a->a,*aval,*ba=b->a,*bval,*ca,*cval;
1285   Mat_MatMatTransMult *abt;
1286   Mat_Product         *product = C->product;
1287 
1288   PetscFunctionBegin;
1289   if (!product) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing product struct");
1290   abt = (Mat_MatMatTransMult *)product->data;
1291   if (!abt) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing product struct");
1292   /* clear old values in C */
1293   if (!c->a) {
1294     ierr      = PetscCalloc1(ci[cm]+1,&ca);CHKERRQ(ierr);
1295     c->a      = ca;
1296     c->free_a = PETSC_TRUE;
1297   } else {
1298     ca =  c->a;
1299     ierr = PetscArrayzero(ca,ci[cm]+1);CHKERRQ(ierr);
1300   }
1301 
1302   if (abt->usecoloring) {
1303     MatTransposeColoring matcoloring = abt->matcoloring;
1304     Mat                  Bt_dense,C_dense = abt->ABt_den;
1305 
1306     /* Get Bt_dense by Apply MatTransposeColoring to B */
1307     Bt_dense = abt->Bt_den;
1308     ierr = MatTransColoringApplySpToDen(matcoloring,B,Bt_dense);CHKERRQ(ierr);
1309 
1310     /* C_dense = A*Bt_dense */
1311     ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,Bt_dense,C_dense);CHKERRQ(ierr);
1312 
1313     /* Recover C from C_dense */
1314     ierr = MatTransColoringApplyDenToSp(matcoloring,C_dense,C);CHKERRQ(ierr);
1315     PetscFunctionReturn(0);
1316   }
1317 
1318   for (i=0; i<cm; i++) {
1319     anzi = ai[i+1] - ai[i];
1320     acol = aj + ai[i];
1321     aval = aa + ai[i];
1322     cnzi = ci[i+1] - ci[i];
1323     ccol = cj + ci[i];
1324     cval = ca + ci[i];
1325     for (j=0; j<cnzi; j++) {
1326       brow = ccol[j];
1327       bnzj = bi[brow+1] - bi[brow];
1328       bcol = bj + bi[brow];
1329       bval = ba + bi[brow];
1330 
1331       /* perform sparse inner-product c(i,j)=A[i,:]*B[j,:]^T */
1332       nexta = 0; nextb = 0;
1333       while (nexta<anzi && nextb<bnzj) {
1334         while (nexta < anzi && acol[nexta] < bcol[nextb]) nexta++;
1335         if (nexta == anzi) break;
1336         while (nextb < bnzj && acol[nexta] > bcol[nextb]) nextb++;
1337         if (nextb == bnzj) break;
1338         if (acol[nexta] == bcol[nextb]) {
1339           cval[j] += aval[nexta]*bval[nextb];
1340           nexta++; nextb++;
1341           flops += 2;
1342         }
1343       }
1344     }
1345   }
1346   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1347   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1348   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
1349   PetscFunctionReturn(0);
1350 }
1351 
1352 PetscErrorCode MatDestroy_SeqAIJ_MatTransMatMult(void *data)
1353 {
1354   PetscErrorCode      ierr;
1355   Mat_MatTransMatMult *atb = (Mat_MatTransMatMult*)data;
1356 
1357   PetscFunctionBegin;
1358   ierr = MatDestroy(&atb->At);CHKERRQ(ierr);
1359   if (atb->destroy) {
1360     ierr = (*atb->destroy)(atb->data);CHKERRQ(ierr);
1361   }
1362   ierr = PetscFree(atb);CHKERRQ(ierr);
1363   PetscFunctionReturn(0);
1364 }
1365 
1366 PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat C)
1367 {
1368   PetscErrorCode ierr;
1369   Mat            At = NULL;
1370   PetscInt       *ati,*atj;
1371   Mat_Product    *product = C->product;
1372   PetscBool      flg,def,square;
1373 
1374   PetscFunctionBegin;
1375   MatCheckProduct(C,4);
1376   square = (PetscBool)(A == B && A->symmetric && A->symmetric_set);
1377   /* outerproduct */
1378   ierr = PetscStrcmp(product->alg,"outerproduct",&flg);CHKERRQ(ierr);
1379   if (flg) {
1380     /* create symbolic At */
1381     if (!square) {
1382       ierr = MatGetSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr);
1383       ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,A->cmap->n,A->rmap->n,ati,atj,NULL,&At);CHKERRQ(ierr);
1384       ierr = MatSetBlockSizes(At,PetscAbs(A->cmap->bs),PetscAbs(B->cmap->bs));CHKERRQ(ierr);
1385       ierr = MatSetType(At,((PetscObject)A)->type_name);CHKERRQ(ierr);
1386     }
1387     /* get symbolic C=At*B */
1388     ierr = MatProductSetAlgorithm(C,"sorted");CHKERRQ(ierr);
1389     ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(square ? A : At,B,fill,C);CHKERRQ(ierr);
1390 
1391     /* clean up */
1392     if (!square) {
1393       ierr = MatDestroy(&At);CHKERRQ(ierr);
1394       ierr = MatRestoreSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr);
1395     }
1396 
1397     C->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ; /* outerproduct */
1398     ierr = MatProductSetAlgorithm(C,"outerproduct");CHKERRQ(ierr);
1399     PetscFunctionReturn(0);
1400   }
1401 
1402   /* matmatmult */
1403   ierr = PetscStrcmp(product->alg,"default",&def);CHKERRQ(ierr);
1404   ierr = PetscStrcmp(product->alg,"at*b",&flg);CHKERRQ(ierr);
1405   if (flg || def) {
1406     Mat_MatTransMatMult *atb;
1407 
1408     if (product->data) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Extra product struct not empty");
1409     ierr = PetscNew(&atb);CHKERRQ(ierr);
1410     if (!square) {
1411       ierr = MatTranspose_SeqAIJ(A,MAT_INITIAL_MATRIX,&At);CHKERRQ(ierr);
1412     }
1413     ierr = MatProductSetAlgorithm(C,"sorted");CHKERRQ(ierr);
1414     ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(square ? A : At,B,fill,C);CHKERRQ(ierr);
1415     ierr = MatProductSetAlgorithm(C,"at*b");CHKERRQ(ierr);
1416     product->data    = atb;
1417     product->destroy = MatDestroy_SeqAIJ_MatTransMatMult;
1418     atb->At          = At;
1419     atb->updateAt    = PETSC_FALSE; /* because At is computed here */
1420 
1421     C->ops->mattransposemultnumeric = NULL; /* see MatProductNumeric_AtB_SeqAIJ_SeqAIJ */
1422     PetscFunctionReturn(0);
1423   }
1424 
1425   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Mat Product Algorithm is not supported");
1426 }
1427 
1428 PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
1429 {
1430   PetscErrorCode ierr;
1431   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data;
1432   PetscInt       am=A->rmap->n,anzi,*ai=a->i,*aj=a->j,*bi=b->i,*bj,bnzi,nextb;
1433   PetscInt       cm=C->rmap->n,*ci=c->i,*cj=c->j,crow,*cjj,i,j,k;
1434   PetscLogDouble flops=0.0;
1435   MatScalar      *aa=a->a,*ba,*ca,*caj;
1436 
1437   PetscFunctionBegin;
1438   if (!c->a) {
1439     ierr = PetscCalloc1(ci[cm]+1,&ca);CHKERRQ(ierr);
1440 
1441     c->a      = ca;
1442     c->free_a = PETSC_TRUE;
1443   } else {
1444     ca   = c->a;
1445     ierr = PetscArrayzero(ca,ci[cm]);CHKERRQ(ierr);
1446   }
1447 
1448   /* compute A^T*B using outer product (A^T)[:,i]*B[i,:] */
1449   for (i=0; i<am; i++) {
1450     bj   = b->j + bi[i];
1451     ba   = b->a + bi[i];
1452     bnzi = bi[i+1] - bi[i];
1453     anzi = ai[i+1] - ai[i];
1454     for (j=0; j<anzi; j++) {
1455       nextb = 0;
1456       crow  = *aj++;
1457       cjj   = cj + ci[crow];
1458       caj   = ca + ci[crow];
1459       /* perform sparse axpy operation.  Note cjj includes bj. */
1460       for (k=0; nextb<bnzi; k++) {
1461         if (cjj[k] == *(bj+nextb)) { /* ccol == bcol */
1462           caj[k] += (*aa)*(*(ba+nextb));
1463           nextb++;
1464         }
1465       }
1466       flops += 2*bnzi;
1467       aa++;
1468     }
1469   }
1470 
1471   /* Assemble the final matrix and clean up */
1472   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1473   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1474   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
1475   PetscFunctionReturn(0);
1476 }
1477 
1478 PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqDense(Mat A,Mat B,PetscReal fill,Mat C)
1479 {
1480   PetscErrorCode ierr;
1481 
1482   PetscFunctionBegin;
1483   ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,0.0,C);CHKERRQ(ierr);
1484   C->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqDense;
1485   PetscFunctionReturn(0);
1486 }
1487 
1488 PETSC_INTERN PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat A,Mat B,Mat C,const PetscBool add)
1489 {
1490   Mat_SeqAIJ        *a=(Mat_SeqAIJ*)A->data;
1491   Mat_SeqDense      *bd=(Mat_SeqDense*)B->data;
1492   Mat_SeqDense      *cd=(Mat_SeqDense*)C->data;
1493   PetscErrorCode    ierr;
1494   PetscScalar       *c,r1,r2,r3,r4,*c1,*c2,*c3,*c4;
1495   const PetscScalar *aa,*b,*b1,*b2,*b3,*b4,*av;
1496   const PetscInt    *aj;
1497   PetscInt          cm=C->rmap->n,cn=B->cmap->n,bm=bd->lda,am=A->rmap->n;
1498   PetscInt          clda=cd->lda;
1499   PetscInt          am4=4*clda,bm4=4*bm,col,i,j,n;
1500 
1501   PetscFunctionBegin;
1502   if (!cm || !cn) PetscFunctionReturn(0);
1503   ierr = MatSeqAIJGetArrayRead(A,&av);CHKERRQ(ierr);
1504   if (add) {
1505     ierr = MatDenseGetArray(C,&c);CHKERRQ(ierr);
1506   } else {
1507     ierr = MatDenseGetArrayWrite(C,&c);CHKERRQ(ierr);
1508   }
1509   ierr = MatDenseGetArrayRead(B,&b);CHKERRQ(ierr);
1510   b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;
1511   c1 = c; c2 = c1 + clda; c3 = c2 + clda; c4 = c3 + clda;
1512   for (col=0; col<(cn/4)*4; col += 4) {  /* over columns of C */
1513     for (i=0; i<am; i++) {        /* over rows of A in those columns */
1514       r1 = r2 = r3 = r4 = 0.0;
1515       n  = a->i[i+1] - a->i[i];
1516       aj = a->j + a->i[i];
1517       aa = av + a->i[i];
1518       for (j=0; j<n; j++) {
1519         const PetscScalar aatmp = aa[j];
1520         const PetscInt    ajtmp = aj[j];
1521         r1 += aatmp*b1[ajtmp];
1522         r2 += aatmp*b2[ajtmp];
1523         r3 += aatmp*b3[ajtmp];
1524         r4 += aatmp*b4[ajtmp];
1525       }
1526       if (add) {
1527         c1[i] += r1;
1528         c2[i] += r2;
1529         c3[i] += r3;
1530         c4[i] += r4;
1531       } else {
1532         c1[i] = r1;
1533         c2[i] = r2;
1534         c3[i] = r3;
1535         c4[i] = r4;
1536       }
1537     }
1538     b1 += bm4; b2 += bm4; b3 += bm4; b4 += bm4;
1539     c1 += am4; c2 += am4; c3 += am4; c4 += am4;
1540   }
1541   /* process remaining columns */
1542   if (col != cn) {
1543     PetscInt rc = cn-col;
1544 
1545     if (rc == 1) {
1546       for (i=0; i<am; i++) {
1547         r1 = 0.0;
1548         n  = a->i[i+1] - a->i[i];
1549         aj = a->j + a->i[i];
1550         aa = av + a->i[i];
1551         for (j=0; j<n; j++) r1 += aa[j]*b1[aj[j]];
1552         if (add) c1[i] += r1;
1553         else c1[i] = r1;
1554       }
1555     } else if (rc == 2) {
1556       for (i=0; i<am; i++) {
1557         r1 = r2 = 0.0;
1558         n  = a->i[i+1] - a->i[i];
1559         aj = a->j + a->i[i];
1560         aa = av + a->i[i];
1561         for (j=0; j<n; j++) {
1562           const PetscScalar aatmp = aa[j];
1563           const PetscInt    ajtmp = aj[j];
1564           r1 += aatmp*b1[ajtmp];
1565           r2 += aatmp*b2[ajtmp];
1566         }
1567         if (add) {
1568           c1[i] += r1;
1569           c2[i] += r2;
1570         } else {
1571           c1[i] = r1;
1572           c2[i] = r2;
1573         }
1574       }
1575     } else {
1576       for (i=0; i<am; i++) {
1577         r1 = r2 = r3 = 0.0;
1578         n  = a->i[i+1] - a->i[i];
1579         aj = a->j + a->i[i];
1580         aa = av + a->i[i];
1581         for (j=0; j<n; j++) {
1582           const PetscScalar aatmp = aa[j];
1583           const PetscInt    ajtmp = aj[j];
1584           r1 += aatmp*b1[ajtmp];
1585           r2 += aatmp*b2[ajtmp];
1586           r3 += aatmp*b3[ajtmp];
1587         }
1588         if (add) {
1589           c1[i] += r1;
1590           c2[i] += r2;
1591           c3[i] += r3;
1592         } else {
1593           c1[i] = r1;
1594           c2[i] = r2;
1595           c3[i] = r3;
1596         }
1597       }
1598     }
1599   }
1600   ierr = PetscLogFlops(cn*(2.0*a->nz));CHKERRQ(ierr);
1601   if (add) {
1602     ierr = MatDenseRestoreArray(C,&c);CHKERRQ(ierr);
1603   } else {
1604     ierr = MatDenseRestoreArrayWrite(C,&c);CHKERRQ(ierr);
1605   }
1606   ierr = MatDenseRestoreArrayRead(B,&b);CHKERRQ(ierr);
1607   ierr = MatSeqAIJRestoreArrayRead(A,&av);CHKERRQ(ierr);
1608   PetscFunctionReturn(0);
1609 }
1610 
1611 PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
1612 {
1613   PetscErrorCode ierr;
1614 
1615   PetscFunctionBegin;
1616   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);
1617   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);
1618   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);
1619 
1620   ierr = MatMatMultNumericAdd_SeqAIJ_SeqDense(A,B,C,PETSC_FALSE);CHKERRQ(ierr);
1621   PetscFunctionReturn(0);
1622 }
1623 
1624 /* ------------------------------------------------------- */
1625 static PetscErrorCode MatProductSetFromOptions_SeqAIJ_SeqDense_AB(Mat C)
1626 {
1627   PetscFunctionBegin;
1628   C->ops->matmultsymbolic = MatMatMultSymbolic_SeqAIJ_SeqDense;
1629   C->ops->productsymbolic = MatProductSymbolic_AB;
1630   PetscFunctionReturn(0);
1631 }
1632 
1633 PETSC_INTERN PetscErrorCode MatTMatTMultSymbolic_SeqAIJ_SeqDense(Mat,Mat,PetscReal,Mat);
1634 
1635 static PetscErrorCode MatProductSetFromOptions_SeqAIJ_SeqDense_AtB(Mat C)
1636 {
1637   PetscFunctionBegin;
1638   C->ops->transposematmultsymbolic = MatTMatTMultSymbolic_SeqAIJ_SeqDense;
1639   C->ops->productsymbolic          = MatProductSymbolic_AtB;
1640   PetscFunctionReturn(0);
1641 }
1642 
1643 static PetscErrorCode MatProductSetFromOptions_SeqAIJ_SeqDense_ABt(Mat C)
1644 {
1645   PetscFunctionBegin;
1646   C->ops->mattransposemultsymbolic = MatTMatTMultSymbolic_SeqAIJ_SeqDense;
1647   C->ops->productsymbolic          = MatProductSymbolic_ABt;
1648   PetscFunctionReturn(0);
1649 }
1650 
1651 PETSC_INTERN PetscErrorCode MatProductSetFromOptions_SeqAIJ_SeqDense(Mat C)
1652 {
1653   PetscErrorCode ierr;
1654   Mat_Product    *product = C->product;
1655 
1656   PetscFunctionBegin;
1657   switch (product->type) {
1658   case MATPRODUCT_AB:
1659     ierr = MatProductSetFromOptions_SeqAIJ_SeqDense_AB(C);CHKERRQ(ierr);
1660     break;
1661   case MATPRODUCT_AtB:
1662     ierr = MatProductSetFromOptions_SeqAIJ_SeqDense_AtB(C);CHKERRQ(ierr);
1663     break;
1664   case MATPRODUCT_ABt:
1665     ierr = MatProductSetFromOptions_SeqAIJ_SeqDense_ABt(C);CHKERRQ(ierr);
1666     break;
1667   default:
1668     break;
1669   }
1670   PetscFunctionReturn(0);
1671 }
1672 /* ------------------------------------------------------- */
1673 static PetscErrorCode MatProductSetFromOptions_SeqXBAIJ_SeqDense_AB(Mat C)
1674 {
1675   PetscErrorCode ierr;
1676   Mat_Product    *product = C->product;
1677   Mat            A = product->A;
1678   PetscBool      baij;
1679 
1680   PetscFunctionBegin;
1681   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQBAIJ,&baij);CHKERRQ(ierr);
1682   if (!baij) { /* A is seqsbaij */
1683     PetscBool sbaij;
1684     ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQSBAIJ,&sbaij);CHKERRQ(ierr);
1685     if (!sbaij) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_WRONGSTATE,"Mat must be either seqbaij or seqsbaij format");
1686 
1687     C->ops->matmultsymbolic = MatMatMultSymbolic_SeqSBAIJ_SeqDense;
1688   } else { /* A is seqbaij */
1689     C->ops->matmultsymbolic = MatMatMultSymbolic_SeqBAIJ_SeqDense;
1690   }
1691 
1692   C->ops->productsymbolic = MatProductSymbolic_AB;
1693   PetscFunctionReturn(0);
1694 }
1695 
1696 PETSC_INTERN PetscErrorCode MatProductSetFromOptions_SeqXBAIJ_SeqDense(Mat C)
1697 {
1698   PetscErrorCode ierr;
1699   Mat_Product    *product = C->product;
1700 
1701   PetscFunctionBegin;
1702   MatCheckProduct(C,1);
1703   if (!product->A) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing A");
1704   if (product->type == MATPRODUCT_AB || (product->type == MATPRODUCT_AtB && product->A->symmetric)) {
1705     ierr = MatProductSetFromOptions_SeqXBAIJ_SeqDense_AB(C);CHKERRQ(ierr);
1706   }
1707   PetscFunctionReturn(0);
1708 }
1709 
1710 /* ------------------------------------------------------- */
1711 static PetscErrorCode MatProductSetFromOptions_SeqDense_SeqAIJ_AB(Mat C)
1712 {
1713   PetscFunctionBegin;
1714   C->ops->matmultsymbolic = MatMatMultSymbolic_SeqDense_SeqAIJ;
1715   C->ops->productsymbolic = MatProductSymbolic_AB;
1716   PetscFunctionReturn(0);
1717 }
1718 
1719 PETSC_INTERN PetscErrorCode MatProductSetFromOptions_SeqDense_SeqAIJ(Mat C)
1720 {
1721   PetscErrorCode ierr;
1722   Mat_Product    *product = C->product;
1723 
1724   PetscFunctionBegin;
1725   if (product->type == MATPRODUCT_AB) {
1726     ierr = MatProductSetFromOptions_SeqDense_SeqAIJ_AB(C);CHKERRQ(ierr);
1727   }
1728   PetscFunctionReturn(0);
1729 }
1730 /* ------------------------------------------------------- */
1731 
1732 PetscErrorCode  MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring coloring,Mat B,Mat Btdense)
1733 {
1734   PetscErrorCode ierr;
1735   Mat_SeqAIJ     *b       = (Mat_SeqAIJ*)B->data;
1736   Mat_SeqDense   *btdense = (Mat_SeqDense*)Btdense->data;
1737   PetscInt       *bi      = b->i,*bj=b->j;
1738   PetscInt       m        = Btdense->rmap->n,n=Btdense->cmap->n,j,k,l,col,anz,*btcol,brow,ncolumns;
1739   MatScalar      *btval,*btval_den,*ba=b->a;
1740   PetscInt       *columns=coloring->columns,*colorforcol=coloring->colorforcol,ncolors=coloring->ncolors;
1741 
1742   PetscFunctionBegin;
1743   btval_den=btdense->v;
1744   ierr     = PetscArrayzero(btval_den,m*n);CHKERRQ(ierr);
1745   for (k=0; k<ncolors; k++) {
1746     ncolumns = coloring->ncolumns[k];
1747     for (l=0; l<ncolumns; l++) { /* insert a row of B to a column of Btdense */
1748       col   = *(columns + colorforcol[k] + l);
1749       btcol = bj + bi[col];
1750       btval = ba + bi[col];
1751       anz   = bi[col+1] - bi[col];
1752       for (j=0; j<anz; j++) {
1753         brow            = btcol[j];
1754         btval_den[brow] = btval[j];
1755       }
1756     }
1757     btval_den += m;
1758   }
1759   PetscFunctionReturn(0);
1760 }
1761 
1762 PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring matcoloring,Mat Cden,Mat Csp)
1763 {
1764   PetscErrorCode    ierr;
1765   Mat_SeqAIJ        *csp = (Mat_SeqAIJ*)Csp->data;
1766   const PetscScalar *ca_den,*ca_den_ptr;
1767   PetscScalar       *ca=csp->a;
1768   PetscInt          k,l,m=Cden->rmap->n,ncolors=matcoloring->ncolors;
1769   PetscInt          brows=matcoloring->brows,*den2sp=matcoloring->den2sp;
1770   PetscInt          nrows,*row,*idx;
1771   PetscInt          *rows=matcoloring->rows,*colorforrow=matcoloring->colorforrow;
1772 
1773   PetscFunctionBegin;
1774   ierr   = MatDenseGetArrayRead(Cden,&ca_den);CHKERRQ(ierr);
1775 
1776   if (brows > 0) {
1777     PetscInt *lstart,row_end,row_start;
1778     lstart = matcoloring->lstart;
1779     ierr = PetscArrayzero(lstart,ncolors);CHKERRQ(ierr);
1780 
1781     row_end = brows;
1782     if (row_end > m) row_end = m;
1783     for (row_start=0; row_start<m; row_start+=brows) { /* loop over row blocks of Csp */
1784       ca_den_ptr = ca_den;
1785       for (k=0; k<ncolors; k++) { /* loop over colors (columns of Cden) */
1786         nrows = matcoloring->nrows[k];
1787         row   = rows  + colorforrow[k];
1788         idx   = den2sp + colorforrow[k];
1789         for (l=lstart[k]; l<nrows; l++) {
1790           if (row[l] >= row_end) {
1791             lstart[k] = l;
1792             break;
1793           } else {
1794             ca[idx[l]] = ca_den_ptr[row[l]];
1795           }
1796         }
1797         ca_den_ptr += m;
1798       }
1799       row_end += brows;
1800       if (row_end > m) row_end = m;
1801     }
1802   } else { /* non-blocked impl: loop over columns of Csp - slow if Csp is large */
1803     ca_den_ptr = ca_den;
1804     for (k=0; k<ncolors; k++) {
1805       nrows = matcoloring->nrows[k];
1806       row   = rows  + colorforrow[k];
1807       idx   = den2sp + colorforrow[k];
1808       for (l=0; l<nrows; l++) {
1809         ca[idx[l]] = ca_den_ptr[row[l]];
1810       }
1811       ca_den_ptr += m;
1812     }
1813   }
1814 
1815   ierr = MatDenseRestoreArrayRead(Cden,&ca_den);CHKERRQ(ierr);
1816 #if defined(PETSC_USE_INFO)
1817   if (matcoloring->brows > 0) {
1818     ierr = PetscInfo1(Csp,"Loop over %D row blocks for den2sp\n",brows);CHKERRQ(ierr);
1819   } else {
1820     ierr = PetscInfo(Csp,"Loop over colors/columns of Cden, inefficient for large sparse matrix product \n");CHKERRQ(ierr);
1821   }
1822 #endif
1823   PetscFunctionReturn(0);
1824 }
1825 
1826 PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat mat,ISColoring iscoloring,MatTransposeColoring c)
1827 {
1828   PetscErrorCode ierr;
1829   PetscInt       i,n,nrows,Nbs,j,k,m,ncols,col,cm;
1830   const PetscInt *is,*ci,*cj,*row_idx;
1831   PetscInt       nis = iscoloring->n,*rowhit,bs = 1;
1832   IS             *isa;
1833   Mat_SeqAIJ     *csp = (Mat_SeqAIJ*)mat->data;
1834   PetscInt       *colorforrow,*rows,*rows_i,*idxhit,*spidx,*den2sp,*den2sp_i;
1835   PetscInt       *colorforcol,*columns,*columns_i,brows;
1836   PetscBool      flg;
1837 
1838   PetscFunctionBegin;
1839   ierr = ISColoringGetIS(iscoloring,PETSC_USE_POINTER,PETSC_IGNORE,&isa);CHKERRQ(ierr);
1840 
1841   /* bs >1 is not being tested yet! */
1842   Nbs       = mat->cmap->N/bs;
1843   c->M      = mat->rmap->N/bs;  /* set total rows, columns and local rows */
1844   c->N      = Nbs;
1845   c->m      = c->M;
1846   c->rstart = 0;
1847   c->brows  = 100;
1848 
1849   c->ncolors = nis;
1850   ierr = PetscMalloc3(nis,&c->ncolumns,nis,&c->nrows,nis+1,&colorforrow);CHKERRQ(ierr);
1851   ierr = PetscMalloc1(csp->nz+1,&rows);CHKERRQ(ierr);
1852   ierr = PetscMalloc1(csp->nz+1,&den2sp);CHKERRQ(ierr);
1853 
1854   brows = c->brows;
1855   ierr = PetscOptionsGetInt(NULL,NULL,"-matden2sp_brows",&brows,&flg);CHKERRQ(ierr);
1856   if (flg) c->brows = brows;
1857   if (brows > 0) {
1858     ierr = PetscMalloc1(nis+1,&c->lstart);CHKERRQ(ierr);
1859   }
1860 
1861   colorforrow[0] = 0;
1862   rows_i         = rows;
1863   den2sp_i       = den2sp;
1864 
1865   ierr = PetscMalloc1(nis+1,&colorforcol);CHKERRQ(ierr);
1866   ierr = PetscMalloc1(Nbs+1,&columns);CHKERRQ(ierr);
1867 
1868   colorforcol[0] = 0;
1869   columns_i      = columns;
1870 
1871   /* get column-wise storage of mat */
1872   ierr = MatGetColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr);
1873 
1874   cm   = c->m;
1875   ierr = PetscMalloc1(cm+1,&rowhit);CHKERRQ(ierr);
1876   ierr = PetscMalloc1(cm+1,&idxhit);CHKERRQ(ierr);
1877   for (i=0; i<nis; i++) { /* loop over color */
1878     ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr);
1879     ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr);
1880 
1881     c->ncolumns[i] = n;
1882     if (n) {
1883       ierr = PetscArraycpy(columns_i,is,n);CHKERRQ(ierr);
1884     }
1885     colorforcol[i+1] = colorforcol[i] + n;
1886     columns_i       += n;
1887 
1888     /* fast, crude version requires O(N*N) work */
1889     ierr = PetscArrayzero(rowhit,cm);CHKERRQ(ierr);
1890 
1891     for (j=0; j<n; j++) { /* loop over columns*/
1892       col     = is[j];
1893       row_idx = cj + ci[col];
1894       m       = ci[col+1] - ci[col];
1895       for (k=0; k<m; k++) { /* loop over columns marking them in rowhit */
1896         idxhit[*row_idx]   = spidx[ci[col] + k];
1897         rowhit[*row_idx++] = col + 1;
1898       }
1899     }
1900     /* count the number of hits */
1901     nrows = 0;
1902     for (j=0; j<cm; j++) {
1903       if (rowhit[j]) nrows++;
1904     }
1905     c->nrows[i]      = nrows;
1906     colorforrow[i+1] = colorforrow[i] + nrows;
1907 
1908     nrows = 0;
1909     for (j=0; j<cm; j++) { /* loop over rows */
1910       if (rowhit[j]) {
1911         rows_i[nrows]   = j;
1912         den2sp_i[nrows] = idxhit[j];
1913         nrows++;
1914       }
1915     }
1916     den2sp_i += nrows;
1917 
1918     ierr    = ISRestoreIndices(isa[i],&is);CHKERRQ(ierr);
1919     rows_i += nrows;
1920   }
1921   ierr = MatRestoreColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr);
1922   ierr = PetscFree(rowhit);CHKERRQ(ierr);
1923   ierr = ISColoringRestoreIS(iscoloring,PETSC_USE_POINTER,&isa);CHKERRQ(ierr);
1924   if (csp->nz != colorforrow[nis]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"csp->nz %d != colorforrow[nis] %d",csp->nz,colorforrow[nis]);
1925 
1926   c->colorforrow = colorforrow;
1927   c->rows        = rows;
1928   c->den2sp      = den2sp;
1929   c->colorforcol = colorforcol;
1930   c->columns     = columns;
1931 
1932   ierr = PetscFree(idxhit);CHKERRQ(ierr);
1933   PetscFunctionReturn(0);
1934 }
1935 
1936 /* --------------------------------------------------------------- */
1937 static PetscErrorCode MatProductNumeric_AtB_SeqAIJ_SeqAIJ(Mat C)
1938 {
1939   PetscErrorCode ierr;
1940   Mat_Product    *product = C->product;
1941   Mat            A=product->A,B=product->B;
1942 
1943   PetscFunctionBegin;
1944   if (C->ops->mattransposemultnumeric) {
1945     /* Alg: "outerproduct" */
1946     ierr = (*C->ops->mattransposemultnumeric)(A,B,C);CHKERRQ(ierr);
1947   } else {
1948     /* Alg: "matmatmult" -- C = At*B */
1949     Mat_MatTransMatMult *atb = (Mat_MatTransMatMult *)product->data;
1950     Mat                 At;
1951 
1952     if (!atb) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing product struct");
1953     At = atb->At;
1954     if (atb->updateAt && At) { /* At is computed in MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ() */
1955       ierr = MatTranspose_SeqAIJ(A,MAT_REUSE_MATRIX,&At);CHKERRQ(ierr);
1956     }
1957     ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ(At ? At : A,B,C);CHKERRQ(ierr);
1958     atb->updateAt = PETSC_TRUE;
1959   }
1960   PetscFunctionReturn(0);
1961 }
1962 
1963 static PetscErrorCode MatProductSymbolic_AtB_SeqAIJ_SeqAIJ(Mat C)
1964 {
1965   PetscErrorCode ierr;
1966   Mat_Product    *product = C->product;
1967   Mat            A=product->A,B=product->B;
1968   PetscReal      fill=product->fill;
1969 
1970   PetscFunctionBegin;
1971   ierr = MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
1972 
1973   C->ops->productnumeric = MatProductNumeric_AtB_SeqAIJ_SeqAIJ;
1974   PetscFunctionReturn(0);
1975 }
1976 
1977 /* --------------------------------------------------------------- */
1978 static PetscErrorCode MatProductSetFromOptions_SeqAIJ_AB(Mat C)
1979 {
1980   PetscErrorCode ierr;
1981   Mat_Product    *product = C->product;
1982   PetscInt       alg = 0; /* default algorithm */
1983   PetscBool      flg = PETSC_FALSE;
1984 #if !defined(PETSC_HAVE_HYPRE)
1985   const char     *algTypes[7] = {"sorted","scalable","scalable_fast","heap","btheap","llcondensed","rowmerge"};
1986   PetscInt       nalg = 7;
1987 #else
1988   const char     *algTypes[8] = {"sorted","scalable","scalable_fast","heap","btheap","llcondensed","rowmerge","hypre"};
1989   PetscInt       nalg = 8;
1990 #endif
1991 
1992   PetscFunctionBegin;
1993   /* Set default algorithm */
1994   ierr = PetscStrcmp(C->product->alg,"default",&flg);CHKERRQ(ierr);
1995   if (flg) {
1996     ierr = MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);CHKERRQ(ierr);
1997   }
1998 
1999   /* Get runtime option */
2000   if (product->api_user) {
2001     ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatMatMult","Mat");CHKERRQ(ierr);
2002     ierr = PetscOptionsEList("-matmatmult_via","Algorithmic approach","MatMatMult",algTypes,nalg,algTypes[0],&alg,&flg);CHKERRQ(ierr);
2003     ierr = PetscOptionsEnd();CHKERRQ(ierr);
2004   } else {
2005     ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatProduct_AB","Mat");CHKERRQ(ierr);
2006     ierr = PetscOptionsEList("-matproduct_ab_via","Algorithmic approach","MatProduct_AB",algTypes,nalg,algTypes[0],&alg,&flg);CHKERRQ(ierr);
2007     ierr = PetscOptionsEnd();CHKERRQ(ierr);
2008   }
2009   if (flg) {
2010     ierr = MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);CHKERRQ(ierr);
2011   }
2012 
2013   C->ops->productsymbolic = MatProductSymbolic_AB;
2014   C->ops->matmultsymbolic = MatMatMultSymbolic_SeqAIJ_SeqAIJ;
2015   PetscFunctionReturn(0);
2016 }
2017 
2018 static PetscErrorCode MatProductSetFromOptions_SeqAIJ_AtB(Mat C)
2019 {
2020   PetscErrorCode ierr;
2021   Mat_Product    *product = C->product;
2022   PetscInt       alg = 0; /* default algorithm */
2023   PetscBool      flg = PETSC_FALSE;
2024   const char     *algTypes[3] = {"default","at*b","outerproduct"};
2025   PetscInt       nalg = 3;
2026 
2027   PetscFunctionBegin;
2028   /* Get runtime option */
2029   if (product->api_user) {
2030     ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatTransposeMatMult","Mat");CHKERRQ(ierr);
2031     ierr = PetscOptionsEList("-mattransposematmult_via","Algorithmic approach","MatTransposeMatMult",algTypes,nalg,algTypes[alg],&alg,&flg);CHKERRQ(ierr);
2032     ierr = PetscOptionsEnd();CHKERRQ(ierr);
2033   } else {
2034     ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatProduct_AtB","Mat");CHKERRQ(ierr);
2035     ierr = PetscOptionsEList("-matproduct_atb_via","Algorithmic approach","MatProduct_AtB",algTypes,nalg,algTypes[alg],&alg,&flg);CHKERRQ(ierr);
2036     ierr = PetscOptionsEnd();CHKERRQ(ierr);
2037   }
2038   if (flg) {
2039     ierr = MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);CHKERRQ(ierr);
2040   }
2041 
2042   C->ops->productsymbolic = MatProductSymbolic_AtB_SeqAIJ_SeqAIJ;
2043   PetscFunctionReturn(0);
2044 }
2045 
2046 static PetscErrorCode MatProductSetFromOptions_SeqAIJ_ABt(Mat C)
2047 {
2048   PetscErrorCode ierr;
2049   Mat_Product    *product = C->product;
2050   PetscInt       alg = 0; /* default algorithm */
2051   PetscBool      flg = PETSC_FALSE;
2052   const char     *algTypes[2] = {"default","color"};
2053   PetscInt       nalg = 2;
2054 
2055   PetscFunctionBegin;
2056   /* Set default algorithm */
2057   ierr = PetscStrcmp(C->product->alg,"default",&flg);CHKERRQ(ierr);
2058   if (!flg) {
2059     alg = 1;
2060     ierr = MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);CHKERRQ(ierr);
2061   }
2062 
2063   /* Get runtime option */
2064   if (product->api_user) {
2065     ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatMatTransposeMult","Mat");CHKERRQ(ierr);
2066     ierr = PetscOptionsEList("-matmattransmult_via","Algorithmic approach","MatMatTransposeMult",algTypes,nalg,algTypes[alg],&alg,&flg);CHKERRQ(ierr);
2067     ierr = PetscOptionsEnd();CHKERRQ(ierr);
2068   } else {
2069     ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatProduct_ABt","Mat");CHKERRQ(ierr);
2070     ierr = PetscOptionsEList("-matproduct_abt_via","Algorithmic approach","MatProduct_ABt",algTypes,nalg,algTypes[alg],&alg,&flg);CHKERRQ(ierr);
2071     ierr = PetscOptionsEnd();CHKERRQ(ierr);
2072   }
2073   if (flg) {
2074     ierr = MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);CHKERRQ(ierr);
2075   }
2076 
2077   C->ops->mattransposemultsymbolic = MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ;
2078   C->ops->productsymbolic          = MatProductSymbolic_ABt;
2079   PetscFunctionReturn(0);
2080 }
2081 
2082 static PetscErrorCode MatProductSetFromOptions_SeqAIJ_PtAP(Mat C)
2083 {
2084   PetscErrorCode ierr;
2085   Mat_Product    *product = C->product;
2086   PetscBool      flg = PETSC_FALSE;
2087   PetscInt       alg = 0; /* default algorithm -- alg=1 should be default!!! */
2088 #if !defined(PETSC_HAVE_HYPRE)
2089   const char      *algTypes[2] = {"scalable","rap"};
2090   PetscInt        nalg = 2;
2091 #else
2092   const char      *algTypes[3] = {"scalable","rap","hypre"};
2093   PetscInt        nalg = 3;
2094 #endif
2095 
2096   PetscFunctionBegin;
2097   /* Set default algorithm */
2098   ierr = PetscStrcmp(product->alg,"default",&flg);CHKERRQ(ierr);
2099   if (flg) {
2100     ierr = MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);CHKERRQ(ierr);
2101   }
2102 
2103   /* Get runtime option */
2104   if (product->api_user) {
2105     ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatPtAP","Mat");CHKERRQ(ierr);
2106     ierr = PetscOptionsEList("-matptap_via","Algorithmic approach","MatPtAP",algTypes,nalg,algTypes[0],&alg,&flg);CHKERRQ(ierr);
2107     ierr = PetscOptionsEnd();CHKERRQ(ierr);
2108   } else {
2109     ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatProduct_PtAP","Mat");CHKERRQ(ierr);
2110     ierr = PetscOptionsEList("-matproduct_ptap_via","Algorithmic approach","MatProduct_PtAP",algTypes,nalg,algTypes[0],&alg,&flg);CHKERRQ(ierr);
2111     ierr = PetscOptionsEnd();CHKERRQ(ierr);
2112   }
2113   if (flg) {
2114     ierr = MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);CHKERRQ(ierr);
2115   }
2116 
2117   C->ops->productsymbolic = MatProductSymbolic_PtAP_SeqAIJ_SeqAIJ;
2118   PetscFunctionReturn(0);
2119 }
2120 
2121 static PetscErrorCode MatProductSetFromOptions_SeqAIJ_RARt(Mat C)
2122 {
2123   PetscErrorCode ierr;
2124   Mat_Product    *product = C->product;
2125   PetscBool      flg = PETSC_FALSE;
2126   PetscInt       alg = 0; /* default algorithm */
2127   const char     *algTypes[3] = {"r*a*rt","r*art","coloring_rart"};
2128   PetscInt        nalg = 3;
2129 
2130   PetscFunctionBegin;
2131   /* Set default algorithm */
2132   ierr = PetscStrcmp(product->alg,"default",&flg);CHKERRQ(ierr);
2133   if (flg) {
2134     ierr = MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);CHKERRQ(ierr);
2135   }
2136 
2137   /* Get runtime option */
2138   if (product->api_user) {
2139     ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatRARt","Mat");CHKERRQ(ierr);
2140     ierr = PetscOptionsEList("-matrart_via","Algorithmic approach","MatRARt",algTypes,nalg,algTypes[0],&alg,&flg);CHKERRQ(ierr);
2141     ierr = PetscOptionsEnd();CHKERRQ(ierr);
2142   } else {
2143     ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatProduct_RARt","Mat");CHKERRQ(ierr);
2144     ierr = PetscOptionsEList("-matproduct_rart_via","Algorithmic approach","MatProduct_RARt",algTypes,nalg,algTypes[0],&alg,&flg);CHKERRQ(ierr);
2145     ierr = PetscOptionsEnd();CHKERRQ(ierr);
2146   }
2147   if (flg) {
2148     ierr = MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);CHKERRQ(ierr);
2149   }
2150 
2151   C->ops->productsymbolic = MatProductSymbolic_RARt_SeqAIJ_SeqAIJ;
2152   PetscFunctionReturn(0);
2153 }
2154 
2155 /* ABC = A*B*C = A*(B*C); ABC's algorithm must be chosen from AB's algorithm */
2156 static PetscErrorCode MatProductSetFromOptions_SeqAIJ_ABC(Mat C)
2157 {
2158   PetscErrorCode ierr;
2159   Mat_Product    *product = C->product;
2160   PetscInt       alg = 0; /* default algorithm */
2161   PetscBool      flg = PETSC_FALSE;
2162   const char     *algTypes[7] = {"sorted","scalable","scalable_fast","heap","btheap","llcondensed","rowmerge"};
2163   PetscInt       nalg = 7;
2164 
2165   PetscFunctionBegin;
2166   /* Set default algorithm */
2167   ierr = PetscStrcmp(product->alg,"default",&flg);CHKERRQ(ierr);
2168   if (flg) {
2169     ierr = MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);CHKERRQ(ierr);
2170   }
2171 
2172   /* Get runtime option */
2173   if (product->api_user) {
2174     ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatMatMatMult","Mat");CHKERRQ(ierr);
2175     ierr = PetscOptionsEList("-matmatmatmult_via","Algorithmic approach","MatMatMatMult",algTypes,nalg,algTypes[alg],&alg,&flg);CHKERRQ(ierr);
2176     ierr = PetscOptionsEnd();CHKERRQ(ierr);
2177   } else {
2178     ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatProduct_ABC","Mat");CHKERRQ(ierr);
2179     ierr = PetscOptionsEList("-matproduct_abc_via","Algorithmic approach","MatProduct_ABC",algTypes,nalg,algTypes[alg],&alg,&flg);CHKERRQ(ierr);
2180     ierr = PetscOptionsEnd();CHKERRQ(ierr);
2181   }
2182   if (flg) {
2183     ierr = MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);CHKERRQ(ierr);
2184   }
2185 
2186   C->ops->matmatmultsymbolic = MatMatMatMultSymbolic_SeqAIJ_SeqAIJ_SeqAIJ;
2187   C->ops->productsymbolic    = MatProductSymbolic_ABC;
2188   PetscFunctionReturn(0);
2189 }
2190 
2191 PetscErrorCode MatProductSetFromOptions_SeqAIJ(Mat C)
2192 {
2193   PetscErrorCode ierr;
2194   Mat_Product    *product = C->product;
2195 
2196   PetscFunctionBegin;
2197   switch (product->type) {
2198   case MATPRODUCT_AB:
2199     ierr = MatProductSetFromOptions_SeqAIJ_AB(C);CHKERRQ(ierr);
2200     break;
2201   case MATPRODUCT_AtB:
2202     ierr = MatProductSetFromOptions_SeqAIJ_AtB(C);CHKERRQ(ierr);
2203     break;
2204   case MATPRODUCT_ABt:
2205     ierr = MatProductSetFromOptions_SeqAIJ_ABt(C);CHKERRQ(ierr);
2206     break;
2207   case MATPRODUCT_PtAP:
2208     ierr = MatProductSetFromOptions_SeqAIJ_PtAP(C);CHKERRQ(ierr);
2209     break;
2210   case MATPRODUCT_RARt:
2211     ierr = MatProductSetFromOptions_SeqAIJ_RARt(C);CHKERRQ(ierr);
2212     break;
2213   case MATPRODUCT_ABC:
2214     ierr = MatProductSetFromOptions_SeqAIJ_ABC(C);CHKERRQ(ierr);
2215     break;
2216   default:
2217     break;
2218   }
2219   PetscFunctionReturn(0);
2220 }
2221