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