xref: /petsc/src/mat/impls/aij/seq/matptap.c (revision bf388a1f1fb7cbb324cb4eba80ecc0a331d4a847)
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=PETSC_NULL,current_space=PETSC_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;
52   PetscInt           i,j,k,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi,nlnk,*lnk;
53   MatScalar          *ca;
54   PetscBT            lnkbt;
55 
56   PetscFunctionBegin;
57   /* Get ij structure of P^T */
58   ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
59   ptJ=ptj;
60 
61   /* Allocate ci array, arrays for fill computation and */
62   /* free space for accumulating nonzero column info */
63   ierr = PetscMalloc((pn+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
64   ci[0] = 0;
65 
66   ierr = PetscMalloc((2*an+1)*sizeof(PetscInt),&ptadenserow);CHKERRQ(ierr);
67   ierr = PetscMemzero(ptadenserow,(2*an+1)*sizeof(PetscInt));CHKERRQ(ierr);
68   ptasparserow = ptadenserow  + an;
69 
70   /* create and initialize a linked list */
71   nlnk = pn+1;
72   ierr = PetscLLCreate(pn,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr);
73 
74   /* Set initial free space to be fill*nnz(A). */
75   /* This should be reasonable if sparsity of PtAP is similar to that of A. */
76   ierr          = PetscFreeSpaceGet((PetscInt)(fill*ai[am]),&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     current_space->array           += cnzi;
117     current_space->local_used      += cnzi;
118     current_space->local_remaining -= cnzi;
119 
120     for (j=0;j<ptanzi;j++) {
121       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(((PetscObject)A)->comm,pn,pn,ci,cj,ca,C);CHKERRQ(ierr);
141   (*C)->rmap->bs = P->cmap->bs;
142   (*C)->cmap->bs = P->cmap->bs;
143 
144   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
145   /* Since these are PETSc arrays, change flags to free them as necessary. */
146   c = (Mat_SeqAIJ *)((*C)->data);
147   c->free_a  = PETSC_TRUE;
148   c->free_ij = PETSC_TRUE;
149   c->nonew   = 0;
150   (*C)->ops->ptapnumeric = MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy;
151 
152   /* Clean up. */
153   ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
154 #if defined(PETSC_USE_INFO)
155   if (ci[pn] != 0) {
156     PetscReal afill = ((PetscReal)ci[pn])/ai[am];
157     if (afill < 1.0) afill = 1.0;
158     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
159     ierr = PetscInfo1((*C),"Use MatPtAP(A,P,MatReuse,%G,&C) for best performance.\n",afill);CHKERRQ(ierr);
160   } else {
161     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
162   }
163 #endif
164   PetscFunctionReturn(0);
165 }
166 
167 #undef __FUNCT__
168 #define __FUNCT__ "MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy"
169 PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy(Mat A,Mat P,Mat C)
170 {
171   PetscErrorCode ierr;
172   Mat_SeqAIJ     *a  = (Mat_SeqAIJ *) A->data;
173   Mat_SeqAIJ     *p  = (Mat_SeqAIJ *) P->data;
174   Mat_SeqAIJ     *c  = (Mat_SeqAIJ *) C->data;
175   PetscInt       *ai=a->i,*aj=a->j,*apj,*apjdense,*pi=p->i,*pj=p->j,*pJ=p->j,*pjj;
176   PetscInt       *ci=c->i,*cj=c->j,*cjj;
177   PetscInt       am=A->rmap->N,cn=C->cmap->N,cm=C->rmap->N;
178   PetscInt       i,j,k,anzi,pnzi,apnzj,nextap,pnzj,prow,crow;
179   MatScalar      *aa=a->a,*apa,*pa=p->a,*pA=p->a,*paj,*ca=c->a,*caj;
180 
181   PetscFunctionBegin;
182   /* Allocate temporary array for storage of one row of A*P (cn: non-scalable) */
183   ierr = PetscMalloc(cn*(sizeof(MatScalar)+sizeof(PetscInt))+c->rmax*sizeof(PetscInt),&apa);CHKERRQ(ierr);
184 
185   apjdense = (PetscInt *)(apa + cn);
186   apj      = apjdense + cn;
187   ierr = PetscMemzero(apa,cn*(sizeof(MatScalar)+sizeof(PetscInt)));CHKERRQ(ierr);
188 
189   /* Clear old values in C */
190   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
191 
192   for (i=0; i<am; i++) {
193     /* Form sparse row of A*P */
194     anzi  = ai[i+1] - ai[i];
195     apnzj = 0;
196     for (j=0; j<anzi; j++) {
197       prow = *aj++;
198       pnzj = pi[prow+1] - pi[prow];
199       pjj  = pj + pi[prow];
200       paj  = pa + pi[prow];
201       for (k=0;k<pnzj;k++) {
202         if (!apjdense[pjj[k]]) {
203           apjdense[pjj[k]] = -1;
204           apj[apnzj++]     = pjj[k];
205         }
206         apa[pjj[k]] += (*aa)*paj[k];
207       }
208       ierr = PetscLogFlops(2.0*pnzj);CHKERRQ(ierr);
209       aa++;
210     }
211 
212     /* Sort the j index array for quick sparse axpy. */
213     /* Note: a array does not need sorting as it is in dense storage locations. */
214     ierr = PetscSortInt(apnzj,apj);CHKERRQ(ierr);
215 
216     /* Compute P^T*A*P using outer product (P^T)[:,j]*(A*P)[j,:]. */
217     pnzi = pi[i+1] - pi[i];
218     for (j=0; j<pnzi; j++) {
219       nextap = 0;
220       crow   = *pJ++;
221       cjj    = cj + ci[crow];
222       caj    = ca + ci[crow];
223       /* Perform sparse axpy operation.  Note cjj includes apj. */
224       for (k=0;nextap<apnzj;k++) {
225 #if defined(PETSC_USE_DEBUG)
226         if (k >= ci[crow+1] - ci[crow]) {
227           SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"k too large k %d, crow %d",k,crow);
228         }
229 #endif
230         if (cjj[k]==apj[nextap]) {
231           caj[k] += (*pA)*apa[apj[nextap++]];
232         }
233       }
234       ierr = PetscLogFlops(2.0*apnzj);CHKERRQ(ierr);
235       pA++;
236     }
237 
238     /* Zero the current row info for A*P */
239     for (j=0; j<apnzj; j++) {
240       apa[apj[j]]      = 0.;
241       apjdense[apj[j]] = 0;
242     }
243   }
244 
245   /* Assemble the final matrix and clean up */
246   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
247   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
248 
249   ierr = PetscFree(apa);CHKERRQ(ierr);
250   PetscFunctionReturn(0);
251 }
252 
253 #undef __FUNCT__
254 #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqAIJ"
255 PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
256 {
257   PetscErrorCode     ierr;
258   Mat_SeqAIJ         *ap,*c;
259   PetscInt           *api,*apj,*ci,pn=P->cmap->N;
260   MatScalar          *ca;
261   Mat_PtAP           *ptap;
262   Mat                Pt,AP;
263   PetscBool          sparse_axpy=PETSC_TRUE;
264 
265   PetscFunctionBegin;
266   ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
267   /* flag 'sparse_axpy' determines which implementations to be used:
268        0: do dense axpy in MatPtAPNumeric() - fastest, but requires storage of struct A*P;
269        1: do two sparse axpy in MatPtAPNumeric() - slowest, does not store structure of A*P. */
270   ierr = PetscOptionsBool("-matptap_scalable","Use sparse axpy but slower MatPtAPNumeric()","",sparse_axpy,&sparse_axpy,PETSC_NULL);CHKERRQ(ierr);
271   ierr = PetscOptionsEnd();CHKERRQ(ierr);
272   if (sparse_axpy){
273     ierr = MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(A,P,fill,C);CHKERRQ(ierr);
274     PetscFunctionReturn(0);
275   }
276 
277   /* Get symbolic Pt = P^T */
278   ierr = MatTransposeSymbolic_SeqAIJ(P,&Pt);CHKERRQ(ierr);
279 
280   /* Get symbolic AP = A*P */
281   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,P,fill,&AP);CHKERRQ(ierr);
282 
283   ap          = (Mat_SeqAIJ*)AP->data;
284   api         = ap->i;
285   apj         = ap->j;
286   ap->free_ij = PETSC_FALSE; /* api and apj are kept in struct ptap, cannot be destroyed with AP */
287 
288   /* Get C = Pt*AP */
289   ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(Pt,AP,fill,C);CHKERRQ(ierr);
290 
291   c  = (Mat_SeqAIJ*)(*C)->data;
292   ci = c->i;
293   ierr = PetscMalloc((ci[pn]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
294   ierr = PetscMemzero(ca,(ci[pn]+1)*sizeof(MatScalar));CHKERRQ(ierr);
295   c->a       = ca;
296   c->free_a  = PETSC_TRUE;
297 
298   /* Create a supporting struct for reuse by MatPtAPNumeric() */
299   ierr = PetscNew(Mat_PtAP,&ptap);CHKERRQ(ierr);
300   c->ptap            = ptap;
301   ptap->destroy      = (*C)->ops->destroy;
302   (*C)->ops->destroy = MatDestroy_SeqAIJ_PtAP;
303 
304   /* Allocate temporary array for storage of one row of A*P */
305   ierr = PetscMalloc((pn+1)*sizeof(PetscScalar),&ptap->apa);CHKERRQ(ierr);
306   ierr = PetscMemzero(ptap->apa,(pn+1)*sizeof(PetscScalar));CHKERRQ(ierr);
307   (*C)->ops->ptapnumeric = MatPtAPNumeric_SeqAIJ_SeqAIJ;
308 
309   ptap->api = api;
310   ptap->apj = apj;
311 
312   /* Clean up. */
313   ierr = MatDestroy(&Pt);CHKERRQ(ierr);
314   ierr = MatDestroy(&AP);CHKERRQ(ierr);
315 #if defined(PETSC_USE_INFO)
316   ierr = PetscInfo2((*C),"given fill %G, use scalable %d\n",fill,sparse_axpy);CHKERRQ(ierr);
317 #endif
318   PetscFunctionReturn(0);
319 }
320 
321 /* #define PROFILE_MatPtAPNumeric */
322 #undef __FUNCT__
323 #define __FUNCT__ "MatPtAPNumeric_SeqAIJ_SeqAIJ"
324 PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat A,Mat P,Mat C)
325 {
326   PetscErrorCode    ierr;
327   Mat_SeqAIJ        *a  = (Mat_SeqAIJ *) A->data;
328   Mat_SeqAIJ        *p  = (Mat_SeqAIJ *) P->data;
329   Mat_SeqAIJ        *c  = (Mat_SeqAIJ *) C->data;
330   const PetscInt    *ai=a->i,*aj=a->j,*pi=p->i,*pj=p->j,*ci=c->i,*cj=c->j;
331   const PetscScalar *aa=a->a,*pa=p->a,*pval;
332   const PetscInt    *apj,*pcol,*cjj;
333   const PetscInt    am=A->rmap->N,cm=C->rmap->N;
334   PetscInt          i,j,k,anz,apnz,pnz,prow,crow,cnz;
335   PetscScalar       *apa,*ca=c->a,*caj,pvalj;
336   Mat_PtAP          *ptap = c->ptap;
337 #if defined(PROFILE_MatPtAPNumeric)
338   PetscLogDouble    t0,tf,time_Cseq0=0.0,time_Cseq1=0.0;
339   PetscInt          flops0=0,flops1=0;
340 #endif
341 
342   PetscFunctionBegin;
343   /* Get temporary array for storage of one row of A*P */
344   apa = ptap->apa;
345 
346   /* Clear old values in C */
347   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
348 
349   for (i=0;i<am;i++) {
350     /* Form sparse row of AP[i,:] = A[i,:]*P */
351 #if defined(PROFILE_MatPtAPNumeric)
352     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
353 #endif
354     anz  = ai[i+1] - ai[i];
355     apnz = 0;
356     for (j=0; j<anz; j++) {
357       prow = aj[j];
358       pnz  = pi[prow+1] - pi[prow];
359       pcol = pj + pi[prow];
360       pval = pa + pi[prow];
361       for (k=0; k<pnz; k++) {
362         apa[pcol[k]] += aa[j]*pval[k];
363       }
364       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
365 #if defined(PROFILE_MatPtAPNumeric)
366       flops0 += 2.0*pnz;
367 #endif
368     }
369     aj += anz; aa += anz;
370 #if defined(PROFILE_MatPtAPNumeric)
371     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
372     time_Cseq0 += tf - t0;
373 #endif
374 
375     /* Compute P^T*A*P using outer product P[i,:]^T*AP[i,:]. */
376 #if defined(PROFILE_MatPtAPNumeric)
377     ierr = PetscGetTime(&t0);CHKERRQ(ierr);
378 #endif
379     apj  = ptap->apj + ptap->api[i];
380     apnz = ptap->api[i+1] - ptap->api[i];
381     pnz  = pi[i+1] - pi[i];
382     pcol = pj + pi[i];
383     pval = pa + pi[i];
384 
385     /* Perform dense axpy */
386     for (j=0; j<pnz; j++) {
387       crow  = pcol[j];
388       cjj   = cj + ci[crow];
389       caj   = ca + ci[crow];
390       pvalj = pval[j];
391       cnz   = ci[crow+1] - ci[crow];
392       for (k=0; k<cnz; k++){
393         caj[k] += pvalj*apa[cjj[k]];
394       }
395       ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
396 #if defined(PROFILE_MatPtAPNumeric)
397       flops1 += 2.0*cnz;
398 #endif
399     }
400 #if defined(PROFILE_MatPtAPNumeric)
401     ierr = PetscGetTime(&tf);CHKERRQ(ierr);
402     time_Cseq1 += tf - t0;
403 #endif
404 
405     /* Zero the current row info for A*P */
406     for (j=0; j<apnz; j++) apa[apj[j]] = 0.0;
407   }
408 
409   /* Assemble the final matrix and clean up */
410   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
411   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
412 #if defined(PROFILE_MatPtAPNumeric)
413   printf("PtAPNumeric_SeqAIJ time %g + %g, flops %d %d\n",time_Cseq0,time_Cseq1,flops0,flops1);
414 #endif
415   PetscFunctionReturn(0);
416 }
417