Actual source code: matmatmult.c

petsc-3.7.0 2016-04-25
Report Typos and Errors
  2: /*
  3:   Defines matrix-matrix product routines for pairs of SeqAIJ matrices
  4:           C = A * B
  5: */

  7: #include <../src/mat/impls/aij/seq/aij.h> /*I "petscmat.h" I*/
  8: #include <../src/mat/utils/freespace.h>
  9: #include <../src/mat/utils/petscheap.h>
 10: #include <petscbt.h>
 11: #include <petsc/private/isimpl.h>
 12: #include <../src/mat/impls/dense/seq/dense.h>

 14: static PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_LLCondensed(Mat,Mat,PetscReal,Mat*);

 18: PETSC_INTERN PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
 19: {
 21:   const char     *algTypes[6] = {"sorted","scalable","scalable_fast","heap","btheap","llcondensed"};
 22:   PetscInt       alg=0; /* set default algorithm */

 25:   if (scall == MAT_INITIAL_MATRIX) {
 26:     PetscObjectOptionsBegin((PetscObject)A);
 27:     PetscOptionsEList("-matmatmult_via","Algorithmic approach","MatMatMult",algTypes,6,algTypes[0],&alg,NULL);
 28:     PetscOptionsEnd();
 29:     PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);
 30:     switch (alg) {
 31:     case 1:
 32:       MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(A,B,fill,C);
 33:       break;
 34:     case 2:
 35:       MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(A,B,fill,C);
 36:       break;
 37:     case 3:
 38:       MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(A,B,fill,C);
 39:       break;
 40:     case 4:
 41:       MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(A,B,fill,C);
 42:       break;
 43:     case 5:
 44:       MatMatMultSymbolic_SeqAIJ_SeqAIJ_LLCondensed(A,B,fill,C);
 45:       break;
 46:     default:
 47:       MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);
 48:      break;
 49:     }
 50:     PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);
 51:   }

 53:   PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);
 54:   (*(*C)->ops->matmultnumeric)(A,B,*C);
 55:   PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);
 56:   return(0);
 57: }

 61: static PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_LLCondensed(Mat A,Mat B,PetscReal fill,Mat *C)
 62: {
 63:   PetscErrorCode     ierr;
 64:   Mat_SeqAIJ         *a =(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
 65:   PetscInt           *ai=a->i,*bi=b->i,*ci,*cj;
 66:   PetscInt           am =A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
 67:   PetscReal          afill;
 68:   PetscInt           i,j,anzi,brow,bnzj,cnzi,*bj,*aj,*lnk,ndouble=0,Crmax;
 69:   PetscTable         ta;
 70:   PetscBT            lnkbt;
 71:   PetscFreeSpaceList free_space=NULL,current_space=NULL;

 74:   /* Get ci and cj */
 75:   /*---------------*/
 76:   /* Allocate ci array, arrays for fill computation and */
 77:   /* free space for accumulating nonzero column info */
 78:   PetscMalloc1(am+2,&ci);
 79:   ci[0] = 0;

 81:   /* create and initialize a linked list */
 82:   Crmax = 2*(b->rmax + (PetscInt)(1.e-2*bn)); /* expected Crmax */;
 83:   if (Crmax > bn) Crmax = bn;
 84:   PetscTableCreate(Crmax,bn,&ta);
 85:   MatRowMergeMax_SeqAIJ(b,bn,ta);
 86:   PetscTableGetCount(ta,&Crmax);
 87:   PetscTableDestroy(&ta);

 89:   PetscLLCondensedCreate(Crmax,bn,&lnk,&lnkbt);

 91:   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
 92:   PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am],bi[bm])),&free_space);

 94:   current_space = free_space;

 96:   /* Determine ci and cj */
 97:   for (i=0; i<am; i++) {
 98:     anzi = ai[i+1] - ai[i];
 99:     aj   = a->j + ai[i];
100:     for (j=0; j<anzi; j++) {
101:       brow = aj[j];
102:       bnzj = bi[brow+1] - bi[brow];
103:       bj   = b->j + bi[brow];
104:       /* add non-zero cols of B into the sorted linked list lnk */
105:       PetscLLCondensedAddSorted(bnzj,bj,lnk,lnkbt);
106:     }
107:     cnzi = lnk[0];

109:     /* If free space is not available, make more free space */
110:     /* Double the amount of total space in the list */
111:     if (current_space->local_remaining<cnzi) {
112:       PetscFreeSpaceGet(PetscIntSumTruncate(cnzi,current_space->total_array_size),&current_space);
113:       ndouble++;
114:     }

116:     /* Copy data into free space, then initialize lnk */
117:     PetscLLCondensedClean(bn,cnzi,current_space->array,lnk,lnkbt);

119:     current_space->array           += cnzi;
120:     current_space->local_used      += cnzi;
121:     current_space->local_remaining -= cnzi;

123:     ci[i+1] = ci[i] + cnzi;
124:   }

126:   /* Column indices are in the list of free space */
127:   /* Allocate space for cj, initialize cj, and */
128:   /* destroy list of free space and other temporary array(s) */
129:   PetscMalloc1(ci[am]+1,&cj);
130:   PetscFreeSpaceContiguous(&free_space,cj);
131:   PetscLLCondensedDestroy(lnk,lnkbt);

133:   /* put together the new symbolic matrix */
134:   MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,C);
135:   MatSetBlockSizesFromMats(*C,A,B);

137:   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
138:   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
139:   c                         = (Mat_SeqAIJ*)((*C)->data);
140:   c->free_a                 = PETSC_FALSE;
141:   c->free_ij                = PETSC_TRUE;
142:   c->nonew                  = 0;
143:   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ; /* fast, needs non-scalable O(bn) array 'abdense' */

145:   /* set MatInfo */
146:   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
147:   if (afill < 1.0) afill = 1.0;
148:   c->maxnz                     = ci[am];
149:   c->nz                        = ci[am];
150:   (*C)->info.mallocs           = ndouble;
151:   (*C)->info.fill_ratio_given  = fill;
152:   (*C)->info.fill_ratio_needed = afill;

154: #if defined(PETSC_USE_INFO)
155:   if (ci[am]) {
156:     PetscInfo3((*C),"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);
157:     PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);
158:   } else {
159:     PetscInfo((*C),"Empty matrix product\n");
160:   }
161: #endif
162:   return(0);
163: }

167: PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
168: {
170:   PetscLogDouble flops=0.0;
171:   Mat_SeqAIJ     *a   = (Mat_SeqAIJ*)A->data;
172:   Mat_SeqAIJ     *b   = (Mat_SeqAIJ*)B->data;
173:   Mat_SeqAIJ     *c   = (Mat_SeqAIJ*)C->data;
174:   PetscInt       *ai  =a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j;
175:   PetscInt       am   =A->rmap->n,cm=C->rmap->n;
176:   PetscInt       i,j,k,anzi,bnzi,cnzi,brow;
177:   PetscScalar    *aa=a->a,*ba=b->a,*baj,*ca,valtmp;
178:   PetscScalar    *ab_dense;

181:   if (!c->a) { /* first call of MatMatMultNumeric_SeqAIJ_SeqAIJ, allocate ca and matmult_abdense */
182:     PetscMalloc1(ci[cm]+1,&ca);
183:     c->a      = ca;
184:     c->free_a = PETSC_TRUE;
185:   } else {
186:     ca        = c->a;
187:   }
188:   if (!c->matmult_abdense) {
189:     PetscCalloc1(B->cmap->N,&ab_dense);
190:     c->matmult_abdense = ab_dense;
191:   } else {
192:     ab_dense = c->matmult_abdense;
193:   }

195:   /* clean old values in C */
196:   PetscMemzero(ca,ci[cm]*sizeof(MatScalar));
197:   /* Traverse A row-wise. */
198:   /* Build the ith row in C by summing over nonzero columns in A, */
199:   /* the rows of B corresponding to nonzeros of A. */
200:   for (i=0; i<am; i++) {
201:     anzi = ai[i+1] - ai[i];
202:     for (j=0; j<anzi; j++) {
203:       brow = aj[j];
204:       bnzi = bi[brow+1] - bi[brow];
205:       bjj  = bj + bi[brow];
206:       baj  = ba + bi[brow];
207:       /* perform dense axpy */
208:       valtmp = aa[j];
209:       for (k=0; k<bnzi; k++) {
210:         ab_dense[bjj[k]] += valtmp*baj[k];
211:       }
212:       flops += 2*bnzi;
213:     }
214:     aj += anzi; aa += anzi;

216:     cnzi = ci[i+1] - ci[i];
217:     for (k=0; k<cnzi; k++) {
218:       ca[k]          += ab_dense[cj[k]];
219:       ab_dense[cj[k]] = 0.0; /* zero ab_dense */
220:     }
221:     flops += cnzi;
222:     cj    += cnzi; ca += cnzi;
223:   }
224:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
225:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
226:   PetscLogFlops(flops);
227:   return(0);
228: }

232: PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(Mat A,Mat B,Mat C)
233: {
235:   PetscLogDouble flops=0.0;
236:   Mat_SeqAIJ     *a   = (Mat_SeqAIJ*)A->data;
237:   Mat_SeqAIJ     *b   = (Mat_SeqAIJ*)B->data;
238:   Mat_SeqAIJ     *c   = (Mat_SeqAIJ*)C->data;
239:   PetscInt       *ai  = a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j;
240:   PetscInt       am   = A->rmap->N,cm=C->rmap->N;
241:   PetscInt       i,j,k,anzi,bnzi,cnzi,brow;
242:   PetscScalar    *aa=a->a,*ba=b->a,*baj,*ca=c->a,valtmp;
243:   PetscInt       nextb;

246:   /* clean old values in C */
247:   PetscMemzero(ca,ci[cm]*sizeof(MatScalar));
248:   /* Traverse A row-wise. */
249:   /* Build the ith row in C by summing over nonzero columns in A, */
250:   /* the rows of B corresponding to nonzeros of A. */
251:   for (i=0; i<am; i++) {
252:     anzi = ai[i+1] - ai[i];
253:     cnzi = ci[i+1] - ci[i];
254:     for (j=0; j<anzi; j++) {
255:       brow = aj[j];
256:       bnzi = bi[brow+1] - bi[brow];
257:       bjj  = bj + bi[brow];
258:       baj  = ba + bi[brow];
259:       /* perform sparse axpy */
260:       valtmp = aa[j];
261:       nextb  = 0;
262:       for (k=0; nextb<bnzi; k++) {
263:         if (cj[k] == bjj[nextb]) { /* ccol == bcol */
264:           ca[k] += valtmp*baj[nextb++];
265:         }
266:       }
267:       flops += 2*bnzi;
268:     }
269:     aj += anzi; aa += anzi;
270:     cj += cnzi; ca += cnzi;
271:   }

273:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
274:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
275:   PetscLogFlops(flops);
276:   return(0);
277: }

281: PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(Mat A,Mat B,PetscReal fill,Mat *C)
282: {
283:   PetscErrorCode     ierr;
284:   Mat_SeqAIJ         *a  = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
285:   PetscInt           *ai = a->i,*bi=b->i,*ci,*cj;
286:   PetscInt           am  = A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
287:   MatScalar          *ca;
288:   PetscReal          afill;
289:   PetscInt           i,j,anzi,brow,bnzj,cnzi,*bj,*aj,*lnk,ndouble=0,Crmax;
290:   PetscTable         ta;
291:   PetscFreeSpaceList free_space=NULL,current_space=NULL;

294:   /* Get ci and cj - same as MatMatMultSymbolic_SeqAIJ_SeqAIJ except using PetscLLxxx_fast() */
295:   /*-----------------------------------------------------------------------------------------*/
296:   /* Allocate arrays for fill computation and free space for accumulating nonzero column */
297:   PetscMalloc1(am+2,&ci);
298:   ci[0] = 0;

300:   /* create and initialize a linked list */
301:   Crmax = 2*(b->rmax + (PetscInt)(1.e-2*bn)); /* expected Crmax */
302:   if (Crmax > bn) Crmax = bn;
303:   PetscTableCreate(Crmax,bn,&ta);
304:   MatRowMergeMax_SeqAIJ(b,bn,ta);
305:   PetscTableGetCount(ta,&Crmax);
306:   PetscTableDestroy(&ta);

308:   PetscLLCondensedCreate_fast(Crmax,&lnk);

310:   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
311:   PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am],bi[bm])),&free_space);
312:   current_space = free_space;

314:   /* Determine ci and cj */
315:   for (i=0; i<am; i++) {
316:     anzi = ai[i+1] - ai[i];
317:     aj   = a->j + ai[i];
318:     for (j=0; j<anzi; j++) {
319:       brow = aj[j];
320:       bnzj = bi[brow+1] - bi[brow];
321:       bj   = b->j + bi[brow];
322:       /* add non-zero cols of B into the sorted linked list lnk */
323:       PetscLLCondensedAddSorted_fast(bnzj,bj,lnk);
324:     }
325:     cnzi = lnk[1];

327:     /* If free space is not available, make more free space */
328:     /* Double the amount of total space in the list */
329:     if (current_space->local_remaining<cnzi) {
330:       PetscFreeSpaceGet(PetscIntSumTruncate(cnzi,current_space->total_array_size),&current_space);
331:       ndouble++;
332:     }

334:     /* Copy data into free space, then initialize lnk */
335:     PetscLLCondensedClean_fast(cnzi,current_space->array,lnk);

337:     current_space->array           += cnzi;
338:     current_space->local_used      += cnzi;
339:     current_space->local_remaining -= cnzi;

341:     ci[i+1] = ci[i] + cnzi;
342:   }

344:   /* Column indices are in the list of free space */
345:   /* Allocate space for cj, initialize cj, and */
346:   /* destroy list of free space and other temporary array(s) */
347:   PetscMalloc1(ci[am]+1,&cj);
348:   PetscFreeSpaceContiguous(&free_space,cj);
349:   PetscLLCondensedDestroy_fast(lnk);

351:   /* Allocate space for ca */
352:   PetscMalloc1(ci[am]+1,&ca);
353:   PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));

355:   /* put together the new symbolic matrix */
356:   MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,ca,C);
357:   MatSetBlockSizesFromMats(*C,A,B);

359:   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
360:   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
361:   c          = (Mat_SeqAIJ*)((*C)->data);
362:   c->free_a  = PETSC_TRUE;
363:   c->free_ij = PETSC_TRUE;
364:   c->nonew   = 0;

366:   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable; /* slower, less memory */

368:   /* set MatInfo */
369:   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
370:   if (afill < 1.0) afill = 1.0;
371:   c->maxnz                     = ci[am];
372:   c->nz                        = ci[am];
373:   (*C)->info.mallocs           = ndouble;
374:   (*C)->info.fill_ratio_given  = fill;
375:   (*C)->info.fill_ratio_needed = afill;

377: #if defined(PETSC_USE_INFO)
378:   if (ci[am]) {
379:     PetscInfo3((*C),"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);
380:     PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);
381:   } else {
382:     PetscInfo((*C),"Empty matrix product\n");
383:   }
384: #endif
385:   return(0);
386: }


391: PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(Mat A,Mat B,PetscReal fill,Mat *C)
392: {
393:   PetscErrorCode     ierr;
394:   Mat_SeqAIJ         *a  = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
395:   PetscInt           *ai = a->i,*bi=b->i,*ci,*cj;
396:   PetscInt           am  = A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
397:   MatScalar          *ca;
398:   PetscReal          afill;
399:   PetscInt           i,j,anzi,brow,bnzj,cnzi,*bj,*aj,*lnk,ndouble=0,Crmax;
400:   PetscTable         ta;
401:   PetscFreeSpaceList free_space=NULL,current_space=NULL;

404:   /* Get ci and cj - same as MatMatMultSymbolic_SeqAIJ_SeqAIJ except using PetscLLxxx_Scalalbe() */
405:   /*---------------------------------------------------------------------------------------------*/
406:   /* Allocate arrays for fill computation and free space for accumulating nonzero column */
407:   PetscMalloc1(am+2,&ci);
408:   ci[0] = 0;

410:   /* create and initialize a linked list */
411:   Crmax = 2*(b->rmax  + (PetscInt)(1.e-2*bn)); /* expected Crmax */
412:   if (Crmax > bn) Crmax = bn;
413:   PetscTableCreate(Crmax,bn,&ta);
414:   MatRowMergeMax_SeqAIJ(b,bn,ta);
415:   PetscTableGetCount(ta,&Crmax);
416:   PetscTableDestroy(&ta);

418:   PetscLLCondensedCreate_Scalable(Crmax,&lnk);

420:   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
421:   PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am],bi[bm])),&free_space);
422:   current_space = free_space;

424:   /* Determine ci and cj */
425:   for (i=0; i<am; i++) {
426:     anzi = ai[i+1] - ai[i];
427:     aj   = a->j + ai[i];
428:     for (j=0; j<anzi; j++) {
429:       brow = aj[j];
430:       bnzj = bi[brow+1] - bi[brow];
431:       bj   = b->j + bi[brow];
432:       /* add non-zero cols of B into the sorted linked list lnk */
433:       PetscLLCondensedAddSorted_Scalable(bnzj,bj,lnk);
434:     }
435:     cnzi = lnk[0];

437:     /* If free space is not available, make more free space */
438:     /* Double the amount of total space in the list */
439:     if (current_space->local_remaining<cnzi) {
440:       PetscFreeSpaceGet(PetscIntSumTruncate(cnzi,current_space->total_array_size),&current_space);
441:       ndouble++;
442:     }

444:     /* Copy data into free space, then initialize lnk */
445:     PetscLLCondensedClean_Scalable(cnzi,current_space->array,lnk);

447:     current_space->array           += cnzi;
448:     current_space->local_used      += cnzi;
449:     current_space->local_remaining -= cnzi;

451:     ci[i+1] = ci[i] + cnzi;
452:   }

454:   /* Column indices are in the list of free space */
455:   /* Allocate space for cj, initialize cj, and */
456:   /* destroy list of free space and other temporary array(s) */
457:   PetscMalloc1(ci[am]+1,&cj);
458:   PetscFreeSpaceContiguous(&free_space,cj);
459:   PetscLLCondensedDestroy_Scalable(lnk);

461:   /* Allocate space for ca */
462:   /*-----------------------*/
463:   PetscMalloc1(ci[am]+1,&ca);
464:   PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));

466:   /* put together the new symbolic matrix */
467:   MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,ca,C);
468:   MatSetBlockSizesFromMats(*C,A,B);

470:   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
471:   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
472:   c          = (Mat_SeqAIJ*)((*C)->data);
473:   c->free_a  = PETSC_TRUE;
474:   c->free_ij = PETSC_TRUE;
475:   c->nonew   = 0;

477:   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable; /* slower, less memory */

479:   /* set MatInfo */
480:   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
481:   if (afill < 1.0) afill = 1.0;
482:   c->maxnz                     = ci[am];
483:   c->nz                        = ci[am];
484:   (*C)->info.mallocs           = ndouble;
485:   (*C)->info.fill_ratio_given  = fill;
486:   (*C)->info.fill_ratio_needed = afill;

488: #if defined(PETSC_USE_INFO)
489:   if (ci[am]) {
490:     PetscInfo3((*C),"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);
491:     PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);
492:   } else {
493:     PetscInfo((*C),"Empty matrix product\n");
494:   }
495: #endif
496:   return(0);
497: }

501: PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(Mat A,Mat B,PetscReal fill,Mat *C)
502: {
503:   PetscErrorCode     ierr;
504:   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
505:   const PetscInt     *ai=a->i,*bi=b->i,*aj=a->j,*bj=b->j;
506:   PetscInt           *ci,*cj,*bb;
507:   PetscInt           am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
508:   PetscReal          afill;
509:   PetscInt           i,j,col,ndouble = 0;
510:   PetscFreeSpaceList free_space=NULL,current_space=NULL;
511:   PetscHeap          h;

514:   /* Get ci and cj - by merging sorted rows using a heap */
515:   /*---------------------------------------------------------------------------------------------*/
516:   /* Allocate arrays for fill computation and free space for accumulating nonzero column */
517:   PetscMalloc1(am+2,&ci);
518:   ci[0] = 0;

520:   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
521:   PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am],bi[bm])),&free_space);
522:   current_space = free_space;

524:   PetscHeapCreate(a->rmax,&h);
525:   PetscMalloc1(a->rmax,&bb);

527:   /* Determine ci and cj */
528:   for (i=0; i<am; i++) {
529:     const PetscInt anzi  = ai[i+1] - ai[i]; /* number of nonzeros in this row of A, this is the number of rows of B that we merge */
530:     const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */
531:     ci[i+1] = ci[i];
532:     /* Populate the min heap */
533:     for (j=0; j<anzi; j++) {
534:       bb[j] = bi[acol[j]];         /* bb points at the start of the row */
535:       if (bb[j] < bi[acol[j]+1]) { /* Add if row is nonempty */
536:         PetscHeapAdd(h,j,bj[bb[j]++]);
537:       }
538:     }
539:     /* Pick off the min element, adding it to free space */
540:     PetscHeapPop(h,&j,&col);
541:     while (j >= 0) {
542:       if (current_space->local_remaining < 1) { /* double the size, but don't exceed 16 MiB */
543:         PetscFreeSpaceGet(PetscMin(PetscIntMultTruncate(2,current_space->total_array_size),16 << 20),&current_space);
544:         ndouble++;
545:       }
546:       *(current_space->array++) = col;
547:       current_space->local_used++;
548:       current_space->local_remaining--;
549:       ci[i+1]++;

551:       /* stash if anything else remains in this row of B */
552:       if (bb[j] < bi[acol[j]+1]) {PetscHeapStash(h,j,bj[bb[j]++]);}
553:       while (1) {               /* pop and stash any other rows of B that also had an entry in this column */
554:         PetscInt j2,col2;
555:         PetscHeapPeek(h,&j2,&col2);
556:         if (col2 != col) break;
557:         PetscHeapPop(h,&j2,&col2);
558:         if (bb[j2] < bi[acol[j2]+1]) {PetscHeapStash(h,j2,bj[bb[j2]++]);}
559:       }
560:       /* Put any stashed elements back into the min heap */
561:       PetscHeapUnstash(h);
562:       PetscHeapPop(h,&j,&col);
563:     }
564:   }
565:   PetscFree(bb);
566:   PetscHeapDestroy(&h);

568:   /* Column indices are in the list of free space */
569:   /* Allocate space for cj, initialize cj, and */
570:   /* destroy list of free space and other temporary array(s) */
571:   PetscMalloc1(ci[am],&cj);
572:   PetscFreeSpaceContiguous(&free_space,cj);

574:   /* put together the new symbolic matrix */
575:   MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,C);
576:   MatSetBlockSizesFromMats(*C,A,B);

578:   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
579:   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
580:   c          = (Mat_SeqAIJ*)((*C)->data);
581:   c->free_a  = PETSC_TRUE;
582:   c->free_ij = PETSC_TRUE;
583:   c->nonew   = 0;

585:   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ;

587:   /* set MatInfo */
588:   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
589:   if (afill < 1.0) afill = 1.0;
590:   c->maxnz                     = ci[am];
591:   c->nz                        = ci[am];
592:   (*C)->info.mallocs           = ndouble;
593:   (*C)->info.fill_ratio_given  = fill;
594:   (*C)->info.fill_ratio_needed = afill;

596: #if defined(PETSC_USE_INFO)
597:   if (ci[am]) {
598:     PetscInfo3((*C),"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);
599:     PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);
600:   } else {
601:     PetscInfo((*C),"Empty matrix product\n");
602:   }
603: #endif
604:   return(0);
605: }

609: PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(Mat A,Mat B,PetscReal fill,Mat *C)
610: {
611:   PetscErrorCode     ierr;
612:   Mat_SeqAIJ         *a  = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
613:   const PetscInt     *ai = a->i,*bi=b->i,*aj=a->j,*bj=b->j;
614:   PetscInt           *ci,*cj,*bb;
615:   PetscInt           am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
616:   PetscReal          afill;
617:   PetscInt           i,j,col,ndouble = 0;
618:   PetscFreeSpaceList free_space=NULL,current_space=NULL;
619:   PetscHeap          h;
620:   PetscBT            bt;

623:   /* Get ci and cj - using a heap for the sorted rows, but use BT so that each index is only added once */
624:   /*---------------------------------------------------------------------------------------------*/
625:   /* Allocate arrays for fill computation and free space for accumulating nonzero column */
626:   PetscMalloc1(am+2,&ci);
627:   ci[0] = 0;

629:   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
630:   PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am],bi[bm])),&free_space);

632:   current_space = free_space;

634:   PetscHeapCreate(a->rmax,&h);
635:   PetscMalloc1(a->rmax,&bb);
636:   PetscBTCreate(bn,&bt);

638:   /* Determine ci and cj */
639:   for (i=0; i<am; i++) {
640:     const PetscInt anzi  = ai[i+1] - ai[i]; /* number of nonzeros in this row of A, this is the number of rows of B that we merge */
641:     const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */
642:     const PetscInt *fptr = current_space->array; /* Save beginning of the row so we can clear the BT later */
643:     ci[i+1] = ci[i];
644:     /* Populate the min heap */
645:     for (j=0; j<anzi; j++) {
646:       PetscInt brow = acol[j];
647:       for (bb[j] = bi[brow]; bb[j] < bi[brow+1]; bb[j]++) {
648:         PetscInt bcol = bj[bb[j]];
649:         if (!PetscBTLookupSet(bt,bcol)) { /* new entry */
650:           PetscHeapAdd(h,j,bcol);
651:           bb[j]++;
652:           break;
653:         }
654:       }
655:     }
656:     /* Pick off the min element, adding it to free space */
657:     PetscHeapPop(h,&j,&col);
658:     while (j >= 0) {
659:       if (current_space->local_remaining < 1) { /* double the size, but don't exceed 16 MiB */
660:         fptr = NULL;                      /* need PetscBTMemzero */
661:         PetscFreeSpaceGet(PetscMin(PetscIntMultTruncate(2,current_space->total_array_size),16 << 20),&current_space);
662:         ndouble++;
663:       }
664:       *(current_space->array++) = col;
665:       current_space->local_used++;
666:       current_space->local_remaining--;
667:       ci[i+1]++;

669:       /* stash if anything else remains in this row of B */
670:       for (; bb[j] < bi[acol[j]+1]; bb[j]++) {
671:         PetscInt bcol = bj[bb[j]];
672:         if (!PetscBTLookupSet(bt,bcol)) { /* new entry */
673:           PetscHeapAdd(h,j,bcol);
674:           bb[j]++;
675:           break;
676:         }
677:       }
678:       PetscHeapPop(h,&j,&col);
679:     }
680:     if (fptr) {                 /* Clear the bits for this row */
681:       for (; fptr<current_space->array; fptr++) {PetscBTClear(bt,*fptr);}
682:     } else {                    /* We reallocated so we don't remember (easily) how to clear only the bits we changed */
683:       PetscBTMemzero(bn,bt);
684:     }
685:   }
686:   PetscFree(bb);
687:   PetscHeapDestroy(&h);
688:   PetscBTDestroy(&bt);

690:   /* Column indices are in the list of free space */
691:   /* Allocate space for cj, initialize cj, and */
692:   /* destroy list of free space and other temporary array(s) */
693:   PetscMalloc1(ci[am],&cj);
694:   PetscFreeSpaceContiguous(&free_space,cj);

696:   /* put together the new symbolic matrix */
697:   MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,C);
698:   MatSetBlockSizesFromMats(*C,A,B);

700:   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
701:   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
702:   c          = (Mat_SeqAIJ*)((*C)->data);
703:   c->free_a  = PETSC_TRUE;
704:   c->free_ij = PETSC_TRUE;
705:   c->nonew   = 0;

707:   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ;

709:   /* set MatInfo */
710:   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
711:   if (afill < 1.0) afill = 1.0;
712:   c->maxnz                     = ci[am];
713:   c->nz                        = ci[am];
714:   (*C)->info.mallocs           = ndouble;
715:   (*C)->info.fill_ratio_given  = fill;
716:   (*C)->info.fill_ratio_needed = afill;

718: #if defined(PETSC_USE_INFO)
719:   if (ci[am]) {
720:     PetscInfo3((*C),"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);
721:     PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);
722:   } else {
723:     PetscInfo((*C),"Empty matrix product\n");
724:   }
725: #endif
726:   return(0);
727: }

731: /* concatenate unique entries and then sort */
732: PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
733: {
734:   PetscErrorCode     ierr;
735:   Mat_SeqAIJ         *a  = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
736:   const PetscInt     *ai = a->i,*bi=b->i,*aj=a->j,*bj=b->j;
737:   PetscInt           *ci,*cj;
738:   PetscInt           am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
739:   PetscReal          afill;
740:   PetscInt           i,j,ndouble = 0;
741:   PetscSegBuffer     seg,segrow;
742:   char               *seen;

745:   PetscMalloc1(am+1,&ci);
746:   ci[0] = 0;

748:   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
749:   PetscSegBufferCreate(sizeof(PetscInt),(PetscInt)(fill*(ai[am]+bi[bm])),&seg);
750:   PetscSegBufferCreate(sizeof(PetscInt),100,&segrow);
751:   PetscMalloc1(bn,&seen);
752:   PetscMemzero(seen,bn*sizeof(char));

754:   /* Determine ci and cj */
755:   for (i=0; i<am; i++) {
756:     const PetscInt anzi  = ai[i+1] - ai[i]; /* number of nonzeros in this row of A, this is the number of rows of B that we merge */
757:     const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */
758:     PetscInt packlen = 0,*PETSC_RESTRICT crow;
759:     /* Pack segrow */
760:     for (j=0; j<anzi; j++) {
761:       PetscInt brow = acol[j],bjstart = bi[brow],bjend = bi[brow+1],k;
762:       for (k=bjstart; k<bjend; k++) {
763:         PetscInt bcol = bj[k];
764:         if (!seen[bcol]) { /* new entry */
765:           PetscInt *PETSC_RESTRICT slot;
766:           PetscSegBufferGetInts(segrow,1,&slot);
767:           *slot = bcol;
768:           seen[bcol] = 1;
769:           packlen++;
770:         }
771:       }
772:     }
773:     PetscSegBufferGetInts(seg,packlen,&crow);
774:     PetscSegBufferExtractTo(segrow,crow);
775:     PetscSortInt(packlen,crow);
776:     ci[i+1] = ci[i] + packlen;
777:     for (j=0; j<packlen; j++) seen[crow[j]] = 0;
778:   }
779:   PetscSegBufferDestroy(&segrow);
780:   PetscFree(seen);

782:   /* Column indices are in the segmented buffer */
783:   PetscSegBufferExtractAlloc(seg,&cj);
784:   PetscSegBufferDestroy(&seg);

786:   /* put together the new symbolic matrix */
787:   MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,C);
788:   MatSetBlockSizesFromMats(*C,A,B);

790:   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
791:   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
792:   c          = (Mat_SeqAIJ*)((*C)->data);
793:   c->free_a  = PETSC_TRUE;
794:   c->free_ij = PETSC_TRUE;
795:   c->nonew   = 0;

797:   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ;

799:   /* set MatInfo */
800:   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
801:   if (afill < 1.0) afill = 1.0;
802:   c->maxnz                     = ci[am];
803:   c->nz                        = ci[am];
804:   (*C)->info.mallocs           = ndouble;
805:   (*C)->info.fill_ratio_given  = fill;
806:   (*C)->info.fill_ratio_needed = afill;

808: #if defined(PETSC_USE_INFO)
809:   if (ci[am]) {
810:     PetscInfo3((*C),"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);
811:     PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);
812:   } else {
813:     PetscInfo((*C),"Empty matrix product\n");
814:   }
815: #endif
816:   return(0);
817: }

819: /* This routine is not used. Should be removed! */
822: PetscErrorCode MatMatTransposeMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
823: {

827:   if (scall == MAT_INITIAL_MATRIX) {
828:     PetscLogEventBegin(MAT_MatTransposeMultSymbolic,A,B,0,0);
829:     MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);
830:     PetscLogEventEnd(MAT_MatTransposeMultSymbolic,A,B,0,0);
831:   }
832:   PetscLogEventBegin(MAT_MatTransposeMultNumeric,A,B,0,0);
833:   MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);
834:   PetscLogEventEnd(MAT_MatTransposeMultNumeric,A,B,0,0);
835:   return(0);
836: }

840: PetscErrorCode MatDestroy_SeqAIJ_MatMatMultTrans(Mat A)
841: {
842:   PetscErrorCode      ierr;
843:   Mat_SeqAIJ          *a=(Mat_SeqAIJ*)A->data;
844:   Mat_MatMatTransMult *abt=a->abt;

847:   (abt->destroy)(A);
848:   MatTransposeColoringDestroy(&abt->matcoloring);
849:   MatDestroy(&abt->Bt_den);
850:   MatDestroy(&abt->ABt_den);
851:   PetscFree(abt);
852:   return(0);
853: }

857: PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
858: {
859:   PetscErrorCode      ierr;
860:   Mat                 Bt;
861:   PetscInt            *bti,*btj;
862:   Mat_MatMatTransMult *abt;
863:   Mat_SeqAIJ          *c;

866:   /* create symbolic Bt */
867:   MatGetSymbolicTranspose_SeqAIJ(B,&bti,&btj);
868:   MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,B->cmap->n,B->rmap->n,bti,btj,NULL,&Bt);
869:   MatSetBlockSizes(Bt,PetscAbs(A->cmap->bs),PetscAbs(B->cmap->bs));

871:   /* get symbolic C=A*Bt */
872:   MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,Bt,fill,C);

874:   /* create a supporting struct for reuse intermidiate dense matrices with matcoloring */
875:   PetscNew(&abt);
876:   c      = (Mat_SeqAIJ*)(*C)->data;
877:   c->abt = abt;

879:   abt->usecoloring = PETSC_FALSE;
880:   abt->destroy     = (*C)->ops->destroy;
881:   (*C)->ops->destroy     = MatDestroy_SeqAIJ_MatMatMultTrans;

883:   PetscOptionsGetBool(((PetscObject)A)->options,NULL,"-matmattransmult_color",&abt->usecoloring,NULL);
884:   if (abt->usecoloring) {
885:     /* Create MatTransposeColoring from symbolic C=A*B^T */
886:     MatTransposeColoring matcoloring;
887:     MatColoring          coloring;
888:     ISColoring           iscoloring;
889:     Mat                  Bt_dense,C_dense;
890:     Mat_SeqAIJ           *c=(Mat_SeqAIJ*)(*C)->data;
891:     /* inode causes memory problem, don't know why */
892:     if (c->inode.use) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MAT_USE_INODES is not supported. Use '-mat_no_inode'");

894:     MatColoringCreate(*C,&coloring);
895:     MatColoringSetDistance(coloring,2);
896:     MatColoringSetType(coloring,MATCOLORINGSL);
897:     MatColoringSetFromOptions(coloring);
898:     MatColoringApply(coloring,&iscoloring);
899:     MatColoringDestroy(&coloring);
900:     MatTransposeColoringCreate(*C,iscoloring,&matcoloring);

902:     abt->matcoloring = matcoloring;

904:     ISColoringDestroy(&iscoloring);

906:     /* Create Bt_dense and C_dense = A*Bt_dense */
907:     MatCreate(PETSC_COMM_SELF,&Bt_dense);
908:     MatSetSizes(Bt_dense,A->cmap->n,matcoloring->ncolors,A->cmap->n,matcoloring->ncolors);
909:     MatSetType(Bt_dense,MATSEQDENSE);
910:     MatSeqDenseSetPreallocation(Bt_dense,NULL);

912:     Bt_dense->assembled = PETSC_TRUE;
913:     abt->Bt_den   = Bt_dense;

915:     MatCreate(PETSC_COMM_SELF,&C_dense);
916:     MatSetSizes(C_dense,A->rmap->n,matcoloring->ncolors,A->rmap->n,matcoloring->ncolors);
917:     MatSetType(C_dense,MATSEQDENSE);
918:     MatSeqDenseSetPreallocation(C_dense,NULL);

920:     Bt_dense->assembled = PETSC_TRUE;
921:     abt->ABt_den  = C_dense;

923: #if defined(PETSC_USE_INFO)
924:     {
925:       Mat_SeqAIJ *c = (Mat_SeqAIJ*)(*C)->data;
926:       PetscInfo7(*C,"Use coloring of C=A*B^T; B^T: %D %D, Bt_dense: %D,%D; Cnz %D / (cm*ncolors %D) = %g\n",B->cmap->n,B->rmap->n,Bt_dense->rmap->n,Bt_dense->cmap->n,c->nz,A->rmap->n*matcoloring->ncolors,(PetscReal)(c->nz)/(A->rmap->n*matcoloring->ncolors));
927:     }
928: #endif
929:   }
930:   /* clean up */
931:   MatDestroy(&Bt);
932:   MatRestoreSymbolicTranspose_SeqAIJ(B,&bti,&btj);
933:   return(0);
934: }

938: PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
939: {
940:   PetscErrorCode      ierr;
941:   Mat_SeqAIJ          *a   =(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data;
942:   PetscInt            *ai  =a->i,*aj=a->j,*bi=b->i,*bj=b->j,anzi,bnzj,nexta,nextb,*acol,*bcol,brow;
943:   PetscInt            cm   =C->rmap->n,*ci=c->i,*cj=c->j,i,j,cnzi,*ccol;
944:   PetscLogDouble      flops=0.0;
945:   MatScalar           *aa  =a->a,*aval,*ba=b->a,*bval,*ca,*cval;
946:   Mat_MatMatTransMult *abt = c->abt;

949:   /* clear old values in C */
950:   if (!c->a) {
951:     PetscMalloc1(ci[cm]+1,&ca);
952:     c->a      = ca;
953:     c->free_a = PETSC_TRUE;
954:   } else {
955:     ca =  c->a;
956:   }
957:   PetscMemzero(ca,ci[cm]*sizeof(MatScalar));

959:   if (abt->usecoloring) {
960:     MatTransposeColoring matcoloring = abt->matcoloring;
961:     Mat                  Bt_dense,C_dense = abt->ABt_den;

963:     /* Get Bt_dense by Apply MatTransposeColoring to B */
964:     Bt_dense = abt->Bt_den;
965:     MatTransColoringApplySpToDen(matcoloring,B,Bt_dense);

967:     /* C_dense = A*Bt_dense */
968:     MatMatMultNumeric_SeqAIJ_SeqDense(A,Bt_dense,C_dense);

970:     /* Recover C from C_dense */
971:     MatTransColoringApplyDenToSp(matcoloring,C_dense,C);
972:     return(0);
973:   }

975:   for (i=0; i<cm; i++) {
976:     anzi = ai[i+1] - ai[i];
977:     acol = aj + ai[i];
978:     aval = aa + ai[i];
979:     cnzi = ci[i+1] - ci[i];
980:     ccol = cj + ci[i];
981:     cval = ca + ci[i];
982:     for (j=0; j<cnzi; j++) {
983:       brow = ccol[j];
984:       bnzj = bi[brow+1] - bi[brow];
985:       bcol = bj + bi[brow];
986:       bval = ba + bi[brow];

988:       /* perform sparse inner-product c(i,j)=A[i,:]*B[j,:]^T */
989:       nexta = 0; nextb = 0;
990:       while (nexta<anzi && nextb<bnzj) {
991:         while (nexta < anzi && acol[nexta] < bcol[nextb]) nexta++;
992:         if (nexta == anzi) break;
993:         while (nextb < bnzj && acol[nexta] > bcol[nextb]) nextb++;
994:         if (nextb == bnzj) break;
995:         if (acol[nexta] == bcol[nextb]) {
996:           cval[j] += aval[nexta]*bval[nextb];
997:           nexta++; nextb++;
998:           flops += 2;
999:         }
1000:       }
1001:     }
1002:   }
1003:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
1004:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
1005:   PetscLogFlops(flops);
1006:   return(0);
1007: }

1011: PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
1012: {

1016:   if (scall == MAT_INITIAL_MATRIX) {
1017:     PetscLogEventBegin(MAT_TransposeMatMultSymbolic,A,B,0,0);
1018:     MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);
1019:     PetscLogEventEnd(MAT_TransposeMatMultSymbolic,A,B,0,0);
1020:   }
1021:   PetscLogEventBegin(MAT_TransposeMatMultNumeric,A,B,0,0);
1022:   MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);
1023:   PetscLogEventEnd(MAT_TransposeMatMultNumeric,A,B,0,0);
1024:   return(0);
1025: }

1029: PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
1030: {
1032:   Mat            At;
1033:   PetscInt       *ati,*atj;

1036:   /* create symbolic At */
1037:   MatGetSymbolicTranspose_SeqAIJ(A,&ati,&atj);
1038:   MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,A->cmap->n,A->rmap->n,ati,atj,NULL,&At);
1039:   MatSetBlockSizes(At,PetscAbs(A->cmap->bs),PetscAbs(B->cmap->bs));

1041:   /* get symbolic C=At*B */
1042:   MatMatMultSymbolic_SeqAIJ_SeqAIJ(At,B,fill,C);

1044:   /* clean up */
1045:   MatDestroy(&At);
1046:   MatRestoreSymbolicTranspose_SeqAIJ(A,&ati,&atj);
1047:   return(0);
1048: }

1052: PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
1053: {
1055:   Mat_SeqAIJ     *a   =(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data;
1056:   PetscInt       am   =A->rmap->n,anzi,*ai=a->i,*aj=a->j,*bi=b->i,*bj,bnzi,nextb;
1057:   PetscInt       cm   =C->rmap->n,*ci=c->i,*cj=c->j,crow,*cjj,i,j,k;
1058:   PetscLogDouble flops=0.0;
1059:   MatScalar      *aa  =a->a,*ba,*ca,*caj;

1062:   if (!c->a) {
1063:     PetscMalloc1(ci[cm]+1,&ca);

1065:     c->a      = ca;
1066:     c->free_a = PETSC_TRUE;
1067:   } else {
1068:     ca = c->a;
1069:   }
1070:   /* clear old values in C */
1071:   PetscMemzero(ca,ci[cm]*sizeof(MatScalar));

1073:   /* compute A^T*B using outer product (A^T)[:,i]*B[i,:] */
1074:   for (i=0; i<am; i++) {
1075:     bj   = b->j + bi[i];
1076:     ba   = b->a + bi[i];
1077:     bnzi = bi[i+1] - bi[i];
1078:     anzi = ai[i+1] - ai[i];
1079:     for (j=0; j<anzi; j++) {
1080:       nextb = 0;
1081:       crow  = *aj++;
1082:       cjj   = cj + ci[crow];
1083:       caj   = ca + ci[crow];
1084:       /* perform sparse axpy operation.  Note cjj includes bj. */
1085:       for (k=0; nextb<bnzi; k++) {
1086:         if (cjj[k] == *(bj+nextb)) { /* ccol == bcol */
1087:           caj[k] += (*aa)*(*(ba+nextb));
1088:           nextb++;
1089:         }
1090:       }
1091:       flops += 2*bnzi;
1092:       aa++;
1093:     }
1094:   }

1096:   /* Assemble the final matrix and clean up */
1097:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
1098:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
1099:   PetscLogFlops(flops);
1100:   return(0);
1101: }

1105: PETSC_INTERN PetscErrorCode MatMatMult_SeqAIJ_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
1106: {

1110:   if (scall == MAT_INITIAL_MATRIX) {
1111:     PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);
1112:     MatMatMultSymbolic_SeqAIJ_SeqDense(A,B,fill,C);
1113:     PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);
1114:   }
1115:   PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);
1116:   MatMatMultNumeric_SeqAIJ_SeqDense(A,B,*C);
1117:   PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);
1118:   return(0);
1119: }

1123: PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
1124: {

1128:   MatMatMultSymbolic_SeqDense_SeqDense(A,B,0.0,C);

1130:   (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqDense;
1131:   return(0);
1132: }

1136: PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
1137: {
1138:   Mat_SeqAIJ        *a=(Mat_SeqAIJ*)A->data;
1139:   Mat_SeqDense      *bd = (Mat_SeqDense*)B->data;
1140:   PetscErrorCode    ierr;
1141:   PetscScalar       *c,*b,r1,r2,r3,r4,*c1,*c2,*c3,*c4,aatmp;
1142:   const PetscScalar *aa,*b1,*b2,*b3,*b4;
1143:   const PetscInt    *aj;
1144:   PetscInt          cm=C->rmap->n,cn=B->cmap->n,bm=bd->lda,am=A->rmap->n;
1145:   PetscInt          am4=4*am,bm4=4*bm,col,i,j,n,ajtmp;

1148:   if (!cm || !cn) return(0);
1149:   if (B->rmap->n != 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,B->rmap->n);
1150:   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);
1151:   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);
1152:   b = bd->v;
1153:   MatDenseGetArray(C,&c);
1154:   b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;
1155:   c1 = c; c2 = c1 + am; c3 = c2 + am; c4 = c3 + am;
1156:   for (col=0; col<cn-4; col += 4) {  /* over columns of C */
1157:     for (i=0; i<am; i++) {        /* over rows of C in those columns */
1158:       r1 = r2 = r3 = r4 = 0.0;
1159:       n  = a->i[i+1] - a->i[i];
1160:       aj = a->j + a->i[i];
1161:       aa = a->a + a->i[i];
1162:       for (j=0; j<n; j++) {
1163:         aatmp = aa[j]; ajtmp = aj[j];
1164:         r1 += aatmp*b1[ajtmp];
1165:         r2 += aatmp*b2[ajtmp];
1166:         r3 += aatmp*b3[ajtmp];
1167:         r4 += aatmp*b4[ajtmp];
1168:       }
1169:       c1[i] = r1;
1170:       c2[i] = r2;
1171:       c3[i] = r3;
1172:       c4[i] = r4;
1173:     }
1174:     b1 += bm4; b2 += bm4; b3 += bm4; b4 += bm4;
1175:     c1 += am4; c2 += am4; c3 += am4; c4 += am4;
1176:   }
1177:   for (; col<cn; col++) {   /* over extra columns of C */
1178:     for (i=0; i<am; i++) {  /* over rows of C in those columns */
1179:       r1 = 0.0;
1180:       n  = a->i[i+1] - a->i[i];
1181:       aj = a->j + a->i[i];
1182:       aa = a->a + a->i[i];
1183:       for (j=0; j<n; j++) {
1184:         r1 += aa[j]*b1[aj[j]];
1185:       }
1186:       c1[i] = r1;
1187:     }
1188:     b1 += bm;
1189:     c1 += am;
1190:   }
1191:   PetscLogFlops(cn*(2.0*a->nz));
1192:   MatDenseRestoreArray(C,&c);
1193:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
1194:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
1195:   return(0);
1196: }

1198: /*
1199:    Note very similar to MatMult_SeqAIJ(), should generate both codes from same base
1200: */
1203: PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
1204: {
1205:   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
1206:   Mat_SeqDense   *bd = (Mat_SeqDense*)B->data;
1208:   PetscScalar    *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4;
1209:   MatScalar      *aa;
1210:   PetscInt       cm  = C->rmap->n, cn=B->cmap->n, bm=bd->lda, col, i,j,n,*aj, am = A->rmap->n,*ii,arm;
1211:   PetscInt       am2 = 2*am, am3 = 3*am,  bm4 = 4*bm,colam,*ridx;

1214:   if (!cm || !cn) return(0);
1215:   b = bd->v;
1216:   MatDenseGetArray(C,&c);
1217:   b1   = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;

1219:   if (a->compressedrow.use) { /* use compressed row format */
1220:     for (col=0; col<cn-4; col += 4) {  /* over columns of C */
1221:       colam = col*am;
1222:       arm   = a->compressedrow.nrows;
1223:       ii    = a->compressedrow.i;
1224:       ridx  = a->compressedrow.rindex;
1225:       for (i=0; i<arm; i++) {        /* over rows of C in those columns */
1226:         r1 = r2 = r3 = r4 = 0.0;
1227:         n  = ii[i+1] - ii[i];
1228:         aj = a->j + ii[i];
1229:         aa = a->a + ii[i];
1230:         for (j=0; j<n; j++) {
1231:           r1 += (*aa)*b1[*aj];
1232:           r2 += (*aa)*b2[*aj];
1233:           r3 += (*aa)*b3[*aj];
1234:           r4 += (*aa++)*b4[*aj++];
1235:         }
1236:         c[colam       + ridx[i]] += r1;
1237:         c[colam + am  + ridx[i]] += r2;
1238:         c[colam + am2 + ridx[i]] += r3;
1239:         c[colam + am3 + ridx[i]] += r4;
1240:       }
1241:       b1 += bm4;
1242:       b2 += bm4;
1243:       b3 += bm4;
1244:       b4 += bm4;
1245:     }
1246:     for (; col<cn; col++) {     /* over extra columns of C */
1247:       colam = col*am;
1248:       arm   = a->compressedrow.nrows;
1249:       ii    = a->compressedrow.i;
1250:       ridx  = a->compressedrow.rindex;
1251:       for (i=0; i<arm; i++) {  /* over rows of C in those columns */
1252:         r1 = 0.0;
1253:         n  = ii[i+1] - ii[i];
1254:         aj = a->j + ii[i];
1255:         aa = a->a + ii[i];

1257:         for (j=0; j<n; j++) {
1258:           r1 += (*aa++)*b1[*aj++];
1259:         }
1260:         c[colam + ridx[i]] += r1;
1261:       }
1262:       b1 += bm;
1263:     }
1264:   } else {
1265:     for (col=0; col<cn-4; col += 4) {  /* over columns of C */
1266:       colam = col*am;
1267:       for (i=0; i<am; i++) {        /* over rows of C in those columns */
1268:         r1 = r2 = r3 = r4 = 0.0;
1269:         n  = a->i[i+1] - a->i[i];
1270:         aj = a->j + a->i[i];
1271:         aa = a->a + a->i[i];
1272:         for (j=0; j<n; j++) {
1273:           r1 += (*aa)*b1[*aj];
1274:           r2 += (*aa)*b2[*aj];
1275:           r3 += (*aa)*b3[*aj];
1276:           r4 += (*aa++)*b4[*aj++];
1277:         }
1278:         c[colam + i]       += r1;
1279:         c[colam + am + i]  += r2;
1280:         c[colam + am2 + i] += r3;
1281:         c[colam + am3 + i] += r4;
1282:       }
1283:       b1 += bm4;
1284:       b2 += bm4;
1285:       b3 += bm4;
1286:       b4 += bm4;
1287:     }
1288:     for (; col<cn; col++) {     /* over extra columns of C */
1289:       colam = col*am;
1290:       for (i=0; i<am; i++) {  /* over rows of C in those columns */
1291:         r1 = 0.0;
1292:         n  = a->i[i+1] - a->i[i];
1293:         aj = a->j + a->i[i];
1294:         aa = a->a + a->i[i];

1296:         for (j=0; j<n; j++) {
1297:           r1 += (*aa++)*b1[*aj++];
1298:         }
1299:         c[colam + i] += r1;
1300:       }
1301:       b1 += bm;
1302:     }
1303:   }
1304:   PetscLogFlops(cn*2.0*a->nz);
1305:   MatDenseRestoreArray(C,&c);
1306:   return(0);
1307: }

1311: PetscErrorCode  MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring coloring,Mat B,Mat Btdense)
1312: {
1314:   Mat_SeqAIJ     *b       = (Mat_SeqAIJ*)B->data;
1315:   Mat_SeqDense   *btdense = (Mat_SeqDense*)Btdense->data;
1316:   PetscInt       *bi      = b->i,*bj=b->j;
1317:   PetscInt       m        = Btdense->rmap->n,n=Btdense->cmap->n,j,k,l,col,anz,*btcol,brow,ncolumns;
1318:   MatScalar      *btval,*btval_den,*ba=b->a;
1319:   PetscInt       *columns=coloring->columns,*colorforcol=coloring->colorforcol,ncolors=coloring->ncolors;

1322:   btval_den=btdense->v;
1323:   PetscMemzero(btval_den,(m*n)*sizeof(MatScalar));
1324:   for (k=0; k<ncolors; k++) {
1325:     ncolumns = coloring->ncolumns[k];
1326:     for (l=0; l<ncolumns; l++) { /* insert a row of B to a column of Btdense */
1327:       col   = *(columns + colorforcol[k] + l);
1328:       btcol = bj + bi[col];
1329:       btval = ba + bi[col];
1330:       anz   = bi[col+1] - bi[col];
1331:       for (j=0; j<anz; j++) {
1332:         brow            = btcol[j];
1333:         btval_den[brow] = btval[j];
1334:       }
1335:     }
1336:     btval_den += m;
1337:   }
1338:   return(0);
1339: }

1343: PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring matcoloring,Mat Cden,Mat Csp)
1344: {
1346:   Mat_SeqAIJ     *csp = (Mat_SeqAIJ*)Csp->data;
1347:   PetscScalar    *ca_den,*ca_den_ptr,*ca=csp->a;
1348:   PetscInt       k,l,m=Cden->rmap->n,ncolors=matcoloring->ncolors;
1349:   PetscInt       brows=matcoloring->brows,*den2sp=matcoloring->den2sp;
1350:   PetscInt       nrows,*row,*idx;
1351:   PetscInt       *rows=matcoloring->rows,*colorforrow=matcoloring->colorforrow;

1354:   MatDenseGetArray(Cden,&ca_den);

1356:   if (brows > 0) {
1357:     PetscInt *lstart,row_end,row_start;
1358:     lstart = matcoloring->lstart;
1359:     PetscMemzero(lstart,ncolors*sizeof(PetscInt));

1361:     row_end = brows;
1362:     if (row_end > m) row_end = m;
1363:     for (row_start=0; row_start<m; row_start+=brows) { /* loop over row blocks of Csp */
1364:       ca_den_ptr = ca_den;
1365:       for (k=0; k<ncolors; k++) { /* loop over colors (columns of Cden) */
1366:         nrows = matcoloring->nrows[k];
1367:         row   = rows  + colorforrow[k];
1368:         idx   = den2sp + colorforrow[k];
1369:         for (l=lstart[k]; l<nrows; l++) {
1370:           if (row[l] >= row_end) {
1371:             lstart[k] = l;
1372:             break;
1373:           } else {
1374:             ca[idx[l]] = ca_den_ptr[row[l]];
1375:           }
1376:         }
1377:         ca_den_ptr += m;
1378:       }
1379:       row_end += brows;
1380:       if (row_end > m) row_end = m;
1381:     }
1382:   } else { /* non-blocked impl: loop over columns of Csp - slow if Csp is large */
1383:     ca_den_ptr = ca_den;
1384:     for (k=0; k<ncolors; k++) {
1385:       nrows = matcoloring->nrows[k];
1386:       row   = rows  + colorforrow[k];
1387:       idx   = den2sp + colorforrow[k];
1388:       for (l=0; l<nrows; l++) {
1389:         ca[idx[l]] = ca_den_ptr[row[l]];
1390:       }
1391:       ca_den_ptr += m;
1392:     }
1393:   }

1395:   MatDenseRestoreArray(Cden,&ca_den);
1396: #if defined(PETSC_USE_INFO)
1397:   if (matcoloring->brows > 0) {
1398:     PetscInfo1(Csp,"Loop over %D row blocks for den2sp\n",brows);
1399:   } else {
1400:     PetscInfo(Csp,"Loop over colors/columns of Cden, inefficient for large sparse matrix product \n");
1401:   }
1402: #endif
1403:   return(0);
1404: }

1408: PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat mat,ISColoring iscoloring,MatTransposeColoring c)
1409: {
1411:   PetscInt       i,n,nrows,Nbs,j,k,m,ncols,col,cm;
1412:   const PetscInt *is,*ci,*cj,*row_idx;
1413:   PetscInt       nis = iscoloring->n,*rowhit,bs = 1;
1414:   IS             *isa;
1415:   Mat_SeqAIJ     *csp = (Mat_SeqAIJ*)mat->data;
1416:   PetscInt       *colorforrow,*rows,*rows_i,*idxhit,*spidx,*den2sp,*den2sp_i;
1417:   PetscInt       *colorforcol,*columns,*columns_i,brows;
1418:   PetscBool      flg;

1421:   ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);

1423:   /* bs >1 is not being tested yet! */
1424:   Nbs       = mat->cmap->N/bs;
1425:   c->M      = mat->rmap->N/bs;  /* set total rows, columns and local rows */
1426:   c->N      = Nbs;
1427:   c->m      = c->M;
1428:   c->rstart = 0;
1429:   c->brows  = 100;

1431:   c->ncolors = nis;
1432:   PetscMalloc3(nis,&c->ncolumns,nis,&c->nrows,nis+1,&colorforrow);
1433:   PetscMalloc1(csp->nz+1,&rows);
1434:   PetscMalloc1(csp->nz+1,&den2sp);

1436:   brows = c->brows;
1437:   PetscOptionsGetInt(NULL,NULL,"-matden2sp_brows",&brows,&flg);
1438:   if (flg) c->brows = brows;
1439:   if (brows > 0) {
1440:     PetscMalloc1(nis+1,&c->lstart);
1441:   }

1443:   colorforrow[0] = 0;
1444:   rows_i         = rows;
1445:   den2sp_i       = den2sp;

1447:   PetscMalloc1(nis+1,&colorforcol);
1448:   PetscMalloc1(Nbs+1,&columns);

1450:   colorforcol[0] = 0;
1451:   columns_i      = columns;

1453:   /* get column-wise storage of mat */
1454:   MatGetColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);

1456:   cm   = c->m;
1457:   PetscMalloc1(cm+1,&rowhit);
1458:   PetscMalloc1(cm+1,&idxhit);
1459:   for (i=0; i<nis; i++) { /* loop over color */
1460:     ISGetLocalSize(isa[i],&n);
1461:     ISGetIndices(isa[i],&is);

1463:     c->ncolumns[i] = n;
1464:     if (n) {
1465:       PetscMemcpy(columns_i,is,n*sizeof(PetscInt));
1466:     }
1467:     colorforcol[i+1] = colorforcol[i] + n;
1468:     columns_i       += n;

1470:     /* fast, crude version requires O(N*N) work */
1471:     PetscMemzero(rowhit,cm*sizeof(PetscInt));

1473:     for (j=0; j<n; j++) { /* loop over columns*/
1474:       col     = is[j];
1475:       row_idx = cj + ci[col];
1476:       m       = ci[col+1] - ci[col];
1477:       for (k=0; k<m; k++) { /* loop over columns marking them in rowhit */
1478:         idxhit[*row_idx]   = spidx[ci[col] + k];
1479:         rowhit[*row_idx++] = col + 1;
1480:       }
1481:     }
1482:     /* count the number of hits */
1483:     nrows = 0;
1484:     for (j=0; j<cm; j++) {
1485:       if (rowhit[j]) nrows++;
1486:     }
1487:     c->nrows[i]      = nrows;
1488:     colorforrow[i+1] = colorforrow[i] + nrows;

1490:     nrows = 0;
1491:     for (j=0; j<cm; j++) { /* loop over rows */
1492:       if (rowhit[j]) {
1493:         rows_i[nrows]   = j;
1494:         den2sp_i[nrows] = idxhit[j];
1495:         nrows++;
1496:       }
1497:     }
1498:     den2sp_i += nrows;

1500:     ISRestoreIndices(isa[i],&is);
1501:     rows_i += nrows;
1502:   }
1503:   MatRestoreColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);
1504:   PetscFree(rowhit);
1505:   ISColoringRestoreIS(iscoloring,&isa);
1506:   if (csp->nz != colorforrow[nis]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"csp->nz %d != colorforrow[nis] %d",csp->nz,colorforrow[nis]);

1508:   c->colorforrow = colorforrow;
1509:   c->rows        = rows;
1510:   c->den2sp      = den2sp;
1511:   c->colorforcol = colorforcol;
1512:   c->columns     = columns;

1514:   PetscFree(idxhit);
1515:   return(0);
1516: }