xref: /petsc/src/mat/impls/aij/seq/matptap.c (revision 715a534632cc6029c3649e8dc2fa74cc65e81325)
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 PETSC_INTERN PetscErrorCode MatPtAP_SeqAIJ_SeqAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
17 {
18   PetscErrorCode      ierr;
19 #if !defined(PETSC_HAVE_HYPRE)
20   const char          *algTypes[2] = {"scalable","rap"};
21   PetscInt            nalg = 2;
22 #else
23   const char          *algTypes[3] = {"scalable","rap","hypre"};
24   PetscInt            nalg = 3;
25 #endif
26   PetscInt            alg = 1; /* set default algorithm */
27   Mat                 Pt;
28   Mat_MatTransMatMult *atb;
29   Mat_SeqAIJ          *c;
30 
31   PetscFunctionBegin;
32   if (scall == MAT_INITIAL_MATRIX) {
33     /*
34      Alg 'scalable' determines which implementations to be used:
35        "rap":      Pt = P^T and C = Pt*A*P
36        "scalable": do outer product and two sparse axpy in MatPtAPNumeric() - might slow, does not store structure of A*P.
37        "hypre":    use boomerAMGBuildCoarseOperator.
38      */
39     ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"MatPtAP","Mat");CHKERRQ(ierr);
40     ierr = PetscOptionsEList("-matptap_via","Algorithmic approach","MatPtAP",algTypes,nalg,algTypes[0],&alg,NULL);CHKERRQ(ierr);
41     ierr = PetscOptionsEnd();CHKERRQ(ierr);
42     switch (alg) {
43     case 1:
44       ierr = PetscNew(&atb);CHKERRQ(ierr);
45       ierr = MatTranspose_SeqAIJ(P,MAT_INITIAL_MATRIX,&Pt);CHKERRQ(ierr);
46       ierr = MatMatMatMult(Pt,A,P,MAT_INITIAL_MATRIX,fill,C);CHKERRQ(ierr);
47 
48       c                      = (Mat_SeqAIJ*)(*C)->data;
49       c->atb                 = atb;
50       atb->At                = Pt;
51       atb->destroy           = (*C)->ops->destroy;
52       (*C)->ops->destroy     = MatDestroy_SeqAIJ_MatTransMatMult;
53       (*C)->ops->ptapnumeric = MatPtAPNumeric_SeqAIJ_SeqAIJ;
54       PetscFunctionReturn(0);
55       break;
56 #if defined(PETSC_HAVE_HYPRE)
57     case 2:
58       ierr = MatPtAPSymbolic_AIJ_AIJ_wHYPRE(A,P,fill,C);CHKERRQ(ierr);
59       break;
60 #endif
61     default:
62       ierr = MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(A,P,fill,C);CHKERRQ(ierr);
63       break;
64     }
65   }
66   ierr = (*(*C)->ops->ptapnumeric)(A,P,*C);CHKERRQ(ierr);
67   PetscFunctionReturn(0);
68 }
69 
70 PetscErrorCode MatDestroy_SeqAIJ_PtAP(Mat A)
71 {
72   PetscErrorCode ierr;
73   Mat_SeqAIJ     *a  = (Mat_SeqAIJ*)A->data;
74   Mat_AP         *ap = a->ap;
75 
76   PetscFunctionBegin;
77   ierr = PetscFree(ap->apa);CHKERRQ(ierr);
78   ierr = PetscFree(ap->api);CHKERRQ(ierr);
79   ierr = PetscFree(ap->apj);CHKERRQ(ierr);
80   ierr = (ap->destroy)(A);CHKERRQ(ierr);
81   ierr = PetscFree(ap);CHKERRQ(ierr);
82   PetscFunctionReturn(0);
83 }
84 
85 PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(Mat A,Mat P,PetscReal fill,Mat *C)
86 {
87   PetscErrorCode     ierr;
88   PetscFreeSpaceList free_space=NULL,current_space=NULL;
89   Mat_SeqAIJ         *a        = (Mat_SeqAIJ*)A->data,*p = (Mat_SeqAIJ*)P->data,*c;
90   PetscInt           *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj;
91   PetscInt           *ci,*cj,*ptadenserow,*ptasparserow,*ptaj,nspacedouble=0;
92   PetscInt           an=A->cmap->N,am=A->rmap->N,pn=P->cmap->N,pm=P->rmap->N;
93   PetscInt           i,j,k,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi,nlnk,*lnk;
94   MatScalar          *ca;
95   PetscBT            lnkbt;
96   PetscReal          afill;
97 
98   PetscFunctionBegin;
99   /* Get ij structure of P^T */
100   ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
101   ptJ  = ptj;
102 
103   /* Allocate ci array, arrays for fill computation and */
104   /* free space for accumulating nonzero column info */
105   ierr  = PetscMalloc1(pn+1,&ci);CHKERRQ(ierr);
106   ci[0] = 0;
107 
108   ierr         = PetscCalloc1(2*an+1,&ptadenserow);CHKERRQ(ierr);
109   ptasparserow = ptadenserow  + an;
110 
111   /* create and initialize a linked list */
112   nlnk = pn+1;
113   ierr = PetscLLCreate(pn,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr);
114 
115   /* Set initial free space to be fill*(nnz(A)+ nnz(P)) */
116   ierr          = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am],pi[pm])),&free_space);CHKERRQ(ierr);
117   current_space = free_space;
118 
119   /* Determine symbolic info for each row of C: */
120   for (i=0; i<pn; i++) {
121     ptnzi  = pti[i+1] - pti[i];
122     ptanzi = 0;
123     /* Determine symbolic row of PtA: */
124     for (j=0; j<ptnzi; j++) {
125       arow = *ptJ++;
126       anzj = ai[arow+1] - ai[arow];
127       ajj  = aj + ai[arow];
128       for (k=0; k<anzj; k++) {
129         if (!ptadenserow[ajj[k]]) {
130           ptadenserow[ajj[k]]    = -1;
131           ptasparserow[ptanzi++] = ajj[k];
132         }
133       }
134     }
135     /* Using symbolic info for row of PtA, determine symbolic info for row of C: */
136     ptaj = ptasparserow;
137     cnzi = 0;
138     for (j=0; j<ptanzi; j++) {
139       prow = *ptaj++;
140       pnzj = pi[prow+1] - pi[prow];
141       pjj  = pj + pi[prow];
142       /* add non-zero cols of P into the sorted linked list lnk */
143       ierr  = PetscLLAddSorted(pnzj,pjj,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr);
144       cnzi += nlnk;
145     }
146 
147     /* If free space is not available, make more free space */
148     /* Double the amount of total space in the list */
149     if (current_space->local_remaining<cnzi) {
150       ierr = PetscFreeSpaceGet(PetscIntSumTruncate(cnzi,current_space->total_array_size),&current_space);CHKERRQ(ierr);
151       nspacedouble++;
152     }
153 
154     /* Copy data into free space, and zero out denserows */
155     ierr = PetscLLClean(pn,pn,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
156 
157     current_space->array           += cnzi;
158     current_space->local_used      += cnzi;
159     current_space->local_remaining -= cnzi;
160 
161     for (j=0; j<ptanzi; j++) ptadenserow[ptasparserow[j]] = 0;
162 
163     /* Aside: Perhaps we should save the pta info for the numerical factorization. */
164     /*        For now, we will recompute what is needed. */
165     ci[i+1] = ci[i] + cnzi;
166   }
167   /* nnz is now stored in ci[ptm], column indices are in the list of free space */
168   /* Allocate space for cj, initialize cj, and */
169   /* destroy list of free space and other temporary array(s) */
170   ierr = PetscMalloc1(ci[pn]+1,&cj);CHKERRQ(ierr);
171   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
172   ierr = PetscFree(ptadenserow);CHKERRQ(ierr);
173   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
174 
175   ierr = PetscCalloc1(ci[pn]+1,&ca);CHKERRQ(ierr);
176 
177   /* put together the new matrix */
178   ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),pn,pn,ci,cj,ca,C);CHKERRQ(ierr);
179   ierr = MatSetBlockSizes(*C,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr);
180 
181   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
182   /* Since these are PETSc arrays, change flags to free them as necessary. */
183   c          = (Mat_SeqAIJ*)((*C)->data);
184   c->free_a  = PETSC_TRUE;
185   c->free_ij = PETSC_TRUE;
186   c->nonew   = 0;
187   (*C)->ops->ptapnumeric = MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy;
188 
189   /* set MatInfo */
190   afill = (PetscReal)ci[pn]/(ai[am]+pi[pm] + 1.e-5);
191   if (afill < 1.0) afill = 1.0;
192   c->maxnz                     = ci[pn];
193   c->nz                        = ci[pn];
194   (*C)->info.mallocs           = nspacedouble;
195   (*C)->info.fill_ratio_given  = fill;
196   (*C)->info.fill_ratio_needed = afill;
197 
198   /* Clean up. */
199   ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
200 #if defined(PETSC_USE_INFO)
201   if (ci[pn] != 0) {
202     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr);
203     ierr = PetscInfo1((*C),"Use MatPtAP(A,P,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr);
204   } else {
205     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
206   }
207 #endif
208   PetscFunctionReturn(0);
209 }
210 
211 PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy(Mat A,Mat P,Mat C)
212 {
213   PetscErrorCode ierr;
214   Mat_SeqAIJ     *a = (Mat_SeqAIJ*) A->data;
215   Mat_SeqAIJ     *p = (Mat_SeqAIJ*) P->data;
216   Mat_SeqAIJ     *c = (Mat_SeqAIJ*) C->data;
217   PetscInt       *ai=a->i,*aj=a->j,*apj,*apjdense,*pi=p->i,*pj=p->j,*pJ=p->j,*pjj;
218   PetscInt       *ci=c->i,*cj=c->j,*cjj;
219   PetscInt       am =A->rmap->N,cn=C->cmap->N,cm=C->rmap->N;
220   PetscInt       i,j,k,anzi,pnzi,apnzj,nextap,pnzj,prow,crow;
221   MatScalar      *aa=a->a,*apa,*pa=p->a,*pA=p->a,*paj,*ca=c->a,*caj;
222 
223   PetscFunctionBegin;
224   /* Allocate temporary array for storage of one row of A*P (cn: non-scalable) */
225   ierr = PetscMalloc3(cn,&apa,cn,&apjdense,cn,&apj);CHKERRQ(ierr);
226   ierr = PetscMemzero(apa,cn*sizeof(MatScalar));CHKERRQ(ierr);
227   ierr = PetscMemzero(apjdense,cn*sizeof(PetscInt));CHKERRQ(ierr);
228 
229   /* Clear old values in C */
230   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
231 
232   for (i=0; i<am; i++) {
233     /* Form sparse row of A*P */
234     anzi  = ai[i+1] - ai[i];
235     apnzj = 0;
236     for (j=0; j<anzi; j++) {
237       prow = *aj++;
238       pnzj = pi[prow+1] - pi[prow];
239       pjj  = pj + pi[prow];
240       paj  = pa + pi[prow];
241       for (k=0; k<pnzj; k++) {
242         if (!apjdense[pjj[k]]) {
243           apjdense[pjj[k]] = -1;
244           apj[apnzj++]     = pjj[k];
245         }
246         apa[pjj[k]] += (*aa)*paj[k];
247       }
248       ierr = PetscLogFlops(2.0*pnzj);CHKERRQ(ierr);
249       aa++;
250     }
251 
252     /* Sort the j index array for quick sparse axpy. */
253     /* Note: a array does not need sorting as it is in dense storage locations. */
254     ierr = PetscSortInt(apnzj,apj);CHKERRQ(ierr);
255 
256     /* Compute P^T*A*P using outer product (P^T)[:,j]*(A*P)[j,:]. */
257     pnzi = pi[i+1] - pi[i];
258     for (j=0; j<pnzi; j++) {
259       nextap = 0;
260       crow   = *pJ++;
261       cjj    = cj + ci[crow];
262       caj    = ca + ci[crow];
263       /* Perform sparse axpy operation.  Note cjj includes apj. */
264       for (k=0; nextap<apnzj; k++) {
265 #if defined(PETSC_USE_DEBUG)
266         if (k >= ci[crow+1] - ci[crow]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"k too large k %d, crow %d",k,crow);
267 #endif
268         if (cjj[k]==apj[nextap]) {
269           caj[k] += (*pA)*apa[apj[nextap++]];
270         }
271       }
272       ierr = PetscLogFlops(2.0*apnzj);CHKERRQ(ierr);
273       pA++;
274     }
275 
276     /* Zero the current row info for A*P */
277     for (j=0; j<apnzj; j++) {
278       apa[apj[j]]      = 0.;
279       apjdense[apj[j]] = 0;
280     }
281   }
282 
283   /* Assemble the final matrix and clean up */
284   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
285   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
286 
287   ierr = PetscFree3(apa,apjdense,apj);CHKERRQ(ierr);
288   PetscFunctionReturn(0);
289 }
290 
291 PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat A,Mat P,Mat C)
292 {
293   PetscErrorCode      ierr;
294   Mat_SeqAIJ          *c = (Mat_SeqAIJ*)C->data;
295   Mat_MatTransMatMult *atb = c->atb;
296   Mat                 Pt = atb->At;
297 
298   PetscFunctionBegin;
299   ierr = MatTranspose_SeqAIJ(P,MAT_REUSE_MATRIX,&Pt);CHKERRQ(ierr);
300   ierr = (C->ops->matmatmultnumeric)(Pt,A,P,C);CHKERRQ(ierr);
301   PetscFunctionReturn(0);
302 }
303