xref: /petsc/src/mat/impls/aij/seq/matptap.c (revision 3ba1676111f5c958fe6c2729b46ca4d523958bb3)
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 
16d71ae5a4SJacob Faibussowitsch PetscErrorCode MatProductSymbolic_PtAP_SeqAIJ_SeqAIJ(Mat C)
17d71ae5a4SJacob Faibussowitsch {
184222ddf1SHong Zhang   Mat_Product        *product = C->product;
194222ddf1SHong Zhang   Mat                 A = product->A, P = product->B;
204222ddf1SHong Zhang   MatProductAlgorithm alg  = product->alg;
214222ddf1SHong Zhang   PetscReal           fill = product->fill;
224222ddf1SHong Zhang   PetscBool           flg;
238fa4b5a6SHong Zhang   Mat                 Pt;
247ba1a0bfSKris Buschelman 
257ba1a0bfSKris Buschelman   PetscFunctionBegin;
264222ddf1SHong Zhang   /* "scalable" */
279566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(alg, "scalable", &flg));
284222ddf1SHong Zhang   if (flg) {
299566063dSJacob Faibussowitsch     PetscCall(MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(A, P, fill, C));
304222ddf1SHong Zhang     C->ops->productnumeric = MatProductNumeric_PtAP;
31*3ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
324222ddf1SHong Zhang   }
334222ddf1SHong Zhang 
344222ddf1SHong Zhang   /* "rap" */
359566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(alg, "rap", &flg));
366718818eSStefano Zampini   if (flg) {
376718818eSStefano Zampini     Mat_MatTransMatMult *atb;
386718818eSStefano Zampini 
399566063dSJacob Faibussowitsch     PetscCall(PetscNew(&atb));
40acd337a6SBarry Smith     PetscCall(MatTranspose(P, MAT_INITIAL_MATRIX, &Pt));
419566063dSJacob Faibussowitsch     PetscCall(MatMatMatMultSymbolic_SeqAIJ_SeqAIJ_SeqAIJ(Pt, A, P, fill, C));
428fa4b5a6SHong Zhang 
438fa4b5a6SHong Zhang     atb->At                = Pt;
446718818eSStefano Zampini     atb->data              = C->product->data;
456718818eSStefano Zampini     atb->destroy           = C->product->destroy;
466718818eSStefano Zampini     C->product->data       = atb;
476718818eSStefano Zampini     C->product->destroy    = MatDestroy_SeqAIJ_MatTransMatMult;
484222ddf1SHong Zhang     C->ops->ptapnumeric    = MatPtAPNumeric_SeqAIJ_SeqAIJ;
494222ddf1SHong Zhang     C->ops->productnumeric = MatProductNumeric_PtAP;
50*3ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
514222ddf1SHong Zhang   }
524222ddf1SHong Zhang 
534222ddf1SHong Zhang   /* hypre */
546f231fbdSstefano_zampini #if defined(PETSC_HAVE_HYPRE)
559566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(alg, "hypre", &flg));
564222ddf1SHong Zhang   if (flg) {
579566063dSJacob Faibussowitsch     PetscCall(MatPtAPSymbolic_AIJ_AIJ_wHYPRE(A, P, fill, C));
58*3ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
594222ddf1SHong Zhang   }
606f231fbdSstefano_zampini #endif
614222ddf1SHong Zhang 
624222ddf1SHong Zhang   SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "MatProductType is not supported");
637ba1a0bfSKris Buschelman }
647ba1a0bfSKris Buschelman 
65d71ae5a4SJacob Faibussowitsch PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(Mat A, Mat P, PetscReal fill, Mat C)
66d71ae5a4SJacob Faibussowitsch {
670298fd71SBarry Smith   PetscFreeSpaceList free_space = NULL, current_space = NULL;
68d20bfe6fSHong Zhang   Mat_SeqAIJ        *a = (Mat_SeqAIJ *)A->data, *p = (Mat_SeqAIJ *)P->data, *c;
6953565b12SHong Zhang   PetscInt          *pti, *ptj, *ptJ, *ai = a->i, *aj = a->j, *ajj, *pi = p->i, *pj = p->j, *pjj;
7053565b12SHong Zhang   PetscInt          *ci, *cj, *ptadenserow, *ptasparserow, *ptaj, nspacedouble = 0;
717d69fa63SHong Zhang   PetscInt           an = A->cmap->N, am = A->rmap->N, pn = P->cmap->N, pm = P->rmap->N;
7253565b12SHong Zhang   PetscInt           i, j, k, ptnzi, arow, anzj, ptanzi, prow, pnzj, cnzi, nlnk, *lnk;
73d20bfe6fSHong Zhang   MatScalar         *ca;
7453565b12SHong Zhang   PetscBT            lnkbt;
757d69fa63SHong Zhang   PetscReal          afill;
76eb9c0419SKris Buschelman 
77eb9c0419SKris Buschelman   PetscFunctionBegin;
7853565b12SHong Zhang   /* Get ij structure of P^T */
799566063dSJacob Faibussowitsch   PetscCall(MatGetSymbolicTranspose_SeqAIJ(P, &pti, &ptj));
80d20bfe6fSHong Zhang   ptJ = ptj;
81eb9c0419SKris Buschelman 
82d20bfe6fSHong Zhang   /* Allocate ci array, arrays for fill computation and */
83d20bfe6fSHong Zhang   /* free space for accumulating nonzero column info */
849566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(pn + 1, &ci));
85d20bfe6fSHong Zhang   ci[0] = 0;
86eb9c0419SKris Buschelman 
879566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(2 * an + 1, &ptadenserow));
88d20bfe6fSHong Zhang   ptasparserow = ptadenserow + an;
8955a3bba9SHong Zhang 
9055a3bba9SHong Zhang   /* create and initialize a linked list */
9155a3bba9SHong Zhang   nlnk = pn + 1;
929566063dSJacob Faibussowitsch   PetscCall(PetscLLCreate(pn, pn, nlnk, lnk, lnkbt));
93eb9c0419SKris Buschelman 
947d69fa63SHong Zhang   /* Set initial free space to be fill*(nnz(A)+ nnz(P)) */
959566063dSJacob Faibussowitsch   PetscCall(PetscFreeSpaceGet(PetscRealIntMultTruncate(fill, PetscIntSumTruncate(ai[am], pi[pm])), &free_space));
96d20bfe6fSHong Zhang   current_space = free_space;
97d20bfe6fSHong Zhang 
98d20bfe6fSHong Zhang   /* Determine symbolic info for each row of C: */
99d20bfe6fSHong Zhang   for (i = 0; i < pn; i++) {
100d20bfe6fSHong Zhang     ptnzi  = pti[i + 1] - pti[i];
101d20bfe6fSHong Zhang     ptanzi = 0;
102d20bfe6fSHong Zhang     /* Determine symbolic row of PtA: */
103d20bfe6fSHong Zhang     for (j = 0; j < ptnzi; j++) {
104d20bfe6fSHong Zhang       arow = *ptJ++;
105d20bfe6fSHong Zhang       anzj = ai[arow + 1] - ai[arow];
106d20bfe6fSHong Zhang       ajj  = aj + ai[arow];
107d20bfe6fSHong Zhang       for (k = 0; k < anzj; k++) {
108d20bfe6fSHong Zhang         if (!ptadenserow[ajj[k]]) {
109d20bfe6fSHong Zhang           ptadenserow[ajj[k]]    = -1;
110d20bfe6fSHong Zhang           ptasparserow[ptanzi++] = ajj[k];
111d20bfe6fSHong Zhang         }
112d20bfe6fSHong Zhang       }
113d20bfe6fSHong Zhang     }
114d20bfe6fSHong Zhang     /* Using symbolic info for row of PtA, determine symbolic info for row of C: */
115d20bfe6fSHong Zhang     ptaj = ptasparserow;
116d20bfe6fSHong Zhang     cnzi = 0;
117d20bfe6fSHong Zhang     for (j = 0; j < ptanzi; j++) {
118d20bfe6fSHong Zhang       prow = *ptaj++;
119d20bfe6fSHong Zhang       pnzj = pi[prow + 1] - pi[prow];
120d20bfe6fSHong Zhang       pjj  = pj + pi[prow];
12155a3bba9SHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
1229566063dSJacob Faibussowitsch       PetscCall(PetscLLAddSorted(pnzj, pjj, pn, &nlnk, lnk, lnkbt));
12355a3bba9SHong Zhang       cnzi += nlnk;
124d20bfe6fSHong Zhang     }
125d20bfe6fSHong Zhang 
126d20bfe6fSHong Zhang     /* If free space is not available, make more free space */
127d20bfe6fSHong Zhang     /* Double the amount of total space in the list */
128d20bfe6fSHong Zhang     if (current_space->local_remaining < cnzi) {
1299566063dSJacob Faibussowitsch       PetscCall(PetscFreeSpaceGet(PetscIntSumTruncate(cnzi, current_space->total_array_size), &current_space));
130f2b054eeSHong Zhang       nspacedouble++;
131d20bfe6fSHong Zhang     }
132d20bfe6fSHong Zhang 
133d20bfe6fSHong Zhang     /* Copy data into free space, and zero out denserows */
1349566063dSJacob Faibussowitsch     PetscCall(PetscLLClean(pn, pn, cnzi, lnk, current_space->array, lnkbt));
1352205254eSKarl Rupp 
136d20bfe6fSHong Zhang     current_space->array += cnzi;
137d20bfe6fSHong Zhang     current_space->local_used += cnzi;
138d20bfe6fSHong Zhang     current_space->local_remaining -= cnzi;
139d20bfe6fSHong Zhang 
1402205254eSKarl Rupp     for (j = 0; j < ptanzi; j++) ptadenserow[ptasparserow[j]] = 0;
1412205254eSKarl Rupp 
142d20bfe6fSHong Zhang     /* Aside: Perhaps we should save the pta info for the numerical factorization. */
143d20bfe6fSHong Zhang     /*        For now, we will recompute what is needed. */
144d20bfe6fSHong Zhang     ci[i + 1] = ci[i] + cnzi;
145d20bfe6fSHong Zhang   }
146d20bfe6fSHong Zhang   /* nnz is now stored in ci[ptm], column indices are in the list of free space */
147d20bfe6fSHong Zhang   /* Allocate space for cj, initialize cj, and */
148d20bfe6fSHong Zhang   /* destroy list of free space and other temporary array(s) */
1499566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(ci[pn] + 1, &cj));
1509566063dSJacob Faibussowitsch   PetscCall(PetscFreeSpaceContiguous(&free_space, cj));
1519566063dSJacob Faibussowitsch   PetscCall(PetscFree(ptadenserow));
1529566063dSJacob Faibussowitsch   PetscCall(PetscLLDestroy(lnk, lnkbt));
153d20bfe6fSHong Zhang 
1549566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(ci[pn] + 1, &ca));
15553565b12SHong Zhang 
15653565b12SHong Zhang   /* put together the new matrix */
1579566063dSJacob Faibussowitsch   PetscCall(MatSetSeqAIJWithArrays_private(PetscObjectComm((PetscObject)A), pn, pn, ci, cj, ca, ((PetscObject)A)->type_name, C));
1589566063dSJacob Faibussowitsch   PetscCall(MatSetBlockSizes(C, PetscAbs(P->cmap->bs), PetscAbs(P->cmap->bs)));
159ae32dd66SHong Zhang 
16053565b12SHong Zhang   /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */
16153565b12SHong Zhang   /* Since these are PETSc arrays, change flags to free them as necessary. */
1624222ddf1SHong Zhang   c          = (Mat_SeqAIJ *)((C)->data);
16353565b12SHong Zhang   c->free_a  = PETSC_TRUE;
16453565b12SHong Zhang   c->free_ij = PETSC_TRUE;
16553565b12SHong Zhang   c->nonew   = 0;
1666718818eSStefano Zampini 
1674222ddf1SHong Zhang   C->ops->ptapnumeric = MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy;
16853565b12SHong Zhang 
1697d69fa63SHong Zhang   /* set MatInfo */
1707d69fa63SHong Zhang   afill = (PetscReal)ci[pn] / (ai[am] + pi[pm] + 1.e-5);
1717d69fa63SHong Zhang   if (afill < 1.0) afill = 1.0;
1724222ddf1SHong Zhang   C->info.mallocs           = nspacedouble;
1734222ddf1SHong Zhang   C->info.fill_ratio_given  = fill;
1744222ddf1SHong Zhang   C->info.fill_ratio_needed = afill;
1757d69fa63SHong Zhang 
17653565b12SHong Zhang   /* Clean up. */
1779566063dSJacob Faibussowitsch   PetscCall(MatRestoreSymbolicTranspose_SeqAIJ(P, &pti, &ptj));
17853565b12SHong Zhang #if defined(PETSC_USE_INFO)
17953565b12SHong Zhang   if (ci[pn] != 0) {
1809566063dSJacob Faibussowitsch     PetscCall(PetscInfo(C, "Reallocs %" PetscInt_FMT "; Fill ratio: given %g needed %g.\n", nspacedouble, (double)fill, (double)afill));
1819566063dSJacob Faibussowitsch     PetscCall(PetscInfo(C, "Use MatPtAP(A,P,MatReuse,%g,&C) for best performance.\n", (double)afill));
18253565b12SHong Zhang   } else {
1839566063dSJacob Faibussowitsch     PetscCall(PetscInfo(C, "Empty matrix product\n"));
18453565b12SHong Zhang   }
185f79958ebSHong Zhang #endif
186*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
18753565b12SHong Zhang }
18853565b12SHong Zhang 
189d71ae5a4SJacob Faibussowitsch PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ_SparseAxpy(Mat A, Mat P, Mat C)
190d71ae5a4SJacob Faibussowitsch {
19153565b12SHong Zhang   Mat_SeqAIJ *a  = (Mat_SeqAIJ *)A->data;
19253565b12SHong Zhang   Mat_SeqAIJ *p  = (Mat_SeqAIJ *)P->data;
19353565b12SHong Zhang   Mat_SeqAIJ *c  = (Mat_SeqAIJ *)C->data;
19453565b12SHong Zhang   PetscInt   *ai = a->i, *aj = a->j, *apj, *apjdense, *pi = p->i, *pj = p->j, *pJ = p->j, *pjj;
19553565b12SHong Zhang   PetscInt   *ci = c->i, *cj = c->j, *cjj;
19653565b12SHong Zhang   PetscInt    am = A->rmap->N, cn = C->cmap->N, cm = C->rmap->N;
19753565b12SHong Zhang   PetscInt    i, j, k, anzi, pnzi, apnzj, nextap, pnzj, prow, crow;
19865e4b4d4SStefano Zampini   MatScalar  *aa, *apa, *pa, *pA, *paj, *ca, *caj;
19953565b12SHong Zhang 
20053565b12SHong Zhang   PetscFunctionBegin;
201ae32dd66SHong Zhang   /* Allocate temporary array for storage of one row of A*P (cn: non-scalable) */
2029566063dSJacob Faibussowitsch   PetscCall(PetscCalloc2(cn, &apa, cn, &apjdense));
2039566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(cn, &apj));
20465e4b4d4SStefano Zampini   /* trigger CPU copies if needed and flag CPU mask for C */
20565e4b4d4SStefano Zampini #if defined(PETSC_HAVE_DEVICE)
20665e4b4d4SStefano Zampini   {
20765e4b4d4SStefano Zampini     const PetscScalar *dummy;
2089566063dSJacob Faibussowitsch     PetscCall(MatSeqAIJGetArrayRead(A, &dummy));
2099566063dSJacob Faibussowitsch     PetscCall(MatSeqAIJRestoreArrayRead(A, &dummy));
2109566063dSJacob Faibussowitsch     PetscCall(MatSeqAIJGetArrayRead(P, &dummy));
2119566063dSJacob Faibussowitsch     PetscCall(MatSeqAIJRestoreArrayRead(P, &dummy));
21265e4b4d4SStefano Zampini     if (C->offloadmask != PETSC_OFFLOAD_UNALLOCATED) C->offloadmask = PETSC_OFFLOAD_CPU;
21365e4b4d4SStefano Zampini   }
21465e4b4d4SStefano Zampini #endif
21565e4b4d4SStefano Zampini   aa = a->a;
21665e4b4d4SStefano Zampini   pa = p->a;
21765e4b4d4SStefano Zampini   pA = p->a;
21865e4b4d4SStefano Zampini   ca = c->a;
21953565b12SHong Zhang 
22053565b12SHong Zhang   /* Clear old values in C */
2219566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(ca, ci[cm]));
22253565b12SHong Zhang 
22353565b12SHong Zhang   for (i = 0; i < am; i++) {
22453565b12SHong Zhang     /* Form sparse row of A*P */
22553565b12SHong Zhang     anzi  = ai[i + 1] - ai[i];
22653565b12SHong Zhang     apnzj = 0;
22753565b12SHong Zhang     for (j = 0; j < anzi; j++) {
22853565b12SHong Zhang       prow = *aj++;
22953565b12SHong Zhang       pnzj = pi[prow + 1] - pi[prow];
23053565b12SHong Zhang       pjj  = pj + pi[prow];
23153565b12SHong Zhang       paj  = pa + pi[prow];
23253565b12SHong Zhang       for (k = 0; k < pnzj; k++) {
23353565b12SHong Zhang         if (!apjdense[pjj[k]]) {
23453565b12SHong Zhang           apjdense[pjj[k]] = -1;
23553565b12SHong Zhang           apj[apnzj++]     = pjj[k];
23653565b12SHong Zhang         }
23753565b12SHong Zhang         apa[pjj[k]] += (*aa) * paj[k];
23853565b12SHong Zhang       }
2399566063dSJacob Faibussowitsch       PetscCall(PetscLogFlops(2.0 * pnzj));
24053565b12SHong Zhang       aa++;
24153565b12SHong Zhang     }
24253565b12SHong Zhang 
24353565b12SHong Zhang     /* Sort the j index array for quick sparse axpy. */
24453565b12SHong Zhang     /* Note: a array does not need sorting as it is in dense storage locations. */
2459566063dSJacob Faibussowitsch     PetscCall(PetscSortInt(apnzj, apj));
24653565b12SHong Zhang 
24753565b12SHong Zhang     /* Compute P^T*A*P using outer product (P^T)[:,j]*(A*P)[j,:]. */
24853565b12SHong Zhang     pnzi = pi[i + 1] - pi[i];
24953565b12SHong Zhang     for (j = 0; j < pnzi; j++) {
25053565b12SHong Zhang       nextap = 0;
25153565b12SHong Zhang       crow   = *pJ++;
25253565b12SHong Zhang       cjj    = cj + ci[crow];
25353565b12SHong Zhang       caj    = ca + ci[crow];
25453565b12SHong Zhang       /* Perform sparse axpy operation.  Note cjj includes apj. */
25553565b12SHong Zhang       for (k = 0; nextap < apnzj; k++) {
2566bdcaf15SBarry 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);
257ad540459SPierre Jolivet         if (cjj[k] == apj[nextap]) caj[k] += (*pA) * apa[apj[nextap++]];
25853565b12SHong Zhang       }
2599566063dSJacob Faibussowitsch       PetscCall(PetscLogFlops(2.0 * apnzj));
26053565b12SHong Zhang       pA++;
26153565b12SHong Zhang     }
26253565b12SHong Zhang 
26353565b12SHong Zhang     /* Zero the current row info for A*P */
26453565b12SHong Zhang     for (j = 0; j < apnzj; j++) {
26553565b12SHong Zhang       apa[apj[j]]      = 0.;
26653565b12SHong Zhang       apjdense[apj[j]] = 0;
26753565b12SHong Zhang     }
26853565b12SHong Zhang   }
26953565b12SHong Zhang 
27053565b12SHong Zhang   /* Assemble the final matrix and clean up */
2719566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
2729566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
273ae32dd66SHong Zhang 
2749566063dSJacob Faibussowitsch   PetscCall(PetscFree2(apa, apjdense));
2759566063dSJacob Faibussowitsch   PetscCall(PetscFree(apj));
276*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
27753565b12SHong Zhang }
27853565b12SHong Zhang 
279d71ae5a4SJacob Faibussowitsch PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat A, Mat P, Mat C)
280d71ae5a4SJacob Faibussowitsch {
2816718818eSStefano Zampini   Mat_MatTransMatMult *atb;
282eb9c0419SKris Buschelman 
283eb9c0419SKris Buschelman   PetscFunctionBegin;
2846718818eSStefano Zampini   MatCheckProduct(C, 3);
2856718818eSStefano Zampini   atb = (Mat_MatTransMatMult *)C->product->data;
28628b400f6SJacob Faibussowitsch   PetscCheck(atb, PetscObjectComm((PetscObject)C), PETSC_ERR_PLIB, "Missing data structure");
287acd337a6SBarry Smith   PetscCall(MatTranspose(P, MAT_REUSE_MATRIX, &atb->At));
28828b400f6SJacob Faibussowitsch   PetscCheck(C->ops->matmultnumeric, PetscObjectComm((PetscObject)C), PETSC_ERR_PLIB, "Missing numeric operation");
2896718818eSStefano Zampini   /* when using rap, MatMatMatMultSymbolic used a different data */
2906718818eSStefano Zampini   if (atb->data) C->product->data = atb->data;
2919566063dSJacob Faibussowitsch   PetscCall((*C->ops->matmatmultnumeric)(atb->At, A, P, C));
2926718818eSStefano Zampini   C->product->data = atb;
293*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
29453565b12SHong Zhang }
295