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