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