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