Actual source code: matmatmult.c
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: }
138: PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_LLCondensed(Mat A,Mat B,PetscReal fill,Mat C)
139: {
140: PetscErrorCode ierr;
141: Mat_SeqAIJ *a =(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
142: PetscInt *ai=a->i,*bi=b->i,*ci,*cj;
143: PetscInt am =A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
144: PetscReal afill;
145: PetscInt i,j,anzi,brow,bnzj,cnzi,*bj,*aj,*lnk,ndouble=0,Crmax;
146: PetscTable ta;
147: PetscBT lnkbt;
148: PetscFreeSpaceList free_space=NULL,current_space=NULL;
151: /* Get ci and cj */
152: /*---------------*/
153: /* Allocate ci array, arrays for fill computation and */
154: /* free space for accumulating nonzero column info */
155: PetscMalloc1(am+2,&ci);
156: ci[0] = 0;
158: /* create and initialize a linked list */
159: PetscTableCreate(bn,bn,&ta);
160: MatRowMergeMax_SeqAIJ(b,bm,ta);
161: PetscTableGetCount(ta,&Crmax);
162: PetscTableDestroy(&ta);
164: PetscLLCondensedCreate(Crmax,bn,&lnk,&lnkbt);
166: /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
167: PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am],bi[bm])),&free_space);
169: current_space = free_space;
171: /* Determine ci and cj */
172: for (i=0; i<am; i++) {
173: anzi = ai[i+1] - ai[i];
174: aj = a->j + ai[i];
175: for (j=0; j<anzi; j++) {
176: brow = aj[j];
177: bnzj = bi[brow+1] - bi[brow];
178: bj = b->j + bi[brow];
179: /* add non-zero cols of B into the sorted linked list lnk */
180: PetscLLCondensedAddSorted(bnzj,bj,lnk,lnkbt);
181: }
182: /* add possible missing diagonal entry */
183: if (C->force_diagonals) {
184: PetscLLCondensedAddSorted(1,&i,lnk,lnkbt);
185: }
186: cnzi = lnk[0];
188: /* If free space is not available, make more free space */
189: /* Double the amount of total space in the list */
190: if (current_space->local_remaining<cnzi) {
191: PetscFreeSpaceGet(PetscIntSumTruncate(cnzi,current_space->total_array_size),¤t_space);
192: ndouble++;
193: }
195: /* Copy data into free space, then initialize lnk */
196: PetscLLCondensedClean(bn,cnzi,current_space->array,lnk,lnkbt);
198: current_space->array += cnzi;
199: current_space->local_used += cnzi;
200: current_space->local_remaining -= cnzi;
202: ci[i+1] = ci[i] + cnzi;
203: }
205: /* Column indices are in the list of free space */
206: /* Allocate space for cj, initialize cj, and */
207: /* destroy list of free space and other temporary array(s) */
208: PetscMalloc1(ci[am]+1,&cj);
209: PetscFreeSpaceContiguous(&free_space,cj);
210: PetscLLCondensedDestroy(lnk,lnkbt);
212: /* put together the new symbolic matrix */
213: MatSetSeqAIJWithArrays_private(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,((PetscObject)A)->type_name,C);
214: MatSetBlockSizesFromMats(C,A,B);
216: /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
217: /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
218: c = (Mat_SeqAIJ*)(C->data);
219: c->free_a = PETSC_FALSE;
220: c->free_ij = PETSC_TRUE;
221: c->nonew = 0;
223: /* fast, needs non-scalable O(bn) array 'abdense' */
224: C->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Sorted;
226: /* set MatInfo */
227: afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
228: if (afill < 1.0) afill = 1.0;
229: c->maxnz = ci[am];
230: c->nz = ci[am];
231: C->info.mallocs = ndouble;
232: C->info.fill_ratio_given = fill;
233: C->info.fill_ratio_needed = afill;
235: #if defined(PETSC_USE_INFO)
236: if (ci[am]) {
237: PetscInfo3(C,"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);
238: PetscInfo1(C,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);
239: } else {
240: PetscInfo(C,"Empty matrix product\n");
241: }
242: #endif
243: return(0);
244: }
246: PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Sorted(Mat A,Mat B,Mat C)
247: {
248: PetscErrorCode ierr;
249: PetscLogDouble flops=0.0;
250: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
251: Mat_SeqAIJ *b = (Mat_SeqAIJ*)B->data;
252: Mat_SeqAIJ *c = (Mat_SeqAIJ*)C->data;
253: PetscInt *ai =a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j;
254: PetscInt am =A->rmap->n,cm=C->rmap->n;
255: PetscInt i,j,k,anzi,bnzi,cnzi,brow;
256: PetscScalar *ca,valtmp;
257: PetscScalar *ab_dense;
258: PetscContainer cab_dense;
259: const PetscScalar *aa,*ba,*baj;
262: MatSeqAIJGetArrayRead(A,&aa);
263: MatSeqAIJGetArrayRead(B,&ba);
264: if (!c->a) { /* first call of MatMatMultNumeric_SeqAIJ_SeqAIJ, allocate ca and matmult_abdense */
265: PetscMalloc1(ci[cm]+1,&ca);
266: c->a = ca;
267: c->free_a = PETSC_TRUE;
268: } else ca = c->a;
270: /* TODO this should be done in the symbolic phase */
271: /* However, this function is so heavily used (sometimes in an hidden way through multnumeric function pointers
272: that is hard to eradicate) */
273: PetscObjectQuery((PetscObject)C,"__PETSc__ab_dense",(PetscObject*)&cab_dense);
274: if (!cab_dense) {
275: PetscMalloc1(B->cmap->N,&ab_dense);
276: PetscContainerCreate(PETSC_COMM_SELF,&cab_dense);
277: PetscContainerSetPointer(cab_dense,ab_dense);
278: PetscContainerSetUserDestroy(cab_dense,PetscContainerUserDestroyDefault);
279: PetscObjectCompose((PetscObject)C,"__PETSc__ab_dense",(PetscObject)cab_dense);
280: PetscObjectDereference((PetscObject)cab_dense);
281: }
282: PetscContainerGetPointer(cab_dense,(void**)&ab_dense);
283: PetscArrayzero(ab_dense,B->cmap->N);
285: /* clean old values in C */
286: PetscArrayzero(ca,ci[cm]);
287: /* Traverse A row-wise. */
288: /* Build the ith row in C by summing over nonzero columns in A, */
289: /* the rows of B corresponding to nonzeros of A. */
290: for (i=0; i<am; i++) {
291: anzi = ai[i+1] - ai[i];
292: for (j=0; j<anzi; j++) {
293: brow = aj[j];
294: bnzi = bi[brow+1] - bi[brow];
295: bjj = bj + bi[brow];
296: baj = ba + bi[brow];
297: /* perform dense axpy */
298: valtmp = aa[j];
299: for (k=0; k<bnzi; k++) {
300: ab_dense[bjj[k]] += valtmp*baj[k];
301: }
302: flops += 2*bnzi;
303: }
304: aj += anzi; aa += anzi;
306: cnzi = ci[i+1] - ci[i];
307: for (k=0; k<cnzi; k++) {
308: ca[k] += ab_dense[cj[k]];
309: ab_dense[cj[k]] = 0.0; /* zero ab_dense */
310: }
311: flops += cnzi;
312: cj += cnzi; ca += cnzi;
313: }
314: #if defined(PETSC_HAVE_DEVICE)
315: if (C->offloadmask != PETSC_OFFLOAD_UNALLOCATED) C->offloadmask = PETSC_OFFLOAD_CPU;
316: #endif
317: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
318: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
319: PetscLogFlops(flops);
320: MatSeqAIJRestoreArrayRead(A,&aa);
321: MatSeqAIJRestoreArrayRead(B,&ba);
322: return(0);
323: }
325: PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(Mat A,Mat B,Mat C)
326: {
327: PetscErrorCode ierr;
328: PetscLogDouble flops=0.0;
329: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
330: Mat_SeqAIJ *b = (Mat_SeqAIJ*)B->data;
331: Mat_SeqAIJ *c = (Mat_SeqAIJ*)C->data;
332: PetscInt *ai = a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j;
333: PetscInt am = A->rmap->N,cm=C->rmap->N;
334: PetscInt i,j,k,anzi,bnzi,cnzi,brow;
335: PetscScalar *ca=c->a,valtmp;
336: const PetscScalar *aa,*ba,*baj;
337: PetscInt nextb;
340: MatSeqAIJGetArrayRead(A,&aa);
341: MatSeqAIJGetArrayRead(B,&ba);
342: if (!ca) { /* first call of MatMatMultNumeric_SeqAIJ_SeqAIJ, allocate ca and matmult_abdense */
343: PetscMalloc1(ci[cm]+1,&ca);
344: c->a = ca;
345: c->free_a = PETSC_TRUE;
346: }
348: /* clean old values in C */
349: PetscArrayzero(ca,ci[cm]);
350: /* Traverse A row-wise. */
351: /* Build the ith row in C by summing over nonzero columns in A, */
352: /* the rows of B corresponding to nonzeros of A. */
353: for (i=0; i<am; i++) {
354: anzi = ai[i+1] - ai[i];
355: cnzi = ci[i+1] - ci[i];
356: for (j=0; j<anzi; j++) {
357: brow = aj[j];
358: bnzi = bi[brow+1] - bi[brow];
359: bjj = bj + bi[brow];
360: baj = ba + bi[brow];
361: /* perform sparse axpy */
362: valtmp = aa[j];
363: nextb = 0;
364: for (k=0; nextb<bnzi; k++) {
365: if (cj[k] == bjj[nextb]) { /* ccol == bcol */
366: ca[k] += valtmp*baj[nextb++];
367: }
368: }
369: flops += 2*bnzi;
370: }
371: aj += anzi; aa += anzi;
372: cj += cnzi; ca += cnzi;
373: }
374: #if defined(PETSC_HAVE_DEVICE)
375: if (C->offloadmask != PETSC_OFFLOAD_UNALLOCATED) C->offloadmask = PETSC_OFFLOAD_CPU;
376: #endif
377: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
378: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
379: PetscLogFlops(flops);
380: MatSeqAIJRestoreArrayRead(A,&aa);
381: MatSeqAIJRestoreArrayRead(B,&ba);
382: return(0);
383: }
385: PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(Mat A,Mat B,PetscReal fill,Mat C)
386: {
387: PetscErrorCode ierr;
388: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
389: PetscInt *ai = a->i,*bi=b->i,*ci,*cj;
390: PetscInt am = A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
391: MatScalar *ca;
392: PetscReal afill;
393: PetscInt i,j,anzi,brow,bnzj,cnzi,*bj,*aj,*lnk,ndouble=0,Crmax;
394: PetscTable ta;
395: PetscFreeSpaceList free_space=NULL,current_space=NULL;
398: /* Get ci and cj - same as MatMatMultSymbolic_SeqAIJ_SeqAIJ except using PetscLLxxx_fast() */
399: /*-----------------------------------------------------------------------------------------*/
400: /* Allocate arrays for fill computation and free space for accumulating nonzero column */
401: PetscMalloc1(am+2,&ci);
402: ci[0] = 0;
404: /* create and initialize a linked list */
405: PetscTableCreate(bn,bn,&ta);
406: MatRowMergeMax_SeqAIJ(b,bm,ta);
407: PetscTableGetCount(ta,&Crmax);
408: PetscTableDestroy(&ta);
410: PetscLLCondensedCreate_fast(Crmax,&lnk);
412: /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
413: PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am],bi[bm])),&free_space);
414: current_space = free_space;
416: /* Determine ci and cj */
417: for (i=0; i<am; i++) {
418: anzi = ai[i+1] - ai[i];
419: aj = a->j + ai[i];
420: for (j=0; j<anzi; j++) {
421: brow = aj[j];
422: bnzj = bi[brow+1] - bi[brow];
423: bj = b->j + bi[brow];
424: /* add non-zero cols of B into the sorted linked list lnk */
425: PetscLLCondensedAddSorted_fast(bnzj,bj,lnk);
426: }
427: /* add possible missing diagonal entry */
428: if (C->force_diagonals) {
429: PetscLLCondensedAddSorted_fast(1,&i,lnk);
430: }
431: cnzi = lnk[1];
433: /* If free space is not available, make more free space */
434: /* Double the amount of total space in the list */
435: if (current_space->local_remaining<cnzi) {
436: PetscFreeSpaceGet(PetscIntSumTruncate(cnzi,current_space->total_array_size),¤t_space);
437: ndouble++;
438: }
440: /* Copy data into free space, then initialize lnk */
441: PetscLLCondensedClean_fast(cnzi,current_space->array,lnk);
443: current_space->array += cnzi;
444: current_space->local_used += cnzi;
445: current_space->local_remaining -= cnzi;
447: ci[i+1] = ci[i] + cnzi;
448: }
450: /* Column indices are in the list of free space */
451: /* Allocate space for cj, initialize cj, and */
452: /* destroy list of free space and other temporary array(s) */
453: PetscMalloc1(ci[am]+1,&cj);
454: PetscFreeSpaceContiguous(&free_space,cj);
455: PetscLLCondensedDestroy_fast(lnk);
457: /* Allocate space for ca */
458: PetscCalloc1(ci[am]+1,&ca);
460: /* put together the new symbolic matrix */
461: MatSetSeqAIJWithArrays_private(PetscObjectComm((PetscObject)A),am,bn,ci,cj,ca,((PetscObject)A)->type_name,C);
462: MatSetBlockSizesFromMats(C,A,B);
464: /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
465: /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
466: c = (Mat_SeqAIJ*)(C->data);
467: c->free_a = PETSC_TRUE;
468: c->free_ij = PETSC_TRUE;
469: c->nonew = 0;
471: /* slower, less memory */
472: C->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable;
474: /* set MatInfo */
475: afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
476: if (afill < 1.0) afill = 1.0;
477: c->maxnz = ci[am];
478: c->nz = ci[am];
479: C->info.mallocs = ndouble;
480: C->info.fill_ratio_given = fill;
481: C->info.fill_ratio_needed = afill;
483: #if defined(PETSC_USE_INFO)
484: if (ci[am]) {
485: PetscInfo3(C,"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);
486: PetscInfo1(C,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);
487: } else {
488: PetscInfo(C,"Empty matrix product\n");
489: }
490: #endif
491: return(0);
492: }
494: PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(Mat A,Mat B,PetscReal fill,Mat C)
495: {
496: PetscErrorCode ierr;
497: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
498: PetscInt *ai = a->i,*bi=b->i,*ci,*cj;
499: PetscInt am = A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
500: MatScalar *ca;
501: PetscReal afill;
502: PetscInt i,j,anzi,brow,bnzj,cnzi,*bj,*aj,*lnk,ndouble=0,Crmax;
503: PetscTable ta;
504: PetscFreeSpaceList free_space=NULL,current_space=NULL;
507: /* Get ci and cj - same as MatMatMultSymbolic_SeqAIJ_SeqAIJ except using PetscLLxxx_Scalalbe() */
508: /*---------------------------------------------------------------------------------------------*/
509: /* Allocate arrays for fill computation and free space for accumulating nonzero column */
510: PetscMalloc1(am+2,&ci);
511: ci[0] = 0;
513: /* create and initialize a linked list */
514: PetscTableCreate(bn,bn,&ta);
515: MatRowMergeMax_SeqAIJ(b,bm,ta);
516: PetscTableGetCount(ta,&Crmax);
517: PetscTableDestroy(&ta);
518: PetscLLCondensedCreate_Scalable(Crmax,&lnk);
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: /* Determine ci and cj */
525: for (i=0; i<am; i++) {
526: anzi = ai[i+1] - ai[i];
527: aj = a->j + ai[i];
528: for (j=0; j<anzi; j++) {
529: brow = aj[j];
530: bnzj = bi[brow+1] - bi[brow];
531: bj = b->j + bi[brow];
532: /* add non-zero cols of B into the sorted linked list lnk */
533: PetscLLCondensedAddSorted_Scalable(bnzj,bj,lnk);
534: }
535: /* add possible missing diagonal entry */
536: if (C->force_diagonals) {
537: PetscLLCondensedAddSorted_Scalable(1,&i,lnk);
538: }
540: cnzi = lnk[0];
542: /* If free space is not available, make more free space */
543: /* Double the amount of total space in the list */
544: if (current_space->local_remaining<cnzi) {
545: PetscFreeSpaceGet(PetscIntSumTruncate(cnzi,current_space->total_array_size),¤t_space);
546: ndouble++;
547: }
549: /* Copy data into free space, then initialize lnk */
550: PetscLLCondensedClean_Scalable(cnzi,current_space->array,lnk);
552: current_space->array += cnzi;
553: current_space->local_used += cnzi;
554: current_space->local_remaining -= cnzi;
556: ci[i+1] = ci[i] + cnzi;
557: }
559: /* Column indices are in the list of free space */
560: /* Allocate space for cj, initialize cj, and */
561: /* destroy list of free space and other temporary array(s) */
562: PetscMalloc1(ci[am]+1,&cj);
563: PetscFreeSpaceContiguous(&free_space,cj);
564: PetscLLCondensedDestroy_Scalable(lnk);
566: /* Allocate space for ca */
567: /*-----------------------*/
568: PetscCalloc1(ci[am]+1,&ca);
570: /* put together the new symbolic matrix */
571: MatSetSeqAIJWithArrays_private(PetscObjectComm((PetscObject)A),am,bn,ci,cj,ca,((PetscObject)A)->type_name,C);
572: MatSetBlockSizesFromMats(C,A,B);
574: /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
575: /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
576: c = (Mat_SeqAIJ*)(C->data);
577: c->free_a = PETSC_TRUE;
578: c->free_ij = PETSC_TRUE;
579: c->nonew = 0;
581: /* slower, less memory */
582: C->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable;
584: /* set MatInfo */
585: afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
586: if (afill < 1.0) afill = 1.0;
587: c->maxnz = ci[am];
588: c->nz = ci[am];
589: C->info.mallocs = ndouble;
590: C->info.fill_ratio_given = fill;
591: C->info.fill_ratio_needed = afill;
593: #if defined(PETSC_USE_INFO)
594: if (ci[am]) {
595: PetscInfo3(C,"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);
596: PetscInfo1(C,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);
597: } else {
598: PetscInfo(C,"Empty matrix product\n");
599: }
600: #endif
601: return(0);
602: }
604: PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(Mat A,Mat B,PetscReal fill,Mat C)
605: {
606: PetscErrorCode ierr;
607: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
608: const PetscInt *ai=a->i,*bi=b->i,*aj=a->j,*bj=b->j;
609: PetscInt *ci,*cj,*bb;
610: PetscInt am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
611: PetscReal afill;
612: PetscInt i,j,col,ndouble = 0;
613: PetscFreeSpaceList free_space=NULL,current_space=NULL;
614: PetscHeap h;
617: /* Get ci and cj - by merging sorted rows using a heap */
618: /*---------------------------------------------------------------------------------------------*/
619: /* Allocate arrays for fill computation and free space for accumulating nonzero column */
620: PetscMalloc1(am+2,&ci);
621: ci[0] = 0;
623: /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
624: PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am],bi[bm])),&free_space);
625: current_space = free_space;
627: PetscHeapCreate(a->rmax,&h);
628: PetscMalloc1(a->rmax,&bb);
630: /* Determine ci and cj */
631: for (i=0; i<am; i++) {
632: 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 */
633: const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */
634: ci[i+1] = ci[i];
635: /* Populate the min heap */
636: for (j=0; j<anzi; j++) {
637: bb[j] = bi[acol[j]]; /* bb points at the start of the row */
638: if (bb[j] < bi[acol[j]+1]) { /* Add if row is nonempty */
639: PetscHeapAdd(h,j,bj[bb[j]++]);
640: }
641: }
642: /* Pick off the min element, adding it to free space */
643: PetscHeapPop(h,&j,&col);
644: while (j >= 0) {
645: if (current_space->local_remaining < 1) { /* double the size, but don't exceed 16 MiB */
646: PetscFreeSpaceGet(PetscMin(PetscIntMultTruncate(2,current_space->total_array_size),16 << 20),¤t_space);
647: ndouble++;
648: }
649: *(current_space->array++) = col;
650: current_space->local_used++;
651: current_space->local_remaining--;
652: ci[i+1]++;
654: /* stash if anything else remains in this row of B */
655: if (bb[j] < bi[acol[j]+1]) {PetscHeapStash(h,j,bj[bb[j]++]);}
656: while (1) { /* pop and stash any other rows of B that also had an entry in this column */
657: PetscInt j2,col2;
658: PetscHeapPeek(h,&j2,&col2);
659: if (col2 != col) break;
660: PetscHeapPop(h,&j2,&col2);
661: if (bb[j2] < bi[acol[j2]+1]) {PetscHeapStash(h,j2,bj[bb[j2]++]);}
662: }
663: /* Put any stashed elements back into the min heap */
664: PetscHeapUnstash(h);
665: PetscHeapPop(h,&j,&col);
666: }
667: }
668: PetscFree(bb);
669: PetscHeapDestroy(&h);
671: /* Column indices are in the list of free space */
672: /* Allocate space for cj, initialize cj, and */
673: /* destroy list of free space and other temporary array(s) */
674: PetscMalloc1(ci[am],&cj);
675: PetscFreeSpaceContiguous(&free_space,cj);
677: /* put together the new symbolic matrix */
678: MatSetSeqAIJWithArrays_private(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,((PetscObject)A)->type_name,C);
679: MatSetBlockSizesFromMats(C,A,B);
681: /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
682: /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
683: c = (Mat_SeqAIJ*)(C->data);
684: c->free_a = PETSC_TRUE;
685: c->free_ij = PETSC_TRUE;
686: c->nonew = 0;
688: C->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Sorted;
690: /* set MatInfo */
691: afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
692: if (afill < 1.0) afill = 1.0;
693: c->maxnz = ci[am];
694: c->nz = ci[am];
695: C->info.mallocs = ndouble;
696: C->info.fill_ratio_given = fill;
697: C->info.fill_ratio_needed = afill;
699: #if defined(PETSC_USE_INFO)
700: if (ci[am]) {
701: PetscInfo3(C,"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);
702: PetscInfo1(C,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);
703: } else {
704: PetscInfo(C,"Empty matrix product\n");
705: }
706: #endif
707: return(0);
708: }
710: PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(Mat A,Mat B,PetscReal fill,Mat C)
711: {
712: PetscErrorCode ierr;
713: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
714: const PetscInt *ai = a->i,*bi=b->i,*aj=a->j,*bj=b->j;
715: PetscInt *ci,*cj,*bb;
716: PetscInt am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
717: PetscReal afill;
718: PetscInt i,j,col,ndouble = 0;
719: PetscFreeSpaceList free_space=NULL,current_space=NULL;
720: PetscHeap h;
721: PetscBT bt;
724: /* Get ci and cj - using a heap for the sorted rows, but use BT so that each index is only added once */
725: /*---------------------------------------------------------------------------------------------*/
726: /* Allocate arrays for fill computation and free space for accumulating nonzero column */
727: PetscMalloc1(am+2,&ci);
728: ci[0] = 0;
730: /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
731: PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am],bi[bm])),&free_space);
733: current_space = free_space;
735: PetscHeapCreate(a->rmax,&h);
736: PetscMalloc1(a->rmax,&bb);
737: PetscBTCreate(bn,&bt);
739: /* Determine ci and cj */
740: for (i=0; i<am; i++) {
741: 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 */
742: const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */
743: const PetscInt *fptr = current_space->array; /* Save beginning of the row so we can clear the BT later */
744: ci[i+1] = ci[i];
745: /* Populate the min heap */
746: for (j=0; j<anzi; j++) {
747: PetscInt brow = acol[j];
748: for (bb[j] = bi[brow]; bb[j] < bi[brow+1]; bb[j]++) {
749: PetscInt bcol = bj[bb[j]];
750: if (!PetscBTLookupSet(bt,bcol)) { /* new entry */
751: PetscHeapAdd(h,j,bcol);
752: bb[j]++;
753: break;
754: }
755: }
756: }
757: /* Pick off the min element, adding it to free space */
758: PetscHeapPop(h,&j,&col);
759: while (j >= 0) {
760: if (current_space->local_remaining < 1) { /* double the size, but don't exceed 16 MiB */
761: fptr = NULL; /* need PetscBTMemzero */
762: PetscFreeSpaceGet(PetscMin(PetscIntMultTruncate(2,current_space->total_array_size),16 << 20),¤t_space);
763: ndouble++;
764: }
765: *(current_space->array++) = col;
766: current_space->local_used++;
767: current_space->local_remaining--;
768: ci[i+1]++;
770: /* stash if anything else remains in this row of B */
771: for (; bb[j] < bi[acol[j]+1]; bb[j]++) {
772: PetscInt bcol = bj[bb[j]];
773: if (!PetscBTLookupSet(bt,bcol)) { /* new entry */
774: PetscHeapAdd(h,j,bcol);
775: bb[j]++;
776: break;
777: }
778: }
779: PetscHeapPop(h,&j,&col);
780: }
781: if (fptr) { /* Clear the bits for this row */
782: for (; fptr<current_space->array; fptr++) {PetscBTClear(bt,*fptr);}
783: } else { /* We reallocated so we don't remember (easily) how to clear only the bits we changed */
784: PetscBTMemzero(bn,bt);
785: }
786: }
787: PetscFree(bb);
788: PetscHeapDestroy(&h);
789: PetscBTDestroy(&bt);
791: /* Column indices are in the list of free space */
792: /* Allocate space for cj, initialize cj, and */
793: /* destroy list of free space and other temporary array(s) */
794: PetscMalloc1(ci[am],&cj);
795: PetscFreeSpaceContiguous(&free_space,cj);
797: /* put together the new symbolic matrix */
798: MatSetSeqAIJWithArrays_private(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,((PetscObject)A)->type_name,C);
799: MatSetBlockSizesFromMats(C,A,B);
801: /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
802: /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
803: c = (Mat_SeqAIJ*)(C->data);
804: c->free_a = PETSC_TRUE;
805: c->free_ij = PETSC_TRUE;
806: c->nonew = 0;
808: C->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Sorted;
810: /* set MatInfo */
811: afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
812: if (afill < 1.0) afill = 1.0;
813: c->maxnz = ci[am];
814: c->nz = ci[am];
815: C->info.mallocs = ndouble;
816: C->info.fill_ratio_given = fill;
817: C->info.fill_ratio_needed = afill;
819: #if defined(PETSC_USE_INFO)
820: if (ci[am]) {
821: PetscInfo3(C,"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);
822: PetscInfo1(C,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);
823: } else {
824: PetscInfo(C,"Empty matrix product\n");
825: }
826: #endif
827: return(0);
828: }
830: PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_RowMerge(Mat A,Mat B,PetscReal fill,Mat C)
831: {
832: PetscErrorCode ierr;
833: Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
834: const PetscInt *ai=a->i,*bi=b->i,*aj=a->j,*bj=b->j,*inputi,*inputj,*inputcol,*inputcol_L1;
835: PetscInt *ci,*cj,*outputj,worki_L1[9],worki_L2[9];
836: PetscInt c_maxmem,a_maxrownnz=0,a_rownnz;
837: const PetscInt workcol[8]={0,1,2,3,4,5,6,7};
838: const PetscInt am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
839: const PetscInt *brow_ptr[8],*brow_end[8];
840: PetscInt window[8];
841: PetscInt window_min,old_window_min,ci_nnz,outputi_nnz=0,L1_nrows,L2_nrows;
842: PetscInt i,k,ndouble=0,L1_rowsleft,rowsleft;
843: PetscReal afill;
844: PetscInt *workj_L1,*workj_L2,*workj_L3;
845: PetscInt L1_nnz,L2_nnz;
847: /* Step 1: Get upper bound on memory required for allocation.
848: Because of the way virtual memory works,
849: only the memory pages that are actually needed will be physically allocated. */
851: PetscMalloc1(am+1,&ci);
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: a_rownnz = 0;
856: for (k=0; k<anzi; ++k) {
857: a_rownnz += bi[acol[k]+1] - bi[acol[k]];
858: if (a_rownnz > bn) {
859: a_rownnz = bn;
860: break;
861: }
862: }
863: a_maxrownnz = PetscMax(a_maxrownnz, a_rownnz);
864: }
865: /* temporary work areas for merging rows */
866: PetscMalloc1(a_maxrownnz*8,&workj_L1);
867: PetscMalloc1(a_maxrownnz*8,&workj_L2);
868: PetscMalloc1(a_maxrownnz,&workj_L3);
870: /* This should be enough for almost all matrices. If not, memory is reallocated later. */
871: c_maxmem = 8*(ai[am]+bi[bm]);
872: /* Step 2: Populate pattern for C */
873: PetscMalloc1(c_maxmem,&cj);
875: ci_nnz = 0;
876: ci[0] = 0;
877: worki_L1[0] = 0;
878: worki_L2[0] = 0;
879: for (i=0; i<am; i++) {
880: 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 */
881: const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */
882: rowsleft = anzi;
883: inputcol_L1 = acol;
884: L2_nnz = 0;
885: L2_nrows = 1; /* Number of rows to be merged on Level 3. output of L3 already exists -> initial value 1 */
886: worki_L2[1] = 0;
887: outputi_nnz = 0;
889: /* If the number of indices in C so far + the max number of columns in the next row > c_maxmem -> allocate more memory */
890: while (ci_nnz+a_maxrownnz > c_maxmem) {
891: c_maxmem *= 2;
892: ndouble++;
893: PetscRealloc(sizeof(PetscInt)*c_maxmem,&cj);
894: }
896: while (rowsleft) {
897: L1_rowsleft = PetscMin(64, rowsleft); /* In the inner loop max 64 rows of B can be merged */
898: L1_nrows = 0;
899: L1_nnz = 0;
900: inputcol = inputcol_L1;
901: inputi = bi;
902: inputj = bj;
904: /* The following macro is used to specialize for small rows in A.
905: This helps with compiler unrolling, improving performance substantially.
906: Input: inputj inputi inputcol bn
907: Output: outputj outputi_nnz */
908: #define MatMatMultSymbolic_RowMergeMacro(ANNZ) \
909: window_min = bn; \
910: outputi_nnz = 0; \
911: for (k=0; k<ANNZ; ++k) { \
912: brow_ptr[k] = inputj + inputi[inputcol[k]]; \
913: brow_end[k] = inputj + inputi[inputcol[k]+1]; \
914: window[k] = (brow_ptr[k] != brow_end[k]) ? *brow_ptr[k] : bn; \
915: window_min = PetscMin(window[k], window_min); \
916: } \
917: while (window_min < bn) { \
918: outputj[outputi_nnz++] = window_min; \
919: /* advance front and compute new minimum */ \
920: old_window_min = window_min; \
921: window_min = bn; \
922: for (k=0; k<ANNZ; ++k) { \
923: if (window[k] == old_window_min) { \
924: brow_ptr[k]++; \
925: window[k] = (brow_ptr[k] != brow_end[k]) ? *brow_ptr[k] : bn; \
926: } \
927: window_min = PetscMin(window[k], window_min); \
928: } \
929: }
931: /************** L E V E L 1 ***************/
932: /* Merge up to 8 rows of B to L1 work array*/
933: while (L1_rowsleft) {
934: outputi_nnz = 0;
935: if (anzi > 8) outputj = workj_L1 + L1_nnz; /* Level 1 rowmerge*/
936: else outputj = cj + ci_nnz; /* Merge directly to C */
938: switch (L1_rowsleft) {
939: case 1: brow_ptr[0] = inputj + inputi[inputcol[0]];
940: brow_end[0] = inputj + inputi[inputcol[0]+1];
941: for (; brow_ptr[0] != brow_end[0]; ++brow_ptr[0]) outputj[outputi_nnz++] = *brow_ptr[0]; /* copy row in b over */
942: inputcol += L1_rowsleft;
943: rowsleft -= L1_rowsleft;
944: L1_rowsleft = 0;
945: break;
946: case 2: MatMatMultSymbolic_RowMergeMacro(2);
947: inputcol += L1_rowsleft;
948: rowsleft -= L1_rowsleft;
949: L1_rowsleft = 0;
950: break;
951: case 3: MatMatMultSymbolic_RowMergeMacro(3);
952: inputcol += L1_rowsleft;
953: rowsleft -= L1_rowsleft;
954: L1_rowsleft = 0;
955: break;
956: case 4: MatMatMultSymbolic_RowMergeMacro(4);
957: inputcol += L1_rowsleft;
958: rowsleft -= L1_rowsleft;
959: L1_rowsleft = 0;
960: break;
961: case 5: MatMatMultSymbolic_RowMergeMacro(5);
962: inputcol += L1_rowsleft;
963: rowsleft -= L1_rowsleft;
964: L1_rowsleft = 0;
965: break;
966: case 6: MatMatMultSymbolic_RowMergeMacro(6);
967: inputcol += L1_rowsleft;
968: rowsleft -= L1_rowsleft;
969: L1_rowsleft = 0;
970: break;
971: case 7: MatMatMultSymbolic_RowMergeMacro(7);
972: inputcol += L1_rowsleft;
973: rowsleft -= L1_rowsleft;
974: L1_rowsleft = 0;
975: break;
976: default: MatMatMultSymbolic_RowMergeMacro(8);
977: inputcol += 8;
978: rowsleft -= 8;
979: L1_rowsleft -= 8;
980: break;
981: }
982: inputcol_L1 = inputcol;
983: L1_nnz += outputi_nnz;
984: worki_L1[++L1_nrows] = L1_nnz;
985: }
987: /********************** L E V E L 2 ************************/
988: /* Merge from L1 work array to either C or to L2 work array */
989: if (anzi > 8) {
990: inputi = worki_L1;
991: inputj = workj_L1;
992: inputcol = workcol;
993: outputi_nnz = 0;
995: if (anzi <= 64) outputj = cj + ci_nnz; /* Merge from L1 work array to C */
996: else outputj = workj_L2 + L2_nnz; /* Merge from L1 work array to L2 work array */
998: switch (L1_nrows) {
999: case 1: brow_ptr[0] = inputj + inputi[inputcol[0]];
1000: brow_end[0] = inputj + inputi[inputcol[0]+1];
1001: for (; brow_ptr[0] != brow_end[0]; ++brow_ptr[0]) outputj[outputi_nnz++] = *brow_ptr[0]; /* copy row in b over */
1002: break;
1003: case 2: MatMatMultSymbolic_RowMergeMacro(2); break;
1004: case 3: MatMatMultSymbolic_RowMergeMacro(3); break;
1005: case 4: MatMatMultSymbolic_RowMergeMacro(4); break;
1006: case 5: MatMatMultSymbolic_RowMergeMacro(5); break;
1007: case 6: MatMatMultSymbolic_RowMergeMacro(6); break;
1008: case 7: MatMatMultSymbolic_RowMergeMacro(7); break;
1009: case 8: MatMatMultSymbolic_RowMergeMacro(8); break;
1010: default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatMatMult logic error: Not merging 1-8 rows from L1 work array!");
1011: }
1012: L2_nnz += outputi_nnz;
1013: worki_L2[++L2_nrows] = L2_nnz;
1015: /************************ L E V E L 3 **********************/
1016: /* Merge from L2 work array to either C or to L2 work array */
1017: if (anzi > 64 && (L2_nrows == 8 || rowsleft == 0)) {
1018: inputi = worki_L2;
1019: inputj = workj_L2;
1020: inputcol = workcol;
1021: outputi_nnz = 0;
1022: if (rowsleft) outputj = workj_L3;
1023: else outputj = cj + ci_nnz;
1024: switch (L2_nrows) {
1025: case 1: brow_ptr[0] = inputj + inputi[inputcol[0]];
1026: brow_end[0] = inputj + inputi[inputcol[0]+1];
1027: for (; brow_ptr[0] != brow_end[0]; ++brow_ptr[0]) outputj[outputi_nnz++] = *brow_ptr[0]; /* copy row in b over */
1028: break;
1029: case 2: MatMatMultSymbolic_RowMergeMacro(2); break;
1030: case 3: MatMatMultSymbolic_RowMergeMacro(3); break;
1031: case 4: MatMatMultSymbolic_RowMergeMacro(4); break;
1032: case 5: MatMatMultSymbolic_RowMergeMacro(5); break;
1033: case 6: MatMatMultSymbolic_RowMergeMacro(6); break;
1034: case 7: MatMatMultSymbolic_RowMergeMacro(7); break;
1035: case 8: MatMatMultSymbolic_RowMergeMacro(8); break;
1036: default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatMatMult logic error: Not merging 1-8 rows from L2 work array!");
1037: }
1038: L2_nrows = 1;
1039: L2_nnz = outputi_nnz;
1040: worki_L2[1] = outputi_nnz;
1041: /* Copy to workj_L2 */
1042: if (rowsleft) {
1043: for (k=0; k<outputi_nnz; ++k) workj_L2[k] = outputj[k];
1044: }
1045: }
1046: }
1047: } /* while (rowsleft) */
1048: #undef MatMatMultSymbolic_RowMergeMacro
1050: /* terminate current row */
1051: ci_nnz += outputi_nnz;
1052: ci[i+1] = ci_nnz;
1053: }
1055: /* Step 3: Create the new symbolic matrix */
1056: MatSetSeqAIJWithArrays_private(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,((PetscObject)A)->type_name,C);
1057: MatSetBlockSizesFromMats(C,A,B);
1059: /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
1060: /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
1061: c = (Mat_SeqAIJ*)(C->data);
1062: c->free_a = PETSC_TRUE;
1063: c->free_ij = PETSC_TRUE;
1064: c->nonew = 0;
1066: C->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Sorted;
1068: /* set MatInfo */
1069: afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
1070: if (afill < 1.0) afill = 1.0;
1071: c->maxnz = ci[am];
1072: c->nz = ci[am];
1073: C->info.mallocs = ndouble;
1074: C->info.fill_ratio_given = fill;
1075: C->info.fill_ratio_needed = afill;
1077: #if defined(PETSC_USE_INFO)
1078: if (ci[am]) {
1079: PetscInfo3(C,"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);
1080: PetscInfo1(C,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);
1081: } else {
1082: PetscInfo(C,"Empty matrix product\n");
1083: }
1084: #endif
1086: /* Step 4: Free temporary work areas */
1087: PetscFree(workj_L1);
1088: PetscFree(workj_L2);
1089: PetscFree(workj_L3);
1090: return(0);
1091: }
1093: /* concatenate unique entries and then sort */
1094: PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Sorted(Mat A,Mat B,PetscReal fill,Mat C)
1095: {
1097: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
1098: const PetscInt *ai = a->i,*bi=b->i,*aj=a->j,*bj=b->j;
1099: PetscInt *ci,*cj,bcol;
1100: PetscInt am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
1101: PetscReal afill;
1102: PetscInt i,j,ndouble = 0;
1103: PetscSegBuffer seg,segrow;
1104: char *seen;
1107: PetscMalloc1(am+1,&ci);
1108: ci[0] = 0;
1110: /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
1111: PetscSegBufferCreate(sizeof(PetscInt),(PetscInt)(fill*(ai[am]+bi[bm])),&seg);
1112: PetscSegBufferCreate(sizeof(PetscInt),100,&segrow);
1113: PetscCalloc1(bn,&seen);
1115: /* Determine ci and cj */
1116: for (i=0; i<am; i++) {
1117: 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 */
1118: const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */
1119: PetscInt packlen = 0,*PETSC_RESTRICT crow;
1121: /* Pack segrow */
1122: for (j=0; j<anzi; j++) {
1123: PetscInt brow = acol[j],bjstart = bi[brow],bjend = bi[brow+1],k;
1124: for (k=bjstart; k<bjend; k++) {
1125: bcol = bj[k];
1126: if (!seen[bcol]) { /* new entry */
1127: PetscInt *PETSC_RESTRICT slot;
1128: PetscSegBufferGetInts(segrow,1,&slot);
1129: *slot = bcol;
1130: seen[bcol] = 1;
1131: packlen++;
1132: }
1133: }
1134: }
1136: /* Check i-th diagonal entry */
1137: if (C->force_diagonals && !seen[i]) {
1138: PetscInt *PETSC_RESTRICT slot;
1139: PetscSegBufferGetInts(segrow,1,&slot);
1140: *slot = i;
1141: seen[i] = 1;
1142: packlen++;
1143: }
1145: PetscSegBufferGetInts(seg,packlen,&crow);
1146: PetscSegBufferExtractTo(segrow,crow);
1147: PetscSortInt(packlen,crow);
1148: ci[i+1] = ci[i] + packlen;
1149: for (j=0; j<packlen; j++) seen[crow[j]] = 0;
1150: }
1151: PetscSegBufferDestroy(&segrow);
1152: PetscFree(seen);
1154: /* Column indices are in the segmented buffer */
1155: PetscSegBufferExtractAlloc(seg,&cj);
1156: PetscSegBufferDestroy(&seg);
1158: /* put together the new symbolic matrix */
1159: MatSetSeqAIJWithArrays_private(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,((PetscObject)A)->type_name,C);
1160: MatSetBlockSizesFromMats(C,A,B);
1162: /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
1163: /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
1164: c = (Mat_SeqAIJ*)(C->data);
1165: c->free_a = PETSC_TRUE;
1166: c->free_ij = PETSC_TRUE;
1167: c->nonew = 0;
1169: C->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Sorted;
1171: /* set MatInfo */
1172: afill = (PetscReal)ci[am]/PetscMax(ai[am]+bi[bm],1) + 1.e-5;
1173: if (afill < 1.0) afill = 1.0;
1174: c->maxnz = ci[am];
1175: c->nz = ci[am];
1176: C->info.mallocs = ndouble;
1177: C->info.fill_ratio_given = fill;
1178: C->info.fill_ratio_needed = afill;
1180: #if defined(PETSC_USE_INFO)
1181: if (ci[am]) {
1182: PetscInfo3(C,"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);
1183: PetscInfo1(C,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);
1184: } else {
1185: PetscInfo(C,"Empty matrix product\n");
1186: }
1187: #endif
1188: return(0);
1189: }
1191: PetscErrorCode MatDestroy_SeqAIJ_MatMatMultTrans(void *data)
1192: {
1193: PetscErrorCode ierr;
1194: Mat_MatMatTransMult *abt=(Mat_MatMatTransMult *)data;
1197: MatTransposeColoringDestroy(&abt->matcoloring);
1198: MatDestroy(&abt->Bt_den);
1199: MatDestroy(&abt->ABt_den);
1200: PetscFree(abt);
1201: return(0);
1202: }
1204: PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat C)
1205: {
1206: PetscErrorCode ierr;
1207: Mat Bt;
1208: PetscInt *bti,*btj;
1209: Mat_MatMatTransMult *abt;
1210: Mat_Product *product = C->product;
1211: char *alg;
1214: if (!product) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing product struct");
1215: if (product->data) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Extra product struct not empty");
1217: /* create symbolic Bt */
1218: MatGetSymbolicTranspose_SeqAIJ(B,&bti,&btj);
1219: MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,B->cmap->n,B->rmap->n,bti,btj,NULL,&Bt);
1220: MatSetBlockSizes(Bt,PetscAbs(A->cmap->bs),PetscAbs(B->cmap->bs));
1221: MatSetType(Bt,((PetscObject)A)->type_name);
1223: /* get symbolic C=A*Bt */
1224: PetscStrallocpy(product->alg,&alg);
1225: MatProductSetAlgorithm(C,"sorted"); /* set algorithm for C = A*Bt */
1226: MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,Bt,fill,C);
1227: MatProductSetAlgorithm(C,alg); /* resume original algorithm for ABt product */
1228: PetscFree(alg);
1230: /* create a supporting struct for reuse intermediate dense matrices with matcoloring */
1231: PetscNew(&abt);
1233: product->data = abt;
1234: product->destroy = MatDestroy_SeqAIJ_MatMatMultTrans;
1236: C->ops->mattransposemultnumeric = MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ;
1238: abt->usecoloring = PETSC_FALSE;
1239: PetscStrcmp(product->alg,"color",&abt->usecoloring);
1240: if (abt->usecoloring) {
1241: /* Create MatTransposeColoring from symbolic C=A*B^T */
1242: MatTransposeColoring matcoloring;
1243: MatColoring coloring;
1244: ISColoring iscoloring;
1245: Mat Bt_dense,C_dense;
1247: /* inode causes memory problem */
1248: MatSetOption(C,MAT_USE_INODES,PETSC_FALSE);
1250: MatColoringCreate(C,&coloring);
1251: MatColoringSetDistance(coloring,2);
1252: MatColoringSetType(coloring,MATCOLORINGSL);
1253: MatColoringSetFromOptions(coloring);
1254: MatColoringApply(coloring,&iscoloring);
1255: MatColoringDestroy(&coloring);
1256: MatTransposeColoringCreate(C,iscoloring,&matcoloring);
1258: abt->matcoloring = matcoloring;
1260: ISColoringDestroy(&iscoloring);
1262: /* Create Bt_dense and C_dense = A*Bt_dense */
1263: MatCreate(PETSC_COMM_SELF,&Bt_dense);
1264: MatSetSizes(Bt_dense,A->cmap->n,matcoloring->ncolors,A->cmap->n,matcoloring->ncolors);
1265: MatSetType(Bt_dense,MATSEQDENSE);
1266: MatSeqDenseSetPreallocation(Bt_dense,NULL);
1268: Bt_dense->assembled = PETSC_TRUE;
1269: abt->Bt_den = Bt_dense;
1271: MatCreate(PETSC_COMM_SELF,&C_dense);
1272: MatSetSizes(C_dense,A->rmap->n,matcoloring->ncolors,A->rmap->n,matcoloring->ncolors);
1273: MatSetType(C_dense,MATSEQDENSE);
1274: MatSeqDenseSetPreallocation(C_dense,NULL);
1276: Bt_dense->assembled = PETSC_TRUE;
1277: abt->ABt_den = C_dense;
1279: #if defined(PETSC_USE_INFO)
1280: {
1281: Mat_SeqAIJ *c = (Mat_SeqAIJ*)C->data;
1282: 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));
1283: }
1284: #endif
1285: }
1286: /* clean up */
1287: MatDestroy(&Bt);
1288: MatRestoreSymbolicTranspose_SeqAIJ(B,&bti,&btj);
1289: return(0);
1290: }
1292: PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
1293: {
1294: PetscErrorCode ierr;
1295: Mat_SeqAIJ *a =(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data;
1296: PetscInt *ai =a->i,*aj=a->j,*bi=b->i,*bj=b->j,anzi,bnzj,nexta,nextb,*acol,*bcol,brow;
1297: PetscInt cm =C->rmap->n,*ci=c->i,*cj=c->j,i,j,cnzi,*ccol;
1298: PetscLogDouble flops=0.0;
1299: MatScalar *aa =a->a,*aval,*ba=b->a,*bval,*ca,*cval;
1300: Mat_MatMatTransMult *abt;
1301: Mat_Product *product = C->product;
1304: if (!product) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing product struct");
1305: abt = (Mat_MatMatTransMult *)product->data;
1306: if (!abt) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing product struct");
1307: /* clear old values in C */
1308: if (!c->a) {
1309: PetscCalloc1(ci[cm]+1,&ca);
1310: c->a = ca;
1311: c->free_a = PETSC_TRUE;
1312: } else {
1313: ca = c->a;
1314: PetscArrayzero(ca,ci[cm]+1);
1315: }
1317: if (abt->usecoloring) {
1318: MatTransposeColoring matcoloring = abt->matcoloring;
1319: Mat Bt_dense,C_dense = abt->ABt_den;
1321: /* Get Bt_dense by Apply MatTransposeColoring to B */
1322: Bt_dense = abt->Bt_den;
1323: MatTransColoringApplySpToDen(matcoloring,B,Bt_dense);
1325: /* C_dense = A*Bt_dense */
1326: MatMatMultNumeric_SeqAIJ_SeqDense(A,Bt_dense,C_dense);
1328: /* Recover C from C_dense */
1329: MatTransColoringApplyDenToSp(matcoloring,C_dense,C);
1330: return(0);
1331: }
1333: for (i=0; i<cm; i++) {
1334: anzi = ai[i+1] - ai[i];
1335: acol = aj + ai[i];
1336: aval = aa + ai[i];
1337: cnzi = ci[i+1] - ci[i];
1338: ccol = cj + ci[i];
1339: cval = ca + ci[i];
1340: for (j=0; j<cnzi; j++) {
1341: brow = ccol[j];
1342: bnzj = bi[brow+1] - bi[brow];
1343: bcol = bj + bi[brow];
1344: bval = ba + bi[brow];
1346: /* perform sparse inner-product c(i,j)=A[i,:]*B[j,:]^T */
1347: nexta = 0; nextb = 0;
1348: while (nexta<anzi && nextb<bnzj) {
1349: while (nexta < anzi && acol[nexta] < bcol[nextb]) nexta++;
1350: if (nexta == anzi) break;
1351: while (nextb < bnzj && acol[nexta] > bcol[nextb]) nextb++;
1352: if (nextb == bnzj) break;
1353: if (acol[nexta] == bcol[nextb]) {
1354: cval[j] += aval[nexta]*bval[nextb];
1355: nexta++; nextb++;
1356: flops += 2;
1357: }
1358: }
1359: }
1360: }
1361: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
1362: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
1363: PetscLogFlops(flops);
1364: return(0);
1365: }
1367: PetscErrorCode MatDestroy_SeqAIJ_MatTransMatMult(void *data)
1368: {
1369: PetscErrorCode ierr;
1370: Mat_MatTransMatMult *atb = (Mat_MatTransMatMult*)data;
1373: MatDestroy(&atb->At);
1374: if (atb->destroy) {
1375: (*atb->destroy)(atb->data);
1376: }
1377: PetscFree(atb);
1378: return(0);
1379: }
1381: PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat C)
1382: {
1384: Mat At = NULL;
1385: PetscInt *ati,*atj;
1386: Mat_Product *product = C->product;
1387: PetscBool flg,def,square;
1390: MatCheckProduct(C,4);
1391: square = (PetscBool)(A == B && A->symmetric && A->symmetric_set);
1392: /* outerproduct */
1393: PetscStrcmp(product->alg,"outerproduct",&flg);
1394: if (flg) {
1395: /* create symbolic At */
1396: if (!square) {
1397: MatGetSymbolicTranspose_SeqAIJ(A,&ati,&atj);
1398: MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,A->cmap->n,A->rmap->n,ati,atj,NULL,&At);
1399: MatSetBlockSizes(At,PetscAbs(A->cmap->bs),PetscAbs(B->cmap->bs));
1400: MatSetType(At,((PetscObject)A)->type_name);
1401: }
1402: /* get symbolic C=At*B */
1403: MatProductSetAlgorithm(C,"sorted");
1404: MatMatMultSymbolic_SeqAIJ_SeqAIJ(square ? A : At,B,fill,C);
1406: /* clean up */
1407: if (!square) {
1408: MatDestroy(&At);
1409: MatRestoreSymbolicTranspose_SeqAIJ(A,&ati,&atj);
1410: }
1412: C->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ; /* outerproduct */
1413: MatProductSetAlgorithm(C,"outerproduct");
1414: return(0);
1415: }
1417: /* matmatmult */
1418: PetscStrcmp(product->alg,"default",&def);
1419: PetscStrcmp(product->alg,"at*b",&flg);
1420: if (flg || def) {
1421: Mat_MatTransMatMult *atb;
1423: if (product->data) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Extra product struct not empty");
1424: PetscNew(&atb);
1425: if (!square) {
1426: MatTranspose_SeqAIJ(A,MAT_INITIAL_MATRIX,&At);
1427: }
1428: MatProductSetAlgorithm(C,"sorted");
1429: MatMatMultSymbolic_SeqAIJ_SeqAIJ(square ? A : At,B,fill,C);
1430: MatProductSetAlgorithm(C,"at*b");
1431: product->data = atb;
1432: product->destroy = MatDestroy_SeqAIJ_MatTransMatMult;
1433: atb->At = At;
1434: atb->updateAt = PETSC_FALSE; /* because At is computed here */
1436: C->ops->mattransposemultnumeric = NULL; /* see MatProductNumeric_AtB_SeqAIJ_SeqAIJ */
1437: return(0);
1438: }
1440: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Mat Product Algorithm is not supported");
1441: }
1443: PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
1444: {
1446: Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data;
1447: PetscInt am=A->rmap->n,anzi,*ai=a->i,*aj=a->j,*bi=b->i,*bj,bnzi,nextb;
1448: PetscInt cm=C->rmap->n,*ci=c->i,*cj=c->j,crow,*cjj,i,j,k;
1449: PetscLogDouble flops=0.0;
1450: MatScalar *aa=a->a,*ba,*ca,*caj;
1453: if (!c->a) {
1454: PetscCalloc1(ci[cm]+1,&ca);
1456: c->a = ca;
1457: c->free_a = PETSC_TRUE;
1458: } else {
1459: ca = c->a;
1460: PetscArrayzero(ca,ci[cm]);
1461: }
1463: /* compute A^T*B using outer product (A^T)[:,i]*B[i,:] */
1464: for (i=0; i<am; i++) {
1465: bj = b->j + bi[i];
1466: ba = b->a + bi[i];
1467: bnzi = bi[i+1] - bi[i];
1468: anzi = ai[i+1] - ai[i];
1469: for (j=0; j<anzi; j++) {
1470: nextb = 0;
1471: crow = *aj++;
1472: cjj = cj + ci[crow];
1473: caj = ca + ci[crow];
1474: /* perform sparse axpy operation. Note cjj includes bj. */
1475: for (k=0; nextb<bnzi; k++) {
1476: if (cjj[k] == *(bj+nextb)) { /* ccol == bcol */
1477: caj[k] += (*aa)*(*(ba+nextb));
1478: nextb++;
1479: }
1480: }
1481: flops += 2*bnzi;
1482: aa++;
1483: }
1484: }
1486: /* Assemble the final matrix and clean up */
1487: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
1488: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
1489: PetscLogFlops(flops);
1490: return(0);
1491: }
1493: PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqDense(Mat A,Mat B,PetscReal fill,Mat C)
1494: {
1498: MatMatMultSymbolic_SeqDense_SeqDense(A,B,0.0,C);
1499: C->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqDense;
1500: return(0);
1501: }
1503: PETSC_INTERN PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat A,Mat B,Mat C,const PetscBool add)
1504: {
1505: Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data;
1506: Mat_SeqDense *bd=(Mat_SeqDense*)B->data;
1507: Mat_SeqDense *cd=(Mat_SeqDense*)C->data;
1508: PetscErrorCode ierr;
1509: PetscScalar *c,r1,r2,r3,r4,*c1,*c2,*c3,*c4;
1510: const PetscScalar *aa,*b,*b1,*b2,*b3,*b4,*av;
1511: const PetscInt *aj;
1512: PetscInt cm=C->rmap->n,cn=B->cmap->n,bm=bd->lda,am=A->rmap->n;
1513: PetscInt clda=cd->lda;
1514: PetscInt am4=4*clda,bm4=4*bm,col,i,j,n;
1517: if (!cm || !cn) return(0);
1518: MatSeqAIJGetArrayRead(A,&av);
1519: if (add) {
1520: MatDenseGetArray(C,&c);
1521: } else {
1522: MatDenseGetArrayWrite(C,&c);
1523: }
1524: MatDenseGetArrayRead(B,&b);
1525: b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;
1526: c1 = c; c2 = c1 + clda; c3 = c2 + clda; c4 = c3 + clda;
1527: for (col=0; col<(cn/4)*4; col += 4) { /* over columns of C */
1528: for (i=0; i<am; i++) { /* over rows of A in those columns */
1529: r1 = r2 = r3 = r4 = 0.0;
1530: n = a->i[i+1] - a->i[i];
1531: aj = a->j + a->i[i];
1532: aa = av + a->i[i];
1533: for (j=0; j<n; j++) {
1534: const PetscScalar aatmp = aa[j];
1535: const PetscInt ajtmp = aj[j];
1536: r1 += aatmp*b1[ajtmp];
1537: r2 += aatmp*b2[ajtmp];
1538: r3 += aatmp*b3[ajtmp];
1539: r4 += aatmp*b4[ajtmp];
1540: }
1541: if (add) {
1542: c1[i] += r1;
1543: c2[i] += r2;
1544: c3[i] += r3;
1545: c4[i] += r4;
1546: } else {
1547: c1[i] = r1;
1548: c2[i] = r2;
1549: c3[i] = r3;
1550: c4[i] = r4;
1551: }
1552: }
1553: b1 += bm4; b2 += bm4; b3 += bm4; b4 += bm4;
1554: c1 += am4; c2 += am4; c3 += am4; c4 += am4;
1555: }
1556: /* process remaining columns */
1557: if (col != cn) {
1558: PetscInt rc = cn-col;
1560: if (rc == 1) {
1561: for (i=0; i<am; i++) {
1562: r1 = 0.0;
1563: n = a->i[i+1] - a->i[i];
1564: aj = a->j + a->i[i];
1565: aa = av + a->i[i];
1566: for (j=0; j<n; j++) r1 += aa[j]*b1[aj[j]];
1567: if (add) c1[i] += r1;
1568: else c1[i] = r1;
1569: }
1570: } else if (rc == 2) {
1571: for (i=0; i<am; i++) {
1572: r1 = r2 = 0.0;
1573: n = a->i[i+1] - a->i[i];
1574: aj = a->j + a->i[i];
1575: aa = av + a->i[i];
1576: for (j=0; j<n; j++) {
1577: const PetscScalar aatmp = aa[j];
1578: const PetscInt ajtmp = aj[j];
1579: r1 += aatmp*b1[ajtmp];
1580: r2 += aatmp*b2[ajtmp];
1581: }
1582: if (add) {
1583: c1[i] += r1;
1584: c2[i] += r2;
1585: } else {
1586: c1[i] = r1;
1587: c2[i] = r2;
1588: }
1589: }
1590: } else {
1591: for (i=0; i<am; i++) {
1592: r1 = r2 = r3 = 0.0;
1593: n = a->i[i+1] - a->i[i];
1594: aj = a->j + a->i[i];
1595: aa = av + a->i[i];
1596: for (j=0; j<n; j++) {
1597: const PetscScalar aatmp = aa[j];
1598: const PetscInt ajtmp = aj[j];
1599: r1 += aatmp*b1[ajtmp];
1600: r2 += aatmp*b2[ajtmp];
1601: r3 += aatmp*b3[ajtmp];
1602: }
1603: if (add) {
1604: c1[i] += r1;
1605: c2[i] += r2;
1606: c3[i] += r3;
1607: } else {
1608: c1[i] = r1;
1609: c2[i] = r2;
1610: c3[i] = r3;
1611: }
1612: }
1613: }
1614: }
1615: PetscLogFlops(cn*(2.0*a->nz));
1616: if (add) {
1617: MatDenseRestoreArray(C,&c);
1618: } else {
1619: MatDenseRestoreArrayWrite(C,&c);
1620: }
1621: MatDenseRestoreArrayRead(B,&b);
1622: MatSeqAIJRestoreArrayRead(A,&av);
1623: return(0);
1624: }
1626: PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
1627: {
1631: 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);
1632: 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);
1633: 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);
1635: MatMatMultNumericAdd_SeqAIJ_SeqDense(A,B,C,PETSC_FALSE);
1636: return(0);
1637: }
1639: /* ------------------------------------------------------- */
1640: static PetscErrorCode MatProductSetFromOptions_SeqAIJ_SeqDense_AB(Mat C)
1641: {
1643: C->ops->matmultsymbolic = MatMatMultSymbolic_SeqAIJ_SeqDense;
1644: C->ops->productsymbolic = MatProductSymbolic_AB;
1645: return(0);
1646: }
1648: PETSC_INTERN PetscErrorCode MatTMatTMultSymbolic_SeqAIJ_SeqDense(Mat,Mat,PetscReal,Mat);
1650: static PetscErrorCode MatProductSetFromOptions_SeqAIJ_SeqDense_AtB(Mat C)
1651: {
1653: C->ops->transposematmultsymbolic = MatTMatTMultSymbolic_SeqAIJ_SeqDense;
1654: C->ops->productsymbolic = MatProductSymbolic_AtB;
1655: return(0);
1656: }
1658: static PetscErrorCode MatProductSetFromOptions_SeqAIJ_SeqDense_ABt(Mat C)
1659: {
1661: C->ops->mattransposemultsymbolic = MatTMatTMultSymbolic_SeqAIJ_SeqDense;
1662: C->ops->productsymbolic = MatProductSymbolic_ABt;
1663: return(0);
1664: }
1666: PETSC_INTERN PetscErrorCode MatProductSetFromOptions_SeqAIJ_SeqDense(Mat C)
1667: {
1669: Mat_Product *product = C->product;
1672: switch (product->type) {
1673: case MATPRODUCT_AB:
1674: MatProductSetFromOptions_SeqAIJ_SeqDense_AB(C);
1675: break;
1676: case MATPRODUCT_AtB:
1677: MatProductSetFromOptions_SeqAIJ_SeqDense_AtB(C);
1678: break;
1679: case MATPRODUCT_ABt:
1680: MatProductSetFromOptions_SeqAIJ_SeqDense_ABt(C);
1681: break;
1682: default:
1683: break;
1684: }
1685: return(0);
1686: }
1687: /* ------------------------------------------------------- */
1688: static PetscErrorCode MatProductSetFromOptions_SeqXBAIJ_SeqDense_AB(Mat C)
1689: {
1691: Mat_Product *product = C->product;
1692: Mat A = product->A;
1693: PetscBool baij;
1696: PetscObjectTypeCompare((PetscObject)A,MATSEQBAIJ,&baij);
1697: if (!baij) { /* A is seqsbaij */
1698: PetscBool sbaij;
1699: PetscObjectTypeCompare((PetscObject)A,MATSEQSBAIJ,&sbaij);
1700: if (!sbaij) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_WRONGSTATE,"Mat must be either seqbaij or seqsbaij format");
1702: C->ops->matmultsymbolic = MatMatMultSymbolic_SeqSBAIJ_SeqDense;
1703: } else { /* A is seqbaij */
1704: C->ops->matmultsymbolic = MatMatMultSymbolic_SeqBAIJ_SeqDense;
1705: }
1707: C->ops->productsymbolic = MatProductSymbolic_AB;
1708: return(0);
1709: }
1711: PETSC_INTERN PetscErrorCode MatProductSetFromOptions_SeqXBAIJ_SeqDense(Mat C)
1712: {
1714: Mat_Product *product = C->product;
1717: MatCheckProduct(C,1);
1718: if (!product->A) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing A");
1719: if (product->type == MATPRODUCT_AB || (product->type == MATPRODUCT_AtB && product->A->symmetric)) {
1720: MatProductSetFromOptions_SeqXBAIJ_SeqDense_AB(C);
1721: }
1722: return(0);
1723: }
1725: /* ------------------------------------------------------- */
1726: static PetscErrorCode MatProductSetFromOptions_SeqDense_SeqAIJ_AB(Mat C)
1727: {
1729: C->ops->matmultsymbolic = MatMatMultSymbolic_SeqDense_SeqAIJ;
1730: C->ops->productsymbolic = MatProductSymbolic_AB;
1731: return(0);
1732: }
1734: PETSC_INTERN PetscErrorCode MatProductSetFromOptions_SeqDense_SeqAIJ(Mat C)
1735: {
1737: Mat_Product *product = C->product;
1740: if (product->type == MATPRODUCT_AB) {
1741: MatProductSetFromOptions_SeqDense_SeqAIJ_AB(C);
1742: }
1743: return(0);
1744: }
1745: /* ------------------------------------------------------- */
1747: PetscErrorCode MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring coloring,Mat B,Mat Btdense)
1748: {
1750: Mat_SeqAIJ *b = (Mat_SeqAIJ*)B->data;
1751: Mat_SeqDense *btdense = (Mat_SeqDense*)Btdense->data;
1752: PetscInt *bi = b->i,*bj=b->j;
1753: PetscInt m = Btdense->rmap->n,n=Btdense->cmap->n,j,k,l,col,anz,*btcol,brow,ncolumns;
1754: MatScalar *btval,*btval_den,*ba=b->a;
1755: PetscInt *columns=coloring->columns,*colorforcol=coloring->colorforcol,ncolors=coloring->ncolors;
1758: btval_den=btdense->v;
1759: PetscArrayzero(btval_den,m*n);
1760: for (k=0; k<ncolors; k++) {
1761: ncolumns = coloring->ncolumns[k];
1762: for (l=0; l<ncolumns; l++) { /* insert a row of B to a column of Btdense */
1763: col = *(columns + colorforcol[k] + l);
1764: btcol = bj + bi[col];
1765: btval = ba + bi[col];
1766: anz = bi[col+1] - bi[col];
1767: for (j=0; j<anz; j++) {
1768: brow = btcol[j];
1769: btval_den[brow] = btval[j];
1770: }
1771: }
1772: btval_den += m;
1773: }
1774: return(0);
1775: }
1777: PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring matcoloring,Mat Cden,Mat Csp)
1778: {
1779: PetscErrorCode ierr;
1780: Mat_SeqAIJ *csp = (Mat_SeqAIJ*)Csp->data;
1781: const PetscScalar *ca_den,*ca_den_ptr;
1782: PetscScalar *ca=csp->a;
1783: PetscInt k,l,m=Cden->rmap->n,ncolors=matcoloring->ncolors;
1784: PetscInt brows=matcoloring->brows,*den2sp=matcoloring->den2sp;
1785: PetscInt nrows,*row,*idx;
1786: PetscInt *rows=matcoloring->rows,*colorforrow=matcoloring->colorforrow;
1789: MatDenseGetArrayRead(Cden,&ca_den);
1791: if (brows > 0) {
1792: PetscInt *lstart,row_end,row_start;
1793: lstart = matcoloring->lstart;
1794: PetscArrayzero(lstart,ncolors);
1796: row_end = brows;
1797: if (row_end > m) row_end = m;
1798: for (row_start=0; row_start<m; row_start+=brows) { /* loop over row blocks of Csp */
1799: ca_den_ptr = ca_den;
1800: for (k=0; k<ncolors; k++) { /* loop over colors (columns of Cden) */
1801: nrows = matcoloring->nrows[k];
1802: row = rows + colorforrow[k];
1803: idx = den2sp + colorforrow[k];
1804: for (l=lstart[k]; l<nrows; l++) {
1805: if (row[l] >= row_end) {
1806: lstart[k] = l;
1807: break;
1808: } else {
1809: ca[idx[l]] = ca_den_ptr[row[l]];
1810: }
1811: }
1812: ca_den_ptr += m;
1813: }
1814: row_end += brows;
1815: if (row_end > m) row_end = m;
1816: }
1817: } else { /* non-blocked impl: loop over columns of Csp - slow if Csp is large */
1818: ca_den_ptr = ca_den;
1819: for (k=0; k<ncolors; k++) {
1820: nrows = matcoloring->nrows[k];
1821: row = rows + colorforrow[k];
1822: idx = den2sp + colorforrow[k];
1823: for (l=0; l<nrows; l++) {
1824: ca[idx[l]] = ca_den_ptr[row[l]];
1825: }
1826: ca_den_ptr += m;
1827: }
1828: }
1830: MatDenseRestoreArrayRead(Cden,&ca_den);
1831: #if defined(PETSC_USE_INFO)
1832: if (matcoloring->brows > 0) {
1833: PetscInfo1(Csp,"Loop over %D row blocks for den2sp\n",brows);
1834: } else {
1835: PetscInfo(Csp,"Loop over colors/columns of Cden, inefficient for large sparse matrix product \n");
1836: }
1837: #endif
1838: return(0);
1839: }
1841: PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat mat,ISColoring iscoloring,MatTransposeColoring c)
1842: {
1844: PetscInt i,n,nrows,Nbs,j,k,m,ncols,col,cm;
1845: const PetscInt *is,*ci,*cj,*row_idx;
1846: PetscInt nis = iscoloring->n,*rowhit,bs = 1;
1847: IS *isa;
1848: Mat_SeqAIJ *csp = (Mat_SeqAIJ*)mat->data;
1849: PetscInt *colorforrow,*rows,*rows_i,*idxhit,*spidx,*den2sp,*den2sp_i;
1850: PetscInt *colorforcol,*columns,*columns_i,brows;
1851: PetscBool flg;
1854: ISColoringGetIS(iscoloring,PETSC_USE_POINTER,PETSC_IGNORE,&isa);
1856: /* bs >1 is not being tested yet! */
1857: Nbs = mat->cmap->N/bs;
1858: c->M = mat->rmap->N/bs; /* set total rows, columns and local rows */
1859: c->N = Nbs;
1860: c->m = c->M;
1861: c->rstart = 0;
1862: c->brows = 100;
1864: c->ncolors = nis;
1865: PetscMalloc3(nis,&c->ncolumns,nis,&c->nrows,nis+1,&colorforrow);
1866: PetscMalloc1(csp->nz+1,&rows);
1867: PetscMalloc1(csp->nz+1,&den2sp);
1869: brows = c->brows;
1870: PetscOptionsGetInt(NULL,NULL,"-matden2sp_brows",&brows,&flg);
1871: if (flg) c->brows = brows;
1872: if (brows > 0) {
1873: PetscMalloc1(nis+1,&c->lstart);
1874: }
1876: colorforrow[0] = 0;
1877: rows_i = rows;
1878: den2sp_i = den2sp;
1880: PetscMalloc1(nis+1,&colorforcol);
1881: PetscMalloc1(Nbs+1,&columns);
1883: colorforcol[0] = 0;
1884: columns_i = columns;
1886: /* get column-wise storage of mat */
1887: MatGetColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);
1889: cm = c->m;
1890: PetscMalloc1(cm+1,&rowhit);
1891: PetscMalloc1(cm+1,&idxhit);
1892: for (i=0; i<nis; i++) { /* loop over color */
1893: ISGetLocalSize(isa[i],&n);
1894: ISGetIndices(isa[i],&is);
1896: c->ncolumns[i] = n;
1897: if (n) {
1898: PetscArraycpy(columns_i,is,n);
1899: }
1900: colorforcol[i+1] = colorforcol[i] + n;
1901: columns_i += n;
1903: /* fast, crude version requires O(N*N) work */
1904: PetscArrayzero(rowhit,cm);
1906: for (j=0; j<n; j++) { /* loop over columns*/
1907: col = is[j];
1908: row_idx = cj + ci[col];
1909: m = ci[col+1] - ci[col];
1910: for (k=0; k<m; k++) { /* loop over columns marking them in rowhit */
1911: idxhit[*row_idx] = spidx[ci[col] + k];
1912: rowhit[*row_idx++] = col + 1;
1913: }
1914: }
1915: /* count the number of hits */
1916: nrows = 0;
1917: for (j=0; j<cm; j++) {
1918: if (rowhit[j]) nrows++;
1919: }
1920: c->nrows[i] = nrows;
1921: colorforrow[i+1] = colorforrow[i] + nrows;
1923: nrows = 0;
1924: for (j=0; j<cm; j++) { /* loop over rows */
1925: if (rowhit[j]) {
1926: rows_i[nrows] = j;
1927: den2sp_i[nrows] = idxhit[j];
1928: nrows++;
1929: }
1930: }
1931: den2sp_i += nrows;
1933: ISRestoreIndices(isa[i],&is);
1934: rows_i += nrows;
1935: }
1936: MatRestoreColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);
1937: PetscFree(rowhit);
1938: ISColoringRestoreIS(iscoloring,PETSC_USE_POINTER,&isa);
1939: if (csp->nz != colorforrow[nis]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"csp->nz %d != colorforrow[nis] %d",csp->nz,colorforrow[nis]);
1941: c->colorforrow = colorforrow;
1942: c->rows = rows;
1943: c->den2sp = den2sp;
1944: c->colorforcol = colorforcol;
1945: c->columns = columns;
1947: PetscFree(idxhit);
1948: return(0);
1949: }
1951: /* --------------------------------------------------------------- */
1952: static PetscErrorCode MatProductNumeric_AtB_SeqAIJ_SeqAIJ(Mat C)
1953: {
1955: Mat_Product *product = C->product;
1956: Mat A=product->A,B=product->B;
1959: if (C->ops->mattransposemultnumeric) {
1960: /* Alg: "outerproduct" */
1961: (*C->ops->mattransposemultnumeric)(A,B,C);
1962: } else {
1963: /* Alg: "matmatmult" -- C = At*B */
1964: Mat_MatTransMatMult *atb = (Mat_MatTransMatMult *)product->data;
1965: Mat At;
1967: if (!atb) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing product struct");
1968: At = atb->At;
1969: if (atb->updateAt && At) { /* At is computed in MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ() */
1970: MatTranspose_SeqAIJ(A,MAT_REUSE_MATRIX,&At);
1971: }
1972: MatMatMultNumeric_SeqAIJ_SeqAIJ(At ? At : A,B,C);
1973: atb->updateAt = PETSC_TRUE;
1974: }
1975: return(0);
1976: }
1978: static PetscErrorCode MatProductSymbolic_AtB_SeqAIJ_SeqAIJ(Mat C)
1979: {
1981: Mat_Product *product = C->product;
1982: Mat A=product->A,B=product->B;
1983: PetscReal fill=product->fill;
1986: MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);
1988: C->ops->productnumeric = MatProductNumeric_AtB_SeqAIJ_SeqAIJ;
1989: return(0);
1990: }
1992: /* --------------------------------------------------------------- */
1993: static PetscErrorCode MatProductSetFromOptions_SeqAIJ_AB(Mat C)
1994: {
1996: Mat_Product *product = C->product;
1997: PetscInt alg = 0; /* default algorithm */
1998: PetscBool flg = PETSC_FALSE;
1999: #if !defined(PETSC_HAVE_HYPRE)
2000: const char *algTypes[7] = {"sorted","scalable","scalable_fast","heap","btheap","llcondensed","rowmerge"};
2001: PetscInt nalg = 7;
2002: #else
2003: const char *algTypes[8] = {"sorted","scalable","scalable_fast","heap","btheap","llcondensed","rowmerge","hypre"};
2004: PetscInt nalg = 8;
2005: #endif
2008: /* Set default algorithm */
2009: PetscStrcmp(C->product->alg,"default",&flg);
2010: if (flg) {
2011: MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);
2012: }
2014: /* Get runtime option */
2015: if (product->api_user) {
2016: PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatMatMult","Mat");
2017: PetscOptionsEList("-matmatmult_via","Algorithmic approach","MatMatMult",algTypes,nalg,algTypes[0],&alg,&flg);
2018: PetscOptionsEnd();
2019: } else {
2020: PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatProduct_AB","Mat");
2021: PetscOptionsEList("-matproduct_ab_via","Algorithmic approach","MatProduct_AB",algTypes,nalg,algTypes[0],&alg,&flg);
2022: PetscOptionsEnd();
2023: }
2024: if (flg) {
2025: MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);
2026: }
2028: C->ops->productsymbolic = MatProductSymbolic_AB;
2029: C->ops->matmultsymbolic = MatMatMultSymbolic_SeqAIJ_SeqAIJ;
2030: return(0);
2031: }
2033: static PetscErrorCode MatProductSetFromOptions_SeqAIJ_AtB(Mat C)
2034: {
2036: Mat_Product *product = C->product;
2037: PetscInt alg = 0; /* default algorithm */
2038: PetscBool flg = PETSC_FALSE;
2039: const char *algTypes[3] = {"default","at*b","outerproduct"};
2040: PetscInt nalg = 3;
2043: /* Get runtime option */
2044: if (product->api_user) {
2045: PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatTransposeMatMult","Mat");
2046: PetscOptionsEList("-mattransposematmult_via","Algorithmic approach","MatTransposeMatMult",algTypes,nalg,algTypes[alg],&alg,&flg);
2047: PetscOptionsEnd();
2048: } else {
2049: PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatProduct_AtB","Mat");
2050: PetscOptionsEList("-matproduct_atb_via","Algorithmic approach","MatProduct_AtB",algTypes,nalg,algTypes[alg],&alg,&flg);
2051: PetscOptionsEnd();
2052: }
2053: if (flg) {
2054: MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);
2055: }
2057: C->ops->productsymbolic = MatProductSymbolic_AtB_SeqAIJ_SeqAIJ;
2058: return(0);
2059: }
2061: static PetscErrorCode MatProductSetFromOptions_SeqAIJ_ABt(Mat C)
2062: {
2064: Mat_Product *product = C->product;
2065: PetscInt alg = 0; /* default algorithm */
2066: PetscBool flg = PETSC_FALSE;
2067: const char *algTypes[2] = {"default","color"};
2068: PetscInt nalg = 2;
2071: /* Set default algorithm */
2072: PetscStrcmp(C->product->alg,"default",&flg);
2073: if (!flg) {
2074: alg = 1;
2075: MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);
2076: }
2078: /* Get runtime option */
2079: if (product->api_user) {
2080: PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatMatTransposeMult","Mat");
2081: PetscOptionsEList("-matmattransmult_via","Algorithmic approach","MatMatTransposeMult",algTypes,nalg,algTypes[alg],&alg,&flg);
2082: PetscOptionsEnd();
2083: } else {
2084: PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatProduct_ABt","Mat");
2085: PetscOptionsEList("-matproduct_abt_via","Algorithmic approach","MatProduct_ABt",algTypes,nalg,algTypes[alg],&alg,&flg);
2086: PetscOptionsEnd();
2087: }
2088: if (flg) {
2089: MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);
2090: }
2092: C->ops->mattransposemultsymbolic = MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ;
2093: C->ops->productsymbolic = MatProductSymbolic_ABt;
2094: return(0);
2095: }
2097: static PetscErrorCode MatProductSetFromOptions_SeqAIJ_PtAP(Mat C)
2098: {
2100: Mat_Product *product = C->product;
2101: PetscBool flg = PETSC_FALSE;
2102: PetscInt alg = 0; /* default algorithm -- alg=1 should be default!!! */
2103: #if !defined(PETSC_HAVE_HYPRE)
2104: const char *algTypes[2] = {"scalable","rap"};
2105: PetscInt nalg = 2;
2106: #else
2107: const char *algTypes[3] = {"scalable","rap","hypre"};
2108: PetscInt nalg = 3;
2109: #endif
2112: /* Set default algorithm */
2113: PetscStrcmp(product->alg,"default",&flg);
2114: if (flg) {
2115: MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);
2116: }
2118: /* Get runtime option */
2119: if (product->api_user) {
2120: PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatPtAP","Mat");
2121: PetscOptionsEList("-matptap_via","Algorithmic approach","MatPtAP",algTypes,nalg,algTypes[0],&alg,&flg);
2122: PetscOptionsEnd();
2123: } else {
2124: PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatProduct_PtAP","Mat");
2125: PetscOptionsEList("-matproduct_ptap_via","Algorithmic approach","MatProduct_PtAP",algTypes,nalg,algTypes[0],&alg,&flg);
2126: PetscOptionsEnd();
2127: }
2128: if (flg) {
2129: MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);
2130: }
2132: C->ops->productsymbolic = MatProductSymbolic_PtAP_SeqAIJ_SeqAIJ;
2133: return(0);
2134: }
2136: static PetscErrorCode MatProductSetFromOptions_SeqAIJ_RARt(Mat C)
2137: {
2139: Mat_Product *product = C->product;
2140: PetscBool flg = PETSC_FALSE;
2141: PetscInt alg = 0; /* default algorithm */
2142: const char *algTypes[3] = {"r*a*rt","r*art","coloring_rart"};
2143: PetscInt nalg = 3;
2146: /* Set default algorithm */
2147: PetscStrcmp(product->alg,"default",&flg);
2148: if (flg) {
2149: MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);
2150: }
2152: /* Get runtime option */
2153: if (product->api_user) {
2154: PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatRARt","Mat");
2155: PetscOptionsEList("-matrart_via","Algorithmic approach","MatRARt",algTypes,nalg,algTypes[0],&alg,&flg);
2156: PetscOptionsEnd();
2157: } else {
2158: PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatProduct_RARt","Mat");
2159: PetscOptionsEList("-matproduct_rart_via","Algorithmic approach","MatProduct_RARt",algTypes,nalg,algTypes[0],&alg,&flg);
2160: PetscOptionsEnd();
2161: }
2162: if (flg) {
2163: MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);
2164: }
2166: C->ops->productsymbolic = MatProductSymbolic_RARt_SeqAIJ_SeqAIJ;
2167: return(0);
2168: }
2170: /* ABC = A*B*C = A*(B*C); ABC's algorithm must be chosen from AB's algorithm */
2171: static PetscErrorCode MatProductSetFromOptions_SeqAIJ_ABC(Mat C)
2172: {
2174: Mat_Product *product = C->product;
2175: PetscInt alg = 0; /* default algorithm */
2176: PetscBool flg = PETSC_FALSE;
2177: const char *algTypes[7] = {"sorted","scalable","scalable_fast","heap","btheap","llcondensed","rowmerge"};
2178: PetscInt nalg = 7;
2181: /* Set default algorithm */
2182: PetscStrcmp(product->alg,"default",&flg);
2183: if (flg) {
2184: MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);
2185: }
2187: /* Get runtime option */
2188: if (product->api_user) {
2189: PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatMatMatMult","Mat");
2190: PetscOptionsEList("-matmatmatmult_via","Algorithmic approach","MatMatMatMult",algTypes,nalg,algTypes[alg],&alg,&flg);
2191: PetscOptionsEnd();
2192: } else {
2193: PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatProduct_ABC","Mat");
2194: PetscOptionsEList("-matproduct_abc_via","Algorithmic approach","MatProduct_ABC",algTypes,nalg,algTypes[alg],&alg,&flg);
2195: PetscOptionsEnd();
2196: }
2197: if (flg) {
2198: MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);
2199: }
2201: C->ops->matmatmultsymbolic = MatMatMatMultSymbolic_SeqAIJ_SeqAIJ_SeqAIJ;
2202: C->ops->productsymbolic = MatProductSymbolic_ABC;
2203: return(0);
2204: }
2206: PetscErrorCode MatProductSetFromOptions_SeqAIJ(Mat C)
2207: {
2209: Mat_Product *product = C->product;
2212: switch (product->type) {
2213: case MATPRODUCT_AB:
2214: MatProductSetFromOptions_SeqAIJ_AB(C);
2215: break;
2216: case MATPRODUCT_AtB:
2217: MatProductSetFromOptions_SeqAIJ_AtB(C);
2218: break;
2219: case MATPRODUCT_ABt:
2220: MatProductSetFromOptions_SeqAIJ_ABt(C);
2221: break;
2222: case MATPRODUCT_PtAP:
2223: MatProductSetFromOptions_SeqAIJ_PtAP(C);
2224: break;
2225: case MATPRODUCT_RARt:
2226: MatProductSetFromOptions_SeqAIJ_RARt(C);
2227: break;
2228: case MATPRODUCT_ABC:
2229: MatProductSetFromOptions_SeqAIJ_ABC(C);
2230: break;
2231: default:
2232: break;
2233: }
2234: return(0);
2235: }