xref: /petsc/src/mat/impls/aij/seq/matrart.c (revision 1ce71dff06bb5e1e97c3afb7c51fa9ca0d212a27)
1807171c5SHong Zhang 
2807171c5SHong Zhang /*
3807171c5SHong Zhang   Defines matrix-matrix product routines for pairs of SeqAIJ matrices
4807171c5SHong Zhang           C = P * A * P^T
5807171c5SHong Zhang */
6807171c5SHong 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 
11807171c5SHong Zhang /*
12807171c5SHong Zhang      MatApplyPAPt_Symbolic_SeqAIJ_SeqAIJ - Forms the symbolic product of two SeqAIJ matrices
13807171c5SHong Zhang            C = P * A * P^T;
14807171c5SHong Zhang 
15807171c5SHong Zhang      Note: C is assumed to be uncreated.
16807171c5SHong Zhang            If this is not the case, Destroy C before calling this routine.
17807171c5SHong Zhang */
18807171c5SHong Zhang #undef __FUNCT__
19807171c5SHong Zhang #define __FUNCT__ "MatApplyPAPt_Symbolic_SeqAIJ_SeqAIJ"
20807171c5SHong Zhang PetscErrorCode MatApplyPAPt_Symbolic_SeqAIJ_SeqAIJ(Mat A,Mat P,Mat *C)
21807171c5SHong Zhang {
22807171c5SHong Zhang   /* Note: This code is virtually identical to that of MatApplyPtAP_SeqAIJ_Symbolic */
23807171c5SHong Zhang   /*        and MatMatMult_SeqAIJ_SeqAIJ_Symbolic.  Perhaps they could be merged nicely. */
24807171c5SHong Zhang   PetscErrorCode     ierr;
25807171c5SHong Zhang   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
26807171c5SHong Zhang   Mat_SeqAIJ         *a=(Mat_SeqAIJ*)A->data,*p=(Mat_SeqAIJ*)P->data,*c;
27807171c5SHong Zhang   PetscInt           *ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pti,*ptj,*ptjj;
28807171c5SHong Zhang   PetscInt           *ci,*cj,*paj,*padenserow,*pasparserow,*denserow,*sparserow;
29807171c5SHong Zhang   PetscInt           an=A->cmap->N,am=A->rmap->N,pn=P->cmap->N,pm=P->rmap->N;
30807171c5SHong Zhang   PetscInt           i,j,k,pnzi,arow,anzj,panzi,ptrow,ptnzj,cnzi;
31807171c5SHong Zhang   MatScalar          *ca;
32807171c5SHong Zhang 
33807171c5SHong Zhang   PetscFunctionBegin;
34807171c5SHong Zhang   /* some error checking which could be moved into interface layer */
35807171c5SHong Zhang   if (pn!=am) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %D != %D",pn,am);
36807171c5SHong Zhang   if (am!=an) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrix 'A' must be square, %D != %D",am, an);
37807171c5SHong Zhang 
38807171c5SHong Zhang   /* Set up timers */
39807171c5SHong Zhang   ierr = PetscLogEventBegin(MAT_Applypapt_symbolic,A,P,0,0);CHKERRQ(ierr);
40807171c5SHong Zhang 
41807171c5SHong Zhang   /* Create ij structure of P^T */
42807171c5SHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
43807171c5SHong Zhang 
44807171c5SHong Zhang   /* Allocate ci array, arrays for fill computation and */
45807171c5SHong Zhang   /* free space for accumulating nonzero column info */
46807171c5SHong Zhang   ierr = PetscMalloc(((pm+1)*1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
47807171c5SHong Zhang   ci[0] = 0;
48807171c5SHong Zhang 
49807171c5SHong Zhang   ierr = PetscMalloc4(an,PetscInt,&padenserow,an,PetscInt,&pasparserow,pm,PetscInt,&denserow,pm,PetscInt,&sparserow);CHKERRQ(ierr);
50807171c5SHong Zhang   ierr = PetscMemzero(padenserow,an*sizeof(PetscInt));CHKERRQ(ierr);
51807171c5SHong Zhang   ierr = PetscMemzero(pasparserow,an*sizeof(PetscInt));CHKERRQ(ierr);
52807171c5SHong Zhang   ierr = PetscMemzero(denserow,pm*sizeof(PetscInt));CHKERRQ(ierr);
53807171c5SHong Zhang   ierr = PetscMemzero(sparserow,pm*sizeof(PetscInt));CHKERRQ(ierr);
54807171c5SHong Zhang 
55807171c5SHong Zhang   /* Set initial free space to be nnz(A) scaled by aspect ratio of Pt. */
56807171c5SHong Zhang   /* This should be reasonable if sparsity of PAPt is similar to that of A. */
57807171c5SHong Zhang   ierr          = PetscFreeSpaceGet((ai[am]/pn)*pm,&free_space);
58807171c5SHong Zhang   current_space = free_space;
59807171c5SHong Zhang 
60807171c5SHong Zhang   /* Determine fill for each row of C: */
61807171c5SHong Zhang   for (i=0;i<pm;i++) {
62807171c5SHong Zhang     pnzi  = pi[i+1] - pi[i];
63807171c5SHong Zhang     panzi = 0;
64807171c5SHong Zhang     /* Get symbolic sparse row of PA: */
65807171c5SHong Zhang     for (j=0;j<pnzi;j++) {
66807171c5SHong Zhang       arow = *pj++;
67807171c5SHong Zhang       anzj = ai[arow+1] - ai[arow];
68807171c5SHong Zhang       ajj  = aj + ai[arow];
69807171c5SHong Zhang       for (k=0;k<anzj;k++) {
70807171c5SHong Zhang         if (!padenserow[ajj[k]]) {
71807171c5SHong Zhang           padenserow[ajj[k]]   = -1;
72807171c5SHong Zhang           pasparserow[panzi++] = ajj[k];
73807171c5SHong Zhang         }
74807171c5SHong Zhang       }
75807171c5SHong Zhang     }
76807171c5SHong Zhang     /* Using symbolic row of PA, determine symbolic row of C: */
77807171c5SHong Zhang     paj    = pasparserow;
78807171c5SHong Zhang     cnzi   = 0;
79807171c5SHong Zhang     for (j=0;j<panzi;j++) {
80807171c5SHong Zhang       ptrow = *paj++;
81807171c5SHong Zhang       ptnzj = pti[ptrow+1] - pti[ptrow];
82807171c5SHong Zhang       ptjj  = ptj + pti[ptrow];
83807171c5SHong Zhang       for (k=0;k<ptnzj;k++) {
84807171c5SHong Zhang         if (!denserow[ptjj[k]]) {
85807171c5SHong Zhang           denserow[ptjj[k]] = -1;
86807171c5SHong Zhang           sparserow[cnzi++] = ptjj[k];
87807171c5SHong Zhang         }
88807171c5SHong Zhang       }
89807171c5SHong Zhang     }
90807171c5SHong Zhang 
91807171c5SHong Zhang     /* sort sparse representation */
92807171c5SHong Zhang     ierr = PetscSortInt(cnzi,sparserow);CHKERRQ(ierr);
93807171c5SHong Zhang 
94807171c5SHong Zhang     /* If free space is not available, make more free space */
95807171c5SHong Zhang     /* Double the amount of total space in the list */
96807171c5SHong Zhang     if (current_space->local_remaining<cnzi) {
97807171c5SHong Zhang       ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,&current_space);CHKERRQ(ierr);
98807171c5SHong Zhang     }
99807171c5SHong Zhang 
100807171c5SHong Zhang     /* Copy data into free space, and zero out dense row */
101807171c5SHong Zhang     ierr = PetscMemcpy(current_space->array,sparserow,cnzi*sizeof(PetscInt));CHKERRQ(ierr);
102807171c5SHong Zhang     current_space->array           += cnzi;
103807171c5SHong Zhang     current_space->local_used      += cnzi;
104807171c5SHong Zhang     current_space->local_remaining -= cnzi;
105807171c5SHong Zhang 
106807171c5SHong Zhang     for (j=0;j<panzi;j++) {
107807171c5SHong Zhang       padenserow[pasparserow[j]] = 0;
108807171c5SHong Zhang     }
109807171c5SHong Zhang     for (j=0;j<cnzi;j++) {
110807171c5SHong Zhang       denserow[sparserow[j]] = 0;
111807171c5SHong Zhang     }
112807171c5SHong Zhang     ci[i+1] = ci[i] + cnzi;
113807171c5SHong Zhang   }
114807171c5SHong Zhang   /* column indices are in the list of free space */
115807171c5SHong Zhang   /* Allocate space for cj, initialize cj, and */
116807171c5SHong Zhang   /* destroy list of free space and other temporary array(s) */
117807171c5SHong Zhang   ierr = PetscMalloc((ci[pm]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
118807171c5SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
119807171c5SHong Zhang   ierr = PetscFree4(padenserow,pasparserow,denserow,sparserow);CHKERRQ(ierr);
120807171c5SHong Zhang 
121807171c5SHong Zhang   /* Allocate space for ca */
122807171c5SHong Zhang   ierr = PetscMalloc((ci[pm]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
123807171c5SHong Zhang   ierr = PetscMemzero(ca,(ci[pm]+1)*sizeof(MatScalar));CHKERRQ(ierr);
124807171c5SHong Zhang 
125807171c5SHong Zhang   /* put together the new matrix */
126807171c5SHong Zhang   ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,pm,pm,ci,cj,ca,C);CHKERRQ(ierr);
127807171c5SHong Zhang 
128807171c5SHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
129807171c5SHong Zhang   /* Since these are PETSc arrays, change flags to free them as necessary. */
130807171c5SHong Zhang   c = (Mat_SeqAIJ *)((*C)->data);
131807171c5SHong Zhang   c->free_a  = PETSC_TRUE;
132807171c5SHong Zhang   c->free_ij = PETSC_TRUE;
133807171c5SHong Zhang   c->nonew   = 0;
134807171c5SHong Zhang 
135807171c5SHong Zhang   /* Clean up. */
136807171c5SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
137807171c5SHong Zhang 
138807171c5SHong Zhang   ierr = PetscLogEventEnd(MAT_Applypapt_symbolic,A,P,0,0);CHKERRQ(ierr);
139807171c5SHong Zhang   PetscFunctionReturn(0);
140807171c5SHong Zhang }
141807171c5SHong Zhang 
142807171c5SHong Zhang /*
143807171c5SHong Zhang      MatApplyPAPt_Numeric_SeqAIJ - Forms the numeric product of two SeqAIJ matrices
144807171c5SHong Zhang            C = P * A * P^T;
145807171c5SHong Zhang      Note: C must have been created by calling MatApplyPAPt_Symbolic_SeqAIJ.
146807171c5SHong Zhang */
147807171c5SHong Zhang #undef __FUNCT__
148807171c5SHong Zhang #define __FUNCT__ "MatApplyPAPt_Numeric_SeqAIJ_SeqAIJ"
149807171c5SHong Zhang PetscErrorCode MatApplyPAPt_Numeric_SeqAIJ_SeqAIJ(Mat A,Mat P,Mat C)
150807171c5SHong Zhang {
151807171c5SHong Zhang   PetscErrorCode ierr;
152807171c5SHong Zhang   PetscInt       flops=0;
153807171c5SHong Zhang   Mat_SeqAIJ     *a  = (Mat_SeqAIJ *) A->data;
154807171c5SHong Zhang   Mat_SeqAIJ     *p  = (Mat_SeqAIJ *) P->data;
155807171c5SHong Zhang   Mat_SeqAIJ     *c  = (Mat_SeqAIJ *) C->data;
156807171c5SHong Zhang   PetscInt       *ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj=p->j,*paj,*pajdense,*ptj;
157807171c5SHong Zhang   PetscInt       *ci=c->i,*cj=c->j;
158807171c5SHong Zhang   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;
159807171c5SHong Zhang   PetscInt       i,j,k,k1,k2,pnzi,anzj,panzj,arow,ptcol,ptnzj,cnzi;
160807171c5SHong Zhang   MatScalar      *aa=a->a,*pa=p->a,*pta=p->a,*ptaj,*paa,*aaj,*ca=c->a,sum;
161807171c5SHong Zhang 
162807171c5SHong Zhang   PetscFunctionBegin;
163807171c5SHong Zhang   /* This error checking should be unnecessary if the symbolic was performed */
164807171c5SHong Zhang   if (pm!=cm) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %D != %D",pm,cm);
165807171c5SHong Zhang   if (pn!=am) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %D != %D",pn,am);
166807171c5SHong Zhang   if (am!=an) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrix 'A' must be square, %D != %D",am, an);
167807171c5SHong Zhang   if (pm!=cn) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %D != %D",pm, cn);
168807171c5SHong Zhang 
169807171c5SHong Zhang   /* Set up timers */
170807171c5SHong Zhang   ierr = PetscLogEventBegin(MAT_Applypapt_numeric,A,P,C,0);CHKERRQ(ierr);
171807171c5SHong Zhang   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
172807171c5SHong Zhang 
173807171c5SHong Zhang   ierr = PetscMalloc3(an,MatScalar,&paa,an,PetscInt,&paj,an,PetscInt,&pajdense);CHKERRQ(ierr);
174807171c5SHong Zhang   ierr = PetscMemzero(paa,an*(sizeof(MatScalar)+2*sizeof(PetscInt)));CHKERRQ(ierr);
175807171c5SHong Zhang 
176807171c5SHong Zhang   for (i=0;i<pm;i++) {
177807171c5SHong Zhang     /* Form sparse row of P*A */
178807171c5SHong Zhang     pnzi  = pi[i+1] - pi[i];
179807171c5SHong Zhang     panzj = 0;
180807171c5SHong Zhang     for (j=0;j<pnzi;j++) {
181807171c5SHong Zhang       arow = *pj++;
182807171c5SHong Zhang       anzj = ai[arow+1] - ai[arow];
183807171c5SHong Zhang       ajj  = aj + ai[arow];
184807171c5SHong Zhang       aaj  = aa + ai[arow];
185807171c5SHong Zhang       for (k=0;k<anzj;k++) {
186807171c5SHong Zhang         if (!pajdense[ajj[k]]) {
187807171c5SHong Zhang           pajdense[ajj[k]] = -1;
188807171c5SHong Zhang           paj[panzj++]     = ajj[k];
189807171c5SHong Zhang         }
190807171c5SHong Zhang         paa[ajj[k]] += (*pa)*aaj[k];
191807171c5SHong Zhang       }
192807171c5SHong Zhang       flops += 2*anzj;
193807171c5SHong Zhang       pa++;
194807171c5SHong Zhang     }
195807171c5SHong Zhang 
196807171c5SHong Zhang     /* Sort the j index array for quick sparse axpy. */
197807171c5SHong Zhang     ierr = PetscSortInt(panzj,paj);CHKERRQ(ierr);
198807171c5SHong Zhang 
199807171c5SHong Zhang     /* Compute P*A*P^T using sparse inner products. */
200807171c5SHong Zhang     /* Take advantage of pre-computed (i,j) of C for locations of non-zeros. */
201807171c5SHong Zhang     cnzi = ci[i+1] - ci[i];
202807171c5SHong Zhang     for (j=0;j<cnzi;j++) {
203807171c5SHong Zhang       /* Form sparse inner product of current row of P*A with (*cj++) col of P^T. */
204807171c5SHong Zhang       ptcol = *cj++;
205807171c5SHong Zhang       ptnzj = pi[ptcol+1] - pi[ptcol];
206807171c5SHong Zhang       ptj   = pjj + pi[ptcol];
207807171c5SHong Zhang       ptaj  = pta + pi[ptcol];
208807171c5SHong Zhang       sum   = 0.;
209807171c5SHong Zhang       k1    = 0;
210807171c5SHong Zhang       k2    = 0;
211807171c5SHong Zhang       while ((k1<panzj) && (k2<ptnzj)) {
212807171c5SHong Zhang         if (paj[k1]==ptj[k2]) {
213807171c5SHong Zhang           sum += paa[paj[k1++]]*ptaj[k2++];
214807171c5SHong Zhang         } else if (paj[k1] < ptj[k2]) {
215807171c5SHong Zhang           k1++;
216807171c5SHong Zhang         } else /* if (paj[k1] > ptj[k2]) */ {
217807171c5SHong Zhang           k2++;
218807171c5SHong Zhang         }
219807171c5SHong Zhang       }
220807171c5SHong Zhang       *ca++ = sum;
221807171c5SHong Zhang     }
222807171c5SHong Zhang 
223807171c5SHong Zhang     /* Zero the current row info for P*A */
224807171c5SHong Zhang     for (j=0;j<panzj;j++) {
225807171c5SHong Zhang       paa[paj[j]]      = 0.;
226807171c5SHong Zhang       pajdense[paj[j]] = 0;
227807171c5SHong Zhang     }
228807171c5SHong Zhang   }
229807171c5SHong Zhang 
230807171c5SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
231807171c5SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
232807171c5SHong Zhang   ierr = PetscFree3(paa,paj,pajdense);CHKERRQ(ierr);
233807171c5SHong Zhang   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
234807171c5SHong Zhang   ierr = PetscLogEventEnd(MAT_Applypapt_numeric,A,P,C,0);CHKERRQ(ierr);
235807171c5SHong Zhang   PetscFunctionReturn(0);
236807171c5SHong Zhang }
237807171c5SHong Zhang 
238807171c5SHong Zhang #undef __FUNCT__
239807171c5SHong Zhang #define __FUNCT__ "MatApplyPAPt_SeqAIJ_SeqAIJ"
240807171c5SHong Zhang PetscErrorCode MatApplyPAPt_SeqAIJ_SeqAIJ(Mat A,Mat P,Mat *C)
241807171c5SHong Zhang {
242807171c5SHong Zhang   PetscErrorCode ierr;
243807171c5SHong Zhang 
244807171c5SHong Zhang   PetscFunctionBegin;
245807171c5SHong Zhang   ierr = PetscLogEventBegin(MAT_Applypapt,A,P,0,0);CHKERRQ(ierr);
246807171c5SHong Zhang   ierr = MatApplyPAPt_Symbolic_SeqAIJ_SeqAIJ(A,P,C);CHKERRQ(ierr);
247807171c5SHong Zhang   ierr = MatApplyPAPt_Numeric_SeqAIJ_SeqAIJ(A,P,*C);CHKERRQ(ierr);
248807171c5SHong Zhang   ierr = PetscLogEventEnd(MAT_Applypapt,A,P,0,0);CHKERRQ(ierr);
249807171c5SHong Zhang   PetscFunctionReturn(0);
250807171c5SHong Zhang }
251807171c5SHong Zhang 
252807171c5SHong Zhang /*--------------------------------------------------*/
253807171c5SHong Zhang /*
254807171c5SHong Zhang   Defines projective product routines where A is a SeqAIJ matrix
255807171c5SHong Zhang           C = R * A * R^T
256807171c5SHong Zhang */
257807171c5SHong Zhang 
258807171c5SHong Zhang #undef __FUNCT__
259807171c5SHong Zhang #define __FUNCT__ "PetscContainerDestroy_Mat_RARt"
260807171c5SHong Zhang PetscErrorCode PetscContainerDestroy_Mat_RARt(void *ptr)
261807171c5SHong Zhang {
262807171c5SHong Zhang   PetscErrorCode ierr;
263807171c5SHong Zhang   Mat_RARt       *rart=(Mat_RARt*)ptr;
264807171c5SHong Zhang 
265807171c5SHong Zhang   PetscFunctionBegin;
266807171c5SHong Zhang   ierr = MatTransposeColoringDestroy(&rart->matcoloring);CHKERRQ(ierr);
267807171c5SHong Zhang   ierr = MatDestroy(&rart->Rt);CHKERRQ(ierr);
268807171c5SHong Zhang   ierr = MatDestroy(&rart->RARt);CHKERRQ(ierr);
269807171c5SHong Zhang   ierr = PetscFree(rart->work);CHKERRQ(ierr);
270807171c5SHong Zhang   ierr = PetscFree(rart);CHKERRQ(ierr);
271807171c5SHong Zhang   PetscFunctionReturn(0);
272807171c5SHong Zhang }
273807171c5SHong Zhang 
274807171c5SHong Zhang #undef __FUNCT__
275807171c5SHong Zhang #define __FUNCT__ "MatDestroy_SeqAIJ_RARt"
276807171c5SHong Zhang PetscErrorCode MatDestroy_SeqAIJ_RARt(Mat A)
277807171c5SHong Zhang {
278807171c5SHong Zhang   PetscErrorCode ierr;
279807171c5SHong Zhang   PetscContainer container;
280807171c5SHong Zhang   Mat_RARt       *rart=PETSC_NULL;
281807171c5SHong Zhang 
282807171c5SHong Zhang   PetscFunctionBegin;
283807171c5SHong Zhang   ierr = PetscObjectQuery((PetscObject)A,"Mat_RARt",(PetscObject *)&container);CHKERRQ(ierr);
284807171c5SHong Zhang   if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit");
285807171c5SHong Zhang   ierr = PetscContainerGetPointer(container,(void **)&rart);CHKERRQ(ierr);
286807171c5SHong Zhang   A->ops->destroy   = rart->destroy;
287807171c5SHong Zhang   if (A->ops->destroy) {
288807171c5SHong Zhang     ierr = (*A->ops->destroy)(A);CHKERRQ(ierr);
289807171c5SHong Zhang   }
290807171c5SHong Zhang   ierr = PetscObjectCompose((PetscObject)A,"Mat_RARt",0);CHKERRQ(ierr);
291807171c5SHong Zhang   PetscFunctionReturn(0);
292807171c5SHong Zhang }
293807171c5SHong Zhang 
294807171c5SHong Zhang #undef __FUNCT__
295807171c5SHong Zhang #define __FUNCT__ "MatRARtSymbolic_SeqAIJ_SeqAIJ"
296807171c5SHong Zhang PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat R,PetscReal fill,Mat *C)
297807171c5SHong Zhang {
298807171c5SHong Zhang   PetscErrorCode      ierr;
299807171c5SHong Zhang   Mat                 P;
300807171c5SHong Zhang   PetscInt            *rti,*rtj;
301807171c5SHong Zhang   Mat_RARt            *rart;
302807171c5SHong Zhang   PetscContainer      container;
303807171c5SHong Zhang   MatTransposeColoring matcoloring;
304807171c5SHong Zhang   ISColoring           iscoloring;
3058d0a38e4SHong Zhang   Mat                  Rt_dense,RARt_dense;
306807171c5SHong Zhang   PetscLogDouble       GColor=0.0,MCCreate=0.0,MDenCreate=0.0,t0,tf,etime=0.0;
307807171c5SHong Zhang   Mat_SeqAIJ           *c;
308807171c5SHong Zhang 
309807171c5SHong Zhang   PetscFunctionBegin;
310807171c5SHong Zhang   ierr = PetscGetTime(&t0);CHKERRQ(ierr);
311807171c5SHong Zhang   /* create symbolic P=Rt */
312807171c5SHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(R,&rti,&rtj);CHKERRQ(ierr);
313807171c5SHong Zhang   ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,R->cmap->n,R->rmap->n,rti,rtj,PETSC_NULL,&P);CHKERRQ(ierr);
314807171c5SHong Zhang 
315807171c5SHong Zhang   /* get symbolic C=Pt*A*P */
316807171c5SHong Zhang   ierr = MatPtAPSymbolic_SeqAIJ_SeqAIJ(A,P,fill,C);CHKERRQ(ierr);
317807171c5SHong Zhang 
318807171c5SHong Zhang   /* create a supporting struct */
319807171c5SHong Zhang   ierr = PetscNew(Mat_RARt,&rart);CHKERRQ(ierr);
320807171c5SHong Zhang 
321807171c5SHong Zhang   /* attach the supporting struct to C */
322807171c5SHong Zhang   ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr);
323807171c5SHong Zhang   ierr = PetscContainerSetPointer(container,rart);CHKERRQ(ierr);
324807171c5SHong Zhang   ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_Mat_RARt);CHKERRQ(ierr);
325807171c5SHong Zhang   ierr = PetscObjectCompose((PetscObject)(*C),"Mat_RARt",(PetscObject)container);CHKERRQ(ierr);
326807171c5SHong Zhang   ierr = PetscContainerDestroy(&container);CHKERRQ(ierr);
327807171c5SHong Zhang   ierr = PetscGetTime(&tf);CHKERRQ(ierr);
328807171c5SHong Zhang   etime += tf - t0;
329807171c5SHong Zhang 
330807171c5SHong Zhang   /* Create MatTransposeColoring from symbolic C=R*A*R^T */
331807171c5SHong Zhang   c=(Mat_SeqAIJ*)(*C)->data;
332807171c5SHong Zhang   ierr = PetscGetTime(&t0);CHKERRQ(ierr);
333807171c5SHong Zhang   ierr = MatGetColoring(*C,MATCOLORINGLF,&iscoloring);CHKERRQ(ierr);
334807171c5SHong Zhang   ierr = PetscGetTime(&tf);CHKERRQ(ierr);
335807171c5SHong Zhang   GColor += tf - t0;
336807171c5SHong Zhang 
337807171c5SHong Zhang   ierr = PetscGetTime(&t0);CHKERRQ(ierr);
338807171c5SHong Zhang   ierr = MatTransposeColoringCreate(*C,iscoloring,&matcoloring);CHKERRQ(ierr);
339807171c5SHong Zhang   rart->matcoloring = matcoloring;
340807171c5SHong Zhang   ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr);
341807171c5SHong Zhang   ierr = PetscGetTime(&tf);CHKERRQ(ierr);
342807171c5SHong Zhang   MCCreate += tf - t0;
343807171c5SHong Zhang 
344807171c5SHong Zhang   ierr = PetscGetTime(&t0);CHKERRQ(ierr);
345807171c5SHong Zhang   /* Create Rt_dense */
346807171c5SHong Zhang   ierr = MatCreate(PETSC_COMM_SELF,&Rt_dense);CHKERRQ(ierr);
347807171c5SHong Zhang   ierr = MatSetSizes(Rt_dense,A->cmap->n,matcoloring->ncolors,A->cmap->n,matcoloring->ncolors);CHKERRQ(ierr);
348807171c5SHong Zhang   ierr = MatSetType(Rt_dense,MATSEQDENSE);CHKERRQ(ierr);
349807171c5SHong Zhang   ierr = MatSeqDenseSetPreallocation(Rt_dense,PETSC_NULL);CHKERRQ(ierr);
350807171c5SHong Zhang   Rt_dense->assembled = PETSC_TRUE;
351807171c5SHong Zhang   rart->Rt            = Rt_dense;
352807171c5SHong Zhang 
353807171c5SHong Zhang   /* Create RARt_dense = R*A*Rt_dense */
354807171c5SHong Zhang   ierr = MatCreate(PETSC_COMM_SELF,&RARt_dense);CHKERRQ(ierr);
355807171c5SHong Zhang   ierr = MatSetSizes(RARt_dense,(*C)->rmap->n,matcoloring->ncolors,(*C)->rmap->n,matcoloring->ncolors);CHKERRQ(ierr);
356807171c5SHong Zhang   ierr = MatSetType(RARt_dense,MATSEQDENSE);CHKERRQ(ierr);
357807171c5SHong Zhang   ierr = MatSeqDenseSetPreallocation(RARt_dense,PETSC_NULL);CHKERRQ(ierr);
358807171c5SHong Zhang   rart->RARt = RARt_dense;
359807171c5SHong Zhang 
360807171c5SHong Zhang   /* Allocate work array to store columns of A*R^T used in MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense() */
361807171c5SHong Zhang   ierr = PetscMalloc(A->rmap->n*4*sizeof(PetscScalar),&rart->work);CHKERRQ(ierr);
362807171c5SHong Zhang 
363807171c5SHong Zhang   ierr = PetscGetTime(&tf);CHKERRQ(ierr);
364807171c5SHong Zhang   MDenCreate += tf - t0;
365807171c5SHong Zhang 
366807171c5SHong Zhang   rart->destroy = (*C)->ops->destroy;
367807171c5SHong Zhang   (*C)->ops->destroy = MatDestroy_SeqAIJ_RARt;
368807171c5SHong Zhang 
369807171c5SHong Zhang   /* clean up */
370807171c5SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(R,&rti,&rtj);CHKERRQ(ierr);
371807171c5SHong Zhang   ierr = MatDestroy(&P);CHKERRQ(ierr);
372807171c5SHong Zhang 
373807171c5SHong Zhang #if defined(PETSC_USE_INFO)
374*1ce71dffSSatish Balay   {
375807171c5SHong Zhang   PetscReal density= (PetscReal)(c->nz)/(RARt_dense->rmap->n*RARt_dense->cmap->n);
376807171c5SHong Zhang   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);
377807171c5SHong Zhang   ierr = PetscInfo5(*C,"Sym = GetColor %g + MColorCreate %g + MDenCreate %g + other %g = %g\n",GColor,MCCreate,MDenCreate,etime,GColor+MCCreate+MDenCreate+etime);
378*1ce71dffSSatish Balay   }
379807171c5SHong Zhang #endif
380807171c5SHong Zhang   PetscFunctionReturn(0);
381807171c5SHong Zhang }
382807171c5SHong Zhang 
383807171c5SHong Zhang /*
384807171c5SHong Zhang  RAB = R * A * B, R and A in seqaij format, B in dense format;
385807171c5SHong Zhang */
386807171c5SHong Zhang #undef __FUNCT__
387807171c5SHong Zhang #define __FUNCT__ "MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense"
388807171c5SHong Zhang PetscErrorCode MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense(Mat R,Mat A,Mat B,Mat RAB,PetscScalar *work)
389807171c5SHong Zhang {
390807171c5SHong Zhang   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*r=(Mat_SeqAIJ*)R->data;
391807171c5SHong Zhang   PetscErrorCode ierr;
392807171c5SHong Zhang   PetscScalar    *b,r1,r2,r3,r4,*b1,*b2,*b3,*b4;
393807171c5SHong Zhang   MatScalar      *aa,*ra;
394807171c5SHong Zhang   PetscInt       cn=B->cmap->n,bm=B->rmap->n,col,i,j,n,*ai=a->i,*aj,am=A->rmap->n;
395807171c5SHong Zhang   PetscInt       am2=2*am,am3=3*am,bm4=4*bm;
396807171c5SHong Zhang   PetscScalar    *d,*c,*c2,*c3,*c4;
397807171c5SHong Zhang   PetscInt       *rj,rm=R->rmap->n,dm=RAB->rmap->n,dn=RAB->cmap->n;
398807171c5SHong Zhang   PetscInt       rm2=2*rm,rm3=3*rm,colrm;
399807171c5SHong Zhang 
400807171c5SHong Zhang   PetscFunctionBegin;
401807171c5SHong Zhang   if (!dm || !dn) PetscFunctionReturn(0);
402807171c5SHong Zhang   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);
403807171c5SHong Zhang   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);
404807171c5SHong Zhang   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);
405807171c5SHong Zhang   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);
406807171c5SHong Zhang 
407807171c5SHong Zhang   ierr = MatGetArray(B,&b);CHKERRQ(ierr);
408807171c5SHong Zhang   ierr = MatGetArray(RAB,&d);CHKERRQ(ierr);
409807171c5SHong Zhang   b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm;
410807171c5SHong Zhang   c = work; c2 = c + am; c3 = c2 + am; c4 = c3 + am;
411807171c5SHong Zhang   for (col=0; col<cn-4; col += 4){  /* over columns of C */
412807171c5SHong Zhang     for (i=0; i<am; i++) {        /* over rows of A in those columns */
413807171c5SHong Zhang       r1 = r2 = r3 = r4 = 0.0;
414807171c5SHong Zhang       n   = ai[i+1] - ai[i];
415807171c5SHong Zhang       aj  = a->j + ai[i];
416807171c5SHong Zhang       aa  = a->a + ai[i];
417807171c5SHong Zhang       for (j=0; j<n; j++) {
418807171c5SHong Zhang         r1 += (*aa)*b1[*aj];
419807171c5SHong Zhang         r2 += (*aa)*b2[*aj];
420807171c5SHong Zhang         r3 += (*aa)*b3[*aj];
421807171c5SHong Zhang         r4 += (*aa++)*b4[*aj++];
422807171c5SHong Zhang       }
423807171c5SHong Zhang       c[i]       = r1;
424807171c5SHong Zhang       c[am  + i] = r2;
425807171c5SHong Zhang       c[am2 + i] = r3;
426807171c5SHong Zhang       c[am3 + i] = r4;
427807171c5SHong Zhang     }
428807171c5SHong Zhang     b1 += bm4;
429807171c5SHong Zhang     b2 += bm4;
430807171c5SHong Zhang     b3 += bm4;
431807171c5SHong Zhang     b4 += bm4;
432807171c5SHong Zhang 
433807171c5SHong Zhang     /* RAB[:,col] = R*C[:,col] */
434807171c5SHong Zhang     colrm = col*rm;
435807171c5SHong Zhang     for (i=0; i<rm; i++) {        /* over rows of R in those columns */
436807171c5SHong Zhang       r1 = r2 = r3 = r4 = 0.0;
437807171c5SHong Zhang       n   = r->i[i+1] - r->i[i];
438807171c5SHong Zhang       rj  = r->j + r->i[i];
439807171c5SHong Zhang       ra  = r->a + r->i[i];
440807171c5SHong Zhang       for (j=0; j<n; j++) {
441807171c5SHong Zhang         r1 += (*ra)*c[*rj];
442807171c5SHong Zhang         r2 += (*ra)*c2[*rj];
443807171c5SHong Zhang         r3 += (*ra)*c3[*rj];
444807171c5SHong Zhang         r4 += (*ra++)*c4[*rj++];
445807171c5SHong Zhang       }
446807171c5SHong Zhang       d[colrm + i]       = r1;
447807171c5SHong Zhang       d[colrm + rm + i]  = r2;
448807171c5SHong Zhang       d[colrm + rm2 + i] = r3;
449807171c5SHong Zhang       d[colrm + rm3 + i] = r4;
450807171c5SHong Zhang     }
451807171c5SHong Zhang   }
452807171c5SHong Zhang   for (;col<cn; col++){     /* over extra columns of C */
453807171c5SHong Zhang     for (i=0; i<am; i++) {  /* over rows of A in those columns */
454807171c5SHong Zhang       r1 = 0.0;
455807171c5SHong Zhang       n   = a->i[i+1] - a->i[i];
456807171c5SHong Zhang       aj  = a->j + a->i[i];
457807171c5SHong Zhang       aa  = a->a + a->i[i];
458807171c5SHong Zhang       for (j=0; j<n; j++) {
459807171c5SHong Zhang         r1 += (*aa++)*b1[*aj++];
460807171c5SHong Zhang       }
461807171c5SHong Zhang       c[i]     = r1;
462807171c5SHong Zhang     }
463807171c5SHong Zhang     b1 += bm;
464807171c5SHong Zhang 
465807171c5SHong Zhang     for (i=0; i<rm; i++) {  /* over rows of R in those columns */
466807171c5SHong Zhang       r1 = 0.0;
467807171c5SHong Zhang       n   = r->i[i+1] - r->i[i];
468807171c5SHong Zhang       rj  = r->j + r->i[i];
469807171c5SHong Zhang       ra  = r->a + r->i[i];
470807171c5SHong Zhang       for (j=0; j<n; j++) {
471807171c5SHong Zhang         r1 += (*ra++)*c[*rj++];
472807171c5SHong Zhang       }
473807171c5SHong Zhang       d[col*rm + i]     = r1;
474807171c5SHong Zhang     }
475807171c5SHong Zhang   }
476807171c5SHong Zhang   ierr = PetscLogFlops(cn*2.0*(a->nz + r->nz));CHKERRQ(ierr);
477807171c5SHong Zhang 
478807171c5SHong Zhang   ierr = MatRestoreArray(B,&b);CHKERRQ(ierr);
479807171c5SHong Zhang   ierr = MatRestoreArray(RAB,&d);CHKERRQ(ierr);
480807171c5SHong Zhang   ierr = MatAssemblyBegin(RAB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
481807171c5SHong Zhang   ierr = MatAssemblyEnd(RAB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
482807171c5SHong Zhang   PetscFunctionReturn(0);
483807171c5SHong Zhang }
484807171c5SHong Zhang 
485807171c5SHong Zhang #undef __FUNCT__
486807171c5SHong Zhang #define __FUNCT__ "MatRARtNumeric_SeqAIJ_SeqAIJ"
487807171c5SHong Zhang PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ(Mat A,Mat R,Mat C)
488807171c5SHong Zhang {
489807171c5SHong Zhang   PetscErrorCode        ierr;
490807171c5SHong Zhang   Mat_RARt              *rart;
491807171c5SHong Zhang   PetscContainer        container;
492807171c5SHong Zhang   MatTransposeColoring  matcoloring;
4938d0a38e4SHong Zhang   Mat                   Rt,RARt;
494807171c5SHong Zhang   PetscLogDouble        Mult_sp_den=0.0,app1=0.0,app2=0.0,t0,tf;
495807171c5SHong Zhang 
496807171c5SHong Zhang   PetscFunctionBegin;
497807171c5SHong Zhang   ierr = PetscObjectQuery((PetscObject)C,"Mat_RARt",(PetscObject *)&container);CHKERRQ(ierr);
498807171c5SHong Zhang   if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit");
499807171c5SHong Zhang   ierr  = PetscContainerGetPointer(container,(void **)&rart);CHKERRQ(ierr);
500807171c5SHong Zhang 
501807171c5SHong Zhang   /* Get dense Rt by Apply MatTransposeColoring to R */
502807171c5SHong Zhang   matcoloring = rart->matcoloring;
503807171c5SHong Zhang   Rt          = rart->Rt;
504807171c5SHong Zhang   ierr = PetscGetTime(&t0);CHKERRQ(ierr);
505807171c5SHong Zhang   ierr = MatTransColoringApplySpToDen(matcoloring,R,Rt);CHKERRQ(ierr);
506807171c5SHong Zhang   ierr = PetscGetTime(&tf);CHKERRQ(ierr);
507807171c5SHong Zhang   app1 += tf - t0;
508807171c5SHong Zhang 
509807171c5SHong Zhang   /* Get dense RARt = R*A*Rt */
510807171c5SHong Zhang   ierr = PetscGetTime(&t0);CHKERRQ(ierr);
511807171c5SHong Zhang   RARt = rart->RARt;
512807171c5SHong Zhang   ierr = MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense(R,A,Rt,RARt,rart->work);CHKERRQ(ierr);
513807171c5SHong Zhang   ierr = PetscGetTime(&tf);CHKERRQ(ierr);
514807171c5SHong Zhang   Mult_sp_den += tf - t0;
515807171c5SHong Zhang 
516807171c5SHong Zhang   /* Recover C from C_dense */
517807171c5SHong Zhang   ierr = PetscGetTime(&t0);CHKERRQ(ierr);
518807171c5SHong Zhang   ierr = MatTransColoringApplyDenToSp(matcoloring,RARt,C);CHKERRQ(ierr);
519807171c5SHong Zhang   ierr = PetscGetTime(&tf);CHKERRQ(ierr);
520807171c5SHong Zhang   app2 += tf - t0;
521807171c5SHong Zhang 
522807171c5SHong Zhang #if defined(PETSC_USE_INFO)
523807171c5SHong Zhang   ierr = PetscInfo4(C,"Num = ColorApp %g + %g + Mult_sp_den %g  = %g\n",app1,app2,Mult_sp_den,app1+app2+Mult_sp_den);
524807171c5SHong Zhang #endif
525807171c5SHong Zhang   PetscFunctionReturn(0);
526807171c5SHong Zhang }
527807171c5SHong Zhang 
528807171c5SHong Zhang EXTERN_C_BEGIN
529807171c5SHong Zhang #undef __FUNCT__
530807171c5SHong Zhang #define __FUNCT__ "MatRARt_SeqAIJ_SeqAIJ"
531807171c5SHong Zhang PetscErrorCode MatRARt_SeqAIJ_SeqAIJ(Mat A,Mat R,MatReuse scall,PetscReal fill,Mat *C)
532807171c5SHong Zhang {
533807171c5SHong Zhang   PetscErrorCode ierr;
534807171c5SHong Zhang 
535807171c5SHong Zhang   PetscFunctionBegin;
536807171c5SHong Zhang   if (scall == MAT_INITIAL_MATRIX){
537807171c5SHong Zhang     ierr = MatRARtSymbolic_SeqAIJ_SeqAIJ(A,R,fill,C);CHKERRQ(ierr);
538807171c5SHong Zhang   }
539807171c5SHong Zhang   ierr = MatRARtNumeric_SeqAIJ_SeqAIJ(A,R,*C);CHKERRQ(ierr);
540807171c5SHong Zhang   PetscFunctionReturn(0);
541807171c5SHong Zhang }
542807171c5SHong Zhang EXTERN_C_END
543