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