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