Actual source code: matptap.c
petsc-3.12.5 2020-03-29
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: PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"MatPtAP","Mat");
40: PetscOptionsEList("-matptap_via","Algorithmic approach","MatPtAP",algTypes,nalg,algTypes[0],&alg,NULL);
41: PetscOptionsEnd();
42: switch (alg) {
43: case 1:
44: PetscNew(&atb);
45: MatTranspose_SeqAIJ(P,MAT_INITIAL_MATRIX,&Pt);
46: MatMatMatMult(Pt,A,P,MAT_INITIAL_MATRIX,fill,C);
48: c = (Mat_SeqAIJ*)(*C)->data;
49: c->atb = atb;
50: atb->At = Pt;
51: atb->destroy = (*C)->ops->destroy;
52: (*C)->ops->destroy = MatDestroy_SeqAIJ_MatTransMatMult;
53: (*C)->ops->ptapnumeric = MatPtAPNumeric_SeqAIJ_SeqAIJ;
54: return(0);
55: break;
56: #if defined(PETSC_HAVE_HYPRE)
57: case 2:
58: PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);
59: MatPtAPSymbolic_AIJ_AIJ_wHYPRE(A,P,fill,C);
60: PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);
61: break;
62: #endif
63: default:
64: PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);
65: MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(A,P,fill,C);
66: PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);
67: break;
68: }
69: }
70: PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);
71: (*(*C)->ops->ptapnumeric)(A,P,*C);
72: PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);
73: return(0);
74: }
76: PetscErrorCode MatDestroy_SeqAIJ_PtAP(Mat A)
77: {
79: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data;
80: Mat_AP *ap = a->ap;
83: PetscFree(ap->apa);
84: PetscFree(ap->api);
85: PetscFree(ap->apj);
86: (ap->destroy)(A);
87: PetscFree(ap);
88: return(0);
89: }
91: PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(Mat A,Mat P,PetscReal fill,Mat *C)
92: {
93: PetscErrorCode ierr;
94: PetscFreeSpaceList free_space=NULL,current_space=NULL;
95: Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*p = (Mat_SeqAIJ*)P->data,*c;
96: PetscInt *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj;
97: PetscInt *ci,*cj,*ptadenserow,*ptasparserow,*ptaj,nspacedouble=0;
98: PetscInt an=A->cmap->N,am=A->rmap->N,pn=P->cmap->N,pm=P->rmap->N;
99: PetscInt i,j,k,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi,nlnk,*lnk;
100: MatScalar *ca;
101: PetscBT lnkbt;
102: PetscReal afill;
105: /* Get ij structure of P^T */
106: MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);
107: ptJ = ptj;
109: /* Allocate ci array, arrays for fill computation and */
110: /* free space for accumulating nonzero column info */
111: PetscMalloc1(pn+1,&ci);
112: ci[0] = 0;
114: PetscCalloc1(2*an+1,&ptadenserow);
115: ptasparserow = ptadenserow + an;
117: /* create and initialize a linked list */
118: nlnk = pn+1;
119: PetscLLCreate(pn,pn,nlnk,lnk,lnkbt);
121: /* Set initial free space to be fill*(nnz(A)+ nnz(P)) */
122: PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am],pi[pm])),&free_space);
123: current_space = free_space;
125: /* Determine symbolic info for each row of C: */
126: for (i=0; i<pn; i++) {
127: ptnzi = pti[i+1] - pti[i];
128: ptanzi = 0;
129: /* Determine symbolic row of PtA: */
130: for (j=0; j<ptnzi; j++) {
131: arow = *ptJ++;
132: anzj = ai[arow+1] - ai[arow];
133: ajj = aj + ai[arow];
134: for (k=0; k<anzj; k++) {
135: if (!ptadenserow[ajj[k]]) {
136: ptadenserow[ajj[k]] = -1;
137: ptasparserow[ptanzi++] = ajj[k];
138: }
139: }
140: }
141: /* Using symbolic info for row of PtA, determine symbolic info for row of C: */
142: ptaj = ptasparserow;
143: cnzi = 0;
144: for (j=0; j<ptanzi; j++) {
145: prow = *ptaj++;
146: pnzj = pi[prow+1] - pi[prow];
147: pjj = pj + pi[prow];
148: /* add non-zero cols of P into the sorted linked list lnk */
149: PetscLLAddSorted(pnzj,pjj,pn,nlnk,lnk,lnkbt);
150: cnzi += nlnk;
151: }
153: /* If free space is not available, make more free space */
154: /* Double the amount of total space in the list */
155: if (current_space->local_remaining<cnzi) {
156: PetscFreeSpaceGet(PetscIntSumTruncate(cnzi,current_space->total_array_size),¤t_space);
157: nspacedouble++;
158: }
160: /* Copy data into free space, and zero out denserows */
161: PetscLLClean(pn,pn,cnzi,lnk,current_space->array,lnkbt);
163: current_space->array += cnzi;
164: current_space->local_used += cnzi;
165: current_space->local_remaining -= cnzi;
167: for (j=0; j<ptanzi; j++) ptadenserow[ptasparserow[j]] = 0;
169: /* Aside: Perhaps we should save the pta info for the numerical factorization. */
170: /* For now, we will recompute what is needed. */
171: ci[i+1] = ci[i] + cnzi;
172: }
173: /* nnz is now stored in ci[ptm], column indices are in the list of free space */
174: /* Allocate space for cj, initialize cj, and */
175: /* destroy list of free space and other temporary array(s) */
176: PetscMalloc1(ci[pn]+1,&cj);
177: PetscFreeSpaceContiguous(&free_space,cj);
178: PetscFree(ptadenserow);
179: PetscLLDestroy(lnk,lnkbt);
181: PetscCalloc1(ci[pn]+1,&ca);
183: /* put together the new matrix */
184: MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),pn,pn,ci,cj,ca,C);
185: MatSetBlockSizes(*C,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));
186: MatSetType(*C,((PetscObject)A)->type_name);
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: PetscCalloc2(cn,&apa,cn,&apjdense);
233: PetscMalloc1(cn,&apj);
235: /* Clear old values in C */
236: PetscArrayzero(ca,ci[cm]);
238: for (i=0; i<am; i++) {
239: /* Form sparse row of A*P */
240: anzi = ai[i+1] - ai[i];
241: apnzj = 0;
242: for (j=0; j<anzi; j++) {
243: prow = *aj++;
244: pnzj = pi[prow+1] - pi[prow];
245: pjj = pj + pi[prow];
246: paj = pa + pi[prow];
247: for (k=0; k<pnzj; k++) {
248: if (!apjdense[pjj[k]]) {
249: apjdense[pjj[k]] = -1;
250: apj[apnzj++] = pjj[k];
251: }
252: apa[pjj[k]] += (*aa)*paj[k];
253: }
254: PetscLogFlops(2.0*pnzj);
255: aa++;
256: }
258: /* Sort the j index array for quick sparse axpy. */
259: /* Note: a array does not need sorting as it is in dense storage locations. */
260: PetscSortInt(apnzj,apj);
262: /* Compute P^T*A*P using outer product (P^T)[:,j]*(A*P)[j,:]. */
263: pnzi = pi[i+1] - pi[i];
264: for (j=0; j<pnzi; j++) {
265: nextap = 0;
266: crow = *pJ++;
267: cjj = cj + ci[crow];
268: caj = ca + ci[crow];
269: /* Perform sparse axpy operation. Note cjj includes apj. */
270: for (k=0; nextap<apnzj; k++) {
271: #if defined(PETSC_USE_DEBUG)
272: if (k >= ci[crow+1] - ci[crow]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"k too large k %d, crow %d",k,crow);
273: #endif
274: if (cjj[k]==apj[nextap]) {
275: caj[k] += (*pA)*apa[apj[nextap++]];
276: }
277: }
278: PetscLogFlops(2.0*apnzj);
279: pA++;
280: }
282: /* Zero the current row info for A*P */
283: for (j=0; j<apnzj; j++) {
284: apa[apj[j]] = 0.;
285: apjdense[apj[j]] = 0;
286: }
287: }
289: /* Assemble the final matrix and clean up */
290: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
291: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
293: PetscFree2(apa,apjdense);
294: PetscFree(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: }