Actual source code: matptap.c


  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:   Mat_Product        *product = C->product;
 19:   Mat                 A = product->A, P = product->B;
 20:   MatProductAlgorithm alg  = product->alg;
 21:   PetscReal           fill = product->fill;
 22:   PetscBool           flg;
 23:   Mat                 Pt;

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

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

 38:     PetscNew(&atb);
 39:     MatTranspose(P, MAT_INITIAL_MATRIX, &Pt);
 40:     MatMatMatMultSymbolic_SeqAIJ_SeqAIJ_SeqAIJ(Pt, A, P, fill, C);

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

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

 61:   SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "MatProductType is not supported");
 62: }

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

 76:   /* Get ij structure of P^T */
 77:   MatGetSymbolicTranspose_SeqAIJ(P, &pti, &ptj);
 78:   ptJ = ptj;

 80:   /* Allocate ci array, arrays for fill computation and */
 81:   /* free space for accumulating nonzero column info */
 82:   PetscMalloc1(pn + 1, &ci);
 83:   ci[0] = 0;

 85:   PetscCalloc1(2 * an + 1, &ptadenserow);
 86:   ptasparserow = ptadenserow + an;

 88:   /* create and initialize a linked list */
 89:   nlnk = pn + 1;
 90:   PetscLLCreate(pn, pn, nlnk, lnk, lnkbt);

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

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

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

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

134:     current_space->array += cnzi;
135:     current_space->local_used += cnzi;
136:     current_space->local_remaining -= cnzi;

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

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

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

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

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

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

167:   /* set MatInfo */
168:   afill = (PetscReal)ci[pn] / (ai[am] + pi[pm] + 1.e-5);
169:   if (afill < 1.0) afill = 1.0;
170:   C->info.mallocs           = nspacedouble;
171:   C->info.fill_ratio_given  = fill;
172:   C->info.fill_ratio_needed = afill;

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

187: PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy(Mat A, Mat P, Mat C)
188: {
189:   Mat_SeqAIJ *a  = (Mat_SeqAIJ *)A->data;
190:   Mat_SeqAIJ *p  = (Mat_SeqAIJ *)P->data;
191:   Mat_SeqAIJ *c  = (Mat_SeqAIJ *)C->data;
192:   PetscInt   *ai = a->i, *aj = a->j, *apj, *apjdense, *pi = p->i, *pj = p->j, *pJ = p->j, *pjj;
193:   PetscInt   *ci = c->i, *cj = c->j, *cjj;
194:   PetscInt    am = A->rmap->N, cn = C->cmap->N, cm = C->rmap->N;
195:   PetscInt    i, j, k, anzi, pnzi, apnzj, nextap, pnzj, prow, crow;
196:   MatScalar  *aa, *apa, *pa, *pA, *paj, *ca, *caj;

198:   /* Allocate temporary array for storage of one row of A*P (cn: non-scalable) */
199:   PetscCalloc2(cn, &apa, cn, &apjdense);
200:   PetscMalloc1(cn, &apj);
201:   /* trigger CPU copies if needed and flag CPU mask for C */
202: #if defined(PETSC_HAVE_DEVICE)
203:   {
204:     const PetscScalar *dummy;
205:     MatSeqAIJGetArrayRead(A, &dummy);
206:     MatSeqAIJRestoreArrayRead(A, &dummy);
207:     MatSeqAIJGetArrayRead(P, &dummy);
208:     MatSeqAIJRestoreArrayRead(P, &dummy);
209:     if (C->offloadmask != PETSC_OFFLOAD_UNALLOCATED) C->offloadmask = PETSC_OFFLOAD_CPU;
210:   }
211: #endif
212:   aa = a->a;
213:   pa = p->a;
214:   pA = p->a;
215:   ca = c->a;

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

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

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

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

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

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

271:   PetscFree2(apa, apjdense);
272:   PetscFree(apj);
273:   return 0;
274: }

276: PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat A, Mat P, Mat C)
277: {
278:   Mat_MatTransMatMult *atb;

280:   MatCheckProduct(C, 3);
281:   atb = (Mat_MatTransMatMult *)C->product->data;
283:   MatTranspose(P, MAT_REUSE_MATRIX, &atb->At);
285:   /* when using rap, MatMatMatMultSymbolic used a different data */
286:   if (atb->data) C->product->data = atb->data;
287:   (*C->ops->matmatmultnumeric)(atb->At, A, P, C);
288:   C->product->data = atb;
289:   return 0;
290: }