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