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