xref: /petsc/src/mat/impls/aij/seq/matrart.c (revision 2c71b3e237ead271e4f3aa1505f92bf476e3413d)
1807171c5SHong Zhang 
23b1b9624SHong Zhang /*
33b1b9624SHong Zhang   Defines projective product routines where A is a SeqAIJ matrix
43b1b9624SHong Zhang           C = R * A * R^T
53b1b9624SHong Zhang */
63b1b9624SHong Zhang 
7807171c5SHong Zhang #include <../src/mat/impls/aij/seq/aij.h>
8807171c5SHong Zhang #include <../src/mat/utils/freespace.h>
9807171c5SHong Zhang #include <../src/mat/impls/dense/seq/dense.h> /*I "petscmat.h" I*/
10807171c5SHong Zhang 
116718818eSStefano Zampini PetscErrorCode MatDestroy_SeqAIJ_RARt(void *data)
12807171c5SHong Zhang {
13807171c5SHong Zhang   PetscErrorCode ierr;
146718818eSStefano Zampini   Mat_RARt       *rart = (Mat_RARt*)data;
15807171c5SHong Zhang 
16807171c5SHong Zhang   PetscFunctionBegin;
17807171c5SHong Zhang   ierr = MatTransposeColoringDestroy(&rart->matcoloring);CHKERRQ(ierr);
18807171c5SHong Zhang   ierr = MatDestroy(&rart->Rt);CHKERRQ(ierr);
19807171c5SHong Zhang   ierr = MatDestroy(&rart->RARt);CHKERRQ(ierr);
203b1b9624SHong Zhang   ierr = MatDestroy(&rart->ARt);CHKERRQ(ierr);
21807171c5SHong Zhang   ierr = PetscFree(rart->work);CHKERRQ(ierr);
226718818eSStefano Zampini   if (rart->destroy) {
236718818eSStefano Zampini     ierr = (*rart->destroy)(rart->data);CHKERRQ(ierr);
24807171c5SHong Zhang   }
25257c235dSHong Zhang   ierr = PetscFree(rart);CHKERRQ(ierr);
26807171c5SHong Zhang   PetscFunctionReturn(0);
27807171c5SHong Zhang }
28807171c5SHong Zhang 
294222ddf1SHong Zhang PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_colorrart(Mat A,Mat R,PetscReal fill,Mat C)
30807171c5SHong Zhang {
31807171c5SHong Zhang   PetscErrorCode       ierr;
32807171c5SHong Zhang   Mat                  P;
33807171c5SHong Zhang   PetscInt             *rti,*rtj;
34807171c5SHong Zhang   Mat_RARt             *rart;
35335efc43SPeter Brune   MatColoring          coloring;
36807171c5SHong Zhang   MatTransposeColoring matcoloring;
37807171c5SHong Zhang   ISColoring           iscoloring;
388d0a38e4SHong Zhang   Mat                  Rt_dense,RARt_dense;
39807171c5SHong Zhang 
40807171c5SHong Zhang   PetscFunctionBegin;
416718818eSStefano Zampini   MatCheckProduct(C,4);
42*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(C->product->data,PetscObjectComm((PetscObject)C),PETSC_ERR_PLIB,"Product data not empty");
43807171c5SHong Zhang   /* create symbolic P=Rt */
44807171c5SHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(R,&rti,&rtj);CHKERRQ(ierr);
450298fd71SBarry Smith   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,R->cmap->n,R->rmap->n,rti,rtj,NULL,&P);CHKERRQ(ierr);
46807171c5SHong Zhang 
47807171c5SHong Zhang   /* get symbolic C=Pt*A*P */
484a1b09b7SHong Zhang   ierr = MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(A,P,fill,C);CHKERRQ(ierr);
494222ddf1SHong Zhang   ierr = MatSetBlockSizes(C,PetscAbs(R->rmap->bs),PetscAbs(R->rmap->bs));CHKERRQ(ierr);
504222ddf1SHong Zhang   C->ops->rartnumeric = MatRARtNumeric_SeqAIJ_SeqAIJ_colorrart;
51807171c5SHong Zhang 
52807171c5SHong Zhang   /* create a supporting struct */
53b00a9115SJed Brown   ierr = PetscNew(&rart);CHKERRQ(ierr);
546718818eSStefano Zampini   C->product->data    = rart;
556718818eSStefano Zampini   C->product->destroy = MatDestroy_SeqAIJ_RARt;
56807171c5SHong Zhang 
5722e94b5dSHong Zhang   /* ------ Use coloring ---------- */
584222ddf1SHong Zhang   /* inode causes memory problem */
594222ddf1SHong Zhang   ierr = MatSetOption(C,MAT_USE_INODES,PETSC_FALSE);CHKERRQ(ierr);
604d478ae7SHong Zhang 
61807171c5SHong Zhang   /* Create MatTransposeColoring from symbolic C=R*A*R^T */
624222ddf1SHong Zhang   ierr = MatColoringCreate(C,&coloring);CHKERRQ(ierr);
63335efc43SPeter Brune   ierr = MatColoringSetDistance(coloring,2);CHKERRQ(ierr);
64335efc43SPeter Brune   ierr = MatColoringSetType(coloring,MATCOLORINGSL);CHKERRQ(ierr);
65335efc43SPeter Brune   ierr = MatColoringSetFromOptions(coloring);CHKERRQ(ierr);
66335efc43SPeter Brune   ierr = MatColoringApply(coloring,&iscoloring);CHKERRQ(ierr);
67335efc43SPeter Brune   ierr = MatColoringDestroy(&coloring);CHKERRQ(ierr);
684222ddf1SHong Zhang   ierr = MatTransposeColoringCreate(C,iscoloring,&matcoloring);CHKERRQ(ierr);
692205254eSKarl Rupp 
70807171c5SHong Zhang   rart->matcoloring = matcoloring;
71807171c5SHong Zhang   ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr);
723b1b9624SHong Zhang 
73807171c5SHong Zhang   /* Create Rt_dense */
74807171c5SHong Zhang   ierr = MatCreate(PETSC_COMM_SELF,&Rt_dense);CHKERRQ(ierr);
75807171c5SHong Zhang   ierr = MatSetSizes(Rt_dense,A->cmap->n,matcoloring->ncolors,A->cmap->n,matcoloring->ncolors);CHKERRQ(ierr);
76807171c5SHong Zhang   ierr = MatSetType(Rt_dense,MATSEQDENSE);CHKERRQ(ierr);
770298fd71SBarry Smith   ierr = MatSeqDenseSetPreallocation(Rt_dense,NULL);CHKERRQ(ierr);
782205254eSKarl Rupp 
79807171c5SHong Zhang   Rt_dense->assembled = PETSC_TRUE;
80807171c5SHong Zhang   rart->Rt            = Rt_dense;
81807171c5SHong Zhang 
82807171c5SHong Zhang   /* Create RARt_dense = R*A*Rt_dense */
83807171c5SHong Zhang   ierr = MatCreate(PETSC_COMM_SELF,&RARt_dense);CHKERRQ(ierr);
844222ddf1SHong Zhang   ierr = MatSetSizes(RARt_dense,C->rmap->n,matcoloring->ncolors,C->rmap->n,matcoloring->ncolors);CHKERRQ(ierr);
85807171c5SHong Zhang   ierr = MatSetType(RARt_dense,MATSEQDENSE);CHKERRQ(ierr);
860298fd71SBarry Smith   ierr = MatSeqDenseSetPreallocation(RARt_dense,NULL);CHKERRQ(ierr);
872205254eSKarl Rupp 
88807171c5SHong Zhang   rart->RARt = RARt_dense;
89807171c5SHong Zhang 
90807171c5SHong Zhang   /* Allocate work array to store columns of A*R^T used in MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense() */
91785e854fSJed Brown   ierr = PetscMalloc1(A->rmap->n*4,&rart->work);CHKERRQ(ierr);
92807171c5SHong Zhang 
93807171c5SHong Zhang   /* clean up */
94807171c5SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(R,&rti,&rtj);CHKERRQ(ierr);
95807171c5SHong Zhang   ierr = MatDestroy(&P);CHKERRQ(ierr);
96807171c5SHong Zhang 
97807171c5SHong Zhang #if defined(PETSC_USE_INFO)
981ce71dffSSatish Balay   {
996718818eSStefano Zampini     Mat_SeqAIJ *c = (Mat_SeqAIJ*)C->data;
100807171c5SHong Zhang     PetscReal density = (PetscReal)(c->nz)/(RARt_dense->rmap->n*RARt_dense->cmap->n);
1014222ddf1SHong Zhang     ierr = PetscInfo(C,"C=R*(A*Rt) via coloring C - use sparse-dense inner products\n");CHKERRQ(ierr);
1027d3de750SJacob Faibussowitsch     ierr = PetscInfo(C,"RARt_den %" PetscInt_FMT " %" PetscInt_FMT "; Rt %" PetscInt_FMT " %" PetscInt_FMT " (RARt->nz %" PetscInt_FMT ")/(m*ncolors)=%g\n",RARt_dense->rmap->n,RARt_dense->cmap->n,R->cmap->n,R->rmap->n,c->nz,density);CHKERRQ(ierr);
1031ce71dffSSatish Balay   }
104807171c5SHong Zhang #endif
105807171c5SHong Zhang   PetscFunctionReturn(0);
106807171c5SHong Zhang }
107807171c5SHong Zhang 
108807171c5SHong Zhang /*
109807171c5SHong Zhang  RAB = R * A * B, R and A in seqaij format, B in dense format;
110807171c5SHong Zhang */
111807171c5SHong Zhang PetscErrorCode MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense(Mat R,Mat A,Mat B,Mat RAB,PetscScalar *work)
112807171c5SHong Zhang {
113807171c5SHong Zhang   Mat_SeqAIJ        *a=(Mat_SeqAIJ*)A->data,*r=(Mat_SeqAIJ*)R->data;
114807171c5SHong Zhang   PetscErrorCode    ierr;
1151683a169SBarry Smith   PetscScalar       r1,r2,r3,r4;
1161683a169SBarry Smith   const PetscScalar *b,*b1,*b2,*b3,*b4;
117807171c5SHong Zhang   MatScalar         *aa,*ra;
118807171c5SHong Zhang   PetscInt          cn =B->cmap->n,bm=B->rmap->n,col,i,j,n,*ai=a->i,*aj,am=A->rmap->n;
119807171c5SHong Zhang   PetscInt          am2=2*am,am3=3*am,bm4=4*bm;
120807171c5SHong Zhang   PetscScalar       *d,*c,*c2,*c3,*c4;
121807171c5SHong Zhang   PetscInt          *rj,rm=R->rmap->n,dm=RAB->rmap->n,dn=RAB->cmap->n;
122807171c5SHong Zhang   PetscInt         rm2=2*rm,rm3=3*rm,colrm;
123807171c5SHong Zhang 
124807171c5SHong Zhang   PetscFunctionBegin;
125807171c5SHong Zhang   if (!dm || !dn) PetscFunctionReturn(0);
126*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(bm != A->cmap->n,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Number columns in A %" PetscInt_FMT " not equal rows in B %" PetscInt_FMT,A->cmap->n,bm);
127*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(am != R->cmap->n,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Number columns in R %" PetscInt_FMT " not equal rows in A %" PetscInt_FMT,R->cmap->n,am);
128*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(R->rmap->n != RAB->rmap->n,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Number rows in RAB %" PetscInt_FMT " not equal rows in R %" PetscInt_FMT,RAB->rmap->n,R->rmap->n);
129*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(B->cmap->n != RAB->cmap->n,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Number columns in RAB %" PetscInt_FMT " not equal columns in B %" PetscInt_FMT,RAB->cmap->n,B->cmap->n);
130807171c5SHong Zhang 
131274010c0SHong Zhang   { /*
132274010c0SHong Zhang      This approach is not as good as original ones (will be removed later), but it reveals that
133c4762a1bSJed Brown      AB_den=A*B takes almost all execution time in R*A*B for src/ksp/ksp/tutorials/ex56.c
134274010c0SHong Zhang      */
135274010c0SHong Zhang     PetscBool via_matmatmult=PETSC_FALSE;
136c5929fdfSBarry Smith     ierr = PetscOptionsGetBool(NULL,NULL,"-matrart_via_matmatmult",&via_matmatmult,NULL);CHKERRQ(ierr);
137274010c0SHong Zhang     if (via_matmatmult) {
1384222ddf1SHong Zhang       Mat AB_den = NULL;
1394222ddf1SHong Zhang       ierr = MatCreate(PetscObjectComm((PetscObject)A),&AB_den);CHKERRQ(ierr);
1404222ddf1SHong Zhang       ierr = MatMatMultSymbolic_SeqAIJ_SeqDense(A,B,0.0,AB_den);CHKERRQ(ierr);
141274010c0SHong Zhang       ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,B,AB_den);CHKERRQ(ierr);
142274010c0SHong Zhang       ierr = MatMatMultNumeric_SeqAIJ_SeqDense(R,AB_den,RAB);CHKERRQ(ierr);
143274010c0SHong Zhang       ierr = MatDestroy(&AB_den);CHKERRQ(ierr);
144274010c0SHong Zhang       PetscFunctionReturn(0);
145274010c0SHong Zhang     }
146274010c0SHong Zhang   }
147274010c0SHong Zhang 
1481683a169SBarry Smith   ierr = MatDenseGetArrayRead(B,&b);CHKERRQ(ierr);
1498c778c55SBarry Smith   ierr = MatDenseGetArray(RAB,&d);CHKERRQ(ierr);
150807171c5SHong Zhang   b1   = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;
151807171c5SHong Zhang   c    = work; c2 = c + am; c3 = c2 + am; c4 = c3 + am;
152807171c5SHong Zhang   for (col=0; col<cn-4; col += 4) {  /* over columns of C */
153807171c5SHong Zhang     for (i=0; i<am; i++) {        /* over rows of A in those columns */
154807171c5SHong Zhang       r1 = r2 = r3 = r4 = 0.0;
155807171c5SHong Zhang       n  = ai[i+1] - ai[i];
156807171c5SHong Zhang       aj = a->j + ai[i];
157807171c5SHong Zhang       aa = a->a + ai[i];
158807171c5SHong Zhang       for (j=0; j<n; j++) {
159807171c5SHong Zhang         r1 += (*aa)*b1[*aj];
160807171c5SHong Zhang         r2 += (*aa)*b2[*aj];
161807171c5SHong Zhang         r3 += (*aa)*b3[*aj];
162807171c5SHong Zhang         r4 += (*aa++)*b4[*aj++];
163807171c5SHong Zhang       }
164807171c5SHong Zhang       c[i]       = r1;
165807171c5SHong Zhang       c[am  + i] = r2;
166807171c5SHong Zhang       c[am2 + i] = r3;
167807171c5SHong Zhang       c[am3 + i] = r4;
168807171c5SHong Zhang     }
169807171c5SHong Zhang     b1 += bm4;
170807171c5SHong Zhang     b2 += bm4;
171807171c5SHong Zhang     b3 += bm4;
172807171c5SHong Zhang     b4 += bm4;
173807171c5SHong Zhang 
174807171c5SHong Zhang     /* RAB[:,col] = R*C[:,col] */
175807171c5SHong Zhang     colrm = col*rm;
176807171c5SHong Zhang     for (i=0; i<rm; i++) {        /* over rows of R in those columns */
177807171c5SHong Zhang       r1 = r2 = r3 = r4 = 0.0;
178807171c5SHong Zhang       n  = r->i[i+1] - r->i[i];
179807171c5SHong Zhang       rj = r->j + r->i[i];
180807171c5SHong Zhang       ra = r->a + r->i[i];
181807171c5SHong Zhang       for (j=0; j<n; j++) {
182807171c5SHong Zhang         r1 += (*ra)*c[*rj];
183807171c5SHong Zhang         r2 += (*ra)*c2[*rj];
184807171c5SHong Zhang         r3 += (*ra)*c3[*rj];
185807171c5SHong Zhang         r4 += (*ra++)*c4[*rj++];
186807171c5SHong Zhang       }
187807171c5SHong Zhang       d[colrm + i]       = r1;
188807171c5SHong Zhang       d[colrm + rm + i]  = r2;
189807171c5SHong Zhang       d[colrm + rm2 + i] = r3;
190807171c5SHong Zhang       d[colrm + rm3 + i] = r4;
191807171c5SHong Zhang     }
192807171c5SHong Zhang   }
193807171c5SHong Zhang   for (; col<cn; col++) {     /* over extra columns of C */
194807171c5SHong Zhang     for (i=0; i<am; i++) {  /* over rows of A in those columns */
195807171c5SHong Zhang       r1 = 0.0;
196807171c5SHong Zhang       n  = a->i[i+1] - a->i[i];
197807171c5SHong Zhang       aj = a->j + a->i[i];
198807171c5SHong Zhang       aa = a->a + a->i[i];
199807171c5SHong Zhang       for (j=0; j<n; j++) {
200807171c5SHong Zhang         r1 += (*aa++)*b1[*aj++];
201807171c5SHong Zhang       }
202807171c5SHong Zhang       c[i] = r1;
203807171c5SHong Zhang     }
204807171c5SHong Zhang     b1 += bm;
205807171c5SHong Zhang 
206807171c5SHong Zhang     for (i=0; i<rm; i++) {  /* over rows of R in those columns */
207807171c5SHong Zhang       r1 = 0.0;
208807171c5SHong Zhang       n  = r->i[i+1] - r->i[i];
209807171c5SHong Zhang       rj = r->j + r->i[i];
210807171c5SHong Zhang       ra = r->a + r->i[i];
211807171c5SHong Zhang       for (j=0; j<n; j++) {
212807171c5SHong Zhang         r1 += (*ra++)*c[*rj++];
213807171c5SHong Zhang       }
214807171c5SHong Zhang       d[col*rm + i] = r1;
215807171c5SHong Zhang     }
216807171c5SHong Zhang   }
217807171c5SHong Zhang   ierr = PetscLogFlops(cn*2.0*(a->nz + r->nz));CHKERRQ(ierr);
218807171c5SHong Zhang 
2191683a169SBarry Smith   ierr = MatDenseRestoreArrayRead(B,&b);CHKERRQ(ierr);
2208c778c55SBarry Smith   ierr = MatDenseRestoreArray(RAB,&d);CHKERRQ(ierr);
221807171c5SHong Zhang   ierr = MatAssemblyBegin(RAB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
222807171c5SHong Zhang   ierr = MatAssemblyEnd(RAB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
223807171c5SHong Zhang   PetscFunctionReturn(0);
224807171c5SHong Zhang }
225807171c5SHong Zhang 
22655bea0ebSHong Zhang PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_colorrart(Mat A,Mat R,Mat C)
227807171c5SHong Zhang {
228807171c5SHong Zhang   PetscErrorCode       ierr;
2296718818eSStefano Zampini   Mat_RARt             *rart;
230807171c5SHong Zhang   MatTransposeColoring matcoloring;
2318d0a38e4SHong Zhang   Mat                  Rt,RARt;
232807171c5SHong Zhang 
233807171c5SHong Zhang   PetscFunctionBegin;
2346718818eSStefano Zampini   MatCheckProduct(C,3);
235*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!C->product->data,PetscObjectComm((PetscObject)C),PETSC_ERR_PLIB,"Product data empty");
2366718818eSStefano Zampini   rart = (Mat_RARt*)C->product->data;
2376718818eSStefano Zampini 
238807171c5SHong Zhang   /* Get dense Rt by Apply MatTransposeColoring to R */
239807171c5SHong Zhang   matcoloring = rart->matcoloring;
240807171c5SHong Zhang   Rt          = rart->Rt;
241807171c5SHong Zhang   ierr  = MatTransColoringApplySpToDen(matcoloring,R,Rt);CHKERRQ(ierr);
242807171c5SHong Zhang 
24350647e95SHong Zhang   /* Get dense RARt = R*A*Rt -- dominates! */
244807171c5SHong Zhang   RARt = rart->RARt;
245807171c5SHong Zhang   ierr = MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense(R,A,Rt,RARt,rart->work);CHKERRQ(ierr);
246807171c5SHong Zhang 
247807171c5SHong Zhang   /* Recover C from C_dense */
248807171c5SHong Zhang   ierr = MatTransColoringApplyDenToSp(matcoloring,RARt,C);CHKERRQ(ierr);
249807171c5SHong Zhang   PetscFunctionReturn(0);
250807171c5SHong Zhang }
251807171c5SHong Zhang 
2524222ddf1SHong Zhang PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_matmattransposemult(Mat A,Mat R,PetscReal fill,Mat C)
253807171c5SHong Zhang {
254807171c5SHong Zhang   PetscErrorCode ierr;
2554222ddf1SHong Zhang   Mat            ARt;
2564fa85fdeSHong Zhang   Mat_RARt       *rart;
2576718818eSStefano Zampini   char           *alg;
2584fa85fdeSHong Zhang 
25995a72cc5SHong Zhang   PetscFunctionBegin;
2606718818eSStefano Zampini   MatCheckProduct(C,4);
261*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(C->product->data,PetscObjectComm((PetscObject)C),PETSC_ERR_PLIB,"Product data not empty");
2624222ddf1SHong Zhang   /* create symbolic ARt = A*R^T  */
2634222ddf1SHong Zhang   ierr = MatProductCreate(A,R,NULL,&ARt);CHKERRQ(ierr);
2644222ddf1SHong Zhang   ierr = MatProductSetType(ARt,MATPRODUCT_ABt);CHKERRQ(ierr);
2654222ddf1SHong Zhang   ierr = MatProductSetAlgorithm(ARt,"sorted");CHKERRQ(ierr);
2664222ddf1SHong Zhang   ierr = MatProductSetFill(ARt,fill);CHKERRQ(ierr);
2674222ddf1SHong Zhang   ierr = MatProductSetFromOptions(ARt);CHKERRQ(ierr);
2684222ddf1SHong Zhang   ierr = MatProductSymbolic(ARt);CHKERRQ(ierr);
2694222ddf1SHong Zhang 
2704222ddf1SHong Zhang   /* compute symbolic C = R*ARt */
2717a3c3d58SStefano Zampini   /* set algorithm for C = R*ARt */
2726718818eSStefano Zampini   ierr = PetscStrallocpy(C->product->alg,&alg);CHKERRQ(ierr);
2737a3c3d58SStefano Zampini   ierr = MatProductSetAlgorithm(C,"sorted");CHKERRQ(ierr);
2744222ddf1SHong Zhang   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(R,ARt,fill,C);CHKERRQ(ierr);
2757a3c3d58SStefano Zampini   /* resume original algorithm for C */
2767a3c3d58SStefano Zampini   ierr = MatProductSetAlgorithm(C,alg);CHKERRQ(ierr);
2776718818eSStefano Zampini   ierr = PetscFree(alg);CHKERRQ(ierr);
2784222ddf1SHong Zhang 
2794222ddf1SHong Zhang   C->ops->rartnumeric = MatRARtNumeric_SeqAIJ_SeqAIJ_matmattransposemult;
28055bea0ebSHong Zhang 
281b00a9115SJed Brown   ierr = PetscNew(&rart);CHKERRQ(ierr);
2823b1b9624SHong Zhang   rart->ARt = ARt;
2836718818eSStefano Zampini   C->product->data    = rart;
2846718818eSStefano Zampini   C->product->destroy = MatDestroy_SeqAIJ_RARt;
2854222ddf1SHong Zhang   ierr = PetscInfo(C,"Use ARt=A*R^T, C=R*ARt via MatMatTransposeMult(). Coloring can be applied to A*R^T.\n");CHKERRQ(ierr);
28655bea0ebSHong Zhang   PetscFunctionReturn(0);
28736104f73SHong Zhang }
28855bea0ebSHong Zhang 
28955bea0ebSHong Zhang PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_matmattransposemult(Mat A,Mat R,Mat C)
29055bea0ebSHong Zhang {
29155bea0ebSHong Zhang   PetscErrorCode ierr;
2926718818eSStefano Zampini   Mat_RARt       *rart;
29355bea0ebSHong Zhang 
29455bea0ebSHong Zhang   PetscFunctionBegin;
2956718818eSStefano Zampini   MatCheckProduct(C,3);
296*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!C->product->data,PetscObjectComm((PetscObject)C),PETSC_ERR_PLIB,"Product data empty");
2976718818eSStefano Zampini   rart = (Mat_RARt*)C->product->data;
2986718818eSStefano Zampini   ierr = MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(A,R,rart->ARt);CHKERRQ(ierr); /* dominate! */
2996718818eSStefano Zampini   ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ(R,rart->ARt,C);CHKERRQ(ierr);
30055bea0ebSHong Zhang   PetscFunctionReturn(0);
301807171c5SHong Zhang }
30255bea0ebSHong Zhang 
3034222ddf1SHong Zhang PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat R,PetscReal fill,Mat C)
30455bea0ebSHong Zhang {
30555bea0ebSHong Zhang   PetscErrorCode ierr;
30695a72cc5SHong Zhang   Mat            Rt;
30755bea0ebSHong Zhang   Mat_RARt       *rart;
30855bea0ebSHong Zhang 
30955bea0ebSHong Zhang   PetscFunctionBegin;
3106718818eSStefano Zampini   MatCheckProduct(C,4);
311*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(C->product->data,PetscObjectComm((PetscObject)C),PETSC_ERR_PLIB,"Product data not empty");
3125008f5a7SHong Zhang   ierr = MatTranspose_SeqAIJ(R,MAT_INITIAL_MATRIX,&Rt);CHKERRQ(ierr);
3135008f5a7SHong Zhang   ierr = MatMatMatMultSymbolic_SeqAIJ_SeqAIJ_SeqAIJ(R,A,Rt,fill,C);CHKERRQ(ierr);
31495a72cc5SHong Zhang 
315b00a9115SJed Brown   ierr = PetscNew(&rart);CHKERRQ(ierr);
3166718818eSStefano Zampini   rart->data = C->product->data;
3176718818eSStefano Zampini   rart->destroy = C->product->destroy;
31895a72cc5SHong Zhang   rart->Rt = Rt;
3196718818eSStefano Zampini   C->product->data    = rart;
3206718818eSStefano Zampini   C->product->destroy = MatDestroy_SeqAIJ_RARt;
3214222ddf1SHong Zhang   C->ops->rartnumeric = MatRARtNumeric_SeqAIJ_SeqAIJ;
3224222ddf1SHong Zhang   ierr = PetscInfo(C,"Use Rt=R^T and C=R*A*Rt via MatMatMatMult() to avoid sparse inner products\n");CHKERRQ(ierr);
32355bea0ebSHong Zhang   PetscFunctionReturn(0);
32455bea0ebSHong Zhang }
32555bea0ebSHong Zhang 
32655bea0ebSHong Zhang PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ(Mat A,Mat R,Mat C)
32755bea0ebSHong Zhang {
32855bea0ebSHong Zhang   PetscErrorCode ierr;
3296718818eSStefano Zampini   Mat_RARt       *rart;
33055bea0ebSHong Zhang 
33155bea0ebSHong Zhang   PetscFunctionBegin;
3326718818eSStefano Zampini   MatCheckProduct(C,3);
333*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!C->product->data,PetscObjectComm((PetscObject)C),PETSC_ERR_PLIB,"Product data empty");
3346718818eSStefano Zampini   rart = (Mat_RARt*)C->product->data;
3356718818eSStefano Zampini   ierr = MatTranspose_SeqAIJ(R,MAT_REUSE_MATRIX,&rart->Rt);CHKERRQ(ierr);
3366718818eSStefano Zampini   /* MatMatMatMultSymbolic used a different data */
3376718818eSStefano Zampini   C->product->data = rart->data;
3386718818eSStefano Zampini   ierr = MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqAIJ(R,A,rart->Rt,C);CHKERRQ(ierr);
3396718818eSStefano Zampini   C->product->data = rart;
34055bea0ebSHong Zhang   PetscFunctionReturn(0);
34155bea0ebSHong Zhang }
34255bea0ebSHong Zhang 
34355bea0ebSHong Zhang PetscErrorCode MatRARt_SeqAIJ_SeqAIJ(Mat A,Mat R,MatReuse scall,PetscReal fill,Mat *C)
34455bea0ebSHong Zhang {
34555bea0ebSHong Zhang   PetscErrorCode ierr;
34655bea0ebSHong Zhang   const char     *algTypes[3] = {"matmatmatmult","matmattransposemult","coloring_rart"};
34755bea0ebSHong Zhang   PetscInt       alg=0; /* set default algorithm */
34855bea0ebSHong Zhang 
34955bea0ebSHong Zhang   PetscFunctionBegin;
35055bea0ebSHong Zhang   if (scall == MAT_INITIAL_MATRIX) {
351715a5346SHong Zhang     ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"MatRARt","Mat");CHKERRQ(ierr);
35255bea0ebSHong Zhang     ierr = PetscOptionsEList("-matrart_via","Algorithmic approach","MatRARt",algTypes,3,algTypes[0],&alg,NULL);CHKERRQ(ierr);
353b56132cbSHong Zhang     ierr = PetscOptionsEnd();CHKERRQ(ierr);
354b56132cbSHong Zhang 
35555bea0ebSHong Zhang     ierr = PetscLogEventBegin(MAT_RARtSymbolic,A,R,0,0);CHKERRQ(ierr);
3564222ddf1SHong Zhang     ierr = MatCreate(PETSC_COMM_SELF,C);CHKERRQ(ierr);
35755bea0ebSHong Zhang     switch (alg) {
35855bea0ebSHong Zhang     case 1:
35955bea0ebSHong Zhang       /* via matmattransposemult: ARt=A*R^T, C=R*ARt - matrix coloring can be applied to A*R^T */
3604222ddf1SHong Zhang       ierr = MatRARtSymbolic_SeqAIJ_SeqAIJ_matmattransposemult(A,R,fill,*C);CHKERRQ(ierr);
36155bea0ebSHong Zhang       break;
36255bea0ebSHong Zhang     case 2:
36355bea0ebSHong Zhang       /* via coloring_rart: apply coloring C = R*A*R^T                          */
3644222ddf1SHong Zhang       ierr = MatRARtSymbolic_SeqAIJ_SeqAIJ_colorrart(A,R,fill,*C);CHKERRQ(ierr);
36555bea0ebSHong Zhang       break;
36655bea0ebSHong Zhang     default:
36755bea0ebSHong Zhang       /* via matmatmatmult: Rt=R^T, C=R*A*Rt - avoid inefficient sparse inner products */
3684222ddf1SHong Zhang       ierr = MatRARtSymbolic_SeqAIJ_SeqAIJ(A,R,fill,*C);CHKERRQ(ierr);
36955bea0ebSHong Zhang       break;
37055bea0ebSHong Zhang     }
3715008f5a7SHong Zhang     ierr = PetscLogEventEnd(MAT_RARtSymbolic,A,R,0,0);CHKERRQ(ierr);
3725008f5a7SHong Zhang   }
3735008f5a7SHong Zhang 
3745008f5a7SHong Zhang   ierr = PetscLogEventBegin(MAT_RARtNumeric,A,R,0,0);CHKERRQ(ierr);
3754222ddf1SHong Zhang   ierr = ((*C)->ops->rartnumeric)(A,R,*C);CHKERRQ(ierr);
3765008f5a7SHong Zhang   ierr = PetscLogEventEnd(MAT_RARtNumeric,A,R,0,0);CHKERRQ(ierr);
377807171c5SHong Zhang   PetscFunctionReturn(0);
378807171c5SHong Zhang }
3794222ddf1SHong Zhang 
3804222ddf1SHong Zhang /* ------------------------------------------------------------- */
3814222ddf1SHong Zhang PetscErrorCode MatProductSymbolic_RARt_SeqAIJ_SeqAIJ(Mat C)
3824222ddf1SHong Zhang {
3834222ddf1SHong Zhang   PetscErrorCode      ierr;
3844222ddf1SHong Zhang   Mat_Product         *product = C->product;
3854222ddf1SHong Zhang   Mat                 A=product->A,R=product->B;
3864222ddf1SHong Zhang   MatProductAlgorithm alg=product->alg;
3874222ddf1SHong Zhang   PetscReal           fill=product->fill;
3884222ddf1SHong Zhang   PetscBool           flg;
3894222ddf1SHong Zhang 
3904222ddf1SHong Zhang   PetscFunctionBegin;
3914222ddf1SHong Zhang   ierr = PetscStrcmp(alg,"r*a*rt",&flg);CHKERRQ(ierr);
3924222ddf1SHong Zhang   if (flg) {
3934222ddf1SHong Zhang     ierr = MatRARtSymbolic_SeqAIJ_SeqAIJ(A,R,fill,C);CHKERRQ(ierr);
3944222ddf1SHong Zhang     goto next;
3954222ddf1SHong Zhang   }
3964222ddf1SHong Zhang 
3974222ddf1SHong Zhang   ierr = PetscStrcmp(alg,"r*art",&flg);CHKERRQ(ierr);
3984222ddf1SHong Zhang   if (flg) {
3994222ddf1SHong Zhang     ierr = MatRARtSymbolic_SeqAIJ_SeqAIJ_matmattransposemult(A,R,fill,C);CHKERRQ(ierr);
4004222ddf1SHong Zhang     goto next;
4014222ddf1SHong Zhang   }
4024222ddf1SHong Zhang 
4034222ddf1SHong Zhang   ierr = PetscStrcmp(alg,"coloring_rart",&flg);CHKERRQ(ierr);
4044222ddf1SHong Zhang   if (flg) {
4054222ddf1SHong Zhang     ierr = MatRARtSymbolic_SeqAIJ_SeqAIJ_colorrart(A,R,fill,C);CHKERRQ(ierr);
4064222ddf1SHong Zhang     goto next;
4074222ddf1SHong Zhang   }
4084222ddf1SHong Zhang 
4094222ddf1SHong Zhang   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatProductAlgorithm is not supported");
4104222ddf1SHong Zhang 
4114222ddf1SHong Zhang next:
4124222ddf1SHong Zhang   C->ops->productnumeric = MatProductNumeric_RARt;
4134222ddf1SHong Zhang   PetscFunctionReturn(0);
4144222ddf1SHong Zhang }
415