xref: /petsc/src/mat/impls/aij/seq/matmatmult.c (revision 0790b1abffaab3bdb1fe3a1da018b28197f55b5e)
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 <../src/mat/impls/dense/seq/dense.h> /*I "petscmat.h" I*/
11 
12 EXTERN_C_BEGIN
13 #undef __FUNCT__
14 #define __FUNCT__ "MatMatMult_SeqAIJ_SeqAIJ"
15 PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
16 {
17   PetscErrorCode ierr;
18 
19   PetscFunctionBegin;
20   if (scall == MAT_INITIAL_MATRIX){
21     /* ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); */
22     ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
23     /* ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);   */
24   }
25   /* ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); */
26   ierr = (*(*C)->ops->matmultnumeric)(A,B,*C);CHKERRQ(ierr);
27   /* ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); */
28   PetscFunctionReturn(0);
29 }
30 EXTERN_C_END
31 
32 #undef __FUNCT__
33 #define __FUNCT__ "MatGetSymbolicMatMatMult_SeqAIJ_SeqAIJ"
34 PetscErrorCode MatGetSymbolicMatMatMult_SeqAIJ_SeqAIJ(PetscInt am,PetscInt *Ai,PetscInt *Aj,PetscInt bm,PetscInt bn,PetscInt *Bi,PetscInt *Bj,PetscReal fill,PetscInt *Ci[],PetscInt *Cj[],PetscInt *nspacedouble)
35 {
36   PetscErrorCode     ierr;
37   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
38   PetscInt           *ai=Ai,*aj=Aj,*bi=Bi,*bj=Bj,*bjj,*ci,*cj;
39   PetscInt           i,j,anzi,brow,bnzj,cnzi,nlnk,*lnk,ndouble=0;
40   PetscBT            lnkbt;
41 
42   PetscFunctionBegin;
43   /* Allocate ci array, arrays for fill computation and */
44   /* free space for accumulating nonzero column info */
45   ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
46   ci[0] = 0;
47 
48   /* create and initialize a linked list */
49   nlnk = bn+1;
50   ierr = PetscLLCreate(bn,bn,nlnk,lnk,lnkbt);CHKERRQ(ierr);
51 
52   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
53   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr);
54   current_space = free_space;
55 
56   /* Determine symbolic info for each row of the product: */
57   for (i=0;i<am;i++) {
58     anzi = ai[i+1] - ai[i];
59     cnzi = 0;
60     j    = anzi;
61     aj   = Aj + ai[i];
62     while (j){/* assume cols are almost in increasing order, starting from its end saves computation */
63       j--;
64       brow = aj[j];
65       bnzj = bi[brow+1] - bi[brow];
66       bjj  = bj + bi[brow];
67       /* add non-zero cols of B into the sorted linked list lnk */
68       ierr = PetscLLAdd(bnzj,bjj,bn,nlnk,lnk,lnkbt);CHKERRQ(ierr);
69       cnzi += nlnk;
70     }
71 
72     /* If free space is not available, make more free space */
73     /* Double the amount of total space in the list */
74     if (current_space->local_remaining<cnzi) {
75       ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,&current_space);CHKERRQ(ierr);
76       ndouble++;
77     }
78 
79     /* Copy data into free space, then initialize lnk */
80     ierr = PetscLLClean(bn,bn,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
81     current_space->array           += cnzi;
82     current_space->local_used      += cnzi;
83     current_space->local_remaining -= cnzi;
84     ci[i+1] = ci[i] + cnzi;
85   }
86 
87   /* Column indices are in the list of free space */
88   /* Allocate space for cj, initialize cj, and */
89   /* destroy list of free space and other temporary array(s) */
90   ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
91   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
92   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
93 
94   *Ci           = ci;
95   *Cj           = cj;
96   *nspacedouble = ndouble;
97   PetscFunctionReturn(0);
98 }
99 
100 #undef __FUNCT__
101 #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ"
102 PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
103 {
104   PetscErrorCode ierr;
105   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
106   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*ci,*cj;
107   PetscInt       am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N,nspacedouble;
108   MatScalar      *ca;
109   PetscReal      afill;
110   PetscBool      dense_axpy; /* false: use sparse axpy; otherwise use dense axpy in MatMatMultNumeric_SeqAIJ_SeqAIJ() */
111 
112   PetscFunctionBegin;
113   /* Get ci and cj */
114   ierr = MatGetSymbolicMatMatMult_SeqAIJ_SeqAIJ(am,ai,aj,bm,bn,bi,bj,fill,&ci,&cj,&nspacedouble);CHKERRQ(ierr);
115 
116   /* Allocate space for ca */
117   ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
118   ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr);
119 
120   /* put together the new symbolic matrix */
121   ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bn,ci,cj,ca,C);CHKERRQ(ierr);
122 
123   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
124   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
125   c = (Mat_SeqAIJ *)((*C)->data);
126   c->free_a   = PETSC_TRUE;
127   c->free_ij  = PETSC_TRUE;
128   c->nonew    = 0;
129   (*C)->ops->matmult = MatMatMult_SeqAIJ_SeqAIJ;
130 
131   /* Determine which MatMatMultNumeric_SeqAIJ_SeqAIJ() to be used */
132   dense_axpy = PETSC_TRUE;
133   ierr = PetscOptionsGetBool(PETSC_NULL,"-matmatmult_denseaxpy",&dense_axpy,PETSC_NULL);CHKERRQ(ierr);
134   if (dense_axpy){
135     ierr = PetscMalloc(bn*sizeof(PetscScalar),&c->matmult_abdense);CHKERRQ(ierr);
136     ierr = PetscMemzero(c->matmult_abdense,dense_axpy*bn*sizeof(PetscScalar));CHKERRQ(ierr);
137     (*C)->ops->matmultnumeric =  MatMatMultNumeric_SeqAIJ_SeqAIJ; /* fast, takes additional dense_axpy*bn*sizeof(PetscScalar) space */
138   } else { /* slower, but use less memory */
139     (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_SparseAxpy; /* slower, less memory */
140   }
141 
142   /* set MatInfo */
143   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
144   if (afill < 1.0) afill = 1.0;
145   c->maxnz                     = ci[am];
146   c->nz                        = ci[am];
147   (*C)->info.mallocs           = nspacedouble;
148   (*C)->info.fill_ratio_given  = fill;
149   (*C)->info.fill_ratio_needed = afill;
150 
151 #if defined(PETSC_USE_INFO)
152   if (ci[am]) {
153     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
154     ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
155   } else {
156     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
157   }
158 #endif
159   PetscFunctionReturn(0);
160 }
161 
162 #undef __FUNCT__
163 #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqAIJ"
164 PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
165 {
166   PetscErrorCode ierr;
167   PetscLogDouble flops=0.0;
168   Mat_SeqAIJ     *a = (Mat_SeqAIJ *)A->data;
169   Mat_SeqAIJ     *b = (Mat_SeqAIJ *)B->data;
170   Mat_SeqAIJ     *c = (Mat_SeqAIJ *)C->data;
171   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j;
172   PetscInt       am=A->rmap->n,cm=C->rmap->n;
173   PetscInt       i,j,k,anzi,bnzi,cnzi,brow;
174   PetscScalar    *aa=a->a,*ba=b->a,*baj,*ca=c->a,valtmp;
175   PetscScalar    *ab_dense=c->matmult_abdense;
176 
177   PetscFunctionBegin;
178   /* clean old values in C */
179   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
180   /* Traverse A row-wise. */
181   /* Build the ith row in C by summing over nonzero columns in A, */
182   /* the rows of B corresponding to nonzeros of A. */
183   for (i=0; i<am; i++) {
184     anzi = ai[i+1] - ai[i];
185     for (j=0; j<anzi; j++) {
186       brow = aj[j];
187       bnzi = bi[brow+1] - bi[brow];
188       bjj  = bj + bi[brow];
189       baj  = ba + bi[brow];
190       /* perform dense axpy */
191       valtmp = aa[j];
192       for (k=0; k<bnzi; k++) {
193         ab_dense[bjj[k]] += valtmp*baj[k];
194       }
195       flops += 2*bnzi;
196     }
197     aj += anzi; aa += anzi;
198 
199     cnzi = ci[i+1] - ci[i];
200     for (k=0; k<cnzi; k++) {
201       ca[k]          += ab_dense[cj[k]];
202       ab_dense[cj[k]] = 0.0; /* zero ab_dense */
203     }
204     flops += cnzi;
205     cj += cnzi; ca += cnzi;
206   }
207   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
208   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
209   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
210   C->ops->matmultnumeric =  MatMatMultNumeric_SeqAIJ_SeqAIJ;
211   PetscFunctionReturn(0);
212 }
213 
214 #undef __FUNCT__
215 #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqAIJ_SparseAxpy"
216 PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_SparseAxpy(Mat A,Mat B,Mat C)
217 {
218   PetscErrorCode ierr;
219   PetscLogDouble flops=0.0;
220   Mat_SeqAIJ     *a = (Mat_SeqAIJ *)A->data;
221   Mat_SeqAIJ     *b = (Mat_SeqAIJ *)B->data;
222   Mat_SeqAIJ     *c = (Mat_SeqAIJ *)C->data;
223   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j;
224   PetscInt       am=A->rmap->N,cm=C->rmap->N;
225   PetscInt       i,j,k,anzi,bnzi,cnzi,brow;
226   PetscScalar    *aa=a->a,*ba=b->a,*baj,*ca=c->a,valtmp;
227   PetscInt       nextb;
228 
229   PetscFunctionBegin;
230   /* clean old values in C */
231   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
232   /* Traverse A row-wise. */
233   /* Build the ith row in C by summing over nonzero columns in A, */
234   /* the rows of B corresponding to nonzeros of A. */
235   for (i=0;i<am;i++) {
236     anzi = ai[i+1] - ai[i];
237     cnzi = ci[i+1] - ci[i];
238     for (j=0;j<anzi;j++) {
239       brow = aj[j];
240       bnzi = bi[brow+1] - bi[brow];
241       bjj  = bj + bi[brow];
242       baj  = ba + bi[brow];
243       /* perform sparse axpy */
244       valtmp = aa[j];
245       nextb  = 0;
246       for (k=0; nextb<bnzi; k++) {
247         if (cj[k] == bjj[nextb]){ /* ccol == bcol */
248           ca[k] += valtmp*baj[nextb++];
249         }
250       }
251       flops += 2*bnzi;
252     }
253     aj += anzi; aa += anzi;
254     cj += cnzi; ca += cnzi;
255   }
256 
257   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
258   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
259   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
260   PetscFunctionReturn(0);
261 }
262 
263 /* This routine is not used. Should be removed! */
264 #undef __FUNCT__
265 #define __FUNCT__ "MatMatTransposeMult_SeqAIJ_SeqAIJ"
266 PetscErrorCode MatMatTransposeMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
267 {
268   PetscErrorCode ierr;
269 
270   PetscFunctionBegin;
271   if (scall == MAT_INITIAL_MATRIX){
272     ierr = MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
273   }
274   ierr = MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr);
275   PetscFunctionReturn(0);
276 }
277 
278 #undef __FUNCT__
279 #define __FUNCT__ "PetscContainerDestroy_Mat_MatMatTransMult"
280 PetscErrorCode PetscContainerDestroy_Mat_MatMatTransMult(void *ptr)
281 {
282   PetscErrorCode      ierr;
283   Mat_MatMatTransMult *multtrans=(Mat_MatMatTransMult*)ptr;
284 
285   PetscFunctionBegin;
286   ierr = MatTransposeColoringDestroy(&multtrans->matcoloring);CHKERRQ(ierr);
287   ierr = MatDestroy(&multtrans->Bt_den);CHKERRQ(ierr);
288   ierr = MatDestroy(&multtrans->ABt_den);CHKERRQ(ierr);
289   ierr = PetscFree(multtrans);CHKERRQ(ierr);
290   PetscFunctionReturn(0);
291 }
292 
293 #undef __FUNCT__
294 #define __FUNCT__ "MatDestroy_SeqAIJ_MatMatMultTrans"
295 PetscErrorCode MatDestroy_SeqAIJ_MatMatMultTrans(Mat A)
296 {
297   PetscErrorCode      ierr;
298   PetscContainer      container;
299   Mat_MatMatTransMult *multtrans=PETSC_NULL;
300 
301   PetscFunctionBegin;
302   ierr = PetscObjectQuery((PetscObject)A,"Mat_MatMatTransMult",(PetscObject *)&container);CHKERRQ(ierr);
303   if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit");
304   ierr = PetscContainerGetPointer(container,(void **)&multtrans);CHKERRQ(ierr);
305   A->ops->destroy   = multtrans->destroy;
306   if (A->ops->destroy) {
307     ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
308   }
309   ierr = PetscObjectCompose((PetscObject)A,"Mat_MatMatTransMult",0);CHKERRQ(ierr);
310   PetscFunctionReturn(0);
311 }
312 
313 #undef __FUNCT__
314 #define __FUNCT__ "MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ"
315 PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
316 {
317   PetscErrorCode      ierr;
318   Mat                 Bt;
319   PetscInt            *bti,*btj;
320   Mat_MatMatTransMult *multtrans;
321   PetscContainer      container;
322   PetscLogDouble      t0,tf,etime2=0.0;
323 
324   PetscFunctionBegin;
325   ierr = PetscGetTime(&t0);CHKERRQ(ierr);
326    /* create symbolic Bt */
327   ierr = MatGetSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr);
328   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,B->cmap->n,B->rmap->n,bti,btj,PETSC_NULL,&Bt);CHKERRQ(ierr);
329 
330   /* get symbolic C=A*Bt */
331   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,Bt,fill,C);CHKERRQ(ierr);
332 
333   /* create a supporting struct for reuse intermidiate dense matrices with matcoloring */
334   ierr = PetscNew(Mat_MatMatTransMult,&multtrans);CHKERRQ(ierr);
335 
336   /* attach the supporting struct to C */
337   ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr);
338   ierr = PetscContainerSetPointer(container,multtrans);CHKERRQ(ierr);
339   ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_Mat_MatMatTransMult);CHKERRQ(ierr);
340   ierr = PetscObjectCompose((PetscObject)(*C),"Mat_MatMatTransMult",(PetscObject)container);CHKERRQ(ierr);
341   ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
342 
343   multtrans->usecoloring = PETSC_FALSE;
344   multtrans->destroy = (*C)->ops->destroy;
345   (*C)->ops->destroy = MatDestroy_SeqAIJ_MatMatMultTrans;
346 
347   ierr = PetscGetTime(&tf);CHKERRQ(ierr);
348   etime2 += tf - t0;
349 
350   ierr = PetscOptionsGetBool(PETSC_NULL,"-matmattransmult_color",&multtrans->usecoloring,PETSC_NULL);CHKERRQ(ierr);
351   if (multtrans->usecoloring){
352     /* Create MatTransposeColoring from symbolic C=A*B^T */
353     MatTransposeColoring matcoloring;
354     ISColoring           iscoloring;
355     Mat                  Bt_dense,C_dense;
356     PetscLogDouble       etime0=0.0,etime01=0.0,etime1=0.0;
357 
358     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
359     ierr = MatGetColoring(*C,MATCOLORINGLF,&iscoloring);CHKERRQ(ierr);
360     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
361     etime0 += tf - t0;
362 
363     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
364     ierr = MatTransposeColoringCreate(*C,iscoloring,&matcoloring);CHKERRQ(ierr);
365     multtrans->matcoloring = matcoloring;
366     ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr);
367     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
368     etime01 += tf - t0;
369 
370     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
371     /* Create Bt_dense and C_dense = A*Bt_dense */
372     ierr = MatCreate(PETSC_COMM_SELF,&Bt_dense);CHKERRQ(ierr);
373     ierr = MatSetSizes(Bt_dense,A->cmap->n,matcoloring->ncolors,A->cmap->n,matcoloring->ncolors);CHKERRQ(ierr);
374     ierr = MatSetType(Bt_dense,MATSEQDENSE);CHKERRQ(ierr);
375     ierr = MatSeqDenseSetPreallocation(Bt_dense,PETSC_NULL);CHKERRQ(ierr);
376     Bt_dense->assembled = PETSC_TRUE;
377     multtrans->Bt_den = Bt_dense;
378 
379     ierr = MatCreate(PETSC_COMM_SELF,&C_dense);CHKERRQ(ierr);
380     ierr = MatSetSizes(C_dense,A->rmap->n,matcoloring->ncolors,A->rmap->n,matcoloring->ncolors);CHKERRQ(ierr);
381     ierr = MatSetType(C_dense,MATSEQDENSE);CHKERRQ(ierr);
382     ierr = MatSeqDenseSetPreallocation(C_dense,PETSC_NULL);CHKERRQ(ierr);
383     Bt_dense->assembled = PETSC_TRUE;
384     multtrans->ABt_den = C_dense;
385     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
386     etime1 += tf - t0;
387 
388 #if defined(PETSC_USE_INFO)
389     {
390     Mat_SeqAIJ *c=(Mat_SeqAIJ*)(*C)->data;
391     ierr = PetscInfo5(*C,"Bt_dense: %D,%D; Cnz %D / (cm*ncolors %D) = %g\n",A->cmap->n,matcoloring->ncolors,c->nz,A->rmap->n*matcoloring->ncolors,(PetscReal)(c->nz)/(A->rmap->n*matcoloring->ncolors));
392     ierr = PetscInfo5(*C,"Sym = GetColor %g + ColorCreate %g + MatDenseCreate %g + non-colorSym %g = %g\n",etime0,etime01,etime1,etime2,etime0+etime01+etime1+etime2);
393     }
394 #endif
395   }
396   /* clean up */
397   ierr = MatDestroy(&Bt);CHKERRQ(ierr);
398   ierr = MatRestoreSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr);
399 
400 
401 
402 #if defined(INEFFICIENT_ALGORITHM)
403   /* The algorithm below computes am*bm sparse inner-product - inefficient! It will be deleted later. */
404   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
405   Mat_SeqAIJ         *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
406   PetscInt           *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*ci,*cj,*acol,*bcol;
407   PetscInt           am=A->rmap->N,bm=B->rmap->N;
408   PetscInt           i,j,anzi,bnzj,cnzi,nlnk,*lnk,nspacedouble=0,ka,kb,index[1];
409   MatScalar          *ca;
410   PetscBT            lnkbt;
411   PetscReal          afill;
412 
413   /* Allocate row pointer array ci  */
414   ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
415   ci[0] = 0;
416 
417   /* Create and initialize a linked list for C columns */
418   nlnk = bm+1;
419   ierr = PetscLLCreate(bm,bm,nlnk,lnk,lnkbt);CHKERRQ(ierr);
420 
421   /* Initial FreeSpace with size fill*(nnz(A)+nnz(B)) */
422   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr);
423   current_space = free_space;
424 
425   /* Determine symbolic info for each row of the product A*B^T: */
426   for (i=0; i<am; i++) {
427     anzi = ai[i+1] - ai[i];
428     cnzi = 0;
429     acol = aj + ai[i];
430     for (j=0; j<bm; j++){
431       bnzj = bi[j+1] - bi[j];
432       bcol= bj + bi[j];
433       /* sparse inner-product c(i,j)=A[i,:]*B[j,:]^T */
434       ka = 0; kb = 0;
435       while (ka < anzi && kb < bnzj){
436         while (acol[ka] < bcol[kb] && ka < anzi) ka++;
437         if (ka == anzi) break;
438         while (acol[ka] > bcol[kb] && kb < bnzj) kb++;
439         if (kb == bnzj) break;
440         if (acol[ka] == bcol[kb]){ /* add nonzero c(i,j) to lnk */
441           index[0] = j;
442           ierr = PetscLLAdd(1,index,bm,nlnk,lnk,lnkbt);CHKERRQ(ierr);
443           cnzi++;
444           break;
445         }
446       }
447     }
448 
449     /* If free space is not available, make more free space */
450     /* Double the amount of total space in the list */
451     if (current_space->local_remaining<cnzi) {
452       ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,&current_space);CHKERRQ(ierr);
453       nspacedouble++;
454     }
455 
456     /* Copy data into free space, then initialize lnk */
457     ierr = PetscLLClean(bm,bm,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
458     current_space->array           += cnzi;
459     current_space->local_used      += cnzi;
460     current_space->local_remaining -= cnzi;
461 
462     ci[i+1] = ci[i] + cnzi;
463   }
464 
465 
466   /* Column indices are in the list of free space.
467      Allocate array cj, copy column indices to cj, and destroy list of free space */
468   ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
469   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
470   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
471 
472   /* Allocate space for ca */
473   ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
474   ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr);
475 
476   /* put together the new symbolic matrix */
477   ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bm,ci,cj,ca,C);CHKERRQ(ierr);
478 
479   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
480   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
481   c = (Mat_SeqAIJ *)((*C)->data);
482   c->free_a   = PETSC_TRUE;
483   c->free_ij  = PETSC_TRUE;
484   c->nonew    = 0;
485 
486   /* set MatInfo */
487   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
488   if (afill < 1.0) afill = 1.0;
489   c->maxnz                     = ci[am];
490   c->nz                        = ci[am];
491   (*C)->info.mallocs           = nspacedouble;
492   (*C)->info.fill_ratio_given  = fill;
493   (*C)->info.fill_ratio_needed = afill;
494 
495 #if defined(PETSC_USE_INFO)
496   if (ci[am]) {
497     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
498     ierr = PetscInfo1((*C),"Use MatMatTransposeMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr);
499   } else {
500     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
501   }
502 #endif
503 #endif
504   PetscFunctionReturn(0);
505 }
506 
507 /* #define USE_ARRAY - for sparse dot product. Slower than !USE_ARRAY */
508 #undef __FUNCT__
509 #define __FUNCT__ "MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ"
510 PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
511 {
512   PetscErrorCode ierr;
513   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data;
514   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,anzi,bnzj,nexta,nextb,*acol,*bcol,brow;
515   PetscInt       cm=C->rmap->n,*ci=c->i,*cj=c->j,i,j,cnzi,*ccol;
516   PetscLogDouble flops=0.0;
517   MatScalar      *aa=a->a,*aval,*ba=b->a,*bval,*ca=c->a,*cval;
518   Mat_MatMatTransMult *multtrans;
519   PetscContainer      container;
520 #if defined(USE_ARRAY)
521   MatScalar      *spdot;
522 #endif
523 
524   PetscFunctionBegin;
525   ierr = PetscObjectQuery((PetscObject)C,"Mat_MatMatTransMult",(PetscObject *)&container);CHKERRQ(ierr);
526   if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit");
527   ierr  = PetscContainerGetPointer(container,(void **)&multtrans);CHKERRQ(ierr);
528   if (multtrans->usecoloring){
529     MatTransposeColoring  matcoloring = multtrans->matcoloring;
530     Mat                   Bt_dense;
531     PetscInt              m,n;
532     PetscLogDouble t0,tf,etime0=0.0,etime1=0.0,etime2=0.0;
533     Mat C_dense = multtrans->ABt_den;
534 
535     Bt_dense = multtrans->Bt_den;
536     ierr = MatGetLocalSize(Bt_dense,&m,&n);CHKERRQ(ierr);
537 
538     /* Get Bt_dense by Apply MatTransposeColoring to B */
539     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
540     ierr = MatTransColoringApplySpToDen(matcoloring,B,Bt_dense);CHKERRQ(ierr);
541     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
542     etime0 += tf - t0;
543 
544     /* C_dense = A*Bt_dense */
545     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
546     ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,Bt_dense,C_dense);CHKERRQ(ierr);
547     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
548     etime2 += tf - t0;
549 
550     /* Recover C from C_dense */
551     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
552     ierr = MatTransColoringApplyDenToSp(matcoloring,C_dense,C);CHKERRQ(ierr);
553     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
554     etime1 += tf - t0;
555 #if defined(PETSC_USE_INFO)
556     ierr = PetscInfo4(C,"Num = ColoringApply: %g %g + Mult_sp_dense %g = %g\n",etime0,etime1,etime2,etime0+etime1+etime2);
557 #endif
558     PetscFunctionReturn(0);
559   }
560 
561 #if defined(USE_ARRAY)
562   /* allocate an array for implementing sparse inner-product */
563   ierr = PetscMalloc((A->cmap->n+1)*sizeof(MatScalar),&spdot);CHKERRQ(ierr);
564   ierr = PetscMemzero(spdot,(A->cmap->n+1)*sizeof(MatScalar));CHKERRQ(ierr);
565 #endif
566 
567   /* clear old values in C */
568   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
569 
570   for (i=0; i<cm; i++) {
571     anzi = ai[i+1] - ai[i];
572     acol = aj + ai[i];
573     aval = aa + ai[i];
574     cnzi = ci[i+1] - ci[i];
575     ccol = cj + ci[i];
576     cval = ca + ci[i];
577     for (j=0; j<cnzi; j++){
578       brow = ccol[j];
579       bnzj = bi[brow+1] - bi[brow];
580       bcol = bj + bi[brow];
581       bval = ba + bi[brow];
582 
583       /* perform sparse inner-product c(i,j)=A[i,:]*B[j,:]^T */
584 #if defined(USE_ARRAY)
585       /* put ba to spdot array */
586       for (nextb=0; nextb<bnzj; nextb++) spdot[bcol[nextb]] = bval[nextb];
587       /* c(i,j)=A[i,:]*B[j,:]^T */
588       for (nexta=0; nexta<anzi; nexta++){
589         cval[j] += spdot[acol[nexta]]*aval[nexta];
590       }
591       /* zero spdot array */
592       for (nextb=0; nextb<bnzj; nextb++) spdot[bcol[nextb]] = 0.0;
593 #else
594       nexta = 0; nextb = 0;
595       while (nexta<anzi && nextb<bnzj){
596         while (acol[nexta] < bcol[nextb] && nexta < anzi) nexta++;
597         if (nexta == anzi) break;
598         while (acol[nexta] > bcol[nextb] && nextb < bnzj) nextb++;
599         if (nextb == bnzj) break;
600         if (acol[nexta] == bcol[nextb]){
601           cval[j] += aval[nexta]*bval[nextb];
602           nexta++; nextb++;
603           flops += 2;
604         }
605       }
606 #endif
607     }
608   }
609   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
610   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
611   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
612 #if defined(USE_ARRAY)
613   ierr = PetscFree(spdot);
614 #endif
615   PetscFunctionReturn(0);
616 }
617 
618 #undef __FUNCT__
619 #define __FUNCT__ "MatTransposeMatMult_SeqAIJ_SeqAIJ"
620 PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) {
621   PetscErrorCode ierr;
622 
623   PetscFunctionBegin;
624   if (scall == MAT_INITIAL_MATRIX){
625     ierr = MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr);
626   }
627   ierr = MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr);
628   PetscFunctionReturn(0);
629 }
630 
631 #undef __FUNCT__
632 #define __FUNCT__ "MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ"
633 PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
634 {
635   PetscErrorCode ierr;
636   Mat            At;
637   PetscInt       *ati,*atj;
638 
639   PetscFunctionBegin;
640   /* create symbolic At */
641   ierr = MatGetSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr);
642   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,A->cmap->n,A->rmap->n,ati,atj,PETSC_NULL,&At);CHKERRQ(ierr);
643 
644   /* get symbolic C=At*B */
645   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(At,B,fill,C);CHKERRQ(ierr);
646 
647   /* clean up */
648   ierr = MatDestroy(&At);CHKERRQ(ierr);
649   ierr = MatRestoreSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr);
650   PetscFunctionReturn(0);
651 }
652 
653 #undef __FUNCT__
654 #define __FUNCT__ "MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ"
655 PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
656 {
657   PetscErrorCode ierr;
658   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data;
659   PetscInt       am=A->rmap->n,anzi,*ai=a->i,*aj=a->j,*bi=b->i,*bj,bnzi,nextb;
660   PetscInt       cm=C->rmap->n,*ci=c->i,*cj=c->j,crow,*cjj,i,j,k;
661   PetscLogDouble flops=0.0;
662   MatScalar      *aa=a->a,*ba,*ca=c->a,*caj;
663 
664   PetscFunctionBegin;
665   /* clear old values in C */
666   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
667 
668   /* compute A^T*B using outer product (A^T)[:,i]*B[i,:] */
669   for (i=0;i<am;i++) {
670     bj   = b->j + bi[i];
671     ba   = b->a + bi[i];
672     bnzi = bi[i+1] - bi[i];
673     anzi = ai[i+1] - ai[i];
674     for (j=0; j<anzi; j++) {
675       nextb = 0;
676       crow  = *aj++;
677       cjj   = cj + ci[crow];
678       caj   = ca + ci[crow];
679       /* perform sparse axpy operation.  Note cjj includes bj. */
680       for (k=0; nextb<bnzi; k++) {
681         if (cjj[k] == *(bj+nextb)) { /* ccol == bcol */
682           caj[k] += (*aa)*(*(ba+nextb));
683           nextb++;
684         }
685       }
686       flops += 2*bnzi;
687       aa++;
688     }
689   }
690 
691   /* Assemble the final matrix and clean up */
692   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
693   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
694   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
695   PetscFunctionReturn(0);
696 }
697 
698 EXTERN_C_BEGIN
699 #undef __FUNCT__
700 #define __FUNCT__ "MatMatMult_SeqAIJ_SeqDense"
701 PetscErrorCode MatMatMult_SeqAIJ_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
702 {
703   PetscErrorCode ierr;
704 
705   PetscFunctionBegin;
706   if (scall == MAT_INITIAL_MATRIX){
707     ierr = MatMatMultSymbolic_SeqAIJ_SeqDense(A,B,fill,C);CHKERRQ(ierr);
708   }
709   ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,B,*C);CHKERRQ(ierr);
710   PetscFunctionReturn(0);
711 }
712 EXTERN_C_END
713 
714 #undef __FUNCT__
715 #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqDense"
716 PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
717 {
718   PetscErrorCode ierr;
719 
720   PetscFunctionBegin;
721   ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,0.0,C);CHKERRQ(ierr);
722   (*C)->ops->matmult = MatMatMult_SeqAIJ_SeqDense;
723   PetscFunctionReturn(0);
724 }
725 
726 #undef __FUNCT__
727 #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqDense"
728 PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
729 {
730   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
731   PetscErrorCode ierr;
732   PetscScalar    *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4;
733   MatScalar      *aa;
734   PetscInt       cm=C->rmap->n, cn=B->cmap->n, bm=B->rmap->n, col, i,j,n,*aj, am = A->rmap->n;
735   PetscInt       am2 = 2*am, am3 = 3*am,  bm4 = 4*bm,colam;
736 
737   PetscFunctionBegin;
738   if (!cm || !cn) PetscFunctionReturn(0);
739   if (bm != A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Number columns in A %D not equal rows in B %D\n",A->cmap->n,bm);
740   if (A->rmap->n != C->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Number rows in C %D not equal rows in A %D\n",C->rmap->n,A->rmap->n);
741   if (B->cmap->n != C->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Number columns in B %D not equal columns in C %D\n",B->cmap->n,C->cmap->n);
742   ierr = MatGetArray(B,&b);CHKERRQ(ierr);
743   ierr = MatGetArray(C,&c);CHKERRQ(ierr);
744   b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;
745   for (col=0; col<cn-4; col += 4){  /* over columns of C */
746     colam = col*am;
747     for (i=0; i<am; i++) {        /* over rows of C in those columns */
748       r1 = r2 = r3 = r4 = 0.0;
749       n   = a->i[i+1] - a->i[i];
750       aj  = a->j + a->i[i];
751       aa  = a->a + a->i[i];
752       for (j=0; j<n; j++) {
753         r1 += (*aa)*b1[*aj];
754         r2 += (*aa)*b2[*aj];
755         r3 += (*aa)*b3[*aj];
756         r4 += (*aa++)*b4[*aj++];
757       }
758       c[colam + i]       = r1;
759       c[colam + am + i]  = r2;
760       c[colam + am2 + i] = r3;
761       c[colam + am3 + i] = r4;
762     }
763     b1 += bm4;
764     b2 += bm4;
765     b3 += bm4;
766     b4 += bm4;
767   }
768   for (;col<cn; col++){     /* over extra columns of C */
769     for (i=0; i<am; i++) {  /* over rows of C in those columns */
770       r1 = 0.0;
771       n   = a->i[i+1] - a->i[i];
772       aj  = a->j + a->i[i];
773       aa  = a->a + a->i[i];
774 
775       for (j=0; j<n; j++) {
776         r1 += (*aa++)*b1[*aj++];
777       }
778       c[col*am + i]     = r1;
779     }
780     b1 += bm;
781   }
782   ierr = PetscLogFlops(cn*(2.0*a->nz));CHKERRQ(ierr);
783   ierr = MatRestoreArray(B,&b);CHKERRQ(ierr);
784   ierr = MatRestoreArray(C,&c);CHKERRQ(ierr);
785   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
786   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
787   PetscFunctionReturn(0);
788 }
789 
790 /*
791    Note very similar to MatMult_SeqAIJ(), should generate both codes from same base
792 */
793 #undef __FUNCT__
794 #define __FUNCT__ "MatMatMultNumericAdd_SeqAIJ_SeqDense"
795 PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
796 {
797   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
798   PetscErrorCode ierr;
799   PetscScalar    *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4;
800   MatScalar      *aa;
801   PetscInt       cm=C->rmap->n, cn=B->cmap->n, bm=B->rmap->n, col, i,j,n,*aj, am = A->rmap->n,*ii,arm;
802   PetscInt       am2 = 2*am, am3 = 3*am,  bm4 = 4*bm,colam,*ridx;
803 
804   PetscFunctionBegin;
805   if (!cm || !cn) PetscFunctionReturn(0);
806   ierr = MatGetArray(B,&b);CHKERRQ(ierr);
807   ierr = MatGetArray(C,&c);CHKERRQ(ierr);
808   b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;
809 
810   if (a->compressedrow.use){ /* use compressed row format */
811     for (col=0; col<cn-4; col += 4){  /* over columns of C */
812       colam = col*am;
813       arm   = a->compressedrow.nrows;
814       ii    = a->compressedrow.i;
815       ridx  = a->compressedrow.rindex;
816       for (i=0; i<arm; i++) {        /* over rows of C in those columns */
817 	r1 = r2 = r3 = r4 = 0.0;
818 	n   = ii[i+1] - ii[i];
819 	aj  = a->j + ii[i];
820 	aa  = a->a + ii[i];
821 	for (j=0; j<n; j++) {
822 	  r1 += (*aa)*b1[*aj];
823 	  r2 += (*aa)*b2[*aj];
824 	  r3 += (*aa)*b3[*aj];
825 	  r4 += (*aa++)*b4[*aj++];
826 	}
827 	c[colam       + ridx[i]] += r1;
828 	c[colam + am  + ridx[i]] += r2;
829 	c[colam + am2 + ridx[i]] += r3;
830 	c[colam + am3 + ridx[i]] += r4;
831       }
832       b1 += bm4;
833       b2 += bm4;
834       b3 += bm4;
835       b4 += bm4;
836     }
837     for (;col<cn; col++){     /* over extra columns of C */
838       colam = col*am;
839       arm   = a->compressedrow.nrows;
840       ii    = a->compressedrow.i;
841       ridx  = a->compressedrow.rindex;
842       for (i=0; i<arm; i++) {  /* over rows of C in those columns */
843 	r1 = 0.0;
844 	n   = ii[i+1] - ii[i];
845 	aj  = a->j + ii[i];
846 	aa  = a->a + ii[i];
847 
848 	for (j=0; j<n; j++) {
849 	  r1 += (*aa++)*b1[*aj++];
850 	}
851 	c[col*am + ridx[i]] += r1;
852       }
853       b1 += bm;
854     }
855   } else {
856     for (col=0; col<cn-4; col += 4){  /* over columns of C */
857       colam = col*am;
858       for (i=0; i<am; i++) {        /* over rows of C in those columns */
859 	r1 = r2 = r3 = r4 = 0.0;
860 	n   = a->i[i+1] - a->i[i];
861 	aj  = a->j + a->i[i];
862 	aa  = a->a + a->i[i];
863 	for (j=0; j<n; j++) {
864 	  r1 += (*aa)*b1[*aj];
865 	  r2 += (*aa)*b2[*aj];
866 	  r3 += (*aa)*b3[*aj];
867 	  r4 += (*aa++)*b4[*aj++];
868 	}
869 	c[colam + i]       += r1;
870 	c[colam + am + i]  += r2;
871 	c[colam + am2 + i] += r3;
872 	c[colam + am3 + i] += r4;
873       }
874       b1 += bm4;
875       b2 += bm4;
876       b3 += bm4;
877       b4 += bm4;
878     }
879     for (;col<cn; col++){     /* over extra columns of C */
880       for (i=0; i<am; i++) {  /* over rows of C in those columns */
881 	r1 = 0.0;
882 	n   = a->i[i+1] - a->i[i];
883 	aj  = a->j + a->i[i];
884 	aa  = a->a + a->i[i];
885 
886 	for (j=0; j<n; j++) {
887 	  r1 += (*aa++)*b1[*aj++];
888 	}
889 	c[col*am + i]     += r1;
890       }
891       b1 += bm;
892     }
893   }
894   ierr = PetscLogFlops(cn*2.0*a->nz);CHKERRQ(ierr);
895   ierr = MatRestoreArray(B,&b);CHKERRQ(ierr);
896   ierr = MatRestoreArray(C,&c);CHKERRQ(ierr);
897   PetscFunctionReturn(0);
898 }
899 
900 #undef __FUNCT__
901 #define __FUNCT__ "MatTransColoringApplySpToDen_SeqAIJ"
902 PetscErrorCode  MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring coloring,Mat B,Mat Btdense)
903 {
904   PetscErrorCode ierr;
905   Mat_SeqAIJ     *b = (Mat_SeqAIJ*)B->data;
906   Mat_SeqDense   *btdense = (Mat_SeqDense*)Btdense->data;
907   PetscInt       *bi=b->i,*bj=b->j;
908   PetscInt       m=Btdense->rmap->n,n=Btdense->cmap->n,j,k,l,col,anz,*btcol,brow,ncolumns;
909   MatScalar      *btval,*btval_den,*ba=b->a;
910   PetscInt       *columns=coloring->columns,*colorforcol=coloring->colorforcol,ncolors=coloring->ncolors;
911 
912   PetscFunctionBegin;
913   btval_den=btdense->v;
914   ierr = PetscMemzero(btval_den,(m*n)*sizeof(MatScalar));CHKERRQ(ierr);
915   for (k=0; k<ncolors; k++) {
916     ncolumns = coloring->ncolumns[k];
917     for (l=0; l<ncolumns; l++) { /* insert a row of B to a column of Btdense */
918       col   = *(columns + colorforcol[k] + l);
919       btcol = bj + bi[col];
920       btval = ba + bi[col];
921       anz   = bi[col+1] - bi[col];
922       for (j=0; j<anz; j++){
923         brow            = btcol[j];
924         btval_den[brow] = btval[j];
925       }
926     }
927     btval_den += m;
928   }
929   PetscFunctionReturn(0);
930 }
931 
932 #undef __FUNCT__
933 #define __FUNCT__ "MatTransColoringApplyDenToSp_SeqAIJ"
934 PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring matcoloring,Mat Cden,Mat Csp)
935 {
936   PetscErrorCode ierr;
937   Mat_SeqAIJ     *csp = (Mat_SeqAIJ*)Csp->data;
938   PetscInt       k,l,*row,*idx,m,ncolors=matcoloring->ncolors,nrows;
939   PetscScalar    *ca_den,*cp_den,*ca=csp->a;
940   PetscInt       *rows=matcoloring->rows,*spidx=matcoloring->columnsforspidx,*colorforrow=matcoloring->colorforrow;
941 
942   PetscFunctionBegin;
943   ierr = MatGetLocalSize(Csp,&m,PETSC_NULL);CHKERRQ(ierr);
944   ierr = MatGetArray(Cden,&ca_den);CHKERRQ(ierr);
945   cp_den = ca_den;
946   for (k=0; k<ncolors; k++) {
947     nrows = matcoloring->nrows[k];
948     row   = rows  + colorforrow[k];
949     idx   = spidx + colorforrow[k];
950     for (l=0; l<nrows; l++){
951       ca[idx[l]] = cp_den[row[l]];
952     }
953     cp_den += m;
954   }
955   ierr = MatRestoreArray(Cden,&ca_den);CHKERRQ(ierr);
956   PetscFunctionReturn(0);
957 }
958 
959 /*
960  MatGetColumnIJ_SeqAIJ_Color() and MatRestoreColumnIJ_SeqAIJ_Color() are customized from
961  MatGetColumnIJ_SeqAIJ() and MatRestoreColumnIJ_SeqAIJ() by adding an output
962  spidx[], index of a->j, to be used for setting 'columnsforspidx' in MatTransposeColoringCreate_SeqAIJ().
963  */
964 #undef __FUNCT__
965 #define __FUNCT__ "MatGetColumnIJ_SeqAIJ_Color"
966 PetscErrorCode MatGetColumnIJ_SeqAIJ_Color(Mat A,PetscInt oshift,PetscBool  symmetric,PetscBool  inodecompressed,PetscInt *nn,PetscInt *ia[],PetscInt *ja[],PetscInt *spidx[],PetscBool  *done)
967 {
968   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data;
969   PetscErrorCode ierr;
970   PetscInt       i,*collengths,*cia,*cja,n = A->cmap->n,m = A->rmap->n;
971   PetscInt       nz = a->i[m],row,*jj,mr,col;
972   PetscInt       *cspidx;
973 
974   PetscFunctionBegin;
975   *nn = n;
976   if (!ia) PetscFunctionReturn(0);
977   if (symmetric) {
978     SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"MatGetColumnIJ_SeqAIJ_Color() not supported for the case symmetric");
979     ierr = MatToSymmetricIJ_SeqAIJ(A->rmap->n,a->i,a->j,0,oshift,ia,ja);CHKERRQ(ierr);
980   } else {
981     ierr = PetscMalloc((n+1)*sizeof(PetscInt),&collengths);CHKERRQ(ierr);
982     ierr = PetscMemzero(collengths,n*sizeof(PetscInt));CHKERRQ(ierr);
983     ierr = PetscMalloc((n+1)*sizeof(PetscInt),&cia);CHKERRQ(ierr);
984     ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&cja);CHKERRQ(ierr);
985     ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&cspidx);CHKERRQ(ierr);
986     jj = a->j;
987     for (i=0; i<nz; i++) {
988       collengths[jj[i]]++;
989     }
990     cia[0] = oshift;
991     for (i=0; i<n; i++) {
992       cia[i+1] = cia[i] + collengths[i];
993     }
994     ierr = PetscMemzero(collengths,n*sizeof(PetscInt));CHKERRQ(ierr);
995     jj   = a->j;
996     for (row=0; row<m; row++) {
997       mr = a->i[row+1] - a->i[row];
998       for (i=0; i<mr; i++) {
999         col = *jj++;
1000         cspidx[cia[col] + collengths[col] - oshift] = a->i[row] + i; /* index of a->j */
1001         cja[cia[col] + collengths[col]++ - oshift] = row + oshift;
1002       }
1003     }
1004     ierr = PetscFree(collengths);CHKERRQ(ierr);
1005     *ia = cia; *ja = cja;
1006     *spidx = cspidx;
1007   }
1008   PetscFunctionReturn(0);
1009 }
1010 
1011 #undef __FUNCT__
1012 #define __FUNCT__ "MatRestoreColumnIJ_SeqAIJ_Color"
1013 PetscErrorCode MatRestoreColumnIJ_SeqAIJ_Color(Mat A,PetscInt oshift,PetscBool  symmetric,PetscBool  inodecompressed,PetscInt *n,PetscInt *ia[],PetscInt *ja[],PetscInt *spidx[],PetscBool  *done)
1014 {
1015   PetscErrorCode ierr;
1016 
1017   PetscFunctionBegin;
1018   if (!ia) PetscFunctionReturn(0);
1019 
1020   ierr = PetscFree(*ia);CHKERRQ(ierr);
1021   ierr = PetscFree(*ja);CHKERRQ(ierr);
1022   ierr = PetscFree(*spidx);CHKERRQ(ierr);
1023   PetscFunctionReturn(0);
1024 }
1025 
1026 #undef __FUNCT__
1027 #define __FUNCT__ "MatTransposeColoringCreate_SeqAIJ"
1028 PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat mat,ISColoring iscoloring,MatTransposeColoring c)
1029 {
1030   PetscErrorCode ierr;
1031   PetscInt       i,n,nrows,N,j,k,m,*row_idx,*ci,*cj,ncols,col,cm;
1032   const PetscInt *is;
1033   PetscInt       nis = iscoloring->n,*rowhit,bs = 1;
1034   IS             *isa;
1035   PetscBool      done;
1036   PetscBool      flg1,flg2;
1037   Mat_SeqAIJ     *csp = (Mat_SeqAIJ*)mat->data;
1038   PetscInt       *colorforrow,*rows,*rows_i,*columnsforspidx,*columnsforspidx_i,*idxhit,*spidx;
1039   PetscInt       *colorforcol,*columns,*columns_i;
1040 
1041   PetscFunctionBegin;
1042   ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr);
1043 
1044   /* this is ugly way to get blocksize but cannot call MatGetBlockSize() because AIJ can have bs > 1 */
1045   ierr = PetscTypeCompare((PetscObject)mat,MATSEQBAIJ,&flg1);CHKERRQ(ierr);
1046   ierr = PetscTypeCompare((PetscObject)mat,MATMPIBAIJ,&flg2);CHKERRQ(ierr);
1047   if (flg1 || flg2) {
1048     ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr);
1049   }
1050 
1051   N         = mat->cmap->N/bs;
1052   c->M      = mat->rmap->N/bs;  /* set total rows, columns and local rows */
1053   c->N      = mat->cmap->N/bs;
1054   c->m      = mat->rmap->N/bs;
1055   c->rstart = 0;
1056 
1057   c->ncolors = nis;
1058   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr);
1059   ierr       = PetscMalloc(nis*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr);
1060   ierr       = PetscMalloc2(csp->nz+1,PetscInt,&rows,csp->nz+1,PetscInt,&columnsforspidx);CHKERRQ(ierr);
1061   ierr       = PetscMalloc((nis+1)*sizeof(PetscInt),&colorforrow);CHKERRQ(ierr);
1062   colorforrow[0]    = 0;
1063   rows_i            = rows;
1064   columnsforspidx_i = columnsforspidx;
1065 
1066   ierr = PetscMalloc((nis+1)*sizeof(PetscInt),&colorforcol);CHKERRQ(ierr);
1067   ierr = PetscMalloc((N+1)*sizeof(PetscInt),&columns);CHKERRQ(ierr);
1068   colorforcol[0] = 0;
1069   columns_i      = columns;
1070 
1071   ierr = MatGetColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,&done);CHKERRQ(ierr); /* column-wise storage of mat */
1072   if (!done) SETERRQ1(((PetscObject)mat)->comm,PETSC_ERR_SUP,"MatGetColumnIJ() not supported for matrix type %s",((PetscObject)mat)->type_name);
1073 
1074   cm = c->m;
1075   ierr = PetscMalloc((cm+1)*sizeof(PetscInt),&rowhit);CHKERRQ(ierr);
1076   ierr = PetscMalloc((cm+1)*sizeof(PetscInt),&idxhit);CHKERRQ(ierr);
1077   for (i=0; i<nis; i++) {
1078     ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr);
1079     ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr);
1080     c->ncolumns[i] = n;
1081     if (n) {
1082       ierr = PetscMemcpy(columns_i,is,n*sizeof(PetscInt));CHKERRQ(ierr);
1083     }
1084     colorforcol[i+1] = colorforcol[i] + n;
1085     columns_i       += n;
1086 
1087     /* fast, crude version requires O(N*N) work */
1088     ierr = PetscMemzero(rowhit,cm*sizeof(PetscInt));CHKERRQ(ierr);
1089 
1090     /* loop over columns*/
1091     for (j=0; j<n; j++) {
1092       col     = is[j];
1093       row_idx = cj + ci[col];
1094       m       = ci[col+1] - ci[col];
1095       /* loop over columns marking them in rowhit */
1096       for (k=0; k<m; k++) {
1097         idxhit[*row_idx]   = spidx[ci[col] + k];
1098         rowhit[*row_idx++] = col + 1;
1099       }
1100     }
1101     /* count the number of hits */
1102     nrows = 0;
1103     for (j=0; j<cm; j++) {
1104       if (rowhit[j]) nrows++;
1105     }
1106     c->nrows[i]      = nrows;
1107     colorforrow[i+1] = colorforrow[i] + nrows;
1108 
1109     nrows       = 0;
1110     for (j=0; j<cm; j++) {
1111       if (rowhit[j]) {
1112         rows_i[nrows]            = j;
1113         columnsforspidx_i[nrows] = idxhit[j];
1114         nrows++;
1115       }
1116     }
1117     ierr = ISRestoreIndices(isa[i],&is);CHKERRQ(ierr);
1118     rows_i += nrows; columnsforspidx_i += nrows;
1119   }
1120   ierr = MatRestoreColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,&done);CHKERRQ(ierr);
1121   ierr = PetscFree(rowhit);CHKERRQ(ierr);
1122   ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr);
1123 #if defined(PETSC_USE_DEBUG)
1124   if (csp->nz != colorforrow[nis]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"csp->nz %d != colorforrow[nis] %d",csp->nz,colorforrow[nis]);
1125 #endif
1126 
1127   c->colorforrow     = colorforrow;
1128   c->rows            = rows;
1129   c->columnsforspidx = columnsforspidx;
1130   c->colorforcol     = colorforcol;
1131   c->columns         = columns;
1132   ierr = PetscFree(idxhit);CHKERRQ(ierr);
1133   PetscFunctionReturn(0);
1134 }
1135