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)); 216718818eSStefano Zampini if (rart->destroy) { 229566063dSJacob Faibussowitsch PetscCall((*rart->destroy)(rart->data)); 23807171c5SHong Zhang } 249566063dSJacob Faibussowitsch PetscCall(PetscFree(rart)); 25807171c5SHong Zhang PetscFunctionReturn(0); 26807171c5SHong Zhang } 27807171c5SHong Zhang 284222ddf1SHong Zhang PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_colorrart(Mat A,Mat R,PetscReal fill,Mat C) 29807171c5SHong Zhang { 30807171c5SHong Zhang Mat P; 31807171c5SHong Zhang PetscInt *rti,*rtj; 32807171c5SHong Zhang Mat_RARt *rart; 33335efc43SPeter Brune MatColoring coloring; 34807171c5SHong Zhang MatTransposeColoring matcoloring; 35807171c5SHong Zhang ISColoring iscoloring; 368d0a38e4SHong Zhang Mat Rt_dense,RARt_dense; 37807171c5SHong Zhang 38807171c5SHong Zhang PetscFunctionBegin; 396718818eSStefano Zampini MatCheckProduct(C,4); 4028b400f6SJacob Faibussowitsch PetscCheck(!C->product->data,PetscObjectComm((PetscObject)C),PETSC_ERR_PLIB,"Product data not empty"); 41807171c5SHong Zhang /* create symbolic P=Rt */ 429566063dSJacob Faibussowitsch PetscCall(MatGetSymbolicTranspose_SeqAIJ(R,&rti,&rtj)); 439566063dSJacob Faibussowitsch PetscCall(MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,R->cmap->n,R->rmap->n,rti,rtj,NULL,&P)); 44807171c5SHong Zhang 45807171c5SHong Zhang /* get symbolic C=Pt*A*P */ 469566063dSJacob Faibussowitsch PetscCall(MatPtAPSymbolic_SeqAIJ_SeqAIJ_SparseAxpy(A,P,fill,C)); 479566063dSJacob Faibussowitsch PetscCall(MatSetBlockSizes(C,PetscAbs(R->rmap->bs),PetscAbs(R->rmap->bs))); 484222ddf1SHong Zhang C->ops->rartnumeric = MatRARtNumeric_SeqAIJ_SeqAIJ_colorrart; 49807171c5SHong Zhang 50807171c5SHong Zhang /* create a supporting struct */ 519566063dSJacob Faibussowitsch PetscCall(PetscNew(&rart)); 526718818eSStefano Zampini C->product->data = rart; 536718818eSStefano Zampini C->product->destroy = MatDestroy_SeqAIJ_RARt; 54807171c5SHong Zhang 5522e94b5dSHong Zhang /* ------ Use coloring ---------- */ 564222ddf1SHong Zhang /* inode causes memory problem */ 579566063dSJacob Faibussowitsch PetscCall(MatSetOption(C,MAT_USE_INODES,PETSC_FALSE)); 584d478ae7SHong Zhang 59807171c5SHong Zhang /* Create MatTransposeColoring from symbolic C=R*A*R^T */ 609566063dSJacob Faibussowitsch PetscCall(MatColoringCreate(C,&coloring)); 619566063dSJacob Faibussowitsch PetscCall(MatColoringSetDistance(coloring,2)); 629566063dSJacob Faibussowitsch PetscCall(MatColoringSetType(coloring,MATCOLORINGSL)); 639566063dSJacob Faibussowitsch PetscCall(MatColoringSetFromOptions(coloring)); 649566063dSJacob Faibussowitsch PetscCall(MatColoringApply(coloring,&iscoloring)); 659566063dSJacob Faibussowitsch PetscCall(MatColoringDestroy(&coloring)); 669566063dSJacob Faibussowitsch PetscCall(MatTransposeColoringCreate(C,iscoloring,&matcoloring)); 672205254eSKarl Rupp 68807171c5SHong Zhang rart->matcoloring = matcoloring; 699566063dSJacob Faibussowitsch PetscCall(ISColoringDestroy(&iscoloring)); 703b1b9624SHong Zhang 71807171c5SHong Zhang /* Create Rt_dense */ 729566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF,&Rt_dense)); 739566063dSJacob Faibussowitsch PetscCall(MatSetSizes(Rt_dense,A->cmap->n,matcoloring->ncolors,A->cmap->n,matcoloring->ncolors)); 749566063dSJacob Faibussowitsch PetscCall(MatSetType(Rt_dense,MATSEQDENSE)); 759566063dSJacob Faibussowitsch PetscCall(MatSeqDenseSetPreallocation(Rt_dense,NULL)); 762205254eSKarl Rupp 77807171c5SHong Zhang Rt_dense->assembled = PETSC_TRUE; 78807171c5SHong Zhang rart->Rt = Rt_dense; 79807171c5SHong Zhang 80807171c5SHong Zhang /* Create RARt_dense = R*A*Rt_dense */ 819566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF,&RARt_dense)); 829566063dSJacob Faibussowitsch PetscCall(MatSetSizes(RARt_dense,C->rmap->n,matcoloring->ncolors,C->rmap->n,matcoloring->ncolors)); 839566063dSJacob Faibussowitsch PetscCall(MatSetType(RARt_dense,MATSEQDENSE)); 849566063dSJacob Faibussowitsch PetscCall(MatSeqDenseSetPreallocation(RARt_dense,NULL)); 852205254eSKarl Rupp 86807171c5SHong Zhang rart->RARt = RARt_dense; 87807171c5SHong Zhang 88807171c5SHong Zhang /* Allocate work array to store columns of A*R^T used in MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense() */ 899566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(A->rmap->n*4,&rart->work)); 90807171c5SHong Zhang 91807171c5SHong Zhang /* clean up */ 929566063dSJacob Faibussowitsch PetscCall(MatRestoreSymbolicTranspose_SeqAIJ(R,&rti,&rtj)); 939566063dSJacob Faibussowitsch PetscCall(MatDestroy(&P)); 94807171c5SHong Zhang 95807171c5SHong Zhang #if defined(PETSC_USE_INFO) 961ce71dffSSatish Balay { 976718818eSStefano Zampini Mat_SeqAIJ *c = (Mat_SeqAIJ*)C->data; 98807171c5SHong Zhang PetscReal density = (PetscReal)(c->nz)/(RARt_dense->rmap->n*RARt_dense->cmap->n); 999566063dSJacob Faibussowitsch PetscCall(PetscInfo(C,"C=R*(A*Rt) via coloring C - use sparse-dense inner products\n")); 1009566063dSJacob 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,density)); 1011ce71dffSSatish Balay } 102807171c5SHong Zhang #endif 103807171c5SHong Zhang PetscFunctionReturn(0); 104807171c5SHong Zhang } 105807171c5SHong Zhang 106807171c5SHong Zhang /* 107807171c5SHong Zhang RAB = R * A * B, R and A in seqaij format, B in dense format; 108807171c5SHong Zhang */ 109807171c5SHong Zhang PetscErrorCode MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense(Mat R,Mat A,Mat B,Mat RAB,PetscScalar *work) 110807171c5SHong Zhang { 111807171c5SHong Zhang Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*r=(Mat_SeqAIJ*)R->data; 1121683a169SBarry Smith PetscScalar r1,r2,r3,r4; 1131683a169SBarry Smith const PetscScalar *b,*b1,*b2,*b3,*b4; 114807171c5SHong Zhang MatScalar *aa,*ra; 115807171c5SHong Zhang PetscInt cn =B->cmap->n,bm=B->rmap->n,col,i,j,n,*ai=a->i,*aj,am=A->rmap->n; 116807171c5SHong Zhang PetscInt am2=2*am,am3=3*am,bm4=4*bm; 117807171c5SHong Zhang PetscScalar *d,*c,*c2,*c3,*c4; 118807171c5SHong Zhang PetscInt *rj,rm=R->rmap->n,dm=RAB->rmap->n,dn=RAB->cmap->n; 119807171c5SHong Zhang PetscInt rm2=2*rm,rm3=3*rm,colrm; 120807171c5SHong Zhang 121807171c5SHong Zhang PetscFunctionBegin; 122807171c5SHong Zhang if (!dm || !dn) PetscFunctionReturn(0); 123*08401ef6SPierre 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); 124*08401ef6SPierre 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); 125*08401ef6SPierre 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); 126*08401ef6SPierre 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); 127807171c5SHong Zhang 128274010c0SHong Zhang { /* 129274010c0SHong Zhang This approach is not as good as original ones (will be removed later), but it reveals that 130c4762a1bSJed Brown AB_den=A*B takes almost all execution time in R*A*B for src/ksp/ksp/tutorials/ex56.c 131274010c0SHong Zhang */ 132274010c0SHong Zhang PetscBool via_matmatmult=PETSC_FALSE; 1339566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-matrart_via_matmatmult",&via_matmatmult,NULL)); 134274010c0SHong Zhang if (via_matmatmult) { 1354222ddf1SHong Zhang Mat AB_den = NULL; 1369566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A),&AB_den)); 1379566063dSJacob Faibussowitsch PetscCall(MatMatMultSymbolic_SeqAIJ_SeqDense(A,B,0.0,AB_den)); 1389566063dSJacob Faibussowitsch PetscCall(MatMatMultNumeric_SeqAIJ_SeqDense(A,B,AB_den)); 1399566063dSJacob Faibussowitsch PetscCall(MatMatMultNumeric_SeqAIJ_SeqDense(R,AB_den,RAB)); 1409566063dSJacob Faibussowitsch PetscCall(MatDestroy(&AB_den)); 141274010c0SHong Zhang PetscFunctionReturn(0); 142274010c0SHong Zhang } 143274010c0SHong Zhang } 144274010c0SHong Zhang 1459566063dSJacob Faibussowitsch PetscCall(MatDenseGetArrayRead(B,&b)); 1469566063dSJacob Faibussowitsch PetscCall(MatDenseGetArray(RAB,&d)); 147807171c5SHong Zhang b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm; 148807171c5SHong Zhang c = work; c2 = c + am; c3 = c2 + am; c4 = c3 + am; 149807171c5SHong Zhang for (col=0; col<cn-4; col += 4) { /* over columns of C */ 150807171c5SHong Zhang for (i=0; i<am; i++) { /* over rows of A in those columns */ 151807171c5SHong Zhang r1 = r2 = r3 = r4 = 0.0; 152807171c5SHong Zhang n = ai[i+1] - ai[i]; 153807171c5SHong Zhang aj = a->j + ai[i]; 154807171c5SHong Zhang aa = a->a + ai[i]; 155807171c5SHong Zhang for (j=0; j<n; j++) { 156807171c5SHong Zhang r1 += (*aa)*b1[*aj]; 157807171c5SHong Zhang r2 += (*aa)*b2[*aj]; 158807171c5SHong Zhang r3 += (*aa)*b3[*aj]; 159807171c5SHong Zhang r4 += (*aa++)*b4[*aj++]; 160807171c5SHong Zhang } 161807171c5SHong Zhang c[i] = r1; 162807171c5SHong Zhang c[am + i] = r2; 163807171c5SHong Zhang c[am2 + i] = r3; 164807171c5SHong Zhang c[am3 + i] = r4; 165807171c5SHong Zhang } 166807171c5SHong Zhang b1 += bm4; 167807171c5SHong Zhang b2 += bm4; 168807171c5SHong Zhang b3 += bm4; 169807171c5SHong Zhang b4 += bm4; 170807171c5SHong Zhang 171807171c5SHong Zhang /* RAB[:,col] = R*C[:,col] */ 172807171c5SHong Zhang colrm = col*rm; 173807171c5SHong Zhang for (i=0; i<rm; i++) { /* over rows of R in those columns */ 174807171c5SHong Zhang r1 = r2 = r3 = r4 = 0.0; 175807171c5SHong Zhang n = r->i[i+1] - r->i[i]; 176807171c5SHong Zhang rj = r->j + r->i[i]; 177807171c5SHong Zhang ra = r->a + r->i[i]; 178807171c5SHong Zhang for (j=0; j<n; j++) { 179807171c5SHong Zhang r1 += (*ra)*c[*rj]; 180807171c5SHong Zhang r2 += (*ra)*c2[*rj]; 181807171c5SHong Zhang r3 += (*ra)*c3[*rj]; 182807171c5SHong Zhang r4 += (*ra++)*c4[*rj++]; 183807171c5SHong Zhang } 184807171c5SHong Zhang d[colrm + i] = r1; 185807171c5SHong Zhang d[colrm + rm + i] = r2; 186807171c5SHong Zhang d[colrm + rm2 + i] = r3; 187807171c5SHong Zhang d[colrm + rm3 + i] = r4; 188807171c5SHong Zhang } 189807171c5SHong Zhang } 190807171c5SHong Zhang for (; col<cn; col++) { /* over extra columns of C */ 191807171c5SHong Zhang for (i=0; i<am; i++) { /* over rows of A in those columns */ 192807171c5SHong Zhang r1 = 0.0; 193807171c5SHong Zhang n = a->i[i+1] - a->i[i]; 194807171c5SHong Zhang aj = a->j + a->i[i]; 195807171c5SHong Zhang aa = a->a + a->i[i]; 196807171c5SHong Zhang for (j=0; j<n; j++) { 197807171c5SHong Zhang r1 += (*aa++)*b1[*aj++]; 198807171c5SHong Zhang } 199807171c5SHong Zhang c[i] = r1; 200807171c5SHong Zhang } 201807171c5SHong Zhang b1 += bm; 202807171c5SHong Zhang 203807171c5SHong Zhang for (i=0; i<rm; i++) { /* over rows of R in those columns */ 204807171c5SHong Zhang r1 = 0.0; 205807171c5SHong Zhang n = r->i[i+1] - r->i[i]; 206807171c5SHong Zhang rj = r->j + r->i[i]; 207807171c5SHong Zhang ra = r->a + r->i[i]; 208807171c5SHong Zhang for (j=0; j<n; j++) { 209807171c5SHong Zhang r1 += (*ra++)*c[*rj++]; 210807171c5SHong Zhang } 211807171c5SHong Zhang d[col*rm + i] = r1; 212807171c5SHong Zhang } 213807171c5SHong Zhang } 2149566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(cn*2.0*(a->nz + r->nz))); 215807171c5SHong Zhang 2169566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArrayRead(B,&b)); 2179566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArray(RAB,&d)); 2189566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(RAB,MAT_FINAL_ASSEMBLY)); 2199566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(RAB,MAT_FINAL_ASSEMBLY)); 220807171c5SHong Zhang PetscFunctionReturn(0); 221807171c5SHong Zhang } 222807171c5SHong Zhang 22355bea0ebSHong Zhang PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_colorrart(Mat A,Mat R,Mat C) 224807171c5SHong Zhang { 2256718818eSStefano Zampini Mat_RARt *rart; 226807171c5SHong Zhang MatTransposeColoring matcoloring; 2278d0a38e4SHong Zhang Mat Rt,RARt; 228807171c5SHong Zhang 229807171c5SHong Zhang PetscFunctionBegin; 2306718818eSStefano Zampini MatCheckProduct(C,3); 23128b400f6SJacob Faibussowitsch PetscCheck(C->product->data,PetscObjectComm((PetscObject)C),PETSC_ERR_PLIB,"Product data empty"); 2326718818eSStefano Zampini rart = (Mat_RARt*)C->product->data; 2336718818eSStefano Zampini 234807171c5SHong Zhang /* Get dense Rt by Apply MatTransposeColoring to R */ 235807171c5SHong Zhang matcoloring = rart->matcoloring; 236807171c5SHong Zhang Rt = rart->Rt; 2379566063dSJacob Faibussowitsch PetscCall(MatTransColoringApplySpToDen(matcoloring,R,Rt)); 238807171c5SHong Zhang 23950647e95SHong Zhang /* Get dense RARt = R*A*Rt -- dominates! */ 240807171c5SHong Zhang RARt = rart->RARt; 2419566063dSJacob Faibussowitsch PetscCall(MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense(R,A,Rt,RARt,rart->work)); 242807171c5SHong Zhang 243807171c5SHong Zhang /* Recover C from C_dense */ 2449566063dSJacob Faibussowitsch PetscCall(MatTransColoringApplyDenToSp(matcoloring,RARt,C)); 245807171c5SHong Zhang PetscFunctionReturn(0); 246807171c5SHong Zhang } 247807171c5SHong Zhang 2484222ddf1SHong Zhang PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_matmattransposemult(Mat A,Mat R,PetscReal fill,Mat C) 249807171c5SHong Zhang { 2504222ddf1SHong Zhang Mat ARt; 2514fa85fdeSHong Zhang Mat_RARt *rart; 2526718818eSStefano Zampini char *alg; 2534fa85fdeSHong Zhang 25495a72cc5SHong Zhang PetscFunctionBegin; 2556718818eSStefano Zampini MatCheckProduct(C,4); 25628b400f6SJacob Faibussowitsch PetscCheck(!C->product->data,PetscObjectComm((PetscObject)C),PETSC_ERR_PLIB,"Product data not empty"); 2574222ddf1SHong Zhang /* create symbolic ARt = A*R^T */ 2589566063dSJacob Faibussowitsch PetscCall(MatProductCreate(A,R,NULL,&ARt)); 2599566063dSJacob Faibussowitsch PetscCall(MatProductSetType(ARt,MATPRODUCT_ABt)); 2609566063dSJacob Faibussowitsch PetscCall(MatProductSetAlgorithm(ARt,"sorted")); 2619566063dSJacob Faibussowitsch PetscCall(MatProductSetFill(ARt,fill)); 2629566063dSJacob Faibussowitsch PetscCall(MatProductSetFromOptions(ARt)); 2639566063dSJacob Faibussowitsch PetscCall(MatProductSymbolic(ARt)); 2644222ddf1SHong Zhang 2654222ddf1SHong Zhang /* compute symbolic C = R*ARt */ 2667a3c3d58SStefano Zampini /* set algorithm for C = R*ARt */ 2679566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(C->product->alg,&alg)); 2689566063dSJacob Faibussowitsch PetscCall(MatProductSetAlgorithm(C,"sorted")); 2699566063dSJacob Faibussowitsch PetscCall(MatMatMultSymbolic_SeqAIJ_SeqAIJ(R,ARt,fill,C)); 2707a3c3d58SStefano Zampini /* resume original algorithm for C */ 2719566063dSJacob Faibussowitsch PetscCall(MatProductSetAlgorithm(C,alg)); 2729566063dSJacob Faibussowitsch PetscCall(PetscFree(alg)); 2734222ddf1SHong Zhang 2744222ddf1SHong Zhang C->ops->rartnumeric = MatRARtNumeric_SeqAIJ_SeqAIJ_matmattransposemult; 27555bea0ebSHong Zhang 2769566063dSJacob Faibussowitsch PetscCall(PetscNew(&rart)); 2773b1b9624SHong Zhang rart->ARt = ARt; 2786718818eSStefano Zampini C->product->data = rart; 2796718818eSStefano Zampini C->product->destroy = MatDestroy_SeqAIJ_RARt; 2809566063dSJacob Faibussowitsch PetscCall(PetscInfo(C,"Use ARt=A*R^T, C=R*ARt via MatMatTransposeMult(). Coloring can be applied to A*R^T.\n")); 28155bea0ebSHong Zhang PetscFunctionReturn(0); 28236104f73SHong Zhang } 28355bea0ebSHong Zhang 28455bea0ebSHong Zhang PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_matmattransposemult(Mat A,Mat R,Mat C) 28555bea0ebSHong Zhang { 2866718818eSStefano Zampini Mat_RARt *rart; 28755bea0ebSHong Zhang 28855bea0ebSHong Zhang PetscFunctionBegin; 2896718818eSStefano Zampini MatCheckProduct(C,3); 29028b400f6SJacob Faibussowitsch PetscCheck(C->product->data,PetscObjectComm((PetscObject)C),PETSC_ERR_PLIB,"Product data empty"); 2916718818eSStefano Zampini rart = (Mat_RARt*)C->product->data; 2929566063dSJacob Faibussowitsch PetscCall(MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(A,R,rart->ARt)); /* dominate! */ 2939566063dSJacob Faibussowitsch PetscCall(MatMatMultNumeric_SeqAIJ_SeqAIJ(R,rart->ARt,C)); 29455bea0ebSHong Zhang PetscFunctionReturn(0); 295807171c5SHong Zhang } 29655bea0ebSHong Zhang 2974222ddf1SHong Zhang PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat R,PetscReal fill,Mat C) 29855bea0ebSHong Zhang { 29995a72cc5SHong Zhang Mat Rt; 30055bea0ebSHong Zhang Mat_RARt *rart; 30155bea0ebSHong Zhang 30255bea0ebSHong Zhang PetscFunctionBegin; 3036718818eSStefano Zampini MatCheckProduct(C,4); 30428b400f6SJacob Faibussowitsch PetscCheck(!C->product->data,PetscObjectComm((PetscObject)C),PETSC_ERR_PLIB,"Product data not empty"); 3059566063dSJacob Faibussowitsch PetscCall(MatTranspose_SeqAIJ(R,MAT_INITIAL_MATRIX,&Rt)); 3069566063dSJacob Faibussowitsch PetscCall(MatMatMatMultSymbolic_SeqAIJ_SeqAIJ_SeqAIJ(R,A,Rt,fill,C)); 30795a72cc5SHong Zhang 3089566063dSJacob Faibussowitsch PetscCall(PetscNew(&rart)); 3096718818eSStefano Zampini rart->data = C->product->data; 3106718818eSStefano Zampini rart->destroy = C->product->destroy; 31195a72cc5SHong Zhang rart->Rt = Rt; 3126718818eSStefano Zampini C->product->data = rart; 3136718818eSStefano Zampini C->product->destroy = MatDestroy_SeqAIJ_RARt; 3144222ddf1SHong Zhang C->ops->rartnumeric = MatRARtNumeric_SeqAIJ_SeqAIJ; 3159566063dSJacob Faibussowitsch PetscCall(PetscInfo(C,"Use Rt=R^T and C=R*A*Rt via MatMatMatMult() to avoid sparse inner products\n")); 31655bea0ebSHong Zhang PetscFunctionReturn(0); 31755bea0ebSHong Zhang } 31855bea0ebSHong Zhang 31955bea0ebSHong Zhang PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ(Mat A,Mat R,Mat C) 32055bea0ebSHong Zhang { 3216718818eSStefano Zampini Mat_RARt *rart; 32255bea0ebSHong Zhang 32355bea0ebSHong Zhang PetscFunctionBegin; 3246718818eSStefano Zampini MatCheckProduct(C,3); 32528b400f6SJacob Faibussowitsch PetscCheck(C->product->data,PetscObjectComm((PetscObject)C),PETSC_ERR_PLIB,"Product data empty"); 3266718818eSStefano Zampini rart = (Mat_RARt*)C->product->data; 3279566063dSJacob Faibussowitsch PetscCall(MatTranspose_SeqAIJ(R,MAT_REUSE_MATRIX,&rart->Rt)); 3286718818eSStefano Zampini /* MatMatMatMultSymbolic used a different data */ 3296718818eSStefano Zampini C->product->data = rart->data; 3309566063dSJacob Faibussowitsch PetscCall(MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqAIJ(R,A,rart->Rt,C)); 3316718818eSStefano Zampini C->product->data = rart; 33255bea0ebSHong Zhang PetscFunctionReturn(0); 33355bea0ebSHong Zhang } 33455bea0ebSHong Zhang 33555bea0ebSHong Zhang PetscErrorCode MatRARt_SeqAIJ_SeqAIJ(Mat A,Mat R,MatReuse scall,PetscReal fill,Mat *C) 33655bea0ebSHong Zhang { 33755bea0ebSHong Zhang PetscErrorCode ierr; 33855bea0ebSHong Zhang const char *algTypes[3] = {"matmatmatmult","matmattransposemult","coloring_rart"}; 33955bea0ebSHong Zhang PetscInt alg=0; /* set default algorithm */ 34055bea0ebSHong Zhang 34155bea0ebSHong Zhang PetscFunctionBegin; 34255bea0ebSHong Zhang if (scall == MAT_INITIAL_MATRIX) { 3439566063dSJacob Faibussowitsch ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"MatRARt","Mat");PetscCall(ierr); 3449566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-matrart_via","Algorithmic approach","MatRARt",algTypes,3,algTypes[0],&alg,NULL)); 3459566063dSJacob Faibussowitsch ierr = PetscOptionsEnd();PetscCall(ierr); 346b56132cbSHong Zhang 3479566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(MAT_RARtSymbolic,A,R,0,0)); 3489566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF,C)); 34955bea0ebSHong Zhang switch (alg) { 35055bea0ebSHong Zhang case 1: 35155bea0ebSHong Zhang /* via matmattransposemult: ARt=A*R^T, C=R*ARt - matrix coloring can be applied to A*R^T */ 3529566063dSJacob Faibussowitsch PetscCall(MatRARtSymbolic_SeqAIJ_SeqAIJ_matmattransposemult(A,R,fill,*C)); 35355bea0ebSHong Zhang break; 35455bea0ebSHong Zhang case 2: 35555bea0ebSHong Zhang /* via coloring_rart: apply coloring C = R*A*R^T */ 3569566063dSJacob Faibussowitsch PetscCall(MatRARtSymbolic_SeqAIJ_SeqAIJ_colorrart(A,R,fill,*C)); 35755bea0ebSHong Zhang break; 35855bea0ebSHong Zhang default: 35955bea0ebSHong Zhang /* via matmatmatmult: Rt=R^T, C=R*A*Rt - avoid inefficient sparse inner products */ 3609566063dSJacob Faibussowitsch PetscCall(MatRARtSymbolic_SeqAIJ_SeqAIJ(A,R,fill,*C)); 36155bea0ebSHong Zhang break; 36255bea0ebSHong Zhang } 3639566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(MAT_RARtSymbolic,A,R,0,0)); 3645008f5a7SHong Zhang } 3655008f5a7SHong Zhang 3669566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(MAT_RARtNumeric,A,R,0,0)); 3679566063dSJacob Faibussowitsch PetscCall(((*C)->ops->rartnumeric)(A,R,*C)); 3689566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(MAT_RARtNumeric,A,R,0,0)); 369807171c5SHong Zhang PetscFunctionReturn(0); 370807171c5SHong Zhang } 3714222ddf1SHong Zhang 3724222ddf1SHong Zhang /* ------------------------------------------------------------- */ 3734222ddf1SHong Zhang PetscErrorCode MatProductSymbolic_RARt_SeqAIJ_SeqAIJ(Mat C) 3744222ddf1SHong Zhang { 3754222ddf1SHong Zhang Mat_Product *product = C->product; 3764222ddf1SHong Zhang Mat A=product->A,R=product->B; 3774222ddf1SHong Zhang MatProductAlgorithm alg=product->alg; 3784222ddf1SHong Zhang PetscReal fill=product->fill; 3794222ddf1SHong Zhang PetscBool flg; 3804222ddf1SHong Zhang 3814222ddf1SHong Zhang PetscFunctionBegin; 3829566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(alg,"r*a*rt",&flg)); 3834222ddf1SHong Zhang if (flg) { 3849566063dSJacob Faibussowitsch PetscCall(MatRARtSymbolic_SeqAIJ_SeqAIJ(A,R,fill,C)); 3854222ddf1SHong Zhang goto next; 3864222ddf1SHong Zhang } 3874222ddf1SHong Zhang 3889566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(alg,"r*art",&flg)); 3894222ddf1SHong Zhang if (flg) { 3909566063dSJacob Faibussowitsch PetscCall(MatRARtSymbolic_SeqAIJ_SeqAIJ_matmattransposemult(A,R,fill,C)); 3914222ddf1SHong Zhang goto next; 3924222ddf1SHong Zhang } 3934222ddf1SHong Zhang 3949566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(alg,"coloring_rart",&flg)); 3954222ddf1SHong Zhang if (flg) { 3969566063dSJacob Faibussowitsch PetscCall(MatRARtSymbolic_SeqAIJ_SeqAIJ_colorrart(A,R,fill,C)); 3974222ddf1SHong Zhang goto next; 3984222ddf1SHong Zhang } 3994222ddf1SHong Zhang 4004222ddf1SHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatProductAlgorithm is not supported"); 4014222ddf1SHong Zhang 4024222ddf1SHong Zhang next: 4034222ddf1SHong Zhang C->ops->productnumeric = MatProductNumeric_RARt; 4044222ddf1SHong Zhang PetscFunctionReturn(0); 4054222ddf1SHong Zhang } 406