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