xref: /petsc/src/mat/impls/aij/seq/matrart.c (revision b56132cb8045a4c06d79b76dd784c967d609210c)
1 
2 /*
3   Defines projective product routines where A is a SeqAIJ matrix
4           C = R * A * R^T
5 */
6 
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*/
10 #include <petsctime.h>
11 
12 /* #define RART_PROFILE */
13 
14 #undef __FUNCT__
15 #define __FUNCT__ "MatDestroy_SeqAIJ_RARt"
16 PetscErrorCode MatDestroy_SeqAIJ_RARt(Mat A)
17 {
18   PetscErrorCode ierr;
19   Mat_SeqAIJ     *a    = (Mat_SeqAIJ*)A->data;
20   Mat_RARt       *rart = a->rart;
21 
22   PetscFunctionBegin;
23   ierr = MatTransposeColoringDestroy(&rart->matcoloring);CHKERRQ(ierr);
24   ierr = MatDestroy(&rart->Rt);CHKERRQ(ierr);
25   ierr = MatDestroy(&rart->RARt);CHKERRQ(ierr);
26   ierr = MatDestroy(&rart->ARt);CHKERRQ(ierr);
27   ierr = PetscFree(rart->work);CHKERRQ(ierr);
28 
29   A->ops->destroy = rart->destroy;
30   if (A->ops->destroy) {
31     ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
32   }
33   ierr = PetscFree(rart);CHKERRQ(ierr);
34   PetscFunctionReturn(0);
35 }
36 
37 #undef __FUNCT__
38 #define __FUNCT__ "MatRARtSymbolic_SeqAIJ_SeqAIJ_colorrart"
39 PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_colorrart(Mat A,Mat R,PetscReal fill,Mat *C)
40 {
41   PetscErrorCode       ierr;
42   Mat                  P;
43   PetscInt             *rti,*rtj;
44   Mat_RARt             *rart;
45   MatTransposeColoring matcoloring;
46   ISColoring           iscoloring;
47   Mat                  Rt_dense,RARt_dense;
48 #if defined(RART_PROFILE)
49   PetscLogDouble       GColor=0.0,MCCreate=0.0,MDenCreate=0.0,t0,tf,etime=0.0;
50 #endif
51   Mat_SeqAIJ           *c;
52 
53   PetscFunctionBegin;
54 #if defined(RART_PROFILE)
55   ierr = PetscTime(&t0);CHKERRQ(ierr);
56 #endif
57   /* create symbolic P=Rt */
58   ierr = MatGetSymbolicTranspose_SeqAIJ(R,&rti,&rtj);CHKERRQ(ierr);
59   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,R->cmap->n,R->rmap->n,rti,rtj,NULL,&P);CHKERRQ(ierr);
60 
61   /* get symbolic C=Pt*A*P */
62   ierr           = MatPtAPSymbolic_SeqAIJ_SeqAIJ(A,P,fill,C);CHKERRQ(ierr);
63   (*C)->rmap->bs = R->rmap->bs;
64   (*C)->cmap->bs = R->rmap->bs;
65   (*C)->ops->rartnumeric = MatRARtNumeric_SeqAIJ_SeqAIJ_colorrart;
66 
67   /* create a supporting struct */
68   ierr = PetscNew(Mat_RARt,&rart);CHKERRQ(ierr);
69   c       = (Mat_SeqAIJ*)(*C)->data;
70   c->rart = rart;
71 #if defined(RART_PROFILE)
72   ierr   = PetscTime(&tf);CHKERRQ(ierr);
73   etime += tf - t0;
74 #endif
75 
76   /* ------ Use coloring ---------- */
77   /* inode causes memory problem, don't know why */
78   if (c->inode.use) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MAT_USE_INODES is not supported. Use '-mat_no_inode'");
79 
80   /* Create MatTransposeColoring from symbolic C=R*A*R^T */
81 #if defined(RART_PROFILE)
82   ierr    = PetscTime(&t0);CHKERRQ(ierr);
83 #endif
84   ierr    = MatGetColoring(*C,MATCOLORINGLF,&iscoloring);CHKERRQ(ierr);
85 #if defined(RART_PROFILE)
86   ierr    = PetscTime(&tf);CHKERRQ(ierr);
87   GColor += tf - t0;
88   ierr = PetscTime(&t0);CHKERRQ(ierr);
89 #endif
90   ierr = MatTransposeColoringCreate(*C,iscoloring,&matcoloring);CHKERRQ(ierr);
91 
92   rart->matcoloring = matcoloring;
93 
94   ierr      = ISColoringDestroy(&iscoloring);CHKERRQ(ierr);
95 #if defined(RART_PROFILE)
96   ierr      = PetscTime(&tf);CHKERRQ(ierr);
97   MCCreate += tf - t0;
98   ierr = PetscTime(&t0);CHKERRQ(ierr);
99 #endif
100 
101   /* Create Rt_dense */
102   ierr = MatCreate(PETSC_COMM_SELF,&Rt_dense);CHKERRQ(ierr);
103   ierr = MatSetSizes(Rt_dense,A->cmap->n,matcoloring->ncolors,A->cmap->n,matcoloring->ncolors);CHKERRQ(ierr);
104   ierr = MatSetType(Rt_dense,MATSEQDENSE);CHKERRQ(ierr);
105   ierr = MatSeqDenseSetPreallocation(Rt_dense,NULL);CHKERRQ(ierr);
106 
107   Rt_dense->assembled = PETSC_TRUE;
108   rart->Rt            = Rt_dense;
109 
110   /* Create RARt_dense = R*A*Rt_dense */
111   ierr = MatCreate(PETSC_COMM_SELF,&RARt_dense);CHKERRQ(ierr);
112   ierr = MatSetSizes(RARt_dense,(*C)->rmap->n,matcoloring->ncolors,(*C)->rmap->n,matcoloring->ncolors);CHKERRQ(ierr);
113   ierr = MatSetType(RARt_dense,MATSEQDENSE);CHKERRQ(ierr);
114   ierr = MatSeqDenseSetPreallocation(RARt_dense,NULL);CHKERRQ(ierr);
115 
116   rart->RARt = RARt_dense;
117 
118   /* Allocate work array to store columns of A*R^T used in MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense() */
119   ierr = PetscMalloc(A->rmap->n*4*sizeof(PetscScalar),&rart->work);CHKERRQ(ierr);
120 
121 #if defined(RART_PROFILE)
122   ierr        = PetscTime(&tf);CHKERRQ(ierr);
123   MDenCreate += tf - t0;
124 #endif
125 
126   rart->destroy      = (*C)->ops->destroy;
127   (*C)->ops->destroy = MatDestroy_SeqAIJ_RARt;
128 
129   /* clean up */
130   ierr = MatRestoreSymbolicTranspose_SeqAIJ(R,&rti,&rtj);CHKERRQ(ierr);
131   ierr = MatDestroy(&P);CHKERRQ(ierr);
132 
133 #if defined(PETSC_USE_INFO)
134   {
135     PetscReal density= (PetscReal)(c->nz)/(RARt_dense->rmap->n*RARt_dense->cmap->n);
136     ierr = PetscInfo6(*C,"RARt_den %D %D; Rt_den %D %D, (RARt->nz %D)/(m*ncolors)=%g\n",RARt_dense->rmap->n,RARt_dense->cmap->n,Rt_dense->rmap->n,Rt_dense->cmap->n,c->nz,density);CHKERRQ(ierr);
137 #if defined(RART_PROFILE)
138     ierr = PetscInfo5(*C,"Sym = GetColor %g + MColorCreate %g + MDenCreate %g + other %g = %g\n",GColor,MCCreate,MDenCreate,etime,GColor+MCCreate+MDenCreate+etime);CHKERRQ(ierr);
139 #endif
140   }
141 #endif
142   PetscFunctionReturn(0);
143 }
144 
145 /*
146  RAB = R * A * B, R and A in seqaij format, B in dense format;
147 */
148 #undef __FUNCT__
149 #define __FUNCT__ "MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense"
150 PetscErrorCode MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense(Mat R,Mat A,Mat B,Mat RAB,PetscScalar *work)
151 {
152   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*r=(Mat_SeqAIJ*)R->data;
153   PetscErrorCode ierr;
154   PetscScalar    *b,r1,r2,r3,r4,*b1,*b2,*b3,*b4;
155   MatScalar      *aa,*ra;
156   PetscInt       cn =B->cmap->n,bm=B->rmap->n,col,i,j,n,*ai=a->i,*aj,am=A->rmap->n;
157   PetscInt       am2=2*am,am3=3*am,bm4=4*bm;
158   PetscScalar    *d,*c,*c2,*c3,*c4;
159   PetscInt       *rj,rm=R->rmap->n,dm=RAB->rmap->n,dn=RAB->cmap->n;
160   PetscInt       rm2=2*rm,rm3=3*rm,colrm;
161 #if defined(RART_PROFILE)
162   PetscLogDouble t0,tf,Mult_sp_den1=0.0,Mult_sp_den2=0.0;
163 #endif
164 
165   PetscFunctionBegin;
166   if (!dm || !dn) PetscFunctionReturn(0);
167   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);
168   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);
169   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);
170   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);
171 
172   { /*
173      This approach is not as good as original ones (will be removed later), but it reveals that
174      AB_den=A*B takes almost all execution time in R*A*B for src/ksp/ksp/examples/tutorials/ex56.c
175      */
176     PetscBool via_matmatmult=PETSC_FALSE;
177     ierr = PetscOptionsGetBool(NULL,"-matrart_via_matmatmult",&via_matmatmult,NULL);CHKERRQ(ierr);
178     if (via_matmatmult) {
179       Mat AB_den;
180       ierr = MatMatMultSymbolic_SeqAIJ_SeqDense(A,B,0.0,&AB_den);CHKERRQ(ierr);
181 
182 #if defined(RART_PROFILE)
183       ierr = PetscTime(&t0);CHKERRQ(ierr);
184 #endif
185       ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,B,AB_den);CHKERRQ(ierr);
186 #if defined(RART_PROFILE)
187       ierr  = PetscTime(&tf);CHKERRQ(ierr);
188       Mult_sp_den1 += tf - t0;
189       ierr = PetscTime(&t0);CHKERRQ(ierr);
190 #endif
191       ierr = MatMatMultNumeric_SeqAIJ_SeqDense(R,AB_den,RAB);CHKERRQ(ierr);
192 #if defined(RART_PROFILE)
193       ierr  = PetscTime(&tf);CHKERRQ(ierr);
194       Mult_sp_den2 += tf - t0;
195       printf(" MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense...time = %g + %g = %g; ncolor: %d\n",Mult_sp_den1, Mult_sp_den2, Mult_sp_den1+ Mult_sp_den2,B->cmap->n);
196 #endif
197       ierr = MatDestroy(&AB_den);CHKERRQ(ierr);
198       PetscFunctionReturn(0);
199     }
200   }
201 
202   ierr = MatDenseGetArray(B,&b);CHKERRQ(ierr);
203   ierr = MatDenseGetArray(RAB,&d);CHKERRQ(ierr);
204   b1   = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;
205   c    = work; c2 = c + am; c3 = c2 + am; c4 = c3 + am;
206   for (col=0; col<cn-4; col += 4) {  /* over columns of C */
207     for (i=0; i<am; i++) {        /* over rows of A in those columns */
208       r1 = r2 = r3 = r4 = 0.0;
209       n  = ai[i+1] - ai[i];
210       aj = a->j + ai[i];
211       aa = a->a + ai[i];
212       for (j=0; j<n; j++) {
213         r1 += (*aa)*b1[*aj];
214         r2 += (*aa)*b2[*aj];
215         r3 += (*aa)*b3[*aj];
216         r4 += (*aa++)*b4[*aj++];
217       }
218       c[i]       = r1;
219       c[am  + i] = r2;
220       c[am2 + i] = r3;
221       c[am3 + i] = r4;
222     }
223     b1 += bm4;
224     b2 += bm4;
225     b3 += bm4;
226     b4 += bm4;
227 
228     /* RAB[:,col] = R*C[:,col] */
229     colrm = col*rm;
230     for (i=0; i<rm; i++) {        /* over rows of R in those columns */
231       r1 = r2 = r3 = r4 = 0.0;
232       n  = r->i[i+1] - r->i[i];
233       rj = r->j + r->i[i];
234       ra = r->a + r->i[i];
235       for (j=0; j<n; j++) {
236         r1 += (*ra)*c[*rj];
237         r2 += (*ra)*c2[*rj];
238         r3 += (*ra)*c3[*rj];
239         r4 += (*ra++)*c4[*rj++];
240       }
241       d[colrm + i]       = r1;
242       d[colrm + rm + i]  = r2;
243       d[colrm + rm2 + i] = r3;
244       d[colrm + rm3 + i] = r4;
245     }
246   }
247   for (; col<cn; col++) {     /* over extra columns of C */
248     for (i=0; i<am; i++) {  /* over rows of A in those columns */
249       r1 = 0.0;
250       n  = a->i[i+1] - a->i[i];
251       aj = a->j + a->i[i];
252       aa = a->a + a->i[i];
253       for (j=0; j<n; j++) {
254         r1 += (*aa++)*b1[*aj++];
255       }
256       c[i] = r1;
257     }
258     b1 += bm;
259 
260     for (i=0; i<rm; i++) {  /* over rows of R in those columns */
261       r1 = 0.0;
262       n  = r->i[i+1] - r->i[i];
263       rj = r->j + r->i[i];
264       ra = r->a + r->i[i];
265       for (j=0; j<n; j++) {
266         r1 += (*ra++)*c[*rj++];
267       }
268       d[col*rm + i] = r1;
269     }
270   }
271   ierr = PetscLogFlops(cn*2.0*(a->nz + r->nz));CHKERRQ(ierr);
272 
273   ierr = MatDenseRestoreArray(B,&b);CHKERRQ(ierr);
274   ierr = MatDenseRestoreArray(RAB,&d);CHKERRQ(ierr);
275   ierr = MatAssemblyBegin(RAB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
276   ierr = MatAssemblyEnd(RAB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
277   PetscFunctionReturn(0);
278 }
279 
280 #undef __FUNCT__
281 #define __FUNCT__ "MatRARtNumeric_SeqAIJ_SeqAIJ_colorrart"
282 PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_colorrart(Mat A,Mat R,Mat C)
283 {
284   PetscErrorCode       ierr;
285   Mat_SeqAIJ           *c = (Mat_SeqAIJ*)C->data;
286   Mat_RARt             *rart=c->rart;
287   MatTransposeColoring matcoloring;
288   Mat                  Rt,RARt;
289 #if defined(RART_PROFILE)
290   PetscLogDouble       Mult_sp_den=0.0,app1=0.0,app2=0.0,t0,tf;
291 #endif
292 
293   PetscFunctionBegin;
294   /* Get dense Rt by Apply MatTransposeColoring to R */
295   matcoloring = rart->matcoloring;
296   Rt          = rart->Rt;
297 
298 #if defined(RART_PROFILE)
299   ierr  = PetscTime(&t0);CHKERRQ(ierr);
300 #endif
301   ierr  = MatTransColoringApplySpToDen(matcoloring,R,Rt);CHKERRQ(ierr);
302 #if defined(RART_PROFILE)
303   ierr  = PetscTime(&tf);CHKERRQ(ierr);
304   app1 += tf - t0;
305 #endif
306 
307   /* Get dense RARt = R*A*Rt -- dominates! */
308 #if defined(RART_PROFILE)
309   ierr = PetscTime(&t0);CHKERRQ(ierr);
310 #endif
311   RARt = rart->RARt;
312   ierr = MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense(R,A,Rt,RARt,rart->work);CHKERRQ(ierr);
313 #if defined(RART_PROFILE)
314   ierr = PetscTime(&tf);CHKERRQ(ierr);
315   Mult_sp_den += tf - t0;
316 #endif
317 
318   /* Recover C from C_dense */
319 #if defined(RART_PROFILE)
320   ierr = PetscTime(&t0);CHKERRQ(ierr);
321 #endif
322   ierr = MatTransColoringApplyDenToSp(matcoloring,RARt,C);CHKERRQ(ierr);
323 #if defined(RART_PROFILE)
324   ierr = PetscTime(&tf);CHKERRQ(ierr);
325   app2 += tf - t0;
326 #endif
327 
328 #if defined(RART_PROFILE)
329   printf("MatRARtNumeric_SeqAIJ_SeqAIJ: ColorApp %g + %g + Mult_sp_den %g  = %g\n",app1,app2,Mult_sp_den,app1+app2+Mult_sp_den);CHKERRQ(ierr);
330 #endif
331 #if defined(PETSC_USE_INFO)
332   ierr = PetscInfo(C,"C=R*(A*Rt) via coloring C - use sparse-dense inner products\n");CHKERRQ(ierr);
333 #endif
334   PetscFunctionReturn(0);
335 }
336 
337 #undef __FUNCT__
338 #define __FUNCT__ "MatRARtSymbolic_SeqAIJ_SeqAIJ_matmattransposemult"
339 PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_matmattransposemult(Mat A,Mat R,PetscReal fill,Mat *C)
340 {
341   PetscErrorCode  ierr;
342   Mat             ARt,RARt;
343   Mat_SeqAIJ     *c;
344   Mat_RARt       *rart;
345 
346   PetscFunctionBegin;
347   /* must use '-mat_no_inode' with '-matmattransmult_color 1' - do not knwo why? */
348   ierr = MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(A,R,fill,&ARt);CHKERRQ(ierr);
349   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(R,ARt,fill,&RARt);CHKERRQ(ierr);
350   *C                     = RARt;
351   RARt->ops->rartnumeric = MatRARtNumeric_SeqAIJ_SeqAIJ_matmattransposemult;
352 
353   ierr = PetscNew(Mat_RARt,&rart);CHKERRQ(ierr);
354   c         = (Mat_SeqAIJ*)(*C)->data;
355   c->rart   = rart;
356   rart->ARt = ARt;
357   rart->destroy      = RARt->ops->destroy;
358   RARt->ops->destroy = MatDestroy_SeqAIJ_RARt;
359   PetscFunctionReturn(0);
360 }
361 
362 #undef __FUNCT__
363 #define __FUNCT__ "MatRARtNumeric_SeqAIJ_SeqAIJ_matmattransposemult"
364 PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_matmattransposemult(Mat A,Mat R,Mat C)
365 {
366   PetscErrorCode  ierr;
367   Mat_SeqAIJ      *c=(Mat_SeqAIJ*)C->data;
368   Mat_RARt        *rart=c->rart;
369   Mat             ARt=rart->ARt;
370 #if defined(RART_PROFILE)
371   PetscLogDouble  t0,t1,t2;
372 #endif
373 
374   PetscFunctionBegin;
375 #if defined(RART_PROFILE)
376     ierr = PetscTime(&t0);CHKERRQ(ierr);
377 #endif
378     ierr = MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(A,R,ARt);CHKERRQ(ierr); /* dominate! */
379 #if defined(RART_PROFILE)
380     ierr = PetscTime(&t1);CHKERRQ(ierr);
381 #endif
382     ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ(R,ARt,C);CHKERRQ(ierr);
383 #if defined(RART_PROFILE)
384     ierr = PetscTime(&t2);CHKERRQ(ierr);
385     printf(" matrart_color_art_num = %g + %g = %g\n",t1-t0,t2-t1,t2-t0);
386 #endif
387 
388 #if defined(PETSC_USE_INFO)
389     ierr = PetscInfo(C,"Use ARt=A*R^T, C=R*ARt via MatMatTransposeMult(). Coloring can be applied to A*R^T.\n");CHKERRQ(ierr);
390 #endif
391   PetscFunctionReturn(0);
392 }
393 
394 #undef __FUNCT__
395 #define __FUNCT__ "MatRARtSymbolic_SeqAIJ_SeqAIJ"
396 PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat R,PetscReal fill,Mat *C)
397 {
398   PetscErrorCode  ierr;
399   Mat             Rt;
400   Mat_SeqAIJ      *c;
401   Mat_RARt        *rart;
402 
403   PetscFunctionBegin;
404   ierr = MatTranspose_SeqAIJ(R,MAT_INITIAL_MATRIX,&Rt);CHKERRQ(ierr);
405   ierr = MatMatMatMultSymbolic_SeqAIJ_SeqAIJ_SeqAIJ(R,A,Rt,fill,C);CHKERRQ(ierr);
406 
407   ierr = PetscNew(Mat_RARt,&rart);CHKERRQ(ierr);
408   rart->Rt = Rt;
409   c        = (Mat_SeqAIJ*)(*C)->data;
410   c->rart  = rart;
411   rart->destroy          = (*C)->ops->destroy;
412   (*C)->ops->destroy     = MatDestroy_SeqAIJ_RARt;
413   (*C)->ops->rartnumeric = MatRARtNumeric_SeqAIJ_SeqAIJ;
414   PetscFunctionReturn(0);
415 }
416 
417 #undef __FUNCT__
418 #define __FUNCT__ "MatRARtNumeric_SeqAIJ_SeqAIJ"
419 PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ(Mat A,Mat R,Mat C)
420 {
421   PetscErrorCode  ierr;
422   Mat_SeqAIJ      *c = (Mat_SeqAIJ*)C->data;
423   Mat_RARt        *rart = c->rart;
424   Mat             Rt = rart->Rt;
425 
426   PetscFunctionBegin;
427   ierr = MatTranspose_SeqAIJ(R,MAT_REUSE_MATRIX,&Rt);CHKERRQ(ierr);
428   ierr = MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqAIJ(R,A,Rt,C);CHKERRQ(ierr);
429 #if defined(PETSC_USE_INFO)
430   ierr = PetscInfo(C,"Use Rt=R^T and C=R*A*Rt via MatMatMatMult() to avoid sparse inner products\n");CHKERRQ(ierr);
431 #endif
432   PetscFunctionReturn(0);
433 }
434 
435 #undef __FUNCT__
436 #define __FUNCT__ "MatRARt_SeqAIJ_SeqAIJ"
437 PetscErrorCode MatRARt_SeqAIJ_SeqAIJ(Mat A,Mat R,MatReuse scall,PetscReal fill,Mat *C)
438 {
439   PetscErrorCode ierr;
440   const char     *algTypes[3] = {"matmatmatmult","matmattransposemult","coloring_rart"};
441   PetscInt       alg=0; /* set default algorithm */
442 
443   PetscFunctionBegin;
444   if (scall == MAT_INITIAL_MATRIX) {
445     ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
446     ierr = PetscOptionsEList("-matrart_via","Algorithmic approach","MatRARt",algTypes,3,algTypes[0],&alg,NULL);CHKERRQ(ierr);
447     ierr = PetscOptionsEnd();CHKERRQ(ierr);
448 
449     ierr = PetscLogEventBegin(MAT_RARtSymbolic,A,R,0,0);CHKERRQ(ierr);
450     switch (alg) {
451     case 1:
452       /* via matmattransposemult: ARt=A*R^T, C=R*ARt - matrix coloring can be applied to A*R^T */
453       ierr = MatRARtSymbolic_SeqAIJ_SeqAIJ_matmattransposemult(A,R,fill,C);CHKERRQ(ierr);
454       break;
455     case 2:
456       /* via coloring_rart: apply coloring C = R*A*R^T                          */
457       ierr = MatRARtSymbolic_SeqAIJ_SeqAIJ_colorrart(A,R,fill,C);CHKERRQ(ierr);
458       break;
459     default:
460       /* via matmatmatmult: Rt=R^T, C=R*A*Rt - avoid inefficient sparse inner products */
461       ierr = MatRARtSymbolic_SeqAIJ_SeqAIJ(A,R,fill,C);CHKERRQ(ierr);
462       break;
463     }
464     ierr = PetscLogEventEnd(MAT_RARtSymbolic,A,R,0,0);CHKERRQ(ierr);
465   }
466 
467   ierr = PetscLogEventBegin(MAT_RARtNumeric,A,R,0,0);CHKERRQ(ierr);
468   ierr = (*(*C)->ops->rartnumeric)(A,R,*C);CHKERRQ(ierr);
469   ierr = PetscLogEventEnd(MAT_RARtNumeric,A,R,0,0);CHKERRQ(ierr);
470   PetscFunctionReturn(0);
471 }
472