Actual source code: matptap.c

petsc-master 2020-10-26
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: PetscErrorCode MatProductSymbolic_PtAP_SeqAIJ_SeqAIJ(Mat C)
 17: {
 18:   PetscErrorCode      ierr;
 19:   Mat_Product         *product = C->product;
 20:   Mat                 A=product->A,P=product->B;
 21:   MatProductAlgorithm alg=product->alg;
 22:   PetscReal           fill=product->fill;
 23:   PetscBool           flg;
 24:   Mat                 Pt;

 27:   /* "scalable" */
 28:   PetscStrcmp(alg,"scalable",&flg);
 29:   if (flg) {
 30:     MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(A,P,fill,C);
 31:     C->ops->productnumeric = MatProductNumeric_PtAP;
 32:     return(0);
 33:   }

 35:   /* "rap" */
 36:   PetscStrcmp(alg,"rap",&flg);
 37:   if (flg) {
 38:     Mat_MatTransMatMult *atb;

 40:     PetscNew(&atb);
 41:     MatTranspose_SeqAIJ(P,MAT_INITIAL_MATRIX,&Pt);
 42:     MatMatMatMultSymbolic_SeqAIJ_SeqAIJ_SeqAIJ(Pt,A,P,fill,C);

 44:     atb->At                = Pt;
 45:     atb->data              = C->product->data;
 46:     atb->destroy           = C->product->destroy;
 47:     C->product->data       = atb;
 48:     C->product->destroy    = MatDestroy_SeqAIJ_MatTransMatMult;
 49:     C->ops->ptapnumeric    = MatPtAPNumeric_SeqAIJ_SeqAIJ;
 50:     C->ops->productnumeric = MatProductNumeric_PtAP;
 51:     return(0);
 52:   }

 54:   /* hypre */
 55: #if defined(PETSC_HAVE_HYPRE)
 56:   PetscStrcmp(alg,"hypre",&flg);
 57:   if (flg) {
 58:     MatPtAPSymbolic_AIJ_AIJ_wHYPRE(A,P,fill,C);
 59:     return(0);
 60:   }
 61: #endif

 63:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatProductType is not supported");
 64:   return(0);
 65: }

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

170:   C->ops->ptapnumeric = MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy;

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

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

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:   PetscCalloc2(cn,&apa,cn,&apjdense);
209:   PetscMalloc1(cn,&apj);

211:   /* Clear old values in C */
212:   PetscArrayzero(ca,ci[cm]);

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

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

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

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

263:   /* Assemble the final matrix and clean up */
264:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
265:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

267:   PetscFree2(apa,apjdense);
268:   PetscFree(apj);
269:   return(0);
270: }

272: PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat A,Mat P,Mat C)
273: {
274:   PetscErrorCode      ierr;
275:   Mat_MatTransMatMult *atb;

278:   MatCheckProduct(C,3);
279:   atb  = (Mat_MatTransMatMult*)C->product->data;
280:   if (!atb) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_PLIB,"Missing data structure");
281:   MatTranspose_SeqAIJ(P,MAT_REUSE_MATRIX,&atb->At);
282:   if (!C->ops->matmultnumeric) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_PLIB,"Missing numeric operation");
283:   /* when using rap, MatMatMatMultSymbolic used a different data */
284:   if (atb->data) C->product->data = atb->data;
285:   (*C->ops->matmatmultnumeric)(atb->At,A,P,C);
286:   C->product->data = atb;
287:   return(0);
288: }