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