Actual source code: matmatmult.c

petsc-master 2020-11-24
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>
  8: #include <../src/mat/utils/freespace.h>
  9: #include <petscbt.h>
 10: #include <petsc/private/isimpl.h>
 11: #include <../src/mat/impls/dense/seq/dense.h>

 13: PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
 14: {

 18:   if (C->ops->matmultnumeric) {
 19:     if (C->ops->matmultnumeric == MatMatMultNumeric_SeqAIJ_SeqAIJ) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Recursive call");
 20:     (*C->ops->matmultnumeric)(A,B,C);
 21:   } else {
 22:     MatMatMultNumeric_SeqAIJ_SeqAIJ_Sorted(A,B,C);
 23:   }
 24:   return(0);
 25: }

 27: /* Modified from MatCreateSeqAIJWithArrays() */
 28: PETSC_INTERN PetscErrorCode MatSetSeqAIJWithArrays_private(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt i[],PetscInt j[],PetscScalar a[],MatType mtype,Mat mat)
 29: {
 31:   PetscInt       ii;
 32:   Mat_SeqAIJ     *aij;
 33:   PetscBool      isseqaij;

 36:   if (m > 0 && i[0]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"i (row indices) must start with 0");
 37:   MatSetSizes(mat,m,n,m,n);

 39:   if (!mtype) {
 40:     PetscObjectBaseTypeCompare((PetscObject)mat,MATSEQAIJ,&isseqaij);
 41:     if (!isseqaij) { MatSetType(mat,MATSEQAIJ); }
 42:   } else {
 43:     MatSetType(mat,mtype);
 44:   }
 45:   MatSeqAIJSetPreallocation_SeqAIJ(mat,MAT_SKIP_ALLOCATION,NULL);
 46:   aij  = (Mat_SeqAIJ*)(mat)->data;
 47:   PetscMalloc1(m,&aij->imax);
 48:   PetscMalloc1(m,&aij->ilen);

 50:   aij->i            = i;
 51:   aij->j            = j;
 52:   aij->a            = a;
 53:   aij->singlemalloc = PETSC_FALSE;
 54:   aij->nonew        = -1; /*this indicates that inserting a new value in the matrix that generates a new nonzero is an error*/
 55:   aij->free_a       = PETSC_FALSE;
 56:   aij->free_ij      = PETSC_FALSE;

 58:   for (ii=0; ii<m; ii++) {
 59:     aij->ilen[ii] = aij->imax[ii] = i[ii+1] - i[ii];
 60:   }

 62:   return(0);
 63: }

 65: PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat C)
 66: {
 67:   PetscErrorCode      ierr;
 68:   Mat_Product         *product = C->product;
 69:   MatProductAlgorithm alg;
 70:   PetscBool           flg;

 73:   if (product) {
 74:     alg = product->alg;
 75:   } else {
 76:     alg = "sorted";
 77:   }
 78:   /* sorted */
 79:   PetscStrcmp(alg,"sorted",&flg);
 80:   if (flg) {
 81:     MatMatMultSymbolic_SeqAIJ_SeqAIJ_Sorted(A,B,fill,C);
 82:     return(0);
 83:   }

 85:   /* scalable */
 86:   PetscStrcmp(alg,"scalable",&flg);
 87:   if (flg) {
 88:     MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(A,B,fill,C);
 89:     return(0);
 90:   }

 92:   /* scalable_fast */
 93:   PetscStrcmp(alg,"scalable_fast",&flg);
 94:   if (flg) {
 95:     MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(A,B,fill,C);
 96:     return(0);
 97:   }

 99:   /* heap */
100:   PetscStrcmp(alg,"heap",&flg);
101:   if (flg) {
102:     MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(A,B,fill,C);
103:     return(0);
104:   }

106:   /* btheap */
107:   PetscStrcmp(alg,"btheap",&flg);
108:   if (flg) {
109:     MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(A,B,fill,C);
110:     return(0);
111:   }

113:   /* llcondensed */
114:   PetscStrcmp(alg,"llcondensed",&flg);
115:   if (flg) {
116:     MatMatMultSymbolic_SeqAIJ_SeqAIJ_LLCondensed(A,B,fill,C);
117:     return(0);
118:   }

120:   /* rowmerge */
121:   PetscStrcmp(alg,"rowmerge",&flg);
122:   if (flg) {
123:     MatMatMultSymbolic_SeqAIJ_SeqAIJ_RowMerge(A,B,fill,C);
124:     return(0);
125:   }

127: #if defined(PETSC_HAVE_HYPRE)
128:   PetscStrcmp(alg,"hypre",&flg);
129:   if (flg) {
130:     MatMatMultSymbolic_AIJ_AIJ_wHYPRE(A,B,fill,C);
131:     return(0);
132:   }
133: #endif

135:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Mat Product Algorithm is not supported");
136:   return(0);
137: }

139: PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_LLCondensed(Mat A,Mat B,PetscReal fill,Mat C)
140: {
141:   PetscErrorCode     ierr;
142:   Mat_SeqAIJ         *a =(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
143:   PetscInt           *ai=a->i,*bi=b->i,*ci,*cj;
144:   PetscInt           am =A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
145:   PetscReal          afill;
146:   PetscInt           i,j,anzi,brow,bnzj,cnzi,*bj,*aj,*lnk,ndouble=0,Crmax;
147:   PetscTable         ta;
148:   PetscBT            lnkbt;
149:   PetscFreeSpaceList free_space=NULL,current_space=NULL;

152:   /* Get ci and cj */
153:   /*---------------*/
154:   /* Allocate ci array, arrays for fill computation and */
155:   /* free space for accumulating nonzero column info */
156:   PetscMalloc1(am+2,&ci);
157:   ci[0] = 0;

159:   /* create and initialize a linked list */
160:   PetscTableCreate(bn,bn,&ta);
161:   MatRowMergeMax_SeqAIJ(b,bm,ta);
162:   PetscTableGetCount(ta,&Crmax);
163:   PetscTableDestroy(&ta);

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

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

170:   current_space = free_space;

172:   /* Determine ci and cj */
173:   for (i=0; i<am; i++) {
174:     anzi = ai[i+1] - ai[i];
175:     aj   = a->j + ai[i];
176:     for (j=0; j<anzi; j++) {
177:       brow = aj[j];
178:       bnzj = bi[brow+1] - bi[brow];
179:       bj   = b->j + bi[brow];
180:       /* add non-zero cols of B into the sorted linked list lnk */
181:       PetscLLCondensedAddSorted(bnzj,bj,lnk,lnkbt);
182:     }
183:     cnzi = lnk[0];

185:     /* If free space is not available, make more free space */
186:     /* Double the amount of total space in the list */
187:     if (current_space->local_remaining<cnzi) {
188:       PetscFreeSpaceGet(PetscIntSumTruncate(cnzi,current_space->total_array_size),&current_space);
189:       ndouble++;
190:     }

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

195:     current_space->array           += cnzi;
196:     current_space->local_used      += cnzi;
197:     current_space->local_remaining -= cnzi;

199:     ci[i+1] = ci[i] + cnzi;
200:   }

202:   /* Column indices are in the list of free space */
203:   /* Allocate space for cj, initialize cj, and */
204:   /* destroy list of free space and other temporary array(s) */
205:   PetscMalloc1(ci[am]+1,&cj);
206:   PetscFreeSpaceContiguous(&free_space,cj);
207:   PetscLLCondensedDestroy(lnk,lnkbt);

209:   /* put together the new symbolic matrix */
210:   MatSetSeqAIJWithArrays_private(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,((PetscObject)A)->type_name,C);
211:   MatSetBlockSizesFromMats(C,A,B);

213:   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
214:   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
215:   c          = (Mat_SeqAIJ*)(C->data);
216:   c->free_a  = PETSC_FALSE;
217:   c->free_ij = PETSC_TRUE;
218:   c->nonew   = 0;

220:   /* fast, needs non-scalable O(bn) array 'abdense' */
221:   C->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Sorted;

223:   /* set MatInfo */
224:   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
225:   if (afill < 1.0) afill = 1.0;
226:   c->maxnz                  = ci[am];
227:   c->nz                     = ci[am];
228:   C->info.mallocs           = ndouble;
229:   C->info.fill_ratio_given  = fill;
230:   C->info.fill_ratio_needed = afill;

232: #if defined(PETSC_USE_INFO)
233:   if (ci[am]) {
234:     PetscInfo3(C,"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);
235:     PetscInfo1(C,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);
236:   } else {
237:     PetscInfo(C,"Empty matrix product\n");
238:   }
239: #endif
240:   return(0);
241: }

243: PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Sorted(Mat A,Mat B,Mat C)
244: {
246:   PetscLogDouble flops=0.0;
247:   Mat_SeqAIJ     *a   = (Mat_SeqAIJ*)A->data;
248:   Mat_SeqAIJ     *b   = (Mat_SeqAIJ*)B->data;
249:   Mat_SeqAIJ     *c   = (Mat_SeqAIJ*)C->data;
250:   PetscInt       *ai  =a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j;
251:   PetscInt       am   =A->rmap->n,cm=C->rmap->n;
252:   PetscInt       i,j,k,anzi,bnzi,cnzi,brow;
253:   PetscScalar    *aa=a->a,*ba=b->a,*baj,*ca,valtmp;
254:   PetscScalar    *ab_dense;
255:   PetscContainer cab_dense;

258:   if (!c->a) { /* first call of MatMatMultNumeric_SeqAIJ_SeqAIJ, allocate ca and matmult_abdense */
259:     PetscMalloc1(ci[cm]+1,&ca);
260:     c->a      = ca;
261:     c->free_a = PETSC_TRUE;
262:   } else ca = c->a;

264:   /* TODO this should be done in the symbolic phase */
265:   /* However, this function is so heavily used (sometimes in an hidden way through multnumeric function pointers
266:      that is hard to eradicate) */
267:   PetscObjectQuery((PetscObject)C,"__PETSc__ab_dense",(PetscObject*)&cab_dense);
268:   if (!cab_dense) {
269:     PetscMalloc1(B->cmap->N,&ab_dense);
270:     PetscContainerCreate(PETSC_COMM_SELF,&cab_dense);
271:     PetscContainerSetPointer(cab_dense,ab_dense);
272:     PetscContainerSetUserDestroy(cab_dense,PetscContainerUserDestroyDefault);
273:     PetscObjectCompose((PetscObject)C,"__PETSc__ab_dense",(PetscObject)cab_dense);
274:     PetscObjectDereference((PetscObject)cab_dense);
275:   }
276:   PetscContainerGetPointer(cab_dense,(void**)&ab_dense);
277:   PetscArrayzero(ab_dense,B->cmap->N);

279:   /* clean old values in C */
280:   PetscArrayzero(ca,ci[cm]);
281:   /* Traverse A row-wise. */
282:   /* Build the ith row in C by summing over nonzero columns in A, */
283:   /* the rows of B corresponding to nonzeros of A. */
284:   for (i=0; i<am; i++) {
285:     anzi = ai[i+1] - ai[i];
286:     for (j=0; j<anzi; j++) {
287:       brow = aj[j];
288:       bnzi = bi[brow+1] - bi[brow];
289:       bjj  = bj + bi[brow];
290:       baj  = ba + bi[brow];
291:       /* perform dense axpy */
292:       valtmp = aa[j];
293:       for (k=0; k<bnzi; k++) {
294:         ab_dense[bjj[k]] += valtmp*baj[k];
295:       }
296:       flops += 2*bnzi;
297:     }
298:     aj += anzi; aa += anzi;

300:     cnzi = ci[i+1] - ci[i];
301:     for (k=0; k<cnzi; k++) {
302:       ca[k]          += ab_dense[cj[k]];
303:       ab_dense[cj[k]] = 0.0; /* zero ab_dense */
304:     }
305:     flops += cnzi;
306:     cj    += cnzi; ca += cnzi;
307:   }
308:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
309:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
310:   PetscLogFlops(flops);
311:   return(0);
312: }

314: PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(Mat A,Mat B,Mat C)
315: {
317:   PetscLogDouble flops=0.0;
318:   Mat_SeqAIJ     *a   = (Mat_SeqAIJ*)A->data;
319:   Mat_SeqAIJ     *b   = (Mat_SeqAIJ*)B->data;
320:   Mat_SeqAIJ     *c   = (Mat_SeqAIJ*)C->data;
321:   PetscInt       *ai  = a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j;
322:   PetscInt       am   = A->rmap->N,cm=C->rmap->N;
323:   PetscInt       i,j,k,anzi,bnzi,cnzi,brow;
324:   PetscScalar    *aa=a->a,*ba=b->a,*baj,*ca=c->a,valtmp;
325:   PetscInt       nextb;

328:   if (!ca) { /* first call of MatMatMultNumeric_SeqAIJ_SeqAIJ, allocate ca and matmult_abdense */
329:     PetscMalloc1(ci[cm]+1,&ca);
330:     c->a      = ca;
331:     c->free_a = PETSC_TRUE;
332:   }

334:   /* clean old values in C */
335:   PetscArrayzero(ca,ci[cm]);
336:   /* Traverse A row-wise. */
337:   /* Build the ith row in C by summing over nonzero columns in A, */
338:   /* the rows of B corresponding to nonzeros of A. */
339:   for (i=0; i<am; i++) {
340:     anzi = ai[i+1] - ai[i];
341:     cnzi = ci[i+1] - ci[i];
342:     for (j=0; j<anzi; j++) {
343:       brow = aj[j];
344:       bnzi = bi[brow+1] - bi[brow];
345:       bjj  = bj + bi[brow];
346:       baj  = ba + bi[brow];
347:       /* perform sparse axpy */
348:       valtmp = aa[j];
349:       nextb  = 0;
350:       for (k=0; nextb<bnzi; k++) {
351:         if (cj[k] == bjj[nextb]) { /* ccol == bcol */
352:           ca[k] += valtmp*baj[nextb++];
353:         }
354:       }
355:       flops += 2*bnzi;
356:     }
357:     aj += anzi; aa += anzi;
358:     cj += cnzi; ca += cnzi;
359:   }
360:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
361:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
362:   PetscLogFlops(flops);
363:   return(0);
364: }

366: PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(Mat A,Mat B,PetscReal fill,Mat C)
367: {
368:   PetscErrorCode     ierr;
369:   Mat_SeqAIJ         *a  = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
370:   PetscInt           *ai = a->i,*bi=b->i,*ci,*cj;
371:   PetscInt           am  = A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
372:   MatScalar          *ca;
373:   PetscReal          afill;
374:   PetscInt           i,j,anzi,brow,bnzj,cnzi,*bj,*aj,*lnk,ndouble=0,Crmax;
375:   PetscTable         ta;
376:   PetscFreeSpaceList free_space=NULL,current_space=NULL;

379:   /* Get ci and cj - same as MatMatMultSymbolic_SeqAIJ_SeqAIJ except using PetscLLxxx_fast() */
380:   /*-----------------------------------------------------------------------------------------*/
381:   /* Allocate arrays for fill computation and free space for accumulating nonzero column */
382:   PetscMalloc1(am+2,&ci);
383:   ci[0] = 0;

385:   /* create and initialize a linked list */
386:   PetscTableCreate(bn,bn,&ta);
387:   MatRowMergeMax_SeqAIJ(b,bm,ta);
388:   PetscTableGetCount(ta,&Crmax);
389:   PetscTableDestroy(&ta);

391:   PetscLLCondensedCreate_fast(Crmax,&lnk);

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

397:   /* Determine ci and cj */
398:   for (i=0; i<am; i++) {
399:     anzi = ai[i+1] - ai[i];
400:     aj   = a->j + ai[i];
401:     for (j=0; j<anzi; j++) {
402:       brow = aj[j];
403:       bnzj = bi[brow+1] - bi[brow];
404:       bj   = b->j + bi[brow];
405:       /* add non-zero cols of B into the sorted linked list lnk */
406:       PetscLLCondensedAddSorted_fast(bnzj,bj,lnk);
407:     }
408:     cnzi = lnk[1];

410:     /* If free space is not available, make more free space */
411:     /* Double the amount of total space in the list */
412:     if (current_space->local_remaining<cnzi) {
413:       PetscFreeSpaceGet(PetscIntSumTruncate(cnzi,current_space->total_array_size),&current_space);
414:       ndouble++;
415:     }

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

420:     current_space->array           += cnzi;
421:     current_space->local_used      += cnzi;
422:     current_space->local_remaining -= cnzi;

424:     ci[i+1] = ci[i] + cnzi;
425:   }

427:   /* Column indices are in the list of free space */
428:   /* Allocate space for cj, initialize cj, and */
429:   /* destroy list of free space and other temporary array(s) */
430:   PetscMalloc1(ci[am]+1,&cj);
431:   PetscFreeSpaceContiguous(&free_space,cj);
432:   PetscLLCondensedDestroy_fast(lnk);

434:   /* Allocate space for ca */
435:   PetscCalloc1(ci[am]+1,&ca);

437:   /* put together the new symbolic matrix */
438:   MatSetSeqAIJWithArrays_private(PetscObjectComm((PetscObject)A),am,bn,ci,cj,ca,((PetscObject)A)->type_name,C);
439:   MatSetBlockSizesFromMats(C,A,B);

441:   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
442:   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
443:   c          = (Mat_SeqAIJ*)(C->data);
444:   c->free_a  = PETSC_TRUE;
445:   c->free_ij = PETSC_TRUE;
446:   c->nonew   = 0;

448:   /* slower, less memory */
449:   C->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable;

451:   /* set MatInfo */
452:   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
453:   if (afill < 1.0) afill = 1.0;
454:   c->maxnz                     = ci[am];
455:   c->nz                        = ci[am];
456:   C->info.mallocs           = ndouble;
457:   C->info.fill_ratio_given  = fill;
458:   C->info.fill_ratio_needed = afill;

460: #if defined(PETSC_USE_INFO)
461:   if (ci[am]) {
462:     PetscInfo3(C,"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);
463:     PetscInfo1(C,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);
464:   } else {
465:     PetscInfo(C,"Empty matrix product\n");
466:   }
467: #endif
468:   return(0);
469: }

471: PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(Mat A,Mat B,PetscReal fill,Mat C)
472: {
473:   PetscErrorCode     ierr;
474:   Mat_SeqAIJ         *a  = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
475:   PetscInt           *ai = a->i,*bi=b->i,*ci,*cj;
476:   PetscInt           am  = A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
477:   MatScalar          *ca;
478:   PetscReal          afill;
479:   PetscInt           i,j,anzi,brow,bnzj,cnzi,*bj,*aj,*lnk,ndouble=0,Crmax;
480:   PetscTable         ta;
481:   PetscFreeSpaceList free_space=NULL,current_space=NULL;

484:   /* Get ci and cj - same as MatMatMultSymbolic_SeqAIJ_SeqAIJ except using PetscLLxxx_Scalalbe() */
485:   /*---------------------------------------------------------------------------------------------*/
486:   /* Allocate arrays for fill computation and free space for accumulating nonzero column */
487:   PetscMalloc1(am+2,&ci);
488:   ci[0] = 0;

490:   /* create and initialize a linked list */
491:   PetscTableCreate(bn,bn,&ta);
492:   MatRowMergeMax_SeqAIJ(b,bm,ta);
493:   PetscTableGetCount(ta,&Crmax);
494:   PetscTableDestroy(&ta);
495:   PetscLLCondensedCreate_Scalable(Crmax,&lnk);

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

501:   /* Determine ci and cj */
502:   for (i=0; i<am; i++) {
503:     anzi = ai[i+1] - ai[i];
504:     aj   = a->j + ai[i];
505:     for (j=0; j<anzi; j++) {
506:       brow = aj[j];
507:       bnzj = bi[brow+1] - bi[brow];
508:       bj   = b->j + bi[brow];
509:       /* add non-zero cols of B into the sorted linked list lnk */
510:       PetscLLCondensedAddSorted_Scalable(bnzj,bj,lnk);
511:     }
512:     cnzi = lnk[0];

514:     /* If free space is not available, make more free space */
515:     /* Double the amount of total space in the list */
516:     if (current_space->local_remaining<cnzi) {
517:       PetscFreeSpaceGet(PetscIntSumTruncate(cnzi,current_space->total_array_size),&current_space);
518:       ndouble++;
519:     }

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

524:     current_space->array           += cnzi;
525:     current_space->local_used      += cnzi;
526:     current_space->local_remaining -= cnzi;

528:     ci[i+1] = ci[i] + cnzi;
529:   }

531:   /* Column indices are in the list of free space */
532:   /* Allocate space for cj, initialize cj, and */
533:   /* destroy list of free space and other temporary array(s) */
534:   PetscMalloc1(ci[am]+1,&cj);
535:   PetscFreeSpaceContiguous(&free_space,cj);
536:   PetscLLCondensedDestroy_Scalable(lnk);

538:   /* Allocate space for ca */
539:   /*-----------------------*/
540:   PetscCalloc1(ci[am]+1,&ca);

542:   /* put together the new symbolic matrix */
543:   MatSetSeqAIJWithArrays_private(PetscObjectComm((PetscObject)A),am,bn,ci,cj,ca,((PetscObject)A)->type_name,C);
544:   MatSetBlockSizesFromMats(C,A,B);

546:   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
547:   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
548:   c          = (Mat_SeqAIJ*)(C->data);
549:   c->free_a  = PETSC_TRUE;
550:   c->free_ij = PETSC_TRUE;
551:   c->nonew   = 0;

553:   /* slower, less memory */
554:   C->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable;

556:   /* set MatInfo */
557:   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
558:   if (afill < 1.0) afill = 1.0;
559:   c->maxnz                     = ci[am];
560:   c->nz                        = ci[am];
561:   C->info.mallocs           = ndouble;
562:   C->info.fill_ratio_given  = fill;
563:   C->info.fill_ratio_needed = afill;

565: #if defined(PETSC_USE_INFO)
566:   if (ci[am]) {
567:     PetscInfo3(C,"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);
568:     PetscInfo1(C,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);
569:   } else {
570:     PetscInfo(C,"Empty matrix product\n");
571:   }
572: #endif
573:   return(0);
574: }

576: PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(Mat A,Mat B,PetscReal fill,Mat C)
577: {
578:   PetscErrorCode     ierr;
579:   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
580:   const PetscInt     *ai=a->i,*bi=b->i,*aj=a->j,*bj=b->j;
581:   PetscInt           *ci,*cj,*bb;
582:   PetscInt           am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
583:   PetscReal          afill;
584:   PetscInt           i,j,col,ndouble = 0;
585:   PetscFreeSpaceList free_space=NULL,current_space=NULL;
586:   PetscHeap          h;

589:   /* Get ci and cj - by merging sorted rows using a heap */
590:   /*---------------------------------------------------------------------------------------------*/
591:   /* Allocate arrays for fill computation and free space for accumulating nonzero column */
592:   PetscMalloc1(am+2,&ci);
593:   ci[0] = 0;

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

599:   PetscHeapCreate(a->rmax,&h);
600:   PetscMalloc1(a->rmax,&bb);

602:   /* Determine ci and cj */
603:   for (i=0; i<am; i++) {
604:     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 */
605:     const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */
606:     ci[i+1] = ci[i];
607:     /* Populate the min heap */
608:     for (j=0; j<anzi; j++) {
609:       bb[j] = bi[acol[j]];         /* bb points at the start of the row */
610:       if (bb[j] < bi[acol[j]+1]) { /* Add if row is nonempty */
611:         PetscHeapAdd(h,j,bj[bb[j]++]);
612:       }
613:     }
614:     /* Pick off the min element, adding it to free space */
615:     PetscHeapPop(h,&j,&col);
616:     while (j >= 0) {
617:       if (current_space->local_remaining < 1) { /* double the size, but don't exceed 16 MiB */
618:         PetscFreeSpaceGet(PetscMin(PetscIntMultTruncate(2,current_space->total_array_size),16 << 20),&current_space);
619:         ndouble++;
620:       }
621:       *(current_space->array++) = col;
622:       current_space->local_used++;
623:       current_space->local_remaining--;
624:       ci[i+1]++;

626:       /* stash if anything else remains in this row of B */
627:       if (bb[j] < bi[acol[j]+1]) {PetscHeapStash(h,j,bj[bb[j]++]);}
628:       while (1) {               /* pop and stash any other rows of B that also had an entry in this column */
629:         PetscInt j2,col2;
630:         PetscHeapPeek(h,&j2,&col2);
631:         if (col2 != col) break;
632:         PetscHeapPop(h,&j2,&col2);
633:         if (bb[j2] < bi[acol[j2]+1]) {PetscHeapStash(h,j2,bj[bb[j2]++]);}
634:       }
635:       /* Put any stashed elements back into the min heap */
636:       PetscHeapUnstash(h);
637:       PetscHeapPop(h,&j,&col);
638:     }
639:   }
640:   PetscFree(bb);
641:   PetscHeapDestroy(&h);

643:   /* Column indices are in the list of free space */
644:   /* Allocate space for cj, initialize cj, and */
645:   /* destroy list of free space and other temporary array(s) */
646:   PetscMalloc1(ci[am],&cj);
647:   PetscFreeSpaceContiguous(&free_space,cj);

649:   /* put together the new symbolic matrix */
650:   MatSetSeqAIJWithArrays_private(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,((PetscObject)A)->type_name,C);
651:   MatSetBlockSizesFromMats(C,A,B);

653:   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
654:   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
655:   c          = (Mat_SeqAIJ*)(C->data);
656:   c->free_a  = PETSC_TRUE;
657:   c->free_ij = PETSC_TRUE;
658:   c->nonew   = 0;

660:   C->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Sorted;

662:   /* set MatInfo */
663:   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
664:   if (afill < 1.0) afill = 1.0;
665:   c->maxnz                     = ci[am];
666:   c->nz                        = ci[am];
667:   C->info.mallocs           = ndouble;
668:   C->info.fill_ratio_given  = fill;
669:   C->info.fill_ratio_needed = afill;

671: #if defined(PETSC_USE_INFO)
672:   if (ci[am]) {
673:     PetscInfo3(C,"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);
674:     PetscInfo1(C,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);
675:   } else {
676:     PetscInfo(C,"Empty matrix product\n");
677:   }
678: #endif
679:   return(0);
680: }

682: PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(Mat A,Mat B,PetscReal fill,Mat C)
683: {
684:   PetscErrorCode     ierr;
685:   Mat_SeqAIJ         *a  = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
686:   const PetscInt     *ai = a->i,*bi=b->i,*aj=a->j,*bj=b->j;
687:   PetscInt           *ci,*cj,*bb;
688:   PetscInt           am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
689:   PetscReal          afill;
690:   PetscInt           i,j,col,ndouble = 0;
691:   PetscFreeSpaceList free_space=NULL,current_space=NULL;
692:   PetscHeap          h;
693:   PetscBT            bt;

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

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

705:   current_space = free_space;

707:   PetscHeapCreate(a->rmax,&h);
708:   PetscMalloc1(a->rmax,&bb);
709:   PetscBTCreate(bn,&bt);

711:   /* Determine ci and cj */
712:   for (i=0; i<am; i++) {
713:     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 */
714:     const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */
715:     const PetscInt *fptr = current_space->array; /* Save beginning of the row so we can clear the BT later */
716:     ci[i+1] = ci[i];
717:     /* Populate the min heap */
718:     for (j=0; j<anzi; j++) {
719:       PetscInt brow = acol[j];
720:       for (bb[j] = bi[brow]; bb[j] < bi[brow+1]; bb[j]++) {
721:         PetscInt bcol = bj[bb[j]];
722:         if (!PetscBTLookupSet(bt,bcol)) { /* new entry */
723:           PetscHeapAdd(h,j,bcol);
724:           bb[j]++;
725:           break;
726:         }
727:       }
728:     }
729:     /* Pick off the min element, adding it to free space */
730:     PetscHeapPop(h,&j,&col);
731:     while (j >= 0) {
732:       if (current_space->local_remaining < 1) { /* double the size, but don't exceed 16 MiB */
733:         fptr = NULL;                      /* need PetscBTMemzero */
734:         PetscFreeSpaceGet(PetscMin(PetscIntMultTruncate(2,current_space->total_array_size),16 << 20),&current_space);
735:         ndouble++;
736:       }
737:       *(current_space->array++) = col;
738:       current_space->local_used++;
739:       current_space->local_remaining--;
740:       ci[i+1]++;

742:       /* stash if anything else remains in this row of B */
743:       for (; bb[j] < bi[acol[j]+1]; bb[j]++) {
744:         PetscInt bcol = bj[bb[j]];
745:         if (!PetscBTLookupSet(bt,bcol)) { /* new entry */
746:           PetscHeapAdd(h,j,bcol);
747:           bb[j]++;
748:           break;
749:         }
750:       }
751:       PetscHeapPop(h,&j,&col);
752:     }
753:     if (fptr) {                 /* Clear the bits for this row */
754:       for (; fptr<current_space->array; fptr++) {PetscBTClear(bt,*fptr);}
755:     } else {                    /* We reallocated so we don't remember (easily) how to clear only the bits we changed */
756:       PetscBTMemzero(bn,bt);
757:     }
758:   }
759:   PetscFree(bb);
760:   PetscHeapDestroy(&h);
761:   PetscBTDestroy(&bt);

763:   /* Column indices are in the list of free space */
764:   /* Allocate space for cj, initialize cj, and */
765:   /* destroy list of free space and other temporary array(s) */
766:   PetscMalloc1(ci[am],&cj);
767:   PetscFreeSpaceContiguous(&free_space,cj);

769:   /* put together the new symbolic matrix */
770:   MatSetSeqAIJWithArrays_private(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,((PetscObject)A)->type_name,C);
771:   MatSetBlockSizesFromMats(C,A,B);

773:   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
774:   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
775:   c          = (Mat_SeqAIJ*)(C->data);
776:   c->free_a  = PETSC_TRUE;
777:   c->free_ij = PETSC_TRUE;
778:   c->nonew   = 0;

780:   C->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Sorted;

782:   /* set MatInfo */
783:   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
784:   if (afill < 1.0) afill = 1.0;
785:   c->maxnz                     = ci[am];
786:   c->nz                        = ci[am];
787:   C->info.mallocs           = ndouble;
788:   C->info.fill_ratio_given  = fill;
789:   C->info.fill_ratio_needed = afill;

791: #if defined(PETSC_USE_INFO)
792:   if (ci[am]) {
793:     PetscInfo3(C,"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);
794:     PetscInfo1(C,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);
795:   } else {
796:     PetscInfo(C,"Empty matrix product\n");
797:   }
798: #endif
799:   return(0);
800: }


803: PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_RowMerge(Mat A,Mat B,PetscReal fill,Mat C)
804: {
805:   PetscErrorCode     ierr;
806:   Mat_SeqAIJ         *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
807:   const PetscInt     *ai=a->i,*bi=b->i,*aj=a->j,*bj=b->j,*inputi,*inputj,*inputcol,*inputcol_L1;
808:   PetscInt           *ci,*cj,*outputj,worki_L1[9],worki_L2[9];
809:   PetscInt           c_maxmem,a_maxrownnz=0,a_rownnz;
810:   const PetscInt     workcol[8]={0,1,2,3,4,5,6,7};
811:   const PetscInt     am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
812:   const PetscInt     *brow_ptr[8],*brow_end[8];
813:   PetscInt           window[8];
814:   PetscInt           window_min,old_window_min,ci_nnz,outputi_nnz=0,L1_nrows,L2_nrows;
815:   PetscInt           i,k,ndouble=0,L1_rowsleft,rowsleft;
816:   PetscReal          afill;
817:   PetscInt           *workj_L1,*workj_L2,*workj_L3;
818:   PetscInt           L1_nnz,L2_nnz;

820:   /* Step 1: Get upper bound on memory required for allocation.
821:              Because of the way virtual memory works,
822:              only the memory pages that are actually needed will be physically allocated. */
824:   PetscMalloc1(am+1,&ci);
825:   for (i=0; i<am; i++) {
826:     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 */
827:     const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */
828:     a_rownnz = 0;
829:     for (k=0; k<anzi; ++k) {
830:       a_rownnz += bi[acol[k]+1] - bi[acol[k]];
831:       if (a_rownnz > bn) {
832:         a_rownnz = bn;
833:         break;
834:       }
835:     }
836:     a_maxrownnz = PetscMax(a_maxrownnz, a_rownnz);
837:   }
838:   /* temporary work areas for merging rows */
839:   PetscMalloc1(a_maxrownnz*8,&workj_L1);
840:   PetscMalloc1(a_maxrownnz*8,&workj_L2);
841:   PetscMalloc1(a_maxrownnz,&workj_L3);

843:   /* This should be enough for almost all matrices. If not, memory is reallocated later. */
844:   c_maxmem = 8*(ai[am]+bi[bm]);
845:   /* Step 2: Populate pattern for C */
846:   PetscMalloc1(c_maxmem,&cj);

848:   ci_nnz       = 0;
849:   ci[0]        = 0;
850:   worki_L1[0]  = 0;
851:   worki_L2[0]  = 0;
852:   for (i=0; i<am; i++) {
853:     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 */
854:     const PetscInt *acol = aj + ai[i];      /* column indices of nonzero entries in this row */
855:     rowsleft             = anzi;
856:     inputcol_L1          = acol;
857:     L2_nnz               = 0;
858:     L2_nrows             = 1;  /* Number of rows to be merged on Level 3. output of L3 already exists -> initial value 1   */
859:     worki_L2[1]          = 0;
860:     outputi_nnz          = 0;

862:     /* If the number of indices in C so far + the max number of columns in the next row > c_maxmem  -> allocate more memory */
863:     while (ci_nnz+a_maxrownnz > c_maxmem) {
864:       c_maxmem *= 2;
865:       ndouble++;
866:       PetscRealloc(sizeof(PetscInt)*c_maxmem,&cj);
867:     }

869:     while (rowsleft) {
870:       L1_rowsleft = PetscMin(64, rowsleft); /* In the inner loop max 64 rows of B can be merged */
871:       L1_nrows    = 0;
872:       L1_nnz      = 0;
873:       inputcol    = inputcol_L1;
874:       inputi      = bi;
875:       inputj      = bj;

877:       /* The following macro is used to specialize for small rows in A.
878:          This helps with compiler unrolling, improving performance substantially.
879:           Input:  inputj   inputi  inputcol  bn
880:           Output: outputj  outputi_nnz                       */
881:        #define MatMatMultSymbolic_RowMergeMacro(ANNZ)                        \
882:          window_min  = bn;                                                   \
883:          outputi_nnz = 0;                                                    \
884:          for (k=0; k<ANNZ; ++k) {                                            \
885:            brow_ptr[k] = inputj + inputi[inputcol[k]];                       \
886:            brow_end[k] = inputj + inputi[inputcol[k]+1];                     \
887:            window[k]   = (brow_ptr[k] != brow_end[k]) ? *brow_ptr[k] : bn;   \
888:            window_min  = PetscMin(window[k], window_min);                    \
889:          }                                                                   \
890:          while (window_min < bn) {                                           \
891:            outputj[outputi_nnz++] = window_min;                              \
892:            /* advance front and compute new minimum */                       \
893:            old_window_min = window_min;                                      \
894:            window_min = bn;                                                  \
895:            for (k=0; k<ANNZ; ++k) {                                          \
896:              if (window[k] == old_window_min) {                              \
897:                brow_ptr[k]++;                                                \
898:                window[k] = (brow_ptr[k] != brow_end[k]) ? *brow_ptr[k] : bn; \
899:              }                                                               \
900:              window_min = PetscMin(window[k], window_min);                   \
901:            }                                                                 \
902:          }

904:       /************** L E V E L  1 ***************/
905:       /* Merge up to 8 rows of B to L1 work array*/
906:       while (L1_rowsleft) {
907:         outputi_nnz = 0;
908:         if (anzi > 8)  outputj = workj_L1 + L1_nnz;     /* Level 1 rowmerge*/
909:         else           outputj = cj + ci_nnz;           /* Merge directly to C */

911:         switch (L1_rowsleft) {
912:         case 1:  brow_ptr[0] = inputj + inputi[inputcol[0]];
913:                  brow_end[0] = inputj + inputi[inputcol[0]+1];
914:                  for (; brow_ptr[0] != brow_end[0]; ++brow_ptr[0]) outputj[outputi_nnz++] = *brow_ptr[0]; /* copy row in b over */
915:                  inputcol    += L1_rowsleft;
916:                  rowsleft    -= L1_rowsleft;
917:                  L1_rowsleft  = 0;
918:                  break;
919:         case 2:  MatMatMultSymbolic_RowMergeMacro(2);
920:                  inputcol    += L1_rowsleft;
921:                  rowsleft    -= L1_rowsleft;
922:                  L1_rowsleft  = 0;
923:                  break;
924:         case 3: MatMatMultSymbolic_RowMergeMacro(3);
925:                  inputcol    += L1_rowsleft;
926:                  rowsleft    -= L1_rowsleft;
927:                  L1_rowsleft  = 0;
928:                  break;
929:         case 4:  MatMatMultSymbolic_RowMergeMacro(4);
930:                  inputcol    += L1_rowsleft;
931:                  rowsleft    -= L1_rowsleft;
932:                  L1_rowsleft  = 0;
933:                  break;
934:         case 5:  MatMatMultSymbolic_RowMergeMacro(5);
935:                  inputcol    += L1_rowsleft;
936:                  rowsleft    -= L1_rowsleft;
937:                  L1_rowsleft  = 0;
938:                  break;
939:         case 6:  MatMatMultSymbolic_RowMergeMacro(6);
940:                  inputcol    += L1_rowsleft;
941:                  rowsleft    -= L1_rowsleft;
942:                  L1_rowsleft  = 0;
943:                  break;
944:         case 7:  MatMatMultSymbolic_RowMergeMacro(7);
945:                  inputcol    += L1_rowsleft;
946:                  rowsleft    -= L1_rowsleft;
947:                  L1_rowsleft  = 0;
948:                  break;
949:         default: MatMatMultSymbolic_RowMergeMacro(8);
950:                  inputcol    += 8;
951:                  rowsleft    -= 8;
952:                  L1_rowsleft -= 8;
953:                  break;
954:         }
955:         inputcol_L1           = inputcol;
956:         L1_nnz               += outputi_nnz;
957:         worki_L1[++L1_nrows]  = L1_nnz;
958:       }

960:       /********************** L E V E L  2 ************************/
961:       /* Merge from L1 work array to either C or to L2 work array */
962:       if (anzi > 8) {
963:         inputi      = worki_L1;
964:         inputj      = workj_L1;
965:         inputcol    = workcol;
966:         outputi_nnz = 0;

968:         if (anzi <= 64) outputj = cj + ci_nnz;        /* Merge from L1 work array to C */
969:         else            outputj = workj_L2 + L2_nnz;  /* Merge from L1 work array to L2 work array */

971:         switch (L1_nrows) {
972:         case 1:  brow_ptr[0] = inputj + inputi[inputcol[0]];
973:                  brow_end[0] = inputj + inputi[inputcol[0]+1];
974:                  for (; brow_ptr[0] != brow_end[0]; ++brow_ptr[0]) outputj[outputi_nnz++] = *brow_ptr[0]; /* copy row in b over */
975:                  break;
976:         case 2:  MatMatMultSymbolic_RowMergeMacro(2); break;
977:         case 3:  MatMatMultSymbolic_RowMergeMacro(3); break;
978:         case 4:  MatMatMultSymbolic_RowMergeMacro(4); break;
979:         case 5:  MatMatMultSymbolic_RowMergeMacro(5); break;
980:         case 6:  MatMatMultSymbolic_RowMergeMacro(6); break;
981:         case 7:  MatMatMultSymbolic_RowMergeMacro(7); break;
982:         case 8:  MatMatMultSymbolic_RowMergeMacro(8); break;
983:         default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatMatMult logic error: Not merging 1-8 rows from L1 work array!");
984:         }
985:         L2_nnz               += outputi_nnz;
986:         worki_L2[++L2_nrows]  = L2_nnz;

988:         /************************ L E V E L  3 **********************/
989:         /* Merge from L2 work array to either C or to L2 work array */
990:         if (anzi > 64 && (L2_nrows == 8 || rowsleft == 0)) {
991:           inputi      = worki_L2;
992:           inputj      = workj_L2;
993:           inputcol    = workcol;
994:           outputi_nnz = 0;
995:           if (rowsleft) outputj = workj_L3;
996:           else          outputj = cj + ci_nnz;
997:           switch (L2_nrows) {
998:           case 1:  brow_ptr[0] = inputj + inputi[inputcol[0]];
999:                    brow_end[0] = inputj + inputi[inputcol[0]+1];
1000:                    for (; brow_ptr[0] != brow_end[0]; ++brow_ptr[0]) outputj[outputi_nnz++] = *brow_ptr[0]; /* copy row in b over */
1001:                    break;
1002:           case 2:  MatMatMultSymbolic_RowMergeMacro(2); break;
1003:           case 3:  MatMatMultSymbolic_RowMergeMacro(3); break;
1004:           case 4:  MatMatMultSymbolic_RowMergeMacro(4); break;
1005:           case 5:  MatMatMultSymbolic_RowMergeMacro(5); break;
1006:           case 6:  MatMatMultSymbolic_RowMergeMacro(6); break;
1007:           case 7:  MatMatMultSymbolic_RowMergeMacro(7); break;
1008:           case 8:  MatMatMultSymbolic_RowMergeMacro(8); break;
1009:           default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatMatMult logic error: Not merging 1-8 rows from L2 work array!");
1010:           }
1011:           L2_nrows    = 1;
1012:           L2_nnz      = outputi_nnz;
1013:           worki_L2[1] = outputi_nnz;
1014:           /* Copy to workj_L2 */
1015:           if (rowsleft) {
1016:             for (k=0; k<outputi_nnz; ++k)  workj_L2[k] = outputj[k];
1017:           }
1018:         }
1019:       }
1020:     }  /* while (rowsleft) */
1021: #undef MatMatMultSymbolic_RowMergeMacro

1023:     /* terminate current row */
1024:     ci_nnz += outputi_nnz;
1025:     ci[i+1] = ci_nnz;
1026:   }

1028:   /* Step 3: Create the new symbolic matrix */
1029:   MatSetSeqAIJWithArrays_private(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,((PetscObject)A)->type_name,C);
1030:   MatSetBlockSizesFromMats(C,A,B);

1032:   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
1033:   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
1034:   c          = (Mat_SeqAIJ*)(C->data);
1035:   c->free_a  = PETSC_TRUE;
1036:   c->free_ij = PETSC_TRUE;
1037:   c->nonew   = 0;

1039:   C->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Sorted;

1041:   /* set MatInfo */
1042:   afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
1043:   if (afill < 1.0) afill = 1.0;
1044:   c->maxnz                     = ci[am];
1045:   c->nz                        = ci[am];
1046:   C->info.mallocs           = ndouble;
1047:   C->info.fill_ratio_given  = fill;
1048:   C->info.fill_ratio_needed = afill;

1050: #if defined(PETSC_USE_INFO)
1051:   if (ci[am]) {
1052:     PetscInfo3(C,"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);
1053:     PetscInfo1(C,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);
1054:   } else {
1055:     PetscInfo(C,"Empty matrix product\n");
1056:   }
1057: #endif

1059:   /* Step 4: Free temporary work areas */
1060:   PetscFree(workj_L1);
1061:   PetscFree(workj_L2);
1062:   PetscFree(workj_L3);
1063:   return(0);
1064: }

1066: /* concatenate unique entries and then sort */
1067: PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Sorted(Mat A,Mat B,PetscReal fill,Mat C)
1068: {
1070:   Mat_SeqAIJ     *a  = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
1071:   const PetscInt *ai = a->i,*bi=b->i,*aj=a->j,*bj=b->j;
1072:   PetscInt       *ci,*cj;
1073:   PetscInt       am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
1074:   PetscReal      afill;
1075:   PetscInt       i,j,ndouble = 0;
1076:   PetscSegBuffer seg,segrow;
1077:   char           *seen;

1080:   PetscMalloc1(am+1,&ci);
1081:   ci[0] = 0;

1083:   /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
1084:   PetscSegBufferCreate(sizeof(PetscInt),(PetscInt)(fill*(ai[am]+bi[bm])),&seg);
1085:   PetscSegBufferCreate(sizeof(PetscInt),100,&segrow);
1086:   PetscCalloc1(bn,&seen);

1088:   /* Determine ci and cj */
1089:   for (i=0; i<am; i++) {
1090:     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 */
1091:     const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */
1092:     PetscInt packlen = 0,*PETSC_RESTRICT crow;
1093:     /* Pack segrow */
1094:     for (j=0; j<anzi; j++) {
1095:       PetscInt brow = acol[j],bjstart = bi[brow],bjend = bi[brow+1],k;
1096:       for (k=bjstart; k<bjend; k++) {
1097:         PetscInt bcol = bj[k];
1098:         if (!seen[bcol]) { /* new entry */
1099:           PetscInt *PETSC_RESTRICT slot;
1100:           PetscSegBufferGetInts(segrow,1,&slot);
1101:           *slot = bcol;
1102:           seen[bcol] = 1;
1103:           packlen++;
1104:         }
1105:       }
1106:     }
1107:     PetscSegBufferGetInts(seg,packlen,&crow);
1108:     PetscSegBufferExtractTo(segrow,crow);
1109:     PetscSortInt(packlen,crow);
1110:     ci[i+1] = ci[i] + packlen;
1111:     for (j=0; j<packlen; j++) seen[crow[j]] = 0;
1112:   }
1113:   PetscSegBufferDestroy(&segrow);
1114:   PetscFree(seen);

1116:   /* Column indices are in the segmented buffer */
1117:   PetscSegBufferExtractAlloc(seg,&cj);
1118:   PetscSegBufferDestroy(&seg);

1120:   /* put together the new symbolic matrix */
1121:   MatSetSeqAIJWithArrays_private(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,((PetscObject)A)->type_name,C);
1122:   MatSetBlockSizesFromMats(C,A,B);

1124:   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
1125:   /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
1126:   c          = (Mat_SeqAIJ*)(C->data);
1127:   c->free_a  = PETSC_TRUE;
1128:   c->free_ij = PETSC_TRUE;
1129:   c->nonew   = 0;

1131:   C->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Sorted;

1133:   /* set MatInfo */
1134:   afill = (PetscReal)ci[am]/PetscMax(ai[am]+bi[bm],1) + 1.e-5;
1135:   if (afill < 1.0) afill = 1.0;
1136:   c->maxnz                  = ci[am];
1137:   c->nz                     = ci[am];
1138:   C->info.mallocs           = ndouble;
1139:   C->info.fill_ratio_given  = fill;
1140:   C->info.fill_ratio_needed = afill;

1142: #if defined(PETSC_USE_INFO)
1143:   if (ci[am]) {
1144:     PetscInfo3(C,"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);
1145:     PetscInfo1(C,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);
1146:   } else {
1147:     PetscInfo(C,"Empty matrix product\n");
1148:   }
1149: #endif
1150:   return(0);
1151: }

1153: PetscErrorCode MatDestroy_SeqAIJ_MatMatMultTrans(void *data)
1154: {
1155:   PetscErrorCode      ierr;
1156:   Mat_MatMatTransMult *abt=(Mat_MatMatTransMult *)data;

1159:   MatTransposeColoringDestroy(&abt->matcoloring);
1160:   MatDestroy(&abt->Bt_den);
1161:   MatDestroy(&abt->ABt_den);
1162:   PetscFree(abt);
1163:   return(0);
1164: }

1166: PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat C)
1167: {
1168:   PetscErrorCode      ierr;
1169:   Mat                 Bt;
1170:   PetscInt            *bti,*btj;
1171:   Mat_MatMatTransMult *abt;
1172:   Mat_Product         *product = C->product;
1173:   char                *alg;

1176:   if (!product) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing product struct");
1177:   if (product->data) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Extra product struct not empty");

1179:   /* create symbolic Bt */
1180:   MatGetSymbolicTranspose_SeqAIJ(B,&bti,&btj);
1181:   MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,B->cmap->n,B->rmap->n,bti,btj,NULL,&Bt);
1182:   MatSetBlockSizes(Bt,PetscAbs(A->cmap->bs),PetscAbs(B->cmap->bs));
1183:   MatSetType(Bt,((PetscObject)A)->type_name);

1185:   /* get symbolic C=A*Bt */
1186:   PetscStrallocpy(product->alg,&alg);
1187:   MatProductSetAlgorithm(C,"sorted"); /* set algorithm for C = A*Bt */
1188:   MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,Bt,fill,C);
1189:   MatProductSetAlgorithm(C,alg); /* resume original algorithm for ABt product */
1190:   PetscFree(alg);

1192:   /* create a supporting struct for reuse intermidiate dense matrices with matcoloring */
1193:   PetscNew(&abt);

1195:   product->data    = abt;
1196:   product->destroy = MatDestroy_SeqAIJ_MatMatMultTrans;

1198:   C->ops->mattransposemultnumeric = MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ;

1200:   abt->usecoloring = PETSC_FALSE;
1201:   PetscStrcmp(product->alg,"color",&abt->usecoloring);
1202:   if (abt->usecoloring) {
1203:     /* Create MatTransposeColoring from symbolic C=A*B^T */
1204:     MatTransposeColoring matcoloring;
1205:     MatColoring          coloring;
1206:     ISColoring           iscoloring;
1207:     Mat                  Bt_dense,C_dense;

1209:     /* inode causes memory problem */
1210:     MatSetOption(C,MAT_USE_INODES,PETSC_FALSE);

1212:     MatColoringCreate(C,&coloring);
1213:     MatColoringSetDistance(coloring,2);
1214:     MatColoringSetType(coloring,MATCOLORINGSL);
1215:     MatColoringSetFromOptions(coloring);
1216:     MatColoringApply(coloring,&iscoloring);
1217:     MatColoringDestroy(&coloring);
1218:     MatTransposeColoringCreate(C,iscoloring,&matcoloring);

1220:     abt->matcoloring = matcoloring;

1222:     ISColoringDestroy(&iscoloring);

1224:     /* Create Bt_dense and C_dense = A*Bt_dense */
1225:     MatCreate(PETSC_COMM_SELF,&Bt_dense);
1226:     MatSetSizes(Bt_dense,A->cmap->n,matcoloring->ncolors,A->cmap->n,matcoloring->ncolors);
1227:     MatSetType(Bt_dense,MATSEQDENSE);
1228:     MatSeqDenseSetPreallocation(Bt_dense,NULL);

1230:     Bt_dense->assembled = PETSC_TRUE;
1231:     abt->Bt_den         = Bt_dense;

1233:     MatCreate(PETSC_COMM_SELF,&C_dense);
1234:     MatSetSizes(C_dense,A->rmap->n,matcoloring->ncolors,A->rmap->n,matcoloring->ncolors);
1235:     MatSetType(C_dense,MATSEQDENSE);
1236:     MatSeqDenseSetPreallocation(C_dense,NULL);

1238:     Bt_dense->assembled = PETSC_TRUE;
1239:     abt->ABt_den  = C_dense;

1241: #if defined(PETSC_USE_INFO)
1242:     {
1243:       Mat_SeqAIJ *c = (Mat_SeqAIJ*)C->data;
1244:       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));
1245:     }
1246: #endif
1247:   }
1248:   /* clean up */
1249:   MatDestroy(&Bt);
1250:   MatRestoreSymbolicTranspose_SeqAIJ(B,&bti,&btj);
1251:   return(0);
1252: }

1254: PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
1255: {
1256:   PetscErrorCode      ierr;
1257:   Mat_SeqAIJ          *a   =(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data;
1258:   PetscInt            *ai  =a->i,*aj=a->j,*bi=b->i,*bj=b->j,anzi,bnzj,nexta,nextb,*acol,*bcol,brow;
1259:   PetscInt            cm   =C->rmap->n,*ci=c->i,*cj=c->j,i,j,cnzi,*ccol;
1260:   PetscLogDouble      flops=0.0;
1261:   MatScalar           *aa  =a->a,*aval,*ba=b->a,*bval,*ca,*cval;
1262:   Mat_MatMatTransMult *abt;
1263:   Mat_Product         *product = C->product;

1266:   if (!product) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing product struct");
1267:   abt = (Mat_MatMatTransMult *)product->data;
1268:   if (!abt) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing product struct");
1269:   /* clear old values in C */
1270:   if (!c->a) {
1271:     PetscCalloc1(ci[cm]+1,&ca);
1272:     c->a      = ca;
1273:     c->free_a = PETSC_TRUE;
1274:   } else {
1275:     ca =  c->a;
1276:     PetscArrayzero(ca,ci[cm]+1);
1277:   }

1279:   if (abt->usecoloring) {
1280:     MatTransposeColoring matcoloring = abt->matcoloring;
1281:     Mat                  Bt_dense,C_dense = abt->ABt_den;

1283:     /* Get Bt_dense by Apply MatTransposeColoring to B */
1284:     Bt_dense = abt->Bt_den;
1285:     MatTransColoringApplySpToDen(matcoloring,B,Bt_dense);

1287:     /* C_dense = A*Bt_dense */
1288:     MatMatMultNumeric_SeqAIJ_SeqDense(A,Bt_dense,C_dense);

1290:     /* Recover C from C_dense */
1291:     MatTransColoringApplyDenToSp(matcoloring,C_dense,C);
1292:     return(0);
1293:   }

1295:   for (i=0; i<cm; i++) {
1296:     anzi = ai[i+1] - ai[i];
1297:     acol = aj + ai[i];
1298:     aval = aa + ai[i];
1299:     cnzi = ci[i+1] - ci[i];
1300:     ccol = cj + ci[i];
1301:     cval = ca + ci[i];
1302:     for (j=0; j<cnzi; j++) {
1303:       brow = ccol[j];
1304:       bnzj = bi[brow+1] - bi[brow];
1305:       bcol = bj + bi[brow];
1306:       bval = ba + bi[brow];

1308:       /* perform sparse inner-product c(i,j)=A[i,:]*B[j,:]^T */
1309:       nexta = 0; nextb = 0;
1310:       while (nexta<anzi && nextb<bnzj) {
1311:         while (nexta < anzi && acol[nexta] < bcol[nextb]) nexta++;
1312:         if (nexta == anzi) break;
1313:         while (nextb < bnzj && acol[nexta] > bcol[nextb]) nextb++;
1314:         if (nextb == bnzj) break;
1315:         if (acol[nexta] == bcol[nextb]) {
1316:           cval[j] += aval[nexta]*bval[nextb];
1317:           nexta++; nextb++;
1318:           flops += 2;
1319:         }
1320:       }
1321:     }
1322:   }
1323:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
1324:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
1325:   PetscLogFlops(flops);
1326:   return(0);
1327: }

1329: PetscErrorCode MatDestroy_SeqAIJ_MatTransMatMult(void *data)
1330: {
1331:   PetscErrorCode      ierr;
1332:   Mat_MatTransMatMult *atb = (Mat_MatTransMatMult*)data;

1335:   MatDestroy(&atb->At);
1336:   if (atb->destroy) {
1337:     (*atb->destroy)(atb->data);
1338:   }
1339:   PetscFree(atb);
1340:   return(0);
1341: }

1343: PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat C)
1344: {
1346:   Mat            At = NULL;
1347:   PetscInt       *ati,*atj;
1348:   Mat_Product    *product = C->product;
1349:   PetscBool      flg,def,square;

1352:   MatCheckProduct(C,4);
1353:   square = (PetscBool)(A == B && A->symmetric && A->symmetric_set);
1354:   /* outerproduct */
1355:   PetscStrcmp(product->alg,"outerproduct",&flg);
1356:   if (flg) {
1357:     /* create symbolic At */
1358:     if (!square) {
1359:       MatGetSymbolicTranspose_SeqAIJ(A,&ati,&atj);
1360:       MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,A->cmap->n,A->rmap->n,ati,atj,NULL,&At);
1361:       MatSetBlockSizes(At,PetscAbs(A->cmap->bs),PetscAbs(B->cmap->bs));
1362:       MatSetType(At,((PetscObject)A)->type_name);
1363:     }
1364:     /* get symbolic C=At*B */
1365:     MatProductSetAlgorithm(C,"sorted");
1366:     MatMatMultSymbolic_SeqAIJ_SeqAIJ(square ? A : At,B,fill,C);

1368:     /* clean up */
1369:     if (!square) {
1370:       MatDestroy(&At);
1371:       MatRestoreSymbolicTranspose_SeqAIJ(A,&ati,&atj);
1372:     }

1374:     C->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ; /* outerproduct */
1375:     MatProductSetAlgorithm(C,"outerproduct");
1376:     return(0);
1377:   }

1379:   /* matmatmult */
1380:   PetscStrcmp(product->alg,"default",&def);
1381:   PetscStrcmp(product->alg,"at*b",&flg);
1382:   if (flg || def) {
1383:     Mat_MatTransMatMult *atb;

1385:     if (product->data) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Extra product struct not empty");
1386:     PetscNew(&atb);
1387:     if (!square) {
1388:       MatTranspose_SeqAIJ(A,MAT_INITIAL_MATRIX,&At);
1389:     }
1390:     MatProductSetAlgorithm(C,"sorted");
1391:     MatMatMultSymbolic_SeqAIJ_SeqAIJ(square ? A : At,B,fill,C);
1392:     MatProductSetAlgorithm(C,"at*b");
1393:     product->data    = atb;
1394:     product->destroy = MatDestroy_SeqAIJ_MatTransMatMult;
1395:     atb->At          = At;
1396:     atb->updateAt    = PETSC_FALSE; /* because At is computed here */

1398:     C->ops->mattransposemultnumeric = NULL; /* see MatProductNumeric_AtB_SeqAIJ_SeqAIJ */
1399:     return(0);
1400:   }

1402:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Mat Product Algorithm is not supported");
1403:   return(0);
1404: }

1406: PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
1407: {
1409:   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data;
1410:   PetscInt       am=A->rmap->n,anzi,*ai=a->i,*aj=a->j,*bi=b->i,*bj,bnzi,nextb;
1411:   PetscInt       cm=C->rmap->n,*ci=c->i,*cj=c->j,crow,*cjj,i,j,k;
1412:   PetscLogDouble flops=0.0;
1413:   MatScalar      *aa=a->a,*ba,*ca,*caj;

1416:   if (!c->a) {
1417:     PetscCalloc1(ci[cm]+1,&ca);

1419:     c->a      = ca;
1420:     c->free_a = PETSC_TRUE;
1421:   } else {
1422:     ca   = c->a;
1423:     PetscArrayzero(ca,ci[cm]);
1424:   }

1426:   /* compute A^T*B using outer product (A^T)[:,i]*B[i,:] */
1427:   for (i=0; i<am; i++) {
1428:     bj   = b->j + bi[i];
1429:     ba   = b->a + bi[i];
1430:     bnzi = bi[i+1] - bi[i];
1431:     anzi = ai[i+1] - ai[i];
1432:     for (j=0; j<anzi; j++) {
1433:       nextb = 0;
1434:       crow  = *aj++;
1435:       cjj   = cj + ci[crow];
1436:       caj   = ca + ci[crow];
1437:       /* perform sparse axpy operation.  Note cjj includes bj. */
1438:       for (k=0; nextb<bnzi; k++) {
1439:         if (cjj[k] == *(bj+nextb)) { /* ccol == bcol */
1440:           caj[k] += (*aa)*(*(ba+nextb));
1441:           nextb++;
1442:         }
1443:       }
1444:       flops += 2*bnzi;
1445:       aa++;
1446:     }
1447:   }

1449:   /* Assemble the final matrix and clean up */
1450:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
1451:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
1452:   PetscLogFlops(flops);
1453:   return(0);
1454: }

1456: PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqDense(Mat A,Mat B,PetscReal fill,Mat C)
1457: {

1461:   MatMatMultSymbolic_SeqDense_SeqDense(A,B,0.0,C);
1462:   C->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqDense;
1463:   return(0);
1464: }

1466: PETSC_INTERN PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat A,Mat B,Mat C,const PetscBool add)
1467: {
1468:   Mat_SeqAIJ        *a=(Mat_SeqAIJ*)A->data;
1469:   Mat_SeqDense      *bd=(Mat_SeqDense*)B->data;
1470:   Mat_SeqDense      *cd=(Mat_SeqDense*)C->data;
1471:   PetscErrorCode    ierr;
1472:   PetscScalar       *c,r1,r2,r3,r4,*c1,*c2,*c3,*c4;
1473:   const PetscScalar *aa,*b,*b1,*b2,*b3,*b4,*av;
1474:   const PetscInt    *aj;
1475:   PetscInt          cm=C->rmap->n,cn=B->cmap->n,bm=bd->lda,am=A->rmap->n;
1476:   PetscInt          clda=cd->lda;
1477:   PetscInt          am4=4*clda,bm4=4*bm,col,i,j,n;

1480:   if (!cm || !cn) return(0);
1481:   MatSeqAIJGetArrayRead(A,&av);
1482:   if (add) {
1483:     MatDenseGetArray(C,&c);
1484:   } else {
1485:     MatDenseGetArrayWrite(C,&c);
1486:   }
1487:   MatDenseGetArrayRead(B,&b);
1488:   b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;
1489:   c1 = c; c2 = c1 + clda; c3 = c2 + clda; c4 = c3 + clda;
1490:   for (col=0; col<(cn/4)*4; col += 4) {  /* over columns of C */
1491:     for (i=0; i<am; i++) {        /* over rows of A in those columns */
1492:       r1 = r2 = r3 = r4 = 0.0;
1493:       n  = a->i[i+1] - a->i[i];
1494:       aj = a->j + a->i[i];
1495:       aa = av + a->i[i];
1496:       for (j=0; j<n; j++) {
1497:         const PetscScalar aatmp = aa[j];
1498:         const PetscInt    ajtmp = aj[j];
1499:         r1 += aatmp*b1[ajtmp];
1500:         r2 += aatmp*b2[ajtmp];
1501:         r3 += aatmp*b3[ajtmp];
1502:         r4 += aatmp*b4[ajtmp];
1503:       }
1504:       if (add) {
1505:         c1[i] += r1;
1506:         c2[i] += r2;
1507:         c3[i] += r3;
1508:         c4[i] += r4;
1509:       } else {
1510:         c1[i] = r1;
1511:         c2[i] = r2;
1512:         c3[i] = r3;
1513:         c4[i] = r4;
1514:       }
1515:     }
1516:     b1 += bm4; b2 += bm4; b3 += bm4; b4 += bm4;
1517:     c1 += am4; c2 += am4; c3 += am4; c4 += am4;
1518:   }
1519:   /* process remaining columns */
1520:   if (col != cn) {
1521:     PetscInt rc = cn-col;

1523:     if (rc == 1) {
1524:       for (i=0; i<am; i++) {
1525:         r1 = 0.0;
1526:         n  = a->i[i+1] - a->i[i];
1527:         aj = a->j + a->i[i];
1528:         aa = av + a->i[i];
1529:         for (j=0; j<n; j++) r1 += aa[j]*b1[aj[j]];
1530:         if (add) c1[i] += r1;
1531:         else c1[i] = r1;
1532:       }
1533:     } else if (rc == 2) {
1534:       for (i=0; i<am; i++) {
1535:         r1 = r2 = 0.0;
1536:         n  = a->i[i+1] - a->i[i];
1537:         aj = a->j + a->i[i];
1538:         aa = av + a->i[i];
1539:         for (j=0; j<n; j++) {
1540:           const PetscScalar aatmp = aa[j];
1541:           const PetscInt    ajtmp = aj[j];
1542:           r1 += aatmp*b1[ajtmp];
1543:           r2 += aatmp*b2[ajtmp];
1544:         }
1545:         if (add) {
1546:           c1[i] += r1;
1547:           c2[i] += r2;
1548:         } else {
1549:           c1[i] = r1;
1550:           c2[i] = r2;
1551:         }
1552:       }
1553:     } else {
1554:       for (i=0; i<am; i++) {
1555:         r1 = r2 = r3 = 0.0;
1556:         n  = a->i[i+1] - a->i[i];
1557:         aj = a->j + a->i[i];
1558:         aa = av + a->i[i];
1559:         for (j=0; j<n; j++) {
1560:           const PetscScalar aatmp = aa[j];
1561:           const PetscInt    ajtmp = aj[j];
1562:           r1 += aatmp*b1[ajtmp];
1563:           r2 += aatmp*b2[ajtmp];
1564:           r3 += aatmp*b3[ajtmp];
1565:         }
1566:         if (add) {
1567:           c1[i] += r1;
1568:           c2[i] += r2;
1569:           c3[i] += r3;
1570:         } else {
1571:           c1[i] = r1;
1572:           c2[i] = r2;
1573:           c3[i] = r3;
1574:         }
1575:       }
1576:     }
1577:   }
1578:   PetscLogFlops(cn*(2.0*a->nz));
1579:   if (add) {
1580:     MatDenseRestoreArray(C,&c);
1581:   } else {
1582:     MatDenseRestoreArrayWrite(C,&c);
1583:   }
1584:   MatDenseRestoreArrayRead(B,&b);
1585:   MatSeqAIJRestoreArrayRead(A,&av);
1586:   return(0);
1587: }

1589: PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
1590: {

1594:   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);
1595:   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);
1596:   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);

1598:   MatMatMultNumericAdd_SeqAIJ_SeqDense(A,B,C,PETSC_FALSE);
1599:   return(0);
1600: }

1602: /* ------------------------------------------------------- */
1603: static PetscErrorCode MatProductSetFromOptions_SeqAIJ_SeqDense_AB(Mat C)
1604: {
1606:   C->ops->matmultsymbolic = MatMatMultSymbolic_SeqAIJ_SeqDense;
1607:   C->ops->productsymbolic = MatProductSymbolic_AB;
1608:   return(0);
1609: }

1611: PETSC_INTERN PetscErrorCode MatTMatTMultSymbolic_SeqAIJ_SeqDense(Mat,Mat,PetscReal,Mat);

1613: static PetscErrorCode MatProductSetFromOptions_SeqAIJ_SeqDense_AtB(Mat C)
1614: {
1616:   C->ops->transposematmultsymbolic = MatTMatTMultSymbolic_SeqAIJ_SeqDense;
1617:   C->ops->productsymbolic          = MatProductSymbolic_AtB;
1618:   return(0);
1619: }

1621: static PetscErrorCode MatProductSetFromOptions_SeqAIJ_SeqDense_ABt(Mat C)
1622: {
1624:   C->ops->mattransposemultsymbolic = MatTMatTMultSymbolic_SeqAIJ_SeqDense;
1625:   C->ops->productsymbolic          = MatProductSymbolic_ABt;
1626:   return(0);
1627: }

1629: PETSC_INTERN PetscErrorCode MatProductSetFromOptions_SeqAIJ_SeqDense(Mat C)
1630: {
1632:   Mat_Product    *product = C->product;

1635:   switch (product->type) {
1636:   case MATPRODUCT_AB:
1637:     MatProductSetFromOptions_SeqAIJ_SeqDense_AB(C);
1638:     break;
1639:   case MATPRODUCT_AtB:
1640:     MatProductSetFromOptions_SeqAIJ_SeqDense_AtB(C);
1641:     break;
1642:   case MATPRODUCT_ABt:
1643:     MatProductSetFromOptions_SeqAIJ_SeqDense_ABt(C);
1644:     break;
1645:   default:
1646:     break;
1647:   }
1648:   return(0);
1649: }
1650: /* ------------------------------------------------------- */
1651: static PetscErrorCode MatProductSetFromOptions_SeqXBAIJ_SeqDense_AB(Mat C)
1652: {
1654:   Mat_Product    *product = C->product;
1655:   Mat            A = product->A;
1656:   PetscBool      baij;

1659:   PetscObjectTypeCompare((PetscObject)A,MATSEQBAIJ,&baij);
1660:   if (!baij) { /* A is seqsbaij */
1661:     PetscBool sbaij;
1662:     PetscObjectTypeCompare((PetscObject)A,MATSEQSBAIJ,&sbaij);
1663:     if (!sbaij) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_WRONGSTATE,"Mat must be either seqbaij or seqsbaij format");

1665:     C->ops->matmultsymbolic = MatMatMultSymbolic_SeqSBAIJ_SeqDense;
1666:   } else { /* A is seqbaij */
1667:     C->ops->matmultsymbolic = MatMatMultSymbolic_SeqBAIJ_SeqDense;
1668:   }

1670:   C->ops->productsymbolic = MatProductSymbolic_AB;
1671:   return(0);
1672: }

1674: PETSC_INTERN PetscErrorCode MatProductSetFromOptions_SeqXBAIJ_SeqDense(Mat C)
1675: {
1677:   Mat_Product    *product = C->product;

1680:   MatCheckProduct(C,1);
1681:   if (!product->A) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing A");
1682:   if (product->type == MATPRODUCT_AB || (product->type == MATPRODUCT_AtB && product->A->symmetric)) {
1683:     MatProductSetFromOptions_SeqXBAIJ_SeqDense_AB(C);
1684:   }
1685:   return(0);
1686: }

1688: /* ------------------------------------------------------- */
1689: static PetscErrorCode MatProductSetFromOptions_SeqDense_SeqAIJ_AB(Mat C)
1690: {
1692:   C->ops->matmultsymbolic = MatMatMultSymbolic_SeqDense_SeqAIJ;
1693:   C->ops->productsymbolic = MatProductSymbolic_AB;
1694:   return(0);
1695: }

1697: PETSC_INTERN PetscErrorCode MatProductSetFromOptions_SeqDense_SeqAIJ(Mat C)
1698: {
1700:   Mat_Product    *product = C->product;

1703:   if (product->type == MATPRODUCT_AB) {
1704:     MatProductSetFromOptions_SeqDense_SeqAIJ_AB(C);
1705:   }
1706:   return(0);
1707: }
1708: /* ------------------------------------------------------- */

1710: PetscErrorCode  MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring coloring,Mat B,Mat Btdense)
1711: {
1713:   Mat_SeqAIJ     *b       = (Mat_SeqAIJ*)B->data;
1714:   Mat_SeqDense   *btdense = (Mat_SeqDense*)Btdense->data;
1715:   PetscInt       *bi      = b->i,*bj=b->j;
1716:   PetscInt       m        = Btdense->rmap->n,n=Btdense->cmap->n,j,k,l,col,anz,*btcol,brow,ncolumns;
1717:   MatScalar      *btval,*btval_den,*ba=b->a;
1718:   PetscInt       *columns=coloring->columns,*colorforcol=coloring->colorforcol,ncolors=coloring->ncolors;

1721:   btval_den=btdense->v;
1722:   PetscArrayzero(btval_den,m*n);
1723:   for (k=0; k<ncolors; k++) {
1724:     ncolumns = coloring->ncolumns[k];
1725:     for (l=0; l<ncolumns; l++) { /* insert a row of B to a column of Btdense */
1726:       col   = *(columns + colorforcol[k] + l);
1727:       btcol = bj + bi[col];
1728:       btval = ba + bi[col];
1729:       anz   = bi[col+1] - bi[col];
1730:       for (j=0; j<anz; j++) {
1731:         brow            = btcol[j];
1732:         btval_den[brow] = btval[j];
1733:       }
1734:     }
1735:     btval_den += m;
1736:   }
1737:   return(0);
1738: }

1740: PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring matcoloring,Mat Cden,Mat Csp)
1741: {
1742:   PetscErrorCode    ierr;
1743:   Mat_SeqAIJ        *csp = (Mat_SeqAIJ*)Csp->data;
1744:   const PetscScalar *ca_den,*ca_den_ptr;
1745:   PetscScalar       *ca=csp->a;
1746:   PetscInt          k,l,m=Cden->rmap->n,ncolors=matcoloring->ncolors;
1747:   PetscInt          brows=matcoloring->brows,*den2sp=matcoloring->den2sp;
1748:   PetscInt          nrows,*row,*idx;
1749:   PetscInt          *rows=matcoloring->rows,*colorforrow=matcoloring->colorforrow;

1752:   MatDenseGetArrayRead(Cden,&ca_den);

1754:   if (brows > 0) {
1755:     PetscInt *lstart,row_end,row_start;
1756:     lstart = matcoloring->lstart;
1757:     PetscArrayzero(lstart,ncolors);

1759:     row_end = brows;
1760:     if (row_end > m) row_end = m;
1761:     for (row_start=0; row_start<m; row_start+=brows) { /* loop over row blocks of Csp */
1762:       ca_den_ptr = ca_den;
1763:       for (k=0; k<ncolors; k++) { /* loop over colors (columns of Cden) */
1764:         nrows = matcoloring->nrows[k];
1765:         row   = rows  + colorforrow[k];
1766:         idx   = den2sp + colorforrow[k];
1767:         for (l=lstart[k]; l<nrows; l++) {
1768:           if (row[l] >= row_end) {
1769:             lstart[k] = l;
1770:             break;
1771:           } else {
1772:             ca[idx[l]] = ca_den_ptr[row[l]];
1773:           }
1774:         }
1775:         ca_den_ptr += m;
1776:       }
1777:       row_end += brows;
1778:       if (row_end > m) row_end = m;
1779:     }
1780:   } else { /* non-blocked impl: loop over columns of Csp - slow if Csp is large */
1781:     ca_den_ptr = ca_den;
1782:     for (k=0; k<ncolors; k++) {
1783:       nrows = matcoloring->nrows[k];
1784:       row   = rows  + colorforrow[k];
1785:       idx   = den2sp + colorforrow[k];
1786:       for (l=0; l<nrows; l++) {
1787:         ca[idx[l]] = ca_den_ptr[row[l]];
1788:       }
1789:       ca_den_ptr += m;
1790:     }
1791:   }

1793:   MatDenseRestoreArrayRead(Cden,&ca_den);
1794: #if defined(PETSC_USE_INFO)
1795:   if (matcoloring->brows > 0) {
1796:     PetscInfo1(Csp,"Loop over %D row blocks for den2sp\n",brows);
1797:   } else {
1798:     PetscInfo(Csp,"Loop over colors/columns of Cden, inefficient for large sparse matrix product \n");
1799:   }
1800: #endif
1801:   return(0);
1802: }

1804: PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat mat,ISColoring iscoloring,MatTransposeColoring c)
1805: {
1807:   PetscInt       i,n,nrows,Nbs,j,k,m,ncols,col,cm;
1808:   const PetscInt *is,*ci,*cj,*row_idx;
1809:   PetscInt       nis = iscoloring->n,*rowhit,bs = 1;
1810:   IS             *isa;
1811:   Mat_SeqAIJ     *csp = (Mat_SeqAIJ*)mat->data;
1812:   PetscInt       *colorforrow,*rows,*rows_i,*idxhit,*spidx,*den2sp,*den2sp_i;
1813:   PetscInt       *colorforcol,*columns,*columns_i,brows;
1814:   PetscBool      flg;

1817:   ISColoringGetIS(iscoloring,PETSC_USE_POINTER,PETSC_IGNORE,&isa);

1819:   /* bs >1 is not being tested yet! */
1820:   Nbs       = mat->cmap->N/bs;
1821:   c->M      = mat->rmap->N/bs;  /* set total rows, columns and local rows */
1822:   c->N      = Nbs;
1823:   c->m      = c->M;
1824:   c->rstart = 0;
1825:   c->brows  = 100;

1827:   c->ncolors = nis;
1828:   PetscMalloc3(nis,&c->ncolumns,nis,&c->nrows,nis+1,&colorforrow);
1829:   PetscMalloc1(csp->nz+1,&rows);
1830:   PetscMalloc1(csp->nz+1,&den2sp);

1832:   brows = c->brows;
1833:   PetscOptionsGetInt(NULL,NULL,"-matden2sp_brows",&brows,&flg);
1834:   if (flg) c->brows = brows;
1835:   if (brows > 0) {
1836:     PetscMalloc1(nis+1,&c->lstart);
1837:   }

1839:   colorforrow[0] = 0;
1840:   rows_i         = rows;
1841:   den2sp_i       = den2sp;

1843:   PetscMalloc1(nis+1,&colorforcol);
1844:   PetscMalloc1(Nbs+1,&columns);

1846:   colorforcol[0] = 0;
1847:   columns_i      = columns;

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

1852:   cm   = c->m;
1853:   PetscMalloc1(cm+1,&rowhit);
1854:   PetscMalloc1(cm+1,&idxhit);
1855:   for (i=0; i<nis; i++) { /* loop over color */
1856:     ISGetLocalSize(isa[i],&n);
1857:     ISGetIndices(isa[i],&is);

1859:     c->ncolumns[i] = n;
1860:     if (n) {
1861:       PetscArraycpy(columns_i,is,n);
1862:     }
1863:     colorforcol[i+1] = colorforcol[i] + n;
1864:     columns_i       += n;

1866:     /* fast, crude version requires O(N*N) work */
1867:     PetscArrayzero(rowhit,cm);

1869:     for (j=0; j<n; j++) { /* loop over columns*/
1870:       col     = is[j];
1871:       row_idx = cj + ci[col];
1872:       m       = ci[col+1] - ci[col];
1873:       for (k=0; k<m; k++) { /* loop over columns marking them in rowhit */
1874:         idxhit[*row_idx]   = spidx[ci[col] + k];
1875:         rowhit[*row_idx++] = col + 1;
1876:       }
1877:     }
1878:     /* count the number of hits */
1879:     nrows = 0;
1880:     for (j=0; j<cm; j++) {
1881:       if (rowhit[j]) nrows++;
1882:     }
1883:     c->nrows[i]      = nrows;
1884:     colorforrow[i+1] = colorforrow[i] + nrows;

1886:     nrows = 0;
1887:     for (j=0; j<cm; j++) { /* loop over rows */
1888:       if (rowhit[j]) {
1889:         rows_i[nrows]   = j;
1890:         den2sp_i[nrows] = idxhit[j];
1891:         nrows++;
1892:       }
1893:     }
1894:     den2sp_i += nrows;

1896:     ISRestoreIndices(isa[i],&is);
1897:     rows_i += nrows;
1898:   }
1899:   MatRestoreColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);
1900:   PetscFree(rowhit);
1901:   ISColoringRestoreIS(iscoloring,PETSC_USE_POINTER,&isa);
1902:   if (csp->nz != colorforrow[nis]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"csp->nz %d != colorforrow[nis] %d",csp->nz,colorforrow[nis]);

1904:   c->colorforrow = colorforrow;
1905:   c->rows        = rows;
1906:   c->den2sp      = den2sp;
1907:   c->colorforcol = colorforcol;
1908:   c->columns     = columns;

1910:   PetscFree(idxhit);
1911:   return(0);
1912: }

1914: /* --------------------------------------------------------------- */
1915: static PetscErrorCode MatProductNumeric_AtB_SeqAIJ_SeqAIJ(Mat C)
1916: {
1918:   Mat_Product    *product = C->product;
1919:   Mat            A=product->A,B=product->B;

1922:   if (C->ops->mattransposemultnumeric) {
1923:     /* Alg: "outerproduct" */
1924:     (*C->ops->mattransposemultnumeric)(A,B,C);
1925:   } else {
1926:     /* Alg: "matmatmult" -- C = At*B */
1927:     Mat_MatTransMatMult *atb = (Mat_MatTransMatMult *)product->data;
1928:     Mat                 At;

1930:     if (!atb) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing product struct");
1931:     At = atb->At;
1932:     if (atb->updateAt && At) { /* At is computed in MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ() */
1933:       MatTranspose_SeqAIJ(A,MAT_REUSE_MATRIX,&At);
1934:     }
1935:     MatMatMultNumeric_SeqAIJ_SeqAIJ(At ? At : A,B,C);
1936:     atb->updateAt = PETSC_TRUE;
1937:   }
1938:   return(0);
1939: }

1941: static PetscErrorCode MatProductSymbolic_AtB_SeqAIJ_SeqAIJ(Mat C)
1942: {
1944:   Mat_Product    *product = C->product;
1945:   Mat            A=product->A,B=product->B;
1946:   PetscReal      fill=product->fill;

1949:   MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);

1951:   C->ops->productnumeric = MatProductNumeric_AtB_SeqAIJ_SeqAIJ;
1952:   return(0);
1953: }

1955: /* --------------------------------------------------------------- */
1956: static PetscErrorCode MatProductSetFromOptions_SeqAIJ_AB(Mat C)
1957: {
1959:   Mat_Product    *product = C->product;
1960:   PetscInt       alg = 0; /* default algorithm */
1961:   PetscBool      flg = PETSC_FALSE;
1962: #if !defined(PETSC_HAVE_HYPRE)
1963:   const char     *algTypes[7] = {"sorted","scalable","scalable_fast","heap","btheap","llcondensed","rowmerge"};
1964:   PetscInt       nalg = 7;
1965: #else
1966:   const char     *algTypes[8] = {"sorted","scalable","scalable_fast","heap","btheap","llcondensed","rowmerge","hypre"};
1967:   PetscInt       nalg = 8;
1968: #endif

1971:   /* Set default algorithm */
1972:   PetscStrcmp(C->product->alg,"default",&flg);
1973:   if (flg) {
1974:     MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);
1975:   }

1977:   /* Get runtime option */
1978:   if (product->api_user) {
1979:     PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatMatMult","Mat");
1980:     PetscOptionsEList("-matmatmult_via","Algorithmic approach","MatMatMult",algTypes,nalg,algTypes[0],&alg,&flg);
1981:     PetscOptionsEnd();
1982:   } else {
1983:     PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatProduct_AB","Mat");
1984:     PetscOptionsEList("-matproduct_ab_via","Algorithmic approach","MatProduct_AB",algTypes,nalg,algTypes[0],&alg,&flg);
1985:     PetscOptionsEnd();
1986:   }
1987:   if (flg) {
1988:     MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);
1989:   }

1991:   C->ops->productsymbolic = MatProductSymbolic_AB;
1992:   C->ops->matmultsymbolic = MatMatMultSymbolic_SeqAIJ_SeqAIJ;
1993:   return(0);
1994: }

1996: static PetscErrorCode MatProductSetFromOptions_SeqAIJ_AtB(Mat C)
1997: {
1999:   Mat_Product    *product = C->product;
2000:   PetscInt       alg = 0; /* default algorithm */
2001:   PetscBool      flg = PETSC_FALSE;
2002:   const char     *algTypes[3] = {"default","at*b","outerproduct"};
2003:   PetscInt       nalg = 3;

2006:   /* Get runtime option */
2007:   if (product->api_user) {
2008:     PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatTransposeMatMult","Mat");
2009:     PetscOptionsEList("-mattransposematmult_via","Algorithmic approach","MatTransposeMatMult",algTypes,nalg,algTypes[alg],&alg,&flg);
2010:     PetscOptionsEnd();
2011:   } else {
2012:     PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatProduct_AtB","Mat");
2013:     PetscOptionsEList("-matproduct_atb_via","Algorithmic approach","MatProduct_AtB",algTypes,nalg,algTypes[alg],&alg,&flg);
2014:     PetscOptionsEnd();
2015:   }
2016:   if (flg) {
2017:     MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);
2018:   }

2020:   C->ops->productsymbolic = MatProductSymbolic_AtB_SeqAIJ_SeqAIJ;
2021:   return(0);
2022: }

2024: static PetscErrorCode MatProductSetFromOptions_SeqAIJ_ABt(Mat C)
2025: {
2027:   Mat_Product    *product = C->product;
2028:   PetscInt       alg = 0; /* default algorithm */
2029:   PetscBool      flg = PETSC_FALSE;
2030:   const char     *algTypes[2] = {"default","color"};
2031:   PetscInt       nalg = 2;

2034:   /* Set default algorithm */
2035:   PetscStrcmp(C->product->alg,"default",&flg);
2036:   if (!flg) {
2037:     alg = 1;
2038:     MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);
2039:   }

2041:   /* Get runtime option */
2042:   if (product->api_user) {
2043:     PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatMatTransposeMult","Mat");
2044:     PetscOptionsEList("-matmattransmult_via","Algorithmic approach","MatMatTransposeMult",algTypes,nalg,algTypes[alg],&alg,&flg);
2045:     PetscOptionsEnd();
2046:   } else {
2047:     PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatProduct_ABt","Mat");
2048:     PetscOptionsEList("-matproduct_abt_via","Algorithmic approach","MatProduct_ABt",algTypes,nalg,algTypes[alg],&alg,&flg);
2049:     PetscOptionsEnd();
2050:   }
2051:   if (flg) {
2052:     MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);
2053:   }

2055:   C->ops->mattransposemultsymbolic = MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ;
2056:   C->ops->productsymbolic          = MatProductSymbolic_ABt;
2057:   return(0);
2058: }

2060: static PetscErrorCode MatProductSetFromOptions_SeqAIJ_PtAP(Mat C)
2061: {
2063:   Mat_Product    *product = C->product;
2064:   PetscBool      flg = PETSC_FALSE;
2065:   PetscInt       alg = 0; /* default algorithm -- alg=1 should be default!!! */
2066: #if !defined(PETSC_HAVE_HYPRE)
2067:   const char      *algTypes[2] = {"scalable","rap"};
2068:   PetscInt        nalg = 2;
2069: #else
2070:   const char      *algTypes[3] = {"scalable","rap","hypre"};
2071:   PetscInt        nalg = 3;
2072: #endif

2075:   /* Set default algorithm */
2076:   PetscStrcmp(product->alg,"default",&flg);
2077:   if (flg) {
2078:     MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);
2079:   }

2081:   /* Get runtime option */
2082:   if (product->api_user) {
2083:     PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatPtAP","Mat");
2084:     PetscOptionsEList("-matptap_via","Algorithmic approach","MatPtAP",algTypes,nalg,algTypes[0],&alg,&flg);
2085:     PetscOptionsEnd();
2086:   } else {
2087:     PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatProduct_PtAP","Mat");
2088:     PetscOptionsEList("-matproduct_ptap_via","Algorithmic approach","MatProduct_PtAP",algTypes,nalg,algTypes[0],&alg,&flg);
2089:     PetscOptionsEnd();
2090:   }
2091:   if (flg) {
2092:     MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);
2093:   }

2095:   C->ops->productsymbolic = MatProductSymbolic_PtAP_SeqAIJ_SeqAIJ;
2096:   return(0);
2097: }

2099: static PetscErrorCode MatProductSetFromOptions_SeqAIJ_RARt(Mat C)
2100: {
2102:   Mat_Product    *product = C->product;
2103:   PetscBool      flg = PETSC_FALSE;
2104:   PetscInt       alg = 0; /* default algorithm */
2105:   const char     *algTypes[3] = {"r*a*rt","r*art","coloring_rart"};
2106:   PetscInt        nalg = 3;

2109:   /* Set default algorithm */
2110:   PetscStrcmp(product->alg,"default",&flg);
2111:   if (flg) {
2112:     MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);
2113:   }

2115:   /* Get runtime option */
2116:   if (product->api_user) {
2117:     PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatRARt","Mat");
2118:     PetscOptionsEList("-matrart_via","Algorithmic approach","MatRARt",algTypes,nalg,algTypes[0],&alg,&flg);
2119:     PetscOptionsEnd();
2120:   } else {
2121:     PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatProduct_RARt","Mat");
2122:     PetscOptionsEList("-matproduct_rart_via","Algorithmic approach","MatProduct_RARt",algTypes,nalg,algTypes[0],&alg,&flg);
2123:     PetscOptionsEnd();
2124:   }
2125:   if (flg) {
2126:     MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);
2127:   }

2129:   C->ops->productsymbolic = MatProductSymbolic_RARt_SeqAIJ_SeqAIJ;
2130:   return(0);
2131: }

2133: /* ABC = A*B*C = A*(B*C); ABC's algorithm must be chosen from AB's algorithm */
2134: static PetscErrorCode MatProductSetFromOptions_SeqAIJ_ABC(Mat C)
2135: {
2137:   Mat_Product    *product = C->product;
2138:   PetscInt       alg = 0; /* default algorithm */
2139:   PetscBool      flg = PETSC_FALSE;
2140:   const char     *algTypes[7] = {"sorted","scalable","scalable_fast","heap","btheap","llcondensed","rowmerge"};
2141:   PetscInt       nalg = 7;

2144:   /* Set default algorithm */
2145:   PetscStrcmp(product->alg,"default",&flg);
2146:   if (flg) {
2147:     MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);
2148:   }

2150:   /* Get runtime option */
2151:   if (product->api_user) {
2152:     PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatMatMatMult","Mat");
2153:     PetscOptionsEList("-matmatmatmult_via","Algorithmic approach","MatMatMatMult",algTypes,nalg,algTypes[alg],&alg,&flg);
2154:     PetscOptionsEnd();
2155:   } else {
2156:     PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatProduct_ABC","Mat");
2157:     PetscOptionsEList("-matproduct_abc_via","Algorithmic approach","MatProduct_ABC",algTypes,nalg,algTypes[alg],&alg,&flg);
2158:     PetscOptionsEnd();
2159:   }
2160:   if (flg) {
2161:     MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);
2162:   }

2164:   C->ops->matmatmultsymbolic = MatMatMatMultSymbolic_SeqAIJ_SeqAIJ_SeqAIJ;
2165:   C->ops->productsymbolic    = MatProductSymbolic_ABC;
2166:   return(0);
2167: }

2169: PetscErrorCode MatProductSetFromOptions_SeqAIJ(Mat C)
2170: {
2172:   Mat_Product    *product = C->product;

2175:   switch (product->type) {
2176:   case MATPRODUCT_AB:
2177:     MatProductSetFromOptions_SeqAIJ_AB(C);
2178:     break;
2179:   case MATPRODUCT_AtB:
2180:     MatProductSetFromOptions_SeqAIJ_AtB(C);
2181:     break;
2182:   case MATPRODUCT_ABt:
2183:     MatProductSetFromOptions_SeqAIJ_ABt(C);
2184:     break;
2185:   case MATPRODUCT_PtAP:
2186:     MatProductSetFromOptions_SeqAIJ_PtAP(C);
2187:     break;
2188:   case MATPRODUCT_RARt:
2189:     MatProductSetFromOptions_SeqAIJ_RARt(C);
2190:     break;
2191:   case MATPRODUCT_ABC:
2192:     MatProductSetFromOptions_SeqAIJ_ABC(C);
2193:     break;
2194:   default:
2195:     break;
2196:   }
2197:   return(0);
2198: }