Actual source code: matptap.c
1: /*
2: Defines projective product routines where A is a SeqAIJ matrix
3: C = P^T * A * P
4: */
6: #include src/mat/impls/aij/seq/aij.h
7: #include src/mat/utils/freespace.h
8: #include petscbt.h
13: PetscErrorCode MatPtAP_SeqAIJ_SeqAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
14: {
18: if (scall == MAT_INITIAL_MATRIX){
19: PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);
20: MatPtAPSymbolic_SeqAIJ_SeqAIJ(A,P,fill,C);
21: PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);
22: }
23: PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);
24: MatPtAPNumeric_SeqAIJ_SeqAIJ(A,P,*C);
25: PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);
26: return(0);
27: }
31: PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
32: {
34: FreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
35: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*p = (Mat_SeqAIJ*)P->data,*c;
36: PetscInt *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj;
37: PetscInt *ci,*cj,*ptadenserow,*ptasparserow,*ptaj;
38: PetscInt an=A->N,am=A->M,pn=P->N,pm=P->M;
39: PetscInt i,j,k,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi,nlnk,*lnk;
40: MatScalar *ca;
41: PetscBT lnkbt;
44: /* Get ij structure of P^T */
45: MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);
46: ptJ=ptj;
48: /* Allocate ci array, arrays for fill computation and */
49: /* free space for accumulating nonzero column info */
50: PetscMalloc((pn+1)*sizeof(PetscInt),&ci);
51: ci[0] = 0;
53: PetscMalloc((2*an+1)*sizeof(PetscInt),&ptadenserow);
54: PetscMemzero(ptadenserow,(2*an+1)*sizeof(PetscInt));
55: ptasparserow = ptadenserow + an;
57: /* create and initialize a linked list */
58: nlnk = pn+1;
59: PetscLLCreate(pn,pn,nlnk,lnk,lnkbt);
61: /* Set initial free space to be nnz(A) scaled by aspect ratio of P. */
62: /* This should be reasonable if sparsity of PtAP is similar to that of A. */
63: GetMoreSpace((ai[am]/pm)*pn,&free_space);
64: current_space = free_space;
66: /* Determine symbolic info for each row of C: */
67: for (i=0;i<pn;i++) {
68: ptnzi = pti[i+1] - pti[i];
69: ptanzi = 0;
70: /* Determine symbolic row of PtA: */
71: for (j=0;j<ptnzi;j++) {
72: arow = *ptJ++;
73: anzj = ai[arow+1] - ai[arow];
74: ajj = aj + ai[arow];
75: for (k=0;k<anzj;k++) {
76: if (!ptadenserow[ajj[k]]) {
77: ptadenserow[ajj[k]] = -1;
78: ptasparserow[ptanzi++] = ajj[k];
79: }
80: }
81: }
82: /* Using symbolic info for row of PtA, determine symbolic info for row of C: */
83: ptaj = ptasparserow;
84: cnzi = 0;
85: for (j=0;j<ptanzi;j++) {
86: prow = *ptaj++;
87: pnzj = pi[prow+1] - pi[prow];
88: pjj = pj + pi[prow];
89: /* add non-zero cols of P into the sorted linked list lnk */
90: PetscLLAdd(pnzj,pjj,pn,nlnk,lnk,lnkbt);
91: cnzi += nlnk;
92: }
93:
94: /* If free space is not available, make more free space */
95: /* Double the amount of total space in the list */
96: if (current_space->local_remaining<cnzi) {
97: GetMoreSpace(current_space->total_array_size,¤t_space);
98: }
100: /* Copy data into free space, and zero out denserows */
101: PetscLLClean(pn,pn,cnzi,lnk,current_space->array,lnkbt);
102: current_space->array += cnzi;
103: current_space->local_used += cnzi;
104: current_space->local_remaining -= cnzi;
105:
106: for (j=0;j<ptanzi;j++) {
107: ptadenserow[ptasparserow[j]] = 0;
108: }
109: /* Aside: Perhaps we should save the pta info for the numerical factorization. */
110: /* For now, we will recompute what is needed. */
111: ci[i+1] = ci[i] + cnzi;
112: }
113: /* nnz is now stored in ci[ptm], column indices are in the list of free space */
114: /* Allocate space for cj, initialize cj, and */
115: /* destroy list of free space and other temporary array(s) */
116: PetscMalloc((ci[pn]+1)*sizeof(PetscInt),&cj);
117: MakeSpaceContiguous(&free_space,cj);
118: PetscFree(ptadenserow);
119: PetscLLDestroy(lnk,lnkbt);
120:
121: /* Allocate space for ca */
122: PetscMalloc((ci[pn]+1)*sizeof(MatScalar),&ca);
123: PetscMemzero(ca,(ci[pn]+1)*sizeof(MatScalar));
124:
125: /* put together the new matrix */
126: MatCreateSeqAIJWithArrays(A->comm,pn,pn,ci,cj,ca,C);
128: /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
129: /* Since these are PETSc arrays, change flags to free them as necessary. */
130: c = (Mat_SeqAIJ *)((*C)->data);
131: c->freedata = PETSC_TRUE;
132: c->nonew = 0;
134: /* Clean up. */
135: MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);
137: return(0);
138: }
140: #include src/mat/impls/maij/maij.h
144: PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqMAIJ(Mat A,Mat PP,Mat *C)
145: {
146: /* This routine requires testing -- I don't think it works. */
148: FreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
149: Mat_SeqMAIJ *pp=(Mat_SeqMAIJ*)PP->data;
150: Mat P=pp->AIJ;
151: Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*p=(Mat_SeqAIJ*)P->data,*c;
152: PetscInt *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj;
153: PetscInt *ci,*cj,*denserow,*sparserow,*ptadenserow,*ptasparserow,*ptaj;
154: PetscInt an=A->N,am=A->M,pn=P->N,pm=P->M,ppdof=pp->dof;
155: PetscInt i,j,k,dof,pdof,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi;
156: MatScalar *ca;
159: /* Start timer */
160: PetscLogEventBegin(MAT_PtAPSymbolic,A,PP,0,0);
162: /* Get ij structure of P^T */
163: MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);
165: /* Allocate ci array, arrays for fill computation and */
166: /* free space for accumulating nonzero column info */
167: PetscMalloc((pn+1)*sizeof(PetscInt),&ci);
168: ci[0] = 0;
170: PetscMalloc((2*pn+2*an+1)*sizeof(PetscInt),&ptadenserow);
171: PetscMemzero(ptadenserow,(2*pn+2*an+1)*sizeof(PetscInt));
172: ptasparserow = ptadenserow + an;
173: denserow = ptasparserow + an;
174: sparserow = denserow + pn;
176: /* Set initial free space to be nnz(A) scaled by aspect ratio of P. */
177: /* This should be reasonable if sparsity of PtAP is similar to that of A. */
178: GetMoreSpace((ai[am]/pm)*pn,&free_space);
179: current_space = free_space;
181: /* Determine symbolic info for each row of C: */
182: for (i=0;i<pn/ppdof;i++) {
183: ptnzi = pti[i+1] - pti[i];
184: ptanzi = 0;
185: ptJ = ptj + pti[i];
186: for (dof=0;dof<ppdof;dof++) {
187: /* Determine symbolic row of PtA: */
188: for (j=0;j<ptnzi;j++) {
189: arow = ptJ[j] + dof;
190: anzj = ai[arow+1] - ai[arow];
191: ajj = aj + ai[arow];
192: for (k=0;k<anzj;k++) {
193: if (!ptadenserow[ajj[k]]) {
194: ptadenserow[ajj[k]] = -1;
195: ptasparserow[ptanzi++] = ajj[k];
196: }
197: }
198: }
199: /* Using symbolic info for row of PtA, determine symbolic info for row of C: */
200: ptaj = ptasparserow;
201: cnzi = 0;
202: for (j=0;j<ptanzi;j++) {
203: pdof = *ptaj%dof;
204: prow = (*ptaj++)/dof;
205: pnzj = pi[prow+1] - pi[prow];
206: pjj = pj + pi[prow];
207: for (k=0;k<pnzj;k++) {
208: if (!denserow[pjj[k]+pdof]) {
209: denserow[pjj[k]+pdof] = -1;
210: sparserow[cnzi++] = pjj[k]+pdof;
211: }
212: }
213: }
215: /* sort sparserow */
216: PetscSortInt(cnzi,sparserow);
217:
218: /* If free space is not available, make more free space */
219: /* Double the amount of total space in the list */
220: if (current_space->local_remaining<cnzi) {
221: GetMoreSpace(current_space->total_array_size,¤t_space);
222: }
224: /* Copy data into free space, and zero out denserows */
225: PetscMemcpy(current_space->array,sparserow,cnzi*sizeof(PetscInt));
226: current_space->array += cnzi;
227: current_space->local_used += cnzi;
228: current_space->local_remaining -= cnzi;
230: for (j=0;j<ptanzi;j++) {
231: ptadenserow[ptasparserow[j]] = 0;
232: }
233: for (j=0;j<cnzi;j++) {
234: denserow[sparserow[j]] = 0;
235: }
236: /* Aside: Perhaps we should save the pta info for the numerical factorization. */
237: /* For now, we will recompute what is needed. */
238: ci[i+1+dof] = ci[i+dof] + cnzi;
239: }
240: }
241: /* nnz is now stored in ci[ptm], column indices are in the list of free space */
242: /* Allocate space for cj, initialize cj, and */
243: /* destroy list of free space and other temporary array(s) */
244: PetscMalloc((ci[pn]+1)*sizeof(PetscInt),&cj);
245: MakeSpaceContiguous(&free_space,cj);
246: PetscFree(ptadenserow);
247:
248: /* Allocate space for ca */
249: PetscMalloc((ci[pn]+1)*sizeof(MatScalar),&ca);
250: PetscMemzero(ca,(ci[pn]+1)*sizeof(MatScalar));
251:
252: /* put together the new matrix */
253: MatCreateSeqAIJWithArrays(A->comm,pn,pn,ci,cj,ca,C);
255: /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
256: /* Since these are PETSc arrays, change flags to free them as necessary. */
257: c = (Mat_SeqAIJ *)((*C)->data);
258: c->freedata = PETSC_TRUE;
259: c->nonew = 0;
261: /* Clean up. */
262: MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);
264: PetscLogEventEnd(MAT_PtAPSymbolic,A,PP,0,0);
265: return(0);
266: }
272: PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat A,Mat P,Mat C)
273: {
275: PetscInt flops=0;
276: Mat_SeqAIJ *a = (Mat_SeqAIJ *) A->data;
277: Mat_SeqAIJ *p = (Mat_SeqAIJ *) P->data;
278: Mat_SeqAIJ *c = (Mat_SeqAIJ *) C->data;
279: PetscInt *ai=a->i,*aj=a->j,*apj,*apjdense,*pi=p->i,*pj=p->j,*pJ=p->j,*pjj;
280: PetscInt *ci=c->i,*cj=c->j,*cjj;
281: PetscInt am=A->M,cn=C->N,cm=C->M;
282: PetscInt i,j,k,anzi,pnzi,apnzj,nextap,pnzj,prow,crow;
283: MatScalar *aa=a->a,*apa,*pa=p->a,*pA=p->a,*paj,*ca=c->a,*caj;
286: /* Allocate temporary array for storage of one row of A*P */
287: PetscMalloc(cn*(sizeof(MatScalar)+2*sizeof(PetscInt)),&apa);
288: PetscMemzero(apa,cn*(sizeof(MatScalar)+2*sizeof(PetscInt)));
290: apj = (PetscInt *)(apa + cn);
291: apjdense = apj + cn;
293: /* Clear old values in C */
294: PetscMemzero(ca,ci[cm]*sizeof(MatScalar));
296: for (i=0;i<am;i++) {
297: /* Form sparse row of A*P */
298: anzi = ai[i+1] - ai[i];
299: apnzj = 0;
300: for (j=0;j<anzi;j++) {
301: prow = *aj++;
302: pnzj = pi[prow+1] - pi[prow];
303: pjj = pj + pi[prow];
304: paj = pa + pi[prow];
305: for (k=0;k<pnzj;k++) {
306: if (!apjdense[pjj[k]]) {
307: apjdense[pjj[k]] = -1;
308: apj[apnzj++] = pjj[k];
309: }
310: apa[pjj[k]] += (*aa)*paj[k];
311: }
312: flops += 2*pnzj;
313: aa++;
314: }
316: /* Sort the j index array for quick sparse axpy. */
317: PetscSortInt(apnzj,apj);
319: /* Compute P^T*A*P using outer product (P^T)[:,j]*(A*P)[j,:]. */
320: pnzi = pi[i+1] - pi[i];
321: for (j=0;j<pnzi;j++) {
322: nextap = 0;
323: crow = *pJ++;
324: cjj = cj + ci[crow];
325: caj = ca + ci[crow];
326: /* Perform sparse axpy operation. Note cjj includes apj. */
327: for (k=0;nextap<apnzj;k++) {
328: if (cjj[k]==apj[nextap]) {
329: caj[k] += (*pA)*apa[apj[nextap++]];
330: }
331: }
332: flops += 2*apnzj;
333: pA++;
334: }
336: /* Zero the current row info for A*P */
337: for (j=0;j<apnzj;j++) {
338: apa[apj[j]] = 0.;
339: apjdense[apj[j]] = 0;
340: }
341: }
343: /* Assemble the final matrix and clean up */
344: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
345: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
346: PetscFree(apa);
347: PetscLogFlops(flops);
349: return(0);
350: }