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