1 2 /* 3 Defines projective product routines where A is a SeqAIJ matrix 4 C = P^T * A * P 5 */ 6 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> 11 12 #if defined(PETSC_HAVE_HYPRE) 13 PETSC_INTERN PetscErrorCode MatPtAPSymbolic_AIJ_AIJ_wHYPRE(Mat, Mat, PetscReal, Mat); 14 #endif 15 16 PetscErrorCode MatProductSymbolic_PtAP_SeqAIJ_SeqAIJ(Mat C) { 17 Mat_Product *product = C->product; 18 Mat A = product->A, P = product->B; 19 MatProductAlgorithm alg = product->alg; 20 PetscReal fill = product->fill; 21 PetscBool flg; 22 Mat Pt; 23 24 PetscFunctionBegin; 25 /* "scalable" */ 26 PetscCall(PetscStrcmp(alg, "scalable", &flg)); 27 if (flg) { 28 PetscCall(MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(A, P, fill, C)); 29 C->ops->productnumeric = MatProductNumeric_PtAP; 30 PetscFunctionReturn(0); 31 } 32 33 /* "rap" */ 34 PetscCall(PetscStrcmp(alg, "rap", &flg)); 35 if (flg) { 36 Mat_MatTransMatMult *atb; 37 38 PetscCall(PetscNew(&atb)); 39 PetscCall(MatTranspose(P, MAT_INITIAL_MATRIX, &Pt)); 40 PetscCall(MatMatMatMultSymbolic_SeqAIJ_SeqAIJ_SeqAIJ(Pt, A, P, fill, C)); 41 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 PetscFunctionReturn(0); 50 } 51 52 /* hypre */ 53 #if defined(PETSC_HAVE_HYPRE) 54 PetscCall(PetscStrcmp(alg, "hypre", &flg)); 55 if (flg) { 56 PetscCall(MatPtAPSymbolic_AIJ_AIJ_wHYPRE(A, P, fill, C)); 57 PetscFunctionReturn(0); 58 } 59 #endif 60 61 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "MatProductType is not supported"); 62 } 63 64 PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(Mat A, Mat P, PetscReal fill, Mat C) { 65 PetscFreeSpaceList free_space = NULL, current_space = NULL; 66 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data, *p = (Mat_SeqAIJ *)P->data, *c; 67 PetscInt *pti, *ptj, *ptJ, *ai = a->i, *aj = a->j, *ajj, *pi = p->i, *pj = p->j, *pjj; 68 PetscInt *ci, *cj, *ptadenserow, *ptasparserow, *ptaj, nspacedouble = 0; 69 PetscInt an = A->cmap->N, am = A->rmap->N, pn = P->cmap->N, pm = P->rmap->N; 70 PetscInt i, j, k, ptnzi, arow, anzj, ptanzi, prow, pnzj, cnzi, nlnk, *lnk; 71 MatScalar *ca; 72 PetscBT lnkbt; 73 PetscReal afill; 74 75 PetscFunctionBegin; 76 /* Get ij structure of P^T */ 77 PetscCall(MatGetSymbolicTranspose_SeqAIJ(P, &pti, &ptj)); 78 ptJ = ptj; 79 80 /* Allocate ci array, arrays for fill computation and */ 81 /* free space for accumulating nonzero column info */ 82 PetscCall(PetscMalloc1(pn + 1, &ci)); 83 ci[0] = 0; 84 85 PetscCall(PetscCalloc1(2 * an + 1, &ptadenserow)); 86 ptasparserow = ptadenserow + an; 87 88 /* create and initialize a linked list */ 89 nlnk = pn + 1; 90 PetscCall(PetscLLCreate(pn, pn, nlnk, lnk, lnkbt)); 91 92 /* Set initial free space to be fill*(nnz(A)+ nnz(P)) */ 93 PetscCall(PetscFreeSpaceGet(PetscRealIntMultTruncate(fill, PetscIntSumTruncate(ai[am], pi[pm])), &free_space)); 94 current_space = free_space; 95 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 PetscCall(PetscLLAddSorted(pnzj, pjj, pn, &nlnk, lnk, lnkbt)); 121 cnzi += nlnk; 122 } 123 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 PetscCall(PetscFreeSpaceGet(PetscIntSumTruncate(cnzi, current_space->total_array_size), ¤t_space)); 128 nspacedouble++; 129 } 130 131 /* Copy data into free space, and zero out denserows */ 132 PetscCall(PetscLLClean(pn, pn, cnzi, lnk, current_space->array, lnkbt)); 133 134 current_space->array += cnzi; 135 current_space->local_used += cnzi; 136 current_space->local_remaining -= cnzi; 137 138 for (j = 0; j < ptanzi; j++) ptadenserow[ptasparserow[j]] = 0; 139 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 PetscCall(PetscMalloc1(ci[pn] + 1, &cj)); 148 PetscCall(PetscFreeSpaceContiguous(&free_space, cj)); 149 PetscCall(PetscFree(ptadenserow)); 150 PetscCall(PetscLLDestroy(lnk, lnkbt)); 151 152 PetscCall(PetscCalloc1(ci[pn] + 1, &ca)); 153 154 /* put together the new matrix */ 155 PetscCall(MatSetSeqAIJWithArrays_private(PetscObjectComm((PetscObject)A), pn, pn, ci, cj, ca, ((PetscObject)A)->type_name, C)); 156 PetscCall(MatSetBlockSizes(C, PetscAbs(P->cmap->bs), PetscAbs(P->cmap->bs))); 157 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; 164 165 C->ops->ptapnumeric = MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy; 166 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; 173 174 /* Clean up. */ 175 PetscCall(MatRestoreSymbolicTranspose_SeqAIJ(P, &pti, &ptj)); 176 #if defined(PETSC_USE_INFO) 177 if (ci[pn] != 0) { 178 PetscCall(PetscInfo(C, "Reallocs %" PetscInt_FMT "; Fill ratio: given %g needed %g.\n", nspacedouble, (double)fill, (double)afill)); 179 PetscCall(PetscInfo(C, "Use MatPtAP(A,P,MatReuse,%g,&C) for best performance.\n", (double)afill)); 180 } else { 181 PetscCall(PetscInfo(C, "Empty matrix product\n")); 182 } 183 #endif 184 PetscFunctionReturn(0); 185 } 186 187 PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy(Mat A, Mat P, Mat C) { 188 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 189 Mat_SeqAIJ *p = (Mat_SeqAIJ *)P->data; 190 Mat_SeqAIJ *c = (Mat_SeqAIJ *)C->data; 191 PetscInt *ai = a->i, *aj = a->j, *apj, *apjdense, *pi = p->i, *pj = p->j, *pJ = p->j, *pjj; 192 PetscInt *ci = c->i, *cj = c->j, *cjj; 193 PetscInt am = A->rmap->N, cn = C->cmap->N, cm = C->rmap->N; 194 PetscInt i, j, k, anzi, pnzi, apnzj, nextap, pnzj, prow, crow; 195 MatScalar *aa, *apa, *pa, *pA, *paj, *ca, *caj; 196 197 PetscFunctionBegin; 198 /* Allocate temporary array for storage of one row of A*P (cn: non-scalable) */ 199 PetscCall(PetscCalloc2(cn, &apa, cn, &apjdense)); 200 PetscCall(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 PetscCall(MatSeqAIJGetArrayRead(A, &dummy)); 206 PetscCall(MatSeqAIJRestoreArrayRead(A, &dummy)); 207 PetscCall(MatSeqAIJGetArrayRead(P, &dummy)); 208 PetscCall(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; 216 217 /* Clear old values in C */ 218 PetscCall(PetscArrayzero(ca, ci[cm])); 219 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 PetscCall(PetscLogFlops(2.0 * pnzj)); 237 aa++; 238 } 239 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 PetscCall(PetscSortInt(apnzj, apj)); 243 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 PetscCall(PetscLogFlops(2.0 * apnzj)); 257 pA++; 258 } 259 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 } 266 267 /* Assemble the final matrix and clean up */ 268 PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY)); 269 PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY)); 270 271 PetscCall(PetscFree2(apa, apjdense)); 272 PetscCall(PetscFree(apj)); 273 PetscFunctionReturn(0); 274 } 275 276 PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat A, Mat P, Mat C) { 277 Mat_MatTransMatMult *atb; 278 279 PetscFunctionBegin; 280 MatCheckProduct(C, 3); 281 atb = (Mat_MatTransMatMult *)C->product->data; 282 PetscCheck(atb, PetscObjectComm((PetscObject)C), PETSC_ERR_PLIB, "Missing data structure"); 283 PetscCall(MatTranspose(P, MAT_REUSE_MATRIX, &atb->At)); 284 PetscCheck(C->ops->matmultnumeric, PetscObjectComm((PetscObject)C), PETSC_ERR_PLIB, "Missing numeric operation"); 285 /* when using rap, MatMatMatMultSymbolic used a different data */ 286 if (atb->data) C->product->data = atb->data; 287 PetscCall((*C->ops->matmatmultnumeric)(atb->At, A, P, C)); 288 C->product->data = atb; 289 PetscFunctionReturn(0); 290 } 291