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