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