xref: /petsc/src/mat/impls/aij/seq/matptap.c (revision 0ec63f53104f4b8d70d10e7c493a96759e264eea)
1 
2 /*
3   Defines projective product routines where A is a SeqAIJ matrix
4           C = P^T * A * P
5 */
6 
7 #include <../src/mat/impls/aij/seq/aij.h>   /*I "petscmat.h" I*/
8 #include <../src/mat/utils/freespace.h>
9 #include <petscbt.h>
10 
11 #undef __FUNCT__
12 #define __FUNCT__ "MatPtAP_SeqAIJ_SeqAIJ"
13 PetscErrorCode MatPtAP_SeqAIJ_SeqAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
14 {
15   PetscErrorCode ierr;
16 
17   PetscFunctionBegin;
18   if (scall == MAT_INITIAL_MATRIX){
19     ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
20     ierr = MatPtAPSymbolic_SeqAIJ_SeqAIJ(A,P,fill,C);CHKERRQ(ierr);
21     ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
22   }
23   ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
24   //ierr = (*(A)->ops->ptapnumeric)(A,P,*C);CHKERRQ(ierr); /* should use matrix *C, not A! */
25   ierr = (*(*C)->ops->ptapnumeric)(A,P,*C);CHKERRQ(ierr);
26   ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
27   PetscFunctionReturn(0);
28 }
29 
30 #undef __FUNCT__
31 #define __FUNCT__ "MatDestroy_SeqAIJ_PtAP"
32 PetscErrorCode MatDestroy_SeqAIJ_PtAP(Mat A)
33 {
34   PetscErrorCode ierr;
35   Mat_SeqAIJ     *a = (Mat_SeqAIJ *)A->data;
36   Mat_PtAP       *ptap = a->ptap;
37 
38   PetscFunctionBegin;
39   /* free ptap, then A */
40   ierr = PetscFree(ptap->apa);CHKERRQ(ierr);
41   ierr = PetscFree(ptap->api);CHKERRQ(ierr);
42   ierr = PetscFree(ptap->apj);CHKERRQ(ierr);
43   ierr = (ptap->destroy)(A);CHKERRQ(ierr);
44   ierr = PetscFree(ptap);CHKERRQ(ierr);
45   PetscFunctionReturn(0);
46 }
47 
48 #undef __FUNCT__
49 #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy"
50 PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(Mat A,Mat P,PetscReal fill,Mat *C)
51 {
52   PetscErrorCode     ierr;
53   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
54   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data,*p = (Mat_SeqAIJ*)P->data,*c;
55   PetscInt           *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj;
56   PetscInt           *ci,*cj,*ptadenserow,*ptasparserow,*ptaj,nspacedouble=0;
57   PetscInt           an=A->cmap->N,am=A->rmap->N,pn=P->cmap->N;
58   PetscInt           i,j,k,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi,nlnk,*lnk;
59   MatScalar          *ca;
60   PetscBT            lnkbt;
61 
62   PetscFunctionBegin;
63   /* Get ij structure of P^T */
64   ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
65   ptJ=ptj;
66 
67   /* Allocate ci array, arrays for fill computation and */
68   /* free space for accumulating nonzero column info */
69   ierr = PetscMalloc((pn+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
70   ci[0] = 0;
71 
72   ierr = PetscMalloc((2*an+1)*sizeof(PetscInt),&ptadenserow);CHKERRQ(ierr);
73   ierr = PetscMemzero(ptadenserow,(2*an+1)*sizeof(PetscInt));CHKERRQ(ierr);
74   ptasparserow = ptadenserow  + an;
75 
76   /* create and initialize a linked list */
77   nlnk = pn+1;
78   ierr = PetscLLCreate(pn,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr);
79 
80   /* Set initial free space to be fill*nnz(A). */
81   /* This should be reasonable if sparsity of PtAP is similar to that of A. */
82   ierr          = PetscFreeSpaceGet((PetscInt)(fill*ai[am]),&free_space);
83   current_space = free_space;
84 
85   /* Determine symbolic info for each row of C: */
86   for (i=0;i<pn;i++) {
87     ptnzi  = pti[i+1] - pti[i];
88     ptanzi = 0;
89     /* Determine symbolic row of PtA: */
90     for (j=0;j<ptnzi;j++) {
91       arow = *ptJ++;
92       anzj = ai[arow+1] - ai[arow];
93       ajj  = aj + ai[arow];
94       for (k=0;k<anzj;k++) {
95         if (!ptadenserow[ajj[k]]) {
96           ptadenserow[ajj[k]]    = -1;
97           ptasparserow[ptanzi++] = ajj[k];
98         }
99       }
100     }
101     /* Using symbolic info for row of PtA, determine symbolic info for row of C: */
102     ptaj = ptasparserow;
103     cnzi   = 0;
104     for (j=0;j<ptanzi;j++) {
105       prow = *ptaj++;
106       pnzj = pi[prow+1] - pi[prow];
107       pjj  = pj + pi[prow];
108       /* add non-zero cols of P into the sorted linked list lnk */
109       ierr = PetscLLAddSorted(pnzj,pjj,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr);
110       cnzi += nlnk;
111     }
112 
113     /* If free space is not available, make more free space */
114     /* Double the amount of total space in the list */
115     if (current_space->local_remaining<cnzi) {
116       ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,&current_space);CHKERRQ(ierr);
117       nspacedouble++;
118     }
119 
120     /* Copy data into free space, and zero out denserows */
121     ierr = PetscLLClean(pn,pn,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
122     current_space->array           += cnzi;
123     current_space->local_used      += cnzi;
124     current_space->local_remaining -= cnzi;
125 
126     for (j=0;j<ptanzi;j++) {
127       ptadenserow[ptasparserow[j]] = 0;
128     }
129     /* Aside: Perhaps we should save the pta info for the numerical factorization. */
130     /*        For now, we will recompute what is needed. */
131     ci[i+1] = ci[i] + cnzi;
132   }
133   /* nnz is now stored in ci[ptm], column indices are in the list of free space */
134   /* Allocate space for cj, initialize cj, and */
135   /* destroy list of free space and other temporary array(s) */
136   ierr = PetscMalloc((ci[pn]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
137   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
138   ierr = PetscFree(ptadenserow);CHKERRQ(ierr);
139   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
140 
141   /* Allocate space for ca */
142   ierr = PetscMalloc((ci[pn]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
143   ierr = PetscMemzero(ca,(ci[pn]+1)*sizeof(MatScalar));CHKERRQ(ierr);
144 
145   /* put together the new matrix */
146   ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,pn,pn,ci,cj,ca,C);CHKERRQ(ierr);
147   (*C)->rmap->bs = P->cmap->bs;
148   (*C)->cmap->bs = P->cmap->bs;
149 
150   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
151   /* Since these are PETSc arrays, change flags to free them as necessary. */
152   c = (Mat_SeqAIJ *)((*C)->data);
153   c->free_a  = PETSC_TRUE;
154   c->free_ij = PETSC_TRUE;
155   c->nonew   = 0;
156   (*C)->ops->ptapnumeric = MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy; /* should use *C->ops until PtAP insterface is updated to double dispatch as MatMatMult() */
157 
158   /* Clean up. */
159   ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
160 #if defined(PETSC_USE_INFO)
161   if (ci[pn] != 0) {
162     PetscReal afill = ((PetscReal)ci[pn])/ai[am];
163     if (afill < 1.0) afill = 1.0;
164     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
165     ierr = PetscInfo1((*C),"Use MatPtAP(A,P,MatReuse,%G,&C) for best performance.\n",afill);CHKERRQ(ierr);
166   } else {
167     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
168   }
169 #endif
170   PetscFunctionReturn(0);
171 }
172 
173 #undef __FUNCT__
174 #define __FUNCT__ "MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy"
175 PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy(Mat A,Mat P,Mat C)
176 {
177   PetscErrorCode ierr;
178   Mat_SeqAIJ     *a  = (Mat_SeqAIJ *) A->data;
179   Mat_SeqAIJ     *p  = (Mat_SeqAIJ *) P->data;
180   Mat_SeqAIJ     *c  = (Mat_SeqAIJ *) C->data;
181   PetscInt       *ai=a->i,*aj=a->j,*apj,*apjdense,*pi=p->i,*pj=p->j,*pJ=p->j,*pjj;
182   PetscInt       *ci=c->i,*cj=c->j,*cjj;
183   PetscInt       am=A->rmap->N,cn=C->cmap->N,cm=C->rmap->N;
184   PetscInt       i,j,k,anzi,pnzi,apnzj,nextap,pnzj,prow,crow;
185   MatScalar      *aa=a->a,*apa,*pa=p->a,*pA=p->a,*paj,*ca=c->a,*caj;
186 
187   PetscFunctionBegin;
188   /* Allocate temporary array for storage of one row of A*P (cn: non-scalable) */
189   ierr = PetscMalloc(cn*(sizeof(MatScalar)+sizeof(PetscInt))+c->rmax*sizeof(PetscInt),&apa);CHKERRQ(ierr);
190 
191   apjdense = (PetscInt *)(apa + cn);
192   apj      = apjdense + cn;
193   ierr = PetscMemzero(apa,cn*(sizeof(MatScalar)+sizeof(PetscInt)));CHKERRQ(ierr);
194 
195   /* Clear old values in C */
196   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
197 
198   for (i=0; i<am; i++) {
199     /* Form sparse row of A*P */
200     anzi  = ai[i+1] - ai[i];
201     apnzj = 0;
202     for (j=0; j<anzi; j++) {
203       prow = *aj++;
204       pnzj = pi[prow+1] - pi[prow];
205       pjj  = pj + pi[prow];
206       paj  = pa + pi[prow];
207       for (k=0;k<pnzj;k++) {
208         if (!apjdense[pjj[k]]) {
209           apjdense[pjj[k]] = -1;
210           apj[apnzj++]     = pjj[k];
211         }
212         apa[pjj[k]] += (*aa)*paj[k];
213       }
214       ierr = PetscLogFlops(2.0*pnzj);CHKERRQ(ierr);
215       aa++;
216     }
217 
218     /* Sort the j index array for quick sparse axpy. */
219     /* Note: a array does not need sorting as it is in dense storage locations. */
220     ierr = PetscSortInt(apnzj,apj);CHKERRQ(ierr);
221 
222     /* Compute P^T*A*P using outer product (P^T)[:,j]*(A*P)[j,:]. */
223     pnzi = pi[i+1] - pi[i];
224     for (j=0; j<pnzi; j++) {
225       nextap = 0;
226       crow   = *pJ++;
227       cjj    = cj + ci[crow];
228       caj    = ca + ci[crow];
229       /* Perform sparse axpy operation.  Note cjj includes apj. */
230       for (k=0;nextap<apnzj;k++) {
231 #if defined(PETSC_USE_DEBUG)
232         if (k >= ci[crow+1] - ci[crow]) {
233           SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"k too large k %d, crow %d",k,crow);
234         }
235 #endif
236         if (cjj[k]==apj[nextap]) {
237           caj[k] += (*pA)*apa[apj[nextap++]];
238         }
239       }
240       ierr = PetscLogFlops(2.0*apnzj);CHKERRQ(ierr);
241       pA++;
242     }
243 
244     /* Zero the current row info for A*P */
245     for (j=0; j<apnzj; j++) {
246       apa[apj[j]]      = 0.;
247       apjdense[apj[j]] = 0;
248     }
249   }
250 
251   /* Assemble the final matrix and clean up */
252   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
253   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
254 
255   ierr = PetscFree(apa);CHKERRQ(ierr);
256   PetscFunctionReturn(0);
257 }
258 
259 #undef __FUNCT__
260 #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqAIJ"
261 PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
262 {
263   PetscErrorCode     ierr;
264   Mat_SeqAIJ         *ap,*c;
265   PetscInt           *api,*apj,*ci,pn=P->cmap->N;
266   MatScalar          *ca;
267   Mat_PtAP           *ptap;
268   Mat                Pt,AP;
269   PetscBool          sparse_axpy=PETSC_TRUE;
270 
271   PetscFunctionBegin;
272   ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
273   /* flag 'sparse_axpy' determines which implementations to be used:
274        0: do dense axpy in MatPtAPNumeric() - fastest, but requires storage of struct A*P;
275        1: do two sparse axpy in MatPtAPNumeric() - slowest, does not store structure of A*P. */
276   ierr = PetscOptionsBool("-matptap_scalable","Use sparse axpy but slower MatPtAPNumeric()","",sparse_axpy,&sparse_axpy,PETSC_NULL);CHKERRQ(ierr);
277   ierr = PetscOptionsEnd();CHKERRQ(ierr);
278   if (sparse_axpy){
279     ierr = MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(A,P,fill,C);CHKERRQ(ierr);
280     PetscFunctionReturn(0);
281   }
282 
283   /* Get symbolic Pt = P^T */
284   ierr = MatTransposeSymbolic_SeqAIJ(P,&Pt);CHKERRQ(ierr);
285 
286   /* Get symbolic AP = A*P */
287   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,P,fill,&AP);CHKERRQ(ierr);
288 
289   ap          = (Mat_SeqAIJ*)AP->data;
290   api         = ap->i;
291   apj         = ap->j;
292   ap->free_ij = PETSC_FALSE; /* api and apj are kept in struct ptap, cannot be destroyed with AP */
293 
294   /* Get C = Pt*AP */
295   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(Pt,AP,fill,C);CHKERRQ(ierr);
296 
297   c  = (Mat_SeqAIJ*)(*C)->data;
298   ci = c->i;
299   ierr = PetscMalloc((ci[pn]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
300   ierr = PetscMemzero(ca,(ci[pn]+1)*sizeof(MatScalar));CHKERRQ(ierr);
301   c->a       = ca;
302   c->free_a  = PETSC_TRUE;
303 
304   /* Create a supporting struct for reuse by MatPtAPNumeric() */
305   ierr = PetscNew(Mat_PtAP,&ptap);CHKERRQ(ierr);
306   c->ptap            = ptap;
307   ptap->destroy      = (*C)->ops->destroy;
308   (*C)->ops->destroy = MatDestroy_SeqAIJ_PtAP;
309 
310   /* Allocate temporary array for storage of one row of A*P */
311   ierr = PetscMalloc((pn+1)*sizeof(PetscScalar),&ptap->apa);CHKERRQ(ierr);
312   ierr = PetscMemzero(ptap->apa,(pn+1)*sizeof(PetscScalar));CHKERRQ(ierr);
313   (*C)->ops->ptapnumeric = MatPtAPNumeric_SeqAIJ_SeqAIJ;
314 
315   ptap->api = api;
316   ptap->apj = apj;
317 
318   /* Clean up. */
319   ierr = MatDestroy(&Pt);CHKERRQ(ierr);
320   ierr = MatDestroy(&AP);CHKERRQ(ierr);
321 #if defined(PETSC_USE_INFO)
322   ierr = PetscInfo2((*C),"given fill %G, use scalable %d\n",fill,sparse_axpy);CHKERRQ(ierr);
323 #endif
324   PetscFunctionReturn(0);
325 }
326 
327 /* #define PROFILE_MatPtAPNumeric */
328 #undef __FUNCT__
329 #define __FUNCT__ "MatPtAPNumeric_SeqAIJ_SeqAIJ"
330 PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat A,Mat P,Mat C)
331 {
332   PetscErrorCode    ierr;
333   Mat_SeqAIJ        *a  = (Mat_SeqAIJ *) A->data;
334   Mat_SeqAIJ        *p  = (Mat_SeqAIJ *) P->data;
335   Mat_SeqAIJ        *c  = (Mat_SeqAIJ *) C->data;
336   const PetscInt    *ai=a->i,*aj=a->j,*pi=p->i,*pj=p->j,*ci=c->i,*cj=c->j;
337   const PetscScalar *aa=a->a,*pa=p->a,*pval;
338   const PetscInt    *apj,*pcol,*cjj;
339   const PetscInt    am=A->rmap->N,cm=C->rmap->N;
340   PetscInt          i,j,k,anz,apnz,pnz,prow,crow,cnz;
341   PetscScalar       *apa,*ca=c->a,*caj,pvalj;
342   Mat_PtAP          *ptap = c->ptap;
343 #if defined(PROFILE_MatPtAPNumeric)
344   PetscLogDouble    t0,tf,time_Cseq0=0.0,time_Cseq1=0.0;
345   PetscInt          flops0=0,flops1=0;
346 #endif
347 
348   PetscFunctionBegin;
349   /* Get temporary array for storage of one row of A*P */
350   apa = ptap->apa;
351 
352   /* Clear old values in C */
353   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
354 
355   for (i=0;i<am;i++) {
356     /* Form sparse row of AP[i,:] = A[i,:]*P */
357 #if defined(PROFILE_MatPtAPNumeric)
358     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
359 #endif
360     anz  = ai[i+1] - ai[i];
361     apnz = 0;
362     for (j=0; j<anz; j++) {
363       prow = aj[j];
364       pnz  = pi[prow+1] - pi[prow];
365       pcol = pj + pi[prow];
366       pval = pa + pi[prow];
367       for (k=0; k<pnz; k++) {
368         apa[pcol[k]] += aa[j]*pval[k];
369       }
370       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
371 #if defined(PROFILE_MatPtAPNumeric)
372       flops0 += 2.0*pnz;
373 #endif
374     }
375     aj += anz; aa += anz;
376 #if defined(PROFILE_MatPtAPNumeric)
377     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
378     time_Cseq0 += tf - t0;
379 #endif
380 
381     /* Compute P^T*A*P using outer product P[i,:]^T*AP[i,:]. */
382 #if defined(PROFILE_MatPtAPNumeric)
383     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
384 #endif
385     apj  = ptap->apj + ptap->api[i];
386     apnz = ptap->api[i+1] - ptap->api[i];
387     pnz  = pi[i+1] - pi[i];
388     pcol = pj + pi[i];
389     pval = pa + pi[i];
390 
391     /* Perform dense axpy */
392     for (j=0; j<pnz; j++) {
393       crow  = pcol[j];
394       cjj   = cj + ci[crow];
395       caj   = ca + ci[crow];
396       pvalj = pval[j];
397       cnz   = ci[crow+1] - ci[crow];
398       for (k=0; k<cnz; k++){
399         caj[k] += pvalj*apa[cjj[k]];
400       }
401       ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
402 #if defined(PROFILE_MatPtAPNumeric)
403       flops1 += 2.0*cnz;
404 #endif
405     }
406 #if defined(PROFILE_MatPtAPNumeric)
407     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
408     time_Cseq1 += tf - t0;
409 #endif
410 
411     /* Zero the current row info for A*P */
412     for (j=0; j<apnz; j++) apa[apj[j]] = 0.0;
413   }
414 
415   /* Assemble the final matrix and clean up */
416   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
417   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
418 #if defined(PROFILE_MatPtAPNumeric)
419   printf("PtAPNumeric_SeqAIJ time %g + %g, flops %d %d\n",time_Cseq0,time_Cseq1,flops0,flops1);
420 #endif
421   PetscFunctionReturn(0);
422 }
423