1807171c5SHong Zhang 23b1b9624SHong Zhang /* 33b1b9624SHong Zhang Defines projective product routines where A is a SeqAIJ matrix 43b1b9624SHong Zhang C = R * A * R^T 53b1b9624SHong Zhang */ 63b1b9624SHong Zhang 7807171c5SHong Zhang #include <../src/mat/impls/aij/seq/aij.h> 8807171c5SHong Zhang #include <../src/mat/utils/freespace.h> 9807171c5SHong Zhang #include <../src/mat/impls/dense/seq/dense.h> /*I "petscmat.h" I*/ 10807171c5SHong Zhang 116718818eSStefano Zampini PetscErrorCode MatDestroy_SeqAIJ_RARt(void *data) 12807171c5SHong Zhang { 136718818eSStefano Zampini Mat_RARt *rart = (Mat_RARt*)data; 14807171c5SHong Zhang 15807171c5SHong Zhang PetscFunctionBegin; 169566063dSJacob Faibussowitsch PetscCall(MatTransposeColoringDestroy(&rart->matcoloring)); 179566063dSJacob Faibussowitsch PetscCall(MatDestroy(&rart->Rt)); 189566063dSJacob Faibussowitsch PetscCall(MatDestroy(&rart->RARt)); 199566063dSJacob Faibussowitsch PetscCall(MatDestroy(&rart->ARt)); 209566063dSJacob Faibussowitsch PetscCall(PetscFree(rart->work)); 211baa6e33SBarry Smith if (rart->destroy) PetscCall((*rart->destroy)(rart->data)); 229566063dSJacob Faibussowitsch PetscCall(PetscFree(rart)); 23807171c5SHong Zhang PetscFunctionReturn(0); 24807171c5SHong Zhang } 25807171c5SHong Zhang 264222ddf1SHong Zhang PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_colorrart(Mat A,Mat R,PetscReal fill,Mat C) 27807171c5SHong Zhang { 28807171c5SHong Zhang Mat P; 29807171c5SHong Zhang PetscInt *rti,*rtj; 30807171c5SHong Zhang Mat_RARt *rart; 31335efc43SPeter Brune MatColoring coloring; 32807171c5SHong Zhang MatTransposeColoring matcoloring; 33807171c5SHong Zhang ISColoring iscoloring; 348d0a38e4SHong Zhang Mat Rt_dense,RARt_dense; 35807171c5SHong Zhang 36807171c5SHong Zhang PetscFunctionBegin; 376718818eSStefano Zampini MatCheckProduct(C,4); 3828b400f6SJacob Faibussowitsch PetscCheck(!C->product->data,PetscObjectComm((PetscObject)C),PETSC_ERR_PLIB,"Product data not empty"); 39807171c5SHong Zhang /* create symbolic P=Rt */ 409566063dSJacob Faibussowitsch PetscCall(MatGetSymbolicTranspose_SeqAIJ(R,&rti,&rtj)); 419566063dSJacob Faibussowitsch PetscCall(MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,R->cmap->n,R->rmap->n,rti,rtj,NULL,&P)); 42807171c5SHong Zhang 43807171c5SHong Zhang /* get symbolic C=Pt*A*P */ 449566063dSJacob Faibussowitsch PetscCall(MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(A,P,fill,C)); 459566063dSJacob Faibussowitsch PetscCall(MatSetBlockSizes(C,PetscAbs(R->rmap->bs),PetscAbs(R->rmap->bs))); 464222ddf1SHong Zhang C->ops->rartnumeric = MatRARtNumeric_SeqAIJ_SeqAIJ_colorrart; 47807171c5SHong Zhang 48807171c5SHong Zhang /* create a supporting struct */ 499566063dSJacob Faibussowitsch PetscCall(PetscNew(&rart)); 506718818eSStefano Zampini C->product->data = rart; 516718818eSStefano Zampini C->product->destroy = MatDestroy_SeqAIJ_RARt; 52807171c5SHong Zhang 5322e94b5dSHong Zhang /* ------ Use coloring ---------- */ 544222ddf1SHong Zhang /* inode causes memory problem */ 559566063dSJacob Faibussowitsch PetscCall(MatSetOption(C,MAT_USE_INODES,PETSC_FALSE)); 564d478ae7SHong Zhang 57807171c5SHong Zhang /* Create MatTransposeColoring from symbolic C=R*A*R^T */ 589566063dSJacob Faibussowitsch PetscCall(MatColoringCreate(C,&coloring)); 599566063dSJacob Faibussowitsch PetscCall(MatColoringSetDistance(coloring,2)); 609566063dSJacob Faibussowitsch PetscCall(MatColoringSetType(coloring,MATCOLORINGSL)); 619566063dSJacob Faibussowitsch PetscCall(MatColoringSetFromOptions(coloring)); 629566063dSJacob Faibussowitsch PetscCall(MatColoringApply(coloring,&iscoloring)); 639566063dSJacob Faibussowitsch PetscCall(MatColoringDestroy(&coloring)); 649566063dSJacob Faibussowitsch PetscCall(MatTransposeColoringCreate(C,iscoloring,&matcoloring)); 652205254eSKarl Rupp 66807171c5SHong Zhang rart->matcoloring = matcoloring; 679566063dSJacob Faibussowitsch PetscCall(ISColoringDestroy(&iscoloring)); 683b1b9624SHong Zhang 69807171c5SHong Zhang /* Create Rt_dense */ 709566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF,&Rt_dense)); 719566063dSJacob Faibussowitsch PetscCall(MatSetSizes(Rt_dense,A->cmap->n,matcoloring->ncolors,A->cmap->n,matcoloring->ncolors)); 729566063dSJacob Faibussowitsch PetscCall(MatSetType(Rt_dense,MATSEQDENSE)); 739566063dSJacob Faibussowitsch PetscCall(MatSeqDenseSetPreallocation(Rt_dense,NULL)); 742205254eSKarl Rupp 75807171c5SHong Zhang Rt_dense->assembled = PETSC_TRUE; 76807171c5SHong Zhang rart->Rt = Rt_dense; 77807171c5SHong Zhang 78807171c5SHong Zhang /* Create RARt_dense = R*A*Rt_dense */ 799566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF,&RARt_dense)); 809566063dSJacob Faibussowitsch PetscCall(MatSetSizes(RARt_dense,C->rmap->n,matcoloring->ncolors,C->rmap->n,matcoloring->ncolors)); 819566063dSJacob Faibussowitsch PetscCall(MatSetType(RARt_dense,MATSEQDENSE)); 829566063dSJacob Faibussowitsch PetscCall(MatSeqDenseSetPreallocation(RARt_dense,NULL)); 832205254eSKarl Rupp 84807171c5SHong Zhang rart->RARt = RARt_dense; 85807171c5SHong Zhang 86807171c5SHong Zhang /* Allocate work array to store columns of A*R^T used in MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense() */ 879566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(A->rmap->n*4,&rart->work)); 88807171c5SHong Zhang 89807171c5SHong Zhang /* clean up */ 909566063dSJacob Faibussowitsch PetscCall(MatRestoreSymbolicTranspose_SeqAIJ(R,&rti,&rtj)); 919566063dSJacob Faibussowitsch PetscCall(MatDestroy(&P)); 92807171c5SHong Zhang 93807171c5SHong Zhang #if defined(PETSC_USE_INFO) 941ce71dffSSatish Balay { 956718818eSStefano Zampini Mat_SeqAIJ *c = (Mat_SeqAIJ*)C->data; 96807171c5SHong Zhang PetscReal density = (PetscReal)(c->nz)/(RARt_dense->rmap->n*RARt_dense->cmap->n); 979566063dSJacob Faibussowitsch PetscCall(PetscInfo(C,"C=R*(A*Rt) via coloring C - use sparse-dense inner products\n")); 9863a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(C,"RARt_den %" PetscInt_FMT " %" PetscInt_FMT "; Rt %" PetscInt_FMT " %" PetscInt_FMT " (RARt->nz %" PetscInt_FMT ")/(m*ncolors)=%g\n",RARt_dense->rmap->n,RARt_dense->cmap->n,R->cmap->n,R->rmap->n,c->nz,(double)density)); 991ce71dffSSatish Balay } 100807171c5SHong Zhang #endif 101807171c5SHong Zhang PetscFunctionReturn(0); 102807171c5SHong Zhang } 103807171c5SHong Zhang 104807171c5SHong Zhang /* 105807171c5SHong Zhang RAB = R * A * B, R and A in seqaij format, B in dense format; 106807171c5SHong Zhang */ 107807171c5SHong Zhang PetscErrorCode MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense(Mat R,Mat A,Mat B,Mat RAB,PetscScalar *work) 108807171c5SHong Zhang { 109807171c5SHong Zhang Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*r=(Mat_SeqAIJ*)R->data; 1101683a169SBarry Smith PetscScalar r1,r2,r3,r4; 1111683a169SBarry Smith const PetscScalar *b,*b1,*b2,*b3,*b4; 112807171c5SHong Zhang MatScalar *aa,*ra; 113807171c5SHong Zhang PetscInt cn =B->cmap->n,bm=B->rmap->n,col,i,j,n,*ai=a->i,*aj,am=A->rmap->n; 114807171c5SHong Zhang PetscInt am2=2*am,am3=3*am,bm4=4*bm; 115807171c5SHong Zhang PetscScalar *d,*c,*c2,*c3,*c4; 116807171c5SHong Zhang PetscInt *rj,rm=R->rmap->n,dm=RAB->rmap->n,dn=RAB->cmap->n; 117807171c5SHong Zhang PetscInt rm2=2*rm,rm3=3*rm,colrm; 118807171c5SHong Zhang 119807171c5SHong Zhang PetscFunctionBegin; 120807171c5SHong Zhang if (!dm || !dn) PetscFunctionReturn(0); 12108401ef6SPierre Jolivet PetscCheck(bm == A->cmap->n,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Number columns in A %" PetscInt_FMT " not equal rows in B %" PetscInt_FMT,A->cmap->n,bm); 12208401ef6SPierre Jolivet PetscCheck(am == R->cmap->n,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Number columns in R %" PetscInt_FMT " not equal rows in A %" PetscInt_FMT,R->cmap->n,am); 12308401ef6SPierre Jolivet PetscCheck(R->rmap->n == RAB->rmap->n,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Number rows in RAB %" PetscInt_FMT " not equal rows in R %" PetscInt_FMT,RAB->rmap->n,R->rmap->n); 12408401ef6SPierre Jolivet PetscCheck(B->cmap->n == RAB->cmap->n,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Number columns in RAB %" PetscInt_FMT " not equal columns in B %" PetscInt_FMT,RAB->cmap->n,B->cmap->n); 125807171c5SHong Zhang 126274010c0SHong Zhang { /* 127274010c0SHong Zhang This approach is not as good as original ones (will be removed later), but it reveals that 128c4762a1bSJed Brown AB_den=A*B takes almost all execution time in R*A*B for src/ksp/ksp/tutorials/ex56.c 129274010c0SHong Zhang */ 130274010c0SHong Zhang PetscBool via_matmatmult=PETSC_FALSE; 1319566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-matrart_via_matmatmult",&via_matmatmult,NULL)); 132274010c0SHong Zhang if (via_matmatmult) { 1334222ddf1SHong Zhang Mat AB_den = NULL; 1349566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A),&AB_den)); 1359566063dSJacob Faibussowitsch PetscCall(MatMatMultSymbolic_SeqAIJ_SeqDense(A,B,0.0,AB_den)); 1369566063dSJacob Faibussowitsch PetscCall(MatMatMultNumeric_SeqAIJ_SeqDense(A,B,AB_den)); 1379566063dSJacob Faibussowitsch PetscCall(MatMatMultNumeric_SeqAIJ_SeqDense(R,AB_den,RAB)); 1389566063dSJacob Faibussowitsch PetscCall(MatDestroy(&AB_den)); 139274010c0SHong Zhang PetscFunctionReturn(0); 140274010c0SHong Zhang } 141274010c0SHong Zhang } 142274010c0SHong Zhang 1439566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayRead(B,&b)); 1449566063dSJacob Faibussowitsch PetscCall(MatDenseGetArray(RAB,&d)); 145807171c5SHong Zhang b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm; 146807171c5SHong Zhang c = work; c2 = c + am; c3 = c2 + am; c4 = c3 + am; 147807171c5SHong Zhang for (col=0; col<cn-4; col += 4) { /* over columns of C */ 148807171c5SHong Zhang for (i=0; i<am; i++) { /* over rows of A in those columns */ 149807171c5SHong Zhang r1 = r2 = r3 = r4 = 0.0; 150807171c5SHong Zhang n = ai[i+1] - ai[i]; 151807171c5SHong Zhang aj = a->j + ai[i]; 152807171c5SHong Zhang aa = a->a + ai[i]; 153807171c5SHong Zhang for (j=0; j<n; j++) { 154807171c5SHong Zhang r1 += (*aa)*b1[*aj]; 155807171c5SHong Zhang r2 += (*aa)*b2[*aj]; 156807171c5SHong Zhang r3 += (*aa)*b3[*aj]; 157807171c5SHong Zhang r4 += (*aa++)*b4[*aj++]; 158807171c5SHong Zhang } 159807171c5SHong Zhang c[i] = r1; 160807171c5SHong Zhang c[am + i] = r2; 161807171c5SHong Zhang c[am2 + i] = r3; 162807171c5SHong Zhang c[am3 + i] = r4; 163807171c5SHong Zhang } 164807171c5SHong Zhang b1 += bm4; 165807171c5SHong Zhang b2 += bm4; 166807171c5SHong Zhang b3 += bm4; 167807171c5SHong Zhang b4 += bm4; 168807171c5SHong Zhang 169807171c5SHong Zhang /* RAB[:,col] = R*C[:,col] */ 170807171c5SHong Zhang colrm = col*rm; 171807171c5SHong Zhang for (i=0; i<rm; i++) { /* over rows of R in those columns */ 172807171c5SHong Zhang r1 = r2 = r3 = r4 = 0.0; 173807171c5SHong Zhang n = r->i[i+1] - r->i[i]; 174807171c5SHong Zhang rj = r->j + r->i[i]; 175807171c5SHong Zhang ra = r->a + r->i[i]; 176807171c5SHong Zhang for (j=0; j<n; j++) { 177807171c5SHong Zhang r1 += (*ra)*c[*rj]; 178807171c5SHong Zhang r2 += (*ra)*c2[*rj]; 179807171c5SHong Zhang r3 += (*ra)*c3[*rj]; 180807171c5SHong Zhang r4 += (*ra++)*c4[*rj++]; 181807171c5SHong Zhang } 182807171c5SHong Zhang d[colrm + i] = r1; 183807171c5SHong Zhang d[colrm + rm + i] = r2; 184807171c5SHong Zhang d[colrm + rm2 + i] = r3; 185807171c5SHong Zhang d[colrm + rm3 + i] = r4; 186807171c5SHong Zhang } 187807171c5SHong Zhang } 188807171c5SHong Zhang for (; col<cn; col++) { /* over extra columns of C */ 189807171c5SHong Zhang for (i=0; i<am; i++) { /* over rows of A in those columns */ 190807171c5SHong Zhang r1 = 0.0; 191807171c5SHong Zhang n = a->i[i+1] - a->i[i]; 192807171c5SHong Zhang aj = a->j + a->i[i]; 193807171c5SHong Zhang aa = a->a + a->i[i]; 194807171c5SHong Zhang for (j=0; j<n; j++) { 195807171c5SHong Zhang r1 += (*aa++)*b1[*aj++]; 196807171c5SHong Zhang } 197807171c5SHong Zhang c[i] = r1; 198807171c5SHong Zhang } 199807171c5SHong Zhang b1 += bm; 200807171c5SHong Zhang 201807171c5SHong Zhang for (i=0; i<rm; i++) { /* over rows of R in those columns */ 202807171c5SHong Zhang r1 = 0.0; 203807171c5SHong Zhang n = r->i[i+1] - r->i[i]; 204807171c5SHong Zhang rj = r->j + r->i[i]; 205807171c5SHong Zhang ra = r->a + r->i[i]; 206807171c5SHong Zhang for (j=0; j<n; j++) { 207807171c5SHong Zhang r1 += (*ra++)*c[*rj++]; 208807171c5SHong Zhang } 209807171c5SHong Zhang d[col*rm + i] = r1; 210807171c5SHong Zhang } 211807171c5SHong Zhang } 2129566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(cn*2.0*(a->nz + r->nz))); 213807171c5SHong Zhang 2149566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayRead(B,&b)); 2159566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArray(RAB,&d)); 2169566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(RAB,MAT_FINAL_ASSEMBLY)); 2179566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(RAB,MAT_FINAL_ASSEMBLY)); 218807171c5SHong Zhang PetscFunctionReturn(0); 219807171c5SHong Zhang } 220807171c5SHong Zhang 22155bea0ebSHong Zhang PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_colorrart(Mat A,Mat R,Mat C) 222807171c5SHong Zhang { 2236718818eSStefano Zampini Mat_RARt *rart; 224807171c5SHong Zhang MatTransposeColoring matcoloring; 2258d0a38e4SHong Zhang Mat Rt,RARt; 226807171c5SHong Zhang 227807171c5SHong Zhang PetscFunctionBegin; 2286718818eSStefano Zampini MatCheckProduct(C,3); 22928b400f6SJacob Faibussowitsch PetscCheck(C->product->data,PetscObjectComm((PetscObject)C),PETSC_ERR_PLIB,"Product data empty"); 2306718818eSStefano Zampini rart = (Mat_RARt*)C->product->data; 2316718818eSStefano Zampini 232807171c5SHong Zhang /* Get dense Rt by Apply MatTransposeColoring to R */ 233807171c5SHong Zhang matcoloring = rart->matcoloring; 234807171c5SHong Zhang Rt = rart->Rt; 2359566063dSJacob Faibussowitsch PetscCall(MatTransColoringApplySpToDen(matcoloring,R,Rt)); 236807171c5SHong Zhang 23750647e95SHong Zhang /* Get dense RARt = R*A*Rt -- dominates! */ 238807171c5SHong Zhang RARt = rart->RARt; 2399566063dSJacob Faibussowitsch PetscCall(MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense(R,A,Rt,RARt,rart->work)); 240807171c5SHong Zhang 241807171c5SHong Zhang /* Recover C from C_dense */ 2429566063dSJacob Faibussowitsch PetscCall(MatTransColoringApplyDenToSp(matcoloring,RARt,C)); 243807171c5SHong Zhang PetscFunctionReturn(0); 244807171c5SHong Zhang } 245807171c5SHong Zhang 2464222ddf1SHong Zhang PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_matmattransposemult(Mat A,Mat R,PetscReal fill,Mat C) 247807171c5SHong Zhang { 2484222ddf1SHong Zhang Mat ARt; 2494fa85fdeSHong Zhang Mat_RARt *rart; 2506718818eSStefano Zampini char *alg; 2514fa85fdeSHong Zhang 25295a72cc5SHong Zhang PetscFunctionBegin; 2536718818eSStefano Zampini MatCheckProduct(C,4); 25428b400f6SJacob Faibussowitsch PetscCheck(!C->product->data,PetscObjectComm((PetscObject)C),PETSC_ERR_PLIB,"Product data not empty"); 2554222ddf1SHong Zhang /* create symbolic ARt = A*R^T */ 2569566063dSJacob Faibussowitsch PetscCall(MatProductCreate(A,R,NULL,&ARt)); 2579566063dSJacob Faibussowitsch PetscCall(MatProductSetType(ARt,MATPRODUCT_ABt)); 2589566063dSJacob Faibussowitsch PetscCall(MatProductSetAlgorithm(ARt,"sorted")); 2599566063dSJacob Faibussowitsch PetscCall(MatProductSetFill(ARt,fill)); 2609566063dSJacob Faibussowitsch PetscCall(MatProductSetFromOptions(ARt)); 2619566063dSJacob Faibussowitsch PetscCall(MatProductSymbolic(ARt)); 2624222ddf1SHong Zhang 2634222ddf1SHong Zhang /* compute symbolic C = R*ARt */ 2647a3c3d58SStefano Zampini /* set algorithm for C = R*ARt */ 2659566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(C->product->alg,&alg)); 2669566063dSJacob Faibussowitsch PetscCall(MatProductSetAlgorithm(C,"sorted")); 2679566063dSJacob Faibussowitsch PetscCall(MatMatMultSymbolic_SeqAIJ_SeqAIJ(R,ARt,fill,C)); 2687a3c3d58SStefano Zampini /* resume original algorithm for C */ 2699566063dSJacob Faibussowitsch PetscCall(MatProductSetAlgorithm(C,alg)); 2709566063dSJacob Faibussowitsch PetscCall(PetscFree(alg)); 2714222ddf1SHong Zhang 2724222ddf1SHong Zhang C->ops->rartnumeric = MatRARtNumeric_SeqAIJ_SeqAIJ_matmattransposemult; 27355bea0ebSHong Zhang 2749566063dSJacob Faibussowitsch PetscCall(PetscNew(&rart)); 2753b1b9624SHong Zhang rart->ARt = ARt; 2766718818eSStefano Zampini C->product->data = rart; 2776718818eSStefano Zampini C->product->destroy = MatDestroy_SeqAIJ_RARt; 2789566063dSJacob Faibussowitsch PetscCall(PetscInfo(C,"Use ARt=A*R^T, C=R*ARt via MatMatTransposeMult(). Coloring can be applied to A*R^T.\n")); 27955bea0ebSHong Zhang PetscFunctionReturn(0); 28036104f73SHong Zhang } 28155bea0ebSHong Zhang 28255bea0ebSHong Zhang PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_matmattransposemult(Mat A,Mat R,Mat C) 28355bea0ebSHong Zhang { 2846718818eSStefano Zampini Mat_RARt *rart; 28555bea0ebSHong Zhang 28655bea0ebSHong Zhang PetscFunctionBegin; 2876718818eSStefano Zampini MatCheckProduct(C,3); 28828b400f6SJacob Faibussowitsch PetscCheck(C->product->data,PetscObjectComm((PetscObject)C),PETSC_ERR_PLIB,"Product data empty"); 2896718818eSStefano Zampini rart = (Mat_RARt*)C->product->data; 2909566063dSJacob Faibussowitsch PetscCall(MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(A,R,rart->ARt)); /* dominate! */ 2919566063dSJacob Faibussowitsch PetscCall(MatMatMultNumeric_SeqAIJ_SeqAIJ(R,rart->ARt,C)); 29255bea0ebSHong Zhang PetscFunctionReturn(0); 293807171c5SHong Zhang } 29455bea0ebSHong Zhang 2954222ddf1SHong Zhang PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat R,PetscReal fill,Mat C) 29655bea0ebSHong Zhang { 29795a72cc5SHong Zhang Mat Rt; 29855bea0ebSHong Zhang Mat_RARt *rart; 29955bea0ebSHong Zhang 30055bea0ebSHong Zhang PetscFunctionBegin; 3016718818eSStefano Zampini MatCheckProduct(C,4); 30228b400f6SJacob Faibussowitsch PetscCheck(!C->product->data,PetscObjectComm((PetscObject)C),PETSC_ERR_PLIB,"Product data not empty"); 303*acd337a6SBarry Smith PetscCall(MatTranspose(R,MAT_INITIAL_MATRIX,&Rt)); 3049566063dSJacob Faibussowitsch PetscCall(MatMatMatMultSymbolic_SeqAIJ_SeqAIJ_SeqAIJ(R,A,Rt,fill,C)); 30595a72cc5SHong Zhang 3069566063dSJacob Faibussowitsch PetscCall(PetscNew(&rart)); 3076718818eSStefano Zampini rart->data = C->product->data; 3086718818eSStefano Zampini rart->destroy = C->product->destroy; 30995a72cc5SHong Zhang rart->Rt = Rt; 3106718818eSStefano Zampini C->product->data = rart; 3116718818eSStefano Zampini C->product->destroy = MatDestroy_SeqAIJ_RARt; 3124222ddf1SHong Zhang C->ops->rartnumeric = MatRARtNumeric_SeqAIJ_SeqAIJ; 3139566063dSJacob Faibussowitsch PetscCall(PetscInfo(C,"Use Rt=R^T and C=R*A*Rt via MatMatMatMult() to avoid sparse inner products\n")); 31455bea0ebSHong Zhang PetscFunctionReturn(0); 31555bea0ebSHong Zhang } 31655bea0ebSHong Zhang 31755bea0ebSHong Zhang PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ(Mat A,Mat R,Mat C) 31855bea0ebSHong Zhang { 3196718818eSStefano Zampini Mat_RARt *rart; 32055bea0ebSHong Zhang 32155bea0ebSHong Zhang PetscFunctionBegin; 3226718818eSStefano Zampini MatCheckProduct(C,3); 32328b400f6SJacob Faibussowitsch PetscCheck(C->product->data,PetscObjectComm((PetscObject)C),PETSC_ERR_PLIB,"Product data empty"); 3246718818eSStefano Zampini rart = (Mat_RARt*)C->product->data; 325*acd337a6SBarry Smith PetscCall(MatTranspose(R,MAT_REUSE_MATRIX,&rart->Rt)); 3266718818eSStefano Zampini /* MatMatMatMultSymbolic used a different data */ 3276718818eSStefano Zampini C->product->data = rart->data; 3289566063dSJacob Faibussowitsch PetscCall(MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqAIJ(R,A,rart->Rt,C)); 3296718818eSStefano Zampini C->product->data = rart; 33055bea0ebSHong Zhang PetscFunctionReturn(0); 33155bea0ebSHong Zhang } 33255bea0ebSHong Zhang 33355bea0ebSHong Zhang PetscErrorCode MatRARt_SeqAIJ_SeqAIJ(Mat A,Mat R,MatReuse scall,PetscReal fill,Mat *C) 33455bea0ebSHong Zhang { 33555bea0ebSHong Zhang const char *algTypes[3] = {"matmatmatmult","matmattransposemult","coloring_rart"}; 33655bea0ebSHong Zhang PetscInt alg=0; /* set default algorithm */ 33755bea0ebSHong Zhang 33855bea0ebSHong Zhang PetscFunctionBegin; 33955bea0ebSHong Zhang if (scall == MAT_INITIAL_MATRIX) { 340d0609cedSBarry Smith PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"MatRARt","Mat"); 3419566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-matrart_via","Algorithmic approach","MatRARt",algTypes,3,algTypes[0],&alg,NULL)); 342d0609cedSBarry Smith PetscOptionsEnd(); 343b56132cbSHong Zhang 3449566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(MAT_RARtSymbolic,A,R,0,0)); 3459566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF,C)); 34655bea0ebSHong Zhang switch (alg) { 34755bea0ebSHong Zhang case 1: 34855bea0ebSHong Zhang /* via matmattransposemult: ARt=A*R^T, C=R*ARt - matrix coloring can be applied to A*R^T */ 3499566063dSJacob Faibussowitsch PetscCall(MatRARtSymbolic_SeqAIJ_SeqAIJ_matmattransposemult(A,R,fill,*C)); 35055bea0ebSHong Zhang break; 35155bea0ebSHong Zhang case 2: 35255bea0ebSHong Zhang /* via coloring_rart: apply coloring C = R*A*R^T */ 3539566063dSJacob Faibussowitsch PetscCall(MatRARtSymbolic_SeqAIJ_SeqAIJ_colorrart(A,R,fill,*C)); 35455bea0ebSHong Zhang break; 35555bea0ebSHong Zhang default: 35655bea0ebSHong Zhang /* via matmatmatmult: Rt=R^T, C=R*A*Rt - avoid inefficient sparse inner products */ 3579566063dSJacob Faibussowitsch PetscCall(MatRARtSymbolic_SeqAIJ_SeqAIJ(A,R,fill,*C)); 35855bea0ebSHong Zhang break; 35955bea0ebSHong Zhang } 3609566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(MAT_RARtSymbolic,A,R,0,0)); 3615008f5a7SHong Zhang } 3625008f5a7SHong Zhang 3639566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(MAT_RARtNumeric,A,R,0,0)); 3649566063dSJacob Faibussowitsch PetscCall(((*C)->ops->rartnumeric)(A,R,*C)); 3659566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(MAT_RARtNumeric,A,R,0,0)); 366807171c5SHong Zhang PetscFunctionReturn(0); 367807171c5SHong Zhang } 3684222ddf1SHong Zhang 3694222ddf1SHong Zhang /* ------------------------------------------------------------- */ 3704222ddf1SHong Zhang PetscErrorCode MatProductSymbolic_RARt_SeqAIJ_SeqAIJ(Mat C) 3714222ddf1SHong Zhang { 3724222ddf1SHong Zhang Mat_Product *product = C->product; 3734222ddf1SHong Zhang Mat A=product->A,R=product->B; 3744222ddf1SHong Zhang MatProductAlgorithm alg=product->alg; 3754222ddf1SHong Zhang PetscReal fill=product->fill; 3764222ddf1SHong Zhang PetscBool flg; 3774222ddf1SHong Zhang 3784222ddf1SHong Zhang PetscFunctionBegin; 3799566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(alg,"r*a*rt",&flg)); 3804222ddf1SHong Zhang if (flg) { 3819566063dSJacob Faibussowitsch PetscCall(MatRARtSymbolic_SeqAIJ_SeqAIJ(A,R,fill,C)); 3824222ddf1SHong Zhang goto next; 3834222ddf1SHong Zhang } 3844222ddf1SHong Zhang 3859566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(alg,"r*art",&flg)); 3864222ddf1SHong Zhang if (flg) { 3879566063dSJacob Faibussowitsch PetscCall(MatRARtSymbolic_SeqAIJ_SeqAIJ_matmattransposemult(A,R,fill,C)); 3884222ddf1SHong Zhang goto next; 3894222ddf1SHong Zhang } 3904222ddf1SHong Zhang 3919566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(alg,"coloring_rart",&flg)); 3924222ddf1SHong Zhang if (flg) { 3939566063dSJacob Faibussowitsch PetscCall(MatRARtSymbolic_SeqAIJ_SeqAIJ_colorrart(A,R,fill,C)); 3944222ddf1SHong Zhang goto next; 3954222ddf1SHong Zhang } 3964222ddf1SHong Zhang 3974222ddf1SHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatProductAlgorithm is not supported"); 3984222ddf1SHong Zhang 3994222ddf1SHong Zhang next: 4004222ddf1SHong Zhang C->ops->productnumeric = MatProductNumeric_RARt; 4014222ddf1SHong Zhang PetscFunctionReturn(0); 4024222ddf1SHong Zhang } 403