1be1d678aSKris Buschelman 2eb9c0419SKris Buschelman /* 342c57489SHong Zhang Defines projective product routines where A is a SeqAIJ matrix 4eb9c0419SKris Buschelman C = P^T * A * P 5eb9c0419SKris Buschelman */ 6eb9c0419SKris Buschelman 7c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h> /*I "petscmat.h" I*/ 8c6db04a5SJed Brown #include <../src/mat/utils/freespace.h> 9c6db04a5SJed Brown #include <petscbt.h> 108563dfccSBarry Smith #include <petsctime.h> 11eb9c0419SKris Buschelman 126f231fbdSstefano_zampini #if defined(PETSC_HAVE_HYPRE) 134222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatPtAPSymbolic_AIJ_AIJ_wHYPRE(Mat, Mat, PetscReal, Mat); 146f231fbdSstefano_zampini #endif 156f231fbdSstefano_zampini 16*9371c9d4SSatish Balay PetscErrorCode MatProductSymbolic_PtAP_SeqAIJ_SeqAIJ(Mat C) { 174222ddf1SHong Zhang Mat_Product *product = C->product; 184222ddf1SHong Zhang Mat A = product->A, P = product->B; 194222ddf1SHong Zhang MatProductAlgorithm alg = product->alg; 204222ddf1SHong Zhang PetscReal fill = product->fill; 214222ddf1SHong Zhang PetscBool flg; 228fa4b5a6SHong Zhang Mat Pt; 237ba1a0bfSKris Buschelman 247ba1a0bfSKris Buschelman PetscFunctionBegin; 254222ddf1SHong Zhang /* "scalable" */ 269566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(alg, "scalable", &flg)); 274222ddf1SHong Zhang if (flg) { 289566063dSJacob Faibussowitsch PetscCall(MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(A, P, fill, C)); 294222ddf1SHong Zhang C->ops->productnumeric = MatProductNumeric_PtAP; 304222ddf1SHong Zhang PetscFunctionReturn(0); 314222ddf1SHong Zhang } 324222ddf1SHong Zhang 334222ddf1SHong Zhang /* "rap" */ 349566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(alg, "rap", &flg)); 356718818eSStefano Zampini if (flg) { 366718818eSStefano Zampini Mat_MatTransMatMult *atb; 376718818eSStefano Zampini 389566063dSJacob Faibussowitsch PetscCall(PetscNew(&atb)); 39acd337a6SBarry Smith PetscCall(MatTranspose(P, MAT_INITIAL_MATRIX, &Pt)); 409566063dSJacob Faibussowitsch PetscCall(MatMatMatMultSymbolic_SeqAIJ_SeqAIJ_SeqAIJ(Pt, A, P, fill, C)); 418fa4b5a6SHong Zhang 428fa4b5a6SHong Zhang atb->At = Pt; 436718818eSStefano Zampini atb->data = C->product->data; 446718818eSStefano Zampini atb->destroy = C->product->destroy; 456718818eSStefano Zampini C->product->data = atb; 466718818eSStefano Zampini C->product->destroy = MatDestroy_SeqAIJ_MatTransMatMult; 474222ddf1SHong Zhang C->ops->ptapnumeric = MatPtAPNumeric_SeqAIJ_SeqAIJ; 484222ddf1SHong Zhang C->ops->productnumeric = MatProductNumeric_PtAP; 498fa4b5a6SHong Zhang PetscFunctionReturn(0); 504222ddf1SHong Zhang } 514222ddf1SHong Zhang 524222ddf1SHong Zhang /* hypre */ 536f231fbdSstefano_zampini #if defined(PETSC_HAVE_HYPRE) 549566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(alg, "hypre", &flg)); 554222ddf1SHong Zhang if (flg) { 569566063dSJacob Faibussowitsch PetscCall(MatPtAPSymbolic_AIJ_AIJ_wHYPRE(A, P, fill, C)); 574222ddf1SHong Zhang PetscFunctionReturn(0); 584222ddf1SHong Zhang } 596f231fbdSstefano_zampini #endif 604222ddf1SHong Zhang 614222ddf1SHong Zhang SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "MatProductType is not supported"); 627ba1a0bfSKris Buschelman } 637ba1a0bfSKris Buschelman 64*9371c9d4SSatish Balay PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(Mat A, Mat P, PetscReal fill, Mat C) { 650298fd71SBarry Smith PetscFreeSpaceList free_space = NULL, current_space = NULL; 66d20bfe6fSHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data, *p = (Mat_SeqAIJ *)P->data, *c; 6753565b12SHong Zhang PetscInt *pti, *ptj, *ptJ, *ai = a->i, *aj = a->j, *ajj, *pi = p->i, *pj = p->j, *pjj; 6853565b12SHong Zhang PetscInt *ci, *cj, *ptadenserow, *ptasparserow, *ptaj, nspacedouble = 0; 697d69fa63SHong Zhang PetscInt an = A->cmap->N, am = A->rmap->N, pn = P->cmap->N, pm = P->rmap->N; 7053565b12SHong Zhang PetscInt i, j, k, ptnzi, arow, anzj, ptanzi, prow, pnzj, cnzi, nlnk, *lnk; 71d20bfe6fSHong Zhang MatScalar *ca; 7253565b12SHong Zhang PetscBT lnkbt; 737d69fa63SHong Zhang PetscReal afill; 74eb9c0419SKris Buschelman 75eb9c0419SKris Buschelman PetscFunctionBegin; 7653565b12SHong Zhang /* Get ij structure of P^T */ 779566063dSJacob Faibussowitsch PetscCall(MatGetSymbolicTranspose_SeqAIJ(P, &pti, &ptj)); 78d20bfe6fSHong Zhang ptJ = ptj; 79eb9c0419SKris Buschelman 80d20bfe6fSHong Zhang /* Allocate ci array, arrays for fill computation and */ 81d20bfe6fSHong Zhang /* free space for accumulating nonzero column info */ 829566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(pn + 1, &ci)); 83d20bfe6fSHong Zhang ci[0] = 0; 84eb9c0419SKris Buschelman 859566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(2 * an + 1, &ptadenserow)); 86d20bfe6fSHong Zhang ptasparserow = ptadenserow + an; 8755a3bba9SHong Zhang 8855a3bba9SHong Zhang /* create and initialize a linked list */ 8955a3bba9SHong Zhang nlnk = pn + 1; 909566063dSJacob Faibussowitsch PetscCall(PetscLLCreate(pn, pn, nlnk, lnk, lnkbt)); 91eb9c0419SKris Buschelman 927d69fa63SHong Zhang /* Set initial free space to be fill*(nnz(A)+ nnz(P)) */ 939566063dSJacob Faibussowitsch PetscCall(PetscFreeSpaceGet(PetscRealIntMultTruncate(fill, PetscIntSumTruncate(ai[am], pi[pm])), &free_space)); 94d20bfe6fSHong Zhang current_space = free_space; 95d20bfe6fSHong Zhang 96d20bfe6fSHong Zhang /* Determine symbolic info for each row of C: */ 97d20bfe6fSHong Zhang for (i = 0; i < pn; i++) { 98d20bfe6fSHong Zhang ptnzi = pti[i + 1] - pti[i]; 99d20bfe6fSHong Zhang ptanzi = 0; 100d20bfe6fSHong Zhang /* Determine symbolic row of PtA: */ 101d20bfe6fSHong Zhang for (j = 0; j < ptnzi; j++) { 102d20bfe6fSHong Zhang arow = *ptJ++; 103d20bfe6fSHong Zhang anzj = ai[arow + 1] - ai[arow]; 104d20bfe6fSHong Zhang ajj = aj + ai[arow]; 105d20bfe6fSHong Zhang for (k = 0; k < anzj; k++) { 106d20bfe6fSHong Zhang if (!ptadenserow[ajj[k]]) { 107d20bfe6fSHong Zhang ptadenserow[ajj[k]] = -1; 108d20bfe6fSHong Zhang ptasparserow[ptanzi++] = ajj[k]; 109d20bfe6fSHong Zhang } 110d20bfe6fSHong Zhang } 111d20bfe6fSHong Zhang } 112d20bfe6fSHong Zhang /* Using symbolic info for row of PtA, determine symbolic info for row of C: */ 113d20bfe6fSHong Zhang ptaj = ptasparserow; 114d20bfe6fSHong Zhang cnzi = 0; 115d20bfe6fSHong Zhang for (j = 0; j < ptanzi; j++) { 116d20bfe6fSHong Zhang prow = *ptaj++; 117d20bfe6fSHong Zhang pnzj = pi[prow + 1] - pi[prow]; 118d20bfe6fSHong Zhang pjj = pj + pi[prow]; 11955a3bba9SHong Zhang /* add non-zero cols of P into the sorted linked list lnk */ 1209566063dSJacob Faibussowitsch PetscCall(PetscLLAddSorted(pnzj, pjj, pn, &nlnk, lnk, lnkbt)); 12155a3bba9SHong Zhang cnzi += nlnk; 122d20bfe6fSHong Zhang } 123d20bfe6fSHong Zhang 124d20bfe6fSHong Zhang /* If free space is not available, make more free space */ 125d20bfe6fSHong Zhang /* Double the amount of total space in the list */ 126d20bfe6fSHong Zhang if (current_space->local_remaining < cnzi) { 1279566063dSJacob Faibussowitsch PetscCall(PetscFreeSpaceGet(PetscIntSumTruncate(cnzi, current_space->total_array_size), ¤t_space)); 128f2b054eeSHong Zhang nspacedouble++; 129d20bfe6fSHong Zhang } 130d20bfe6fSHong Zhang 131d20bfe6fSHong Zhang /* Copy data into free space, and zero out denserows */ 1329566063dSJacob Faibussowitsch PetscCall(PetscLLClean(pn, pn, cnzi, lnk, current_space->array, lnkbt)); 1332205254eSKarl Rupp 134d20bfe6fSHong Zhang current_space->array += cnzi; 135d20bfe6fSHong Zhang current_space->local_used += cnzi; 136d20bfe6fSHong Zhang current_space->local_remaining -= cnzi; 137d20bfe6fSHong Zhang 1382205254eSKarl Rupp for (j = 0; j < ptanzi; j++) ptadenserow[ptasparserow[j]] = 0; 1392205254eSKarl Rupp 140d20bfe6fSHong Zhang /* Aside: Perhaps we should save the pta info for the numerical factorization. */ 141d20bfe6fSHong Zhang /* For now, we will recompute what is needed. */ 142d20bfe6fSHong Zhang ci[i + 1] = ci[i] + cnzi; 143d20bfe6fSHong Zhang } 144d20bfe6fSHong Zhang /* nnz is now stored in ci[ptm], column indices are in the list of free space */ 145d20bfe6fSHong Zhang /* Allocate space for cj, initialize cj, and */ 146d20bfe6fSHong Zhang /* destroy list of free space and other temporary array(s) */ 1479566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ci[pn] + 1, &cj)); 1489566063dSJacob Faibussowitsch PetscCall(PetscFreeSpaceContiguous(&free_space, cj)); 1499566063dSJacob Faibussowitsch PetscCall(PetscFree(ptadenserow)); 1509566063dSJacob Faibussowitsch PetscCall(PetscLLDestroy(lnk, lnkbt)); 151d20bfe6fSHong Zhang 1529566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(ci[pn] + 1, &ca)); 15353565b12SHong Zhang 15453565b12SHong Zhang /* put together the new matrix */ 1559566063dSJacob Faibussowitsch PetscCall(MatSetSeqAIJWithArrays_private(PetscObjectComm((PetscObject)A), pn, pn, ci, cj, ca, ((PetscObject)A)->type_name, C)); 1569566063dSJacob Faibussowitsch PetscCall(MatSetBlockSizes(C, PetscAbs(P->cmap->bs), PetscAbs(P->cmap->bs))); 157ae32dd66SHong Zhang 15853565b12SHong Zhang /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 15953565b12SHong Zhang /* Since these are PETSc arrays, change flags to free them as necessary. */ 1604222ddf1SHong Zhang c = (Mat_SeqAIJ *)((C)->data); 16153565b12SHong Zhang c->free_a = PETSC_TRUE; 16253565b12SHong Zhang c->free_ij = PETSC_TRUE; 16353565b12SHong Zhang c->nonew = 0; 1646718818eSStefano Zampini 1654222ddf1SHong Zhang C->ops->ptapnumeric = MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy; 16653565b12SHong Zhang 1677d69fa63SHong Zhang /* set MatInfo */ 1687d69fa63SHong Zhang afill = (PetscReal)ci[pn] / (ai[am] + pi[pm] + 1.e-5); 1697d69fa63SHong Zhang if (afill < 1.0) afill = 1.0; 1704222ddf1SHong Zhang C->info.mallocs = nspacedouble; 1714222ddf1SHong Zhang C->info.fill_ratio_given = fill; 1724222ddf1SHong Zhang C->info.fill_ratio_needed = afill; 1737d69fa63SHong Zhang 17453565b12SHong Zhang /* Clean up. */ 1759566063dSJacob Faibussowitsch PetscCall(MatRestoreSymbolicTranspose_SeqAIJ(P, &pti, &ptj)); 17653565b12SHong Zhang #if defined(PETSC_USE_INFO) 17753565b12SHong Zhang if (ci[pn] != 0) { 1789566063dSJacob Faibussowitsch PetscCall(PetscInfo(C, "Reallocs %" PetscInt_FMT "; Fill ratio: given %g needed %g.\n", nspacedouble, (double)fill, (double)afill)); 1799566063dSJacob Faibussowitsch PetscCall(PetscInfo(C, "Use MatPtAP(A,P,MatReuse,%g,&C) for best performance.\n", (double)afill)); 18053565b12SHong Zhang } else { 1819566063dSJacob Faibussowitsch PetscCall(PetscInfo(C, "Empty matrix product\n")); 18253565b12SHong Zhang } 183f79958ebSHong Zhang #endif 18453565b12SHong Zhang PetscFunctionReturn(0); 18553565b12SHong Zhang } 18653565b12SHong Zhang 187*9371c9d4SSatish Balay PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy(Mat A, Mat P, Mat C) { 18853565b12SHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 18953565b12SHong Zhang Mat_SeqAIJ *p = (Mat_SeqAIJ *)P->data; 19053565b12SHong Zhang Mat_SeqAIJ *c = (Mat_SeqAIJ *)C->data; 19153565b12SHong Zhang PetscInt *ai = a->i, *aj = a->j, *apj, *apjdense, *pi = p->i, *pj = p->j, *pJ = p->j, *pjj; 19253565b12SHong Zhang PetscInt *ci = c->i, *cj = c->j, *cjj; 19353565b12SHong Zhang PetscInt am = A->rmap->N, cn = C->cmap->N, cm = C->rmap->N; 19453565b12SHong Zhang PetscInt i, j, k, anzi, pnzi, apnzj, nextap, pnzj, prow, crow; 19565e4b4d4SStefano Zampini MatScalar *aa, *apa, *pa, *pA, *paj, *ca, *caj; 19653565b12SHong Zhang 19753565b12SHong Zhang PetscFunctionBegin; 198ae32dd66SHong Zhang /* Allocate temporary array for storage of one row of A*P (cn: non-scalable) */ 1999566063dSJacob Faibussowitsch PetscCall(PetscCalloc2(cn, &apa, cn, &apjdense)); 2009566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(cn, &apj)); 20165e4b4d4SStefano Zampini /* trigger CPU copies if needed and flag CPU mask for C */ 20265e4b4d4SStefano Zampini #if defined(PETSC_HAVE_DEVICE) 20365e4b4d4SStefano Zampini { 20465e4b4d4SStefano Zampini const PetscScalar *dummy; 2059566063dSJacob Faibussowitsch PetscCall(MatSeqAIJGetArrayRead(A, &dummy)); 2069566063dSJacob Faibussowitsch PetscCall(MatSeqAIJRestoreArrayRead(A, &dummy)); 2079566063dSJacob Faibussowitsch PetscCall(MatSeqAIJGetArrayRead(P, &dummy)); 2089566063dSJacob Faibussowitsch PetscCall(MatSeqAIJRestoreArrayRead(P, &dummy)); 20965e4b4d4SStefano Zampini if (C->offloadmask != PETSC_OFFLOAD_UNALLOCATED) C->offloadmask = PETSC_OFFLOAD_CPU; 21065e4b4d4SStefano Zampini } 21165e4b4d4SStefano Zampini #endif 21265e4b4d4SStefano Zampini aa = a->a; 21365e4b4d4SStefano Zampini pa = p->a; 21465e4b4d4SStefano Zampini pA = p->a; 21565e4b4d4SStefano Zampini ca = c->a; 21653565b12SHong Zhang 21753565b12SHong Zhang /* Clear old values in C */ 2189566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(ca, ci[cm])); 21953565b12SHong Zhang 22053565b12SHong Zhang for (i = 0; i < am; i++) { 22153565b12SHong Zhang /* Form sparse row of A*P */ 22253565b12SHong Zhang anzi = ai[i + 1] - ai[i]; 22353565b12SHong Zhang apnzj = 0; 22453565b12SHong Zhang for (j = 0; j < anzi; j++) { 22553565b12SHong Zhang prow = *aj++; 22653565b12SHong Zhang pnzj = pi[prow + 1] - pi[prow]; 22753565b12SHong Zhang pjj = pj + pi[prow]; 22853565b12SHong Zhang paj = pa + pi[prow]; 22953565b12SHong Zhang for (k = 0; k < pnzj; k++) { 23053565b12SHong Zhang if (!apjdense[pjj[k]]) { 23153565b12SHong Zhang apjdense[pjj[k]] = -1; 23253565b12SHong Zhang apj[apnzj++] = pjj[k]; 23353565b12SHong Zhang } 23453565b12SHong Zhang apa[pjj[k]] += (*aa) * paj[k]; 23553565b12SHong Zhang } 2369566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(2.0 * pnzj)); 23753565b12SHong Zhang aa++; 23853565b12SHong Zhang } 23953565b12SHong Zhang 24053565b12SHong Zhang /* Sort the j index array for quick sparse axpy. */ 24153565b12SHong Zhang /* Note: a array does not need sorting as it is in dense storage locations. */ 2429566063dSJacob Faibussowitsch PetscCall(PetscSortInt(apnzj, apj)); 24353565b12SHong Zhang 24453565b12SHong Zhang /* Compute P^T*A*P using outer product (P^T)[:,j]*(A*P)[j,:]. */ 24553565b12SHong Zhang pnzi = pi[i + 1] - pi[i]; 24653565b12SHong Zhang for (j = 0; j < pnzi; j++) { 24753565b12SHong Zhang nextap = 0; 24853565b12SHong Zhang crow = *pJ++; 24953565b12SHong Zhang cjj = cj + ci[crow]; 25053565b12SHong Zhang caj = ca + ci[crow]; 25153565b12SHong Zhang /* Perform sparse axpy operation. Note cjj includes apj. */ 25253565b12SHong Zhang for (k = 0; nextap < apnzj; k++) { 2536bdcaf15SBarry Smith 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*9371c9d4SSatish Balay if (cjj[k] == apj[nextap]) { caj[k] += (*pA) * apa[apj[nextap++]]; } 25553565b12SHong Zhang } 2569566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(2.0 * apnzj)); 25753565b12SHong Zhang pA++; 25853565b12SHong Zhang } 25953565b12SHong Zhang 26053565b12SHong Zhang /* Zero the current row info for A*P */ 26153565b12SHong Zhang for (j = 0; j < apnzj; j++) { 26253565b12SHong Zhang apa[apj[j]] = 0.; 26353565b12SHong Zhang apjdense[apj[j]] = 0; 26453565b12SHong Zhang } 26553565b12SHong Zhang } 26653565b12SHong Zhang 26753565b12SHong Zhang /* Assemble the final matrix and clean up */ 2689566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY)); 2699566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY)); 270ae32dd66SHong Zhang 2719566063dSJacob Faibussowitsch PetscCall(PetscFree2(apa, apjdense)); 2729566063dSJacob Faibussowitsch PetscCall(PetscFree(apj)); 27353565b12SHong Zhang PetscFunctionReturn(0); 27453565b12SHong Zhang } 27553565b12SHong Zhang 276*9371c9d4SSatish Balay PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat A, Mat P, Mat C) { 2776718818eSStefano Zampini Mat_MatTransMatMult *atb; 278eb9c0419SKris Buschelman 279eb9c0419SKris Buschelman PetscFunctionBegin; 2806718818eSStefano Zampini MatCheckProduct(C, 3); 2816718818eSStefano Zampini atb = (Mat_MatTransMatMult *)C->product->data; 28228b400f6SJacob Faibussowitsch PetscCheck(atb, PetscObjectComm((PetscObject)C), PETSC_ERR_PLIB, "Missing data structure"); 283acd337a6SBarry Smith PetscCall(MatTranspose(P, MAT_REUSE_MATRIX, &atb->At)); 28428b400f6SJacob Faibussowitsch PetscCheck(C->ops->matmultnumeric, PetscObjectComm((PetscObject)C), PETSC_ERR_PLIB, "Missing numeric operation"); 2856718818eSStefano Zampini /* when using rap, MatMatMatMultSymbolic used a different data */ 2866718818eSStefano Zampini if (atb->data) C->product->data = atb->data; 2879566063dSJacob Faibussowitsch PetscCall((*C->ops->matmatmultnumeric)(atb->At, A, P, C)); 2886718818eSStefano Zampini C->product->data = atb; 28953565b12SHong Zhang PetscFunctionReturn(0); 29053565b12SHong Zhang } 291