Actual source code: matrart.c
petsc-dev 2014-02-02
2: /*
3: Defines projective product routines where A is a SeqAIJ matrix
4: C = R * A * R^T
5: */
7: #include <../src/mat/impls/aij/seq/aij.h>
8: #include <../src/mat/utils/freespace.h>
9: #include <../src/mat/impls/dense/seq/dense.h> /*I "petscmat.h" I*/
13: PetscErrorCode MatDestroy_SeqAIJ_RARt(Mat A)
14: {
16: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
17: Mat_RARt *rart = a->rart;
20: MatTransposeColoringDestroy(&rart->matcoloring);
21: MatDestroy(&rart->Rt);
22: MatDestroy(&rart->RARt);
23: MatDestroy(&rart->ARt);
24: PetscFree(rart->work);
26: A->ops->destroy = rart->destroy;
27: if (A->ops->destroy) {
28: (*A->ops->destroy)(A);
29: }
30: PetscFree(rart);
31: return(0);
32: }
36: PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_colorrart(Mat A,Mat R,PetscReal fill,Mat *C)
37: {
38: PetscErrorCode ierr;
39: Mat P;
40: PetscInt *rti,*rtj;
41: Mat_RARt *rart;
42: MatColoring coloring;
43: MatTransposeColoring matcoloring;
44: ISColoring iscoloring;
45: Mat Rt_dense,RARt_dense;
46: Mat_SeqAIJ *c;
49: /* create symbolic P=Rt */
50: MatGetSymbolicTranspose_SeqAIJ(R,&rti,&rtj);
51: MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,R->cmap->n,R->rmap->n,rti,rtj,NULL,&P);
53: /* get symbolic C=Pt*A*P */
54: MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(A,P,fill,C);
55: (*C)->rmap->bs = R->rmap->bs;
56: (*C)->cmap->bs = R->rmap->bs;
57: (*C)->ops->rartnumeric = MatRARtNumeric_SeqAIJ_SeqAIJ_colorrart;
59: /* create a supporting struct */
60: PetscNew(&rart);
61: c = (Mat_SeqAIJ*)(*C)->data;
62: c->rart = rart;
64: /* ------ Use coloring ---------- */
65: /* inode causes memory problem, don't know why */
66: if (c->inode.use) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MAT_USE_INODES is not supported. Use '-mat_no_inode'");
68: /* Create MatTransposeColoring from symbolic C=R*A*R^T */
69: MatColoringCreate(*C,&coloring);
70: MatColoringSetDistance(coloring,2);
71: MatColoringSetType(coloring,MATCOLORINGSL);
72: MatColoringSetFromOptions(coloring);
73: MatColoringApply(coloring,&iscoloring);
74: MatColoringDestroy(&coloring);
75: MatTransposeColoringCreate(*C,iscoloring,&matcoloring);
77: rart->matcoloring = matcoloring;
78: ISColoringDestroy(&iscoloring);
80: /* Create Rt_dense */
81: MatCreate(PETSC_COMM_SELF,&Rt_dense);
82: MatSetSizes(Rt_dense,A->cmap->n,matcoloring->ncolors,A->cmap->n,matcoloring->ncolors);
83: MatSetType(Rt_dense,MATSEQDENSE);
84: MatSeqDenseSetPreallocation(Rt_dense,NULL);
86: Rt_dense->assembled = PETSC_TRUE;
87: rart->Rt = Rt_dense;
89: /* Create RARt_dense = R*A*Rt_dense */
90: MatCreate(PETSC_COMM_SELF,&RARt_dense);
91: MatSetSizes(RARt_dense,(*C)->rmap->n,matcoloring->ncolors,(*C)->rmap->n,matcoloring->ncolors);
92: MatSetType(RARt_dense,MATSEQDENSE);
93: MatSeqDenseSetPreallocation(RARt_dense,NULL);
95: rart->RARt = RARt_dense;
97: /* Allocate work array to store columns of A*R^T used in MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense() */
98: PetscMalloc1(A->rmap->n*4,&rart->work);
100: rart->destroy = (*C)->ops->destroy;
101: (*C)->ops->destroy = MatDestroy_SeqAIJ_RARt;
103: /* clean up */
104: MatRestoreSymbolicTranspose_SeqAIJ(R,&rti,&rtj);
105: MatDestroy(&P);
107: #if defined(PETSC_USE_INFO)
108: {
109: PetscReal density= (PetscReal)(c->nz)/(RARt_dense->rmap->n*RARt_dense->cmap->n);
110: PetscInfo(*C,"C=R*(A*Rt) via coloring C - use sparse-dense inner products\n");
111: PetscInfo6(*C,"RARt_den %D %D; Rt %D %D (RARt->nz %D)/(m*ncolors)=%g\n",RARt_dense->rmap->n,RARt_dense->cmap->n,R->cmap->n,R->rmap->n,c->nz,density);
112: }
113: #endif
114: return(0);
115: }
117: /*
118: RAB = R * A * B, R and A in seqaij format, B in dense format;
119: */
122: PetscErrorCode MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense(Mat R,Mat A,Mat B,Mat RAB,PetscScalar *work)
123: {
124: Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*r=(Mat_SeqAIJ*)R->data;
126: PetscScalar *b,r1,r2,r3,r4,*b1,*b2,*b3,*b4;
127: MatScalar *aa,*ra;
128: PetscInt cn =B->cmap->n,bm=B->rmap->n,col,i,j,n,*ai=a->i,*aj,am=A->rmap->n;
129: PetscInt am2=2*am,am3=3*am,bm4=4*bm;
130: PetscScalar *d,*c,*c2,*c3,*c4;
131: PetscInt *rj,rm=R->rmap->n,dm=RAB->rmap->n,dn=RAB->cmap->n;
132: PetscInt rm2=2*rm,rm3=3*rm,colrm;
135: if (!dm || !dn) return(0);
136: 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);
137: if (am != R->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Number columns in R %D not equal rows in A %D\n",R->cmap->n,am);
138: if (R->rmap->n != RAB->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Number rows in RAB %D not equal rows in R %D\n",RAB->rmap->n,R->rmap->n);
139: if (B->cmap->n != RAB->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Number columns in RAB %D not equal columns in B %D\n",RAB->cmap->n,B->cmap->n);
141: { /*
142: This approach is not as good as original ones (will be removed later), but it reveals that
143: AB_den=A*B takes almost all execution time in R*A*B for src/ksp/ksp/examples/tutorials/ex56.c
144: */
145: PetscBool via_matmatmult=PETSC_FALSE;
146: PetscOptionsGetBool(NULL,"-matrart_via_matmatmult",&via_matmatmult,NULL);
147: if (via_matmatmult) {
148: Mat AB_den;
149: MatMatMultSymbolic_SeqAIJ_SeqDense(A,B,0.0,&AB_den);
150: MatMatMultNumeric_SeqAIJ_SeqDense(A,B,AB_den);
151: MatMatMultNumeric_SeqAIJ_SeqDense(R,AB_den,RAB);
152: MatDestroy(&AB_den);
153: return(0);
154: }
155: }
157: MatDenseGetArray(B,&b);
158: MatDenseGetArray(RAB,&d);
159: b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;
160: c = work; c2 = c + am; c3 = c2 + am; c4 = c3 + am;
161: for (col=0; col<cn-4; col += 4) { /* over columns of C */
162: for (i=0; i<am; i++) { /* over rows of A in those columns */
163: r1 = r2 = r3 = r4 = 0.0;
164: n = ai[i+1] - ai[i];
165: aj = a->j + ai[i];
166: aa = a->a + ai[i];
167: for (j=0; j<n; j++) {
168: r1 += (*aa)*b1[*aj];
169: r2 += (*aa)*b2[*aj];
170: r3 += (*aa)*b3[*aj];
171: r4 += (*aa++)*b4[*aj++];
172: }
173: c[i] = r1;
174: c[am + i] = r2;
175: c[am2 + i] = r3;
176: c[am3 + i] = r4;
177: }
178: b1 += bm4;
179: b2 += bm4;
180: b3 += bm4;
181: b4 += bm4;
183: /* RAB[:,col] = R*C[:,col] */
184: colrm = col*rm;
185: for (i=0; i<rm; i++) { /* over rows of R in those columns */
186: r1 = r2 = r3 = r4 = 0.0;
187: n = r->i[i+1] - r->i[i];
188: rj = r->j + r->i[i];
189: ra = r->a + r->i[i];
190: for (j=0; j<n; j++) {
191: r1 += (*ra)*c[*rj];
192: r2 += (*ra)*c2[*rj];
193: r3 += (*ra)*c3[*rj];
194: r4 += (*ra++)*c4[*rj++];
195: }
196: d[colrm + i] = r1;
197: d[colrm + rm + i] = r2;
198: d[colrm + rm2 + i] = r3;
199: d[colrm + rm3 + i] = r4;
200: }
201: }
202: for (; col<cn; col++) { /* over extra columns of C */
203: for (i=0; i<am; i++) { /* over rows of A in those columns */
204: r1 = 0.0;
205: n = a->i[i+1] - a->i[i];
206: aj = a->j + a->i[i];
207: aa = a->a + a->i[i];
208: for (j=0; j<n; j++) {
209: r1 += (*aa++)*b1[*aj++];
210: }
211: c[i] = r1;
212: }
213: b1 += bm;
215: for (i=0; i<rm; i++) { /* over rows of R in those columns */
216: r1 = 0.0;
217: n = r->i[i+1] - r->i[i];
218: rj = r->j + r->i[i];
219: ra = r->a + r->i[i];
220: for (j=0; j<n; j++) {
221: r1 += (*ra++)*c[*rj++];
222: }
223: d[col*rm + i] = r1;
224: }
225: }
226: PetscLogFlops(cn*2.0*(a->nz + r->nz));
228: MatDenseRestoreArray(B,&b);
229: MatDenseRestoreArray(RAB,&d);
230: MatAssemblyBegin(RAB,MAT_FINAL_ASSEMBLY);
231: MatAssemblyEnd(RAB,MAT_FINAL_ASSEMBLY);
232: return(0);
233: }
237: PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_colorrart(Mat A,Mat R,Mat C)
238: {
239: PetscErrorCode ierr;
240: Mat_SeqAIJ *c = (Mat_SeqAIJ*)C->data;
241: Mat_RARt *rart=c->rart;
242: MatTransposeColoring matcoloring;
243: Mat Rt,RARt;
246: /* Get dense Rt by Apply MatTransposeColoring to R */
247: matcoloring = rart->matcoloring;
248: Rt = rart->Rt;
249: MatTransColoringApplySpToDen(matcoloring,R,Rt);
251: /* Get dense RARt = R*A*Rt -- dominates! */
252: RARt = rart->RARt;
253: MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense(R,A,Rt,RARt,rart->work);
255: /* Recover C from C_dense */
256: MatTransColoringApplyDenToSp(matcoloring,RARt,C);
257: return(0);
258: }
262: PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_matmattransposemult(Mat A,Mat R,PetscReal fill,Mat *C)
263: {
264: PetscErrorCode ierr;
265: Mat ARt,RARt;
266: Mat_SeqAIJ *c;
267: Mat_RARt *rart;
270: /* must use '-mat_no_inode' with '-matmattransmult_color 1' - do not knwo why? */
271: MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(A,R,fill,&ARt);
272: MatMatMultSymbolic_SeqAIJ_SeqAIJ(R,ARt,fill,&RARt);
273: *C = RARt;
274: RARt->ops->rartnumeric = MatRARtNumeric_SeqAIJ_SeqAIJ_matmattransposemult;
276: PetscNew(&rart);
277: c = (Mat_SeqAIJ*)(*C)->data;
278: c->rart = rart;
279: rart->ARt = ARt;
280: rart->destroy = RARt->ops->destroy;
281: RARt->ops->destroy = MatDestroy_SeqAIJ_RARt;
282: #if defined(PETSC_USE_INFO)
283: PetscInfo(*C,"Use ARt=A*R^T, C=R*ARt via MatMatTransposeMult(). Coloring can be applied to A*R^T.\n");
284: #endif
285: return(0);
286: }
290: PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_matmattransposemult(Mat A,Mat R,Mat C)
291: {
292: PetscErrorCode ierr;
293: Mat_SeqAIJ *c=(Mat_SeqAIJ*)C->data;
294: Mat_RARt *rart=c->rart;
295: Mat ARt=rart->ARt;
296:
298: MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(A,R,ARt); /* dominate! */
299: MatMatMultNumeric_SeqAIJ_SeqAIJ(R,ARt,C);
300: return(0);
301: }
305: PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat R,PetscReal fill,Mat *C)
306: {
307: PetscErrorCode ierr;
308: Mat Rt;
309: Mat_SeqAIJ *c;
310: Mat_RARt *rart;
313: MatTranspose_SeqAIJ(R,MAT_INITIAL_MATRIX,&Rt);
314: MatMatMatMultSymbolic_SeqAIJ_SeqAIJ_SeqAIJ(R,A,Rt,fill,C);
316: PetscNew(&rart);
317: rart->Rt = Rt;
318: c = (Mat_SeqAIJ*)(*C)->data;
319: c->rart = rart;
320: rart->destroy = (*C)->ops->destroy;
321: (*C)->ops->destroy = MatDestroy_SeqAIJ_RARt;
322: (*C)->ops->rartnumeric = MatRARtNumeric_SeqAIJ_SeqAIJ;
323: #if defined(PETSC_USE_INFO)
324: PetscInfo(*C,"Use Rt=R^T and C=R*A*Rt via MatMatMatMult() to avoid sparse inner products\n");
325: #endif
326: return(0);
327: }
331: PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ(Mat A,Mat R,Mat C)
332: {
333: PetscErrorCode ierr;
334: Mat_SeqAIJ *c = (Mat_SeqAIJ*)C->data;
335: Mat_RARt *rart = c->rart;
336: Mat Rt = rart->Rt;
339: MatTranspose_SeqAIJ(R,MAT_REUSE_MATRIX,&Rt);
340: MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqAIJ(R,A,Rt,C);
341: return(0);
342: }
346: PetscErrorCode MatRARt_SeqAIJ_SeqAIJ(Mat A,Mat R,MatReuse scall,PetscReal fill,Mat *C)
347: {
349: const char *algTypes[3] = {"matmatmatmult","matmattransposemult","coloring_rart"};
350: PetscInt alg=0; /* set default algorithm */
351:
353: if (scall == MAT_INITIAL_MATRIX) {
354: PetscObjectOptionsBegin((PetscObject)A);
355: PetscOptionsEList("-matrart_via","Algorithmic approach","MatRARt",algTypes,3,algTypes[0],&alg,NULL);
356: PetscOptionsEnd();
358: PetscLogEventBegin(MAT_RARtSymbolic,A,R,0,0);
359: switch (alg) {
360: case 1:
361: /* via matmattransposemult: ARt=A*R^T, C=R*ARt - matrix coloring can be applied to A*R^T */
362: MatRARtSymbolic_SeqAIJ_SeqAIJ_matmattransposemult(A,R,fill,C);
363: break;
364: case 2:
365: /* via coloring_rart: apply coloring C = R*A*R^T */
366: MatRARtSymbolic_SeqAIJ_SeqAIJ_colorrart(A,R,fill,C);
367: break;
368: default:
369: /* via matmatmatmult: Rt=R^T, C=R*A*Rt - avoid inefficient sparse inner products */
370: MatRARtSymbolic_SeqAIJ_SeqAIJ(A,R,fill,C);
371: break;
372: }
373: PetscLogEventEnd(MAT_RARtSymbolic,A,R,0,0);
374: }
376: PetscLogEventBegin(MAT_RARtNumeric,A,R,0,0);
377: (*(*C)->ops->rartnumeric)(A,R,*C);
378: PetscLogEventEnd(MAT_RARtNumeric,A,R,0,0);
379: return(0);
380: }