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