Actual source code: matmatmult.c
petsc-3.9.4 2018-09-11
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: static PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_LLCondensed(Mat,Mat,PetscReal,Mat*);
15: #if defined(PETSC_HAVE_HYPRE)
16: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_AIJ_AIJ_wHYPRE(Mat,Mat,PetscReal,Mat*);
17: #endif
19: PETSC_INTERN PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
20: {
22: #if !defined(PETSC_HAVE_HYPRE)
23: const char *algTypes[6] = {"sorted","scalable","scalable_fast","heap","btheap","llcondensed"};
24: PetscInt nalg = 6;
25: #else
26: const char *algTypes[7] = {"sorted","scalable","scalable_fast","heap","btheap","llcondensed","hypre"};
27: PetscInt nalg = 7;
28: #endif
29: PetscInt alg = 0; /* set default algorithm */
32: if (scall == MAT_INITIAL_MATRIX) {
33: PetscObjectOptionsBegin((PetscObject)A);
34: PetscOptionsObject->alreadyprinted = PETSC_FALSE; /* a hack to ensure the option shows in '-help' */
35: PetscOptionsEList("-matmatmult_via","Algorithmic approach","MatMatMult",algTypes,nalg,algTypes[0],&alg,NULL);
36: PetscOptionsEnd();
37: PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);
38: switch (alg) {
39: case 1:
40: MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(A,B,fill,C);
41: break;
42: case 2:
43: MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(A,B,fill,C);
44: break;
45: case 3:
46: MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(A,B,fill,C);
47: break;
48: case 4:
49: MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(A,B,fill,C);
50: break;
51: case 5:
52: MatMatMultSymbolic_SeqAIJ_SeqAIJ_LLCondensed(A,B,fill,C);
53: break;
54: #if defined(PETSC_HAVE_HYPRE)
55: case 6:
56: MatMatMultSymbolic_AIJ_AIJ_wHYPRE(A,B,fill,C);
57: break;
58: #endif
59: default:
60: MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);
61: break;
62: }
63: PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);
64: }
66: PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);
67: (*(*C)->ops->matmultnumeric)(A,B,*C);
68: PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);
69: return(0);
70: }
72: static PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_LLCondensed(Mat A,Mat B,PetscReal fill,Mat *C)
73: {
74: PetscErrorCode ierr;
75: Mat_SeqAIJ *a =(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
76: PetscInt *ai=a->i,*bi=b->i,*ci,*cj;
77: PetscInt am =A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
78: PetscReal afill;
79: PetscInt i,j,anzi,brow,bnzj,cnzi,*bj,*aj,*lnk,ndouble=0,Crmax;
80: PetscTable ta;
81: PetscBT lnkbt;
82: PetscFreeSpaceList free_space=NULL,current_space=NULL;
85: /* Get ci and cj */
86: /*---------------*/
87: /* Allocate ci array, arrays for fill computation and */
88: /* free space for accumulating nonzero column info */
89: PetscMalloc1(am+2,&ci);
90: ci[0] = 0;
92: /* create and initialize a linked list */
93: PetscTableCreate(bn,bn,&ta);
94: MatRowMergeMax_SeqAIJ(b,bm,ta);
95: PetscTableGetCount(ta,&Crmax);
96: PetscTableDestroy(&ta);
98: PetscLLCondensedCreate(Crmax,bn,&lnk,&lnkbt);
100: /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
101: PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am],bi[bm])),&free_space);
103: current_space = free_space;
105: /* Determine ci and cj */
106: for (i=0; i<am; i++) {
107: anzi = ai[i+1] - ai[i];
108: aj = a->j + ai[i];
109: for (j=0; j<anzi; j++) {
110: brow = aj[j];
111: bnzj = bi[brow+1] - bi[brow];
112: bj = b->j + bi[brow];
113: /* add non-zero cols of B into the sorted linked list lnk */
114: PetscLLCondensedAddSorted(bnzj,bj,lnk,lnkbt);
115: }
116: cnzi = lnk[0];
118: /* If free space is not available, make more free space */
119: /* Double the amount of total space in the list */
120: if (current_space->local_remaining<cnzi) {
121: PetscFreeSpaceGet(PetscIntSumTruncate(cnzi,current_space->total_array_size),¤t_space);
122: ndouble++;
123: }
125: /* Copy data into free space, then initialize lnk */
126: PetscLLCondensedClean(bn,cnzi,current_space->array,lnk,lnkbt);
128: current_space->array += cnzi;
129: current_space->local_used += cnzi;
130: current_space->local_remaining -= cnzi;
132: ci[i+1] = ci[i] + cnzi;
133: }
135: /* Column indices are in the list of free space */
136: /* Allocate space for cj, initialize cj, and */
137: /* destroy list of free space and other temporary array(s) */
138: PetscMalloc1(ci[am]+1,&cj);
139: PetscFreeSpaceContiguous(&free_space,cj);
140: PetscLLCondensedDestroy(lnk,lnkbt);
142: /* put together the new symbolic matrix */
143: MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,C);
144: MatSetBlockSizesFromMats(*C,A,B);
146: /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
147: /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
148: c = (Mat_SeqAIJ*)((*C)->data);
149: c->free_a = PETSC_FALSE;
150: c->free_ij = PETSC_TRUE;
151: c->nonew = 0;
152: (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ; /* fast, needs non-scalable O(bn) array 'abdense' */
154: /* set MatInfo */
155: afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
156: if (afill < 1.0) afill = 1.0;
157: c->maxnz = ci[am];
158: c->nz = ci[am];
159: (*C)->info.mallocs = ndouble;
160: (*C)->info.fill_ratio_given = fill;
161: (*C)->info.fill_ratio_needed = afill;
163: #if defined(PETSC_USE_INFO)
164: if (ci[am]) {
165: PetscInfo3((*C),"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);
166: PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);
167: } else {
168: PetscInfo((*C),"Empty matrix product\n");
169: }
170: #endif
171: return(0);
172: }
174: PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
175: {
177: PetscLogDouble flops=0.0;
178: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
179: Mat_SeqAIJ *b = (Mat_SeqAIJ*)B->data;
180: Mat_SeqAIJ *c = (Mat_SeqAIJ*)C->data;
181: PetscInt *ai =a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j;
182: PetscInt am =A->rmap->n,cm=C->rmap->n;
183: PetscInt i,j,k,anzi,bnzi,cnzi,brow;
184: PetscScalar *aa=a->a,*ba=b->a,*baj,*ca,valtmp;
185: PetscScalar *ab_dense;
188: if (!c->a) { /* first call of MatMatMultNumeric_SeqAIJ_SeqAIJ, allocate ca and matmult_abdense */
189: PetscMalloc1(ci[cm]+1,&ca);
190: c->a = ca;
191: c->free_a = PETSC_TRUE;
192: } else {
193: ca = c->a;
194: }
195: if (!c->matmult_abdense) {
196: PetscCalloc1(B->cmap->N,&ab_dense);
197: c->matmult_abdense = ab_dense;
198: } else {
199: ab_dense = c->matmult_abdense;
200: }
202: /* clean old values in C */
203: PetscMemzero(ca,ci[cm]*sizeof(MatScalar));
204: /* Traverse A row-wise. */
205: /* Build the ith row in C by summing over nonzero columns in A, */
206: /* the rows of B corresponding to nonzeros of A. */
207: for (i=0; i<am; i++) {
208: anzi = ai[i+1] - ai[i];
209: for (j=0; j<anzi; j++) {
210: brow = aj[j];
211: bnzi = bi[brow+1] - bi[brow];
212: bjj = bj + bi[brow];
213: baj = ba + bi[brow];
214: /* perform dense axpy */
215: valtmp = aa[j];
216: for (k=0; k<bnzi; k++) {
217: ab_dense[bjj[k]] += valtmp*baj[k];
218: }
219: flops += 2*bnzi;
220: }
221: aj += anzi; aa += anzi;
223: cnzi = ci[i+1] - ci[i];
224: for (k=0; k<cnzi; k++) {
225: ca[k] += ab_dense[cj[k]];
226: ab_dense[cj[k]] = 0.0; /* zero ab_dense */
227: }
228: flops += cnzi;
229: cj += cnzi; ca += cnzi;
230: }
231: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
232: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
233: PetscLogFlops(flops);
234: return(0);
235: }
237: PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(Mat A,Mat B,Mat C)
238: {
240: PetscLogDouble flops=0.0;
241: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
242: Mat_SeqAIJ *b = (Mat_SeqAIJ*)B->data;
243: Mat_SeqAIJ *c = (Mat_SeqAIJ*)C->data;
244: PetscInt *ai = a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j;
245: PetscInt am = A->rmap->N,cm=C->rmap->N;
246: PetscInt i,j,k,anzi,bnzi,cnzi,brow;
247: PetscScalar *aa=a->a,*ba=b->a,*baj,*ca=c->a,valtmp;
248: PetscInt nextb;
251: if (!ca) { /* first call of MatMatMultNumeric_SeqAIJ_SeqAIJ, allocate ca and matmult_abdense */
252: PetscMalloc1(ci[cm]+1,&ca);
253: c->a = ca;
254: c->free_a = PETSC_TRUE;
255: }
257: /* clean old values in C */
258: PetscMemzero(ca,ci[cm]*sizeof(MatScalar));
259: /* Traverse A row-wise. */
260: /* Build the ith row in C by summing over nonzero columns in A, */
261: /* the rows of B corresponding to nonzeros of A. */
262: for (i=0; i<am; i++) {
263: anzi = ai[i+1] - ai[i];
264: cnzi = ci[i+1] - ci[i];
265: for (j=0; j<anzi; j++) {
266: brow = aj[j];
267: bnzi = bi[brow+1] - bi[brow];
268: bjj = bj + bi[brow];
269: baj = ba + bi[brow];
270: /* perform sparse axpy */
271: valtmp = aa[j];
272: nextb = 0;
273: for (k=0; nextb<bnzi; k++) {
274: if (cj[k] == bjj[nextb]) { /* ccol == bcol */
275: ca[k] += valtmp*baj[nextb++];
276: }
277: }
278: flops += 2*bnzi;
279: }
280: aj += anzi; aa += anzi;
281: cj += cnzi; ca += cnzi;
282: }
284: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
285: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
286: PetscLogFlops(flops);
287: return(0);
288: }
290: PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(Mat A,Mat B,PetscReal fill,Mat *C)
291: {
292: PetscErrorCode ierr;
293: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
294: PetscInt *ai = a->i,*bi=b->i,*ci,*cj;
295: PetscInt am = A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
296: MatScalar *ca;
297: PetscReal afill;
298: PetscInt i,j,anzi,brow,bnzj,cnzi,*bj,*aj,*lnk,ndouble=0,Crmax;
299: PetscTable ta;
300: PetscFreeSpaceList free_space=NULL,current_space=NULL;
303: /* Get ci and cj - same as MatMatMultSymbolic_SeqAIJ_SeqAIJ except using PetscLLxxx_fast() */
304: /*-----------------------------------------------------------------------------------------*/
305: /* Allocate arrays for fill computation and free space for accumulating nonzero column */
306: PetscMalloc1(am+2,&ci);
307: ci[0] = 0;
309: /* create and initialize a linked list */
310: PetscTableCreate(bn,bn,&ta);
311: MatRowMergeMax_SeqAIJ(b,bm,ta);
312: PetscTableGetCount(ta,&Crmax);
313: PetscTableDestroy(&ta);
315: PetscLLCondensedCreate_fast(Crmax,&lnk);
317: /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
318: PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am],bi[bm])),&free_space);
319: current_space = free_space;
321: /* Determine ci and cj */
322: for (i=0; i<am; i++) {
323: anzi = ai[i+1] - ai[i];
324: aj = a->j + ai[i];
325: for (j=0; j<anzi; j++) {
326: brow = aj[j];
327: bnzj = bi[brow+1] - bi[brow];
328: bj = b->j + bi[brow];
329: /* add non-zero cols of B into the sorted linked list lnk */
330: PetscLLCondensedAddSorted_fast(bnzj,bj,lnk);
331: }
332: cnzi = lnk[1];
334: /* If free space is not available, make more free space */
335: /* Double the amount of total space in the list */
336: if (current_space->local_remaining<cnzi) {
337: PetscFreeSpaceGet(PetscIntSumTruncate(cnzi,current_space->total_array_size),¤t_space);
338: ndouble++;
339: }
341: /* Copy data into free space, then initialize lnk */
342: PetscLLCondensedClean_fast(cnzi,current_space->array,lnk);
344: current_space->array += cnzi;
345: current_space->local_used += cnzi;
346: current_space->local_remaining -= cnzi;
348: ci[i+1] = ci[i] + cnzi;
349: }
351: /* Column indices are in the list of free space */
352: /* Allocate space for cj, initialize cj, and */
353: /* destroy list of free space and other temporary array(s) */
354: PetscMalloc1(ci[am]+1,&cj);
355: PetscFreeSpaceContiguous(&free_space,cj);
356: PetscLLCondensedDestroy_fast(lnk);
358: /* Allocate space for ca */
359: PetscMalloc1(ci[am]+1,&ca);
360: PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));
362: /* put together the new symbolic matrix */
363: MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,ca,C);
364: MatSetBlockSizesFromMats(*C,A,B);
366: /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
367: /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
368: c = (Mat_SeqAIJ*)((*C)->data);
369: c->free_a = PETSC_TRUE;
370: c->free_ij = PETSC_TRUE;
371: c->nonew = 0;
373: (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable; /* slower, less memory */
375: /* set MatInfo */
376: afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
377: if (afill < 1.0) afill = 1.0;
378: c->maxnz = ci[am];
379: c->nz = ci[am];
380: (*C)->info.mallocs = ndouble;
381: (*C)->info.fill_ratio_given = fill;
382: (*C)->info.fill_ratio_needed = afill;
384: #if defined(PETSC_USE_INFO)
385: if (ci[am]) {
386: PetscInfo3((*C),"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);
387: PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);
388: } else {
389: PetscInfo((*C),"Empty matrix product\n");
390: }
391: #endif
392: return(0);
393: }
396: PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(Mat A,Mat B,PetscReal fill,Mat *C)
397: {
398: PetscErrorCode ierr;
399: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
400: PetscInt *ai = a->i,*bi=b->i,*ci,*cj;
401: PetscInt am = A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
402: MatScalar *ca;
403: PetscReal afill;
404: PetscInt i,j,anzi,brow,bnzj,cnzi,*bj,*aj,*lnk,ndouble=0,Crmax;
405: PetscTable ta;
406: PetscFreeSpaceList free_space=NULL,current_space=NULL;
409: /* Get ci and cj - same as MatMatMultSymbolic_SeqAIJ_SeqAIJ except using PetscLLxxx_Scalalbe() */
410: /*---------------------------------------------------------------------------------------------*/
411: /* Allocate arrays for fill computation and free space for accumulating nonzero column */
412: PetscMalloc1(am+2,&ci);
413: ci[0] = 0;
415: /* create and initialize a linked list */
416: PetscTableCreate(bn,bn,&ta);
417: MatRowMergeMax_SeqAIJ(b,bm,ta);
418: PetscTableGetCount(ta,&Crmax);
419: PetscTableDestroy(&ta);
420: PetscLLCondensedCreate_Scalable(Crmax,&lnk);
422: /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
423: PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am],bi[bm])),&free_space);
424: current_space = free_space;
426: /* Determine ci and cj */
427: for (i=0; i<am; i++) {
428: anzi = ai[i+1] - ai[i];
429: aj = a->j + ai[i];
430: for (j=0; j<anzi; j++) {
431: brow = aj[j];
432: bnzj = bi[brow+1] - bi[brow];
433: bj = b->j + bi[brow];
434: /* add non-zero cols of B into the sorted linked list lnk */
435: PetscLLCondensedAddSorted_Scalable(bnzj,bj,lnk);
436: }
437: cnzi = lnk[0];
439: /* If free space is not available, make more free space */
440: /* Double the amount of total space in the list */
441: if (current_space->local_remaining<cnzi) {
442: PetscFreeSpaceGet(PetscIntSumTruncate(cnzi,current_space->total_array_size),¤t_space);
443: ndouble++;
444: }
446: /* Copy data into free space, then initialize lnk */
447: PetscLLCondensedClean_Scalable(cnzi,current_space->array,lnk);
449: current_space->array += cnzi;
450: current_space->local_used += cnzi;
451: current_space->local_remaining -= cnzi;
453: ci[i+1] = ci[i] + cnzi;
454: }
456: /* Column indices are in the list of free space */
457: /* Allocate space for cj, initialize cj, and */
458: /* destroy list of free space and other temporary array(s) */
459: PetscMalloc1(ci[am]+1,&cj);
460: PetscFreeSpaceContiguous(&free_space,cj);
461: PetscLLCondensedDestroy_Scalable(lnk);
463: /* Allocate space for ca */
464: /*-----------------------*/
465: PetscMalloc1(ci[am]+1,&ca);
466: PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));
468: /* put together the new symbolic matrix */
469: MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,ca,C);
470: MatSetBlockSizesFromMats(*C,A,B);
472: /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
473: /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
474: c = (Mat_SeqAIJ*)((*C)->data);
475: c->free_a = PETSC_TRUE;
476: c->free_ij = PETSC_TRUE;
477: c->nonew = 0;
479: (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable; /* slower, less memory */
481: /* set MatInfo */
482: afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
483: if (afill < 1.0) afill = 1.0;
484: c->maxnz = ci[am];
485: c->nz = ci[am];
486: (*C)->info.mallocs = ndouble;
487: (*C)->info.fill_ratio_given = fill;
488: (*C)->info.fill_ratio_needed = afill;
490: #if defined(PETSC_USE_INFO)
491: if (ci[am]) {
492: PetscInfo3((*C),"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);
493: PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);
494: } else {
495: PetscInfo((*C),"Empty matrix product\n");
496: }
497: #endif
498: return(0);
499: }
501: PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(Mat A,Mat B,PetscReal fill,Mat *C)
502: {
503: PetscErrorCode ierr;
504: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
505: const PetscInt *ai=a->i,*bi=b->i,*aj=a->j,*bj=b->j;
506: PetscInt *ci,*cj,*bb;
507: PetscInt am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
508: PetscReal afill;
509: PetscInt i,j,col,ndouble = 0;
510: PetscFreeSpaceList free_space=NULL,current_space=NULL;
511: PetscHeap h;
514: /* Get ci and cj - by merging sorted rows using a heap */
515: /*---------------------------------------------------------------------------------------------*/
516: /* Allocate arrays for fill computation and free space for accumulating nonzero column */
517: PetscMalloc1(am+2,&ci);
518: ci[0] = 0;
520: /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
521: PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am],bi[bm])),&free_space);
522: current_space = free_space;
524: PetscHeapCreate(a->rmax,&h);
525: PetscMalloc1(a->rmax,&bb);
527: /* Determine ci and cj */
528: for (i=0; i<am; i++) {
529: const PetscInt anzi = ai[i+1] - ai[i]; /* number of nonzeros in this row of A, this is the number of rows of B that we merge */
530: const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */
531: ci[i+1] = ci[i];
532: /* Populate the min heap */
533: for (j=0; j<anzi; j++) {
534: bb[j] = bi[acol[j]]; /* bb points at the start of the row */
535: if (bb[j] < bi[acol[j]+1]) { /* Add if row is nonempty */
536: PetscHeapAdd(h,j,bj[bb[j]++]);
537: }
538: }
539: /* Pick off the min element, adding it to free space */
540: PetscHeapPop(h,&j,&col);
541: while (j >= 0) {
542: if (current_space->local_remaining < 1) { /* double the size, but don't exceed 16 MiB */
543: PetscFreeSpaceGet(PetscMin(PetscIntMultTruncate(2,current_space->total_array_size),16 << 20),¤t_space);
544: ndouble++;
545: }
546: *(current_space->array++) = col;
547: current_space->local_used++;
548: current_space->local_remaining--;
549: ci[i+1]++;
551: /* stash if anything else remains in this row of B */
552: if (bb[j] < bi[acol[j]+1]) {PetscHeapStash(h,j,bj[bb[j]++]);}
553: while (1) { /* pop and stash any other rows of B that also had an entry in this column */
554: PetscInt j2,col2;
555: PetscHeapPeek(h,&j2,&col2);
556: if (col2 != col) break;
557: PetscHeapPop(h,&j2,&col2);
558: if (bb[j2] < bi[acol[j2]+1]) {PetscHeapStash(h,j2,bj[bb[j2]++]);}
559: }
560: /* Put any stashed elements back into the min heap */
561: PetscHeapUnstash(h);
562: PetscHeapPop(h,&j,&col);
563: }
564: }
565: PetscFree(bb);
566: PetscHeapDestroy(&h);
568: /* Column indices are in the list of free space */
569: /* Allocate space for cj, initialize cj, and */
570: /* destroy list of free space and other temporary array(s) */
571: PetscMalloc1(ci[am],&cj);
572: PetscFreeSpaceContiguous(&free_space,cj);
574: /* put together the new symbolic matrix */
575: MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,C);
576: MatSetBlockSizesFromMats(*C,A,B);
578: /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
579: /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
580: c = (Mat_SeqAIJ*)((*C)->data);
581: c->free_a = PETSC_TRUE;
582: c->free_ij = PETSC_TRUE;
583: c->nonew = 0;
585: (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ;
587: /* set MatInfo */
588: afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
589: if (afill < 1.0) afill = 1.0;
590: c->maxnz = ci[am];
591: c->nz = ci[am];
592: (*C)->info.mallocs = ndouble;
593: (*C)->info.fill_ratio_given = fill;
594: (*C)->info.fill_ratio_needed = afill;
596: #if defined(PETSC_USE_INFO)
597: if (ci[am]) {
598: PetscInfo3((*C),"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);
599: PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);
600: } else {
601: PetscInfo((*C),"Empty matrix product\n");
602: }
603: #endif
604: return(0);
605: }
607: PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(Mat A,Mat B,PetscReal fill,Mat *C)
608: {
609: PetscErrorCode ierr;
610: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
611: const PetscInt *ai = a->i,*bi=b->i,*aj=a->j,*bj=b->j;
612: PetscInt *ci,*cj,*bb;
613: PetscInt am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
614: PetscReal afill;
615: PetscInt i,j,col,ndouble = 0;
616: PetscFreeSpaceList free_space=NULL,current_space=NULL;
617: PetscHeap h;
618: PetscBT bt;
621: /* Get ci and cj - using a heap for the sorted rows, but use BT so that each index is only added once */
622: /*---------------------------------------------------------------------------------------------*/
623: /* Allocate arrays for fill computation and free space for accumulating nonzero column */
624: PetscMalloc1(am+2,&ci);
625: ci[0] = 0;
627: /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
628: PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am],bi[bm])),&free_space);
630: current_space = free_space;
632: PetscHeapCreate(a->rmax,&h);
633: PetscMalloc1(a->rmax,&bb);
634: PetscBTCreate(bn,&bt);
636: /* Determine ci and cj */
637: for (i=0; i<am; i++) {
638: 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 */
639: const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */
640: const PetscInt *fptr = current_space->array; /* Save beginning of the row so we can clear the BT later */
641: ci[i+1] = ci[i];
642: /* Populate the min heap */
643: for (j=0; j<anzi; j++) {
644: PetscInt brow = acol[j];
645: for (bb[j] = bi[brow]; bb[j] < bi[brow+1]; bb[j]++) {
646: PetscInt bcol = bj[bb[j]];
647: if (!PetscBTLookupSet(bt,bcol)) { /* new entry */
648: PetscHeapAdd(h,j,bcol);
649: bb[j]++;
650: break;
651: }
652: }
653: }
654: /* Pick off the min element, adding it to free space */
655: PetscHeapPop(h,&j,&col);
656: while (j >= 0) {
657: if (current_space->local_remaining < 1) { /* double the size, but don't exceed 16 MiB */
658: fptr = NULL; /* need PetscBTMemzero */
659: PetscFreeSpaceGet(PetscMin(PetscIntMultTruncate(2,current_space->total_array_size),16 << 20),¤t_space);
660: ndouble++;
661: }
662: *(current_space->array++) = col;
663: current_space->local_used++;
664: current_space->local_remaining--;
665: ci[i+1]++;
667: /* stash if anything else remains in this row of B */
668: for (; bb[j] < bi[acol[j]+1]; bb[j]++) {
669: PetscInt bcol = bj[bb[j]];
670: if (!PetscBTLookupSet(bt,bcol)) { /* new entry */
671: PetscHeapAdd(h,j,bcol);
672: bb[j]++;
673: break;
674: }
675: }
676: PetscHeapPop(h,&j,&col);
677: }
678: if (fptr) { /* Clear the bits for this row */
679: for (; fptr<current_space->array; fptr++) {PetscBTClear(bt,*fptr);}
680: } else { /* We reallocated so we don't remember (easily) how to clear only the bits we changed */
681: PetscBTMemzero(bn,bt);
682: }
683: }
684: PetscFree(bb);
685: PetscHeapDestroy(&h);
686: PetscBTDestroy(&bt);
688: /* Column indices are in the list of free space */
689: /* Allocate space for cj, initialize cj, and */
690: /* destroy list of free space and other temporary array(s) */
691: PetscMalloc1(ci[am],&cj);
692: PetscFreeSpaceContiguous(&free_space,cj);
694: /* put together the new symbolic matrix */
695: MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,C);
696: MatSetBlockSizesFromMats(*C,A,B);
698: /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
699: /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
700: c = (Mat_SeqAIJ*)((*C)->data);
701: c->free_a = PETSC_TRUE;
702: c->free_ij = PETSC_TRUE;
703: c->nonew = 0;
705: (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ;
707: /* set MatInfo */
708: afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
709: if (afill < 1.0) afill = 1.0;
710: c->maxnz = ci[am];
711: c->nz = ci[am];
712: (*C)->info.mallocs = ndouble;
713: (*C)->info.fill_ratio_given = fill;
714: (*C)->info.fill_ratio_needed = afill;
716: #if defined(PETSC_USE_INFO)
717: if (ci[am]) {
718: PetscInfo3((*C),"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);
719: PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);
720: } else {
721: PetscInfo((*C),"Empty matrix product\n");
722: }
723: #endif
724: return(0);
725: }
727: /* concatenate unique entries and then sort */
728: PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
729: {
730: PetscErrorCode ierr;
731: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
732: const PetscInt *ai = a->i,*bi=b->i,*aj=a->j,*bj=b->j;
733: PetscInt *ci,*cj;
734: PetscInt am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
735: PetscReal afill;
736: PetscInt i,j,ndouble = 0;
737: PetscSegBuffer seg,segrow;
738: char *seen;
741: PetscMalloc1(am+1,&ci);
742: ci[0] = 0;
744: /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
745: PetscSegBufferCreate(sizeof(PetscInt),(PetscInt)(fill*(ai[am]+bi[bm])),&seg);
746: PetscSegBufferCreate(sizeof(PetscInt),100,&segrow);
747: PetscMalloc1(bn,&seen);
748: PetscMemzero(seen,bn*sizeof(char));
750: /* Determine ci and cj */
751: for (i=0; i<am; i++) {
752: 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 */
753: const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */
754: PetscInt packlen = 0,*PETSC_RESTRICT crow;
755: /* Pack segrow */
756: for (j=0; j<anzi; j++) {
757: PetscInt brow = acol[j],bjstart = bi[brow],bjend = bi[brow+1],k;
758: for (k=bjstart; k<bjend; k++) {
759: PetscInt bcol = bj[k];
760: if (!seen[bcol]) { /* new entry */
761: PetscInt *PETSC_RESTRICT slot;
762: PetscSegBufferGetInts(segrow,1,&slot);
763: *slot = bcol;
764: seen[bcol] = 1;
765: packlen++;
766: }
767: }
768: }
769: PetscSegBufferGetInts(seg,packlen,&crow);
770: PetscSegBufferExtractTo(segrow,crow);
771: PetscSortInt(packlen,crow);
772: ci[i+1] = ci[i] + packlen;
773: for (j=0; j<packlen; j++) seen[crow[j]] = 0;
774: }
775: PetscSegBufferDestroy(&segrow);
776: PetscFree(seen);
778: /* Column indices are in the segmented buffer */
779: PetscSegBufferExtractAlloc(seg,&cj);
780: PetscSegBufferDestroy(&seg);
782: /* put together the new symbolic matrix */
783: MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,C);
784: MatSetBlockSizesFromMats(*C,A,B);
786: /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
787: /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
788: c = (Mat_SeqAIJ*)((*C)->data);
789: c->free_a = PETSC_TRUE;
790: c->free_ij = PETSC_TRUE;
791: c->nonew = 0;
793: (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ;
795: /* set MatInfo */
796: afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
797: if (afill < 1.0) afill = 1.0;
798: c->maxnz = ci[am];
799: c->nz = ci[am];
800: (*C)->info.mallocs = ndouble;
801: (*C)->info.fill_ratio_given = fill;
802: (*C)->info.fill_ratio_needed = afill;
804: #if defined(PETSC_USE_INFO)
805: if (ci[am]) {
806: PetscInfo3((*C),"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);
807: PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);
808: } else {
809: PetscInfo((*C),"Empty matrix product\n");
810: }
811: #endif
812: return(0);
813: }
815: /* This routine is not used. Should be removed! */
816: PetscErrorCode MatMatTransposeMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
817: {
821: if (scall == MAT_INITIAL_MATRIX) {
822: PetscLogEventBegin(MAT_MatTransposeMultSymbolic,A,B,0,0);
823: MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);
824: PetscLogEventEnd(MAT_MatTransposeMultSymbolic,A,B,0,0);
825: }
826: PetscLogEventBegin(MAT_MatTransposeMultNumeric,A,B,0,0);
827: MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);
828: PetscLogEventEnd(MAT_MatTransposeMultNumeric,A,B,0,0);
829: return(0);
830: }
832: PetscErrorCode MatDestroy_SeqAIJ_MatMatMultTrans(Mat A)
833: {
834: PetscErrorCode ierr;
835: Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data;
836: Mat_MatMatTransMult *abt=a->abt;
839: (abt->destroy)(A);
840: MatTransposeColoringDestroy(&abt->matcoloring);
841: MatDestroy(&abt->Bt_den);
842: MatDestroy(&abt->ABt_den);
843: PetscFree(abt);
844: return(0);
845: }
847: PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
848: {
849: PetscErrorCode ierr;
850: Mat Bt;
851: PetscInt *bti,*btj;
852: Mat_MatMatTransMult *abt;
853: Mat_SeqAIJ *c;
856: /* create symbolic Bt */
857: MatGetSymbolicTranspose_SeqAIJ(B,&bti,&btj);
858: MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,B->cmap->n,B->rmap->n,bti,btj,NULL,&Bt);
859: MatSetBlockSizes(Bt,PetscAbs(A->cmap->bs),PetscAbs(B->cmap->bs));
861: /* get symbolic C=A*Bt */
862: MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,Bt,fill,C);
864: /* create a supporting struct for reuse intermidiate dense matrices with matcoloring */
865: PetscNew(&abt);
866: c = (Mat_SeqAIJ*)(*C)->data;
867: c->abt = abt;
869: abt->usecoloring = PETSC_FALSE;
870: abt->destroy = (*C)->ops->destroy;
871: (*C)->ops->destroy = MatDestroy_SeqAIJ_MatMatMultTrans;
873: PetscOptionsGetBool(((PetscObject)A)->options,NULL,"-matmattransmult_color",&abt->usecoloring,NULL);
874: if (abt->usecoloring) {
875: /* Create MatTransposeColoring from symbolic C=A*B^T */
876: MatTransposeColoring matcoloring;
877: MatColoring coloring;
878: ISColoring iscoloring;
879: Mat Bt_dense,C_dense;
880: Mat_SeqAIJ *c=(Mat_SeqAIJ*)(*C)->data;
881: /* inode causes memory problem, don't know why */
882: if (c->inode.use) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MAT_USE_INODES is not supported. Use '-mat_no_inode'");
884: MatColoringCreate(*C,&coloring);
885: MatColoringSetDistance(coloring,2);
886: MatColoringSetType(coloring,MATCOLORINGSL);
887: MatColoringSetFromOptions(coloring);
888: MatColoringApply(coloring,&iscoloring);
889: MatColoringDestroy(&coloring);
890: MatTransposeColoringCreate(*C,iscoloring,&matcoloring);
892: abt->matcoloring = matcoloring;
894: ISColoringDestroy(&iscoloring);
896: /* Create Bt_dense and C_dense = A*Bt_dense */
897: MatCreate(PETSC_COMM_SELF,&Bt_dense);
898: MatSetSizes(Bt_dense,A->cmap->n,matcoloring->ncolors,A->cmap->n,matcoloring->ncolors);
899: MatSetType(Bt_dense,MATSEQDENSE);
900: MatSeqDenseSetPreallocation(Bt_dense,NULL);
902: Bt_dense->assembled = PETSC_TRUE;
903: abt->Bt_den = Bt_dense;
905: MatCreate(PETSC_COMM_SELF,&C_dense);
906: MatSetSizes(C_dense,A->rmap->n,matcoloring->ncolors,A->rmap->n,matcoloring->ncolors);
907: MatSetType(C_dense,MATSEQDENSE);
908: MatSeqDenseSetPreallocation(C_dense,NULL);
910: Bt_dense->assembled = PETSC_TRUE;
911: abt->ABt_den = C_dense;
913: #if defined(PETSC_USE_INFO)
914: {
915: Mat_SeqAIJ *c = (Mat_SeqAIJ*)(*C)->data;
916: 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));
917: }
918: #endif
919: }
920: /* clean up */
921: MatDestroy(&Bt);
922: MatRestoreSymbolicTranspose_SeqAIJ(B,&bti,&btj);
923: return(0);
924: }
926: PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
927: {
928: PetscErrorCode ierr;
929: Mat_SeqAIJ *a =(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data;
930: PetscInt *ai =a->i,*aj=a->j,*bi=b->i,*bj=b->j,anzi,bnzj,nexta,nextb,*acol,*bcol,brow;
931: PetscInt cm =C->rmap->n,*ci=c->i,*cj=c->j,i,j,cnzi,*ccol;
932: PetscLogDouble flops=0.0;
933: MatScalar *aa =a->a,*aval,*ba=b->a,*bval,*ca,*cval;
934: Mat_MatMatTransMult *abt = c->abt;
937: /* clear old values in C */
938: if (!c->a) {
939: PetscMalloc1(ci[cm]+1,&ca);
940: c->a = ca;
941: c->free_a = PETSC_TRUE;
942: } else {
943: ca = c->a;
944: }
945: PetscMemzero(ca,ci[cm]*sizeof(MatScalar));
947: if (abt->usecoloring) {
948: MatTransposeColoring matcoloring = abt->matcoloring;
949: Mat Bt_dense,C_dense = abt->ABt_den;
951: /* Get Bt_dense by Apply MatTransposeColoring to B */
952: Bt_dense = abt->Bt_den;
953: MatTransColoringApplySpToDen(matcoloring,B,Bt_dense);
955: /* C_dense = A*Bt_dense */
956: MatMatMultNumeric_SeqAIJ_SeqDense(A,Bt_dense,C_dense);
958: /* Recover C from C_dense */
959: MatTransColoringApplyDenToSp(matcoloring,C_dense,C);
960: return(0);
961: }
963: for (i=0; i<cm; i++) {
964: anzi = ai[i+1] - ai[i];
965: acol = aj + ai[i];
966: aval = aa + ai[i];
967: cnzi = ci[i+1] - ci[i];
968: ccol = cj + ci[i];
969: cval = ca + ci[i];
970: for (j=0; j<cnzi; j++) {
971: brow = ccol[j];
972: bnzj = bi[brow+1] - bi[brow];
973: bcol = bj + bi[brow];
974: bval = ba + bi[brow];
976: /* perform sparse inner-product c(i,j)=A[i,:]*B[j,:]^T */
977: nexta = 0; nextb = 0;
978: while (nexta<anzi && nextb<bnzj) {
979: while (nexta < anzi && acol[nexta] < bcol[nextb]) nexta++;
980: if (nexta == anzi) break;
981: while (nextb < bnzj && acol[nexta] > bcol[nextb]) nextb++;
982: if (nextb == bnzj) break;
983: if (acol[nexta] == bcol[nextb]) {
984: cval[j] += aval[nexta]*bval[nextb];
985: nexta++; nextb++;
986: flops += 2;
987: }
988: }
989: }
990: }
991: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
992: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
993: PetscLogFlops(flops);
994: return(0);
995: }
997: PetscErrorCode MatDestroy_SeqAIJ_MatTransMatMult(Mat A)
998: {
999: PetscErrorCode ierr;
1000: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
1001: Mat_MatTransMatMult *atb = a->atb;
1004: MatDestroy(&atb->At);
1005: (atb->destroy)(A);
1006: PetscFree(atb);
1007: return(0);
1008: }
1010: PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
1011: {
1012: PetscErrorCode ierr;
1013: const char *algTypes[2] = {"matmatmult","outerproduct"};
1014: PetscInt alg=0; /* set default algorithm */
1015: Mat At;
1016: Mat_MatTransMatMult *atb;
1017: Mat_SeqAIJ *c;
1020: if (scall == MAT_INITIAL_MATRIX) {
1021: PetscObjectOptionsBegin((PetscObject)A);
1022: PetscOptionsObject->alreadyprinted = PETSC_FALSE; /* a hack to ensure the option shows in '-help' */
1023: PetscOptionsEList("-mattransposematmult_via","Algorithmic approach","MatTransposeMatMult",algTypes,2,algTypes[0],&alg,NULL);
1024: PetscOptionsEnd();
1026: switch (alg) {
1027: case 1:
1028: MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);
1029: break;
1030: default:
1031: PetscNew(&atb);
1032: MatTranspose_SeqAIJ(A,MAT_INITIAL_MATRIX,&At);
1033: MatMatMult_SeqAIJ_SeqAIJ(At,B,MAT_INITIAL_MATRIX,fill,C);
1035: c = (Mat_SeqAIJ*)(*C)->data;
1036: c->atb = atb;
1037: atb->At = At;
1038: atb->destroy = (*C)->ops->destroy;
1039: (*C)->ops->destroy = MatDestroy_SeqAIJ_MatTransMatMult;
1041: break;
1042: }
1043: }
1044: if (alg) {
1045: (*(*C)->ops->mattransposemultnumeric)(A,B,*C);
1046: } else if (!alg && scall == MAT_REUSE_MATRIX) {
1047: c = (Mat_SeqAIJ*)(*C)->data;
1048: atb = c->atb;
1049: At = atb->At;
1050: MatTranspose_SeqAIJ(A,MAT_REUSE_MATRIX,&At);
1051: MatMatMult_SeqAIJ_SeqAIJ(At,B,MAT_REUSE_MATRIX,fill,C);
1052: }
1053: return(0);
1054: }
1056: PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
1057: {
1059: Mat At;
1060: PetscInt *ati,*atj;
1063: /* create symbolic At */
1064: MatGetSymbolicTranspose_SeqAIJ(A,&ati,&atj);
1065: MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,A->cmap->n,A->rmap->n,ati,atj,NULL,&At);
1066: MatSetBlockSizes(At,PetscAbs(A->cmap->bs),PetscAbs(B->cmap->bs));
1068: /* get symbolic C=At*B */
1069: MatMatMultSymbolic_SeqAIJ_SeqAIJ(At,B,fill,C);
1071: /* clean up */
1072: MatDestroy(&At);
1073: MatRestoreSymbolicTranspose_SeqAIJ(A,&ati,&atj);
1075: (*C)->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ;
1076: return(0);
1077: }
1079: PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
1080: {
1082: Mat_SeqAIJ *a =(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data;
1083: PetscInt am =A->rmap->n,anzi,*ai=a->i,*aj=a->j,*bi=b->i,*bj,bnzi,nextb;
1084: PetscInt cm =C->rmap->n,*ci=c->i,*cj=c->j,crow,*cjj,i,j,k;
1085: PetscLogDouble flops=0.0;
1086: MatScalar *aa =a->a,*ba,*ca,*caj;
1089: if (!c->a) {
1090: PetscMalloc1(ci[cm]+1,&ca);
1092: c->a = ca;
1093: c->free_a = PETSC_TRUE;
1094: } else {
1095: ca = c->a;
1096: }
1097: /* clear old values in C */
1098: PetscMemzero(ca,ci[cm]*sizeof(MatScalar));
1100: /* compute A^T*B using outer product (A^T)[:,i]*B[i,:] */
1101: for (i=0; i<am; i++) {
1102: bj = b->j + bi[i];
1103: ba = b->a + bi[i];
1104: bnzi = bi[i+1] - bi[i];
1105: anzi = ai[i+1] - ai[i];
1106: for (j=0; j<anzi; j++) {
1107: nextb = 0;
1108: crow = *aj++;
1109: cjj = cj + ci[crow];
1110: caj = ca + ci[crow];
1111: /* perform sparse axpy operation. Note cjj includes bj. */
1112: for (k=0; nextb<bnzi; k++) {
1113: if (cjj[k] == *(bj+nextb)) { /* ccol == bcol */
1114: caj[k] += (*aa)*(*(ba+nextb));
1115: nextb++;
1116: }
1117: }
1118: flops += 2*bnzi;
1119: aa++;
1120: }
1121: }
1123: /* Assemble the final matrix and clean up */
1124: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
1125: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
1126: PetscLogFlops(flops);
1127: return(0);
1128: }
1130: PETSC_INTERN PetscErrorCode MatMatMult_SeqAIJ_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
1131: {
1135: if (scall == MAT_INITIAL_MATRIX) {
1136: PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);
1137: MatMatMultSymbolic_SeqAIJ_SeqDense(A,B,fill,C);
1138: PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);
1139: }
1140: PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);
1141: MatMatMultNumeric_SeqAIJ_SeqDense(A,B,*C);
1142: PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);
1143: return(0);
1144: }
1146: PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
1147: {
1151: MatMatMultSymbolic_SeqDense_SeqDense(A,B,0.0,C);
1153: (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqDense;
1154: return(0);
1155: }
1157: PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
1158: {
1159: Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data;
1160: Mat_SeqDense *bd = (Mat_SeqDense*)B->data;
1161: PetscErrorCode ierr;
1162: PetscScalar *c,*b,r1,r2,r3,r4,*c1,*c2,*c3,*c4,aatmp;
1163: const PetscScalar *aa,*b1,*b2,*b3,*b4;
1164: const PetscInt *aj;
1165: PetscInt cm=C->rmap->n,cn=B->cmap->n,bm=bd->lda,am=A->rmap->n;
1166: PetscInt am4=4*am,bm4=4*bm,col,i,j,n,ajtmp;
1169: if (!cm || !cn) return(0);
1170: 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);
1171: 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);
1172: 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);
1173: b = bd->v;
1174: MatDenseGetArray(C,&c);
1175: b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;
1176: c1 = c; c2 = c1 + am; c3 = c2 + am; c4 = c3 + am;
1177: for (col=0; col<cn-4; col += 4) { /* over columns of C */
1178: for (i=0; i<am; i++) { /* over rows of C in those columns */
1179: r1 = r2 = r3 = r4 = 0.0;
1180: n = a->i[i+1] - a->i[i];
1181: aj = a->j + a->i[i];
1182: aa = a->a + a->i[i];
1183: for (j=0; j<n; j++) {
1184: aatmp = aa[j]; ajtmp = aj[j];
1185: r1 += aatmp*b1[ajtmp];
1186: r2 += aatmp*b2[ajtmp];
1187: r3 += aatmp*b3[ajtmp];
1188: r4 += aatmp*b4[ajtmp];
1189: }
1190: c1[i] = r1;
1191: c2[i] = r2;
1192: c3[i] = r3;
1193: c4[i] = r4;
1194: }
1195: b1 += bm4; b2 += bm4; b3 += bm4; b4 += bm4;
1196: c1 += am4; c2 += am4; c3 += am4; c4 += am4;
1197: }
1198: for (; col<cn; col++) { /* over extra columns of C */
1199: for (i=0; i<am; i++) { /* over rows of C in those columns */
1200: r1 = 0.0;
1201: n = a->i[i+1] - a->i[i];
1202: aj = a->j + a->i[i];
1203: aa = a->a + a->i[i];
1204: for (j=0; j<n; j++) {
1205: r1 += aa[j]*b1[aj[j]];
1206: }
1207: c1[i] = r1;
1208: }
1209: b1 += bm;
1210: c1 += am;
1211: }
1212: PetscLogFlops(cn*(2.0*a->nz));
1213: MatDenseRestoreArray(C,&c);
1214: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
1215: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
1216: return(0);
1217: }
1219: /*
1220: Note very similar to MatMult_SeqAIJ(), should generate both codes from same base
1221: */
1222: PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
1223: {
1224: Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data;
1225: Mat_SeqDense *bd = (Mat_SeqDense*)B->data;
1227: PetscScalar *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4;
1228: MatScalar *aa;
1229: PetscInt cm = C->rmap->n, cn=B->cmap->n, bm=bd->lda, col, i,j,n,*aj, am = A->rmap->n,*ii,arm;
1230: PetscInt am2 = 2*am, am3 = 3*am, bm4 = 4*bm,colam,*ridx;
1233: if (!cm || !cn) return(0);
1234: b = bd->v;
1235: MatDenseGetArray(C,&c);
1236: b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;
1238: if (a->compressedrow.use) { /* use compressed row format */
1239: for (col=0; col<cn-4; col += 4) { /* over columns of C */
1240: colam = col*am;
1241: arm = a->compressedrow.nrows;
1242: ii = a->compressedrow.i;
1243: ridx = a->compressedrow.rindex;
1244: for (i=0; i<arm; i++) { /* over rows of C in those columns */
1245: r1 = r2 = r3 = r4 = 0.0;
1246: n = ii[i+1] - ii[i];
1247: aj = a->j + ii[i];
1248: aa = a->a + ii[i];
1249: for (j=0; j<n; j++) {
1250: r1 += (*aa)*b1[*aj];
1251: r2 += (*aa)*b2[*aj];
1252: r3 += (*aa)*b3[*aj];
1253: r4 += (*aa++)*b4[*aj++];
1254: }
1255: c[colam + ridx[i]] += r1;
1256: c[colam + am + ridx[i]] += r2;
1257: c[colam + am2 + ridx[i]] += r3;
1258: c[colam + am3 + ridx[i]] += r4;
1259: }
1260: b1 += bm4;
1261: b2 += bm4;
1262: b3 += bm4;
1263: b4 += bm4;
1264: }
1265: for (; col<cn; col++) { /* over extra columns of C */
1266: colam = col*am;
1267: arm = a->compressedrow.nrows;
1268: ii = a->compressedrow.i;
1269: ridx = a->compressedrow.rindex;
1270: for (i=0; i<arm; i++) { /* over rows of C in those columns */
1271: r1 = 0.0;
1272: n = ii[i+1] - ii[i];
1273: aj = a->j + ii[i];
1274: aa = a->a + ii[i];
1276: for (j=0; j<n; j++) {
1277: r1 += (*aa++)*b1[*aj++];
1278: }
1279: c[colam + ridx[i]] += r1;
1280: }
1281: b1 += bm;
1282: }
1283: } else {
1284: for (col=0; col<cn-4; col += 4) { /* over columns of C */
1285: colam = col*am;
1286: for (i=0; i<am; i++) { /* over rows of C in those columns */
1287: r1 = r2 = r3 = r4 = 0.0;
1288: n = a->i[i+1] - a->i[i];
1289: aj = a->j + a->i[i];
1290: aa = a->a + a->i[i];
1291: for (j=0; j<n; j++) {
1292: r1 += (*aa)*b1[*aj];
1293: r2 += (*aa)*b2[*aj];
1294: r3 += (*aa)*b3[*aj];
1295: r4 += (*aa++)*b4[*aj++];
1296: }
1297: c[colam + i] += r1;
1298: c[colam + am + i] += r2;
1299: c[colam + am2 + i] += r3;
1300: c[colam + am3 + i] += r4;
1301: }
1302: b1 += bm4;
1303: b2 += bm4;
1304: b3 += bm4;
1305: b4 += bm4;
1306: }
1307: for (; col<cn; col++) { /* over extra columns of C */
1308: colam = col*am;
1309: for (i=0; i<am; i++) { /* over rows of C in those columns */
1310: r1 = 0.0;
1311: n = a->i[i+1] - a->i[i];
1312: aj = a->j + a->i[i];
1313: aa = a->a + a->i[i];
1315: for (j=0; j<n; j++) {
1316: r1 += (*aa++)*b1[*aj++];
1317: }
1318: c[colam + i] += r1;
1319: }
1320: b1 += bm;
1321: }
1322: }
1323: PetscLogFlops(cn*2.0*a->nz);
1324: MatDenseRestoreArray(C,&c);
1325: return(0);
1326: }
1328: PetscErrorCode MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring coloring,Mat B,Mat Btdense)
1329: {
1331: Mat_SeqAIJ *b = (Mat_SeqAIJ*)B->data;
1332: Mat_SeqDense *btdense = (Mat_SeqDense*)Btdense->data;
1333: PetscInt *bi = b->i,*bj=b->j;
1334: PetscInt m = Btdense->rmap->n,n=Btdense->cmap->n,j,k,l,col,anz,*btcol,brow,ncolumns;
1335: MatScalar *btval,*btval_den,*ba=b->a;
1336: PetscInt *columns=coloring->columns,*colorforcol=coloring->colorforcol,ncolors=coloring->ncolors;
1339: btval_den=btdense->v;
1340: PetscMemzero(btval_den,(m*n)*sizeof(MatScalar));
1341: for (k=0; k<ncolors; k++) {
1342: ncolumns = coloring->ncolumns[k];
1343: for (l=0; l<ncolumns; l++) { /* insert a row of B to a column of Btdense */
1344: col = *(columns + colorforcol[k] + l);
1345: btcol = bj + bi[col];
1346: btval = ba + bi[col];
1347: anz = bi[col+1] - bi[col];
1348: for (j=0; j<anz; j++) {
1349: brow = btcol[j];
1350: btval_den[brow] = btval[j];
1351: }
1352: }
1353: btval_den += m;
1354: }
1355: return(0);
1356: }
1358: PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring matcoloring,Mat Cden,Mat Csp)
1359: {
1361: Mat_SeqAIJ *csp = (Mat_SeqAIJ*)Csp->data;
1362: PetscScalar *ca_den,*ca_den_ptr,*ca=csp->a;
1363: PetscInt k,l,m=Cden->rmap->n,ncolors=matcoloring->ncolors;
1364: PetscInt brows=matcoloring->brows,*den2sp=matcoloring->den2sp;
1365: PetscInt nrows,*row,*idx;
1366: PetscInt *rows=matcoloring->rows,*colorforrow=matcoloring->colorforrow;
1369: MatDenseGetArray(Cden,&ca_den);
1371: if (brows > 0) {
1372: PetscInt *lstart,row_end,row_start;
1373: lstart = matcoloring->lstart;
1374: PetscMemzero(lstart,ncolors*sizeof(PetscInt));
1376: row_end = brows;
1377: if (row_end > m) row_end = m;
1378: for (row_start=0; row_start<m; row_start+=brows) { /* loop over row blocks of Csp */
1379: ca_den_ptr = ca_den;
1380: for (k=0; k<ncolors; k++) { /* loop over colors (columns of Cden) */
1381: nrows = matcoloring->nrows[k];
1382: row = rows + colorforrow[k];
1383: idx = den2sp + colorforrow[k];
1384: for (l=lstart[k]; l<nrows; l++) {
1385: if (row[l] >= row_end) {
1386: lstart[k] = l;
1387: break;
1388: } else {
1389: ca[idx[l]] = ca_den_ptr[row[l]];
1390: }
1391: }
1392: ca_den_ptr += m;
1393: }
1394: row_end += brows;
1395: if (row_end > m) row_end = m;
1396: }
1397: } else { /* non-blocked impl: loop over columns of Csp - slow if Csp is large */
1398: ca_den_ptr = ca_den;
1399: for (k=0; k<ncolors; k++) {
1400: nrows = matcoloring->nrows[k];
1401: row = rows + colorforrow[k];
1402: idx = den2sp + colorforrow[k];
1403: for (l=0; l<nrows; l++) {
1404: ca[idx[l]] = ca_den_ptr[row[l]];
1405: }
1406: ca_den_ptr += m;
1407: }
1408: }
1410: MatDenseRestoreArray(Cden,&ca_den);
1411: #if defined(PETSC_USE_INFO)
1412: if (matcoloring->brows > 0) {
1413: PetscInfo1(Csp,"Loop over %D row blocks for den2sp\n",brows);
1414: } else {
1415: PetscInfo(Csp,"Loop over colors/columns of Cden, inefficient for large sparse matrix product \n");
1416: }
1417: #endif
1418: return(0);
1419: }
1421: PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat mat,ISColoring iscoloring,MatTransposeColoring c)
1422: {
1424: PetscInt i,n,nrows,Nbs,j,k,m,ncols,col,cm;
1425: const PetscInt *is,*ci,*cj,*row_idx;
1426: PetscInt nis = iscoloring->n,*rowhit,bs = 1;
1427: IS *isa;
1428: Mat_SeqAIJ *csp = (Mat_SeqAIJ*)mat->data;
1429: PetscInt *colorforrow,*rows,*rows_i,*idxhit,*spidx,*den2sp,*den2sp_i;
1430: PetscInt *colorforcol,*columns,*columns_i,brows;
1431: PetscBool flg;
1434: ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);
1436: /* bs >1 is not being tested yet! */
1437: Nbs = mat->cmap->N/bs;
1438: c->M = mat->rmap->N/bs; /* set total rows, columns and local rows */
1439: c->N = Nbs;
1440: c->m = c->M;
1441: c->rstart = 0;
1442: c->brows = 100;
1444: c->ncolors = nis;
1445: PetscMalloc3(nis,&c->ncolumns,nis,&c->nrows,nis+1,&colorforrow);
1446: PetscMalloc1(csp->nz+1,&rows);
1447: PetscMalloc1(csp->nz+1,&den2sp);
1449: brows = c->brows;
1450: PetscOptionsGetInt(NULL,NULL,"-matden2sp_brows",&brows,&flg);
1451: if (flg) c->brows = brows;
1452: if (brows > 0) {
1453: PetscMalloc1(nis+1,&c->lstart);
1454: }
1456: colorforrow[0] = 0;
1457: rows_i = rows;
1458: den2sp_i = den2sp;
1460: PetscMalloc1(nis+1,&colorforcol);
1461: PetscMalloc1(Nbs+1,&columns);
1463: colorforcol[0] = 0;
1464: columns_i = columns;
1466: /* get column-wise storage of mat */
1467: MatGetColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);
1469: cm = c->m;
1470: PetscMalloc1(cm+1,&rowhit);
1471: PetscMalloc1(cm+1,&idxhit);
1472: for (i=0; i<nis; i++) { /* loop over color */
1473: ISGetLocalSize(isa[i],&n);
1474: ISGetIndices(isa[i],&is);
1476: c->ncolumns[i] = n;
1477: if (n) {
1478: PetscMemcpy(columns_i,is,n*sizeof(PetscInt));
1479: }
1480: colorforcol[i+1] = colorforcol[i] + n;
1481: columns_i += n;
1483: /* fast, crude version requires O(N*N) work */
1484: PetscMemzero(rowhit,cm*sizeof(PetscInt));
1486: for (j=0; j<n; j++) { /* loop over columns*/
1487: col = is[j];
1488: row_idx = cj + ci[col];
1489: m = ci[col+1] - ci[col];
1490: for (k=0; k<m; k++) { /* loop over columns marking them in rowhit */
1491: idxhit[*row_idx] = spidx[ci[col] + k];
1492: rowhit[*row_idx++] = col + 1;
1493: }
1494: }
1495: /* count the number of hits */
1496: nrows = 0;
1497: for (j=0; j<cm; j++) {
1498: if (rowhit[j]) nrows++;
1499: }
1500: c->nrows[i] = nrows;
1501: colorforrow[i+1] = colorforrow[i] + nrows;
1503: nrows = 0;
1504: for (j=0; j<cm; j++) { /* loop over rows */
1505: if (rowhit[j]) {
1506: rows_i[nrows] = j;
1507: den2sp_i[nrows] = idxhit[j];
1508: nrows++;
1509: }
1510: }
1511: den2sp_i += nrows;
1513: ISRestoreIndices(isa[i],&is);
1514: rows_i += nrows;
1515: }
1516: MatRestoreColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);
1517: PetscFree(rowhit);
1518: ISColoringRestoreIS(iscoloring,&isa);
1519: if (csp->nz != colorforrow[nis]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"csp->nz %d != colorforrow[nis] %d",csp->nz,colorforrow[nis]);
1521: c->colorforrow = colorforrow;
1522: c->rows = rows;
1523: c->den2sp = den2sp;
1524: c->colorforcol = colorforcol;
1525: c->columns = columns;
1527: PetscFree(idxhit);
1528: return(0);
1529: }