Actual source code: matmatmult.c
petsc-3.5.4 2015-05-23
2: /*
3: Defines matrix-matrix product routines for pairs of SeqAIJ matrices
4: C = A * B
5: */
7: #include <../src/mat/impls/aij/seq/aij.h> /*I "petscmat.h" I*/
8: #include <../src/mat/utils/freespace.h>
9: #include <../src/mat/utils/petscheap.h>
10: #include <petscbt.h>
11: #include <../src/mat/impls/dense/seq/dense.h>
13: static PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_LLCondensed(Mat,Mat,PetscReal,Mat*);
17: PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
18: {
20: const char *algTypes[6] = {"sorted","scalable","scalable_fast","heap","btheap","llcondensed"};
21: PetscInt alg=0; /* set default algorithm */
24: if (scall == MAT_INITIAL_MATRIX) {
25: PetscObjectOptionsBegin((PetscObject)A);
26: PetscOptionsEList("-matmatmult_via","Algorithmic approach","MatMatMult",algTypes,6,algTypes[0],&alg,NULL);
27: PetscOptionsEnd();
28: PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);
29: switch (alg) {
30: case 1:
31: MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(A,B,fill,C);
32: break;
33: case 2:
34: MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(A,B,fill,C);
35: break;
36: case 3:
37: MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(A,B,fill,C);
38: break;
39: case 4:
40: MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(A,B,fill,C);
41: break;
42: case 5:
43: MatMatMultSymbolic_SeqAIJ_SeqAIJ_LLCondensed(A,B,fill,C);
44: break;
45: default:
46: MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);
47: break;
48: }
49: PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);
50: }
52: PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);
53: (*(*C)->ops->matmultnumeric)(A,B,*C);
54: PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);
55: return(0);
56: }
60: static PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_LLCondensed(Mat A,Mat B,PetscReal fill,Mat *C)
61: {
62: PetscErrorCode ierr;
63: Mat_SeqAIJ *a =(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
64: PetscInt *ai=a->i,*bi=b->i,*ci,*cj;
65: PetscInt am =A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
66: PetscReal afill;
67: PetscInt i,j,anzi,brow,bnzj,cnzi,*bj,*aj,nlnk_max,*lnk,ndouble=0;
68: PetscBT lnkbt;
69: PetscFreeSpaceList free_space=NULL,current_space=NULL;
72: /* Get ci and cj */
73: /*---------------*/
74: /* Allocate ci array, arrays for fill computation and */
75: /* free space for accumulating nonzero column info */
76: PetscMalloc1(((am+1)+1),&ci);
77: ci[0] = 0;
79: /* create and initialize a linked list */
80: nlnk_max = a->rmax*b->rmax;
81: if (!nlnk_max || nlnk_max > bn) nlnk_max = bn;
82: PetscLLCondensedCreate(nlnk_max,bn,&lnk,&lnkbt);
84: /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
85: PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);
87: current_space = free_space;
89: /* Determine ci and cj */
90: for (i=0; i<am; i++) {
91: anzi = ai[i+1] - ai[i];
92: aj = a->j + ai[i];
93: for (j=0; j<anzi; j++) {
94: brow = aj[j];
95: bnzj = bi[brow+1] - bi[brow];
96: bj = b->j + bi[brow];
97: /* add non-zero cols of B into the sorted linked list lnk */
98: PetscLLCondensedAddSorted(bnzj,bj,lnk,lnkbt);
99: }
100: cnzi = lnk[0];
102: /* If free space is not available, make more free space */
103: /* Double the amount of total space in the list */
104: if (current_space->local_remaining<cnzi) {
105: PetscFreeSpaceGet(cnzi+current_space->total_array_size,¤t_space);
106: ndouble++;
107: }
109: /* Copy data into free space, then initialize lnk */
110: PetscLLCondensedClean(bn,cnzi,current_space->array,lnk,lnkbt);
112: current_space->array += cnzi;
113: current_space->local_used += cnzi;
114: current_space->local_remaining -= cnzi;
116: ci[i+1] = ci[i] + cnzi;
117: }
119: /* Column indices are in the list of free space */
120: /* Allocate space for cj, initialize cj, and */
121: /* destroy list of free space and other temporary array(s) */
122: PetscMalloc1((ci[am]+1),&cj);
123: PetscFreeSpaceContiguous(&free_space,cj);
124: PetscLLCondensedDestroy(lnk,lnkbt);
126: /* put together the new symbolic matrix */
127: MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,C);
128: MatSetBlockSizesFromMats(*C,A,B);
130: /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
131: /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
132: c = (Mat_SeqAIJ*)((*C)->data);
133: c->free_a = PETSC_FALSE;
134: c->free_ij = PETSC_TRUE;
135: c->nonew = 0;
136: (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ; /* fast, needs non-scalable O(bn) array 'abdense' */
138: /* set MatInfo */
139: afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
140: if (afill < 1.0) afill = 1.0;
141: c->maxnz = ci[am];
142: c->nz = ci[am];
143: (*C)->info.mallocs = ndouble;
144: (*C)->info.fill_ratio_given = fill;
145: (*C)->info.fill_ratio_needed = afill;
147: #if defined(PETSC_USE_INFO)
148: if (ci[am]) {
149: PetscInfo3((*C),"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);
150: PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);
151: } else {
152: PetscInfo((*C),"Empty matrix product\n");
153: }
154: #endif
155: return(0);
156: }
160: PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
161: {
163: PetscLogDouble flops=0.0;
164: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
165: Mat_SeqAIJ *b = (Mat_SeqAIJ*)B->data;
166: Mat_SeqAIJ *c = (Mat_SeqAIJ*)C->data;
167: PetscInt *ai =a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j;
168: PetscInt am =A->rmap->n,cm=C->rmap->n;
169: PetscInt i,j,k,anzi,bnzi,cnzi,brow;
170: PetscScalar *aa=a->a,*ba=b->a,*baj,*ca,valtmp;
171: PetscScalar *ab_dense;
174: if (!c->a) { /* first call of MatMatMultNumeric_SeqAIJ_SeqAIJ, allocate ca and matmult_abdense */
175: PetscMalloc1((ci[cm]+1),&ca);
176: c->a = ca;
177: c->free_a = PETSC_TRUE;
178: } else {
179: ca = c->a;
180: }
181: if (!c->matmult_abdense) {
182: PetscCalloc1(B->cmap->N,&ab_dense);
183: c->matmult_abdense = ab_dense;
184: } else {
185: ab_dense = c->matmult_abdense;
186: }
188: /* clean old values in C */
189: PetscMemzero(ca,ci[cm]*sizeof(MatScalar));
190: /* Traverse A row-wise. */
191: /* Build the ith row in C by summing over nonzero columns in A, */
192: /* the rows of B corresponding to nonzeros of A. */
193: for (i=0; i<am; i++) {
194: anzi = ai[i+1] - ai[i];
195: for (j=0; j<anzi; j++) {
196: brow = aj[j];
197: bnzi = bi[brow+1] - bi[brow];
198: bjj = bj + bi[brow];
199: baj = ba + bi[brow];
200: /* perform dense axpy */
201: valtmp = aa[j];
202: for (k=0; k<bnzi; k++) {
203: ab_dense[bjj[k]] += valtmp*baj[k];
204: }
205: flops += 2*bnzi;
206: }
207: aj += anzi; aa += anzi;
209: cnzi = ci[i+1] - ci[i];
210: for (k=0; k<cnzi; k++) {
211: ca[k] += ab_dense[cj[k]];
212: ab_dense[cj[k]] = 0.0; /* zero ab_dense */
213: }
214: flops += cnzi;
215: cj += cnzi; ca += cnzi;
216: }
217: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
218: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
219: PetscLogFlops(flops);
220: return(0);
221: }
225: PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(Mat A,Mat B,Mat C)
226: {
228: PetscLogDouble flops=0.0;
229: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
230: Mat_SeqAIJ *b = (Mat_SeqAIJ*)B->data;
231: Mat_SeqAIJ *c = (Mat_SeqAIJ*)C->data;
232: PetscInt *ai = a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j;
233: PetscInt am = A->rmap->N,cm=C->rmap->N;
234: PetscInt i,j,k,anzi,bnzi,cnzi,brow;
235: PetscScalar *aa=a->a,*ba=b->a,*baj,*ca=c->a,valtmp;
236: PetscInt nextb;
239: /* clean old values in C */
240: PetscMemzero(ca,ci[cm]*sizeof(MatScalar));
241: /* Traverse A row-wise. */
242: /* Build the ith row in C by summing over nonzero columns in A, */
243: /* the rows of B corresponding to nonzeros of A. */
244: for (i=0; i<am; i++) {
245: anzi = ai[i+1] - ai[i];
246: cnzi = ci[i+1] - ci[i];
247: for (j=0; j<anzi; j++) {
248: brow = aj[j];
249: bnzi = bi[brow+1] - bi[brow];
250: bjj = bj + bi[brow];
251: baj = ba + bi[brow];
252: /* perform sparse axpy */
253: valtmp = aa[j];
254: nextb = 0;
255: for (k=0; nextb<bnzi; k++) {
256: if (cj[k] == bjj[nextb]) { /* ccol == bcol */
257: ca[k] += valtmp*baj[nextb++];
258: }
259: }
260: flops += 2*bnzi;
261: }
262: aj += anzi; aa += anzi;
263: cj += cnzi; ca += cnzi;
264: }
266: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
267: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
268: PetscLogFlops(flops);
269: return(0);
270: }
274: PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(Mat A,Mat B,PetscReal fill,Mat *C)
275: {
276: PetscErrorCode ierr;
277: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
278: PetscInt *ai = a->i,*bi=b->i,*ci,*cj;
279: PetscInt am = A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
280: MatScalar *ca;
281: PetscReal afill;
282: PetscInt i,j,anzi,brow,bnzj,cnzi,*bj,*aj,nlnk_max,*lnk,ndouble=0;
283: PetscFreeSpaceList free_space=NULL,current_space=NULL;
286: /* Get ci and cj - same as MatMatMultSymbolic_SeqAIJ_SeqAIJ except using PetscLLxxx_fast() */
287: /*-----------------------------------------------------------------------------------------*/
288: /* Allocate arrays for fill computation and free space for accumulating nonzero column */
289: PetscMalloc1(((am+1)+1),&ci);
290: ci[0] = 0;
292: /* create and initialize a linked list */
293: nlnk_max = a->rmax*b->rmax;
294: if (!nlnk_max || nlnk_max > bn) nlnk_max = bn; /* in case rmax is not defined for A or B */
295: PetscLLCondensedCreate_fast(nlnk_max,&lnk);
297: /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
298: PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);
299: current_space = free_space;
301: /* Determine ci and cj */
302: for (i=0; i<am; i++) {
303: anzi = ai[i+1] - ai[i];
304: aj = a->j + ai[i];
305: for (j=0; j<anzi; j++) {
306: brow = aj[j];
307: bnzj = bi[brow+1] - bi[brow];
308: bj = b->j + bi[brow];
309: /* add non-zero cols of B into the sorted linked list lnk */
310: PetscLLCondensedAddSorted_fast(bnzj,bj,lnk);
311: }
312: cnzi = lnk[1];
314: /* If free space is not available, make more free space */
315: /* Double the amount of total space in the list */
316: if (current_space->local_remaining<cnzi) {
317: PetscFreeSpaceGet(cnzi+current_space->total_array_size,¤t_space);
318: ndouble++;
319: }
321: /* Copy data into free space, then initialize lnk */
322: PetscLLCondensedClean_fast(cnzi,current_space->array,lnk);
324: current_space->array += cnzi;
325: current_space->local_used += cnzi;
326: current_space->local_remaining -= cnzi;
328: ci[i+1] = ci[i] + cnzi;
329: }
331: /* Column indices are in the list of free space */
332: /* Allocate space for cj, initialize cj, and */
333: /* destroy list of free space and other temporary array(s) */
334: PetscMalloc1((ci[am]+1),&cj);
335: PetscFreeSpaceContiguous(&free_space,cj);
336: PetscLLCondensedDestroy_fast(lnk);
338: /* Allocate space for ca */
339: PetscMalloc1((ci[am]+1),&ca);
340: PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));
342: /* put together the new symbolic matrix */
343: MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,ca,C);
344: MatSetBlockSizesFromMats(*C,A,B);
346: /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
347: /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
348: c = (Mat_SeqAIJ*)((*C)->data);
349: c->free_a = PETSC_TRUE;
350: c->free_ij = PETSC_TRUE;
351: c->nonew = 0;
353: (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable; /* slower, less memory */
355: /* set MatInfo */
356: afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
357: if (afill < 1.0) afill = 1.0;
358: c->maxnz = ci[am];
359: c->nz = ci[am];
360: (*C)->info.mallocs = ndouble;
361: (*C)->info.fill_ratio_given = fill;
362: (*C)->info.fill_ratio_needed = afill;
364: #if defined(PETSC_USE_INFO)
365: if (ci[am]) {
366: PetscInfo3((*C),"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);
367: PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);
368: } else {
369: PetscInfo((*C),"Empty matrix product\n");
370: }
371: #endif
372: return(0);
373: }
378: PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(Mat A,Mat B,PetscReal fill,Mat *C)
379: {
380: PetscErrorCode ierr;
381: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
382: PetscInt *ai = a->i,*bi=b->i,*ci,*cj;
383: PetscInt am = A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
384: MatScalar *ca;
385: PetscReal afill;
386: PetscInt i,j,anzi,brow,bnzj,cnzi,*bj,*aj,nlnk_max,*lnk,ndouble=0;
387: PetscFreeSpaceList free_space=NULL,current_space=NULL;
390: /* Get ci and cj - same as MatMatMultSymbolic_SeqAIJ_SeqAIJ except using PetscLLxxx_Scalalbe() */
391: /*---------------------------------------------------------------------------------------------*/
392: /* Allocate arrays for fill computation and free space for accumulating nonzero column */
393: PetscMalloc1(((am+1)+1),&ci);
394: ci[0] = 0;
396: /* create and initialize a linked list */
397: nlnk_max = a->rmax*b->rmax;
398: if (!nlnk_max || nlnk_max > bn) nlnk_max = bn; /* in case rmax is not defined for A or B */
399: PetscLLCondensedCreate_Scalable(nlnk_max,&lnk);
401: /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
402: PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);
403: current_space = free_space;
405: /* Determine ci and cj */
406: for (i=0; i<am; i++) {
407: anzi = ai[i+1] - ai[i];
408: aj = a->j + ai[i];
409: for (j=0; j<anzi; j++) {
410: brow = aj[j];
411: bnzj = bi[brow+1] - bi[brow];
412: bj = b->j + bi[brow];
413: /* add non-zero cols of B into the sorted linked list lnk */
414: PetscLLCondensedAddSorted_Scalable(bnzj,bj,lnk);
415: }
416: cnzi = lnk[0];
418: /* If free space is not available, make more free space */
419: /* Double the amount of total space in the list */
420: if (current_space->local_remaining<cnzi) {
421: PetscFreeSpaceGet(cnzi+current_space->total_array_size,¤t_space);
422: ndouble++;
423: }
425: /* Copy data into free space, then initialize lnk */
426: PetscLLCondensedClean_Scalable(cnzi,current_space->array,lnk);
428: current_space->array += cnzi;
429: current_space->local_used += cnzi;
430: current_space->local_remaining -= cnzi;
432: ci[i+1] = ci[i] + cnzi;
433: }
435: /* Column indices are in the list of free space */
436: /* Allocate space for cj, initialize cj, and */
437: /* destroy list of free space and other temporary array(s) */
438: PetscMalloc1((ci[am]+1),&cj);
439: PetscFreeSpaceContiguous(&free_space,cj);
440: PetscLLCondensedDestroy_Scalable(lnk);
442: /* Allocate space for ca */
443: /*-----------------------*/
444: PetscMalloc1((ci[am]+1),&ca);
445: PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));
447: /* put together the new symbolic matrix */
448: MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,ca,C);
449: MatSetBlockSizesFromMats(*C,A,B);
451: /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
452: /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
453: c = (Mat_SeqAIJ*)((*C)->data);
454: c->free_a = PETSC_TRUE;
455: c->free_ij = PETSC_TRUE;
456: c->nonew = 0;
458: (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable; /* slower, less memory */
460: /* set MatInfo */
461: afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
462: if (afill < 1.0) afill = 1.0;
463: c->maxnz = ci[am];
464: c->nz = ci[am];
465: (*C)->info.mallocs = ndouble;
466: (*C)->info.fill_ratio_given = fill;
467: (*C)->info.fill_ratio_needed = afill;
469: #if defined(PETSC_USE_INFO)
470: if (ci[am]) {
471: PetscInfo3((*C),"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);
472: PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);
473: } else {
474: PetscInfo((*C),"Empty matrix product\n");
475: }
476: #endif
477: return(0);
478: }
482: PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(Mat A,Mat B,PetscReal fill,Mat *C)
483: {
484: PetscErrorCode ierr;
485: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
486: const PetscInt *ai=a->i,*bi=b->i,*aj=a->j,*bj=b->j;
487: PetscInt *ci,*cj,*bb;
488: PetscInt am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
489: PetscReal afill;
490: PetscInt i,j,col,ndouble = 0;
491: PetscFreeSpaceList free_space=NULL,current_space=NULL;
492: PetscHeap h;
495: /* Get ci and cj - by merging sorted rows using a heap */
496: /*---------------------------------------------------------------------------------------------*/
497: /* Allocate arrays for fill computation and free space for accumulating nonzero column */
498: PetscMalloc1(((am+1)+1),&ci);
499: ci[0] = 0;
501: /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
502: PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);
503: current_space = free_space;
505: PetscHeapCreate(a->rmax,&h);
506: PetscMalloc1(a->rmax,&bb);
508: /* Determine ci and cj */
509: for (i=0; i<am; i++) {
510: 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 */
511: const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */
512: ci[i+1] = ci[i];
513: /* Populate the min heap */
514: for (j=0; j<anzi; j++) {
515: bb[j] = bi[acol[j]]; /* bb points at the start of the row */
516: if (bb[j] < bi[acol[j]+1]) { /* Add if row is nonempty */
517: PetscHeapAdd(h,j,bj[bb[j]++]);
518: }
519: }
520: /* Pick off the min element, adding it to free space */
521: PetscHeapPop(h,&j,&col);
522: while (j >= 0) {
523: if (current_space->local_remaining < 1) { /* double the size, but don't exceed 16 MiB */
524: PetscFreeSpaceGet(PetscMin(2*current_space->total_array_size,16 << 20),¤t_space);
525: ndouble++;
526: }
527: *(current_space->array++) = col;
528: current_space->local_used++;
529: current_space->local_remaining--;
530: ci[i+1]++;
532: /* stash if anything else remains in this row of B */
533: if (bb[j] < bi[acol[j]+1]) {PetscHeapStash(h,j,bj[bb[j]++]);}
534: while (1) { /* pop and stash any other rows of B that also had an entry in this column */
535: PetscInt j2,col2;
536: PetscHeapPeek(h,&j2,&col2);
537: if (col2 != col) break;
538: PetscHeapPop(h,&j2,&col2);
539: if (bb[j2] < bi[acol[j2]+1]) {PetscHeapStash(h,j2,bj[bb[j2]++]);}
540: }
541: /* Put any stashed elements back into the min heap */
542: PetscHeapUnstash(h);
543: PetscHeapPop(h,&j,&col);
544: }
545: }
546: PetscFree(bb);
547: PetscHeapDestroy(&h);
549: /* Column indices are in the list of free space */
550: /* Allocate space for cj, initialize cj, and */
551: /* destroy list of free space and other temporary array(s) */
552: PetscMalloc1(ci[am],&cj);
553: PetscFreeSpaceContiguous(&free_space,cj);
555: /* put together the new symbolic matrix */
556: MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,C);
557: MatSetBlockSizesFromMats(*C,A,B);
559: /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
560: /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
561: c = (Mat_SeqAIJ*)((*C)->data);
562: c->free_a = PETSC_TRUE;
563: c->free_ij = PETSC_TRUE;
564: c->nonew = 0;
566: (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ;
568: /* set MatInfo */
569: afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
570: if (afill < 1.0) afill = 1.0;
571: c->maxnz = ci[am];
572: c->nz = ci[am];
573: (*C)->info.mallocs = ndouble;
574: (*C)->info.fill_ratio_given = fill;
575: (*C)->info.fill_ratio_needed = afill;
577: #if defined(PETSC_USE_INFO)
578: if (ci[am]) {
579: PetscInfo3((*C),"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);
580: PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);
581: } else {
582: PetscInfo((*C),"Empty matrix product\n");
583: }
584: #endif
585: return(0);
586: }
590: PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(Mat A,Mat B,PetscReal fill,Mat *C)
591: {
592: PetscErrorCode ierr;
593: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
594: const PetscInt *ai = a->i,*bi=b->i,*aj=a->j,*bj=b->j;
595: PetscInt *ci,*cj,*bb;
596: PetscInt am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
597: PetscReal afill;
598: PetscInt i,j,col,ndouble = 0;
599: PetscFreeSpaceList free_space=NULL,current_space=NULL;
600: PetscHeap h;
601: PetscBT bt;
604: /* Get ci and cj - using a heap for the sorted rows, but use BT so that each index is only added once */
605: /*---------------------------------------------------------------------------------------------*/
606: /* Allocate arrays for fill computation and free space for accumulating nonzero column */
607: PetscMalloc1(((am+1)+1),&ci);
608: ci[0] = 0;
610: /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
611: PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);
613: current_space = free_space;
615: PetscHeapCreate(a->rmax,&h);
616: PetscMalloc1(a->rmax,&bb);
617: PetscBTCreate(bn,&bt);
619: /* Determine ci and cj */
620: for (i=0; i<am; i++) {
621: 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 */
622: const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */
623: const PetscInt *fptr = current_space->array; /* Save beginning of the row so we can clear the BT later */
624: ci[i+1] = ci[i];
625: /* Populate the min heap */
626: for (j=0; j<anzi; j++) {
627: PetscInt brow = acol[j];
628: for (bb[j] = bi[brow]; bb[j] < bi[brow+1]; bb[j]++) {
629: PetscInt bcol = bj[bb[j]];
630: if (!PetscBTLookupSet(bt,bcol)) { /* new entry */
631: PetscHeapAdd(h,j,bcol);
632: bb[j]++;
633: break;
634: }
635: }
636: }
637: /* Pick off the min element, adding it to free space */
638: PetscHeapPop(h,&j,&col);
639: while (j >= 0) {
640: if (current_space->local_remaining < 1) { /* double the size, but don't exceed 16 MiB */
641: fptr = NULL; /* need PetscBTMemzero */
642: PetscFreeSpaceGet(PetscMin(2*current_space->total_array_size,16 << 20),¤t_space);
643: ndouble++;
644: }
645: *(current_space->array++) = col;
646: current_space->local_used++;
647: current_space->local_remaining--;
648: ci[i+1]++;
650: /* stash if anything else remains in this row of B */
651: for (; bb[j] < bi[acol[j]+1]; bb[j]++) {
652: PetscInt bcol = bj[bb[j]];
653: if (!PetscBTLookupSet(bt,bcol)) { /* new entry */
654: PetscHeapAdd(h,j,bcol);
655: bb[j]++;
656: break;
657: }
658: }
659: PetscHeapPop(h,&j,&col);
660: }
661: if (fptr) { /* Clear the bits for this row */
662: for (; fptr<current_space->array; fptr++) {PetscBTClear(bt,*fptr);}
663: } else { /* We reallocated so we don't remember (easily) how to clear only the bits we changed */
664: PetscBTMemzero(bn,bt);
665: }
666: }
667: PetscFree(bb);
668: PetscHeapDestroy(&h);
669: PetscBTDestroy(&bt);
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: MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,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;
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: }
712: /* concatenate unique entries and then sort */
713: PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
714: {
715: PetscErrorCode ierr;
716: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c;
717: const PetscInt *ai = a->i,*bi=b->i,*aj=a->j,*bj=b->j;
718: PetscInt *ci,*cj;
719: PetscInt am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N;
720: PetscReal afill;
721: PetscInt i,j,ndouble = 0;
722: PetscSegBuffer seg,segrow;
723: char *seen;
726: PetscMalloc1((am+1),&ci);
727: ci[0] = 0;
729: /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */
730: PetscSegBufferCreate(sizeof(PetscInt),(PetscInt)(fill*(ai[am]+bi[bm])),&seg);
731: PetscSegBufferCreate(sizeof(PetscInt),100,&segrow);
732: PetscMalloc1(bn,&seen);
733: PetscMemzero(seen,bn*sizeof(char));
735: /* Determine ci and cj */
736: for (i=0; i<am; i++) {
737: 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 */
738: const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */
739: PetscInt packlen = 0,*PETSC_RESTRICT crow;
740: /* Pack segrow */
741: for (j=0; j<anzi; j++) {
742: PetscInt brow = acol[j],bjstart = bi[brow],bjend = bi[brow+1],k;
743: for (k=bjstart; k<bjend; k++) {
744: PetscInt bcol = bj[k];
745: if (!seen[bcol]) { /* new entry */
746: PetscInt *PETSC_RESTRICT slot;
747: PetscSegBufferGetInts(segrow,1,&slot);
748: *slot = bcol;
749: seen[bcol] = 1;
750: packlen++;
751: }
752: }
753: }
754: PetscSegBufferGetInts(seg,packlen,&crow);
755: PetscSegBufferExtractTo(segrow,crow);
756: PetscSortInt(packlen,crow);
757: ci[i+1] = ci[i] + packlen;
758: for (j=0; j<packlen; j++) seen[crow[j]] = 0;
759: }
760: PetscSegBufferDestroy(&segrow);
761: PetscFree(seen);
763: /* Column indices are in the segmented buffer */
764: PetscSegBufferExtractAlloc(seg,&cj);
765: PetscSegBufferDestroy(&seg);
767: /* put together the new symbolic matrix */
768: MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,C);
769: MatSetBlockSizesFromMats(*C,A,B);
771: /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
772: /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */
773: c = (Mat_SeqAIJ*)((*C)->data);
774: c->free_a = PETSC_TRUE;
775: c->free_ij = PETSC_TRUE;
776: c->nonew = 0;
778: (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ;
780: /* set MatInfo */
781: afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5;
782: if (afill < 1.0) afill = 1.0;
783: c->maxnz = ci[am];
784: c->nz = ci[am];
785: (*C)->info.mallocs = ndouble;
786: (*C)->info.fill_ratio_given = fill;
787: (*C)->info.fill_ratio_needed = afill;
789: #if defined(PETSC_USE_INFO)
790: if (ci[am]) {
791: PetscInfo3((*C),"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);
792: PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);
793: } else {
794: PetscInfo((*C),"Empty matrix product\n");
795: }
796: #endif
797: return(0);
798: }
800: /* This routine is not used. Should be removed! */
803: PetscErrorCode MatMatTransposeMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
804: {
808: if (scall == MAT_INITIAL_MATRIX) {
809: PetscLogEventBegin(MAT_MatTransposeMultSymbolic,A,B,0,0);
810: MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);
811: PetscLogEventEnd(MAT_MatTransposeMultSymbolic,A,B,0,0);
812: }
813: PetscLogEventBegin(MAT_MatTransposeMultNumeric,A,B,0,0);
814: MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);
815: PetscLogEventEnd(MAT_MatTransposeMultNumeric,A,B,0,0);
816: return(0);
817: }
821: PetscErrorCode MatDestroy_SeqAIJ_MatMatMultTrans(Mat A)
822: {
823: PetscErrorCode ierr;
824: Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data;
825: Mat_MatMatTransMult *abt=a->abt;
828: (abt->destroy)(A);
829: MatTransposeColoringDestroy(&abt->matcoloring);
830: MatDestroy(&abt->Bt_den);
831: MatDestroy(&abt->ABt_den);
832: PetscFree(abt);
833: return(0);
834: }
838: PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
839: {
840: PetscErrorCode ierr;
841: Mat Bt;
842: PetscInt *bti,*btj;
843: Mat_MatMatTransMult *abt;
844: Mat_SeqAIJ *c;
847: /* create symbolic Bt */
848: MatGetSymbolicTranspose_SeqAIJ(B,&bti,&btj);
849: MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,B->cmap->n,B->rmap->n,bti,btj,NULL,&Bt);
850: MatSetBlockSizes(Bt,PetscAbs(A->cmap->bs),PetscAbs(B->cmap->bs));
852: /* get symbolic C=A*Bt */
853: MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,Bt,fill,C);
855: /* create a supporting struct for reuse intermidiate dense matrices with matcoloring */
856: PetscNew(&abt);
857: c = (Mat_SeqAIJ*)(*C)->data;
858: c->abt = abt;
860: abt->usecoloring = PETSC_FALSE;
861: abt->destroy = (*C)->ops->destroy;
862: (*C)->ops->destroy = MatDestroy_SeqAIJ_MatMatMultTrans;
864: PetscOptionsGetBool(NULL,"-matmattransmult_color",&abt->usecoloring,NULL);
865: if (abt->usecoloring) {
866: /* Create MatTransposeColoring from symbolic C=A*B^T */
867: MatTransposeColoring matcoloring;
868: MatColoring coloring;
869: ISColoring iscoloring;
870: Mat Bt_dense,C_dense;
871: Mat_SeqAIJ *c=(Mat_SeqAIJ*)(*C)->data;
872: /* inode causes memory problem, don't know why */
873: if (c->inode.use) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MAT_USE_INODES is not supported. Use '-mat_no_inode'");
875: MatColoringCreate(*C,&coloring);
876: MatColoringSetDistance(coloring,2);
877: MatColoringSetType(coloring,MATCOLORINGSL);
878: MatColoringSetFromOptions(coloring);
879: MatColoringApply(coloring,&iscoloring);
880: MatColoringDestroy(&coloring);
881: MatTransposeColoringCreate(*C,iscoloring,&matcoloring);
883: abt->matcoloring = matcoloring;
885: ISColoringDestroy(&iscoloring);
887: /* Create Bt_dense and C_dense = A*Bt_dense */
888: MatCreate(PETSC_COMM_SELF,&Bt_dense);
889: MatSetSizes(Bt_dense,A->cmap->n,matcoloring->ncolors,A->cmap->n,matcoloring->ncolors);
890: MatSetType(Bt_dense,MATSEQDENSE);
891: MatSeqDenseSetPreallocation(Bt_dense,NULL);
893: Bt_dense->assembled = PETSC_TRUE;
894: abt->Bt_den = Bt_dense;
896: MatCreate(PETSC_COMM_SELF,&C_dense);
897: MatSetSizes(C_dense,A->rmap->n,matcoloring->ncolors,A->rmap->n,matcoloring->ncolors);
898: MatSetType(C_dense,MATSEQDENSE);
899: MatSeqDenseSetPreallocation(C_dense,NULL);
901: Bt_dense->assembled = PETSC_TRUE;
902: abt->ABt_den = C_dense;
904: #if defined(PETSC_USE_INFO)
905: {
906: Mat_SeqAIJ *c = (Mat_SeqAIJ*)(*C)->data;
907: 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));
908: }
909: #endif
910: }
911: /* clean up */
912: MatDestroy(&Bt);
913: MatRestoreSymbolicTranspose_SeqAIJ(B,&bti,&btj);
914: return(0);
915: }
919: PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
920: {
921: PetscErrorCode ierr;
922: Mat_SeqAIJ *a =(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data;
923: PetscInt *ai =a->i,*aj=a->j,*bi=b->i,*bj=b->j,anzi,bnzj,nexta,nextb,*acol,*bcol,brow;
924: PetscInt cm =C->rmap->n,*ci=c->i,*cj=c->j,i,j,cnzi,*ccol;
925: PetscLogDouble flops=0.0;
926: MatScalar *aa =a->a,*aval,*ba=b->a,*bval,*ca,*cval;
927: Mat_MatMatTransMult *abt = c->abt;
930: /* clear old values in C */
931: if (!c->a) {
932: PetscMalloc1((ci[cm]+1),&ca);
933: c->a = ca;
934: c->free_a = PETSC_TRUE;
935: } else {
936: ca = c->a;
937: }
938: PetscMemzero(ca,ci[cm]*sizeof(MatScalar));
940: if (abt->usecoloring) {
941: MatTransposeColoring matcoloring = abt->matcoloring;
942: Mat Bt_dense,C_dense = abt->ABt_den;
944: /* Get Bt_dense by Apply MatTransposeColoring to B */
945: Bt_dense = abt->Bt_den;
946: MatTransColoringApplySpToDen(matcoloring,B,Bt_dense);
948: /* C_dense = A*Bt_dense */
949: MatMatMultNumeric_SeqAIJ_SeqDense(A,Bt_dense,C_dense);
951: /* Recover C from C_dense */
952: MatTransColoringApplyDenToSp(matcoloring,C_dense,C);
953: return(0);
954: }
956: for (i=0; i<cm; i++) {
957: anzi = ai[i+1] - ai[i];
958: acol = aj + ai[i];
959: aval = aa + ai[i];
960: cnzi = ci[i+1] - ci[i];
961: ccol = cj + ci[i];
962: cval = ca + ci[i];
963: for (j=0; j<cnzi; j++) {
964: brow = ccol[j];
965: bnzj = bi[brow+1] - bi[brow];
966: bcol = bj + bi[brow];
967: bval = ba + bi[brow];
969: /* perform sparse inner-product c(i,j)=A[i,:]*B[j,:]^T */
970: nexta = 0; nextb = 0;
971: while (nexta<anzi && nextb<bnzj) {
972: while (nexta < anzi && acol[nexta] < bcol[nextb]) nexta++;
973: if (nexta == anzi) break;
974: while (nextb < bnzj && acol[nexta] > bcol[nextb]) nextb++;
975: if (nextb == bnzj) break;
976: if (acol[nexta] == bcol[nextb]) {
977: cval[j] += aval[nexta]*bval[nextb];
978: nexta++; nextb++;
979: flops += 2;
980: }
981: }
982: }
983: }
984: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
985: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
986: PetscLogFlops(flops);
987: return(0);
988: }
992: PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
993: {
997: if (scall == MAT_INITIAL_MATRIX) {
998: PetscLogEventBegin(MAT_TransposeMatMultSymbolic,A,B,0,0);
999: MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);
1000: PetscLogEventEnd(MAT_TransposeMatMultSymbolic,A,B,0,0);
1001: }
1002: PetscLogEventBegin(MAT_TransposeMatMultNumeric,A,B,0,0);
1003: MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);
1004: PetscLogEventEnd(MAT_TransposeMatMultNumeric,A,B,0,0);
1005: return(0);
1006: }
1010: PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C)
1011: {
1013: Mat At;
1014: PetscInt *ati,*atj;
1017: /* create symbolic At */
1018: MatGetSymbolicTranspose_SeqAIJ(A,&ati,&atj);
1019: MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,A->cmap->n,A->rmap->n,ati,atj,NULL,&At);
1020: MatSetBlockSizes(At,PetscAbs(A->cmap->bs),PetscAbs(B->cmap->bs));
1022: /* get symbolic C=At*B */
1023: MatMatMultSymbolic_SeqAIJ_SeqAIJ(At,B,fill,C);
1025: /* clean up */
1026: MatDestroy(&At);
1027: MatRestoreSymbolicTranspose_SeqAIJ(A,&ati,&atj);
1028: return(0);
1029: }
1033: PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C)
1034: {
1036: Mat_SeqAIJ *a =(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data;
1037: PetscInt am =A->rmap->n,anzi,*ai=a->i,*aj=a->j,*bi=b->i,*bj,bnzi,nextb;
1038: PetscInt cm =C->rmap->n,*ci=c->i,*cj=c->j,crow,*cjj,i,j,k;
1039: PetscLogDouble flops=0.0;
1040: MatScalar *aa =a->a,*ba,*ca,*caj;
1043: if (!c->a) {
1044: PetscMalloc1((ci[cm]+1),&ca);
1046: c->a = ca;
1047: c->free_a = PETSC_TRUE;
1048: } else {
1049: ca = c->a;
1050: }
1051: /* clear old values in C */
1052: PetscMemzero(ca,ci[cm]*sizeof(MatScalar));
1054: /* compute A^T*B using outer product (A^T)[:,i]*B[i,:] */
1055: for (i=0; i<am; i++) {
1056: bj = b->j + bi[i];
1057: ba = b->a + bi[i];
1058: bnzi = bi[i+1] - bi[i];
1059: anzi = ai[i+1] - ai[i];
1060: for (j=0; j<anzi; j++) {
1061: nextb = 0;
1062: crow = *aj++;
1063: cjj = cj + ci[crow];
1064: caj = ca + ci[crow];
1065: /* perform sparse axpy operation. Note cjj includes bj. */
1066: for (k=0; nextb<bnzi; k++) {
1067: if (cjj[k] == *(bj+nextb)) { /* ccol == bcol */
1068: caj[k] += (*aa)*(*(ba+nextb));
1069: nextb++;
1070: }
1071: }
1072: flops += 2*bnzi;
1073: aa++;
1074: }
1075: }
1077: /* Assemble the final matrix and clean up */
1078: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
1079: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
1080: PetscLogFlops(flops);
1081: return(0);
1082: }
1086: PetscErrorCode MatMatMult_SeqAIJ_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
1087: {
1091: if (scall == MAT_INITIAL_MATRIX) {
1092: PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);
1093: MatMatMultSymbolic_SeqAIJ_SeqDense(A,B,fill,C);
1094: PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);
1095: }
1096: PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);
1097: MatMatMultNumeric_SeqAIJ_SeqDense(A,B,*C);
1098: PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);
1099: return(0);
1100: }
1104: PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C)
1105: {
1109: MatMatMultSymbolic_SeqDense_SeqDense(A,B,0.0,C);
1111: (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqDense;
1112: return(0);
1113: }
1117: PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
1118: {
1119: Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data;
1120: PetscErrorCode ierr;
1121: PetscScalar *c,*b,r1,r2,r3,r4,*c1,*c2,*c3,*c4,aatmp;
1122: const PetscScalar *aa,*b1,*b2,*b3,*b4;
1123: const PetscInt *aj;
1124: PetscInt cm=C->rmap->n,cn=B->cmap->n,bm=B->rmap->n,am=A->rmap->n;
1125: PetscInt am4=4*am,bm4=4*bm,col,i,j,n,ajtmp;
1128: if (!cm || !cn) return(0);
1129: if (bm != A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Number columns in A %D not equal rows in B %D\n",A->cmap->n,bm);
1130: 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);
1131: 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);
1132: MatDenseGetArray(B,&b);
1133: MatDenseGetArray(C,&c);
1134: b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;
1135: c1 = c; c2 = c1 + am; c3 = c2 + am; c4 = c3 + am;
1136: for (col=0; col<cn-4; col += 4) { /* over columns of C */
1137: for (i=0; i<am; i++) { /* over rows of C in those columns */
1138: r1 = r2 = r3 = r4 = 0.0;
1139: n = a->i[i+1] - a->i[i];
1140: aj = a->j + a->i[i];
1141: aa = a->a + a->i[i];
1142: for (j=0; j<n; j++) {
1143: aatmp = aa[j]; ajtmp = aj[j];
1144: r1 += aatmp*b1[ajtmp];
1145: r2 += aatmp*b2[ajtmp];
1146: r3 += aatmp*b3[ajtmp];
1147: r4 += aatmp*b4[ajtmp];
1148: }
1149: c1[i] = r1;
1150: c2[i] = r2;
1151: c3[i] = r3;
1152: c4[i] = r4;
1153: }
1154: b1 += bm4; b2 += bm4; b3 += bm4; b4 += bm4;
1155: c1 += am4; c2 += am4; c3 += am4; c4 += am4;
1156: }
1157: for (; col<cn; col++) { /* over extra columns of C */
1158: for (i=0; i<am; i++) { /* over rows of C in those columns */
1159: r1 = 0.0;
1160: n = a->i[i+1] - a->i[i];
1161: aj = a->j + a->i[i];
1162: aa = a->a + a->i[i];
1163: for (j=0; j<n; j++) {
1164: r1 += aa[j]*b1[aj[j]];
1165: }
1166: c1[i] = r1;
1167: }
1168: b1 += bm;
1169: c1 += am;
1170: }
1171: PetscLogFlops(cn*(2.0*a->nz));
1172: MatDenseRestoreArray(B,&b);
1173: MatDenseRestoreArray(C,&c);
1174: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
1175: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
1176: return(0);
1177: }
1179: /*
1180: Note very similar to MatMult_SeqAIJ(), should generate both codes from same base
1181: */
1184: PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat A,Mat B,Mat C)
1185: {
1186: Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data;
1188: PetscScalar *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4;
1189: MatScalar *aa;
1190: PetscInt cm = C->rmap->n, cn=B->cmap->n, bm=B->rmap->n, col, i,j,n,*aj, am = A->rmap->n,*ii,arm;
1191: PetscInt am2 = 2*am, am3 = 3*am, bm4 = 4*bm,colam,*ridx;
1194: if (!cm || !cn) return(0);
1195: MatDenseGetArray(B,&b);
1196: MatDenseGetArray(C,&c);
1197: b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;
1199: if (a->compressedrow.use) { /* use compressed row format */
1200: for (col=0; col<cn-4; col += 4) { /* over columns of C */
1201: colam = col*am;
1202: arm = a->compressedrow.nrows;
1203: ii = a->compressedrow.i;
1204: ridx = a->compressedrow.rindex;
1205: for (i=0; i<arm; i++) { /* over rows of C in those columns */
1206: r1 = r2 = r3 = r4 = 0.0;
1207: n = ii[i+1] - ii[i];
1208: aj = a->j + ii[i];
1209: aa = a->a + ii[i];
1210: for (j=0; j<n; j++) {
1211: r1 += (*aa)*b1[*aj];
1212: r2 += (*aa)*b2[*aj];
1213: r3 += (*aa)*b3[*aj];
1214: r4 += (*aa++)*b4[*aj++];
1215: }
1216: c[colam + ridx[i]] += r1;
1217: c[colam + am + ridx[i]] += r2;
1218: c[colam + am2 + ridx[i]] += r3;
1219: c[colam + am3 + ridx[i]] += r4;
1220: }
1221: b1 += bm4;
1222: b2 += bm4;
1223: b3 += bm4;
1224: b4 += bm4;
1225: }
1226: for (; col<cn; col++) { /* over extra columns of C */
1227: colam = col*am;
1228: arm = a->compressedrow.nrows;
1229: ii = a->compressedrow.i;
1230: ridx = a->compressedrow.rindex;
1231: for (i=0; i<arm; i++) { /* over rows of C in those columns */
1232: r1 = 0.0;
1233: n = ii[i+1] - ii[i];
1234: aj = a->j + ii[i];
1235: aa = a->a + ii[i];
1237: for (j=0; j<n; j++) {
1238: r1 += (*aa++)*b1[*aj++];
1239: }
1240: c[colam + ridx[i]] += r1;
1241: }
1242: b1 += bm;
1243: }
1244: } else {
1245: for (col=0; col<cn-4; col += 4) { /* over columns of C */
1246: colam = col*am;
1247: for (i=0; i<am; i++) { /* over rows of C in those columns */
1248: r1 = r2 = r3 = r4 = 0.0;
1249: n = a->i[i+1] - a->i[i];
1250: aj = a->j + a->i[i];
1251: aa = a->a + a->i[i];
1252: for (j=0; j<n; j++) {
1253: r1 += (*aa)*b1[*aj];
1254: r2 += (*aa)*b2[*aj];
1255: r3 += (*aa)*b3[*aj];
1256: r4 += (*aa++)*b4[*aj++];
1257: }
1258: c[colam + i] += r1;
1259: c[colam + am + i] += r2;
1260: c[colam + am2 + i] += r3;
1261: c[colam + am3 + i] += r4;
1262: }
1263: b1 += bm4;
1264: b2 += bm4;
1265: b3 += bm4;
1266: b4 += bm4;
1267: }
1268: for (; col<cn; col++) { /* over extra columns of C */
1269: colam = col*am;
1270: for (i=0; i<am; i++) { /* over rows of C in those columns */
1271: r1 = 0.0;
1272: n = a->i[i+1] - a->i[i];
1273: aj = a->j + a->i[i];
1274: aa = a->a + a->i[i];
1276: for (j=0; j<n; j++) {
1277: r1 += (*aa++)*b1[*aj++];
1278: }
1279: c[colam + i] += r1;
1280: }
1281: b1 += bm;
1282: }
1283: }
1284: PetscLogFlops(cn*2.0*a->nz);
1285: MatDenseRestoreArray(B,&b);
1286: MatDenseRestoreArray(C,&c);
1287: return(0);
1288: }
1292: PetscErrorCode MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring coloring,Mat B,Mat Btdense)
1293: {
1295: Mat_SeqAIJ *b = (Mat_SeqAIJ*)B->data;
1296: Mat_SeqDense *btdense = (Mat_SeqDense*)Btdense->data;
1297: PetscInt *bi = b->i,*bj=b->j;
1298: PetscInt m = Btdense->rmap->n,n=Btdense->cmap->n,j,k,l,col,anz,*btcol,brow,ncolumns;
1299: MatScalar *btval,*btval_den,*ba=b->a;
1300: PetscInt *columns=coloring->columns,*colorforcol=coloring->colorforcol,ncolors=coloring->ncolors;
1303: btval_den=btdense->v;
1304: PetscMemzero(btval_den,(m*n)*sizeof(MatScalar));
1305: for (k=0; k<ncolors; k++) {
1306: ncolumns = coloring->ncolumns[k];
1307: for (l=0; l<ncolumns; l++) { /* insert a row of B to a column of Btdense */
1308: col = *(columns + colorforcol[k] + l);
1309: btcol = bj + bi[col];
1310: btval = ba + bi[col];
1311: anz = bi[col+1] - bi[col];
1312: for (j=0; j<anz; j++) {
1313: brow = btcol[j];
1314: btval_den[brow] = btval[j];
1315: }
1316: }
1317: btval_den += m;
1318: }
1319: return(0);
1320: }
1324: PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring matcoloring,Mat Cden,Mat Csp)
1325: {
1327: Mat_SeqAIJ *csp = (Mat_SeqAIJ*)Csp->data;
1328: PetscScalar *ca_den,*ca_den_ptr,*ca=csp->a;
1329: PetscInt k,l,m=Cden->rmap->n,ncolors=matcoloring->ncolors;
1330: PetscInt brows=matcoloring->brows,*den2sp=matcoloring->den2sp;
1331: PetscInt nrows,*row,*idx;
1332: PetscInt *rows=matcoloring->rows,*colorforrow=matcoloring->colorforrow;
1335: MatDenseGetArray(Cden,&ca_den);
1337: if (brows > 0) {
1338: PetscInt *lstart,row_end,row_start;
1339: lstart = matcoloring->lstart;
1340: PetscMemzero(lstart,ncolors*sizeof(PetscInt));
1342: row_end = brows;
1343: if (row_end > m) row_end = m;
1344: for (row_start=0; row_start<m; row_start+=brows) { /* loop over row blocks of Csp */
1345: ca_den_ptr = ca_den;
1346: for (k=0; k<ncolors; k++) { /* loop over colors (columns of Cden) */
1347: nrows = matcoloring->nrows[k];
1348: row = rows + colorforrow[k];
1349: idx = den2sp + colorforrow[k];
1350: for (l=lstart[k]; l<nrows; l++) {
1351: if (row[l] >= row_end) {
1352: lstart[k] = l;
1353: break;
1354: } else {
1355: ca[idx[l]] = ca_den_ptr[row[l]];
1356: }
1357: }
1358: ca_den_ptr += m;
1359: }
1360: row_end += brows;
1361: if (row_end > m) row_end = m;
1362: }
1363: } else { /* non-blocked impl: loop over columns of Csp - slow if Csp is large */
1364: ca_den_ptr = ca_den;
1365: for (k=0; k<ncolors; k++) {
1366: nrows = matcoloring->nrows[k];
1367: row = rows + colorforrow[k];
1368: idx = den2sp + colorforrow[k];
1369: for (l=0; l<nrows; l++) {
1370: ca[idx[l]] = ca_den_ptr[row[l]];
1371: }
1372: ca_den_ptr += m;
1373: }
1374: }
1376: MatDenseRestoreArray(Cden,&ca_den);
1377: #if defined(PETSC_USE_INFO)
1378: if (matcoloring->brows > 0) {
1379: PetscInfo1(Csp,"Loop over %D row blocks for den2sp\n",brows);
1380: } else {
1381: PetscInfo(Csp,"Loop over colors/columns of Cden, inefficient for large sparse matrix product \n");
1382: }
1383: #endif
1384: return(0);
1385: }
1389: PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat mat,ISColoring iscoloring,MatTransposeColoring c)
1390: {
1392: PetscInt i,n,nrows,Nbs,j,k,m,ncols,col,cm;
1393: const PetscInt *is,*ci,*cj,*row_idx;
1394: PetscInt nis = iscoloring->n,*rowhit,bs = 1;
1395: IS *isa;
1396: Mat_SeqAIJ *csp = (Mat_SeqAIJ*)mat->data;
1397: PetscInt *colorforrow,*rows,*rows_i,*idxhit,*spidx,*den2sp,*den2sp_i;
1398: PetscInt *colorforcol,*columns,*columns_i,brows;
1399: PetscBool flg;
1402: ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);
1404: /* bs >1 is not being tested yet! */
1405: Nbs = mat->cmap->N/bs;
1406: c->M = mat->rmap->N/bs; /* set total rows, columns and local rows */
1407: c->N = Nbs;
1408: c->m = c->M;
1409: c->rstart = 0;
1410: c->brows = 100;
1412: c->ncolors = nis;
1413: PetscMalloc3(nis,&c->ncolumns,nis,&c->nrows,nis+1,&colorforrow);
1414: PetscMalloc1((csp->nz+1),&rows);
1415: PetscMalloc1((csp->nz+1),&den2sp);
1417: brows = c->brows;
1418: PetscOptionsGetInt(NULL,"-matden2sp_brows",&brows,&flg);
1419: if (flg) c->brows = brows;
1420: if (brows > 0) {
1421: PetscMalloc1((nis+1),&c->lstart);
1422: }
1424: colorforrow[0] = 0;
1425: rows_i = rows;
1426: den2sp_i = den2sp;
1428: PetscMalloc1((nis+1),&colorforcol);
1429: PetscMalloc1((Nbs+1),&columns);
1431: colorforcol[0] = 0;
1432: columns_i = columns;
1434: /* get column-wise storage of mat */
1435: MatGetColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);
1437: cm = c->m;
1438: PetscMalloc1((cm+1),&rowhit);
1439: PetscMalloc1((cm+1),&idxhit);
1440: for (i=0; i<nis; i++) { /* loop over color */
1441: ISGetLocalSize(isa[i],&n);
1442: ISGetIndices(isa[i],&is);
1444: c->ncolumns[i] = n;
1445: if (n) {
1446: PetscMemcpy(columns_i,is,n*sizeof(PetscInt));
1447: }
1448: colorforcol[i+1] = colorforcol[i] + n;
1449: columns_i += n;
1451: /* fast, crude version requires O(N*N) work */
1452: PetscMemzero(rowhit,cm*sizeof(PetscInt));
1454: for (j=0; j<n; j++) { /* loop over columns*/
1455: col = is[j];
1456: row_idx = cj + ci[col];
1457: m = ci[col+1] - ci[col];
1458: for (k=0; k<m; k++) { /* loop over columns marking them in rowhit */
1459: idxhit[*row_idx] = spidx[ci[col] + k];
1460: rowhit[*row_idx++] = col + 1;
1461: }
1462: }
1463: /* count the number of hits */
1464: nrows = 0;
1465: for (j=0; j<cm; j++) {
1466: if (rowhit[j]) nrows++;
1467: }
1468: c->nrows[i] = nrows;
1469: colorforrow[i+1] = colorforrow[i] + nrows;
1471: nrows = 0;
1472: for (j=0; j<cm; j++) { /* loop over rows */
1473: if (rowhit[j]) {
1474: rows_i[nrows] = j;
1475: den2sp_i[nrows] = idxhit[j];
1476: nrows++;
1477: }
1478: }
1479: den2sp_i += nrows;
1481: ISRestoreIndices(isa[i],&is);
1482: rows_i += nrows;
1483: }
1484: MatRestoreColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);
1485: PetscFree(rowhit);
1486: ISColoringRestoreIS(iscoloring,&isa);
1487: if (csp->nz != colorforrow[nis]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"csp->nz %d != colorforrow[nis] %d",csp->nz,colorforrow[nis]);
1489: c->colorforrow = colorforrow;
1490: c->rows = rows;
1491: c->den2sp = den2sp;
1492: c->colorforcol = colorforcol;
1493: c->columns = columns;
1495: PetscFree(idxhit);
1496: return(0);
1497: }