xref: /petsc/src/mat/impls/aij/seq/matptap.c (revision d5c9c0c4eebc2f2a01a1bd0c86fca87e2acd2a03)
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 PetscErrorCode MatProductSymbolic_PtAP_SeqAIJ_SeqAIJ(Mat C)
17 {
18   PetscErrorCode      ierr;
19   Mat_Product         *product = C->product;
20   Mat                 A=product->A,P=product->B;
21   MatProductAlgorithm alg=product->alg;
22   PetscReal           fill=product->fill;
23   PetscBool           flg;
24   Mat                 Pt;
25   Mat_MatTransMatMult *atb;
26   Mat_SeqAIJ          *c;
27 
28   PetscFunctionBegin;
29   /* "scalable" */
30   ierr = PetscStrcmp(alg,"scalable",&flg);CHKERRQ(ierr);
31   if (flg) {
32     ierr = MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(A,P,fill,C);CHKERRQ(ierr);
33     C->ops->productnumeric = MatProductNumeric_PtAP;
34     PetscFunctionReturn(0);
35   }
36 
37   /* "rap" */
38   ierr = PetscStrcmp(alg,"rap",&flg);CHKERRQ(ierr);
39   if (flg) { /* Set default algorithm */
40     ierr = PetscNew(&atb);CHKERRQ(ierr);
41     ierr = MatTranspose_SeqAIJ(P,MAT_INITIAL_MATRIX,&Pt);CHKERRQ(ierr);
42     ierr = MatMatMatMultSymbolic_SeqAIJ_SeqAIJ_SeqAIJ(Pt,A,P,fill,C);CHKERRQ(ierr);
43 
44     c                      = (Mat_SeqAIJ*)C->data;
45     c->atb                 = atb;
46     atb->At                = Pt;
47     atb->destroy           = C->ops->destroy;
48     C->ops->destroy        = MatDestroy_SeqAIJ_MatTransMatMult;
49     C->ops->ptapnumeric    = MatPtAPNumeric_SeqAIJ_SeqAIJ;
50     C->ops->productnumeric = MatProductNumeric_PtAP;
51     PetscFunctionReturn(0);
52   }
53 
54   /* hypre */
55 #if defined(PETSC_HAVE_HYPRE)
56   ierr = PetscStrcmp(alg,"hypre",&flg);CHKERRQ(ierr);
57   if (flg) {
58     ierr = MatPtAPSymbolic_AIJ_AIJ_wHYPRE(A,P,fill,C);CHKERRQ(ierr);
59     PetscFunctionReturn(0);
60   }
61 #endif
62 
63   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatProductType is not supported");
64   PetscFunctionReturn(0);
65 }
66 
67 PetscErrorCode MatDestroy_SeqAIJ_PtAP(Mat A)
68 {
69   PetscErrorCode ierr;
70   Mat_SeqAIJ     *a  = (Mat_SeqAIJ*)A->data;
71   Mat_AP         *ap = a->ap;
72 
73   PetscFunctionBegin;
74   ierr = PetscFree(ap->apa);CHKERRQ(ierr);
75   ierr = PetscFree(ap->api);CHKERRQ(ierr);
76   ierr = PetscFree(ap->apj);CHKERRQ(ierr);
77   ierr = (ap->destroy)(A);CHKERRQ(ierr);
78   ierr = PetscFree(ap);CHKERRQ(ierr);
79   PetscFunctionReturn(0);
80 }
81 
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 = MatSetSeqAIJWithArrays_private(PetscObjectComm((PetscObject)A),pn,pn,ci,cj,ca,((PetscObject)A)->type_name,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 PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy(Mat A,Mat P,Mat C)
209 {
210   PetscErrorCode ierr;
211   Mat_SeqAIJ     *a = (Mat_SeqAIJ*) A->data;
212   Mat_SeqAIJ     *p = (Mat_SeqAIJ*) P->data;
213   Mat_SeqAIJ     *c = (Mat_SeqAIJ*) C->data;
214   PetscInt       *ai=a->i,*aj=a->j,*apj,*apjdense,*pi=p->i,*pj=p->j,*pJ=p->j,*pjj;
215   PetscInt       *ci=c->i,*cj=c->j,*cjj;
216   PetscInt       am =A->rmap->N,cn=C->cmap->N,cm=C->rmap->N;
217   PetscInt       i,j,k,anzi,pnzi,apnzj,nextap,pnzj,prow,crow;
218   MatScalar      *aa=a->a,*apa,*pa=p->a,*pA=p->a,*paj,*ca=c->a,*caj;
219 
220   PetscFunctionBegin;
221   /* Allocate temporary array for storage of one row of A*P (cn: non-scalable) */
222   ierr = PetscCalloc2(cn,&apa,cn,&apjdense);CHKERRQ(ierr);
223   ierr = PetscMalloc1(cn,&apj);CHKERRQ(ierr);
224 
225   /* Clear old values in C */
226   ierr = PetscArrayzero(ca,ci[cm]);CHKERRQ(ierr);
227 
228   for (i=0; i<am; i++) {
229     /* Form sparse row of A*P */
230     anzi  = ai[i+1] - ai[i];
231     apnzj = 0;
232     for (j=0; j<anzi; j++) {
233       prow = *aj++;
234       pnzj = pi[prow+1] - pi[prow];
235       pjj  = pj + pi[prow];
236       paj  = pa + pi[prow];
237       for (k=0; k<pnzj; k++) {
238         if (!apjdense[pjj[k]]) {
239           apjdense[pjj[k]] = -1;
240           apj[apnzj++]     = pjj[k];
241         }
242         apa[pjj[k]] += (*aa)*paj[k];
243       }
244       ierr = PetscLogFlops(2.0*pnzj);CHKERRQ(ierr);
245       aa++;
246     }
247 
248     /* Sort the j index array for quick sparse axpy. */
249     /* Note: a array does not need sorting as it is in dense storage locations. */
250     ierr = PetscSortInt(apnzj,apj);CHKERRQ(ierr);
251 
252     /* Compute P^T*A*P using outer product (P^T)[:,j]*(A*P)[j,:]. */
253     pnzi = pi[i+1] - pi[i];
254     for (j=0; j<pnzi; j++) {
255       nextap = 0;
256       crow   = *pJ++;
257       cjj    = cj + ci[crow];
258       caj    = ca + ci[crow];
259       /* Perform sparse axpy operation.  Note cjj includes apj. */
260       for (k=0; nextap<apnzj; k++) {
261         if (PetscUnlikelyDebug(k >= ci[crow+1] - ci[crow])) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"k too large k %d, crow %d",k,crow);
262         if (cjj[k]==apj[nextap]) {
263           caj[k] += (*pA)*apa[apj[nextap++]];
264         }
265       }
266       ierr = PetscLogFlops(2.0*apnzj);CHKERRQ(ierr);
267       pA++;
268     }
269 
270     /* Zero the current row info for A*P */
271     for (j=0; j<apnzj; j++) {
272       apa[apj[j]]      = 0.;
273       apjdense[apj[j]] = 0;
274     }
275   }
276 
277   /* Assemble the final matrix and clean up */
278   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
279   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
280 
281   ierr = PetscFree2(apa,apjdense);CHKERRQ(ierr);
282   ierr = PetscFree(apj);CHKERRQ(ierr);
283   PetscFunctionReturn(0);
284 }
285 
286 PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat A,Mat P,Mat C)
287 {
288   PetscErrorCode      ierr;
289   Mat_SeqAIJ          *c = (Mat_SeqAIJ*)C->data;
290   Mat_MatTransMatMult *atb = c->atb;
291   Mat                 Pt = atb->At;
292 
293   PetscFunctionBegin;
294   ierr = MatTranspose_SeqAIJ(P,MAT_REUSE_MATRIX,&Pt);CHKERRQ(ierr);
295   ierr = (C->ops->matmatmultnumeric)(Pt,A,P,C);CHKERRQ(ierr);
296   PetscFunctionReturn(0);
297 }
298