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*/ 108563dfccSBarry Smith #include <petsctime.h> 11807171c5SHong Zhang 123b1b9624SHong Zhang /* #define RART_PROFILE */ 13807171c5SHong Zhang 14807171c5SHong Zhang #undef __FUNCT__ 15257c235dSHong Zhang #define __FUNCT__ "MatDestroy_SeqAIJ_RARt" 16257c235dSHong Zhang PetscErrorCode MatDestroy_SeqAIJ_RARt(Mat A) 17807171c5SHong Zhang { 18807171c5SHong Zhang PetscErrorCode ierr; 19257c235dSHong Zhang Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 20257c235dSHong Zhang Mat_RARt *rart = a->rart; 21807171c5SHong Zhang 22807171c5SHong Zhang PetscFunctionBegin; 23807171c5SHong Zhang ierr = MatTransposeColoringDestroy(&rart->matcoloring);CHKERRQ(ierr); 24807171c5SHong Zhang ierr = MatDestroy(&rart->Rt);CHKERRQ(ierr); 25807171c5SHong Zhang ierr = MatDestroy(&rart->RARt);CHKERRQ(ierr); 263b1b9624SHong Zhang ierr = MatDestroy(&rart->ARt);CHKERRQ(ierr); 27807171c5SHong Zhang ierr = PetscFree(rart->work);CHKERRQ(ierr); 28807171c5SHong Zhang 29807171c5SHong Zhang A->ops->destroy = rart->destroy; 30807171c5SHong Zhang if (A->ops->destroy) { 31807171c5SHong Zhang ierr = (*A->ops->destroy)(A);CHKERRQ(ierr); 32807171c5SHong Zhang } 33257c235dSHong Zhang ierr = PetscFree(rart);CHKERRQ(ierr); 34807171c5SHong Zhang PetscFunctionReturn(0); 35807171c5SHong Zhang } 36807171c5SHong Zhang 37807171c5SHong Zhang #undef __FUNCT__ 38*55bea0ebSHong Zhang #define __FUNCT__ "MatRARtSymbolic_SeqAIJ_SeqAIJ_colorrart" 39*55bea0ebSHong Zhang PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_colorrart(Mat A,Mat R,PetscReal fill,Mat *C) 40807171c5SHong Zhang { 41807171c5SHong Zhang PetscErrorCode ierr; 42807171c5SHong Zhang Mat P; 43807171c5SHong Zhang PetscInt *rti,*rtj; 44807171c5SHong Zhang Mat_RARt *rart; 45807171c5SHong Zhang MatTransposeColoring matcoloring; 46807171c5SHong Zhang ISColoring iscoloring; 478d0a38e4SHong Zhang Mat Rt_dense,RARt_dense; 483b1b9624SHong Zhang #if defined(RART_PROFILE) 49807171c5SHong Zhang PetscLogDouble GColor=0.0,MCCreate=0.0,MDenCreate=0.0,t0,tf,etime=0.0; 503b1b9624SHong Zhang #endif 51807171c5SHong Zhang Mat_SeqAIJ *c; 52807171c5SHong Zhang 53807171c5SHong Zhang PetscFunctionBegin; 543b1b9624SHong Zhang #if defined(RART_PROFILE) 558563dfccSBarry Smith ierr = PetscTime(&t0);CHKERRQ(ierr); 563b1b9624SHong Zhang #endif 57807171c5SHong Zhang /* create symbolic P=Rt */ 58807171c5SHong Zhang ierr = MatGetSymbolicTranspose_SeqAIJ(R,&rti,&rtj);CHKERRQ(ierr); 590298fd71SBarry Smith ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,R->cmap->n,R->rmap->n,rti,rtj,NULL,&P);CHKERRQ(ierr); 60807171c5SHong Zhang 61807171c5SHong Zhang /* get symbolic C=Pt*A*P */ 62807171c5SHong Zhang ierr = MatPtAPSymbolic_SeqAIJ_SeqAIJ(A,P,fill,C);CHKERRQ(ierr); 63fd8e3dafSHong Zhang (*C)->rmap->bs = R->rmap->bs; 64fd8e3dafSHong Zhang (*C)->cmap->bs = R->rmap->bs; 65*55bea0ebSHong Zhang (*C)->ops->rartnumeric = MatRARtNumeric_SeqAIJ_SeqAIJ_colorrart; 66807171c5SHong Zhang 67807171c5SHong Zhang /* create a supporting struct */ 68807171c5SHong Zhang ierr = PetscNew(Mat_RARt,&rart);CHKERRQ(ierr); 69257c235dSHong Zhang c = (Mat_SeqAIJ*)(*C)->data; 70257c235dSHong Zhang c->rart = rart; 713b1b9624SHong Zhang #if defined(RART_PROFILE) 728563dfccSBarry Smith ierr = PetscTime(&tf);CHKERRQ(ierr); 73807171c5SHong Zhang etime += tf - t0; 743b1b9624SHong Zhang #endif 75807171c5SHong Zhang 7622e94b5dSHong Zhang /* ------ Use coloring ---------- */ 77807171c5SHong Zhang /* Create MatTransposeColoring from symbolic C=R*A*R^T */ 783b1b9624SHong Zhang #if defined(RART_PROFILE) 798563dfccSBarry Smith ierr = PetscTime(&t0);CHKERRQ(ierr); 803b1b9624SHong Zhang #endif 81807171c5SHong Zhang ierr = MatGetColoring(*C,MATCOLORINGLF,&iscoloring);CHKERRQ(ierr); 823b1b9624SHong Zhang #if defined(RART_PROFILE) 838563dfccSBarry Smith ierr = PetscTime(&tf);CHKERRQ(ierr); 84807171c5SHong Zhang GColor += tf - t0; 858563dfccSBarry Smith ierr = PetscTime(&t0);CHKERRQ(ierr); 863b1b9624SHong Zhang #endif 87807171c5SHong Zhang ierr = MatTransposeColoringCreate(*C,iscoloring,&matcoloring);CHKERRQ(ierr); 882205254eSKarl Rupp 89807171c5SHong Zhang rart->matcoloring = matcoloring; 902205254eSKarl Rupp 91807171c5SHong Zhang ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr); 923b1b9624SHong Zhang #if defined(RART_PROFILE) 938563dfccSBarry Smith ierr = PetscTime(&tf);CHKERRQ(ierr); 94807171c5SHong Zhang MCCreate += tf - t0; 958563dfccSBarry Smith ierr = PetscTime(&t0);CHKERRQ(ierr); 963b1b9624SHong Zhang #endif 973b1b9624SHong Zhang 98807171c5SHong Zhang /* Create Rt_dense */ 99807171c5SHong Zhang ierr = MatCreate(PETSC_COMM_SELF,&Rt_dense);CHKERRQ(ierr); 100807171c5SHong Zhang ierr = MatSetSizes(Rt_dense,A->cmap->n,matcoloring->ncolors,A->cmap->n,matcoloring->ncolors);CHKERRQ(ierr); 101807171c5SHong Zhang ierr = MatSetType(Rt_dense,MATSEQDENSE);CHKERRQ(ierr); 1020298fd71SBarry Smith ierr = MatSeqDenseSetPreallocation(Rt_dense,NULL);CHKERRQ(ierr); 1032205254eSKarl Rupp 104807171c5SHong Zhang Rt_dense->assembled = PETSC_TRUE; 105807171c5SHong Zhang rart->Rt = Rt_dense; 106807171c5SHong Zhang 107807171c5SHong Zhang /* Create RARt_dense = R*A*Rt_dense */ 108807171c5SHong Zhang ierr = MatCreate(PETSC_COMM_SELF,&RARt_dense);CHKERRQ(ierr); 109807171c5SHong Zhang ierr = MatSetSizes(RARt_dense,(*C)->rmap->n,matcoloring->ncolors,(*C)->rmap->n,matcoloring->ncolors);CHKERRQ(ierr); 110807171c5SHong Zhang ierr = MatSetType(RARt_dense,MATSEQDENSE);CHKERRQ(ierr); 1110298fd71SBarry Smith ierr = MatSeqDenseSetPreallocation(RARt_dense,NULL);CHKERRQ(ierr); 1122205254eSKarl Rupp 113807171c5SHong Zhang rart->RARt = RARt_dense; 114807171c5SHong Zhang 115807171c5SHong Zhang /* Allocate work array to store columns of A*R^T used in MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense() */ 116807171c5SHong Zhang ierr = PetscMalloc(A->rmap->n*4*sizeof(PetscScalar),&rart->work);CHKERRQ(ierr); 117807171c5SHong Zhang 1183b1b9624SHong Zhang #if defined(RART_PROFILE) 1198563dfccSBarry Smith ierr = PetscTime(&tf);CHKERRQ(ierr); 120807171c5SHong Zhang MDenCreate += tf - t0; 1213b1b9624SHong Zhang #endif 122807171c5SHong Zhang 123807171c5SHong Zhang rart->destroy = (*C)->ops->destroy; 124807171c5SHong Zhang (*C)->ops->destroy = MatDestroy_SeqAIJ_RARt; 125807171c5SHong Zhang 126807171c5SHong Zhang /* clean up */ 127807171c5SHong Zhang ierr = MatRestoreSymbolicTranspose_SeqAIJ(R,&rti,&rtj);CHKERRQ(ierr); 128807171c5SHong Zhang ierr = MatDestroy(&P);CHKERRQ(ierr); 129807171c5SHong Zhang 130807171c5SHong Zhang #if defined(PETSC_USE_INFO) 1311ce71dffSSatish Balay { 132807171c5SHong Zhang PetscReal density= (PetscReal)(c->nz)/(RARt_dense->rmap->n*RARt_dense->cmap->n); 133a2ea699eSBarry Smith ierr = PetscInfo6(*C,"RARt_den %D %D; Rt_den %D %D, (RARt->nz %D)/(m*ncolors)=%g\n",RARt_dense->rmap->n,RARt_dense->cmap->n,Rt_dense->rmap->n,Rt_dense->cmap->n,c->nz,density);CHKERRQ(ierr); 1343b1b9624SHong Zhang #if defined(RART_PROFILE) 135a2ea699eSBarry Smith ierr = PetscInfo5(*C,"Sym = GetColor %g + MColorCreate %g + MDenCreate %g + other %g = %g\n",GColor,MCCreate,MDenCreate,etime,GColor+MCCreate+MDenCreate+etime);CHKERRQ(ierr); 1363b1b9624SHong Zhang #endif 1371ce71dffSSatish Balay } 138807171c5SHong Zhang #endif 139807171c5SHong Zhang PetscFunctionReturn(0); 140807171c5SHong Zhang } 141807171c5SHong Zhang 142807171c5SHong Zhang /* 143807171c5SHong Zhang RAB = R * A * B, R and A in seqaij format, B in dense format; 144807171c5SHong Zhang */ 145807171c5SHong Zhang #undef __FUNCT__ 146807171c5SHong Zhang #define __FUNCT__ "MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense" 147807171c5SHong Zhang PetscErrorCode MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense(Mat R,Mat A,Mat B,Mat RAB,PetscScalar *work) 148807171c5SHong Zhang { 149807171c5SHong Zhang Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*r=(Mat_SeqAIJ*)R->data; 150807171c5SHong Zhang PetscErrorCode ierr; 151807171c5SHong Zhang PetscScalar *b,r1,r2,r3,r4,*b1,*b2,*b3,*b4; 152807171c5SHong Zhang MatScalar *aa,*ra; 153807171c5SHong Zhang PetscInt cn =B->cmap->n,bm=B->rmap->n,col,i,j,n,*ai=a->i,*aj,am=A->rmap->n; 154807171c5SHong Zhang PetscInt am2=2*am,am3=3*am,bm4=4*bm; 155807171c5SHong Zhang PetscScalar *d,*c,*c2,*c3,*c4; 156807171c5SHong Zhang PetscInt *rj,rm=R->rmap->n,dm=RAB->rmap->n,dn=RAB->cmap->n; 157807171c5SHong Zhang PetscInt rm2=2*rm,rm3=3*rm,colrm; 1583b1b9624SHong Zhang #if defined(RART_PROFILE) 159274010c0SHong Zhang PetscLogDouble t0,tf,Mult_sp_den1=0.0,Mult_sp_den2=0.0; 1603b1b9624SHong Zhang #endif 161807171c5SHong Zhang 162807171c5SHong Zhang PetscFunctionBegin; 163807171c5SHong Zhang if (!dm || !dn) PetscFunctionReturn(0); 164807171c5SHong Zhang if (bm != A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Number columns in A %D not equal rows in B %D\n",A->cmap->n,bm); 165807171c5SHong Zhang if (am != R->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Number columns in R %D not equal rows in A %D\n",R->cmap->n,am); 166807171c5SHong Zhang if (R->rmap->n != RAB->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Number rows in RAB %D not equal rows in R %D\n",RAB->rmap->n,R->rmap->n); 167807171c5SHong Zhang if (B->cmap->n != RAB->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Number columns in RAB %D not equal columns in B %D\n",RAB->cmap->n,B->cmap->n); 168807171c5SHong Zhang 169274010c0SHong Zhang { /* 170274010c0SHong Zhang This approach is not as good as original ones (will be removed later), but it reveals that 171274010c0SHong Zhang AB_den=A*B takes almost all execution time in R*A*B for src/ksp/ksp/examples/tutorials/ex56.c 172274010c0SHong Zhang */ 173274010c0SHong Zhang PetscBool via_matmatmult=PETSC_FALSE; 174274010c0SHong Zhang ierr = PetscOptionsGetBool(NULL,"-matrart_via_matmatmult",&via_matmatmult,NULL);CHKERRQ(ierr); 175274010c0SHong Zhang if (via_matmatmult) { 176274010c0SHong Zhang Mat AB_den; 177274010c0SHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqDense(A,B,0.0,&AB_den);CHKERRQ(ierr); 178274010c0SHong Zhang 1793b1b9624SHong Zhang #if defined(RART_PROFILE) 18026be7272SHong Zhang ierr = PetscTime(&t0);CHKERRQ(ierr); 1813b1b9624SHong Zhang #endif 182274010c0SHong Zhang ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,B,AB_den);CHKERRQ(ierr); 1833b1b9624SHong Zhang #if defined(RART_PROFILE) 18426be7272SHong Zhang ierr = PetscTime(&tf);CHKERRQ(ierr); 185274010c0SHong Zhang Mult_sp_den1 += tf - t0; 18626be7272SHong Zhang ierr = PetscTime(&t0);CHKERRQ(ierr); 1873b1b9624SHong Zhang #endif 188274010c0SHong Zhang ierr = MatMatMultNumeric_SeqAIJ_SeqDense(R,AB_den,RAB);CHKERRQ(ierr); 1893b1b9624SHong Zhang #if defined(RART_PROFILE) 19026be7272SHong Zhang ierr = PetscTime(&tf);CHKERRQ(ierr); 191274010c0SHong Zhang Mult_sp_den2 += tf - t0; 192274010c0SHong Zhang printf(" MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense...time = %g + %g = %g; ncolor: %d\n",Mult_sp_den1, Mult_sp_den2, Mult_sp_den1+ Mult_sp_den2,B->cmap->n); 1933b1b9624SHong Zhang #endif 194274010c0SHong Zhang ierr = MatDestroy(&AB_den);CHKERRQ(ierr); 195274010c0SHong Zhang PetscFunctionReturn(0); 196274010c0SHong Zhang } 197274010c0SHong Zhang } 198274010c0SHong Zhang 1998c778c55SBarry Smith ierr = MatDenseGetArray(B,&b);CHKERRQ(ierr); 2008c778c55SBarry Smith ierr = MatDenseGetArray(RAB,&d);CHKERRQ(ierr); 201807171c5SHong Zhang b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm; 202807171c5SHong Zhang c = work; c2 = c + am; c3 = c2 + am; c4 = c3 + am; 203807171c5SHong Zhang for (col=0; col<cn-4; col += 4) { /* over columns of C */ 204807171c5SHong Zhang for (i=0; i<am; i++) { /* over rows of A in those columns */ 205807171c5SHong Zhang r1 = r2 = r3 = r4 = 0.0; 206807171c5SHong Zhang n = ai[i+1] - ai[i]; 207807171c5SHong Zhang aj = a->j + ai[i]; 208807171c5SHong Zhang aa = a->a + ai[i]; 209807171c5SHong Zhang for (j=0; j<n; j++) { 210807171c5SHong Zhang r1 += (*aa)*b1[*aj]; 211807171c5SHong Zhang r2 += (*aa)*b2[*aj]; 212807171c5SHong Zhang r3 += (*aa)*b3[*aj]; 213807171c5SHong Zhang r4 += (*aa++)*b4[*aj++]; 214807171c5SHong Zhang } 215807171c5SHong Zhang c[i] = r1; 216807171c5SHong Zhang c[am + i] = r2; 217807171c5SHong Zhang c[am2 + i] = r3; 218807171c5SHong Zhang c[am3 + i] = r4; 219807171c5SHong Zhang } 220807171c5SHong Zhang b1 += bm4; 221807171c5SHong Zhang b2 += bm4; 222807171c5SHong Zhang b3 += bm4; 223807171c5SHong Zhang b4 += bm4; 224807171c5SHong Zhang 225807171c5SHong Zhang /* RAB[:,col] = R*C[:,col] */ 226807171c5SHong Zhang colrm = col*rm; 227807171c5SHong Zhang for (i=0; i<rm; i++) { /* over rows of R in those columns */ 228807171c5SHong Zhang r1 = r2 = r3 = r4 = 0.0; 229807171c5SHong Zhang n = r->i[i+1] - r->i[i]; 230807171c5SHong Zhang rj = r->j + r->i[i]; 231807171c5SHong Zhang ra = r->a + r->i[i]; 232807171c5SHong Zhang for (j=0; j<n; j++) { 233807171c5SHong Zhang r1 += (*ra)*c[*rj]; 234807171c5SHong Zhang r2 += (*ra)*c2[*rj]; 235807171c5SHong Zhang r3 += (*ra)*c3[*rj]; 236807171c5SHong Zhang r4 += (*ra++)*c4[*rj++]; 237807171c5SHong Zhang } 238807171c5SHong Zhang d[colrm + i] = r1; 239807171c5SHong Zhang d[colrm + rm + i] = r2; 240807171c5SHong Zhang d[colrm + rm2 + i] = r3; 241807171c5SHong Zhang d[colrm + rm3 + i] = r4; 242807171c5SHong Zhang } 243807171c5SHong Zhang } 244807171c5SHong Zhang for (; col<cn; col++) { /* over extra columns of C */ 245807171c5SHong Zhang for (i=0; i<am; i++) { /* over rows of A in those columns */ 246807171c5SHong Zhang r1 = 0.0; 247807171c5SHong Zhang n = a->i[i+1] - a->i[i]; 248807171c5SHong Zhang aj = a->j + a->i[i]; 249807171c5SHong Zhang aa = a->a + a->i[i]; 250807171c5SHong Zhang for (j=0; j<n; j++) { 251807171c5SHong Zhang r1 += (*aa++)*b1[*aj++]; 252807171c5SHong Zhang } 253807171c5SHong Zhang c[i] = r1; 254807171c5SHong Zhang } 255807171c5SHong Zhang b1 += bm; 256807171c5SHong Zhang 257807171c5SHong Zhang for (i=0; i<rm; i++) { /* over rows of R in those columns */ 258807171c5SHong Zhang r1 = 0.0; 259807171c5SHong Zhang n = r->i[i+1] - r->i[i]; 260807171c5SHong Zhang rj = r->j + r->i[i]; 261807171c5SHong Zhang ra = r->a + r->i[i]; 262807171c5SHong Zhang for (j=0; j<n; j++) { 263807171c5SHong Zhang r1 += (*ra++)*c[*rj++]; 264807171c5SHong Zhang } 265807171c5SHong Zhang d[col*rm + i] = r1; 266807171c5SHong Zhang } 267807171c5SHong Zhang } 268807171c5SHong Zhang ierr = PetscLogFlops(cn*2.0*(a->nz + r->nz));CHKERRQ(ierr); 269807171c5SHong Zhang 2708c778c55SBarry Smith ierr = MatDenseRestoreArray(B,&b);CHKERRQ(ierr); 2718c778c55SBarry Smith ierr = MatDenseRestoreArray(RAB,&d);CHKERRQ(ierr); 272807171c5SHong Zhang ierr = MatAssemblyBegin(RAB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 273807171c5SHong Zhang ierr = MatAssemblyEnd(RAB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 274807171c5SHong Zhang PetscFunctionReturn(0); 275807171c5SHong Zhang } 276807171c5SHong Zhang 277807171c5SHong Zhang #undef __FUNCT__ 278*55bea0ebSHong Zhang #define __FUNCT__ "MatRARtNumeric_SeqAIJ_SeqAIJ_colorrart" 279*55bea0ebSHong Zhang PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_colorrart(Mat A,Mat R,Mat C) 280807171c5SHong Zhang { 281807171c5SHong Zhang PetscErrorCode ierr; 282257c235dSHong Zhang Mat_SeqAIJ *c = (Mat_SeqAIJ*)C->data; 283257c235dSHong Zhang Mat_RARt *rart=c->rart; 284807171c5SHong Zhang MatTransposeColoring matcoloring; 2858d0a38e4SHong Zhang Mat Rt,RARt; 2863b1b9624SHong Zhang #if defined(RART_PROFILE) 287807171c5SHong Zhang PetscLogDouble Mult_sp_den=0.0,app1=0.0,app2=0.0,t0,tf; 2883b1b9624SHong Zhang #endif 289807171c5SHong Zhang 290807171c5SHong Zhang PetscFunctionBegin; 291807171c5SHong Zhang /* Get dense Rt by Apply MatTransposeColoring to R */ 292807171c5SHong Zhang matcoloring = rart->matcoloring; 293807171c5SHong Zhang Rt = rart->Rt; 2942205254eSKarl Rupp 2953b1b9624SHong Zhang #if defined(RART_PROFILE) 2968563dfccSBarry Smith ierr = PetscTime(&t0);CHKERRQ(ierr); 2973b1b9624SHong Zhang #endif 298807171c5SHong Zhang ierr = MatTransColoringApplySpToDen(matcoloring,R,Rt);CHKERRQ(ierr); 2993b1b9624SHong Zhang #if defined(RART_PROFILE) 3008563dfccSBarry Smith ierr = PetscTime(&tf);CHKERRQ(ierr); 301807171c5SHong Zhang app1 += tf - t0; 3023b1b9624SHong Zhang #endif 303807171c5SHong Zhang 30450647e95SHong Zhang /* Get dense RARt = R*A*Rt -- dominates! */ 3053b1b9624SHong Zhang #if defined(RART_PROFILE) 3068563dfccSBarry Smith ierr = PetscTime(&t0);CHKERRQ(ierr); 3073b1b9624SHong Zhang #endif 308807171c5SHong Zhang RARt = rart->RARt; 309807171c5SHong Zhang ierr = MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqDense(R,A,Rt,RARt,rart->work);CHKERRQ(ierr); 3103b1b9624SHong Zhang #if defined(RART_PROFILE) 3118563dfccSBarry Smith ierr = PetscTime(&tf);CHKERRQ(ierr); 312807171c5SHong Zhang Mult_sp_den += tf - t0; 3133b1b9624SHong Zhang #endif 314807171c5SHong Zhang 315807171c5SHong Zhang /* Recover C from C_dense */ 3163b1b9624SHong Zhang #if defined(RART_PROFILE) 3178563dfccSBarry Smith ierr = PetscTime(&t0);CHKERRQ(ierr); 3183b1b9624SHong Zhang #endif 319807171c5SHong Zhang ierr = MatTransColoringApplyDenToSp(matcoloring,RARt,C);CHKERRQ(ierr); 3203b1b9624SHong Zhang #if defined(RART_PROFILE) 3218563dfccSBarry Smith ierr = PetscTime(&tf);CHKERRQ(ierr); 322807171c5SHong Zhang app2 += tf - t0; 3233b1b9624SHong Zhang #endif 324807171c5SHong Zhang 3253b1b9624SHong Zhang #if defined(RART_PROFILE) 3263b1b9624SHong Zhang printf("MatRARtNumeric_SeqAIJ_SeqAIJ: ColorApp %g + %g + Mult_sp_den %g = %g\n",app1,app2,Mult_sp_den,app1+app2+Mult_sp_den);CHKERRQ(ierr); 327807171c5SHong Zhang #endif 328807171c5SHong Zhang PetscFunctionReturn(0); 329807171c5SHong Zhang } 330807171c5SHong Zhang 331807171c5SHong Zhang #undef __FUNCT__ 332*55bea0ebSHong Zhang #define __FUNCT__ "MatRARtSymbolic_SeqAIJ_SeqAIJ_matmattransposemult" 333*55bea0ebSHong Zhang PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ_matmattransposemult(Mat A,Mat R,PetscReal fill,Mat *C) 334807171c5SHong Zhang { 335807171c5SHong Zhang PetscErrorCode ierr; 336*55bea0ebSHong Zhang Mat ARt,RARt; 3374fa85fdeSHong Zhang Mat_SeqAIJ *c; 3384fa85fdeSHong Zhang Mat_RARt *rart; 3394fa85fdeSHong Zhang 34095a72cc5SHong Zhang PetscFunctionBegin; 341b3b4fc5aSHong Zhang /* must use '-mat_no_inode' with '-matmattransmult_color 1' - do not knwo why? */ 34236104f73SHong Zhang ierr = MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(A,R,fill,&ARt);CHKERRQ(ierr); 34336104f73SHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(R,ARt,fill,&RARt);CHKERRQ(ierr); 3443b1b9624SHong Zhang *C = RARt; 345*55bea0ebSHong Zhang RARt->ops->rartnumeric = MatRARtNumeric_SeqAIJ_SeqAIJ_matmattransposemult; 346*55bea0ebSHong Zhang 3473b1b9624SHong Zhang ierr = PetscNew(Mat_RARt,&rart);CHKERRQ(ierr); 3483b1b9624SHong Zhang c = (Mat_SeqAIJ*)(*C)->data; 3493b1b9624SHong Zhang c->rart = rart; 3503b1b9624SHong Zhang rart->ARt = ARt; 3513b1b9624SHong Zhang rart->destroy = RARt->ops->destroy; 3523b1b9624SHong Zhang RARt->ops->destroy = MatDestroy_SeqAIJ_RARt; 353*55bea0ebSHong Zhang PetscFunctionReturn(0); 35436104f73SHong Zhang } 355*55bea0ebSHong Zhang 356*55bea0ebSHong Zhang #undef __FUNCT__ 357*55bea0ebSHong Zhang #define __FUNCT__ "MatRARtNumeric_SeqAIJ_SeqAIJ_matmattransposemult" 358*55bea0ebSHong Zhang PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ_matmattransposemult(Mat A,Mat R,Mat C) 359*55bea0ebSHong Zhang { 360*55bea0ebSHong Zhang PetscErrorCode ierr; 361*55bea0ebSHong Zhang Mat_SeqAIJ *c=(Mat_SeqAIJ*)C->data; 362*55bea0ebSHong Zhang Mat_RARt *rart=c->rart; 363*55bea0ebSHong Zhang Mat ARt=rart->ARt; 364*55bea0ebSHong Zhang #if defined(RART_PROFILE) 365*55bea0ebSHong Zhang PetscLogDouble t0,t1,t2; 366*55bea0ebSHong Zhang #endif 367*55bea0ebSHong Zhang 368*55bea0ebSHong Zhang PetscFunctionBegin; 3693b1b9624SHong Zhang #if defined(RART_PROFILE) 37026be7272SHong Zhang ierr = PetscTime(&t0);CHKERRQ(ierr); 3713b1b9624SHong Zhang #endif 372b3b4fc5aSHong Zhang ierr = MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(A,R,ARt);CHKERRQ(ierr); /* dominate! */ 3733b1b9624SHong Zhang #if defined(RART_PROFILE) 37426be7272SHong Zhang ierr = PetscTime(&t1);CHKERRQ(ierr); 3753b1b9624SHong Zhang #endif 376*55bea0ebSHong Zhang ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ(R,ARt,C);CHKERRQ(ierr); 3773b1b9624SHong Zhang #if defined(RART_PROFILE) 37826be7272SHong Zhang ierr = PetscTime(&t2);CHKERRQ(ierr); 379b3b4fc5aSHong Zhang printf(" matrart_color_art_num = %g + %g = %g\n",t1-t0,t2-t1,t2-t0); 3803b1b9624SHong Zhang #endif 381*55bea0ebSHong Zhang 38295a72cc5SHong Zhang #if defined(PETSC_USE_INFO) 383*55bea0ebSHong Zhang ierr = PetscInfo(C,"Use ARt=A*R^T, C=R*ARt via MatMatTransposeMult(). Coloring can be applied to A*R^T.\n");CHKERRQ(ierr); 38495a72cc5SHong Zhang #endif 385*55bea0ebSHong Zhang PetscFunctionReturn(0); 386807171c5SHong Zhang } 387*55bea0ebSHong Zhang 388*55bea0ebSHong Zhang #undef __FUNCT__ 389*55bea0ebSHong Zhang #define __FUNCT__ "MatRARtSymbolic_SeqAIJ_SeqAIJ" 390*55bea0ebSHong Zhang PetscErrorCode MatRARtSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat R,PetscReal fill,Mat *C) 391*55bea0ebSHong Zhang { 392*55bea0ebSHong Zhang PetscErrorCode ierr; 39395a72cc5SHong Zhang Mat Rt; 394*55bea0ebSHong Zhang Mat_SeqAIJ *c; 395*55bea0ebSHong Zhang Mat_RARt *rart; 396*55bea0ebSHong Zhang 397*55bea0ebSHong Zhang PetscFunctionBegin; 3985008f5a7SHong Zhang ierr = MatTranspose_SeqAIJ(R,MAT_INITIAL_MATRIX,&Rt);CHKERRQ(ierr); 3995008f5a7SHong Zhang ierr = MatMatMatMultSymbolic_SeqAIJ_SeqAIJ_SeqAIJ(R,A,Rt,fill,C);CHKERRQ(ierr); 40095a72cc5SHong Zhang 40195a72cc5SHong Zhang ierr = PetscNew(Mat_RARt,&rart);CHKERRQ(ierr); 40295a72cc5SHong Zhang rart->Rt = Rt; 40395a72cc5SHong Zhang c = (Mat_SeqAIJ*)(*C)->data; 40495a72cc5SHong Zhang c->rart = rart; 40595a72cc5SHong Zhang rart->destroy = (*C)->ops->destroy; 40695a72cc5SHong Zhang (*C)->ops->destroy = MatDestroy_SeqAIJ_RARt; 407*55bea0ebSHong Zhang (*C)->ops->rartnumeric = MatRARtNumeric_SeqAIJ_SeqAIJ; 408*55bea0ebSHong Zhang PetscFunctionReturn(0); 409*55bea0ebSHong Zhang } 410*55bea0ebSHong Zhang 411*55bea0ebSHong Zhang #undef __FUNCT__ 412*55bea0ebSHong Zhang #define __FUNCT__ "MatRARtNumeric_SeqAIJ_SeqAIJ" 413*55bea0ebSHong Zhang PetscErrorCode MatRARtNumeric_SeqAIJ_SeqAIJ(Mat A,Mat R,Mat C) 414*55bea0ebSHong Zhang { 415*55bea0ebSHong Zhang PetscErrorCode ierr; 416*55bea0ebSHong Zhang Mat_SeqAIJ *c = (Mat_SeqAIJ*)C->data; 417*55bea0ebSHong Zhang Mat_RARt *rart = c->rart; 418*55bea0ebSHong Zhang Mat Rt = rart->Rt; 419*55bea0ebSHong Zhang 420*55bea0ebSHong Zhang PetscFunctionBegin; 421*55bea0ebSHong Zhang ierr = MatTranspose_SeqAIJ(R,MAT_REUSE_MATRIX,&Rt);CHKERRQ(ierr); 422*55bea0ebSHong Zhang ierr = MatMatMatMultNumeric_SeqAIJ_SeqAIJ_SeqAIJ(R,A,Rt,C);CHKERRQ(ierr); 423*55bea0ebSHong Zhang #if defined(PETSC_USE_INFO) 424*55bea0ebSHong Zhang ierr = PetscInfo(C,"Use Rt=R^T and C=R*A*Rt via MatMatMatMult() to avoid sparse inner products\n");CHKERRQ(ierr); 425*55bea0ebSHong Zhang #endif 426*55bea0ebSHong Zhang PetscFunctionReturn(0); 427*55bea0ebSHong Zhang } 428*55bea0ebSHong Zhang 429*55bea0ebSHong Zhang #undef __FUNCT__ 430*55bea0ebSHong Zhang #define __FUNCT__ "MatRARt_SeqAIJ_SeqAIJ" 431*55bea0ebSHong Zhang PetscErrorCode MatRARt_SeqAIJ_SeqAIJ(Mat A,Mat R,MatReuse scall,PetscReal fill,Mat *C) 432*55bea0ebSHong Zhang { 433*55bea0ebSHong Zhang PetscErrorCode ierr; 434*55bea0ebSHong Zhang const char *algTypes[3] = {"matmatmatmult","matmattransposemult","coloring_rart"}; 435*55bea0ebSHong Zhang PetscInt alg=0; /* set default algorithm */ 436*55bea0ebSHong Zhang 437*55bea0ebSHong Zhang PetscFunctionBegin; 438*55bea0ebSHong Zhang if (scall == MAT_INITIAL_MATRIX) { 439*55bea0ebSHong Zhang /* ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr); -- prefix ? */ 440*55bea0ebSHong Zhang ierr = PetscOptionsEList("-matrart_via","Algorithmic approach","MatRARt",algTypes,3,algTypes[0],&alg,NULL);CHKERRQ(ierr); 441*55bea0ebSHong Zhang /* ierr = PetscOptionsEnd();CHKERRQ(ierr); */ 442*55bea0ebSHong Zhang ierr = PetscLogEventBegin(MAT_RARtSymbolic,A,R,0,0);CHKERRQ(ierr); 443*55bea0ebSHong Zhang switch (alg) { 444*55bea0ebSHong Zhang case 1: 445*55bea0ebSHong Zhang /* via matmattransposemult: ARt=A*R^T, C=R*ARt - matrix coloring can be applied to A*R^T */ 446*55bea0ebSHong Zhang ierr = MatRARtSymbolic_SeqAIJ_SeqAIJ_matmattransposemult(A,R,fill,C);CHKERRQ(ierr); 447*55bea0ebSHong Zhang break; 448*55bea0ebSHong Zhang case 2: 449*55bea0ebSHong Zhang /* via coloring_rart: apply coloring C = R*A*R^T */ 450*55bea0ebSHong Zhang ierr = MatRARtSymbolic_SeqAIJ_SeqAIJ_colorrart(A,R,fill,C);CHKERRQ(ierr); 451*55bea0ebSHong Zhang break; 452*55bea0ebSHong Zhang default: 453*55bea0ebSHong Zhang /* via matmatmatmult: Rt=R^T, C=R*A*Rt - avoid inefficient sparse inner products */ 454*55bea0ebSHong Zhang ierr = MatRARtSymbolic_SeqAIJ_SeqAIJ(A,R,fill,C);CHKERRQ(ierr); 455*55bea0ebSHong Zhang break; 456*55bea0ebSHong Zhang } 4575008f5a7SHong Zhang ierr = PetscLogEventEnd(MAT_RARtSymbolic,A,R,0,0);CHKERRQ(ierr); 4585008f5a7SHong Zhang } 4595008f5a7SHong Zhang 4605008f5a7SHong Zhang ierr = PetscLogEventBegin(MAT_RARtNumeric,A,R,0,0);CHKERRQ(ierr); 461*55bea0ebSHong Zhang ierr = (*(*C)->ops->rartnumeric)(A,R,*C);CHKERRQ(ierr); 4625008f5a7SHong Zhang ierr = PetscLogEventEnd(MAT_RARtNumeric,A,R,0,0);CHKERRQ(ierr); 463807171c5SHong Zhang PetscFunctionReturn(0); 464807171c5SHong Zhang } 465