1eb9c0419SKris Buschelman /*
242c57489SHong Zhang Defines projective product routines where A is a SeqAIJ matrix
3eb9c0419SKris Buschelman C = P^T * A * P
4eb9c0419SKris Buschelman */
5eb9c0419SKris Buschelman
6c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h> /*I "petscmat.h" I*/
7c6db04a5SJed Brown #include <../src/mat/utils/freespace.h>
8c6db04a5SJed Brown #include <petscbt.h>
98563dfccSBarry Smith #include <petsctime.h>
10eb9c0419SKris Buschelman
116f231fbdSstefano_zampini #if defined(PETSC_HAVE_HYPRE)
124222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatPtAPSymbolic_AIJ_AIJ_wHYPRE(Mat, Mat, PetscReal, Mat);
136f231fbdSstefano_zampini #endif
146f231fbdSstefano_zampini
MatProductSymbolic_PtAP_SeqAIJ_SeqAIJ(Mat C)15d71ae5a4SJacob Faibussowitsch PetscErrorCode MatProductSymbolic_PtAP_SeqAIJ_SeqAIJ(Mat C)
16d71ae5a4SJacob Faibussowitsch {
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;
303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
314222ddf1SHong Zhang }
324222ddf1SHong Zhang
334222ddf1SHong Zhang /* "rap" */
349566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(alg, "rap", &flg));
356718818eSStefano Zampini if (flg) {
36cc1eb50dSBarry Smith MatProductCtx_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;
46cc1eb50dSBarry Smith C->product->destroy = MatProductCtxDestroy_SeqAIJ_MatTransMatMult;
474222ddf1SHong Zhang C->ops->ptapnumeric = MatPtAPNumeric_SeqAIJ_SeqAIJ;
484222ddf1SHong Zhang C->ops->productnumeric = MatProductNumeric_PtAP;
493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
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));
573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
584222ddf1SHong Zhang }
596f231fbdSstefano_zampini #endif
604222ddf1SHong Zhang
614222ddf1SHong Zhang SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "MatProductType is not supported");
627ba1a0bfSKris Buschelman }
637ba1a0bfSKris Buschelman
MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(Mat A,Mat P,PetscReal fill,Mat C)64d71ae5a4SJacob Faibussowitsch PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(Mat A, Mat P, PetscReal fill, Mat C)
65d71ae5a4SJacob Faibussowitsch {
660298fd71SBarry Smith PetscFreeSpaceList free_space = NULL, current_space = NULL;
67d20bfe6fSHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data, *p = (Mat_SeqAIJ *)P->data, *c;
6853565b12SHong Zhang PetscInt *pti, *ptj, *ptJ, *ai = a->i, *aj = a->j, *ajj, *pi = p->i, *pj = p->j, *pjj;
6953565b12SHong Zhang PetscInt *ci, *cj, *ptadenserow, *ptasparserow, *ptaj, nspacedouble = 0;
707d69fa63SHong Zhang PetscInt an = A->cmap->N, am = A->rmap->N, pn = P->cmap->N, pm = P->rmap->N;
7153565b12SHong Zhang PetscInt i, j, k, ptnzi, arow, anzj, ptanzi, prow, pnzj, cnzi, nlnk, *lnk;
72d20bfe6fSHong Zhang MatScalar *ca;
7353565b12SHong Zhang PetscBT lnkbt;
747d69fa63SHong Zhang PetscReal afill;
75eb9c0419SKris Buschelman
76eb9c0419SKris Buschelman PetscFunctionBegin;
7753565b12SHong Zhang /* Get ij structure of P^T */
789566063dSJacob Faibussowitsch PetscCall(MatGetSymbolicTranspose_SeqAIJ(P, &pti, &ptj));
79d20bfe6fSHong Zhang ptJ = ptj;
80eb9c0419SKris Buschelman
81d20bfe6fSHong Zhang /* Allocate ci array, arrays for fill computation and */
82d20bfe6fSHong Zhang /* free space for accumulating nonzero column info */
839566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(pn + 1, &ci));
84d20bfe6fSHong Zhang ci[0] = 0;
85eb9c0419SKris Buschelman
869566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(2 * an + 1, &ptadenserow));
87d20bfe6fSHong Zhang ptasparserow = ptadenserow + an;
8855a3bba9SHong Zhang
8955a3bba9SHong Zhang /* create and initialize a linked list */
9055a3bba9SHong Zhang nlnk = pn + 1;
919566063dSJacob Faibussowitsch PetscCall(PetscLLCreate(pn, pn, nlnk, lnk, lnkbt));
92eb9c0419SKris Buschelman
937d69fa63SHong Zhang /* Set initial free space to be fill*(nnz(A)+ nnz(P)) */
949566063dSJacob Faibussowitsch PetscCall(PetscFreeSpaceGet(PetscRealIntMultTruncate(fill, PetscIntSumTruncate(ai[am], pi[pm])), &free_space));
95d20bfe6fSHong Zhang current_space = free_space;
96d20bfe6fSHong Zhang
97d20bfe6fSHong Zhang /* Determine symbolic info for each row of C: */
98d20bfe6fSHong Zhang for (i = 0; i < pn; i++) {
99d20bfe6fSHong Zhang ptnzi = pti[i + 1] - pti[i];
100d20bfe6fSHong Zhang ptanzi = 0;
101d20bfe6fSHong Zhang /* Determine symbolic row of PtA: */
102d20bfe6fSHong Zhang for (j = 0; j < ptnzi; j++) {
103d20bfe6fSHong Zhang arow = *ptJ++;
104d20bfe6fSHong Zhang anzj = ai[arow + 1] - ai[arow];
105d20bfe6fSHong Zhang ajj = aj + ai[arow];
106d20bfe6fSHong Zhang for (k = 0; k < anzj; k++) {
107d20bfe6fSHong Zhang if (!ptadenserow[ajj[k]]) {
108d20bfe6fSHong Zhang ptadenserow[ajj[k]] = -1;
109d20bfe6fSHong Zhang ptasparserow[ptanzi++] = ajj[k];
110d20bfe6fSHong Zhang }
111d20bfe6fSHong Zhang }
112d20bfe6fSHong Zhang }
113d20bfe6fSHong Zhang /* Using symbolic info for row of PtA, determine symbolic info for row of C: */
114d20bfe6fSHong Zhang ptaj = ptasparserow;
115d20bfe6fSHong Zhang cnzi = 0;
116d20bfe6fSHong Zhang for (j = 0; j < ptanzi; j++) {
117d20bfe6fSHong Zhang prow = *ptaj++;
118d20bfe6fSHong Zhang pnzj = pi[prow + 1] - pi[prow];
119d20bfe6fSHong Zhang pjj = pj + pi[prow];
12055a3bba9SHong Zhang /* add non-zero cols of P into the sorted linked list lnk */
1219566063dSJacob Faibussowitsch PetscCall(PetscLLAddSorted(pnzj, pjj, pn, &nlnk, lnk, lnkbt));
12255a3bba9SHong Zhang cnzi += nlnk;
123d20bfe6fSHong Zhang }
124d20bfe6fSHong Zhang
125d20bfe6fSHong Zhang /* If free space is not available, make more free space */
126d20bfe6fSHong Zhang /* Double the amount of total space in the list */
127d20bfe6fSHong Zhang if (current_space->local_remaining < cnzi) {
1289566063dSJacob Faibussowitsch PetscCall(PetscFreeSpaceGet(PetscIntSumTruncate(cnzi, current_space->total_array_size), ¤t_space));
129f2b054eeSHong Zhang nspacedouble++;
130d20bfe6fSHong Zhang }
131d20bfe6fSHong Zhang
132d20bfe6fSHong Zhang /* Copy data into free space, and zero out denserows */
1339566063dSJacob Faibussowitsch PetscCall(PetscLLClean(pn, pn, cnzi, lnk, current_space->array, lnkbt));
1342205254eSKarl Rupp
135d20bfe6fSHong Zhang current_space->array += cnzi;
136d20bfe6fSHong Zhang current_space->local_used += cnzi;
137d20bfe6fSHong Zhang current_space->local_remaining -= cnzi;
138d20bfe6fSHong Zhang
1392205254eSKarl Rupp for (j = 0; j < ptanzi; j++) ptadenserow[ptasparserow[j]] = 0;
1402205254eSKarl Rupp
141d20bfe6fSHong Zhang /* Aside: Perhaps we should save the pta info for the numerical factorization. */
142d20bfe6fSHong Zhang /* For now, we will recompute what is needed. */
143d20bfe6fSHong Zhang ci[i + 1] = ci[i] + cnzi;
144d20bfe6fSHong Zhang }
145d20bfe6fSHong Zhang /* nnz is now stored in ci[ptm], column indices are in the list of free space */
146d20bfe6fSHong Zhang /* Allocate space for cj, initialize cj, and */
147d20bfe6fSHong Zhang /* destroy list of free space and other temporary array(s) */
14884648c2dSPierre Jolivet PetscCall(PetscMalloc1(ci[pn], &cj));
1499566063dSJacob Faibussowitsch PetscCall(PetscFreeSpaceContiguous(&free_space, cj));
1509566063dSJacob Faibussowitsch PetscCall(PetscFree(ptadenserow));
1519566063dSJacob Faibussowitsch PetscCall(PetscLLDestroy(lnk, lnkbt));
152d20bfe6fSHong Zhang
15384648c2dSPierre Jolivet PetscCall(PetscCalloc1(ci[pn], &ca));
15453565b12SHong Zhang
15553565b12SHong Zhang /* put together the new matrix */
1569566063dSJacob Faibussowitsch PetscCall(MatSetSeqAIJWithArrays_private(PetscObjectComm((PetscObject)A), pn, pn, ci, cj, ca, ((PetscObject)A)->type_name, C));
15758b7e2c1SStefano Zampini PetscCall(MatSetBlockSizes(C, P->cmap->bs, P->cmap->bs));
158ae32dd66SHong Zhang
15953565b12SHong Zhang /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
16053565b12SHong Zhang /* Since these are PETSc arrays, change flags to free them as necessary. */
1614222ddf1SHong Zhang c = (Mat_SeqAIJ *)((C)->data);
16253565b12SHong Zhang c->free_a = PETSC_TRUE;
16353565b12SHong Zhang c->free_ij = PETSC_TRUE;
16453565b12SHong Zhang c->nonew = 0;
1656718818eSStefano Zampini
1664222ddf1SHong Zhang C->ops->ptapnumeric = MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy;
16753565b12SHong Zhang
1687d69fa63SHong Zhang /* set MatInfo */
1697d69fa63SHong Zhang afill = (PetscReal)ci[pn] / (ai[am] + pi[pm] + 1.e-5);
1707d69fa63SHong Zhang if (afill < 1.0) afill = 1.0;
1714222ddf1SHong Zhang C->info.mallocs = nspacedouble;
1724222ddf1SHong Zhang C->info.fill_ratio_given = fill;
1734222ddf1SHong Zhang C->info.fill_ratio_needed = afill;
1747d69fa63SHong Zhang
17553565b12SHong Zhang /* Clean up. */
1769566063dSJacob Faibussowitsch PetscCall(MatRestoreSymbolicTranspose_SeqAIJ(P, &pti, &ptj));
17753565b12SHong Zhang #if defined(PETSC_USE_INFO)
17853565b12SHong Zhang if (ci[pn] != 0) {
1799566063dSJacob Faibussowitsch PetscCall(PetscInfo(C, "Reallocs %" PetscInt_FMT "; Fill ratio: given %g needed %g.\n", nspacedouble, (double)fill, (double)afill));
1809566063dSJacob Faibussowitsch PetscCall(PetscInfo(C, "Use MatPtAP(A,P,MatReuse,%g,&C) for best performance.\n", (double)afill));
18153565b12SHong Zhang } else {
1829566063dSJacob Faibussowitsch PetscCall(PetscInfo(C, "Empty matrix product\n"));
18353565b12SHong Zhang }
184f79958ebSHong Zhang #endif
1853ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
18653565b12SHong Zhang }
18753565b12SHong Zhang
MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy(Mat A,Mat P,Mat C)188d71ae5a4SJacob Faibussowitsch PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy(Mat A, Mat P, Mat C)
189d71ae5a4SJacob Faibussowitsch {
19053565b12SHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data;
19153565b12SHong Zhang Mat_SeqAIJ *p = (Mat_SeqAIJ *)P->data;
19253565b12SHong Zhang Mat_SeqAIJ *c = (Mat_SeqAIJ *)C->data;
19353565b12SHong Zhang PetscInt *ai = a->i, *aj = a->j, *apj, *apjdense, *pi = p->i, *pj = p->j, *pJ = p->j, *pjj;
19453565b12SHong Zhang PetscInt *ci = c->i, *cj = c->j, *cjj;
19553565b12SHong Zhang PetscInt am = A->rmap->N, cn = C->cmap->N, cm = C->rmap->N;
19653565b12SHong Zhang PetscInt i, j, k, anzi, pnzi, apnzj, nextap, pnzj, prow, crow;
19765e4b4d4SStefano Zampini MatScalar *aa, *apa, *pa, *pA, *paj, *ca, *caj;
19853565b12SHong Zhang
19953565b12SHong Zhang PetscFunctionBegin;
200ae32dd66SHong Zhang /* Allocate temporary array for storage of one row of A*P (cn: non-scalable) */
2019566063dSJacob Faibussowitsch PetscCall(PetscCalloc2(cn, &apa, cn, &apjdense));
2029566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(cn, &apj));
20365e4b4d4SStefano Zampini /* trigger CPU copies if needed and flag CPU mask for C */
20465e4b4d4SStefano Zampini #if defined(PETSC_HAVE_DEVICE)
20565e4b4d4SStefano Zampini {
206*2a8381b2SBarry Smith const PetscScalar *unused;
207*2a8381b2SBarry Smith PetscCall(MatSeqAIJGetArrayRead(A, &unused));
208*2a8381b2SBarry Smith PetscCall(MatSeqAIJRestoreArrayRead(A, &unused));
209*2a8381b2SBarry Smith PetscCall(MatSeqAIJGetArrayRead(P, &unused));
210*2a8381b2SBarry Smith PetscCall(MatSeqAIJRestoreArrayRead(P, &unused));
21165e4b4d4SStefano Zampini if (C->offloadmask != PETSC_OFFLOAD_UNALLOCATED) C->offloadmask = PETSC_OFFLOAD_CPU;
21265e4b4d4SStefano Zampini }
21365e4b4d4SStefano Zampini #endif
21465e4b4d4SStefano Zampini aa = a->a;
21565e4b4d4SStefano Zampini pa = p->a;
21665e4b4d4SStefano Zampini pA = p->a;
21765e4b4d4SStefano Zampini ca = c->a;
21853565b12SHong Zhang
21953565b12SHong Zhang /* Clear old values in C */
2209566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(ca, ci[cm]));
22153565b12SHong Zhang
22253565b12SHong Zhang for (i = 0; i < am; i++) {
22353565b12SHong Zhang /* Form sparse row of A*P */
22453565b12SHong Zhang anzi = ai[i + 1] - ai[i];
22553565b12SHong Zhang apnzj = 0;
22653565b12SHong Zhang for (j = 0; j < anzi; j++) {
22753565b12SHong Zhang prow = *aj++;
22853565b12SHong Zhang pnzj = pi[prow + 1] - pi[prow];
22953565b12SHong Zhang pjj = pj + pi[prow];
23053565b12SHong Zhang paj = pa + pi[prow];
23153565b12SHong Zhang for (k = 0; k < pnzj; k++) {
23253565b12SHong Zhang if (!apjdense[pjj[k]]) {
23353565b12SHong Zhang apjdense[pjj[k]] = -1;
23453565b12SHong Zhang apj[apnzj++] = pjj[k];
23553565b12SHong Zhang }
23653565b12SHong Zhang apa[pjj[k]] += (*aa) * paj[k];
23753565b12SHong Zhang }
2389566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(2.0 * pnzj));
23953565b12SHong Zhang aa++;
24053565b12SHong Zhang }
24153565b12SHong Zhang
24253565b12SHong Zhang /* Sort the j index array for quick sparse axpy. */
24353565b12SHong Zhang /* Note: a array does not need sorting as it is in dense storage locations. */
2449566063dSJacob Faibussowitsch PetscCall(PetscSortInt(apnzj, apj));
24553565b12SHong Zhang
24653565b12SHong Zhang /* Compute P^T*A*P using outer product (P^T)[:,j]*(A*P)[j,:]. */
24753565b12SHong Zhang pnzi = pi[i + 1] - pi[i];
24853565b12SHong Zhang for (j = 0; j < pnzi; j++) {
24953565b12SHong Zhang nextap = 0;
25053565b12SHong Zhang crow = *pJ++;
25153565b12SHong Zhang cjj = cj + ci[crow];
25253565b12SHong Zhang caj = ca + ci[crow];
25353565b12SHong Zhang /* Perform sparse axpy operation. Note cjj includes apj. */
25453565b12SHong Zhang for (k = 0; nextap < apnzj; k++) {
2556bdcaf15SBarry 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);
256ad540459SPierre Jolivet if (cjj[k] == apj[nextap]) caj[k] += (*pA) * apa[apj[nextap++]];
25753565b12SHong Zhang }
2589566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(2.0 * apnzj));
25953565b12SHong Zhang pA++;
26053565b12SHong Zhang }
26153565b12SHong Zhang
26253565b12SHong Zhang /* Zero the current row info for A*P */
26353565b12SHong Zhang for (j = 0; j < apnzj; j++) {
26453565b12SHong Zhang apa[apj[j]] = 0.;
26553565b12SHong Zhang apjdense[apj[j]] = 0;
26653565b12SHong Zhang }
26753565b12SHong Zhang }
26853565b12SHong Zhang
26953565b12SHong Zhang /* Assemble the final matrix and clean up */
2709566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
2719566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
272ae32dd66SHong Zhang
2739566063dSJacob Faibussowitsch PetscCall(PetscFree2(apa, apjdense));
2749566063dSJacob Faibussowitsch PetscCall(PetscFree(apj));
2753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
27653565b12SHong Zhang }
27753565b12SHong Zhang
MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat A,Mat P,Mat C)278d71ae5a4SJacob Faibussowitsch PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat A, Mat P, Mat C)
279d71ae5a4SJacob Faibussowitsch {
280cc1eb50dSBarry Smith MatProductCtx_MatTransMatMult *atb;
281eb9c0419SKris Buschelman
282eb9c0419SKris Buschelman PetscFunctionBegin;
2836718818eSStefano Zampini MatCheckProduct(C, 3);
284cc1eb50dSBarry Smith atb = (MatProductCtx_MatTransMatMult *)C->product->data;
28528b400f6SJacob Faibussowitsch PetscCheck(atb, PetscObjectComm((PetscObject)C), PETSC_ERR_PLIB, "Missing data structure");
286acd337a6SBarry Smith PetscCall(MatTranspose(P, MAT_REUSE_MATRIX, &atb->At));
28728b400f6SJacob Faibussowitsch PetscCheck(C->ops->matmultnumeric, PetscObjectComm((PetscObject)C), PETSC_ERR_PLIB, "Missing numeric operation");
2886718818eSStefano Zampini /* when using rap, MatMatMatMultSymbolic used a different data */
2896718818eSStefano Zampini if (atb->data) C->product->data = atb->data;
2909566063dSJacob Faibussowitsch PetscCall((*C->ops->matmatmultnumeric)(atb->At, A, P, C));
2916718818eSStefano Zampini C->product->data = atb;
2923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
29353565b12SHong Zhang }
294