Actual source code: matptap.c

petsc-3.9.4 2018-09-11
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: {
 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:       MatPtAPSymbolic_AIJ_AIJ_wHYPRE(A,P,fill,C);
 60:       break;
 61: #endif
 62:     default:
 63:       MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(A,P,fill,C);
 64:       break;
 65:     }
 66:   }
 67:   (*(*C)->ops->ptapnumeric)(A,P,*C);
 68:   return(0);
 69: }

 71: PetscErrorCode MatDestroy_SeqAIJ_PtAP(Mat A)
 72: {
 74:   Mat_SeqAIJ     *a    = (Mat_SeqAIJ*)A->data;
 75:   Mat_PtAP       *ptap = a->ptap;

 78:   PetscFree(ptap->apa);
 79:   PetscFree(ptap->api);
 80:   PetscFree(ptap->apj);
 81:   (ptap->destroy)(A);
 82:   PetscFree(ptap);
 83:   return(0);
 84: }

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

100:   /* Get ij structure of P^T */
101:   MatGetSymbolicTranspose_SeqAIJ(P,&pti,&ptj);
102:   ptJ  = ptj;

104:   /* Allocate ci array, arrays for fill computation and */
105:   /* free space for accumulating nonzero column info */
106:   PetscMalloc1(pn+1,&ci);
107:   ci[0] = 0;

109:   PetscCalloc1(2*an+1,&ptadenserow);
110:   ptasparserow = ptadenserow  + an;

112:   /* create and initialize a linked list */
113:   nlnk = pn+1;
114:   PetscLLCreate(pn,pn,nlnk,lnk,lnkbt);

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

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

148:     /* If free space is not available, make more free space */
149:     /* Double the amount of total space in the list */
150:     if (current_space->local_remaining<cnzi) {
151:       PetscFreeSpaceGet(PetscIntSumTruncate(cnzi,current_space->total_array_size),&current_space);
152:       nspacedouble++;
153:     }

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

158:     current_space->array           += cnzi;
159:     current_space->local_used      += cnzi;
160:     current_space->local_remaining -= cnzi;

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

164:     /* Aside: Perhaps we should save the pta info for the numerical factorization. */
165:     /*        For now, we will recompute what is needed. */
166:     ci[i+1] = ci[i] + cnzi;
167:   }
168:   /* nnz is now stored in ci[ptm], column indices are in the list of free space */
169:   /* Allocate space for cj, initialize cj, and */
170:   /* destroy list of free space and other temporary array(s) */
171:   PetscMalloc1(ci[pn]+1,&cj);
172:   PetscFreeSpaceContiguous(&free_space,cj);
173:   PetscFree(ptadenserow);
174:   PetscLLDestroy(lnk,lnkbt);

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

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

182:   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
183:   /* Since these are PETSc arrays, change flags to free them as necessary. */
184:   c          = (Mat_SeqAIJ*)((*C)->data);
185:   c->free_a  = PETSC_TRUE;
186:   c->free_ij = PETSC_TRUE;
187:   c->nonew   = 0;
188:   (*C)->ops->ptapnumeric = MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy;

190:   /* set MatInfo */
191:   afill = (PetscReal)ci[pn]/(ai[am]+pi[pm] + 1.e-5);
192:   if (afill < 1.0) afill = 1.0;
193:   c->maxnz                     = ci[pn];
194:   c->nz                        = ci[pn];
195:   (*C)->info.mallocs           = nspacedouble;
196:   (*C)->info.fill_ratio_given  = fill;
197:   (*C)->info.fill_ratio_needed = afill;

199:   /* Clean up. */
200:   MatRestoreSymbolicTranspose_SeqAIJ(P,&pti,&ptj);
201: #if defined(PETSC_USE_INFO)
202:   if (ci[pn] != 0) {
203:     PetscInfo3((*C),"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);
204:     PetscInfo1((*C),"Use MatPtAP(A,P,MatReuse,%g,&C) for best performance.\n",(double)afill);
205:   } else {
206:     PetscInfo((*C),"Empty matrix product\n");
207:   }
208: #endif
209:   return(0);
210: }

212: PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy(Mat A,Mat P,Mat C)
213: {
215:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*) A->data;
216:   Mat_SeqAIJ     *p = (Mat_SeqAIJ*) P->data;
217:   Mat_SeqAIJ     *c = (Mat_SeqAIJ*) C->data;
218:   PetscInt       *ai=a->i,*aj=a->j,*apj,*apjdense,*pi=p->i,*pj=p->j,*pJ=p->j,*pjj;
219:   PetscInt       *ci=c->i,*cj=c->j,*cjj;
220:   PetscInt       am =A->rmap->N,cn=C->cmap->N,cm=C->rmap->N;
221:   PetscInt       i,j,k,anzi,pnzi,apnzj,nextap,pnzj,prow,crow;
222:   MatScalar      *aa=a->a,*apa,*pa=p->a,*pA=p->a,*paj,*ca=c->a,*caj;

225:   /* Allocate temporary array for storage of one row of A*P (cn: non-scalable) */
226:   PetscMalloc3(cn,&apa,cn,&apjdense,cn,&apj);
227:   PetscMemzero(apa,cn*sizeof(MatScalar));
228:   PetscMemzero(apjdense,cn*sizeof(PetscInt));

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

233:   for (i=0; i<am; i++) {
234:     /* Form sparse row of A*P */
235:     anzi  = ai[i+1] - ai[i];
236:     apnzj = 0;
237:     for (j=0; j<anzi; j++) {
238:       prow = *aj++;
239:       pnzj = pi[prow+1] - pi[prow];
240:       pjj  = pj + pi[prow];
241:       paj  = pa + pi[prow];
242:       for (k=0; k<pnzj; k++) {
243:         if (!apjdense[pjj[k]]) {
244:           apjdense[pjj[k]] = -1;
245:           apj[apnzj++]     = pjj[k];
246:         }
247:         apa[pjj[k]] += (*aa)*paj[k];
248:       }
249:       PetscLogFlops(2.0*pnzj);
250:       aa++;
251:     }

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

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

277:     /* Zero the current row info for A*P */
278:     for (j=0; j<apnzj; j++) {
279:       apa[apj[j]]      = 0.;
280:       apjdense[apj[j]] = 0;
281:     }
282:   }

284:   /* Assemble the final matrix and clean up */
285:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
286:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

288:   PetscFree3(apa,apjdense,apj);
289:   return(0);
290: }

292: PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat A,Mat P,Mat C)
293: {
294:   PetscErrorCode      ierr;
295:   Mat_SeqAIJ          *c = (Mat_SeqAIJ*)C->data;
296:   Mat_MatTransMatMult *atb = c->atb;
297:   Mat                 Pt = atb->At;

300:   MatTranspose_SeqAIJ(P,MAT_REUSE_MATRIX,&Pt);
301:   (C->ops->matmatmultnumeric)(Pt,A,P,C);
302:   return(0);
303: }