Actual source code: matptap.c

petsc-3.6.1 2015-08-06
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>   /*I "petscmat.h" I*/
  8: #include <../src/mat/utils/freespace.h>
  9: #include <petscbt.h>
 10: #include <petsctime.h>

 14: PetscErrorCode MatPtAP_SeqAIJ_SeqAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
 15: {
 17:   const char     *algTypes[2] = {"scalable","nonscalable"};
 18:   PetscInt       alg=0; /* set default algorithm */

 21:   if (scall == MAT_INITIAL_MATRIX) {
 22:     /* 
 23:      Alg 'scalable' determines which implementations to be used:
 24:        "nonscalable": do dense axpy in MatPtAPNumeric() - fastest, but requires storage of struct A*P;
 25:        "scalable":    do two sparse axpy in MatPtAPNumeric() - might slow, does not store structure of A*P. 
 26:      */
 27:     PetscObjectOptionsBegin((PetscObject)A);
 28:     PetscOptionsEList("-matptap_via","Algorithmic approach","MatPtAP",algTypes,2,algTypes[0],&alg,NULL);
 29:     PetscOptionsEnd();
 30:     PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);
 31:     switch (alg) {
 32:     case 1:
 33:       MatPtAPSymbolic_SeqAIJ_SeqAIJ_DenseAxpy(A,P,fill,C);
 34:       break;
 35:     default:
 36:       MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(A,P,fill,C);
 37:       break;
 38:     }
 39:     PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);
 40:   }
 41:   PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);
 42:   (*(*C)->ops->ptapnumeric)(A,P,*C);
 43:   PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);
 44:   return(0);
 45: }

 49: PetscErrorCode MatDestroy_SeqAIJ_PtAP(Mat A)
 50: {
 52:   Mat_SeqAIJ     *a    = (Mat_SeqAIJ*)A->data;
 53:   Mat_PtAP       *ptap = a->ptap;

 56:   PetscFree(ptap->apa);
 57:   PetscFree(ptap->api);
 58:   PetscFree(ptap->apj);
 59:   (ptap->destroy)(A);
 60:   PetscFree(ptap);
 61:   return(0);
 62: }

 66: PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(Mat A,Mat P,PetscReal fill,Mat *C)
 67: {
 68:   PetscErrorCode     ierr;
 69:   PetscFreeSpaceList free_space=NULL,current_space=NULL;
 70:   Mat_SeqAIJ         *a        = (Mat_SeqAIJ*)A->data,*p = (Mat_SeqAIJ*)P->data,*c;
 71:   PetscInt           *pti,*ptj,*ptJ,*ai=a->i,*aj=a->j,*ajj,*pi=p->i,*pj=p->j,*pjj;
 72:   PetscInt           *ci,*cj,*ptadenserow,*ptasparserow,*ptaj,nspacedouble=0;
 73:   PetscInt           an=A->cmap->N,am=A->rmap->N,pn=P->cmap->N,pm=P->rmap->N;
 74:   PetscInt           i,j,k,ptnzi,arow,anzj,ptanzi,prow,pnzj,cnzi,nlnk,*lnk;
 75:   MatScalar          *ca;
 76:   PetscBT            lnkbt;
 77:   PetscReal          afill;

 80:   /* Get ij structure of P^T */
 81:   MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);
 82:   ptJ  = ptj;

 84:   /* Allocate ci array, arrays for fill computation and */
 85:   /* free space for accumulating nonzero column info */
 86:   PetscMalloc1(pn+1,&ci);
 87:   ci[0] = 0;

 89:   PetscCalloc1(2*an+1,&ptadenserow);
 90:   ptasparserow = ptadenserow  + an;

 92:   /* create and initialize a linked list */
 93:   nlnk = pn+1;
 94:   PetscLLCreate(pn,pn,nlnk,lnk,lnkbt);

 96:   /* Set initial free space to be fill*(nnz(A)+ nnz(P)) */
 97:   PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+pi[pm])),&free_space);
 98:   current_space = free_space;

100:   /* Determine symbolic info for each row of C: */
101:   for (i=0; i<pn; i++) {
102:     ptnzi  = pti[i+1] - pti[i];
103:     ptanzi = 0;
104:     /* Determine symbolic row of PtA: */
105:     for (j=0; j<ptnzi; j++) {
106:       arow = *ptJ++;
107:       anzj = ai[arow+1] - ai[arow];
108:       ajj  = aj + ai[arow];
109:       for (k=0; k<anzj; k++) {
110:         if (!ptadenserow[ajj[k]]) {
111:           ptadenserow[ajj[k]]    = -1;
112:           ptasparserow[ptanzi++] = ajj[k];
113:         }
114:       }
115:     }
116:     /* Using symbolic info for row of PtA, determine symbolic info for row of C: */
117:     ptaj = ptasparserow;
118:     cnzi = 0;
119:     for (j=0; j<ptanzi; j++) {
120:       prow = *ptaj++;
121:       pnzj = pi[prow+1] - pi[prow];
122:       pjj  = pj + pi[prow];
123:       /* add non-zero cols of P into the sorted linked list lnk */
124:       PetscLLAddSorted(pnzj,pjj,pn,nlnk,lnk,lnkbt);
125:       cnzi += nlnk;
126:     }

128:     /* If free space is not available, make more free space */
129:     /* Double the amount of total space in the list */
130:     if (current_space->local_remaining<cnzi) {
131:       PetscFreeSpaceGet(cnzi+current_space->total_array_size,&current_space);
132:       nspacedouble++;
133:     }

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

138:     current_space->array           += cnzi;
139:     current_space->local_used      += cnzi;
140:     current_space->local_remaining -= cnzi;

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

144:     /* Aside: Perhaps we should save the pta info for the numerical factorization. */
145:     /*        For now, we will recompute what is needed. */
146:     ci[i+1] = ci[i] + cnzi;
147:   }
148:   /* nnz is now stored in ci[ptm], column indices are in the list of free space */
149:   /* Allocate space for cj, initialize cj, and */
150:   /* destroy list of free space and other temporary array(s) */
151:   PetscMalloc1(ci[pn]+1,&cj);
152:   PetscFreeSpaceContiguous(&free_space,cj);
153:   PetscFree(ptadenserow);
154:   PetscLLDestroy(lnk,lnkbt);

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

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

162:   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
163:   /* Since these are PETSc arrays, change flags to free them as necessary. */
164:   c          = (Mat_SeqAIJ*)((*C)->data);
165:   c->free_a  = PETSC_TRUE;
166:   c->free_ij = PETSC_TRUE;
167:   c->nonew   = 0;
168:   (*C)->ops->ptapnumeric = MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy;

170:   /* set MatInfo */
171:   afill = (PetscReal)ci[pn]/(ai[am]+pi[pm] + 1.e-5);
172:   if (afill < 1.0) afill = 1.0;
173:   c->maxnz                     = ci[pn];
174:   c->nz                        = ci[pn];
175:   (*C)->info.mallocs           = nspacedouble;
176:   (*C)->info.fill_ratio_given  = fill;
177:   (*C)->info.fill_ratio_needed = afill;

179:   /* Clean up. */
180:   MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);
181: #if defined(PETSC_USE_INFO)
182:   if (ci[pn] != 0) {
183:     PetscInfo3((*C),"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);
184:     PetscInfo1((*C),"Use MatPtAP(A,P,MatReuse,%g,&C) for best performance.\n",(double)afill);
185:   } else {
186:     PetscInfo((*C),"Empty matrix product\n");
187:   }
188: #endif
189:   return(0);
190: }

194: PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy(Mat A,Mat P,Mat C)
195: {
197:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*) A->data;
198:   Mat_SeqAIJ     *p = (Mat_SeqAIJ*) P->data;
199:   Mat_SeqAIJ     *c = (Mat_SeqAIJ*) C->data;
200:   PetscInt       *ai=a->i,*aj=a->j,*apj,*apjdense,*pi=p->i,*pj=p->j,*pJ=p->j,*pjj;
201:   PetscInt       *ci=c->i,*cj=c->j,*cjj;
202:   PetscInt       am =A->rmap->N,cn=C->cmap->N,cm=C->rmap->N;
203:   PetscInt       i,j,k,anzi,pnzi,apnzj,nextap,pnzj,prow,crow;
204:   MatScalar      *aa=a->a,*apa,*pa=p->a,*pA=p->a,*paj,*ca=c->a,*caj;

207:   /* Allocate temporary array for storage of one row of A*P (cn: non-scalable) */
208:   PetscMalloc3(cn,&apa,cn,&apjdense,c->rmax,&apj);
209:   PetscMemzero(apa,cn*sizeof(MatScalar));
210:   PetscMemzero(apjdense,cn*sizeof(PetscInt));

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

215:   for (i=0; i<am; i++) {
216:     /* Form sparse row of A*P */
217:     anzi  = ai[i+1] - ai[i];
218:     apnzj = 0;
219:     for (j=0; j<anzi; j++) {
220:       prow = *aj++;
221:       pnzj = pi[prow+1] - pi[prow];
222:       pjj  = pj + pi[prow];
223:       paj  = pa + pi[prow];
224:       for (k=0; k<pnzj; k++) {
225:         if (!apjdense[pjj[k]]) {
226:           apjdense[pjj[k]] = -1;
227:           apj[apnzj++]     = pjj[k];
228:         }
229:         apa[pjj[k]] += (*aa)*paj[k];
230:       }
231:       PetscLogFlops(2.0*pnzj);
232:       aa++;
233:     }

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

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

259:     /* Zero the current row info for A*P */
260:     for (j=0; j<apnzj; j++) {
261:       apa[apj[j]]      = 0.;
262:       apjdense[apj[j]] = 0;
263:     }
264:   }

266:   /* Assemble the final matrix and clean up */
267:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
268:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

270:   PetscFree3(apa,apjdense,apj);
271:   return(0);
272: }

276: PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_DenseAxpy(Mat A,Mat P,PetscReal fill,Mat *C)
277: {
279:   Mat_SeqAIJ     *ap,*c;
280:   PetscInt       *api,*apj,*ci,pn=P->cmap->N;
281:   MatScalar      *ca;
282:   Mat_PtAP       *ptap;
283:   Mat            Pt,AP;

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

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

292:   ap          = (Mat_SeqAIJ*)AP->data;
293:   api         = ap->i;
294:   apj         = ap->j;
295:   ap->free_ij = PETSC_FALSE; /* api and apj are kept in struct ptap, cannot be destroyed with AP */

297:   /* Get C = Pt*AP */
298:   MatMatMultSymbolic_SeqAIJ_SeqAIJ(Pt,AP,fill,C);

300:   c         = (Mat_SeqAIJ*)(*C)->data;
301:   ci        = c->i;
302:   PetscCalloc1(ci[pn]+1,&ca);
303:   c->a      = ca;
304:   c->free_a = PETSC_TRUE;

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

309:   c->ptap            = ptap;
310:   ptap->destroy      = (*C)->ops->destroy;
311:   (*C)->ops->destroy = MatDestroy_SeqAIJ_PtAP;

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

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

318:   ptap->api = api;
319:   ptap->apj = apj;

321:   /* Clean up. */
322:   MatDestroy(&Pt);
323:   MatDestroy(&AP);
324: #if defined(PETSC_USE_INFO)
325:   PetscInfo1((*C),"given fill %g\n",(double)fill);
326: #endif
327:   return(0);
328: }

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

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

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

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

382:     time_Cseq0 += tf - t0;
383: #endif

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

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

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

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