Actual source code: matptap.c

petsc-3.8.4 2018-03-24
Report Typos and Errors

  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: {
 19: #if !defined(PETSC_HAVE_HYPRE)
 20:   const char     *algTypes[2] = {"scalable","nonscalable"};
 21:   PetscInt       nalg = 2;
 22: #else
 23:   const char     *algTypes[3] = {"scalable","nonscalable","hypre"};
 24:   PetscInt       nalg = 3;
 25: #endif
 26:   PetscInt       alg = 0; /* set default algorithm */

 29:   if (scall == MAT_INITIAL_MATRIX) {
 30:     /*
 31:      Alg 'scalable' determines which implementations to be used:
 32:        "nonscalable": do dense axpy in MatPtAPNumeric() - fastest, but requires storage of struct A*P;
 33:        "scalable":    do two sparse axpy in MatPtAPNumeric() - might slow, does not store structure of A*P.
 34:        "hypre":    use boomerAMGBuildCoarseOperator.
 35:      */
 36:     PetscObjectOptionsBegin((PetscObject)A);
 37:     PetscOptionsObject->alreadyprinted = PETSC_FALSE; /* a hack to ensure the option shows in '-help' */
 38:     PetscOptionsEList("-matptap_via","Algorithmic approach","MatPtAP",algTypes,nalg,algTypes[0],&alg,NULL);
 39:     PetscOptionsEnd();
 40:     PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);
 41:     switch (alg) {
 42:     case 1:
 43:       MatPtAPSymbolic_SeqAIJ_SeqAIJ_DenseAxpy(A,P,fill,C);
 44:       break;
 45: #if defined(PETSC_HAVE_HYPRE)
 46:     case 2:
 47:       MatPtAPSymbolic_AIJ_AIJ_wHYPRE(A,P,fill,C);
 48:       break;
 49: #endif
 50:     default:
 51:       MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(A,P,fill,C);
 52:       break;
 53:     }
 54:     PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);
 55:   }
 56:   PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);
 57:   (*(*C)->ops->ptapnumeric)(A,P,*C);
 58:   PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);
 59:   return(0);
 60: }

 62: PetscErrorCode MatDestroy_SeqAIJ_PtAP(Mat A)
 63: {
 65:   Mat_SeqAIJ     *a    = (Mat_SeqAIJ*)A->data;
 66:   Mat_PtAP       *ptap = a->ptap;

 69:   PetscFree(ptap->apa);
 70:   PetscFree(ptap->api);
 71:   PetscFree(ptap->apj);
 72:   (ptap->destroy)(A);
 73:   PetscFree(ptap);
 74:   return(0);
 75: }

 77: PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(Mat A,Mat P,PetscReal fill,Mat *C)
 78: {
 79:   PetscErrorCode     ierr;
 80:   PetscFreeSpaceList free_space=NULL,current_space=NULL;
 81:   Mat_SeqAIJ         *a        = (Mat_SeqAIJ*)A->data,*p = (Mat_SeqAIJ*)P->data,*c;
 82:   PetscInt           *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj;
 83:   PetscInt           *ci,*cj,*ptadenserow,*ptasparserow,*ptaj,nspacedouble=0;
 84:   PetscInt           an=A->cmap->N,am=A->rmap->N,pn=P->cmap->N,pm=P->rmap->N;
 85:   PetscInt           i,j,k,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi,nlnk,*lnk;
 86:   MatScalar          *ca;
 87:   PetscBT            lnkbt;
 88:   PetscReal          afill;

 91:   /* Get ij structure of P^T */
 92:   MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);
 93:   ptJ  = ptj;

 95:   /* Allocate ci array, arrays for fill computation and */
 96:   /* free space for accumulating nonzero column info */
 97:   PetscMalloc1(pn+1,&ci);
 98:   ci[0] = 0;

100:   PetscCalloc1(2*an+1,&ptadenserow);
101:   ptasparserow = ptadenserow  + an;

103:   /* create and initialize a linked list */
104:   nlnk = pn+1;
105:   PetscLLCreate(pn,pn,nlnk,lnk,lnkbt);

107:   /* Set initial free space to be fill*(nnz(A)+ nnz(P)) */
108:   PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am],pi[pm])),&free_space);
109:   current_space = free_space;

111:   /* Determine symbolic info for each row of C: */
112:   for (i=0; i<pn; i++) {
113:     ptnzi  = pti[i+1] - pti[i];
114:     ptanzi = 0;
115:     /* Determine symbolic row of PtA: */
116:     for (j=0; j<ptnzi; j++) {
117:       arow = *ptJ++;
118:       anzj = ai[arow+1] - ai[arow];
119:       ajj  = aj + ai[arow];
120:       for (k=0; k<anzj; k++) {
121:         if (!ptadenserow[ajj[k]]) {
122:           ptadenserow[ajj[k]]    = -1;
123:           ptasparserow[ptanzi++] = ajj[k];
124:         }
125:       }
126:     }
127:     /* Using symbolic info for row of PtA, determine symbolic info for row of C: */
128:     ptaj = ptasparserow;
129:     cnzi = 0;
130:     for (j=0; j<ptanzi; j++) {
131:       prow = *ptaj++;
132:       pnzj = pi[prow+1] - pi[prow];
133:       pjj  = pj + pi[prow];
134:       /* add non-zero cols of P into the sorted linked list lnk */
135:       PetscLLAddSorted(pnzj,pjj,pn,nlnk,lnk,lnkbt);
136:       cnzi += nlnk;
137:     }

139:     /* If free space is not available, make more free space */
140:     /* Double the amount of total space in the list */
141:     if (current_space->local_remaining<cnzi) {
142:       PetscFreeSpaceGet(PetscIntSumTruncate(cnzi,current_space->total_array_size),&current_space);
143:       nspacedouble++;
144:     }

146:     /* Copy data into free space, and zero out denserows */
147:     PetscLLClean(pn,pn,cnzi,lnk,current_space->array,lnkbt);

149:     current_space->array           += cnzi;
150:     current_space->local_used      += cnzi;
151:     current_space->local_remaining -= cnzi;

153:     for (j=0; j<ptanzi; j++) ptadenserow[ptasparserow[j]] = 0;

155:     /* Aside: Perhaps we should save the pta info for the numerical factorization. */
156:     /*        For now, we will recompute what is needed. */
157:     ci[i+1] = ci[i] + cnzi;
158:   }
159:   /* nnz is now stored in ci[ptm], column indices are in the list of free space */
160:   /* Allocate space for cj, initialize cj, and */
161:   /* destroy list of free space and other temporary array(s) */
162:   PetscMalloc1(ci[pn]+1,&cj);
163:   PetscFreeSpaceContiguous(&free_space,cj);
164:   PetscFree(ptadenserow);
165:   PetscLLDestroy(lnk,lnkbt);

167:   PetscCalloc1(ci[pn]+1,&ca);

169:   /* put together the new matrix */
170:   MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),pn,pn,ci,cj,ca,C);
171:   MatSetBlockSizes(*C,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));

173:   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
174:   /* Since these are PETSc arrays, change flags to free them as necessary. */
175:   c          = (Mat_SeqAIJ*)((*C)->data);
176:   c->free_a  = PETSC_TRUE;
177:   c->free_ij = PETSC_TRUE;
178:   c->nonew   = 0;
179:   (*C)->ops->ptapnumeric = MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy;

181:   /* set MatInfo */
182:   afill = (PetscReal)ci[pn]/(ai[am]+pi[pm] + 1.e-5);
183:   if (afill < 1.0) afill = 1.0;
184:   c->maxnz                     = ci[pn];
185:   c->nz                        = ci[pn];
186:   (*C)->info.mallocs           = nspacedouble;
187:   (*C)->info.fill_ratio_given  = fill;
188:   (*C)->info.fill_ratio_needed = afill;

190:   /* Clean up. */
191:   MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);
192: #if defined(PETSC_USE_INFO)
193:   if (ci[pn] != 0) {
194:     PetscInfo3((*C),"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);
195:     PetscInfo1((*C),"Use MatPtAP(A,P,MatReuse,%g,&C) for best performance.\n",(double)afill);
196:   } else {
197:     PetscInfo((*C),"Empty matrix product\n");
198:   }
199: #endif
200:   return(0);
201: }

203: PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy(Mat A,Mat P,Mat C)
204: {
206:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*) A->data;
207:   Mat_SeqAIJ     *p = (Mat_SeqAIJ*) P->data;
208:   Mat_SeqAIJ     *c = (Mat_SeqAIJ*) C->data;
209:   PetscInt       *ai=a->i,*aj=a->j,*apj,*apjdense,*pi=p->i,*pj=p->j,*pJ=p->j,*pjj;
210:   PetscInt       *ci=c->i,*cj=c->j,*cjj;
211:   PetscInt       am =A->rmap->N,cn=C->cmap->N,cm=C->rmap->N;
212:   PetscInt       i,j,k,anzi,pnzi,apnzj,nextap,pnzj,prow,crow;
213:   MatScalar      *aa=a->a,*apa,*pa=p->a,*pA=p->a,*paj,*ca=c->a,*caj;

216:   /* Allocate temporary array for storage of one row of A*P (cn: non-scalable) */
217:   PetscMalloc3(cn,&apa,cn,&apjdense,cn,&apj);
218:   PetscMemzero(apa,cn*sizeof(MatScalar));
219:   PetscMemzero(apjdense,cn*sizeof(PetscInt));

221:   /* Clear old values in C */
222:   PetscMemzero(ca,ci[cm]*sizeof(MatScalar));

224:   for (i=0; i<am; i++) {
225:     /* Form sparse row of A*P */
226:     anzi  = ai[i+1] - ai[i];
227:     apnzj = 0;
228:     for (j=0; j<anzi; j++) {
229:       prow = *aj++;
230:       pnzj = pi[prow+1] - pi[prow];
231:       pjj  = pj + pi[prow];
232:       paj  = pa + pi[prow];
233:       for (k=0; k<pnzj; k++) {
234:         if (!apjdense[pjj[k]]) {
235:           apjdense[pjj[k]] = -1;
236:           apj[apnzj++]     = pjj[k];
237:         }
238:         apa[pjj[k]] += (*aa)*paj[k];
239:       }
240:       PetscLogFlops(2.0*pnzj);
241:       aa++;
242:     }

244:     /* Sort the j index array for quick sparse axpy. */
245:     /* Note: a array does not need sorting as it is in dense storage locations. */
246:     PetscSortInt(apnzj,apj);

248:     /* Compute P^T*A*P using outer product (P^T)[:,j]*(A*P)[j,:]. */
249:     pnzi = pi[i+1] - pi[i];
250:     for (j=0; j<pnzi; j++) {
251:       nextap = 0;
252:       crow   = *pJ++;
253:       cjj    = cj + ci[crow];
254:       caj    = ca + ci[crow];
255:       /* Perform sparse axpy operation.  Note cjj includes apj. */
256:       for (k=0; nextap<apnzj; k++) {
257: #if defined(PETSC_USE_DEBUG)
258:         if (k >= ci[crow+1] - ci[crow]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"k too large k %d, crow %d",k,crow);
259: #endif
260:         if (cjj[k]==apj[nextap]) {
261:           caj[k] += (*pA)*apa[apj[nextap++]];
262:         }
263:       }
264:       PetscLogFlops(2.0*apnzj);
265:       pA++;
266:     }

268:     /* Zero the current row info for A*P */
269:     for (j=0; j<apnzj; j++) {
270:       apa[apj[j]]      = 0.;
271:       apjdense[apj[j]] = 0;
272:     }
273:   }

275:   /* Assemble the final matrix and clean up */
276:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
277:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

279:   PetscFree3(apa,apjdense,apj);
280:   return(0);
281: }

283: PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_DenseAxpy(Mat A,Mat P,PetscReal fill,Mat *C)
284: {
286:   Mat_SeqAIJ     *ap,*c;
287:   PetscInt       *api,*apj,*ci,pn=P->cmap->N;
288:   MatScalar      *ca;
289:   Mat_PtAP       *ptap;
290:   Mat            Pt,AP;

293:   /* Get symbolic Pt = P^T */
294:   MatTransposeSymbolic_SeqAIJ(P,&Pt);

296:   /* Get symbolic AP = A*P */
297:   MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,P,fill,&AP);

299:   ap          = (Mat_SeqAIJ*)AP->data;
300:   api         = ap->i;
301:   apj         = ap->j;
302:   ap->free_ij = PETSC_FALSE; /* api and apj are kept in struct ptap, cannot be destroyed with AP */

304:   /* Get C = Pt*AP */
305:   MatMatMultSymbolic_SeqAIJ_SeqAIJ(Pt,AP,fill,C);

307:   c         = (Mat_SeqAIJ*)(*C)->data;
308:   ci        = c->i;
309:   PetscCalloc1(ci[pn]+1,&ca);
310:   c->a      = ca;
311:   c->free_a = PETSC_TRUE;

313:   /* Create a supporting struct for reuse by MatPtAPNumeric() */
314:   PetscNew(&ptap);

316:   c->ptap            = ptap;
317:   ptap->destroy      = (*C)->ops->destroy;
318:   (*C)->ops->destroy = MatDestroy_SeqAIJ_PtAP;

320:   /* Allocate temporary array for storage of one row of A*P */
321:   PetscCalloc1(pn+1,&ptap->apa);

323:   (*C)->ops->ptapnumeric = MatPtAPNumeric_SeqAIJ_SeqAIJ;

325:   ptap->api = api;
326:   ptap->apj = apj;

328:   /* Clean up. */
329:   MatDestroy(&Pt);
330:   MatDestroy(&AP);
331: #if defined(PETSC_USE_INFO)
332:   PetscInfo1((*C),"given fill %g\n",(double)fill);
333: #endif
334:   return(0);
335: }

337: /* #define PROFILE_MatPtAPNumeric */
338: PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat A,Mat P,Mat C)
339: {
340:   PetscErrorCode    ierr;
341:   Mat_SeqAIJ        *a = (Mat_SeqAIJ*) A->data;
342:   Mat_SeqAIJ        *p = (Mat_SeqAIJ*) P->data;
343:   Mat_SeqAIJ        *c = (Mat_SeqAIJ*) C->data;
344:   const PetscInt    *ai=a->i,*aj=a->j,*pi=p->i,*pj=p->j,*ci=c->i,*cj=c->j;
345:   const PetscScalar *aa=a->a,*pa=p->a,*pval;
346:   const PetscInt    *apj,*pcol,*cjj;
347:   const PetscInt    am=A->rmap->N,cm=C->rmap->N;
348:   PetscInt          i,j,k,anz,apnz,pnz,prow,crow,cnz;
349:   PetscScalar       *apa,*ca=c->a,*caj,pvalj;
350:   Mat_PtAP          *ptap = c->ptap;
351: #if defined(PROFILE_MatPtAPNumeric)
352:   PetscLogDouble t0,tf,time_Cseq0=0.0,time_Cseq1=0.0;
353:   PetscInt       flops0=0,flops1=0;
354: #endif

357:   /* Get temporary array for storage of one row of A*P */
358:   apa = ptap->apa;

360:   /* Clear old values in C */
361:   PetscMemzero(ca,ci[cm]*sizeof(MatScalar));

363:   for (i=0; i<am; i++) {
364:     /* Form sparse row of AP[i,:] = A[i,:]*P */
365: #if defined(PROFILE_MatPtAPNumeric)
366:     PetscTime(&t0);
367: #endif
368:     anz  = ai[i+1] - ai[i];
369:     for (j=0; j<anz; j++) {
370:       prow = aj[j];
371:       pnz  = pi[prow+1] - pi[prow];
372:       pcol = pj + pi[prow];
373:       pval = pa + pi[prow];
374:       for (k=0; k<pnz; k++) {
375:         apa[pcol[k]] += aa[j]*pval[k];
376:       }
377:       PetscLogFlops(2.0*pnz);
378: #if defined(PROFILE_MatPtAPNumeric)
379:       flops0 += 2.0*pnz;
380: #endif
381:     }
382:     aj += anz; aa += anz;
383: #if defined(PROFILE_MatPtAPNumeric)
384:     PetscTime(&tf);

386:     time_Cseq0 += tf - t0;
387: #endif

389:     /* Compute P^T*A*P using outer product P[i,:]^T*AP[i,:]. */
390: #if defined(PROFILE_MatPtAPNumeric)
391:     PetscTime(&t0);
392: #endif
393:     apj  = ptap->apj + ptap->api[i];
394:     apnz = ptap->api[i+1] - ptap->api[i];
395:     pnz  = pi[i+1] - pi[i];
396:     pcol = pj + pi[i];
397:     pval = pa + pi[i];

399:     /* Perform dense axpy */
400:     for (j=0; j<pnz; j++) {
401:       crow  = pcol[j];
402:       cjj   = cj + ci[crow];
403:       caj   = ca + ci[crow];
404:       pvalj = pval[j];
405:       cnz   = ci[crow+1] - ci[crow];
406:       for (k=0; k<cnz; k++) caj[k] += pvalj*apa[cjj[k]];
407:       PetscLogFlops(2.0*cnz);
408: #if defined(PROFILE_MatPtAPNumeric)
409:       flops1 += 2.0*cnz;
410: #endif
411:     }
412: #if defined(PROFILE_MatPtAPNumeric)
413:     PetscTime(&tf);
414:     time_Cseq1 += tf - t0;
415: #endif

417:     /* Zero the current row info for A*P */
418:     for (j=0; j<apnz; j++) apa[apj[j]] = 0.0;
419:   }

421:   /* Assemble the final matrix and clean up */
422:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
423:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
424: #if defined(PROFILE_MatPtAPNumeric)
425:   printf("PtAPNumeric_SeqAIJ time %g + %g, flops %d %d\n",time_Cseq0,time_Cseq1,flops0,flops1);
426: #endif
427:   return(0);
428: }