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