xref: /petsc/src/mat/impls/aij/seq/matrart.c (revision 6273346d76c04da314646f01a57ee4ec0081c8d9) !
1 
2 /*
3   Defines matrix-matrix product routines for pairs of SeqAIJ matrices
4           C = P * A * P^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 
11 /*
12      MatApplyPAPt_Symbolic_SeqAIJ_SeqAIJ - Forms the symbolic product of two SeqAIJ matrices
13            C = P * A * P^T;
14 
15      Note: C is assumed to be uncreated.
16            If this is not the case, Destroy C before calling this routine.
17 */
18 #undef __FUNCT__
19 #define __FUNCT__ "MatApplyPAPt_Symbolic_SeqAIJ_SeqAIJ"
20 PetscErrorCode MatApplyPAPt_Symbolic_SeqAIJ_SeqAIJ(Mat A,Mat P,Mat *C)
21 {
22   /* Note: This code is virtually identical to that of MatApplyPtAP_SeqAIJ_Symbolic */
23   /*        and MatMatMult_SeqAIJ_SeqAIJ_Symbolic.  Perhaps they could be merged nicely. */
24   PetscErrorCode     ierr;
25   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
26   Mat_SeqAIJ         *a=(Mat_SeqAIJ*)A->data,*p=(Mat_SeqAIJ*)P->data,*c;
27   PetscInt           *ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pti,*ptj,*ptjj;
28   PetscInt           *ci,*cj,*paj,*padenserow,*pasparserow,*denserow,*sparserow;
29   PetscInt           an=A->cmap->N,am=A->rmap->N,pn=P->cmap->N,pm=P->rmap->N;
30   PetscInt           i,j,k,pnzi,arow,anzj,panzi,ptrow,ptnzj,cnzi;
31   MatScalar          *ca;
32 
33   PetscFunctionBegin;
34   /* some error checking which could be moved into interface layer */
35   if (pn!=am) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %D != %D",pn,am);
36   if (am!=an) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrix 'A' must be square, %D != %D",am, an);
37 
38   /* Set up timers */
39   ierr = PetscLogEventBegin(MAT_Applypapt_symbolic,A,P,0,0);CHKERRQ(ierr);
40 
41   /* Create ij structure of P^T */
42   ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
43 
44   /* Allocate ci array, arrays for fill computation and */
45   /* free space for accumulating nonzero column info */
46   ierr = PetscMalloc(((pm+1)*1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
47   ci[0] = 0;
48 
49   ierr = PetscMalloc4(an,PetscInt,&padenserow,an,PetscInt,&pasparserow,pm,PetscInt,&denserow,pm,PetscInt,&sparserow);CHKERRQ(ierr);
50   ierr = PetscMemzero(padenserow,an*sizeof(PetscInt));CHKERRQ(ierr);
51   ierr = PetscMemzero(pasparserow,an*sizeof(PetscInt));CHKERRQ(ierr);
52   ierr = PetscMemzero(denserow,pm*sizeof(PetscInt));CHKERRQ(ierr);
53   ierr = PetscMemzero(sparserow,pm*sizeof(PetscInt));CHKERRQ(ierr);
54 
55   /* Set initial free space to be nnz(A) scaled by aspect ratio of Pt. */
56   /* This should be reasonable if sparsity of PAPt is similar to that of A. */
57   ierr          = PetscFreeSpaceGet((ai[am]/pn)*pm,&free_space);
58   current_space = free_space;
59 
60   /* Determine fill for each row of C: */
61   for (i=0;i<pm;i++) {
62     pnzi  = pi[i+1] - pi[i];
63     panzi = 0;
64     /* Get symbolic sparse row of PA: */
65     for (j=0;j<pnzi;j++) {
66       arow = *pj++;
67       anzj = ai[arow+1] - ai[arow];
68       ajj  = aj + ai[arow];
69       for (k=0;k<anzj;k++) {
70         if (!padenserow[ajj[k]]) {
71           padenserow[ajj[k]]   = -1;
72           pasparserow[panzi++] = ajj[k];
73         }
74       }
75     }
76     /* Using symbolic row of PA, determine symbolic row of C: */
77     paj    = pasparserow;
78     cnzi   = 0;
79     for (j=0;j<panzi;j++) {
80       ptrow = *paj++;
81       ptnzj = pti[ptrow+1] - pti[ptrow];
82       ptjj  = ptj + pti[ptrow];
83       for (k=0;k<ptnzj;k++) {
84         if (!denserow[ptjj[k]]) {
85           denserow[ptjj[k]] = -1;
86           sparserow[cnzi++] = ptjj[k];
87         }
88       }
89     }
90 
91     /* sort sparse representation */
92     ierr = PetscSortInt(cnzi,sparserow);CHKERRQ(ierr);
93 
94     /* If free space is not available, make more free space */
95     /* Double the amount of total space in the list */
96     if (current_space->local_remaining<cnzi) {
97       ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,&current_space);CHKERRQ(ierr);
98     }
99 
100     /* Copy data into free space, and zero out dense row */
101     ierr = PetscMemcpy(current_space->array,sparserow,cnzi*sizeof(PetscInt));CHKERRQ(ierr);
102     current_space->array           += cnzi;
103     current_space->local_used      += cnzi;
104     current_space->local_remaining -= cnzi;
105 
106     for (j=0;j<panzi;j++) {
107       padenserow[pasparserow[j]] = 0;
108     }
109     for (j=0;j<cnzi;j++) {
110       denserow[sparserow[j]] = 0;
111     }
112     ci[i+1] = ci[i] + cnzi;
113   }
114   /* column indices are in the list of free space */
115   /* Allocate space for cj, initialize cj, and */
116   /* destroy list of free space and other temporary array(s) */
117   ierr = PetscMalloc((ci[pm]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
118   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
119   ierr = PetscFree4(padenserow,pasparserow,denserow,sparserow);CHKERRQ(ierr);
120 
121   /* Allocate space for ca */
122   ierr = PetscMalloc((ci[pm]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
123   ierr = PetscMemzero(ca,(ci[pm]+1)*sizeof(MatScalar));CHKERRQ(ierr);
124 
125   /* put together the new matrix */
126   ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,pm,pm,ci,cj,ca,C);CHKERRQ(ierr);
127 
128   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
129   /* Since these are PETSc arrays, change flags to free them as necessary. */
130   c = (Mat_SeqAIJ *)((*C)->data);
131   c->free_a  = PETSC_TRUE;
132   c->free_ij = PETSC_TRUE;
133   c->nonew   = 0;
134 
135   /* Clean up. */
136   ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
137 
138   ierr = PetscLogEventEnd(MAT_Applypapt_symbolic,A,P,0,0);CHKERRQ(ierr);
139   PetscFunctionReturn(0);
140 }
141 
142 /*
143      MatApplyPAPt_Numeric_SeqAIJ - Forms the numeric product of two SeqAIJ matrices
144            C = P * A * P^T;
145      Note: C must have been created by calling MatApplyPAPt_Symbolic_SeqAIJ.
146 */
147 #undef __FUNCT__
148 #define __FUNCT__ "MatApplyPAPt_Numeric_SeqAIJ_SeqAIJ"
149 PetscErrorCode MatApplyPAPt_Numeric_SeqAIJ_SeqAIJ(Mat A,Mat P,Mat C)
150 {
151   PetscErrorCode ierr;
152   PetscInt       flops=0;
153   Mat_SeqAIJ     *a  = (Mat_SeqAIJ *) A->data;
154   Mat_SeqAIJ     *p  = (Mat_SeqAIJ *) P->data;
155   Mat_SeqAIJ     *c  = (Mat_SeqAIJ *) C->data;
156   PetscInt       *ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj=p->j,*paj,*pajdense,*ptj;
157   PetscInt       *ci=c->i,*cj=c->j;
158   PetscInt       an=A->cmap->N,am=A->rmap->N,pn=P->cmap->N,pm=P->rmap->N,cn=C->cmap->N,cm=C->rmap->N;
159   PetscInt       i,j,k,k1,k2,pnzi,anzj,panzj,arow,ptcol,ptnzj,cnzi;
160   MatScalar      *aa=a->a,*pa=p->a,*pta=p->a,*ptaj,*paa,*aaj,*ca=c->a,sum;
161 
162   PetscFunctionBegin;
163   /* This error checking should be unnecessary if the symbolic was performed */
164   if (pm!=cm) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %D != %D",pm,cm);
165   if (pn!=am) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %D != %D",pn,am);
166   if (am!=an) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrix 'A' must be square, %D != %D",am, an);
167   if (pm!=cn) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %D != %D",pm, cn);
168 
169   /* Set up timers */
170   ierr = PetscLogEventBegin(MAT_Applypapt_numeric,A,P,C,0);CHKERRQ(ierr);
171   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
172 
173   ierr = PetscMalloc3(an,MatScalar,&paa,an,PetscInt,&paj,an,PetscInt,&pajdense);CHKERRQ(ierr);
174   ierr = PetscMemzero(paa,an*(sizeof(MatScalar)+2*sizeof(PetscInt)));CHKERRQ(ierr);
175 
176   for (i=0;i<pm;i++) {
177     /* Form sparse row of P*A */
178     pnzi  = pi[i+1] - pi[i];
179     panzj = 0;
180     for (j=0;j<pnzi;j++) {
181       arow = *pj++;
182       anzj = ai[arow+1] - ai[arow];
183       ajj  = aj + ai[arow];
184       aaj  = aa + ai[arow];
185       for (k=0;k<anzj;k++) {
186         if (!pajdense[ajj[k]]) {
187           pajdense[ajj[k]] = -1;
188           paj[panzj++]     = ajj[k];
189         }
190         paa[ajj[k]] += (*pa)*aaj[k];
191       }
192       flops += 2*anzj;
193       pa++;
194     }
195 
196     /* Sort the j index array for quick sparse axpy. */
197     ierr = PetscSortInt(panzj,paj);CHKERRQ(ierr);
198 
199     /* Compute P*A*P^T using sparse inner products. */
200     /* Take advantage of pre-computed (i,j) of C for locations of non-zeros. */
201     cnzi = ci[i+1] - ci[i];
202     for (j=0;j<cnzi;j++) {
203       /* Form sparse inner product of current row of P*A with (*cj++) col of P^T. */
204       ptcol = *cj++;
205       ptnzj = pi[ptcol+1] - pi[ptcol];
206       ptj   = pjj + pi[ptcol];
207       ptaj  = pta + pi[ptcol];
208       sum   = 0.;
209       k1    = 0;
210       k2    = 0;
211       while ((k1<panzj) && (k2<ptnzj)) {
212         if (paj[k1]==ptj[k2]) {
213           sum += paa[paj[k1++]]*ptaj[k2++];
214         } else if (paj[k1] < ptj[k2]) {
215           k1++;
216         } else /* if (paj[k1] > ptj[k2]) */ {
217           k2++;
218         }
219       }
220       *ca++ = sum;
221     }
222 
223     /* Zero the current row info for P*A */
224     for (j=0;j<panzj;j++) {
225       paa[paj[j]]      = 0.;
226       pajdense[paj[j]] = 0;
227     }
228   }
229 
230   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
231   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
232   ierr = PetscFree3(paa,paj,pajdense);CHKERRQ(ierr);
233   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
234   ierr = PetscLogEventEnd(MAT_Applypapt_numeric,A,P,C,0);CHKERRQ(ierr);
235   PetscFunctionReturn(0);
236 }
237 
238 #undef __FUNCT__
239 #define __FUNCT__ "MatApplyPAPt_SeqAIJ_SeqAIJ"
240 PetscErrorCode MatApplyPAPt_SeqAIJ_SeqAIJ(Mat A,Mat P,Mat *C)
241 {
242   PetscErrorCode ierr;
243 
244   PetscFunctionBegin;
245   ierr = PetscLogEventBegin(MAT_Applypapt,A,P,0,0);CHKERRQ(ierr);
246   ierr = MatApplyPAPt_Symbolic_SeqAIJ_SeqAIJ(A,P,C);CHKERRQ(ierr);
247   ierr = MatApplyPAPt_Numeric_SeqAIJ_SeqAIJ(A,P,*C);CHKERRQ(ierr);
248   ierr = PetscLogEventEnd(MAT_Applypapt,A,P,0,0);CHKERRQ(ierr);
249   PetscFunctionReturn(0);
250 }
251 
252 /*--------------------------------------------------*/
253 /*
254   Defines projective product routines where A is a SeqAIJ matrix
255           C = R * A * R^T
256 */
257 
258 #undef __FUNCT__
259 #define __FUNCT__ "PetscContainerDestroy_Mat_RARt"
260 PetscErrorCode PetscContainerDestroy_Mat_RARt(void *ptr)
261 {
262   PetscErrorCode ierr;
263   Mat_RARt       *rart=(Mat_RARt*)ptr;
264 
265   PetscFunctionBegin;
266   ierr = MatTransposeColoringDestroy(&rart->matcoloring);CHKERRQ(ierr);
267   ierr = MatDestroy(&rart->Rt);CHKERRQ(ierr);
268   ierr = MatDestroy(&rart->RARt);CHKERRQ(ierr);
269   ierr = PetscFree(rart->work);CHKERRQ(ierr);
270   ierr = PetscFree(rart);CHKERRQ(ierr);
271   PetscFunctionReturn(0);
272 }
273 
274 #undef __FUNCT__
275 #define __FUNCT__ "MatDestroy_SeqAIJ_RARt"
276 PetscErrorCode MatDestroy_SeqAIJ_RARt(Mat A)
277 {
278   PetscErrorCode ierr;
279   PetscContainer container;
280   Mat_RARt       *rart=PETSC_NULL;
281 
282   PetscFunctionBegin;
283   ierr = PetscObjectQuery((PetscObject)A,"Mat_RARt",(PetscObject *)&container);CHKERRQ(ierr);
284   if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit");
285   ierr = PetscContainerGetPointer(container,(void **)&rart);CHKERRQ(ierr);
286   A->ops->destroy   = rart->destroy;
287   if (A->ops->destroy) {
288     ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
289   }
290   ierr = PetscObjectCompose((PetscObject)A,"Mat_RARt",0);CHKERRQ(ierr);
291   PetscFunctionReturn(0);
292 }
293 
294 #undef __FUNCT__
295 #define __FUNCT__ "MatRARtSymbolic_SeqAIJ_SeqAIJ"
296 PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat R,PetscReal fill,Mat *C)
297 {
298   PetscErrorCode      ierr;
299   Mat                 P;
300   PetscInt            *rti,*rtj;
301   Mat_RARt            *rart;
302   PetscContainer      container;
303   MatTransposeColoring matcoloring;
304   ISColoring           iscoloring;
305   Mat                  Rt_dense,RARt_dense;
306   PetscLogDouble       GColor=0.0,MCCreate=0.0,MDenCreate=0.0,t0,tf,etime=0.0;
307   Mat_SeqAIJ           *c;
308 
309   PetscFunctionBegin;
310   ierr = PetscGetTime(&t0);CHKERRQ(ierr);
311   /* create symbolic P=Rt */
312   ierr = MatGetSymbolicTranspose_SeqAIJ(R,&rti,&rtj);CHKERRQ(ierr);
313   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,R->cmap->n,R->rmap->n,rti,rtj,PETSC_NULL,&P);CHKERRQ(ierr);
314 
315   /* get symbolic C=Pt*A*P */
316   ierr = MatPtAPSymbolic_SeqAIJ_SeqAIJ(A,P,fill,C);CHKERRQ(ierr);
317 
318   /* create a supporting struct */
319   ierr = PetscNew(Mat_RARt,&rart);CHKERRQ(ierr);
320 
321   /* attach the supporting struct to C */
322   ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr);
323   ierr = PetscContainerSetPointer(container,rart);CHKERRQ(ierr);
324   ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_Mat_RARt);CHKERRQ(ierr);
325   ierr = PetscObjectCompose((PetscObject)(*C),"Mat_RARt",(PetscObject)container);CHKERRQ(ierr);
326   ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
327   ierr = PetscGetTime(&tf);CHKERRQ(ierr);
328   etime += tf - t0;
329 
330   /* Create MatTransposeColoring from symbolic C=R*A*R^T */
331   c=(Mat_SeqAIJ*)(*C)->data;
332   ierr = PetscGetTime(&t0);CHKERRQ(ierr);
333   ierr = MatGetColoring(*C,MATCOLORINGLF,&iscoloring);CHKERRQ(ierr);
334   ierr = PetscGetTime(&tf);CHKERRQ(ierr);
335   GColor += tf - t0;
336 
337   ierr = PetscGetTime(&t0);CHKERRQ(ierr);
338   ierr = MatTransposeColoringCreate(*C,iscoloring,&matcoloring);CHKERRQ(ierr);
339   rart->matcoloring = matcoloring;
340   ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr);
341   ierr = PetscGetTime(&tf);CHKERRQ(ierr);
342   MCCreate += tf - t0;
343 
344   ierr = PetscGetTime(&t0);CHKERRQ(ierr);
345   /* Create Rt_dense */
346   ierr = MatCreate(PETSC_COMM_SELF,&Rt_dense);CHKERRQ(ierr);
347   ierr = MatSetSizes(Rt_dense,A->cmap->n,matcoloring->ncolors,A->cmap->n,matcoloring->ncolors);CHKERRQ(ierr);
348   ierr = MatSetType(Rt_dense,MATSEQDENSE);CHKERRQ(ierr);
349   ierr = MatSeqDenseSetPreallocation(Rt_dense,PETSC_NULL);CHKERRQ(ierr);
350   Rt_dense->assembled = PETSC_TRUE;
351   rart->Rt            = Rt_dense;
352 
353   /* Create RARt_dense = R*A*Rt_dense */
354   ierr = MatCreate(PETSC_COMM_SELF,&RARt_dense);CHKERRQ(ierr);
355   ierr = MatSetSizes(RARt_dense,(*C)->rmap->n,matcoloring->ncolors,(*C)->rmap->n,matcoloring->ncolors);CHKERRQ(ierr);
356   ierr = MatSetType(RARt_dense,MATSEQDENSE);CHKERRQ(ierr);
357   ierr = MatSeqDenseSetPreallocation(RARt_dense,PETSC_NULL);CHKERRQ(ierr);
358   rart->RARt = RARt_dense;
359 
360   /* Allocate work array to store columns of A*R^T used in MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense() */
361   ierr = PetscMalloc(A->rmap->n*4*sizeof(PetscScalar),&rart->work);CHKERRQ(ierr);
362 
363   ierr = PetscGetTime(&tf);CHKERRQ(ierr);
364   MDenCreate += tf - t0;
365 
366   rart->destroy = (*C)->ops->destroy;
367   (*C)->ops->destroy = MatDestroy_SeqAIJ_RARt;
368 
369   /* clean up */
370   ierr = MatRestoreSymbolicTranspose_SeqAIJ(R,&rti,&rtj);CHKERRQ(ierr);
371   ierr = MatDestroy(&P);CHKERRQ(ierr);
372 
373 #if defined(PETSC_USE_INFO)
374   {
375   PetscReal density= (PetscReal)(c->nz)/(RARt_dense->rmap->n*RARt_dense->cmap->n);
376   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);
377   ierr = PetscInfo5(*C,"Sym = GetColor %g + MColorCreate %g + MDenCreate %g + other %g = %g\n",GColor,MCCreate,MDenCreate,etime,GColor+MCCreate+MDenCreate+etime);
378   }
379 #endif
380   PetscFunctionReturn(0);
381 }
382 
383 /*
384  RAB = R * A * B, R and A in seqaij format, B in dense format;
385 */
386 #undef __FUNCT__
387 #define __FUNCT__ "MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense"
388 PetscErrorCode MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense(Mat R,Mat A,Mat B,Mat RAB,PetscScalar *work)
389 {
390   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*r=(Mat_SeqAIJ*)R->data;
391   PetscErrorCode ierr;
392   PetscScalar    *b,r1,r2,r3,r4,*b1,*b2,*b3,*b4;
393   MatScalar      *aa,*ra;
394   PetscInt       cn=B->cmap->n,bm=B->rmap->n,col,i,j,n,*ai=a->i,*aj,am=A->rmap->n;
395   PetscInt       am2=2*am,am3=3*am,bm4=4*bm;
396   PetscScalar    *d,*c,*c2,*c3,*c4;
397   PetscInt       *rj,rm=R->rmap->n,dm=RAB->rmap->n,dn=RAB->cmap->n;
398   PetscInt       rm2=2*rm,rm3=3*rm,colrm;
399 
400   PetscFunctionBegin;
401   if (!dm || !dn) PetscFunctionReturn(0);
402   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);
403   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);
404   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);
405   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);
406 
407   ierr = MatGetArray(B,&b);CHKERRQ(ierr);
408   ierr = MatGetArray(RAB,&d);CHKERRQ(ierr);
409   b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;
410   c = work; c2 = c + am; c3 = c2 + am; c4 = c3 + am;
411   for (col=0; col<cn-4; col += 4){  /* over columns of C */
412     for (i=0; i<am; i++) {        /* over rows of A in those columns */
413       r1 = r2 = r3 = r4 = 0.0;
414       n   = ai[i+1] - ai[i];
415       aj  = a->j + ai[i];
416       aa  = a->a + ai[i];
417       for (j=0; j<n; j++) {
418         r1 += (*aa)*b1[*aj];
419         r2 += (*aa)*b2[*aj];
420         r3 += (*aa)*b3[*aj];
421         r4 += (*aa++)*b4[*aj++];
422       }
423       c[i]       = r1;
424       c[am  + i] = r2;
425       c[am2 + i] = r3;
426       c[am3 + i] = r4;
427     }
428     b1 += bm4;
429     b2 += bm4;
430     b3 += bm4;
431     b4 += bm4;
432 
433     /* RAB[:,col] = R*C[:,col] */
434     colrm = col*rm;
435     for (i=0; i<rm; i++) {        /* over rows of R in those columns */
436       r1 = r2 = r3 = r4 = 0.0;
437       n   = r->i[i+1] - r->i[i];
438       rj  = r->j + r->i[i];
439       ra  = r->a + r->i[i];
440       for (j=0; j<n; j++) {
441         r1 += (*ra)*c[*rj];
442         r2 += (*ra)*c2[*rj];
443         r3 += (*ra)*c3[*rj];
444         r4 += (*ra++)*c4[*rj++];
445       }
446       d[colrm + i]       = r1;
447       d[colrm + rm + i]  = r2;
448       d[colrm + rm2 + i] = r3;
449       d[colrm + rm3 + i] = r4;
450     }
451   }
452   for (;col<cn; col++){     /* over extra columns of C */
453     for (i=0; i<am; i++) {  /* over rows of A in those columns */
454       r1 = 0.0;
455       n   = a->i[i+1] - a->i[i];
456       aj  = a->j + a->i[i];
457       aa  = a->a + a->i[i];
458       for (j=0; j<n; j++) {
459         r1 += (*aa++)*b1[*aj++];
460       }
461       c[i]     = r1;
462     }
463     b1 += bm;
464 
465     for (i=0; i<rm; i++) {  /* over rows of R in those columns */
466       r1 = 0.0;
467       n   = r->i[i+1] - r->i[i];
468       rj  = r->j + r->i[i];
469       ra  = r->a + r->i[i];
470       for (j=0; j<n; j++) {
471         r1 += (*ra++)*c[*rj++];
472       }
473       d[col*rm + i]     = r1;
474     }
475   }
476   ierr = PetscLogFlops(cn*2.0*(a->nz + r->nz));CHKERRQ(ierr);
477 
478   ierr = MatRestoreArray(B,&b);CHKERRQ(ierr);
479   ierr = MatRestoreArray(RAB,&d);CHKERRQ(ierr);
480   ierr = MatAssemblyBegin(RAB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
481   ierr = MatAssemblyEnd(RAB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
482   PetscFunctionReturn(0);
483 }
484 
485 #undef __FUNCT__
486 #define __FUNCT__ "MatRARtNumeric_SeqAIJ_SeqAIJ"
487 PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ(Mat A,Mat R,Mat C)
488 {
489   PetscErrorCode        ierr;
490   Mat_RARt              *rart;
491   PetscContainer        container;
492   MatTransposeColoring  matcoloring;
493   Mat                   Rt,RARt;
494   PetscLogDouble        Mult_sp_den=0.0,app1=0.0,app2=0.0,t0,tf;
495 
496   PetscFunctionBegin;
497   ierr = PetscObjectQuery((PetscObject)C,"Mat_RARt",(PetscObject *)&container);CHKERRQ(ierr);
498   if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit");
499   ierr  = PetscContainerGetPointer(container,(void **)&rart);CHKERRQ(ierr);
500 
501   /* Get dense Rt by Apply MatTransposeColoring to R */
502   matcoloring = rart->matcoloring;
503   Rt          = rart->Rt;
504   ierr = PetscGetTime(&t0);CHKERRQ(ierr);
505   ierr = MatTransColoringApplySpToDen(matcoloring,R,Rt);CHKERRQ(ierr);
506   ierr = PetscGetTime(&tf);CHKERRQ(ierr);
507   app1 += tf - t0;
508 
509   /* Get dense RARt = R*A*Rt */
510   ierr = PetscGetTime(&t0);CHKERRQ(ierr);
511   RARt = rart->RARt;
512   ierr = MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense(R,A,Rt,RARt,rart->work);CHKERRQ(ierr);
513   ierr = PetscGetTime(&tf);CHKERRQ(ierr);
514   Mult_sp_den += tf - t0;
515 
516   /* Recover C from C_dense */
517   ierr = PetscGetTime(&t0);CHKERRQ(ierr);
518   ierr = MatTransColoringApplyDenToSp(matcoloring,RARt,C);CHKERRQ(ierr);
519   ierr = PetscGetTime(&tf);CHKERRQ(ierr);
520   app2 += tf - t0;
521 
522 #if defined(PETSC_USE_INFO)
523   ierr = PetscInfo4(C,"Num = ColorApp %g + %g + Mult_sp_den %g  = %g\n",app1,app2,Mult_sp_den,app1+app2+Mult_sp_den);
524 #endif
525   PetscFunctionReturn(0);
526 }
527 
528 EXTERN_C_BEGIN
529 #undef __FUNCT__
530 #define __FUNCT__ "MatRARt_SeqAIJ_SeqAIJ"
531 PetscErrorCode MatRARt_SeqAIJ_SeqAIJ(Mat A,Mat R,MatReuse scall,PetscReal fill,Mat *C)
532 {
533   PetscErrorCode ierr;
534 
535   PetscFunctionBegin;
536   if (scall == MAT_INITIAL_MATRIX){
537     ierr = MatRARtSymbolic_SeqAIJ_SeqAIJ(A,R,fill,C);CHKERRQ(ierr);
538   }
539   ierr = MatRARtNumeric_SeqAIJ_SeqAIJ(A,R,*C);CHKERRQ(ierr);
540   PetscFunctionReturn(0);
541 }
542 EXTERN_C_END
543