xref: /petsc/src/mat/impls/aij/seq/matptap.c (revision 2f613bf53f46f9356e00a2ca2bd69453be72fc31)
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 
26   PetscFunctionBegin;
27   /* "scalable" */
28   ierr = PetscStrcmp(alg,"scalable",&flg);CHKERRQ(ierr);
29   if (flg) {
30     ierr = MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(A,P,fill,C);CHKERRQ(ierr);
31     C->ops->productnumeric = MatProductNumeric_PtAP;
32     PetscFunctionReturn(0);
33   }
34 
35   /* "rap" */
36   ierr = PetscStrcmp(alg,"rap",&flg);CHKERRQ(ierr);
37   if (flg) {
38     Mat_MatTransMatMult *atb;
39 
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     atb->At                = Pt;
45     atb->data              = C->product->data;
46     atb->destroy           = C->product->destroy;
47     C->product->data       = atb;
48     C->product->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 }
65 
66 PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(Mat A,Mat P,PetscReal fill,Mat C)
67 {
68   PetscErrorCode     ierr;
69   PetscFreeSpaceList free_space=NULL,current_space=NULL;
70   Mat_SeqAIJ         *a        = (Mat_SeqAIJ*)A->data,*p = (Mat_SeqAIJ*)P->data,*c;
71   PetscInt           *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj;
72   PetscInt           *ci,*cj,*ptadenserow,*ptasparserow,*ptaj,nspacedouble=0;
73   PetscInt           an=A->cmap->N,am=A->rmap->N,pn=P->cmap->N,pm=P->rmap->N;
74   PetscInt           i,j,k,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi,nlnk,*lnk;
75   MatScalar          *ca;
76   PetscBT            lnkbt;
77   PetscReal          afill;
78 
79   PetscFunctionBegin;
80   /* Get ij structure of P^T */
81   ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
82   ptJ  = ptj;
83 
84   /* Allocate ci array, arrays for fill computation and */
85   /* free space for accumulating nonzero column info */
86   ierr  = PetscMalloc1(pn+1,&ci);CHKERRQ(ierr);
87   ci[0] = 0;
88 
89   ierr         = PetscCalloc1(2*an+1,&ptadenserow);CHKERRQ(ierr);
90   ptasparserow = ptadenserow  + an;
91 
92   /* create and initialize a linked list */
93   nlnk = pn+1;
94   ierr = PetscLLCreate(pn,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr);
95 
96   /* Set initial free space to be fill*(nnz(A)+ nnz(P)) */
97   ierr          = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am],pi[pm])),&free_space);CHKERRQ(ierr);
98   current_space = free_space;
99 
100   /* Determine symbolic info for each row of C: */
101   for (i=0; i<pn; i++) {
102     ptnzi  = pti[i+1] - pti[i];
103     ptanzi = 0;
104     /* Determine symbolic row of PtA: */
105     for (j=0; j<ptnzi; j++) {
106       arow = *ptJ++;
107       anzj = ai[arow+1] - ai[arow];
108       ajj  = aj + ai[arow];
109       for (k=0; k<anzj; k++) {
110         if (!ptadenserow[ajj[k]]) {
111           ptadenserow[ajj[k]]    = -1;
112           ptasparserow[ptanzi++] = ajj[k];
113         }
114       }
115     }
116     /* Using symbolic info for row of PtA, determine symbolic info for row of C: */
117     ptaj = ptasparserow;
118     cnzi = 0;
119     for (j=0; j<ptanzi; j++) {
120       prow = *ptaj++;
121       pnzj = pi[prow+1] - pi[prow];
122       pjj  = pj + pi[prow];
123       /* add non-zero cols of P into the sorted linked list lnk */
124       ierr  = PetscLLAddSorted(pnzj,pjj,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr);
125       cnzi += nlnk;
126     }
127 
128     /* If free space is not available, make more free space */
129     /* Double the amount of total space in the list */
130     if (current_space->local_remaining<cnzi) {
131       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(cnzi,current_space->total_array_size),&current_space);CHKERRQ(ierr);
132       nspacedouble++;
133     }
134 
135     /* Copy data into free space, and zero out denserows */
136     ierr = PetscLLClean(pn,pn,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
137 
138     current_space->array           += cnzi;
139     current_space->local_used      += cnzi;
140     current_space->local_remaining -= cnzi;
141 
142     for (j=0; j<ptanzi; j++) ptadenserow[ptasparserow[j]] = 0;
143 
144     /* Aside: Perhaps we should save the pta info for the numerical factorization. */
145     /*        For now, we will recompute what is needed. */
146     ci[i+1] = ci[i] + cnzi;
147   }
148   /* nnz is now stored in ci[ptm], column indices are in the list of free space */
149   /* Allocate space for cj, initialize cj, and */
150   /* destroy list of free space and other temporary array(s) */
151   ierr = PetscMalloc1(ci[pn]+1,&cj);CHKERRQ(ierr);
152   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
153   ierr = PetscFree(ptadenserow);CHKERRQ(ierr);
154   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
155 
156   ierr = PetscCalloc1(ci[pn]+1,&ca);CHKERRQ(ierr);
157 
158   /* put together the new matrix */
159   ierr = MatSetSeqAIJWithArrays_private(PetscObjectComm((PetscObject)A),pn,pn,ci,cj,ca,((PetscObject)A)->type_name,C);CHKERRQ(ierr);
160   ierr = MatSetBlockSizes(C,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr);
161 
162   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
163   /* Since these are PETSc arrays, change flags to free them as necessary. */
164   c          = (Mat_SeqAIJ*)((C)->data);
165   c->free_a  = PETSC_TRUE;
166   c->free_ij = PETSC_TRUE;
167   c->nonew   = 0;
168 
169   C->ops->ptapnumeric = MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy;
170 
171   /* set MatInfo */
172   afill = (PetscReal)ci[pn]/(ai[am]+pi[pm] + 1.e-5);
173   if (afill < 1.0) afill = 1.0;
174   c->maxnz                     = ci[pn];
175   c->nz                        = ci[pn];
176   C->info.mallocs           = nspacedouble;
177   C->info.fill_ratio_given  = fill;
178   C->info.fill_ratio_needed = afill;
179 
180   /* Clean up. */
181   ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
182 #if defined(PETSC_USE_INFO)
183   if (ci[pn] != 0) {
184     ierr = PetscInfo3(C,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
185     ierr = PetscInfo1(C,"Use MatPtAP(A,P,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr);
186   } else {
187     ierr = PetscInfo(C,"Empty matrix product\n");CHKERRQ(ierr);
188   }
189 #endif
190   PetscFunctionReturn(0);
191 }
192 
193 PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy(Mat A,Mat P,Mat C)
194 {
195   PetscErrorCode ierr;
196   Mat_SeqAIJ     *a = (Mat_SeqAIJ*) A->data;
197   Mat_SeqAIJ     *p = (Mat_SeqAIJ*) P->data;
198   Mat_SeqAIJ     *c = (Mat_SeqAIJ*) C->data;
199   PetscInt       *ai=a->i,*aj=a->j,*apj,*apjdense,*pi=p->i,*pj=p->j,*pJ=p->j,*pjj;
200   PetscInt       *ci=c->i,*cj=c->j,*cjj;
201   PetscInt       am =A->rmap->N,cn=C->cmap->N,cm=C->rmap->N;
202   PetscInt       i,j,k,anzi,pnzi,apnzj,nextap,pnzj,prow,crow;
203   MatScalar      *aa,*apa,*pa,*pA,*paj,*ca,*caj;
204 
205   PetscFunctionBegin;
206   /* Allocate temporary array for storage of one row of A*P (cn: non-scalable) */
207   ierr = PetscCalloc2(cn,&apa,cn,&apjdense);CHKERRQ(ierr);
208   ierr = PetscMalloc1(cn,&apj);CHKERRQ(ierr);
209   /* trigger CPU copies if needed and flag CPU mask for C */
210 #if defined(PETSC_HAVE_DEVICE)
211   {
212     const PetscScalar *dummy;
213     ierr = MatSeqAIJGetArrayRead(A,&dummy);CHKERRQ(ierr);
214     ierr = MatSeqAIJRestoreArrayRead(A,&dummy);CHKERRQ(ierr);
215     ierr = MatSeqAIJGetArrayRead(P,&dummy);CHKERRQ(ierr);
216     ierr = MatSeqAIJRestoreArrayRead(P,&dummy);CHKERRQ(ierr);
217     if (C->offloadmask != PETSC_OFFLOAD_UNALLOCATED) C->offloadmask = PETSC_OFFLOAD_CPU;
218   }
219 #endif
220   aa = a->a;
221   pa = p->a;
222   pA = p->a;
223   ca = c->a;
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_MatTransMatMult *atb;
290 
291   PetscFunctionBegin;
292   MatCheckProduct(C,3);
293   atb  = (Mat_MatTransMatMult*)C->product->data;
294   if (!atb) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_PLIB,"Missing data structure");
295   ierr = MatTranspose_SeqAIJ(P,MAT_REUSE_MATRIX,&atb->At);CHKERRQ(ierr);
296   if (!C->ops->matmultnumeric) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_PLIB,"Missing numeric operation");
297   /* when using rap, MatMatMatMultSymbolic used a different data */
298   if (atb->data) C->product->data = atb->data;
299   ierr = (*C->ops->matmatmultnumeric)(atb->At,A,P,C);CHKERRQ(ierr);
300   C->product->data = atb;
301   PetscFunctionReturn(0);
302 }
303