xref: /petsc/src/mat/impls/aij/seq/matptap.c (revision f3f8564e68dacfa3f5888fed8abe34f553e62bf2)
1 #define PETSCMAT_DLL
2 
3 /*
4   Defines projective product routines where A is a SeqAIJ matrix
5           C = P^T * A * P
6 */
7 
8 #include "src/mat/impls/aij/seq/aij.h"   /*I "petscmat.h" I*/
9 #include "src/mat/utils/freespace.h"
10 #include "petscbt.h"
11 
12 #undef __FUNCT__
13 #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ"
14 PetscErrorCode MatPtAPSymbolic_SeqAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
15 {
16   PetscErrorCode ierr;
17 
18   PetscFunctionBegin;
19   if (!P->ops->ptapsymbolic_seqaij) {
20     SETERRQ2(PETSC_ERR_SUP,"Not implemented for A=%s and P=%s",A->type_name,P->type_name);
21   }
22   ierr = (*P->ops->ptapsymbolic_seqaij)(A,P,fill,C);CHKERRQ(ierr);
23   PetscFunctionReturn(0);
24 }
25 
26 #undef __FUNCT__
27 #define __FUNCT__ "MatPtAPNumeric_SeqAIJ"
28 PetscErrorCode MatPtAPNumeric_SeqAIJ(Mat A,Mat P,Mat C)
29 {
30   PetscErrorCode ierr;
31 
32   PetscFunctionBegin;
33   if (!P->ops->ptapnumeric_seqaij) {
34     SETERRQ2(PETSC_ERR_SUP,"Not implemented for A=%s and P=%s",A->type_name,P->type_name);
35   }
36   ierr = (*P->ops->ptapnumeric_seqaij)(A,P,C);CHKERRQ(ierr);
37   PetscFunctionReturn(0);
38 }
39 
40 #undef __FUNCT__
41 #define __FUNCT__ "MatPtAPSymbolic_SeqAIJ_SeqAIJ"
42 PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
43 {
44   PetscErrorCode     ierr;
45   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
46   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data,*p = (Mat_SeqAIJ*)P->data,*c;
47   PetscInt           *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj;
48   PetscInt           *ci,*cj,*ptadenserow,*ptasparserow,*ptaj,nspacedouble=0;
49   PetscInt           an=A->cmap.N,am=A->rmap.N,pn=P->cmap.N,pm=P->rmap.N;
50   PetscInt           i,j,k,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi,nlnk,*lnk;
51   MatScalar          *ca;
52   PetscBT            lnkbt;
53 
54   PetscFunctionBegin;
55   /* Get ij structure of P^T */
56   ierr = MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
57   ptJ=ptj;
58 
59   /* Allocate ci array, arrays for fill computation and */
60   /* free space for accumulating nonzero column info */
61   ierr = PetscMalloc((pn+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr);
62   ci[0] = 0;
63 
64   ierr = PetscMalloc((2*an+1)*sizeof(PetscInt),&ptadenserow);CHKERRQ(ierr);
65   ierr = PetscMemzero(ptadenserow,(2*an+1)*sizeof(PetscInt));CHKERRQ(ierr);
66   ptasparserow = ptadenserow  + an;
67 
68   /* create and initialize a linked list */
69   nlnk = pn+1;
70   ierr = PetscLLCreate(pn,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr);
71 
72   /* Set initial free space to be fill*nnz(A). */
73   /* This should be reasonable if sparsity of PtAP is similar to that of A. */
74   ierr          = PetscFreeSpaceGet((PetscInt)(fill*ai[am]),&free_space);
75   current_space = free_space;
76 
77   /* Determine symbolic info for each row of C: */
78   for (i=0;i<pn;i++) {
79     ptnzi  = pti[i+1] - pti[i];
80     ptanzi = 0;
81     /* Determine symbolic row of PtA: */
82     for (j=0;j<ptnzi;j++) {
83       arow = *ptJ++;
84       anzj = ai[arow+1] - ai[arow];
85       ajj  = aj + ai[arow];
86       for (k=0;k<anzj;k++) {
87         if (!ptadenserow[ajj[k]]) {
88           ptadenserow[ajj[k]]    = -1;
89           ptasparserow[ptanzi++] = ajj[k];
90         }
91       }
92     }
93     /* Using symbolic info for row of PtA, determine symbolic info for row of C: */
94     ptaj = ptasparserow;
95     cnzi   = 0;
96     for (j=0;j<ptanzi;j++) {
97       prow = *ptaj++;
98       pnzj = pi[prow+1] - pi[prow];
99       pjj  = pj + pi[prow];
100       /* add non-zero cols of P into the sorted linked list lnk */
101       ierr = PetscLLAdd(pnzj,pjj,pn,nlnk,lnk,lnkbt);CHKERRQ(ierr);
102       cnzi += nlnk;
103     }
104 
105     /* If free space is not available, make more free space */
106     /* Double the amount of total space in the list */
107     if (current_space->local_remaining<cnzi) {
108       ierr = PetscFreeSpaceGet(current_space->total_array_size,&current_space);CHKERRQ(ierr);
109       nspacedouble++;
110     }
111 
112     /* Copy data into free space, and zero out denserows */
113     ierr = PetscLLClean(pn,pn,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
114     current_space->array           += cnzi;
115     current_space->local_used      += cnzi;
116     current_space->local_remaining -= cnzi;
117 
118     for (j=0;j<ptanzi;j++) {
119       ptadenserow[ptasparserow[j]] = 0;
120     }
121     /* Aside: Perhaps we should save the pta info for the numerical factorization. */
122     /*        For now, we will recompute what is needed. */
123     ci[i+1] = ci[i] + cnzi;
124   }
125   /* nnz is now stored in ci[ptm], column indices are in the list of free space */
126   /* Allocate space for cj, initialize cj, and */
127   /* destroy list of free space and other temporary array(s) */
128   ierr = PetscMalloc((ci[pn]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr);
129   ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr);
130   ierr = PetscFree(ptadenserow);CHKERRQ(ierr);
131   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
132 
133   /* Allocate space for ca */
134   ierr = PetscMalloc((ci[pn]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr);
135   ierr = PetscMemzero(ca,(ci[pn]+1)*sizeof(MatScalar));CHKERRQ(ierr);
136 
137   /* put together the new matrix */
138   ierr = MatCreateSeqAIJWithArrays(A->comm,pn,pn,ci,cj,ca,C);CHKERRQ(ierr);
139 
140   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
141   /* Since these are PETSc arrays, change flags to free them as necessary. */
142   c = (Mat_SeqAIJ *)((*C)->data);
143   c->free_a  = PETSC_TRUE;
144   c->free_ij = PETSC_TRUE;
145   c->nonew   = 0;
146 
147   /* Clean up. */
148   ierr = MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);CHKERRQ(ierr);
149 #if defined(PETSC_USE_INFO)
150   if (ci[pn] != 0) {
151     PetscReal afill = ((PetscReal)ci[pn])/ai[am];
152     if (afill < 1.0) afill = 1.0;
153     ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
154     ierr = PetscInfo1((*C),"Use MatPtAP(A,P,MatReuse,%G,&C) for best performance.\n",afill);CHKERRQ(ierr);
155   } else {
156     ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr);
157   }
158 #endif
159   PetscFunctionReturn(0);
160 }
161 
162 #undef __FUNCT__
163 #define __FUNCT__ "MatPtAPNumeric_SeqAIJ_SeqAIJ"
164 PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat A,Mat P,Mat C)
165 {
166   PetscErrorCode ierr;
167   PetscInt       flops=0;
168   Mat_SeqAIJ     *a  = (Mat_SeqAIJ *) A->data;
169   Mat_SeqAIJ     *p  = (Mat_SeqAIJ *) P->data;
170   Mat_SeqAIJ     *c  = (Mat_SeqAIJ *) C->data;
171   PetscInt       *ai=a->i,*aj=a->j,*apj,*apjdense,*pi=p->i,*pj=p->j,*pJ=p->j,*pjj;
172   PetscInt       *ci=c->i,*cj=c->j,*cjj;
173   PetscInt       am=A->rmap.N,cn=C->cmap.N,cm=C->rmap.N;
174   PetscInt       i,j,k,anzi,pnzi,apnzj,nextap,pnzj,prow,crow;
175   MatScalar      *aa=a->a,*apa,*pa=p->a,*pA=p->a,*paj,*ca=c->a,*caj;
176 
177   PetscFunctionBegin;
178   /* Allocate temporary array for storage of one row of A*P */
179   ierr = PetscMalloc(cn*(sizeof(MatScalar)+2*sizeof(PetscInt)),&apa);CHKERRQ(ierr);
180   ierr = PetscMemzero(apa,cn*(sizeof(MatScalar)+2*sizeof(PetscInt)));CHKERRQ(ierr);
181 
182   apj      = (PetscInt *)(apa + cn);
183   apjdense = apj + cn;
184 
185   /* Clear old values in C */
186   ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr);
187 
188   for (i=0;i<am;i++) {
189     /* Form sparse row of A*P */
190     anzi  = ai[i+1] - ai[i];
191     apnzj = 0;
192     for (j=0;j<anzi;j++) {
193       prow = *aj++;
194       pnzj = pi[prow+1] - pi[prow];
195       pjj  = pj + pi[prow];
196       paj  = pa + pi[prow];
197       for (k=0;k<pnzj;k++) {
198         if (!apjdense[pjj[k]]) {
199           apjdense[pjj[k]] = -1;
200           apj[apnzj++]     = pjj[k];
201         }
202         apa[pjj[k]] += (*aa)*paj[k];
203       }
204       flops += 2*pnzj;
205       aa++;
206     }
207 
208     /* Sort the j index array for quick sparse axpy. */
209     /* Note: a array does not need sorting as it is in dense storage locations. */
210     ierr = PetscSortInt(apnzj,apj);CHKERRQ(ierr);
211 
212     /* Compute P^T*A*P using outer product (P^T)[:,j]*(A*P)[j,:]. */
213     pnzi = pi[i+1] - pi[i];
214     for (j=0;j<pnzi;j++) {
215       nextap = 0;
216       crow   = *pJ++;
217       cjj    = cj + ci[crow];
218       caj    = ca + ci[crow];
219       /* Perform sparse axpy operation.  Note cjj includes apj. */
220       for (k=0;nextap<apnzj;k++) {
221 #if defined(PETSC_USE_DEBUG)
222         if (k >= ci[crow+1] - ci[crow]) {
223           SETERRQ2(PETSC_ERR_PLIB,"k too large k %d, crow %d",k,crow);
224         }
225 #endif
226         if (cjj[k]==apj[nextap]) {
227           caj[k] += (*pA)*apa[apj[nextap++]];
228         }
229       }
230       flops += 2*apnzj;
231       pA++;
232     }
233 
234     /* Zero the current row info for A*P */
235     for (j=0;j<apnzj;j++) {
236       apa[apj[j]]      = 0.;
237       apjdense[apj[j]] = 0;
238     }
239   }
240 
241   /* Assemble the final matrix and clean up */
242   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
243   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
244   ierr = PetscFree(apa);CHKERRQ(ierr);
245   ierr = PetscLogFlops(flops);CHKERRQ(ierr);
246 
247   PetscFunctionReturn(0);
248 }
249