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